LifeV
QRKeast.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  @file
29  @brief File for the definition of the quadrature rule defined in
30  the article [Keast85]
31 
32  @author Samuel Quinodoz <samuel.quinodoz@epfl.ch>
33  @date 12-2011
34 
35  @contributor Samuel Quinodoz <samuel.quinodoz@epfl.ch>
36  @mantainer Samuel Quinodoz <samuel.quinodoz@epfl.ch>
37  */
38 
39 
40 #ifndef QR_KEAST_HPP
41 #define QR_KEAST_HPP
42 
43 #include <lifev/core/LifeV.hpp>
44 
45 #include <lifev/core/array/VectorSmall.hpp>
46 
47 #include <lifev/core/mesh/ElementShapes.hpp>
48 
49 #include <lifev/core/fem/QuadraturePoint.hpp>
50 
51 #include <iostream>
52 #include <fstream>
53 
54 namespace LifeV
55 {
56 
57 
58 //! QRKeast - A set of quadrature for tetrahedra
59 /*!
60  We concentrate on the lowest degree, so negative weights
61  are possible!
62 */
63 
64 template< UInt degree >
65 class QRKeast
66 {
67 public:
68 
69 private:
70  QRKeast();
71  ~QRKeast();
72  QRKeast ( const QRKeast&);
73 };
74 
75 
76 /* Specializations */
77 template <>
78 class QRKeast<1>
79 {
80 public:
81 
82  QRKeast() {};
83  ~QRKeast() {};
84 
85  static const Real& weight (const UInt& iPt)
86  {
87  return M_weights[iPt];
88  }
89 
90  static const Real& quadPointCoor ( const UInt iPt, const UInt iCoor)
91  {
92  return M_points[iPt][iCoor];
93  }
94 
95  static VectorSmall<3> quadPointCoor ( const UInt iPt)
96  {
97  return VectorSmall<3> (M_points[iPt][0], M_points[iPt][1], M_points[iPt][2]);
98  }
99 
100  static QuadraturePoint quadPoint (const UInt iPt)
101  {
102  return QuadraturePoint (M_points[iPt][0], M_points[iPt][1], M_points[iPt][2], M_weights[iPt]);
103  }
104 
105  static UInt nbQuadPt()
106  {
107  return 1;
108  }
109 
110  static ReferenceShapes shape()
111  {
112  return TETRA;
113  }
114 
115  static UInt dimension()
116  {
117  return 3;
118  }
119 
120 private:
121 
122  static const Real M_points[1][3];
123  static const Real M_weights[1];
124 };
125 
126 template <>
127 class QRKeast<4>
128 {
129 public:
130 
131  QRKeast() {};
132  ~QRKeast() {};
133 
134  static const Real& weight (const UInt& iPt)
135  {
136  return M_weights[iPt];
137  }
138 
139  static const Real& quadPointCoor ( const UInt iPt, const UInt iCoor)
140  {
141  return M_points[iPt][iCoor];
142  }
143 
144  static VectorSmall<3> quadPointCoor ( const UInt iPt)
145  {
146  return VectorSmall<3> (M_points[iPt][0], M_points[iPt][1], M_points[iPt][2]);
147  }
148 
149  static QuadraturePoint quadPoint (const UInt iPt)
150  {
151  return QuadraturePoint (M_points[iPt][0], M_points[iPt][1], M_points[iPt][2], M_weights[iPt]);
152  }
153 
154  static UInt nbQuadPt()
155  {
156  return 11;
157  }
158 
159  static ReferenceShapes shape()
160  {
161  return TETRA;
162  }
163 
164  static UInt dimension()
165  {
166  return 3;
167  }
168 
169 private:
170 
171  static const Real M_points[11][3];
172  static const Real M_weights[11];
173 };
174 
175 template <>
176 class QRKeast<6>
177 {
178 public:
179 
180  QRKeast() {};
181  ~QRKeast() {};
182 
183  static const Real& weight (const UInt& iPt)
184  {
185  return M_weights[iPt];
186  }
187 
188  static const Real& quadPointCoor ( const UInt iPt, const UInt iCoor)
189  {
190  return M_points[iPt][iCoor];
191  }
192 
193  static VectorSmall<3> quadPointCoor ( const UInt iPt)
194  {
195  return VectorSmall<3> (M_points[iPt][0], M_points[iPt][1], M_points[iPt][2]);
196  }
197 
198  static QuadraturePoint quadPoint (const UInt iPt)
199  {
200  return QuadraturePoint (M_points[iPt][0], M_points[iPt][1], M_points[iPt][2], M_weights[iPt]);
201  }
202 
203  static UInt nbQuadPt()
204  {
205  return 24;
206  }
207 
208  static ReferenceShapes shape()
209  {
210  return TETRA;
211  }
212 
213  static UInt dimension()
214  {
215  return 3;
216  }
217 
218 private:
219 
220  static const Real M_points[24][3];
221  static const Real M_weights[24];
222 };
223 
224 
225 template <>
226 class QRKeast<7>
227 {
228 public:
229 
230  QRKeast() {};
231  ~QRKeast() {};
232 
233  static const Real& weight (const UInt& iPt)
234  {
235  return M_weights[iPt];
236  }
237 
238  static const Real& quadPointCoor ( const UInt iPt, const UInt iCoor)
239  {
240  return M_points[iPt][iCoor];
241  }
242 
243  static VectorSmall<3> quadPointCoor ( const UInt iPt)
244  {
245  return VectorSmall<3> (M_points[iPt][0], M_points[iPt][1], M_points[iPt][2]);
246  }
247 
248  static QuadraturePoint quadPoint (const UInt iPt)
249  {
250  return QuadraturePoint (M_points[iPt][0], M_points[iPt][1], M_points[iPt][2], M_weights[iPt]);
251  }
252 
253  static UInt nbQuadPt()
254  {
255  return 31;
256  }
257 
258  static ReferenceShapes shape()
259  {
260  return TETRA;
261  }
262 
263  static UInt dimension()
264  {
265  return 3;
266  }
267 
268 private:
269 
270  static const Real M_points[31][3];
271  static const Real M_weights[31];
272 };
273 
274 
275 }
276 #endif
QuadraturePoint - Simple container for a point of a quadrature rule.
QuadraturePoint(Real x, Real y, Real z, Real weight)
Full constructor for 3D.
double Real
Generic real data.
Definition: LifeV.hpp:175
VectorSmall(Real const &x, Real const &y, Real const &z)
Full constructor with all components explicitly initialized.
uint32_type UInt
generic unsigned integer (used mainly for addressing)
Definition: LifeV.hpp:191