v0.14.0
Loading...
Searching...
No Matches
EdgeElementForcesAndSourcesCore.hpp
Go to the documentation of this file.
1/** \file EdgeElementForcesAndSourcesCore.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 __EDGEELEMENTFORCESANDSURCESCORE_HPP__
13#define __EDGEELEMENTFORCESANDSURCESCORE_HPP__
14
15using namespace boost::numeric;
16
17namespace MoFEM {
18
19struct FaceElementForcesAndSourcesCoreOnSide;
20
21/** \brief Edge finite element
22 * \ingroup mofem_forces_and_sources_edge_element
23 *
24 * User is implementing own operator at Gauss points level, by own object
25 * derived from EdgeElementForcesAndSourcesCoreL::UserDataOperator. Arbitrary
26 * number of operator added pushing objects to rowOpPtrVector and
27 * rowColOpPtrVector.
28 *
29 */
64
65/** \brief default operator for EDGE element
66 \ingroup mofem_forces_and_sources_edge_element
67 */
70
71 using ForcesAndSourcesCore::UserDataOperator::UserDataOperator;
72
73 /** \brief get element connectivity
74 */
75 inline const EntityHandle *getConn();
76
77 /**
78 * \brief get edge length
79 */
80 inline double getLength();
81
82 /**
83 * \brief get edge direction
84 */
85 inline VectorDouble &getDirection();
86
87 /**
88 * \brief get edge normal
89 * NOTE: it should be used only in 2D analysis
90 */
91 inline auto getFTensor1Normal();
92
93 /**
94 * @brief get ftensor1 edge normal
95 *
96 * @param vec vector in third direction
97 * @return auto
98 */
99 inline auto getFTensor1Normal(const FTensor::Tensor1<double, 3> &vec);
100
101 inline auto
103
104 /** \brief get normal at integration points
105
106 Example:
107 \code
108 double nrm2;
109 FTensor::Index<'i',3> i;
110 auto t_normal = getFTensor1NormalsAtGaussPts();
111 for(int gg = gg!=data.getN().size1();gg++) {
112 nrm2 = sqrt(t_normal(i)*t_normal(i));
113 ++t_normal;
114 }
115 \endcode
116
117 */
118 inline auto getFTensor1NormalsAtGaussPts();
119
120 /**
121 * \brief get edge node coordinates
122 */
123 inline VectorDouble &getCoords();
124
125 /**
126 * \brief get tangent vector to edge curve at integration points
127 */
129
133
134 /**
135 * \brief get pointer to this finite element
136 */
138
140
141 /**
142 * \brief get get coords at gauss points
143
144 \code
145 FTensor::Index<'i',3> i;
146 auto t_center;
147 auto t_coords = getFTensor1Coords();
148 t_center(i) = 0;
149 for(int nn = 0;nn!=2;nn++) {
150 t_center(i) += t_coords(i);
151 ++t_coords;
152 }
153 t_center(i) /= 2;
154 \endcode
155
156 */
158
159 /**
160 * @deprecated use getFTensor1Coords
161 */
164 return getFTensor1Coords();
165 }
166
167 /**
168 * @brief Get tangent vector at integration points
169 *
170 * @return FTensor::Tensor1<FTensor::PackPtr<double *, 3>, DIM>
171 */
172 template <int DIM = 3>
175
177 loopSideFaces(const string fe_name,
179
180protected:
182};
183
184const EntityHandle *
185EdgeElementForcesAndSourcesCore::UserDataOperator::getConn() {
186 return static_cast<EdgeElementForcesAndSourcesCore *>(ptrFE)->cOnn;
187}
188
189double EdgeElementForcesAndSourcesCore::UserDataOperator::getLength() {
190 return static_cast<EdgeElementForcesAndSourcesCore *>(ptrFE)->lEngth;
191}
192
194EdgeElementForcesAndSourcesCore::UserDataOperator::getDirection() {
195 return static_cast<EdgeElementForcesAndSourcesCore *>(ptrFE)->dIrection;
196}
197
199EdgeElementForcesAndSourcesCore::UserDataOperator::getCoords() {
200 return static_cast<EdgeElementForcesAndSourcesCore *>(ptrFE)->cOords;
201}
202
204EdgeElementForcesAndSourcesCore::UserDataOperator::getTangentAtGaussPts() {
205 return static_cast<EdgeElementForcesAndSourcesCore *>(ptrFE)
207}
208
210EdgeElementForcesAndSourcesCore::UserDataOperator::getEdgeFE() {
211 return static_cast<EdgeElementForcesAndSourcesCore *>(ptrFE);
212}
213
215EdgeElementForcesAndSourcesCore::UserDataOperator::getFTensor1Direction() {
216 return FTensor::Tensor1<double, 3>(getDirection()[0], getDirection()[1],
217 getDirection()[2]);
218}
219
221EdgeElementForcesAndSourcesCore::UserDataOperator::getFTensor1Coords() {
222 double *ptr = &*getCoords().data().begin();
224 &ptr[2]);
225}
226
227template <>
229EdgeElementForcesAndSourcesCore::UserDataOperator::getFTensor1TangentAtGaussPts<
230 3>() {
231 double *ptr = &*getTangentAtGaussPts().data().begin();
233 &ptr[2]);
234}
235
236template <>
238EdgeElementForcesAndSourcesCore::UserDataOperator::getFTensor1TangentAtGaussPts<
239 2>() {
240 double *ptr = &*getTangentAtGaussPts().data().begin();
241 return FTensor::Tensor1<FTensor::PackPtr<double *, 3>, 2>(ptr, &ptr[1]);
242}
243
244auto EdgeElementForcesAndSourcesCore::UserDataOperator::getFTensor1Normal(
245 const FTensor::Tensor1<double, 3> &vec) {
246 FTensor::Index<'i', 3> i;
247 FTensor::Index<'j', 3> j;
248 FTensor::Index<'k', 3> k;
250 auto t_dir = getFTensor1Direction();
251 t_normal(i) = FTensor::levi_civita(i, j, k) * t_dir(j) * vec(k);
252 return t_normal;
253}
254
255auto EdgeElementForcesAndSourcesCore::UserDataOperator::
256 getFTensor1Normal() {
257 return getFTensor1Normal(tFaceOrientation);
258}
259
260auto EdgeElementForcesAndSourcesCore::UserDataOperator::
261 getFTensor1NormalsAtGaussPts(const FTensor::Tensor1<double, 3> &vec) {
262 FTensor::Index<'i', 3> i;
263 FTensor::Index<'j', 3> j;
264 FTensor::Index<'k', 3> k;
265
266 auto &normals_at_gauss_pts =
268 normals_at_gauss_pts.resize(3, getGaussPts().size2(), false);
269 auto t_normal = getFTensor1FromMat<3>(normals_at_gauss_pts);
270 auto t_dir = getFTensor1TangentAtGaussPts<3>();
271 for(auto gg = 0; gg!=getGaussPts().size2(); ++gg) {
272 t_normal(i) = FTensor::levi_civita(i, j, k) * t_dir(j) * vec(k);
273 ++t_normal;
274 ++t_dir;
275 }
276
277 return getFTensor1FromMat<3>(normals_at_gauss_pts);
278}
279
280auto EdgeElementForcesAndSourcesCore::UserDataOperator::
281 getFTensor1NormalsAtGaussPts() {
282 return getFTensor1NormalsAtGaussPts(tFaceOrientation);
283}
284
285/**
286 * @deprecated use EdgeElementForcesAndSourcesCore
287 */
290
291} // namespace MoFEM
292
293#endif //__EDGEELEMENTFORCESANDSURCESCORE_HPP__
294
295/**
296 * \defgroup mofem_forces_and_sources_edge_element Edge Element
297 *
298 * \brief Implementation of edge element.
299 *
300 * \ingroup mofem_forces_and_sources
301 */
T data[Tensor_Dim]
#define DEPRECATED
Definition definitions.h:17
FTensor::Index< 'i', SPACE_DIM > i
FTensor::Index< 'j', 3 > j
FTensor::Index< 'k', 3 > k
constexpr std::enable_if<(Dim0<=2 &&Dim1<=2), Tensor2_Expr< Levi_Civita< T >, T, Dim0, Dim1, i, j > >::typ levi_civita)(const Index< i, Dim0 > &, const Index< j, Dim1 > &)
levi_civita functions to make for easy adhoc use
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
implementation of Data Operators for Forces and Sources
Definition Common.hpp:10
DEPRECATED typedef EdgeElementForcesAndSourcesCore EdgeElementForcesAndSourcesCoreBase
Deprecated interface functions.
DEPRECATED FTensor::Tensor1< FTensor::PackPtr< double *, 3 >, 3 > getTensor1Coords()
auto getFTensor1Normal()
get edge normal NOTE: it should be used only in 2D analysis
MatrixDouble & getTangentAtGaussPts()
get tangent vector to edge curve at integration points
MoFEMErrorCode loopSideFaces(const string fe_name, FaceElementForcesAndSourcesCoreOnSide &fe_side)
const EdgeElementForcesAndSourcesCore * getEdgeFE()
get pointer to this finite element
FTensor::Tensor1< FTensor::PackPtr< double *, 3 >, DIM > getFTensor1TangentAtGaussPts()
Get tangent vector at integration points.
FTensor::Tensor1< FTensor::PackPtr< double *, 3 >, 3 > getFTensor1Coords()
get get coords at gauss points
static FTensor::Tensor1< double, 3 > tFaceOrientation
MoFEMErrorCode operator()()
function is run for every finite element
Base face element used to integrate on skeleton.
structure to get information form mofem into EntitiesFieldData