LifeV
QuadratureBoundary.hpp
Go to the documentation of this file.
1 //@HEADER
2 /*
3 *******************************************************************************
4 
5  Copyright (C) 2004, 2005, 2007 EPFL, Politecnico di Milano, INRIA
6  Copyright (C) 2010 EPFL, Politecnico di Milano, Emory University
7 
8  This file is part of LifeV.
9 
10  LifeV is free software; you can redistribute it and/or modify
11  it under the terms of the GNU Lesser General Public License as published by
12  the Free Software Foundation, either version 3 of the License, or
13  (at your option) any later version.
14 
15  LifeV is distributed in the hope that it will be useful,
16  but WITHOUT ANY WARRANTY; without even the implied warranty of
17  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
18  Lesser General Public License for more details.
19 
20  You should have received a copy of the GNU Lesser General Public License
21  along with LifeV. If not, see <http://www.gnu.org/licenses/>.
22 
23 *******************************************************************************
24 */
25 //@HEADER
26 
27 /*!
28 
29  */
30 
31 #ifndef QUADRATUREBOUNDARY_H
32 #define QUADRATUREBOUNDARY_H 1
33 
34 #include <lifev/core/LifeV.hpp>
35 
36 #include <lifev/core/fem/QuadratureRule.hpp>
37 
38 #include <vector>
39 
40 namespace LifeV
41 {
42 
43 //! QuadratureBoundary - Short description of the class
44 
45 class QuadratureBoundary
46 {
47 public:
48 
49  QuadratureBoundary() {};
50 
51  QuadratureBoundary ( const QuadratureBoundary& qrbd )
52  : M_qrs (qrbd.M_qrs)
53  {};
54 
55 
56  void setQuadrature (const UInt i, const QuadratureRule& qr)
57  {
58  while (i >= M_qrs.size() )
59  {
60  M_qrs.push_back (QuadratureRule() );
61  }
62 
63  M_qrs[i] = qr;
64  }
65 
66  QuadratureRule qr (const UInt i) const
67  {
68  ASSERT (i < M_qrs.size(), "Invalid quadrature number");
69  return M_qrs[i];
70  }
71 
72 
73  void showMe() const
74  {
75  for (UInt i (0); i < M_qrs.size(); ++i)
76  {
77  std::cout << " ## QR " << i << " ## " << std::endl;
78  M_qrs[i].showMe();
79  }
80  }
81 private:
82 
83  std::vector< QuadratureRule > M_qrs;
84 
85 };
86 
87 inline
88 QuadratureBoundary
89 buildTetraBDQR (const QuadratureRule& my_qr)
90 {
91  QuadratureBoundary qrbd;
92 
93  // Face 0
94  qrbd.setQuadrature (0, my_qr);
95 
96  // Face 1
97  QuadratureRule qf1 ("none", TRIANGLE, 3, 0, 0);
98  for (UInt iq (0); iq < my_qr.nbQuadPt(); ++iq)
99  {
100  Real x (my_qr.quadPointCoor (iq, 0) );
101  Real y (my_qr.quadPointCoor (iq, 1) );
102  Real w (my_qr.weight (iq) );
103 
104  qf1.addPoint (QuadraturePoint (x, 0, y, w) );
105  }
106  qrbd.setQuadrature (1, qf1);
107 
108  // Face 2
109  QuadratureRule qf2 ("none", TRIANGLE, 3, 0, 0);
110  for (UInt iq (0); iq < my_qr.nbQuadPt(); ++iq)
111  {
112  Real x (my_qr.quadPointCoor (iq, 0) );
113  Real y (my_qr.quadPointCoor (iq, 1) );
114  Real w (my_qr.weight (iq) );
115 
116  qf2.addPoint (QuadraturePoint (x, y, 1 - x - y, w) );
117  }
118  qrbd.setQuadrature (2, qf2);
119 
120  // Face 3
121  QuadratureRule qf3 ("none", TRIANGLE, 3, 0, 0);
122  for (UInt iq (0); iq < my_qr.nbQuadPt(); ++iq)
123  {
124  Real x (my_qr.quadPointCoor (iq, 0) );
125  Real y (my_qr.quadPointCoor (iq, 1) );
126  Real w (my_qr.weight (iq) );
127 
128  qf3.addPoint (QuadraturePoint (0, y, x, w) );
129  }
130  qrbd.setQuadrature (3, qf3);
131 
132  return qrbd;
133 
134 }
135 
136 } // Namespace LifeV
137 
138 #endif /* QUADRATUREBOUNDARY_H */
#define ASSERT(X, A)
Definition: LifeAssert.hpp:90