v0.8.23
Public Member Functions | List of all members
MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator Struct Reference

default operator for TET element More...

#include <src/finite_elements/VolumeElementForcesAndSourcesCore.hpp>

Inheritance diagram for MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator:
[legend]
Collaboration diagram for MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator:
[legend]

Public Member Functions

 UserDataOperator (const FieldSpace space)
 
 UserDataOperator (const std::string &field_name, const char type)
 
 UserDataOperator (const std::string &row_field_name, const std::string &col_field_name, const char type, const bool symm=true)
 
int getNumNodes ()
 get element number of nodes More...
 
const EntityHandlegetConn ()
 get element connectivity More...
 
double getVolume ()
 element volume (linear geometry) More...
 
FTensor::Tensor2< double *, 3, 3 > & getJac ()
 get element Jacobian More...
 
FTensor::Tensor2< double *, 3, 3 > & getInvJac ()
 get element inverse Jacobian More...
 
double getMeasure ()
 get measure of element More...
 
VectorDoublegetCoords ()
 nodal coordinates More...
 
MatrixDoublegetCoordsAtGaussPts ()
 Gauss points and weight, matrix (nb. of points x 3) More...
 
MatrixDoublegetHoCoordsAtGaussPts ()
 coordinate at Gauss points (if hierarchical approximation of element geometry) More...
 
MatrixDoublegetHoGaussPtsJac ()
 
MatrixDoublegetHoGaussPtsInvJac ()
 
VectorDoublegetHoGaussPtsDetJac ()
 
auto getFTenosr0HoMeasure ()
 
auto getFTensor1CoordsAtGaussPts ()
 Get coordinates at integration points assuming linear geometry. More...
 
auto getFTensor1HoCoordsAtGaussPts ()
 Get coordinates at integration points taking geometry from field. More...
 
auto getFTensor2HoGaussPtsJac ()
 
auto getFTensor2HoGaussPtsInvJac ()
 
const VolumeElementForcesAndSourcesCoregetVolumeFE ()
 return pointer to Generic Volume Finite Element object More...
 
MoFEMErrorCode getDivergenceOfHDivBaseFunctions (int side, EntityType type, DataForcesAndSourcesCore::EntData &data, int gg, VectorDouble &div)
 Get divergence of base functions at integration point. More...
 
MoFEMErrorCode getCurlOfHCurlBaseFunctions (int side, EntityType type, DataForcesAndSourcesCore::EntData &data, int gg, MatrixDouble &curl)
 Get curl of base functions at integration point. More...
 
DEPRECATED auto getTensor1CoordsAtGaussPts ()
 
DEPRECATED auto getTensor1HoCoordsAtGaussPts ()
 
- Public Member Functions inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
 UserDataOperator (const FieldSpace space, const char type=OPLAST, const bool symm=true)
 
 UserDataOperator (const std::string &field_name, const char type, const bool symm=true)
 
 UserDataOperator (const std::string &row_field_name, const std::string &col_field_name, const char type, const bool symm=true)
 
virtual ~UserDataOperator ()
 
boost::shared_ptr< const NumeredEntFiniteElementgetNumeredEntFiniteElementPtr () const
 Return raw pointer to NumeredEntFiniteElement. More...
 
EntityHandle getFEEntityHandle () const
 Return finite element entity handle. More...
 
MoFEMErrorCode getProblemRowIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 Get row indices. More...
 
MoFEMErrorCode getProblemColIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 Get col indices. More...
 
virtual MoFEMErrorCode setPtrFE (ForcesAndSourcesCore *ptr)
 
const FEMethodgetFEMethod () const
 Return raw pointer to Finite Element Method object. More...
 
int getOpType () const
 Get operator types. More...
 
void setOpType (const OpType type)
 Set operator type. More...
 
void addOpType (const OpType type)
 Add operator type. More...
 
int getNinTheLoop () const
 get number of finite element in the loop More...
 
int getLoopSize () const
 get size of elements in the loop More...
 
