v0.10.0
Public Member Functions | Protected 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 EntityHandle * getConn ()
 get element connectivity More...
 
double getVolume () const
 element volume (linear geometry) 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 () const
 get measure of element 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 ()
 
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...
 
const PetscData::SwitchesgetDataCtx () const
 
const KspMethod::KSPContext getKSPCtx () const
 
const SnesMethod::SNESContext getSNESCtx () const
 
const TSMethod::TSContext getTSCtx () const
 
Vec getKSPf () const
 
Mat getKSPA () const
 
Mat getKSPB () const
 
Vec getSNESf () const
 
Vec getSNESx () const
 
Mat getSNESA () const
 
Mat getSNESB () const
 
DEPRECATED Vec getSnesF () const
 
DEPRECATED Vec getSnesX () const
 
DEPRECATED Mat getSnesA () const
 
DEPRECATED Mat getSnesB () const
 
Vec getTSu () const
 
Vec getTSu_t () const
 
Vec getTSu_tt () const
 
Vec getTSf () const
 
Mat getTSA () const
 
Mat getTSB () const
 
int getTSstep () const
 
double getTStime () const
 
double getTSa () const
 
double getTSaa () 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)
 
virtual ~DataOperator ()=default
 
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)
 
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 error_if_no_base=false)
 
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...
 

Protected Member Functions

MoFEMErrorCode setPtrFE (ForcesAndSourcesCore *ptr)
 
- Protected Member Functions inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
ForcesAndSourcesCoregetPtrFE () const
 
ForcesAndSourcesCoregetSidePtrFE () const
 

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...
 
std::array< bool, MBMAXTYPE > doEntities
 If true operator is executed for entity. More...
 
booldoVertices
 \deprectaed If false skip vertices More...
 
booldoEdges
 \deprectaed If false skip edges More...
 
booldoQuads
 \deprectaed More...
 
booldoTris
 \deprectaed More...
 
booldoTets
 \deprectaed More...
 
booldoPrisms
 \deprectaed More...
 
- Protected Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
ForcesAndSourcesCoreptrFE
 

Detailed Description

default operator for TET element

Examples
ElasticityMixedFormulation.hpp, EshelbianPlasticity.cpp, HookeElement.cpp, HookeElement.hpp, MagneticElement.hpp, OpPostProcElastic.hpp, PlasticOps.hpp, Remodeling.hpp, SmallStrainPlasticity.hpp, UnsaturatedFlow.hpp, build_large_problem.cpp, contact.cpp, dynamic_elastic.cpp, eigen_elastic.cpp, elasticity.cpp, field_evaluator.cpp, heat_method.cpp, helmholtz.cpp, lorentz_force.cpp, plastic.cpp, and test_cache_on_entities.cpp.

Definition at line 45 of file VolumeElementForcesAndSourcesCore.hpp.

Member Function Documentation

◆ getConn()

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

get element connectivity

Definition at line 398 of file VolumeElementForcesAndSourcesCore.hpp.

398  {
399  return static_cast<VolumeElementForcesAndSourcesCoreBase *>(ptrFE)->conn;
400 }
VolumeElementForcesAndSourcesCoreBase(Interface &m_field, const EntityType type=MBTET)

◆ getCoords()

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

nodal coordinates

Definition at line 431 of file VolumeElementForcesAndSourcesCore.hpp.

◆ 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 436 of file VolumeElementForcesAndSourcesCore.hpp.

◆ 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
}
#define CHKERR
Inline error check.
Definition: definitions.h:604
ublas::matrix< double, ublas::row_major, DoubleAllocator > MatrixDouble
Definition: Types.hpp:76
MoFEMErrorCode getCurlOfHCurlBaseFunctions(int side, EntityType type, DataForcesAndSourcesCore::EntData &data, int gg, MatrixDouble &curl)
Get curl of base functions at integration point.

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 379 of file VolumeElementForcesAndSourcesCore.cpp.

