v0.8.4
Public Member Functions | List of all members
MoFEM::FaceElementForcesAndSourcesCore::UserDataOperator Struct Reference

default operator for TRI element More...

#include <src/finite_elements/FaceElementForcesAndSourcesCore.hpp>

Inheritance diagram for MoFEM::FaceElementForcesAndSourcesCore::UserDataOperator:
[legend]
Collaboration diagram for MoFEM::FaceElementForcesAndSourcesCore::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)
 
double getArea ()
 get area of face More...
 
double getMeasure ()
 get measure of element More...
 
VectorDoublegetNormal ()
 get triangle normal More...
 
VectorDoublegetTangent1 ()
 get triangle tangent 1 More...
 
VectorDoublegetTangent2 ()
 get triangle tangent 2 More...
 
FTensor::Tensor1< double *, 3 > getTensor1Normal ()
 get normal as tensor More...
 
FTensor::Tensor1< double *, 3 > getTensor1Tangent1 ()
 get tangentOne as tensor More...
 
FTensor::Tensor1< double *, 3 > getTensor2Tangent1 ()
 get tangentTwo as tensor More...
 
int getNumNodes ()
 get element number of nodes More...
 
const EntityHandlegetConn ()
 get element connectivity More...
 
VectorDoublegetCoords ()
 get triangle coordinates More...
 
FTensor::Tensor1< double *, 3 > getTensor1Coords ()
 get get coords at gauss points More...
 
MatrixDoublegetGaussPts ()
 get matrix of integration (Gauss) points on Face Element where columns 0,1 are x,y coordinates respectively and column 2 is a value of weight for example getGaussPts()(1,13) returns y coordinate of 13th Gauss point on particular face element More...
 
FTensor::Tensor0< double * > getFTensor0IntegrationWeight ()
 
MatrixDoublegetCoordsAtGaussPts ()
 Gauss points and weight, matrix (nb. of points x 3) More...
 
FTensor::Tensor1< double *, 3 > getTensor1CoordsAtGaussPts ()
 get coordinates at Gauss pts. More...
 
MatrixDoublegetHoCoordsAtGaussPts ()
 coordinate at Gauss points (if hierarchical approximation of element geometry) More...
 
FTensor::Tensor1< double *, 3 > getTensor1HoCoordsAtGaussPts ()
 get coordinates at Gauss pts (takes in account ho approx. of geometry) More...
 
MatrixDoublegetNormalsAtGaussPt ()
 if higher order geometry return normals at Gauss pts. More...
 
ublas::matrix_row< MatrixDoublegetNormalsAtGaussPt (const int gg)
 if higher order geometry return normals at Gauss pts. More...
 
MatrixDoublegetTangent1AtGaussPt ()
 if higher order geometry return tangent vector to triangle at Gauss pts. More...
 
MatrixDoublegetTangent2AtGaussPt ()
 if higher order geometry return tangent vector to triangle at Gauss pts. More...
 
FTensor::Tensor1< double *, 3 > getTensor1NormalsAtGaussPt ()
 get normal at integration points More...
 
FTensor::Tensor1< double *, 3 > getTensor1Tangent1AtGaussPt ()
 get tangent 1 at integration points More...
 
FTensor::Tensor1< double *, 3 > getTensor1Tangent2AtGaussPt ()
 get tangent 2 at integration points More...
 
DEPRECATED const FaceElementForcesAndSourcesCoregetFaceElementForcesAndSourcesCore ()
 
DEPRECATED const FaceElementForcesAndSourcesCoregetTriFE ()
 
const FaceElementForcesAndSourcesCoregetFaceFE ()
 return pointer to Generic Triangle Finite Element object More...
 
MoFEMErrorCode loopSideVolumes (const string &fe_name, VolumeElementForcesAndSourcesCoreOnSide &method)
 
- 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...
 
DEPRECATED MoFEMErrorCode getPorblemRowIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 
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)
 
DEPRECATED MoFEMErrorCode getPorblemColIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 
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
 
- 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 right 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 left 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 TRI element

Examples:
hello_world.cpp, MagneticElement.hpp, simple_elasticity.cpp, simple_interface.cpp, and UnsaturatedFlow.hpp.

Definition at line 86 of file FaceElementForcesAndSourcesCore.hpp.

Constructor & Destructor Documentation

◆ UserDataOperator() [1/3]

MoFEM::FaceElementForcesAndSourcesCore::UserDataOperator::UserDataOperator ( const FieldSpace  space)

Definition at line 88 of file FaceElementForcesAndSourcesCore.hpp.