const std::string & getFEName () const
 Get name of the element. More...
 
Vec getSnesF () const
 
Vec getSnesX () const
 
Mat getSnesA () const
 
Mat getSnesB () const
 
Vec getTSu () const
 
Vec getTSu_t () const
 
Vec getTSf () const
 
Mat getTSA () const
 
Mat getTSB () const
 
int getTSstep () const
 
double getTStime () const
 
double getTSa () const
 
MatrixDoublegetGaussPts ()
 matrix of integration (Gauss) points for Volume Element More...
 
auto getFTensor0IntegrationWeight ()
 Get integration weights. More...
 
DEPRECATED MoFEMErrorCode getPorblemRowIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 
DEPRECATED MoFEMErrorCode getPorblemColIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 
- Public Member Functions inherited from MoFEM::DataOperator
 DataOperator (const bool symm=true, const bool do_vertices=true, const bool do_edges=true, const bool do_quads=true, const bool do_tris=true, const bool do_tets=true, const bool do_prisms=true)
 
virtual ~DataOperator ()
 
virtual MoFEMErrorCode doWork (int row_side, int col_side, EntityType row_type, EntityType col_type, DataForcesAndSourcesCore::EntData &row_data, DataForcesAndSourcesCore::EntData &col_data)
 Operator for bi-linear form, usually to calculate values on left hand side. More...
 
virtual MoFEMErrorCode opLhs (DataForcesAndSourcesCore &row_data, DataForcesAndSourcesCore &col_data, bool symm=true)
 
virtual MoFEMErrorCode opLhs (DataForcesAndSourcesCore &row_data, DataForcesAndSourcesCore &col_data)
 
virtual MoFEMErrorCode doWork (int side, EntityType type, DataForcesAndSourcesCore::EntData &data)
 Operator for linear form, usually to calculate values on right hand side. More...
 
virtual MoFEMErrorCode opRhs (DataForcesAndSourcesCore &data, const bool do_vertices, const bool do_edges, const bool do_quads, const bool do_tris, const bool do_tets, const bool do_prisms, const bool error_if_no_base=true)
 
virtual MoFEMErrorCode opRhs (DataForcesAndSourcesCore &data, const bool error_if_no_base=true)
 
bool getSymm () const
 Get if operator uses symmetry of DOFs or not. More...
 
void setSymm ()
 set if operator is executed taking in account symmetry More...
 
void unSetSymm ()
 unset if operator is executed for non symmetric problem More...
 

Additional Inherited Members

- Public Types inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
enum  OpType { OPROW = 1 << 0, OPCOL = 1 << 1, OPROWCOL = 1 << 2, OPLAST = 1 << 3 }
 Controls loop over entities on element. More...
 
- Public Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
char opType
 
std::string rowFieldName
 
std::string colFieldName
 
FieldSpace sPace
 
- Public Attributes inherited from MoFEM::DataOperator
bool sYmm
 If true assume that matrix is symmetric structure. More...
 
bool doVertices
 If false skip vertices. More...
 
bool doEdges
 If false skip edges. More...
 
bool doQuads
 
bool doTris
 
bool doTets
 
bool doPrisms
 
- Protected Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
ForcesAndSourcesCoreptrFE
 

Detailed Description

default operator for TET element

Examples
elasticity.cpp, ElasticityMixedFormulation.hpp, forces_and_sources_testing_users_base.cpp, hello_world.cpp, HookeElement.cpp, HookeElement.hpp, MagneticElement.hpp, Remodeling.hpp, simple_elasticity.cpp, simple_interface.cpp, SmallStrainPlasticity.hpp, testing_jacobian_of_hook_scaled_with_density_element.cpp, and UnsaturatedFlow.hpp.

Definition at line 78 of file VolumeElementForcesAndSourcesCore.hpp.

Constructor & Destructor Documentation

◆ UserDataOperator() [1/3]

MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator::UserDataOperator ( const FieldSpace  space)
Examples
ElasticityMixedFormulation.hpp, and SmallStrainPlasticity.hpp.

