v0.15.0
Loading...
Searching...
No Matches
VolumeElementForcesAndSourcesCore.hpp
Go to the documentation of this file.
1/** \file VolumeElementForcesAndSourcesCore.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_HPP__
12#define __VOLUMEELEMENTFORCESANDSOURCESCORE_HPP__
13
14using namespace boost::numeric;
15
16namespace MoFEM {
17
18/** \brief Volume finite element base
19 \ingroup mofem_forces_and_sources_volume_element
20
21 User is implementing own operator at Gauss point level, by class
22 derived from VolumeElementForcesAndSourcesCore::UserDataOperator. Arbitrary
23 number of operator can be added by pushing objects to OpPtrVector
24
25 */
27
29
31 const EntityType type);
32
33 std::string meshPositionsFieldName; ///< \deprecated DO NOT USE!
34
35 /** \brief default operator for TET element
36 * \ingroup mofem_forces_and_sources_volume_element
37 */
38 struct UserDataOperator;
39
40 enum Switches {
41 };
42
44
45protected:
46
47 // Note that functions below could be overloaded by user to change default
48 // behavior of the element.
49
50 /**
51 * \brief Set integration points
52 * @return Error code
53 */
55
56 /**
57 * \brief Calculate element volume and Jacobian
58 *
59 * Note that at that point is assumed that geometry is exclusively defined by
60 * corner nodes.
61 *
62 * @return Error code
63 */
65
66 /**
67 * \brief Calculate coordinate at integration points
68 * @return Error code
69 */
71
72 /**
73 * \brief Determine approximation space and order of base functions
74 * @return Error code
75 */
77
78 /**
79 * \brief Transform base functions based on geometric element Jacobian.
80 *
81 * This function apply transformation to base functions and its derivatives.
82 * For example when base functions for H-div are present the
83 * Piola-Transformation is applied to base functions and their derivatives.
84 *
85 * @return Error code
86 */
88
89
93
98
99 double &vOlume;
100
105
106 friend class UserDataOperator;
107};
108
111
112 using ForcesAndSourcesCore::UserDataOperator::UserDataOperator;
113
114 /** \brief get element number of nodes
115 */
116 inline int getNumNodes();
117
118 /** \brief get element connectivity
119 */
120 inline const EntityHandle *getConn();
121
122 /** \brief element volume (linear geometry)
123 */
124 inline double getVolume() const;
125
126 /** \brief element volume (linear geometry)
127 */
128 inline double &getVolume();
129
130 /**
131 * \brief get element Jacobian
132 */
134
135 /**
136 * \brief get element inverse Jacobian
137 */
139
140 /** \brief nodal coordinates
141 */
142 inline VectorDouble &getCoords();
143
144 /** \brief return pointer to Generic Volume Finite Element object
145 */
147
148protected:
150};
151
152int VolumeElementForcesAndSourcesCore::UserDataOperator::getNumNodes() {
153 return static_cast<VolumeElementForcesAndSourcesCore *>(ptrFE)->num_nodes;
154}
155
156const EntityHandle *
157VolumeElementForcesAndSourcesCore::UserDataOperator::getConn() {
158 return static_cast<VolumeElementForcesAndSourcesCore *>(ptrFE)->conn;
159}
160
161double
162VolumeElementForcesAndSourcesCore::UserDataOperator::getVolume() const {
163 return static_cast<VolumeElementForcesAndSourcesCore *>(ptrFE)->vOlume;
164}
165
166double &VolumeElementForcesAndSourcesCore::UserDataOperator::getVolume() {
167 return static_cast<VolumeElementForcesAndSourcesCore *>(ptrFE)->vOlume;
168}
169
171VolumeElementForcesAndSourcesCore::UserDataOperator::getJac() {
172 return static_cast<VolumeElementForcesAndSourcesCore *>(ptrFE)->tJac;
173}
174
176VolumeElementForcesAndSourcesCore::UserDataOperator::getInvJac() {
177 return static_cast<VolumeElementForcesAndSourcesCore *>(ptrFE)->tInvJac;
178}
179
181VolumeElementForcesAndSourcesCore::UserDataOperator::getCoords() {
182 return static_cast<VolumeElementForcesAndSourcesCore *>(ptrFE)->coords;
183}
184
186VolumeElementForcesAndSourcesCore::UserDataOperator::getVolumeFE() const {
187 return static_cast<VolumeElementForcesAndSourcesCore *>(ptrFE);
188}
189
190/**
191 * @deprecated use VolumeElementForcesAndSourcesCore
192 *
193 */
196
197/**
198 * @deprecated obsolete template. it will be removed in next versions
199 */
200template <int SWITCH>
203
204 using VolumeElementForcesAndSourcesCore::VolumeElementForcesAndSourcesCore;
206};
207
208} // namespace MoFEM
209
210#endif //__VOLUMEELEMENTFORCESANDSOURCESCORE_HPP__
211
212/**
213 * \defgroup mofem_forces_and_sources_volume_element Volume Element
214 * \brief Implementation of general volume element.
215 *
216 * \ingroup mofem_forces_and_sources
217 **/
#define DEPRECATED
Definition definitions.h:17
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
implementation of Data Operators for Forces and Sources
Definition Common.hpp:10
DEPRECATED typedef VolumeElementForcesAndSourcesCore VolumeElementForcesAndSourcesCoreBase
Deprecated interface functions.
structure to get information form mofem into EntitiesFieldData
apply contravariant (Piola) transfer to Hdiv space
apply covariant transfer to Hcurl space
Transform local reference derivatives of shape function to global derivatives.
brief Transform local reference derivatives of shape function to global derivatives
VolumeElementForcesAndSourcesCore * getVolumeFE() const
return pointer to Generic Volume Finite Element object
FTensor::Tensor2< double *, 3, 3 > & getInvJac()
get element inverse Jacobian
FTensor::Tensor2< double *, 3, 3 > & getJac()
get element Jacobian
virtual MoFEMErrorCode transformBaseFunctions()
Transform base functions based on geometric element Jacobian.
virtual MoFEMErrorCode setIntegrationPts()
Set integration points.
MoFEMErrorCode operator()()
function is run for every finite element
virtual MoFEMErrorCode getSpaceBaseAndOrderOnElement()
Determine approximation space and order of base functions.
virtual MoFEMErrorCode calculateVolumeAndJacobian()
Calculate element volume and Jacobian.
virtual MoFEMErrorCode calculateCoordinatesAtGaussPts()
Calculate coordinate at integration points.