ForcesAndSourcesCore::UserDataOperator UserDataOperator

◆ UserDataOperator() [2/3]

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

Definition at line 91 of file FaceElementForcesAndSourcesCore.hpp.

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

◆ UserDataOperator() [3/3]

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

Definition at line 94 of file FaceElementForcesAndSourcesCore.hpp.

97  : ForcesAndSourcesCore::UserDataOperator(row_field_name, col_field_name,
98  type, symm){};
ForcesAndSourcesCore::UserDataOperator UserDataOperator

Member Function Documentation

◆ getArea()

double MoFEM::FaceElementForcesAndSourcesCore::UserDataOperator::getArea ( )

get area of face

Returns
area of face

Definition at line 104 of file FaceElementForcesAndSourcesCore.hpp.

◆ getConn()

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

get element connectivity

Definition at line 161 of file FaceElementForcesAndSourcesCore.hpp.

◆ getCoords()

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

◆ getCoordsAtGaussPts()

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

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

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

Definition at line 211 of file FaceElementForcesAndSourcesCore.hpp.

◆ getFaceElementForcesAndSourcesCore()

DEPRECATED const FaceElementForcesAndSourcesCore* MoFEM::FaceElementForcesAndSourcesCore::UserDataOperator::getFaceElementForcesAndSourcesCore ( )
Deprecated:
use getFaceFE instead

Definition at line 331 of file FaceElementForcesAndSourcesCore.hpp.

331  {
332  return getFaceFE();
333  }
const FaceElementForcesAndSourcesCore * getFaceFE()
return pointer to Generic Triangle Finite Element object

◆ getFaceFE()

const FaceElementForcesAndSourcesCore* MoFEM::FaceElementForcesAndSourcesCore::UserDataOperator::getFaceFE ( )

return pointer to Generic Triangle Finite Element object

Definition at line 343 of file FaceElementForcesAndSourcesCore.hpp.

◆ getFTensor0IntegrationWeight()

FTensor::Tensor0<double *> MoFEM::FaceElementForcesAndSourcesCore::UserDataOperator::getFTensor0IntegrationWeight ( )

Definition at line 201 of file FaceElementForcesAndSourcesCore.hpp.

201  {
202  return FTensor::Tensor0<double *>(&(getGaussPts()(2, 0)), 1);
203  }
MatrixDouble & getGaussPts()
get matrix of integration (Gauss) points on Face Element where columns 0,1 are x,y coordinates respec...

◆ getGaussPts()

MatrixDouble& MoFEM::FaceElementForcesAndSourcesCore::UserDataOperator::getGaussPts ( )

get matrix of integration (Gauss) points on Face Element where columns 0,1 are x,y coordinates respectively and column 2 is a value of weight for example getGaussPts()(1,13) returns y coordinate of 13th Gauss point on particular face element

Definition at line 197 of file FaceElementForcesAndSourcesCore.hpp.

◆ getHoCoordsAtGaussPts()

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

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

Note: returned matrix has size 0 in rows and columns if no HO approximation of geometry is available.

Definition at line 230 of file FaceElementForcesAndSourcesCore.hpp.

◆ getMeasure()

double MoFEM::FaceElementForcesAndSourcesCore::UserDataOperator::getMeasure ( )

get measure of element

Returns
area of face

Definition at line 112 of file FaceElementForcesAndSourcesCore.hpp.

◆ getNormal()

VectorDouble& MoFEM::FaceElementForcesAndSourcesCore::UserDataOperator::getNormal ( )

◆ getNormalsAtGaussPt() [1/2]

MatrixDouble& MoFEM::FaceElementForcesAndSourcesCore::UserDataOperator::getNormalsAtGaussPt ( )

if higher order geometry return normals at Gauss pts.

Note: returned matrix has size 0 in rows and columns if no HO approximation of geometry is available.

Definition at line 253 of file FaceElementForcesAndSourcesCore.hpp.

◆ getNormalsAtGaussPt() [2/2]

ublas::matrix_row<MatrixDouble> MoFEM::FaceElementForcesAndSourcesCore::UserDataOperator::getNormalsAtGaussPt ( const int  gg)

if higher order geometry return normals at Gauss pts.

Parameters
gggauss point number

Definition at line 262 of file FaceElementForcesAndSourcesCore.hpp.

262  {
263  return ublas::matrix_row<MatrixDouble>(
264  static_cast<FaceElementForcesAndSourcesCore *>(ptrFE)
266  gg);
267  }

◆ getNumNodes()