Definition at line 80 of file VolumeElementForcesAndSourcesCore.hpp.

ForcesAndSourcesCore::UserDataOperator UserDataOperator

◆ UserDataOperator() [2/3]

MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator::UserDataOperator ( const std::string &  field_name,
const char  type 
)

Definition at line 83 of file VolumeElementForcesAndSourcesCore.hpp.

84  : ForcesAndSourcesCore::UserDataOperator(field_name, type) {}
ForcesAndSourcesCore::UserDataOperator UserDataOperator

◆ UserDataOperator() [3/3]

MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator::UserDataOperator ( const std::string &  row_field_name,
const std::string &  col_field_name,
const char  type,
const bool  symm = true 
)

Definition at line 86 of file VolumeElementForcesAndSourcesCore.hpp.

89  : ForcesAndSourcesCore::UserDataOperator(row_field_name, col_field_name,
90  type, symm) {}
ForcesAndSourcesCore::UserDataOperator UserDataOperator

Member Function Documentation

◆ getConn()

const EntityHandle* MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator::getConn ( )

get element connectivity

Definition at line 100 of file VolumeElementForcesAndSourcesCore.hpp.

100  {
101  return static_cast<VolumeElementForcesAndSourcesCore *>(ptrFE)->conn;
102  }

◆ getCoords()

VectorDouble& MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator::getCoords ( )

nodal coordinates

Definition at line 132 of file VolumeElementForcesAndSourcesCore.hpp.

132  {
133  return static_cast<VolumeElementForcesAndSourcesCore *>(ptrFE)->coords;
134  }

◆ getCoordsAtGaussPts()

MatrixDouble& MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator::getCoordsAtGaussPts ( )

Gauss points and weight, matrix (nb. of points x 3)

Column 0-2 integration points coordinate x, y and z, respectively. At rows are integration points.

Examples
UnsaturatedFlow.hpp.

Definition at line 142 of file VolumeElementForcesAndSourcesCore.hpp.

142  {
143  return static_cast<VolumeElementForcesAndSourcesCore *>(ptrFE)
144  ->coordsAtGaussPts;
145  }

◆ getCurlOfHCurlBaseFunctions()

MoFEMErrorCode UserDataOperator::getCurlOfHCurlBaseFunctions ( int  side,
EntityType  type,
DataForcesAndSourcesCore::EntData data,
int  gg,
MatrixDouble curl 
)

Get curl of base functions at integration point.

\[ \nabla \times \mathbf{\phi} \]

Works only for H-curl and H-div space.

How to use it:

MatrixDouble curl_mat(data.getFieldData().size(),3);
for(int gg = 0;gg<nb_gauss_pts;gg++) {
CHKERR getCurlOfHCurlBaseFunctions(side,type,data,gg,curl_mat);
&curl_mat(0,0), &curl_mat(0,1), &curl_mat(0,2)
);
// do something with curl
}

where curl_mat is matrix which number of rows is equal to nb. of base functions. Number of columns is 3, since we work in 3d here. Rows represents curl of base functions.

Parameters
sideside (local) number of entity on element
typetype of entity
datadata structure
gggauss pts
curlcurl matrix, nb. of rows is equal to number of base functions, columns are curl of base vector
Returns
error code
Examples
MagneticElement.hpp.

Definition at line 401 of file VolumeElementForcesAndSourcesCore.cpp.

