Directory: | ./ |
---|---|
File: | FiniteElement/curBaseFiniteElement.hpp |
Date: | 2024-04-14 07:32:34 |
Exec | Total | Coverage | |
---|---|---|---|
Lines: | 50 | 53 | 94.3% |
Branches: | 7 | 14 | 50.0% |
Line | Branch | Exec | Source |
---|---|---|---|
1 | // ______ ______ _ _ _____ ______ | ||
2 | // | ____| ____| | (_)/ ____| | ____| | ||
3 | // | |__ | |__ | | _| (___ ___| |__ | ||
4 | // | __| | __| | | | |\___ \ / __| __| | ||
5 | // | | | |____| |____| |____) | (__| |____ | ||
6 | // |_| |______|______|_|_____/ \___|______| | ||
7 | // Finite Elements for Life Sciences and Engineering | ||
8 | // | ||
9 | // License: LGL2.1 License | ||
10 | // FELiScE default license: LICENSE in root folder | ||
11 | // | ||
12 | // Main authors: J-F. Gerbeau and V. Martin | ||
13 | // | ||
14 | |||
15 | #ifndef _CUR_BASE_FINITE_ELEMENT_H | ||
16 | #define _CUR_BASE_FINITE_ELEMENT_H | ||
17 | |||
18 | // System includes | ||
19 | #include <vector> | ||
20 | #include <atomic> | ||
21 | |||
22 | // External includes | ||
23 | |||
24 | // Project includes | ||
25 | #include "geoElement.hpp" | ||
26 | #include "refElement.hpp" | ||
27 | #include "quadratureRule.hpp" | ||
28 | #include "Core/felisce.hpp" | ||
29 | #include "Core/shared_pointers.hpp" | ||
30 | |||
31 | namespace felisce | ||
32 | { | ||
33 | ///@name felisce Globals | ||
34 | ///@{ | ||
35 | |||
36 | ///@} | ||
37 | ///@name Type Definitions | ||
38 | ///@{ | ||
39 | |||
40 | ///@} | ||
41 | ///@name Enum's | ||
42 | ///@{ | ||
43 | |||
44 | ///@} | ||
45 | ///@name Functions | ||
46 | ///@{ | ||
47 | |||
48 | ///@} | ||
49 | ///@name felisce Classes | ||
50 | ///@{ | ||
51 | |||
52 | /** | ||
53 | * @class CurBaseFiniteElement | ||
54 | * @authors J-F. Gerbeau and V. Martin | ||
55 | * @brief Base class for the current and curvilinear finite element | ||
56 | */ | ||
57 | class CurBaseFiniteElement | ||
58 | { | ||
59 | public: | ||
60 | ///@name Type Definitions | ||
61 | ///@{ | ||
62 | |||
63 | /// Pointer definition of Point | ||
64 | FELISCE_CLASS_POINTER_DEFINITION(CurBaseFiniteElement); | ||
65 | |||
66 | ///@} | ||
67 | ///@name Life Cycle | ||
68 | ///@{ | ||
69 | |||
70 | /// constructors (vector of DegreeOfExactness) | ||
71 | CurBaseFiniteElement( | ||
72 | const RefElement& refEle, | ||
73 | const GeoElement& geoEle, | ||
74 | const std::vector<DegreeOfExactness>& degOfExactness | ||
75 | ); | ||
76 | |||
77 | /// constructors | ||
78 | CurBaseFiniteElement( | ||
79 | const RefElement& refEle, | ||
80 | const GeoElement& geoEle, | ||
81 | const DegreeOfExactness& degOfExactness | ||
82 | ); | ||
83 | |||
84 | /// Destructor | ||
85 | virtual ~CurBaseFiniteElement(); | ||
86 | |||
87 | /// Copy constructor | ||
88 | 16384 | CurBaseFiniteElement(const CurBaseFiniteElement& rOther) | |
89 | 16384 | : phi(rOther.phi), | |
90 |
1/2✓ Branch 1 taken 16384 times.
✗ Branch 2 not taken.
|
16384 | weightMeas(rOther.weightMeas), |
91 |
1/2✓ Branch 1 taken 16384 times.
✗ Branch 2 not taken.
|
16384 | currentQuadPoint(rOther.currentQuadPoint), |
92 | 16384 | m_refEle(rOther.m_refEle), | |
93 | 16384 | m_geoEle(rOther.m_geoEle), | |
94 | 16384 | m_quadratureRule(rOther.m_quadratureRule), | |
95 | 16384 | m_numDof(rOther.m_numDof), | |
96 | 16384 | m_numQuadraturePoint(rOther.m_numQuadraturePoint), | |
97 | 16384 | m_numRefCoor(rOther.m_numRefCoor), | |
98 | 16384 | m_numCoor(rOther.m_numCoor), | |
99 | 16384 | m_numPoint(rOther.m_numPoint), | |
100 | 16384 | m_hasPoint(rOther.m_hasPoint), | |
101 | 16384 | m_hasMeas(rOther.m_hasMeas), | |
102 | 16384 | m_hasQuadPtCoor(rOther.m_hasQuadPtCoor), | |
103 | 16384 | m_hasOriginalQuadPoint(rOther.m_hasOriginalQuadPoint), | |
104 | 16384 | m_currentId(rOther.m_currentId), | |
105 |
1/2✓ Branch 1 taken 16384 times.
✗ Branch 2 not taken.
|
16384 | m_dPhiRef(rOther.m_dPhiRef), |
106 |
1/2✓ Branch 1 taken 16384 times.
✗ Branch 2 not taken.
|
16384 | m_phiGeo(rOther.m_phiGeo), |
107 |
1/2✓ Branch 1 taken 16384 times.
✗ Branch 2 not taken.
|
16384 | m_dPhiGeo(rOther.m_dPhiGeo), |
108 |
1/2✓ Branch 1 taken 16384 times.
✗ Branch 2 not taken.
|
16384 | m_meas(rOther.m_meas), |
109 |
1/2✓ Branch 2 taken 16384 times.
✗ Branch 3 not taken.
|
32768 | m_point(rOther.m_point) |
110 | 16384 | {} | |
111 | |||
112 | ///@} | ||
113 | ///@name Operators | ||
114 | ///@{ | ||
115 | |||
116 | ///@} | ||
117 | ///@name Operations | ||
118 | ///@{ | ||
119 | |||
120 | // Do nothing: this method is just there to provide a non pure virtual/non inline method to avoid emitting | ||
121 | // weak tables in every translation unit. | ||
122 | virtual void placeholder() const; | ||
123 | |||
124 | |||
125 | |||
126 | /// minimal-minimal update called by all updates: set the points only | ||
127 | void updatePoint(const std::vector<Point*>& point); | ||
128 | |||
129 | void updatePoint(const UBlasMatrix& point, const std::vector<int>& ipt); | ||
130 | /*! | ||
131 | compute the coordinate (x,y,z)= F(xi,eta,zeta), (F: geo mappping) | ||
132 | where (xi,eta,zeta) are the coor in the ref element | ||
133 | and (x,y,z) the coor in the current element | ||
134 | (in 2D, z=0 and zeta is discarded) | ||
135 | */ | ||
136 | void coorMap( Point& pt, const Point& refpt ) const; | ||
137 | |||
138 | /// Return the measure of the current element | ||
139 | double measure() const; | ||
140 | |||
141 | /// Return the diameter of the current element | ||
142 | double diameter() const; | ||
143 | |||
144 | /// measure of the segment from id1 and id2 | ||
145 | double measOfSegment(felInt id1, felInt id2) const; | ||
146 | //----------- | ||
147 | // updates: | ||
148 | //----------- | ||
149 | /// minimal update: we just identify the id of the current element | ||
150 | virtual void update(const int id,const std::vector<Point*>& point); | ||
151 | |||
152 | /// compute the arrays meas, weightDet, jacobian on the current element | ||
153 | 46 | virtual void updateMeas(const int /*id*/,const std::vector<Point*>& /*point*/) {}; | |
154 | |||
155 | /// compute the images of the quadrature points through the geometrical mapping | ||
156 | /// note: the Boundary quad points are computed separately in CurrFEWithBd | ||
157 | void computeCurrentQuadraturePoint(); | ||
158 | |||
159 | /// map point in reference element | ||
160 | void mapPointInReferenceElement(const int id, const std::vector<Point*>& ptElem, const Point& inPoint, Point& outPoint); | ||
161 | |||
162 | ///@} | ||
163 | ///@name Access | ||
164 | ///@{ | ||
165 | |||
166 | inline UBlasVector& jacobian() { | ||
167 | return m_meas; | ||
168 | } | ||
169 | |||
170 | 225516 | inline const RefElement& refEle() const { | |
171 | 225516 | return m_refEle; | |
172 | } | ||
173 | |||
174 | 1947335 | inline const GeoElement& geoEle() const { | |
175 | 1947335 | return m_geoEle; | |
176 | } | ||
177 | |||
178 | 1660524 | inline const QuadratureRule& quadratureRule() const { | |
179 | 1660524 | return *m_quadratureRule; | |
180 | } | ||
181 | |||
182 | 10573979651 | inline int numDof() const { | |
183 | 10573979651 | return m_numDof; | |
184 | } | ||
185 | |||
186 | 6671075245 | inline int numQuadraturePoint() const { | |
187 | 6671075245 | return m_numQuadraturePoint; | |
188 | } | ||
189 | |||
190 | 1527016360 | inline int numCoor() const { | |
191 | 1527016360 | return m_numCoor; | |
192 | } | ||
193 | |||
194 | 6575723 | inline int numRefCoor() const { | |
195 | 6575723 | return m_numRefCoor; | |
196 | } | ||
197 | |||
198 | 3940428 | inline int numPoint() const { | |
199 | 3940428 | return m_numPoint; | |
200 | } | ||
201 | |||
202 | /// return true if the point has been updated | ||
203 | 2059890 | inline bool hasPoint() const { | |
204 | 2059890 | return m_hasPoint; | |
205 | } | ||
206 | |||
207 | /// return true if the measure has been updated | ||
208 | ✗ | inline bool hasMeas() const { | |
209 | ✗ | return m_hasMeas; | |
210 | } | ||
211 | |||
212 | /// return true if the coordinates of the quadrature points have been updated | ||
213 | 554928 | inline bool hasQuadPtCoor() const { | |
214 | 554928 | return m_hasQuadPtCoor; | |
215 | } | ||
216 | |||
217 | /// return true if the basis function are computed at the original quadrature point | ||
218 | 270288 | inline bool hasOriginalQuadPoint() const { | |
219 | 270288 | return m_hasOriginalQuadPoint; | |
220 | } | ||
221 | |||
222 | /// return true id of the element | ||
223 | inline int id() const { | ||
224 | return m_currentId; | ||
225 | } | ||
226 | |||
227 | /// access to dPhiGeo is necessary | ||
228 | 1543099 | const UBlasMatrix& dPhiGeo(int ig) const { | |
229 | 1543099 | return m_dPhiGeo[ig]; | |
230 | } | ||
231 | |||
232 | /// access to phiGeo is necessary | ||
233 | 1543099 | const UBlasVector& phiGeo(int ig) const { | |
234 | 1543099 | return m_phiGeo[ig]; | |
235 | } | ||
236 | |||
237 | ///@} | ||
238 | ///@name Public member Variables | ||
239 | ///@{ | ||
240 | |||
241 | /// phi[ig](idof): (in the current FE (and equal to phiRef)) | ||
242 | std::vector<UBlasVector> phi; | ||
243 | |||
244 | UBlasVector weightMeas; | ||
245 | |||
246 | /// images of the quadrature points on the current element | ||
247 | std::vector<Point> currentQuadPoint; | ||
248 | |||
249 | ///@} | ||
250 | ///@name Inquiry | ||
251 | ///@{ | ||
252 | |||
253 | ///@} | ||
254 | ///@name Input and output | ||
255 | ///@{ | ||
256 | |||
257 | void print(int verbose,std::ostream& c=std::cout) const; | ||
258 | |||
259 | void m_printBase(int verbose,std::ostream& c=std::cout) const; | ||
260 | |||
261 | ///@} | ||
262 | protected: | ||
263 | ///@name Protected static Member Variables | ||
264 | ///@{ | ||
265 | |||
266 | ///@} | ||
267 | ///@name Protected member Variables | ||
268 | ///@{ | ||
269 | |||
270 | // TODO: Define flags!! | ||
271 | // TODO: Move definitions to geometries | ||
272 | const RefElement& m_refEle; | ||
273 | const GeoElement& m_geoEle; | ||
274 | const QuadratureRule* m_quadratureRule; | ||
275 | int m_numDof; | ||
276 | int m_numQuadraturePoint; | ||
277 | int m_numRefCoor; /// Space dimension in the reference element | ||
278 | int m_numCoor; /// Space dimension (=numRefCoor (CurrentFE) or numRefCoor+1 (CurvilinearFE)) | ||
279 | int m_numPoint; | ||
280 | bool m_hasPoint = false; | ||
281 | bool m_hasMeas = false; | ||
282 | bool m_hasQuadPtCoor = false; | ||
283 | bool m_hasOriginalQuadPoint = true; | ||
284 | int m_currentId = 0; | ||
285 | |||
286 | /// m_dPhiRef[ig](icoor,jdof) = derivative with respect to x_icoor of the jdof-th reference basis function at integration point ig | ||
287 | std::vector<UBlasMatrix> m_dPhiRef; | ||
288 | |||
289 | /// Geometric reference functions m_phiGeo[ig](ipoint) : value of the ipoint-th function at integration point ig | ||
290 | std::vector<UBlasVector> m_phiGeo; | ||
291 | |||
292 | /// Derivatives of the geometric reference functions m_dPhiGeo[ig](icoor,jpoint) : value of the jpoint-th function at integration point ig | ||
293 | std::vector<UBlasMatrix> m_dPhiGeo; | ||
294 | |||
295 | /// m_meas(ig): determinant of the jacobian at the integration point ig (in the current FE) | ||
296 | UBlasVector m_meas; | ||
297 | |||
298 | /// Coordinates of the points defining the element: m_point(ipoint,jcoor) (in the current FE) | ||
299 | UBlasMatrix m_point; | ||
300 | |||
301 | ///@} | ||
302 | ///@name Protected Operators | ||
303 | ///@{ | ||
304 | |||
305 | ///@} | ||
306 | ///@name Protected Operations | ||
307 | ///@{ | ||
308 | |||
309 | /// compute the integration weights | ||
310 | ✗ | virtual void computeMeas() {}; | |
311 | |||
312 | /// update basis function with the new quadrature point | ||
313 | void updateBasisWithQuadPoint(const std::vector<Point>& point); | ||
314 | |||
315 | ///@} | ||
316 | ///@name Protected Access | ||
317 | ///@{ | ||
318 | |||
319 | ///@} | ||
320 | ///@name Protected Inquiry | ||
321 | ///@{ | ||
322 | |||
323 | ///@} | ||
324 | ///@name Protected LifeCycle | ||
325 | ///@{ | ||
326 | |||
327 | ///@} | ||
328 | private: | ||
329 | ///@name Private static Member Variables | ||
330 | ///@{ | ||
331 | |||
332 | ///@} | ||
333 | ///@name Private member Variables | ||
334 | ///@{ | ||
335 | |||
336 | ///@} | ||
337 | ///@name Private Operators | ||
338 | ///@{ | ||
339 | |||
340 | ///@} | ||
341 | ///@name Private Operations | ||
342 | ///@{ | ||
343 | |||
344 | ///@} | ||
345 | ///@name Private Access | ||
346 | ///@{ | ||
347 | |||
348 | ///@} | ||
349 | ///@name Private Inquiry | ||
350 | ///@{ | ||
351 | |||
352 | ///@} | ||
353 | ///@name Private LifeCycle | ||
354 | ///@{ | ||
355 | |||
356 | ///@} | ||
357 | |||
358 | }; // class curBaseFiniteElement | ||
359 | |||
360 | ///@} | ||
361 | ///@name Type Definitions | ||
362 | ///@{ | ||
363 | |||
364 | ///@} | ||
365 | } /* namespace felisce.*/ | ||
366 | |||
367 | #endif /* _CUR_BASE_FINITE_ELEMENT_H defined */ | ||
368 |