int MoFEM::FaceElementForcesAndSourcesCore::UserDataOperator::getNumNodes ( )

◆ getTangent1()

VectorDouble& MoFEM::FaceElementForcesAndSourcesCore::UserDataOperator::getTangent1 ( )

◆ getTangent1AtGaussPt()

MatrixDouble& MoFEM::FaceElementForcesAndSourcesCore::UserDataOperator::getTangent1AtGaussPt ( )

if higher order geometry return tangent vector to triangle at Gauss pts.

Note: returned matrix has size 0 in rows and columns if no HO approximation of geometry is avaliable.

Definition at line 276 of file FaceElementForcesAndSourcesCore.hpp.

◆ getTangent2()

VectorDouble& MoFEM::FaceElementForcesAndSourcesCore::UserDataOperator::getTangent2 ( )

◆ getTangent2AtGaussPt()

MatrixDouble& MoFEM::FaceElementForcesAndSourcesCore::UserDataOperator::getTangent2AtGaussPt ( )

if higher order geometry return tangent vector to triangle at Gauss pts.

Note: returned matrix has size 0 in rows and columns if no HO approximation of geometry is avaliable.

Definition at line 288 of file FaceElementForcesAndSourcesCore.hpp.

◆ getTensor1Coords()

FTensor::Tensor1<double *, 3> MoFEM::FaceElementForcesAndSourcesCore::UserDataOperator::getTensor1Coords ( )

get get coords at gauss points

t_center(i) = 0;
for(int nn = 0;nn!=3;nn++) {
t_center(i) += t_coords(i);
++t_coords;
}
t_center(i) /= 3;

Definition at line 187 of file FaceElementForcesAndSourcesCore.hpp.

187  {
188  double *ptr = &*getCoords().data().begin();
189  return FTensor::Tensor1<double *, 3>(ptr, &ptr[1], &ptr[2], 3);
190  }

◆ getTensor1CoordsAtGaussPts()

FTensor::Tensor1<double *, 3> MoFEM::FaceElementForcesAndSourcesCore::UserDataOperator::getTensor1CoordsAtGaussPts ( )

get coordinates at Gauss pts.

Definition at line 218 of file FaceElementForcesAndSourcesCore.hpp.

218  {
219  double *ptr = &*getCoordsAtGaussPts().data().begin();
220  return FTensor::Tensor1<double *, 3>(ptr, &ptr[1], &ptr[2], 3);
221  }
MatrixDouble & getCoordsAtGaussPts()
Gauss points and weight, matrix (nb. of points x 3)

◆ getTensor1HoCoordsAtGaussPts()

FTensor::Tensor1<double *, 3> MoFEM::FaceElementForcesAndSourcesCore::UserDataOperator::getTensor1HoCoordsAtGaussPts ( )

get coordinates at Gauss pts (takes in account ho approx. of geometry)

Definition at line 238 of file FaceElementForcesAndSourcesCore.hpp.

238  {
239  if (getHoCoordsAtGaussPts().size1() == 0 &&
240  getHoCoordsAtGaussPts().size2() != 3) {
241  return getTensor1Coords();
242  }
243  double *ptr = &*getHoCoordsAtGaussPts().data().begin();
244  return FTensor::Tensor1<double *, 3>(ptr, &ptr[1], &ptr[2], 3);
245  }
MatrixDouble & getHoCoordsAtGaussPts()
coordinate at Gauss points (if hierarchical approximation of element geometry)
FTensor::Tensor1< double *, 3 > getTensor1Coords()
get get coords at gauss points

◆ getTensor1Normal()

FTensor::Tensor1<double *, 3> MoFEM::FaceElementForcesAndSourcesCore::UserDataOperator::getTensor1Normal ( )

get normal as tensor

Definition at line 134 of file FaceElementForcesAndSourcesCore.hpp.

134  {
135  double *ptr = &*getNormal().data().begin();
136  return FTensor::Tensor1<double *, 3>(ptr, &ptr[1], &ptr[2]);
137  }

◆ getTensor1NormalsAtGaussPt()

FTensor::Tensor1<double *, 3> MoFEM::FaceElementForcesAndSourcesCore::UserDataOperator::getTensor1NormalsAtGaussPt ( )

get normal at integration points

Example:

double nrm2;
for(int gg = gg!=data.getN().size1();gg++) {
nrm2 = sqrt(t_normal(i)*t_normal(i));
++t_normal;
}

Definition at line 307 of file FaceElementForcesAndSourcesCore.hpp.

