v0.14.0
Loading...
Searching...
No Matches
FatPrismElementForcesAndSourcesCore.hpp
Go to the documentation of this file.
1/** \file ForcesAndSourcesCore.hpp
2
3 \brief Implementation of elements on entities.
4
5 Those element are inherited by user to implement specific implementation of
6 particular problem.
7
8*/
9
10
11
12#ifndef __FATPRISMELEMENTFORCESANDSURCESCORE_HPP__
13#define __FATPRISMELEMENTFORCESANDSURCESCORE_HPP__
14
15using namespace boost::numeric;
16
17namespace MoFEM {
18
19/** \brief FatPrism finite element
20 \ingroup mofem_forces_and_sources_prism_element
21
22 User is implementing own operator at Gauss points level, by own object
23 derived from FatPrismElementForcesAndSourcesCoreL::UserDataOperator. Arbitrary
24 number of operator added pushing objects to rowOpPtrVector and
25 rowColOpPtrVector.
26
27 \todo Need to implement operators that will make this element work as Volume
28 element
29
30 */
33
35
36 virtual int getRuleTrianglesOnly(int order) { return 2 * order; };
37 virtual int getRuleThroughThickness(int order) { return 2 * order; };
38
39 virtual MoFEMErrorCode setGaussPtsTrianglesOnly(int order_triangles_only) {
41 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
43 }
44
45 virtual MoFEMErrorCode setGaussPtsThroughThickness(int order_thickness) {
47 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
49 }
50
51 /** \brief default operator for Flat Prism element
52 * \ingroup mofem_forces_and_sources_prism_element
53 */
56
57 using VolumeElementForcesAndSourcesCore::UserDataOperator::UserDataOperator;
58
59 /** \brief get face aRea
60 \param dd if dd == 0 it is for face F3 if dd == 1 is for face F4
61 */
62 inline double getArea(const int dd);
63
64 inline double getAreaF3();
65 inline double getAreaF4();
66
67 /** \brief get triangle normal
68
69 Normal has 6 elements, first 3 are for face F3 another three for face F4
70
71 */
72 inline VectorDouble &getNormal();
73
75
77
78 /** \brief get Gauss pts. in the prism
79 */
80 inline MatrixDouble &getGaussPts();
81
82 /** \brief get Gauss pts. on triangles
83 */
85
86 /** \brief get Gauss pts. through thickness
87 */
89
90 /** \brief get coordinates at Gauss pts.
91
92 Matrix has size (nb integration points)x(coordinates on F3 and F4 = 6),
93 i.e. coordinates on face F3 and F4
94
95 */
97
98 /** \brief get coordinates at Gauss pts.
99
100 Matrix has size (nb integration points)x(coordinates on F3 and F4 = 6),
101 i.e. coordinates on face F3 and F4
102
103 */
105
106 /** \brief coordinate at Gauss points on face 3 (if hierarchical
107 * approximation of element geometry)
108 */
110
111 /** \brief coordinate at Gauss points on face 4 (if hierarchical
112 * approximation of element geometry)
113 */
115
116 /** \brief if higher order geometry return normals at face F3 at Gauss pts.
117 *
118 * Face 3 is top face in canonical triangle numeration, see \cite
119 * tautges2010canonical
120 *
121 */
123
124 /** \brief if higher order geometry return normals at face F4 at Gauss pts.
125 *
126 * Face 4 is top face in canonical triangle numeration, see \cite
127 * tautges2010canonical
128 *
129 */
131
132 /** \brief if higher order geometry return normals at Gauss pts.
133 *
134 * Face 3 is top face in canonical triangle numeration, see \cite
135 * tautges2010canonical
136 *
137 * \param gg gauss point number
138 */
139 inline ublas::matrix_row<MatrixDouble> getNormalsAtGaussPtsF3(const int gg);
140
141 /** \brief if higher order geometry return normals at Gauss pts.
142 *
143 * Face 3 is top face in canonical triangle numeration, see \cite
144 * tautges2010canonical
145 *
146 * \param gg gauss point number
147 */
148 inline ublas::matrix_row<MatrixDouble> getNormalsAtGaussPtsF4(const int gg);
149
150 /** \brief if higher order geometry return tangent vector to triangle at
151 * Gauss pts.
152 */
154
155 /** \brief if higher order geometry return tangent vector to triangle at
156 * Gauss pts.
157 */
159
160 /** \brief if higher order geometry return tangent vector to triangle at
161 * Gauss pts.
162 */
164
165 /** \brief if higher order geometry return tangent vector to triangle at
166 * Gauss pts.
167 */
169
171
173
174 /** \brief return pointer to fat prism finite element
175 */
177
178 protected:
180
181 };
182
184
185protected:
186 double aRea[2];
188
192
195
205
206 friend class UserDataOperator;
207};
208
209inline double
210FatPrismElementForcesAndSourcesCore::UserDataOperator::getArea(const int dd) {
211 return static_cast<FatPrismElementForcesAndSourcesCore *>(ptrFE)->aRea[0];
212}
213
214inline double
215FatPrismElementForcesAndSourcesCore::UserDataOperator::getAreaF3() {
216 return static_cast<FatPrismElementForcesAndSourcesCore *>(ptrFE)->aRea[0];
217}
218inline double
219FatPrismElementForcesAndSourcesCore::UserDataOperator::getAreaF4() {
220 return static_cast<FatPrismElementForcesAndSourcesCore *>(ptrFE)->aRea[1];
221}
222
223inline VectorDouble &
224FatPrismElementForcesAndSourcesCore::UserDataOperator::getNormal() {
225 return static_cast<FatPrismElementForcesAndSourcesCore *>(ptrFE)->normal;
226}
227
228inline VectorAdaptor
229FatPrismElementForcesAndSourcesCore::UserDataOperator::getNormalF3() {
230 double *data =
231 &(static_cast<FatPrismElementForcesAndSourcesCore *>(ptrFE)->normal[0]);
232 return VectorAdaptor(3, ublas::shallow_array_adaptor<double>(3, data));
233}
234
235inline VectorAdaptor
236FatPrismElementForcesAndSourcesCore::UserDataOperator::getNormalF4() {
237 double *data =
238 &(static_cast<FatPrismElementForcesAndSourcesCore *>(ptrFE)->normal[3]);
239 return VectorAdaptor(3, ublas::shallow_array_adaptor<double>(3, data));
240}
241
242inline MatrixDouble &
243FatPrismElementForcesAndSourcesCore::UserDataOperator::getGaussPts() {
244 return static_cast<FatPrismElementForcesAndSourcesCore *>(ptrFE)->gaussPts;
245}
246
247inline MatrixDouble &FatPrismElementForcesAndSourcesCore::UserDataOperator::
248 getGaussPtsTrianglesOnly() {
249 return static_cast<FatPrismElementForcesAndSourcesCore *>(ptrFE)
251}
252
253inline MatrixDouble &FatPrismElementForcesAndSourcesCore::UserDataOperator::
254 getGaussPtsThroughThickness() {
255 return static_cast<FatPrismElementForcesAndSourcesCore *>(ptrFE)
257}
258
259inline MatrixDouble &
260FatPrismElementForcesAndSourcesCore::UserDataOperator::getCoordsAtGaussPts() {
261 return static_cast<FatPrismElementForcesAndSourcesCore *>(ptrFE)
263}
264
265inline MatrixDouble &FatPrismElementForcesAndSourcesCore::UserDataOperator::
266 getCoordsAtGaussPtsTrianglesOnly() {
267 return static_cast<FatPrismElementForcesAndSourcesCore *>(ptrFE)
269}
270
271inline MatrixDouble &FatPrismElementForcesAndSourcesCore::UserDataOperator::
272 getHOCoordsAtGaussPtsF3() {
273 return static_cast<FatPrismElementForcesAndSourcesCore *>(ptrFE)
275}
276
277inline MatrixDouble &FatPrismElementForcesAndSourcesCore::UserDataOperator::
278 getHOCoordsAtGaussPtsF4() {
279 return static_cast<FatPrismElementForcesAndSourcesCore *>(ptrFE)
281}
282
283inline MatrixDouble &
284FatPrismElementForcesAndSourcesCore::UserDataOperator::getNormalsAtGaussPtsF3() {
285 return static_cast<FatPrismElementForcesAndSourcesCore *>(ptrFE)
287}
288
289inline MatrixDouble &
290FatPrismElementForcesAndSourcesCore::UserDataOperator::getNormalsAtGaussPtsF4() {
291 return static_cast<FatPrismElementForcesAndSourcesCore *>(ptrFE)
293}
294
295inline ublas::matrix_row<MatrixDouble>
296FatPrismElementForcesAndSourcesCore::UserDataOperator::getNormalsAtGaussPtsF3(
297 const int gg) {
298 return ublas::matrix_row<MatrixDouble>(
299 static_cast<FatPrismElementForcesAndSourcesCore *>(ptrFE)
301 gg);
302}
303
304inline ublas::matrix_row<MatrixDouble>
305FatPrismElementForcesAndSourcesCore::UserDataOperator::getNormalsAtGaussPtsF4(
306 const int gg) {
307 return ublas::matrix_row<MatrixDouble>(
308 static_cast<FatPrismElementForcesAndSourcesCore *>(ptrFE)
310 gg);
311}
312
313inline MatrixDouble &FatPrismElementForcesAndSourcesCore::UserDataOperator::
314 getTangent1AtGaussPtF3() {
315 return static_cast<FatPrismElementForcesAndSourcesCore *>(ptrFE)
317}
318
319inline MatrixDouble &FatPrismElementForcesAndSourcesCore::UserDataOperator::
320 getTangent2AtGaussPtF3() {
321 return static_cast<FatPrismElementForcesAndSourcesCore *>(ptrFE)
323}
324
325inline MatrixDouble &FatPrismElementForcesAndSourcesCore::UserDataOperator::
326 getTangent1AtGaussPtF4() {
327 return static_cast<FatPrismElementForcesAndSourcesCore *>(ptrFE)
329}
330
331inline MatrixDouble &FatPrismElementForcesAndSourcesCore::UserDataOperator::
332 getTangent2AtGaussPtF4() {
333 return static_cast<FatPrismElementForcesAndSourcesCore *>(ptrFE)
335}
336
337inline EntitiesFieldData &FatPrismElementForcesAndSourcesCore::
338 UserDataOperator::getTrianglesOnlyDataStructure() {
339 return static_cast<FatPrismElementForcesAndSourcesCore *>(ptrFE)
341}
342
343inline EntitiesFieldData &FatPrismElementForcesAndSourcesCore::
344 UserDataOperator::getTroughThicknessDataStructure() {
345 return static_cast<FatPrismElementForcesAndSourcesCore *>(ptrFE)
347}
348
350FatPrismElementForcesAndSourcesCore::UserDataOperator::getPrismFE() {
351 return static_cast<FatPrismElementForcesAndSourcesCore *>(ptrFE);
352}
353
354} // namespace MoFEM
355
356#endif //__FATPRISMELEMENTFORCESANDSURCESCORE_HPP__
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
@ MOFEM_NOT_IMPLEMENTED
Definition definitions.h:32
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
constexpr int order
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
VectorShallowArrayAdaptor< double > VectorAdaptor
Definition Types.hpp:115
implementation of Data Operators for Forces and Sources
Definition Common.hpp:10
Deprecated interface functions.
data structure for finite element entity
MatrixDouble & getTangent2AtGaussPtF3()
if higher order geometry return tangent vector to triangle at Gauss pts.
MatrixDouble & getHOCoordsAtGaussPtsF4()
coordinate at Gauss points on face 4 (if hierarchical approximation of element geometry)
MatrixDouble & getTangent2AtGaussPtF4()
if higher order geometry return tangent vector to triangle at Gauss pts.
MatrixDouble & getNormalsAtGaussPtsF4()
if higher order geometry return normals at face F4 at Gauss pts.
MatrixDouble & getNormalsAtGaussPtsF3()
if higher order geometry return normals at face F3 at Gauss pts.
const FatPrismElementForcesAndSourcesCore * getPrismFE()
return pointer to fat prism finite element
MatrixDouble & getHOCoordsAtGaussPtsF3()
coordinate at Gauss points on face 3 (if hierarchical approximation of element geometry)
MatrixDouble & getTangent1AtGaussPtF4()
if higher order geometry return tangent vector to triangle at Gauss pts.
MatrixDouble & getTangent1AtGaussPtF3()
if higher order geometry return tangent vector to triangle at Gauss pts.
virtual MoFEMErrorCode setGaussPtsTrianglesOnly(int order_triangles_only)
virtual MoFEMErrorCode setGaussPtsThroughThickness(int order_thickness)
MoFEMErrorCode operator()()
function is run for every finite element
structure to get information form mofem into EntitiesFieldData
MatrixDouble coordsAtGaussPts
coordinated at gauss points
MatrixDouble gaussPts
Matrix of integration points.
calculate normals at Gauss points of triangle element