403  {
405 
406  int nb_dofs = data.getFieldData().size();
407  if (nb_dofs == 0)
409 
410  if (data.getSpace() != HDIV && data.getSpace() != HCURL) {
412  "This function should be used for primarily for HCURL"
413  " but will work with HDIV used but is used with %s",
414  FieldSpaceNames[data.getSpace()]);
415  }
416 
417  if ((unsigned int)nb_dofs != data.getDiffN().size2() / 9) {
419  "Data insistency, wrong number of dofs = %s "
420  "%d != %d/9",
421  FieldSpaceNames[data.getSpace()], nb_dofs,
422  data.getDiffN().size2());
423  }
424 
425  curl.resize(nb_dofs, 3, false);
427  &curl(0, 0), &curl(0, 1), &curl(0, 2));
428  const double *grad_ptr = &data.getDiffN()(gg, 0);
429 
431  grad_ptr, &grad_ptr[HVEC0_1], &grad_ptr[HVEC0_2], &grad_ptr[HVEC1_0],
432  &grad_ptr[HVEC1_1], &grad_ptr[HVEC1_2], &grad_ptr[HVEC2_0],
433  &grad_ptr[HVEC2_1], &grad_ptr[HVEC2_2]);
434  for (int dd = 0; dd != nb_dofs; ++dd) {
435  t_curl(0) = t_grad_base(2, 1) - t_grad_base(1, 2);
436  t_curl(1) = t_grad_base(0, 2) - t_grad_base(2, 0);
437  t_curl(2) = t_grad_base(1, 0) - t_grad_base(0, 1);
438  ++t_curl;
439  ++t_grad_base;
440  }
441 
443 }
field with continuous normal traction
Definition: definitions.h:172
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:476
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:507
const Tensor2_symmetric_Expr< const ddTensor0< T, Dim, i, j >, typename promote< T, double >::V, Dim, i, j > dd(const Tensor0< T * > &a, const Index< i, Dim > index1, const Index< j, Dim > index2, const Tensor1< int, Dim > &d_ijk, const Tensor1< double, Dim > &d_xyz)
Definition: ddTensor0.hpp:33
static const char *const FieldSpaceNames[]
Definition: definitions.h:177
field with continuous tangents
Definition: definitions.h:171
const VolumeElementForcesAndSourcesCore * getVolumeFE()
return pointer to Generic Volume Finite Element object
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:406
virtual MPI_Comm & get_comm() const =0

◆ getDivergenceOfHDivBaseFunctions()

MoFEMErrorCode UserDataOperator::getDivergenceOfHDivBaseFunctions ( int  side,
EntityType  type,
DataForcesAndSourcesCore::EntData data,
int  gg,
VectorDouble div 
)

Get divergence of base functions at integration point.

Works only for H-div space.

How to use it:

VectorDouble div_vec(data.getFieldData().size());
for(int gg = 0;gg<nb_gauss_pts;gg++) {
CHKERR getDivergenceOfHDivBaseFunctions(side,type,data,gg,div_vec);
// do something with vec_div
}

where vec_div has size of nb. of dofs, at each element represents divergence of base functions.

Parameters
sideside (local) number of entity on element
typetype of entity
datadata structure
gggauss pts
divdivergence vector, size of vector is equal to number of base functions
Returns
error code
Examples
UnsaturatedFlow.hpp.

Definition at line 362 of file VolumeElementForcesAndSourcesCore.cpp.

364  {
366 
367  int nb_dofs = data.getFieldData().size();
368  if (nb_dofs == 0)
370 
371  if (data.getSpace() != HDIV && data.getSpace() != HCURL) {
373  "This function should be used for HDIV used but is used with %s",
374  FieldSpaceNames[data.getSpace()]);
375  }
376 
377  if ((unsigned int)nb_dofs != data.getDiffN().size2() / 9) {
379  "Data inositency, wrong number of dofs = %s "
380  "%d != %d/9",
381  FieldSpaceNames[data.getSpace()], nb_dofs,
382  data.getDiffN().size2());
383  }
384 
385  div.resize(nb_dofs, false);
386 
387  FTensor::Tensor0<double *> t_div(&*div.data().begin());
388  const double *grad_ptr = &data.getDiffN()(gg, 0);
390  grad_ptr, &grad_ptr[HVEC1_1], &grad_ptr[HVEC2_2]);
391  for (int dd = 0; dd < nb_dofs; dd++) {
392  t_div = t_grad_base(0) + t_grad_base(1) + t_grad_base(2);
393  ++t_div;
394  ++t_grad_base;
395  }
396 
398 }
field with continuous normal traction
Definition: definitions.h:172
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:476
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:507
const Tensor2_symmetric_Expr< const ddTensor0< T, Dim, i, j >, typename promote< T, double >::V, Dim, i, j > dd(const Tensor0< T * > &a, const Index< i, Dim > index1, const Index< j, Dim > index2, const Tensor1< int, Dim > &d_ijk, const Tensor1< double, Dim > &d_xyz)
Definition: ddTensor0.hpp:33
static const char *const FieldSpaceNames[]
Definition: definitions.h:177
field with continuous tangents
Definition: definitions.h:171
const VolumeElementForcesAndSourcesCore * getVolumeFE()
return pointer to Generic Volume Finite Element object
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:406
virtual MPI_Comm & get_comm() const =0

