v0.9.0
Public Member Functions | List of all members
MoFEM::VolumeElementForcesAndSourcesCoreBase::UserDataOperator Struct Reference

default operator for TET element More...

#include <src/finite_elements/VolumeElementForcesAndSourcesCore.hpp>

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

Public Member Functions

int getNumNodes ()
 get element number of nodes More...
 
const EntityHandlegetConn ()
 get element connectivity More...
 
double getVolume () const
 element volume (linear geometry) More...
 
doublegetVolume ()
 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 () const
 get measure of element More...
 
doublegetMeasure ()
 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 ()
 
VolumeElementForcesAndSourcesCoreBasegetVolumeFE () const
 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...
 
- 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)
 
boost::shared_ptr< const NumeredEntFiniteElementgetNumeredEntFiniteElementPtr () const
 Return raw pointer to NumeredEntFiniteElement. More...
 
EntityHandle getFEEntityHandle () const
 Return finite element entity handle. More...
 
boost::weak_ptr< SideNumber > getSideNumberPtr (const int side_number, const EntityType type)
 Get the side number pointer. More...
 
EntityHandle getSideEntity (const int side_number, const EntityType type)
 Get the side entity. More...
 
int getNumberOfNodesOnElement ()
 Get the number of nodes on finite element. 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...
 
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
 
- 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 Member Functions inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
MoFEMErrorCode setPtrFE (ForcesAndSourcesCore *ptr)
 
ForcesAndSourcesCoregetPtrFE () const
 
ForcesAndSourcesCoregetSidePtrFE () const
 
- 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, hcurl_curl_operator.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 48 of file VolumeElementForcesAndSourcesCore.hpp.

Member Function Documentation

◆ getConn()

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

get element connectivity

Definition at line 392 of file VolumeElementForcesAndSourcesCore.hpp.

392  {
393  return static_cast<VolumeElementForcesAndSourcesCoreBase *>(ptrFE)->conn;
394 }

◆ getCoords()

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

nodal coordinates

Definition at line 425 of file VolumeElementForcesAndSourcesCore.hpp.

425  {
426  return static_cast<VolumeElementForcesAndSourcesCoreBase *>(ptrFE)->coords;
427 }

◆ getCoordsAtGaussPts()

MatrixDouble & MoFEM::VolumeElementForcesAndSourcesCoreBase::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 430 of file VolumeElementForcesAndSourcesCore.hpp.

430  {
431  return static_cast<VolumeElementForcesAndSourcesCoreBase *>(ptrFE)
432  ->coordsAtGaussPts;
433 }

◆ getCurlOfHCurlBaseFunctions()

MoFEMErrorCode MoFEM::VolumeElementForcesAndSourcesCoreBase::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 368 of file VolumeElementForcesAndSourcesCore.cpp.

370  {
372 
373  int nb_dofs = data.getFieldData().size();
374  if (nb_dofs == 0)
376 
377  if (data.getSpace() != HDIV && data.getSpace() != HCURL) {
379  "This function should be used for primarily for HCURL"
380  " but will work with HDIV used but is used with %s",
381  FieldSpaceNames[data.getSpace()]);
382  }
383 
384  if ((unsigned int)nb_dofs != data.getDiffN().size2() / 9) {
386  "Data insistency, wrong number of dofs = %s "
387  "%d != %d/9",
388  FieldSpaceNames[data.getSpace()], nb_dofs,
389  data.getDiffN().size2());
390  }
391 
392  curl.resize(nb_dofs, 3, false);
394  &curl(0, 0), &curl(0, 1), &curl(0, 2));
395  const double *grad_ptr = &data.getDiffN()(gg, 0);
396 
398  grad_ptr, &grad_ptr[HVEC0_1], &grad_ptr[HVEC0_2], &grad_ptr[HVEC1_0],
399  &grad_ptr[HVEC1_1], &grad_ptr[HVEC1_2], &grad_ptr[HVEC2_0],
400  &grad_ptr[HVEC2_1], &grad_ptr[HVEC2_2]);
401  for (int dd = 0; dd != nb_dofs; ++dd) {
402  t_curl(0) = t_grad_base(2, 1) - t_grad_base(1, 2);
403  t_curl(1) = t_grad_base(0, 2) - t_grad_base(2, 0);
404  t_curl(2) = t_grad_base(1, 0) - t_grad_base(0, 1);
405  ++t_curl;
406  ++t_grad_base;
407  }
408 
410 }
field with continuous normal traction
Definition: definitions.h:173
VolumeElementForcesAndSourcesCoreBase * getVolumeFE() const
return pointer to Generic Volume Finite Element object
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:508
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:178
field with continuous tangents
Definition: definitions.h:172
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:407
virtual MPI_Comm & get_comm() const =0