382  {
384 
385  int nb_dofs = data.getFieldData().size();
386  if (nb_dofs == 0)
388 
389  if (data.getSpace() != HDIV && data.getSpace() != HCURL) {
391  "This function should be used for primarily for HCURL"
392  " but will work with HDIV used but is used with %s",
393  FieldSpaceNames[data.getSpace()]);
394  }
395 
396  if ((unsigned int)nb_dofs != data.getDiffN().size2() / 9) {
398  "Data insistency, wrong number of dofs = %s "
399  "%d != %d/9",
400  FieldSpaceNames[data.getSpace()], nb_dofs,
401  data.getDiffN().size2());
402  }
403 
404  curl.resize(nb_dofs, 3, false);
406  &curl(0, 0), &curl(0, 1), &curl(0, 2));
407  const double *grad_ptr = &data.getDiffN()(gg, 0);
408 
410  grad_ptr, &grad_ptr[HVEC0_1], &grad_ptr[HVEC0_2], &grad_ptr[HVEC1_0],
411  &grad_ptr[HVEC1_1], &grad_ptr[HVEC1_2], &grad_ptr[HVEC2_0],
412  &grad_ptr[HVEC2_1], &grad_ptr[HVEC2_2]);
413  for (int dd = 0; dd != nb_dofs; ++dd) {
414  t_curl(0) = t_grad_base(2, 1) - t_grad_base(1, 2);
415  t_curl(1) = t_grad_base(0, 2) - t_grad_base(2, 0);
416  t_curl(2) = t_grad_base(1, 0) - t_grad_base(0, 1);
417  ++t_curl;
418  ++t_grad_base;
419  }
420 
422 }
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:516
@ HCURL
field with continuous tangents
Definition: definitions.h:178
@ HDIV
field with continuous normal traction
Definition: definitions.h:179
static const char *const FieldSpaceNames[]
Definition: definitions.h:184
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:415
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:123
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:485
@ HVEC1_1
Definition: definitions.h:265
@ HVEC0_1
Definition: definitions.h:264
@ HVEC1_0
Definition: definitions.h:262
@ HVEC2_1
Definition: definitions.h:266
@ HVEC1_2
Definition: definitions.h:268
@ HVEC2_2
Definition: definitions.h:269
@ HVEC2_0
Definition: definitions.h:263
@ HVEC0_2
Definition: definitions.h:267
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
virtual MPI_Comm & get_comm() const =0
VolumeElementForcesAndSourcesCoreBase * getVolumeFE() const
return pointer to Generic Volume Finite Element object

◆ 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++) {
// do something with vec_div
}
ublas::vector< double, DoubleAllocator > VectorDouble
Definition: Types.hpp:74
MoFEMErrorCode getDivergenceOfHDivBaseFunctions(int side, EntityType type, DataForcesAndSourcesCore::EntData &data, int gg, VectorDouble &div)
Get divergence of base functions at integration point.

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 340 of file VolumeElementForcesAndSourcesCore.cpp.

343  {
345 
346  int nb_dofs = data.getFieldData().size();
347  if (nb_dofs == 0)
349 
350  if (data.getSpace() != HDIV && data.getSpace() != HCURL) {
352  "This function should be used for HDIV used but is used with %s",
353  FieldSpaceNames[data.getSpace()]);
354  }
355 
356  if ((unsigned int)nb_dofs != data.getDiffN().size2() / 9) {
358  "Data inositency, wrong number of dofs = %s "
359  "%d != %d/9",
360  FieldSpaceNames[data.getSpace()], nb_dofs,
361  data.getDiffN().size2());
362  }
363 
364  div.resize(nb_dofs, false);
365 
366  FTensor::Tensor0<double *> t_div(&*div.data().begin());
367  const double *grad_ptr = &data.getDiffN()(gg, 0);
369  grad_ptr, &grad_ptr[HVEC1_1], &grad_ptr[HVEC2_2]);
370  for (int dd = 0; dd < nb_dofs; dd++) {
371  t_div = t_grad_base(0) + t_grad_base(1) + t_grad_base(2);
372  ++t_div;
373  ++t_grad_base;
374  }
375 
377 }

◆ 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;
}
auto getFTensor1CoordsAtGaussPts()
Get coordinates at integration points assuming linear geometry.
Examples
UnsaturatedFlow.hpp, and prism_polynomial_approximation.cpp.

Definition at line 471 of file VolumeElementForcesAndSourcesCore.hpp.

472  {
474  &getCoordsAtGaussPts()(0, 0), &getCoordsAtGaussPts()(0, 1),
475  &getCoordsAtGaussPts()(0, 2));
476 }
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;
}
auto getFTensor1HoCoordsAtGaussPts()
Get coordinates at integration points taking geometry from field.

