v0.13.2
Loading...
Searching...
No Matches
VolumeElementForcesAndSourcesCoreOnSide.hpp
Go to the documentation of this file.
1/** \file VolumeElementForcesAndSourcesCoreOnSide.hpp
2 \brief Volume element.
3
4 Those element are inherited by user to implement specific implementation of
5 particular problem.
6
7*/
8
9
10
11#ifndef __VOLUMEELEMENTFORCESANDSOURCESCORE_ONSIDE_HPP__
12#define __VOLUMEELEMENTFORCESANDSOURCESCORE_ONSIDE_HPP__
13
14using namespace boost::numeric;
15
16namespace MoFEM {
17
18/**
19 * \brief Base volume element used to integrate on skeleton
20 * \ingroup mofem_forces_and_sources_volume_element
21 */
24
26
27 /**
28 * @brief Get the face nodes mapped on volume element
29 *
30 * \todo That this is not general, e.g., for quad number of nodes is 4.
31 *
32 * @return const std::array<int, 4>&
33 */
34 inline const std::array<int, 4> &getFaceConnMap() const;
35
36 /**
37 * @brief Get face nodes maped on volume
38 *
39 * \todo That this is not general, e.g., for prism or hex, size of fixed array
40 * is wrong.
41 *
42 * @return const sdt::array<int, 4>&
43 */
44 inline const std::array<int, 8> &getTetConnMap() const;
45
46 /**
47 * @brief Get node on volume opposite to volume element
48 *
49 * \todo That this is not general, e.g., for prism or hex, opoosite node is
50 * not unique.
51 *
52 * @return int
53 */
54 inline int getOppositeNode() const;
55
56 /**
57 * @brief Sense face on volume
58 * @deprecated use getSkeletonSense()
59 *
60 * @return int
61 */
62 DEPRECATED inline int getFaceSense() const;
63
64 /**
65 * @brief Sense face on volume
66 *
67 * @return int
68 */
69 inline int getSkeletonSense() const;
70
71 /**
72 * @brief Face number on the volume
73 *
74 * @return int
75 */
76 inline int getFaceSideNumber() const;
77
78 /** \brief default operator for TET element
79 * \ingroup mofem_forces_and_sources_volume_element
80 */
81 struct UserDataOperator;
82
83 int getRule(int order);
85
86private:
87 int faceSense; ///< Sense of face, could be 1 or -1
88 int faceSideNumber; ///< Face side number
89 std::array<int, 4> faceConnMap;
90 std::array<int, 8> tetConnMap;
92};
93
94/** \brief default operator for TET element
95 * \ingroup mofem_forces_and_sources_volume_element
96 */
99
100 using VolumeElementForcesAndSourcesCore::UserDataOperator::UserDataOperator;
101
103
105
106 /**
107 * \brief get face sense in respect to volume
108 * \deprecated use getSkeletonSense
109 * @return error code
110 */
111 DEPRECATED inline int getFaceSense() const;
112
113 /**
114 * @brief Get the skeleton sense
115 *
116 * calls getFaceSense()
117 *
118 * @return int
119 */
120 inline int getSkeletonSense() const;
121
122 /**
123 * \brief get face side number in respect to volume
124 * @return error code
125 */
126 inline int getFaceSideNumber() const;
127
128 inline bool getEdgeFace(const int ee) const;
129
130 /**
131 * get face normal on side which is this element
132 * @return face normal
133 */
134 inline VectorDouble &getNormal();
135
136 /** \brief get triangle tangent 1
137 */
138 inline VectorDouble &getTangent1();
139
140 /** \brief get triangle tangent 2
141 */
142 inline VectorDouble &getTangent2();
143
144 /** \brief get normal as tensor
145 */
146 inline auto getFTensor1Normal();
147
148 /** \brief get tangentOne as tensor
149 */
150 inline auto getFTensor1Tangent1();
151
152 /** \brief get tangentTwo as tensor
153 */
154 inline auto getFTensor1Tangent2();
155
156 /** \brief if higher order geometry return normals at Gauss pts.
157
158 Note: returned matrix has size 0 in rows and columns if no HO approximation
159 of geometry is available.
160
161 */
163
164 /** \brief if higher order geometry return normals at Gauss pts.
165 *
166 * \param gg gauss point number
167 */
168 inline ublas::matrix_row<MatrixDouble> getNormalsAtGaussPts(const int gg);
169
170 /** \brief get normal at integration points
171
172 Example:
173 \code
174 double nrm2;
175 FTensor::Index<'i',3> i;
176 auto t_normal = getFTensor1NormalsAtGaussPts();
177 for(int gg = gg!=data.getN().size1();gg++) {
178 nrm2 = sqrt(t_normal(i)*t_normal(i));
179 ++t_normal;
180 }
181 \endcode
182
183 */
184 inline auto getFTensor1NormalsAtGaussPts();
185
186 /** \brief get face coordinates at Gauss pts.
187
188 \note Coordinates should be the same what function getCoordsAtGaussPts
189 on tets is returning. If both coordinates are different it is error, or you
190 do something very unusual.
191
192 */
194
195protected:
197};
198
199const std::array<int, 4> &
201 return faceConnMap;
202}
203
204const std::array<int, 8> &
206 return tetConnMap;
207}
208
210 return oppositeNode;
211}
212
214 return getSkeletonSense();
215}
216
218 return faceSense;
219}
220
222 return faceSideNumber;
223}
224
226VolumeElementForcesAndSourcesCoreOnSide::UserDataOperator::getFaceFE()
227 const {
228 return static_cast<FaceElementForcesAndSourcesCore *>(
230}
231
233VolumeElementForcesAndSourcesCoreOnSide::UserDataOperator::getNormal() {
234 return getFaceFE()->nOrmal;
235}
236
238VolumeElementForcesAndSourcesCoreOnSide::UserDataOperator::getTangent1() {
239 return getFaceFE()->tangentOne;
240}
241
243VolumeElementForcesAndSourcesCoreOnSide::UserDataOperator::getTangent2() {
244 return getFaceFE()->tangentTwo;
245}
246
248VolumeElementForcesAndSourcesCoreOnSide::UserDataOperator::getVolumeFE()
249 const {
250 return static_cast<VolumeElementForcesAndSourcesCoreOnSide *>(ptrFE);
251}
252
253int VolumeElementForcesAndSourcesCoreOnSide::UserDataOperator::
254 getFaceSense() const {
255 return getSkeletonSense();
256}
257
258int VolumeElementForcesAndSourcesCoreOnSide::UserDataOperator::
259 getSkeletonSense() const {
260 return getVolumeFE()->faceSense;
261}
262
263auto VolumeElementForcesAndSourcesCoreOnSide::UserDataOperator::
264 getFTensor1Normal() {
265 double *ptr = &*getNormal().data().begin();
266 return FTensor::Tensor1<double *, 3>(ptr, &ptr[1], &ptr[2]);
267}
268
269auto VolumeElementForcesAndSourcesCoreOnSide::UserDataOperator::
270 getFTensor1Tangent1() {
271 double *ptr = &*getTangent1().data().begin();
272 return FTensor::Tensor1<double *, 3>(ptr, &ptr[1], &ptr[2]);
273}
274
275auto VolumeElementForcesAndSourcesCoreOnSide::UserDataOperator::
276 getFTensor1Tangent2() {
277 double *ptr = &*getTangent2().data().begin();
278 return FTensor::Tensor1<double *, 3>(ptr, &ptr[1], &ptr[2]);
279}
280
281inline auto VolumeElementForcesAndSourcesCoreOnSide::UserDataOperator::
282 getFTensor1NormalsAtGaussPts() {
283 double *ptr = &*getNormalsAtGaussPts().data().begin();
285 &ptr[2]);
286}
287
288int VolumeElementForcesAndSourcesCoreOnSide::UserDataOperator::
289 getFaceSideNumber() const {
290 return getVolumeFE()->faceSideNumber;
291}
292
293MatrixDouble &VolumeElementForcesAndSourcesCoreOnSide::UserDataOperator::
294 getNormalsAtGaussPts() {
295 return getFaceFE()->normalsAtGaussPts;
296}
297
298MatrixDouble &VolumeElementForcesAndSourcesCoreOnSide::UserDataOperator::
299 getFaceCoordsAtGaussPts() {
300 return getFaceFE()->coordsAtGaussPts;
301}
302
303bool VolumeElementForcesAndSourcesCoreOnSide::UserDataOperator::getEdgeFace(
304 const int ee) const {
305 constexpr bool edges_on_faces[6][4] = {{true, false, false, true}, // e0
306 {false, true, false, true}, // e1
307 {false, false, true, true}, // e2
308 {true, false, true, false}, // e3
309 {true, true, false, false}, // e4
310 {false, true, true, false}};
311 return edges_on_faces[ee][getFaceSideNumber()];
312}
313
314ublas::matrix_row<MatrixDouble> VolumeElementForcesAndSourcesCoreOnSide::
315 UserDataOperator::getNormalsAtGaussPts(const int gg) {
316 return ublas::matrix_row<MatrixDouble>(getNormalsAtGaussPts(), gg);
317}
318
319/**
320 * @deprecated use VolumeElementForcesAndSourcesCore
321 *
322 */
325
326/**
327 * @deprecated do not use this template, it is obsolete
328 */
329template <int SWITCH>
333 VolumeElementForcesAndSourcesCoreOnSide;
336};
337
338
339
340} // namespace MoFEM
341
342#endif //__VOLUMEELEMENTFORCESANDSOURCESCORE_ONSIDE_HPP__
#define DEPRECATED
Definition: definitions.h:17
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
DEPRECATED typedef VolumeElementForcesAndSourcesCoreOnSide VolumeElementForcesAndSourcesCoreOnSideBase
structure to get information form mofem into EntitiesFieldData
ForcesAndSourcesCore * sidePtrFE
Element to integrate on the sides.
VolumeElementForcesAndSourcesCore(Interface &m_field, const EntityType type=MBTET)
MatrixDouble & getNormalsAtGaussPts()
if higher order geometry return normals at Gauss pts.
Base volume element used to integrate on skeleton.
const std::array< int, 8 > & getTetConnMap() const
Get face nodes maped on volume.
const std::array< int, 4 > & getFaceConnMap() const
Get the face nodes mapped on volume element.
int getOppositeNode() const
Get node on volume opposite to volume element.