◆ getDivergenceOfHDivBaseFunctions()

MoFEMErrorCode MoFEM::VolumeElementForcesAndSourcesCoreBase::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 329 of file VolumeElementForcesAndSourcesCore.cpp.

331  {
333 
334  int nb_dofs = data.getFieldData().size();
335  if (nb_dofs == 0)
337 
338  if (data.getSpace() != HDIV && data.getSpace() != HCURL) {
340  "This function should be used for HDIV used but is used with %s",
341  FieldSpaceNames[data.getSpace()]);
342  }
343 
344  if ((unsigned int)nb_dofs != data.getDiffN().size2() / 9) {
346  "Data inositency, wrong number of dofs = %s "
347  "%d != %d/9",
348  FieldSpaceNames[data.getSpace()], nb_dofs,
349  data.getDiffN().size2());
350  }
351 
352  div.resize(nb_dofs, false);
353 
354  FTensor::Tensor0<double *> t_div(&*div.data().begin());
355  const double *grad_ptr = &data.getDiffN()(gg, 0);
357  grad_ptr, &grad_ptr[HVEC1_1], &grad_ptr[HVEC2_2]);
358  for (int dd = 0; dd < nb_dofs; dd++) {
359  t_div = t_grad_base(0) + t_grad_base(1) + t_grad_base(2);
360  ++t_div;
361  ++t_grad_base;
362  }
363 
365 }
field with continuous normal traction
Definition: definitions.h:173
VolumeElementForcesAndSourcesCoreBase * getVolumeFE() const
return pointer to Generic Volume Finite Element object
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:508
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:178
field with continuous tangents
Definition: definitions.h:172
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:407
virtual MPI_Comm & get_comm() const =0

◆ getFTenosr0HoMeasure()

auto MoFEM::VolumeElementForcesAndSourcesCoreBase::UserDataOperator::getFTenosr0HoMeasure ( )

◆ getFTensor1CoordsAtGaussPts()

auto MoFEM::VolumeElementForcesAndSourcesCoreBase::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
prism_polynomial_approximation.cpp, and UnsaturatedFlow.hpp.

Definition at line 466 of file VolumeElementForcesAndSourcesCore.hpp.

466  {
468  &getCoordsAtGaussPts()(0, 0), &getCoordsAtGaussPts()(0, 1),
469  &getCoordsAtGaussPts()(0, 2));
470 }
MatrixDouble & getCoordsAtGaussPts()
Gauss points and weight, matrix (nb. of points x 3)

◆ getFTensor1HoCoordsAtGaussPts()

auto MoFEM::VolumeElementForcesAndSourcesCoreBase::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 473 of file VolumeElementForcesAndSourcesCore.hpp.

473  {
474  double *ptr = &*getHoCoordsAtGaussPts().data().begin();
475  return FTensor::Tensor1<FTensor::PackPtr<double *, 3>, 3>(ptr, ptr + 1,
476  ptr + 2);
477 }
MatrixDouble & getHoCoordsAtGaussPts()
coordinate at Gauss points (if hierarchical approximation of element geometry)

◆ getFTensor2HoGaussPtsInvJac()

auto MoFEM::VolumeElementForcesAndSourcesCoreBase::UserDataOperator::getFTensor2HoGaussPtsInvJac ( )

Definition at line 488 of file VolumeElementForcesAndSourcesCore.hpp.

488  {
489  double *ptr = &*getHoGaussPtsInvJac().data().begin();
491  ptr, ptr + 1, ptr + 2, ptr + 3, ptr + 4, ptr + 5, ptr + 6, ptr + 7,
492  ptr + 8);
493 }

◆ getFTensor2HoGaussPtsJac()

