v0.14.0
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 
14 using namespace boost::numeric;
15 
16 namespace 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  const EntityType type = MBTET);
30 
31 
32  std::string meshPositionsFieldName; ///< \deprecated DO NOT USE!
33 
34  /** \brief default operator for TET element
35  * \ingroup mofem_forces_and_sources_volume_element
36  */
37  struct UserDataOperator;
38 
39  enum Switches {
40  };
41 
42  MoFEMErrorCode operator()();
43 
44 protected:
45 
46  // Note that functions below could be overloaded by user to change default
47  // behavior of the element.
48 
49  /**
50  * \brief Set integration points
51  * @return Error code
52  */
53  virtual MoFEMErrorCode setIntegrationPts();
54 
55  /**
56  * \brief Calculate element volume and Jacobian
57  *
58  * Note that at that point is assumed that geometry is exclusively defined by
59  * corner nodes.
60  *
61  * @return Error code
62  */
63  virtual MoFEMErrorCode calculateVolumeAndJacobian();
64 
65  /**
66  * \brief Calculate coordinate at integration points
67  * @return Error code
68  */
69  virtual MoFEMErrorCode calculateCoordinatesAtGaussPts();
70 
71  /**
72  * \brief Determine approximation space and order of base functions
73  * @return Error code
74  */
75  virtual MoFEMErrorCode getSpaceBaseAndOrderOnElement();
76 
77  /**
78  * \brief Transform base functions based on geometric element Jacobian.
79  *
80  * This function apply transformation to base functions and its derivatives.
81  * For example when base functions for H-div are present the
82  * Piola-Transformarion is applied to base functions and their derivatives.
83  *
84  * @return Error code
85  */
86  virtual MoFEMErrorCode transformBaseFunctions();
87 
88 
92 
97 
98  double &vOlume;
99 
104 
105  friend class UserDataOperator;
106 };
107 
110 
112 
113  /** \brief get element number of nodes
114  */
115  inline int getNumNodes();
116 
117  /** \brief get element connectivity
118  */
119  inline const EntityHandle *getConn();
120 
121  /** \brief element volume (linear geometry)
122  */
123  inline double getVolume() const;
124 
125  /** \brief element volume (linear geometry)
126  */
127  inline double &getVolume();
128 
129  /**
130  * \brief get element Jacobian
131  */
132  inline FTensor::Tensor2<double *, 3, 3> &getJac();
133 
134  /**
135  * \brief get element inverse Jacobian
136  */
137  inline FTensor::Tensor2<double *, 3, 3> &getInvJac();
138 
139  /** \brief nodal coordinates
140  */
141  inline VectorDouble &getCoords();
142 
143  /** \brief return pointer to Generic Volume Finite Element object
144  */
145  inline VolumeElementForcesAndSourcesCore *getVolumeFE() const;
146 
147 protected:
148  MoFEMErrorCode setPtrFE(ForcesAndSourcesCore *ptr);
149 };
150 
151 int VolumeElementForcesAndSourcesCore::UserDataOperator::getNumNodes() {
152  return static_cast<VolumeElementForcesAndSourcesCore *>(ptrFE)->num_nodes;
153 }
154 
155 const EntityHandle *
156 VolumeElementForcesAndSourcesCore::UserDataOperator::getConn() {
157  return static_cast<VolumeElementForcesAndSourcesCore *>(ptrFE)->conn;
158 }
159 
160 double
161 VolumeElementForcesAndSourcesCore::UserDataOperator::getVolume() const {
162  return static_cast<VolumeElementForcesAndSourcesCore *>(ptrFE)->vOlume;
163 }
164 
165 double &VolumeElementForcesAndSourcesCore::UserDataOperator::getVolume() {
166  return static_cast<VolumeElementForcesAndSourcesCore *>(ptrFE)->vOlume;
167 }
168 
170 VolumeElementForcesAndSourcesCore::UserDataOperator::getJac() {
171  return static_cast<VolumeElementForcesAndSourcesCore *>(ptrFE)->tJac;
172 }
173 
175 VolumeElementForcesAndSourcesCore::UserDataOperator::getInvJac() {
176  return static_cast<VolumeElementForcesAndSourcesCore *>(ptrFE)->tInvJac;
177 }
178 
179 VectorDouble &
180 VolumeElementForcesAndSourcesCore::UserDataOperator::getCoords() {
181  return static_cast<VolumeElementForcesAndSourcesCore *>(ptrFE)->coords;
182 }
183 
185 VolumeElementForcesAndSourcesCore::UserDataOperator::getVolumeFE() const {
186  return static_cast<VolumeElementForcesAndSourcesCore *>(ptrFE);
187 }
188 
189 /**
190  * @deprecated use VolumeElementForcesAndSourcesCore
191  *
192  */
195 
196 /**
197  * @deprecated obsolete template. it will be removed in next versions
198  */
199 template <int SWITCH>
202 
203  using VolumeElementForcesAndSourcesCore::VolumeElementForcesAndSourcesCore;
205 };
206 
207 } // namespace MoFEM
208 
209 #endif //__VOLUMEELEMENTFORCESANDSOURCESCORE_HPP__
210 
211 /**
212  * \defgroup mofem_forces_and_sources_volume_element Volume Element
213  * \brief Implementation of general volume element.
214  *
215  * \ingroup mofem_forces_and_sources
216  **/
MoFEM::VolumeElementForcesAndSourcesCore::meshPositionsFieldName
std::string meshPositionsFieldName
Definition: VolumeElementForcesAndSourcesCore.hpp:32
DEPRECATED
#define DEPRECATED
Definition: definitions.h:17
EntityHandle
MoFEM::VolumeElementForcesAndSourcesCoreSwitch
Definition: VolumeElementForcesAndSourcesCore.hpp:200
MoFEM::OpSetInvJacH1
Transform local reference derivatives of shape function to global derivatives.
Definition: DataOperators.hpp:131
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
MoFEM::VolumeElementForcesAndSourcesCore::opCovariantPiolaTransform
OpSetCovariantPiolaTransform opCovariantPiolaTransform
Definition: VolumeElementForcesAndSourcesCore.hpp:95
MoFEM::VolumeElementForcesAndSourcesCore::opSetInvJacHdivAndHcurl
OpSetInvJacHdivAndHcurl opSetInvJacHdivAndHcurl
Definition: VolumeElementForcesAndSourcesCore.hpp:96
FTensor::Tensor2< double *, 3, 3 >
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
MoFEM::VolumeElementForcesAndSourcesCore::coords
VectorDouble coords
Definition: VolumeElementForcesAndSourcesCore.hpp:89
MoFEM::OpSetCovariantPiolaTransform
apply covariant transfer to Hcurl space
Definition: DataOperators.hpp:218
MoFEM::VolumeElementForcesAndSourcesCore::num_nodes
int num_nodes
Definition: VolumeElementForcesAndSourcesCore.hpp:100
MoFEM
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
MoFEM::OpSetContravariantPiolaTransform
apply contravariant (Piola) transfer to Hdiv space
Definition: DataOperators.hpp:183
MoFEM::ForcesAndSourcesCore::UserDataOperator
Definition: ForcesAndSourcesCore.hpp:549
convert.type
type
Definition: convert.py:64
MoFEM::OpSetInvJacHdivAndHcurl
brief Transform local reference derivatives of shape function to global derivatives
Definition: DataOperators.hpp:153
MoFEM::VolumeElementForcesAndSourcesCore::conn
const EntityHandle * conn
Definition: VolumeElementForcesAndSourcesCore.hpp:101
MoFEM::VolumeElementForcesAndSourcesCore::opSetInvJacH1
OpSetInvJacH1 opSetInvJacH1
Definition: VolumeElementForcesAndSourcesCore.hpp:93
MoFEM::VolumeElementForcesAndSourcesCore::tJac
FTensor::Tensor2< double *, 3, 3 > tJac
Definition: VolumeElementForcesAndSourcesCore.hpp:102
MoFEM::VolumeElementForcesAndSourcesCore
Volume finite element base.
Definition: VolumeElementForcesAndSourcesCore.hpp:26
MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator
Definition: VolumeElementForcesAndSourcesCore.hpp:108
MoFEM::VolumeElementForcesAndSourcesCore::Switches
Switches
Definition: VolumeElementForcesAndSourcesCore.hpp:39
MoFEM::ForcesAndSourcesCore
structure to get information form mofem into EntitiesFieldData
Definition: ForcesAndSourcesCore.hpp:22
MoFEM::VolumeElementForcesAndSourcesCore::vOlume
double & vOlume
Definition: VolumeElementForcesAndSourcesCore.hpp:98
MoFEM::VolumeElementForcesAndSourcesCore::invJac
MatrixDouble3by3 invJac
Definition: VolumeElementForcesAndSourcesCore.hpp:91
MoFEM::VolumeElementForcesAndSourcesCore::jAc
MatrixDouble3by3 jAc
Definition: VolumeElementForcesAndSourcesCore.hpp:90
MoFEM::VolumeElementForcesAndSourcesCoreBase
DEPRECATED typedef VolumeElementForcesAndSourcesCore VolumeElementForcesAndSourcesCoreBase
Definition: VolumeElementForcesAndSourcesCore.hpp:194
UserDataOperator
ForcesAndSourcesCore::UserDataOperator UserDataOperator
Definition: HookeElement.hpp:75
MatrixBoundedArray< double, 9 >
UBlasVector< double >
MoFEM::VolumeElementForcesAndSourcesCore::opContravariantPiolaTransform
OpSetContravariantPiolaTransform opContravariantPiolaTransform
Definition: VolumeElementForcesAndSourcesCore.hpp:94
MoFEM::VolumeElementForcesAndSourcesCore::tInvJac
FTensor::Tensor2< double *, 3, 3 > tInvJac
Definition: VolumeElementForcesAndSourcesCore.hpp:103