Definition at line 478 of file VolumeElementForcesAndSourcesCore.hpp.

479  {
480  double *ptr = &*getHoCoordsAtGaussPts().data().begin();
481  return FTensor::Tensor1<FTensor::PackPtr<double *, 3>, 3>(ptr, ptr + 1,
482  ptr + 2);
483 }
MatrixDouble & getHoCoordsAtGaussPts()
coordinate at Gauss points (if hierarchical approximation of element geometry)

◆ getFTensor2HoGaussPtsInvJac()

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

Definition at line 493 of file VolumeElementForcesAndSourcesCore.hpp.

494  {
495  double *ptr = &*getHoGaussPtsInvJac().data().begin();
497  ptr, ptr + 1, ptr + 2, ptr + 3, ptr + 4, ptr + 5, ptr + 6, ptr + 7,
498  ptr + 8);
499 }

◆ getFTensor2HoGaussPtsJac()

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

Definition at line 485 of file VolumeElementForcesAndSourcesCore.hpp.

486  {
487  double *ptr = &*getHoGaussPtsJac().data().begin();
489  ptr, ptr + 1, ptr + 2, ptr + 3, ptr + 4, ptr + 5, ptr + 6, ptr + 7,
490  ptr + 8);
491 }

◆ getHoCoordsAtGaussPts()

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

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

Definition at line 441 of file VolumeElementForcesAndSourcesCore.hpp.

◆ getHoGaussPtsDetJac()

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

◆ getHoGaussPtsInvJac()

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

◆ getHoGaussPtsJac()

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

◆ getInvJac()

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

get element inverse Jacobian

Definition at line 417 of file VolumeElementForcesAndSourcesCore.hpp.

417  {
418  return static_cast<VolumeElementForcesAndSourcesCoreBase *>(ptrFE)->tInvJac;
419 }

◆ getJac()

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

get element Jacobian

Definition at line 412 of file VolumeElementForcesAndSourcesCore.hpp.

412  {
413  return static_cast<VolumeElementForcesAndSourcesCoreBase *>(ptrFE)->tJac;
414 }

◆ getMeasure() [1/2]

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

get measure of element

Returns
volume

Definition at line 426 of file VolumeElementForcesAndSourcesCore.hpp.

426  {
427  return getVolume();
428 }

◆ getMeasure() [2/2]

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

get measure of element

Returns
volume
Examples
PlasticOps.hpp, and heat_method.cpp.

Definition at line 422 of file VolumeElementForcesAndSourcesCore.hpp.

422  {
423  return getVolume();
424 }

◆ getNumNodes()

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

get element number of nodes

Definition at line 393 of file VolumeElementForcesAndSourcesCore.hpp.

◆ getVolume() [1/2]

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

element volume (linear geometry)

Definition at line 407 of file VolumeElementForcesAndSourcesCore.hpp.

◆ getVolume() [2/2]

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

element volume (linear geometry)

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

Definition at line 403 of file VolumeElementForcesAndSourcesCore.hpp.

403  {
404  return static_cast<VolumeElementForcesAndSourcesCoreBase *>(ptrFE)->vOlume;
405 }

◆ getVolumeFE()

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

return pointer to Generic Volume Finite Element object

Definition at line 502 of file VolumeElementForcesAndSourcesCore.hpp.

502  {
503  return static_cast<VolumeElementForcesAndSourcesCoreBase *>(ptrFE);
504 }

◆ setPtrFE()

MoFEMErrorCode MoFEM::VolumeElementForcesAndSourcesCoreBase::UserDataOperator::setPtrFE ( ForcesAndSourcesCore ptr)
protectedvirtual

Reimplemented from MoFEM::ForcesAndSourcesCore::UserDataOperator.

Reimplemented in MoFEM::VolumeElementForcesAndSourcesCoreOnSideBase::UserDataOperator.

Definition at line 425 of file VolumeElementForcesAndSourcesCore.cpp.

426  {
428  if(!(ptrFE = dynamic_cast<VolumeElementForcesAndSourcesCoreBase *>(ptr)))
429  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
430  "User operator and finite element do not work together");
432 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:509

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