307  {
308  double *ptr = &*getNormalsAtGaussPt().data().begin();
309  return FTensor::Tensor1<double *, 3>(ptr, &ptr[1], &ptr[2], 3);
310  }
MatrixDouble & getNormalsAtGaussPt()
if higher order geometry return normals at Gauss pts.

◆ getTensor1Tangent1()

FTensor::Tensor1<double *, 3> MoFEM::FaceElementForcesAndSourcesCore::UserDataOperator::getTensor1Tangent1 ( )

get tangentOne as tensor

Definition at line 141 of file FaceElementForcesAndSourcesCore.hpp.

141  {
142  double *ptr = &*getTangent1().data().begin();
143  return FTensor::Tensor1<double *, 3>(ptr, &ptr[1], &ptr[2]);
144  }

◆ getTensor1Tangent1AtGaussPt()

FTensor::Tensor1<double *, 3> MoFEM::FaceElementForcesAndSourcesCore::UserDataOperator::getTensor1Tangent1AtGaussPt ( )

get tangent 1 at integration points

Definition at line 315 of file FaceElementForcesAndSourcesCore.hpp.

315  {
316  double *ptr = &*getTangent1AtGaussPt().data().begin();
317  return FTensor::Tensor1<double *, 3>(ptr, &ptr[1], &ptr[2], 3);
318  }
MatrixDouble & getTangent1AtGaussPt()
if higher order geometry return tangent vector to triangle at Gauss pts.

◆ getTensor1Tangent2AtGaussPt()

FTensor::Tensor1<double *, 3> MoFEM::FaceElementForcesAndSourcesCore::UserDataOperator::getTensor1Tangent2AtGaussPt ( )

get tangent 2 at integration points

Definition at line 323 of file FaceElementForcesAndSourcesCore.hpp.

323  {
324  double *ptr = &*getTangent2AtGaussPt().data().begin();
325  return FTensor::Tensor1<double *, 3>(ptr, &ptr[1], &ptr[2], 3);
326  }
MatrixDouble & getTangent2AtGaussPt()
if higher order geometry return tangent vector to triangle at Gauss pts.

◆ getTensor2Tangent1()

FTensor::Tensor1<double *, 3> MoFEM::FaceElementForcesAndSourcesCore::UserDataOperator::getTensor2Tangent1 ( )

get tangentTwo as tensor

Definition at line 148 of file FaceElementForcesAndSourcesCore.hpp.

148  {
149  double *ptr = &*getTangent2().data().begin();
150  return FTensor::Tensor1<double *, 3>(ptr, &ptr[1], &ptr[2]);
151  }

◆ getTriFE()

DEPRECATED const FaceElementForcesAndSourcesCore* MoFEM::FaceElementForcesAndSourcesCore::UserDataOperator::getTriFE ( )
Deprecated:
use getFaceFE instead

Definition at line 337 of file FaceElementForcesAndSourcesCore.hpp.

337  {
338  return getFaceFE();
339  }
const FaceElementForcesAndSourcesCore * getFaceFE()
return pointer to Generic Triangle Finite Element object

◆ loopSideVolumes()

MoFEMErrorCode UserDataOperator::loopSideVolumes ( const string &  fe_name,
VolumeElementForcesAndSourcesCoreOnSide method 
)

User call this function to loop over elements on the side of face. This function calls MoFEM::VolumeElementForcesAndSourcesCoreOnSide with is operator to do calculations.

Parameters
fe_nameName of the element
methodFinite element object
Returns
error code

Definition at line 80 of file FaceElementForcesAndSourcesCore.cpp.