◆ getFTenosr0HoMeasure()

auto MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator::getFTenosr0HoMeasure ( )

◆ getFTensor1CoordsAtGaussPts()

auto MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator::getFTensor1CoordsAtGaussPts ( )

Get coordinates at integration points assuming linear geometry.

auto t_coords = getFTensor1CoordsAtGaussPts();
for(int gg = 0;gg!=nb_int_ptrs;gg++) {
// do something
++t_coords;
}
Examples
UnsaturatedFlow.hpp.

Definition at line 187 of file VolumeElementForcesAndSourcesCore.hpp.

187  {
189  &getCoordsAtGaussPts()(0, 0), &getCoordsAtGaussPts()(0, 1),
190  &getCoordsAtGaussPts()(0, 2));
191  }
MatrixDouble & getCoordsAtGaussPts()
Gauss points and weight, matrix (nb. of points x 3)

◆ getFTensor1HoCoordsAtGaussPts()

auto MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator::getFTensor1HoCoordsAtGaussPts ( )

Get coordinates at integration points taking geometry from field.

This is HO geometry given by arbitrary order polynomial

for(int gg = 0;gg!=nb_int_ptrs;gg++) {
// do something
++t_coords;
}

Definition at line 206 of file VolumeElementForcesAndSourcesCore.hpp.

206  {
207  double *ptr = &*getHoCoordsAtGaussPts().data().begin();
208  return FTensor::Tensor1<FTensor::PackPtr<double *, 3>, 3>(ptr, ptr + 1,
209  ptr + 2);
210  }
MatrixDouble & getHoCoordsAtGaussPts()
coordinate at Gauss points (if hierarchical approximation of element geometry)

◆ getFTensor2HoGaussPtsInvJac()

auto MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator::getFTensor2HoGaussPtsInvJac ( )

Definition at line 219 of file VolumeElementForcesAndSourcesCore.hpp.

219  {
220  double *ptr = &*getHoGaussPtsInvJac().data().begin();
222  ptr, ptr + 1, ptr + 2, ptr + 3, ptr + 4, ptr + 5, ptr + 6, ptr + 7,
223  ptr + 8);
224  }

◆ getFTensor2HoGaussPtsJac()

auto MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator::getFTensor2HoGaussPtsJac ( )

Definition at line 212 of file VolumeElementForcesAndSourcesCore.hpp.

212  {
213  double *ptr = &*getHoGaussPtsJac().data().begin();
215  ptr, ptr + 1, ptr + 2, ptr + 3, ptr + 4, ptr + 5, ptr + 6, ptr + 7,
216  ptr + 8);
217  }

◆ getHoCoordsAtGaussPts()

MatrixDouble& MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator::getHoCoordsAtGaussPts ( )

coordinate at Gauss points (if hierarchical approximation of element geometry)

Definition at line 150 of file VolumeElementForcesAndSourcesCore.hpp.

150  {
151  return static_cast<VolumeElementForcesAndSourcesCore *>(ptrFE)
152  ->hoCoordsAtGaussPts;
153  }

◆ getHoGaussPtsDetJac()