auto MoFEM::VolumeElementForcesAndSourcesCoreBase::UserDataOperator::getFTensor2HoGaussPtsJac ( )

Definition at line 480 of file VolumeElementForcesAndSourcesCore.hpp.

480  {
481  double *ptr = &*getHoGaussPtsJac().data().begin();
483  ptr, ptr + 1, ptr + 2, ptr + 3, ptr + 4, ptr + 5, ptr + 6, ptr + 7,
484  ptr + 8);
485 }

◆ getHoCoordsAtGaussPts()

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

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

Definition at line 436 of file VolumeElementForcesAndSourcesCore.hpp.

436  {
437  return static_cast<VolumeElementForcesAndSourcesCoreBase *>(ptrFE)
438  ->hoCoordsAtGaussPts;
439 }

◆ getHoGaussPtsDetJac()

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

Definition at line 454 of file VolumeElementForcesAndSourcesCore.hpp.

454  {
455  return static_cast<VolumeElementForcesAndSourcesCoreBase *>(ptrFE)
456  ->hoGaussPtsDetJac;
457 }

◆ getHoGaussPtsInvJac()

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

Definition at line 448 of file VolumeElementForcesAndSourcesCore.hpp.

448  {
449  return static_cast<VolumeElementForcesAndSourcesCoreBase *>(ptrFE)
450  ->hoGaussPtsInvJac;
451 }

◆ getHoGaussPtsJac()

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

Definition at line 442 of file VolumeElementForcesAndSourcesCore.hpp.

442  {
443  return static_cast<VolumeElementForcesAndSourcesCoreBase *>(ptrFE)
444  ->hoGaussPtsJac;
445 }

◆ getInvJac()

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

get element inverse Jacobian

Definition at line 411 of file VolumeElementForcesAndSourcesCore.hpp.

411  {
412  return static_cast<VolumeElementForcesAndSourcesCoreBase *>(ptrFE)->tInvJac;
413 }

◆ getJac()

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

get element Jacobian

Definition at line 406 of file VolumeElementForcesAndSourcesCore.hpp.

406  {
407  return static_cast<VolumeElementForcesAndSourcesCoreBase *>(ptrFE)->tJac;
408 }

◆ getMeasure() [1/2]

double MoFEM::VolumeElementForcesAndSourcesCoreBase::UserDataOperator::getMeasure ( ) const

get measure of element

Returns
volume

Definition at line 416 of file VolumeElementForcesAndSourcesCore.hpp.

416  {
417  return getVolume();
418 }

◆ getMeasure() [2/2]

double & MoFEM::VolumeElementForcesAndSourcesCoreBase::UserDataOperator::getMeasure ( )

get measure of element

Returns
volume

Definition at line 420 of file VolumeElementForcesAndSourcesCore.hpp.

420  {
421  return getVolume();
422 }

◆ getNumNodes()

int MoFEM::VolumeElementForcesAndSourcesCoreBase::UserDataOperator::getNumNodes ( )

get element number of nodes

Definition at line 387 of file VolumeElementForcesAndSourcesCore.hpp.

387  {
388  return static_cast<VolumeElementForcesAndSourcesCoreBase *>(ptrFE)->num_nodes;
389 }

◆ getVolume() [1/2]

double MoFEM::VolumeElementForcesAndSourcesCoreBase::UserDataOperator::getVolume ( ) const

element volume (linear geometry)

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

Definition at line 397 of file VolumeElementForcesAndSourcesCore.hpp.

397  {
398  return static_cast<VolumeElementForcesAndSourcesCoreBase *>(ptrFE)->vOlume;
399 }

◆ getVolume() [2/2]

double & MoFEM::VolumeElementForcesAndSourcesCoreBase::UserDataOperator::getVolume ( )

element volume (linear geometry)

Definition at line 401 of file VolumeElementForcesAndSourcesCore.hpp.

401  {
402  return static_cast<VolumeElementForcesAndSourcesCoreBase *>(ptrFE)->vOlume;
403 }

◆ getVolumeFE()

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

return pointer to Generic Volume Finite Element object

Definition at line 496 of file VolumeElementForcesAndSourcesCore.hpp.

496  {
497  return static_cast<VolumeElementForcesAndSourcesCoreBase *>(ptrFE);
498 }

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