81  {
82 
84 
85  const EntityHandle ent = getNumeredEntFiniteElementPtr()->getEnt();
86  const Problem *problem_ptr = getFEMethod()->problemPtr;
87  Range adjacent_volumes;
88  ierr = getFaceFE()->mField.getInterface<BitRefManager>()->getAdjacenciesAny(
89  ent, 3, adjacent_volumes);
90  CHKERRG(ierr);
91  typedef NumeredEntFiniteElement_multiIndex::index<
92  Composite_Name_And_Ent_mi_tag>::type FEByComposite;
93  FEByComposite &numered_fe = (const_cast<NumeredEntFiniteElement_multiIndex &>(
94  problem_ptr->numeredFiniteElements))
95  .get<Composite_Name_And_Ent_mi_tag>();
96 
97  method.feName = fe_name;
98 
99  ierr = method.setFaceFEPtr(getFaceFE());
100  CHKERRG(ierr);
101  ierr = method.copyBasicMethod(*getFEMethod());
102  CHKERRG(ierr);
103  ierr = method.copyKsp(*getFEMethod());
104  CHKERRG(ierr);
105  ierr = method.copySnes(*getFEMethod());
106  CHKERRG(ierr);
107  ierr = method.copyTs(*getFEMethod());
108  CHKERRG(ierr);
109 
110  try {
111  ierr = method.preProcess();
112  CHKERRG(ierr);
113  } catch (const std::exception &ex) {
114  std::ostringstream ss;
115  ss << "throw in method: " << ex.what() << " at line " << __LINE__
116  << " in file " << __FILE__ << std::endl;
117  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, ss.str().c_str());
118  }
119 
120  int nn = 0;
121  method.loopSize = adjacent_volumes.size();
122  for (Range::iterator vit = adjacent_volumes.begin();
123  vit != adjacent_volumes.end(); vit++) {
124  FEByComposite::iterator miit =
125  numered_fe.find(boost::make_tuple(fe_name, *vit));
126  if (miit != numered_fe.end()) {
127  // cerr << **miit << endl;
128  // cerr << &(**miit) << endl;
129  // cerr << (*miit)->getEnt() << endl;
130  method.nInTheLoop = nn++;
131  method.numeredEntFiniteElementPtr = *miit;
132  method.dataPtr = (*miit)->sPtr->data_dofs;
133  method.rowPtr = (*miit)->rows_dofs;
134  method.colPtr = (*miit)->cols_dofs;
135 
136  try {
137  ierr = method();
138  CHKERRG(ierr);
139  } catch (const std::exception &ex) {
140  std::ostringstream ss;
141  ss << "throw in method: " << ex.what() << " at line " << __LINE__
142  << " in file " << __FILE__ << std::endl;
143  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, ss.str().c_str());
144  }
145  }
146  }
147 
148  try {
149  ierr = method.postProcess();
150  CHKERRG(ierr);
151  } catch (const std::exception &ex) {
152  std::ostringstream ss;
153  ss << "throw in method: " << ex.what() << " at line " << __LINE__
154  << " in file " << __FILE__ << std::endl;
155  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, ss.str().c_str());
156  }
157 
159 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:522
boost::shared_ptr< const NumeredEntFiniteElement > getNumeredEntFiniteElementPtr() const
Return raw pointer to NumeredEntFiniteElement.
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:565
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:528
const Problem * problemPtr
const FEMethod * getFEMethod() const
Return raw pointer to Finite Element Method object.
multi_index_container< boost::shared_ptr< NumeredEntFiniteElement >, indexed_by< ordered_unique< tag< Unique_mi_tag >, const_mem_fun< NumeredEntFiniteElement::interface_type_EntFiniteElement, UId, &NumeredEntFiniteElement::getGlobalUniqueId > >, ordered_non_unique< tag< Part_mi_tag >, member< NumeredEntFiniteElement, unsigned int, &NumeredEntFiniteElement::part > >, ordered_non_unique< tag< FiniteElement_name_mi_tag >, const_mem_fun< NumeredEntFiniteElement::interface_type_FiniteElement, boost::string_ref, &NumeredEntFiniteElement::getNameRef > >, ordered_non_unique< tag< Ent_mi_tag >, const_mem_fun< NumeredEntFiniteElement::interface_type_EntFiniteElement, EntityHandle, &NumeredEntFiniteElement::getEnt > >, ordered_non_unique< tag< Composite_Name_And_Ent_mi_tag >, composite_key< NumeredEntFiniteElement, const_mem_fun< NumeredEntFiniteElement::interface_type_FiniteElement, boost::string_ref, &NumeredEntFiniteElement::getNameRef >, const_mem_fun< NumeredEntFiniteElement::interface_type_EntFiniteElement, EntityHandle, &NumeredEntFiniteElement::getEnt > > >, ordered_non_unique< tag< Composite_Name_And_Part_mi_tag >, composite_key< NumeredEntFiniteElement, const_mem_fun< NumeredEntFiniteElement::interface_type_FiniteElement, boost::string_ref, &NumeredEntFiniteElement::getNameRef >, member< NumeredEntFiniteElement, unsigned int, &NumeredEntFiniteElement::part > > > > > NumeredEntFiniteElement_multiIndex
MultiIndex for entities for NumeredEntFiniteElement.
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Common.hpp:80
MoFEMErrorCode getInterface(const MOFEMuuid &uuid, IFACE *&iface) const
Get interface by uuid and return reference to pointer of interface.
const FaceElementForcesAndSourcesCore * getFaceFE()
return pointer to Generic Triangle Finite Element object

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