VectorDouble& MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator::getHoGaussPtsDetJac ( )
Examples
HookeElement.cpp, HookeElement.hpp, and MagneticElement.hpp.

Definition at line 165 of file VolumeElementForcesAndSourcesCore.hpp.

165  {
166  return static_cast<VolumeElementForcesAndSourcesCore *>(ptrFE)
167  ->hoGaussPtsDetJac;
168  }

◆ getHoGaussPtsInvJac()

MatrixDouble& MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator::getHoGaussPtsInvJac ( )

Definition at line 160 of file VolumeElementForcesAndSourcesCore.hpp.

160  {
161  return static_cast<VolumeElementForcesAndSourcesCore *>(ptrFE)
162  ->hoGaussPtsInvJac;
163  }

◆ getHoGaussPtsJac()

MatrixDouble& MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator::getHoGaussPtsJac ( )

Definition at line 155 of file VolumeElementForcesAndSourcesCore.hpp.

155  {
156  return static_cast<VolumeElementForcesAndSourcesCore *>(ptrFE)
157  ->hoGaussPtsJac;
158  }

◆ getInvJac()

FTensor::Tensor2<double *, 3, 3>& MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator::getInvJac ( )

get element inverse Jacobian

Definition at line 120 of file VolumeElementForcesAndSourcesCore.hpp.

120  {
121  return static_cast<VolumeElementForcesAndSourcesCore *>(ptrFE)->tInvJac;
122  }

◆ getJac()

FTensor::Tensor2<double *, 3, 3>& MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator::getJac ( )

get element Jacobian

Definition at line 113 of file VolumeElementForcesAndSourcesCore.hpp.

113  {
114  return static_cast<VolumeElementForcesAndSourcesCore *>(ptrFE)->tJac;
115  }

◆ getMeasure()

double MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator::getMeasure ( )

get measure of element

Returns
area of face

Definition at line 128 of file VolumeElementForcesAndSourcesCore.hpp.

128 { return getVolume(); }

◆ getNumNodes()

int MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator::getNumNodes ( )

get element number of nodes

Definition at line 94 of file VolumeElementForcesAndSourcesCore.hpp.

94  {
95  return static_cast<VolumeElementForcesAndSourcesCore *>(ptrFE)->num_nodes;
96  }

◆ getTensor1CoordsAtGaussPts()

DEPRECATED auto MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator::getTensor1CoordsAtGaussPts ( )
Deprecated:
use getFTensor1CoordsAtGaussPts

Definition at line 299 of file VolumeElementForcesAndSourcesCore.hpp.

299  {
301  }
auto getFTensor1CoordsAtGaussPts()
Get coordinates at integration points assuming linear geometry.

◆ getTensor1HoCoordsAtGaussPts()

DEPRECATED auto MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator::getTensor1HoCoordsAtGaussPts ( )
Deprecated:
use getFTensor1HoCoordsAtGaussPts

Definition at line 304 of file VolumeElementForcesAndSourcesCore.hpp.

304  {
306  }
auto getFTensor1HoCoordsAtGaussPts()
Get coordinates at integration points taking geometry from field.

◆ getVolume()

double MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator::getVolume ( )

element volume (linear geometry)

Examples
HookeElement.cpp, HookeElement.hpp, MagneticElement.hpp, and UnsaturatedFlow.hpp.

Definition at line 106 of file VolumeElementForcesAndSourcesCore.hpp.

106  {
107  return static_cast<VolumeElementForcesAndSourcesCore *>(ptrFE)->vOlume;
108  }

◆ getVolumeFE()

const VolumeElementForcesAndSourcesCore* MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator::getVolumeFE ( )

return pointer to Generic Volume Finite Element object

Definition at line 228 of file VolumeElementForcesAndSourcesCore.hpp.

228  {
229  return static_cast<VolumeElementForcesAndSourcesCore *>(ptrFE);
230  }

The documentation for this struct was generated from the following files: