v0.14.0
Functions | Variables
forces_and_sources_getting_mult_H1_H1_atom_test.cpp File Reference

Atom test verifying forces and sources operator on H1 approx. space. More...

#include <MoFEM.hpp>

Go to the source code of this file.

Functions

template<typename T >
void zero_entries (T &t)
 
int main (int argc, char *argv[])
 

Variables

static char help [] = "...\n\n"
 
constexpr double eps = 1e-6
 

Detailed Description

Atom test verifying forces and sources operator on H1 approx. space.

Definition in file forces_and_sources_getting_mult_H1_H1_atom_test.cpp.

Function Documentation

◆ main()

int main ( int  argc,
char *  argv[] 
)

Definition at line 25 of file forces_and_sources_getting_mult_H1_H1_atom_test.cpp.

25  {
26 
27  MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
28 
29  try {
30 
31  moab::Core mb_instance;
32  moab::Interface &moab = mb_instance;
33  int rank;
34  MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
35 
36  PetscBool flg = PETSC_TRUE;
37  char mesh_file_name[255];
38 #if PETSC_VERSION_GE(3, 6, 4)
39  CHKERR PetscOptionsGetString(PETSC_NULL, "", "-my_file", mesh_file_name,
40  255, &flg);
41 #else
42  CHKERR PetscOptionsGetString(PETSC_NULL, PETSC_NULL, "-my_file",
43  mesh_file_name, 255, &flg);
44 #endif
45  if (flg != PETSC_TRUE) {
46  SETERRQ(PETSC_COMM_SELF, 1, "*** ERROR -my_file (MESH FILE NEEDED)");
47  }
48 
49  // Create MoFEM database
50  MoFEM::Core core(moab);
51  MoFEM::Interface &m_field = core;
52 
53  const char *option;
54  option = "";
55  ;
56  CHKERR moab.load_file(mesh_file_name, 0, option);
57 
58  // set entitities bit level
59  BitRefLevel bit_level0;
60  bit_level0.set(0);
61  EntityHandle meshset_level0;
62  CHKERR moab.create_meshset(MESHSET_SET, meshset_level0);
63  CHKERRG(rval);
64  CHKERR m_field.getInterface<BitRefManager>()->setBitRefLevelByDim(
65  0, 3, bit_level0);
66 
67  // Fields
68  CHKERR m_field.add_field("FIELD1", H1, AINSWORTH_LEGENDRE_BASE, 1);
69  CHKERR m_field.add_field("FIELD2", H1, AINSWORTH_LEGENDRE_BASE, 3);
70 
71  // FE
72  CHKERR m_field.add_finite_element("TEST_FE");
73 
74  // Define rows/cols and element data
75  CHKERR m_field.modify_finite_element_add_field_row("TEST_FE", "FIELD1");
76  CHKERR m_field.modify_finite_element_add_field_col("TEST_FE", "FIELD2");
77  CHKERR m_field.modify_finite_element_add_field_data("TEST_FE", "FIELD1");
78  CHKERR m_field.modify_finite_element_add_field_data("TEST_FE", "FIELD2");
79 
80  // Problem
81  CHKERR m_field.add_problem("TEST_PROBLEM");
82 
83  // set finite elements for problem
84  CHKERR m_field.modify_problem_add_finite_element("TEST_PROBLEM", "TEST_FE");
85  // set refinement level for problem
86  CHKERR m_field.modify_problem_ref_level_add_bit("TEST_PROBLEM", bit_level0);
87 
88  // meshset consisting all entities in mesh
89  EntityHandle root_set = moab.get_root_set();
90  // add entities to field
91  CHKERR m_field.add_ents_to_field_by_type(root_set, MBTET, "FIELD1");
92  CHKERR m_field.add_ents_to_field_by_type(root_set, MBTET, "FIELD2");
93  // add entities to finite element
94  CHKERR m_field.add_ents_to_finite_element_by_type(root_set, MBTET,
95  "TEST_FE");
96 
97  // set app. order
98  // see Hierarchic Finite Element Bases on Unstructured Tetrahedral Meshes
99  // (Mark Ainsworth & Joe Coyle)
100  int order = 5;
101  CHKERR m_field.set_field_order(root_set, MBTET, "FIELD1", order);
102  CHKERR m_field.set_field_order(root_set, MBTRI, "FIELD1", order);
103  CHKERR m_field.set_field_order(root_set, MBEDGE, "FIELD1", order);
104  CHKERR m_field.set_field_order(root_set, MBVERTEX, "FIELD1", 1);
105  CHKERR m_field.set_field_order(root_set, MBTET, "FIELD2", order);
106  CHKERR m_field.set_field_order(root_set, MBTRI, "FIELD2", order);
107  CHKERR m_field.set_field_order(root_set, MBEDGE, "FIELD2", order);
108  CHKERR m_field.set_field_order(root_set, MBVERTEX, "FIELD2", 1);
109 
110  /****/
111  // build database
112  // build field
113  CHKERR m_field.build_fields();
114  // build finite elemnts
115  CHKERR m_field.build_finite_elements();
116  // build adjacencies
117  CHKERR m_field.build_adjacencies(bit_level0);
118  // build problem
119  //
120  ProblemsManager *prb_mng_ptr;
121  CHKERR m_field.getInterface(prb_mng_ptr);
122  // const Problem_multiIndex *problems_ptr;
123  CHKERR prb_mng_ptr->buildProblem("TEST_PROBLEM", true);
124 
125  /****/
126  // mesh partitioning
127  // partition
128  CHKERR prb_mng_ptr->partitionSimpleProblem("TEST_PROBLEM");
129  CHKERR prb_mng_ptr->partitionFiniteElements("TEST_PROBLEM");
130  // what are ghost nodes, see Petsc Manual
131  CHKERR prb_mng_ptr->partitionGhostDofs("TEST_PROBLEM");
132 
133  struct ForcesAndSourcesCore_TestFE : public ForcesAndSourcesCore {
134 
135  typedef tee_device<std::ostream, std::ofstream> TeeDevice;
136  typedef stream<TeeDevice> TeeStream;
137 
138  struct my_mult_H1_H1 : public DataOperator {
139 
140  std::ofstream ofs;
141  TeeDevice my_tee;
142  TeeStream my_split;
143 
144  my_mult_H1_H1()
145  : ofs("forces_and_sources_getting_mult_H1_H1_atom_test.txt"),
146  my_tee(std::cout, ofs), my_split(my_tee){};
147 
148  ~my_mult_H1_H1() { my_split.close(); }
149 
150  ublas::matrix<FieldData> NN;
151 
152  MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type,
153  EntityType col_type,
154  EntitiesFieldData::EntData &row_data,
155  EntitiesFieldData::EntData &col_data) {
157 
158  row_data.getBase() = AINSWORTH_LEGENDRE_BASE;
159  col_data.getBase() = AINSWORTH_LEGENDRE_BASE;
160  int nb_row_dofs = row_data.getN().size2();
161  int nb_col_dofs = col_data.getN().size2();
162 
163  my_split << row_side << " " << col_side << " " << row_type << " "
164  << col_type << std::endl;
165  my_split << "nb_row_dofs " << nb_row_dofs << " nb_col_dofs "
166  << nb_col_dofs << std::endl;
167  NN.resize(nb_row_dofs, nb_col_dofs);
168 
169  zero_entries(row_data.getN().data());
170  zero_entries(row_data.getDiffN().data());
171 
172  my_split << std::setprecision(3);
173  my_split << std::fixed;
174  my_split << row_data.getN() << std::endl;
175  my_split << col_data.getN() << std::endl;
176 
177  for (unsigned int gg = 0; gg < row_data.getN().size1(); gg++) {
178 
179  bzero(&*NN.data().begin(),
180  nb_row_dofs * nb_col_dofs * sizeof(FieldData));
181 
182  cblas_dger(CblasRowMajor, nb_row_dofs, nb_col_dofs, 1,
183  &row_data.getN()(gg, 0), 1, &col_data.getN()(gg, 0), 1,
184  &*NN.data().begin(), nb_col_dofs);
185 
186  my_split << "gg " << gg << " : ";
187  my_split << std::setprecision(3);
188  my_split << std::fixed;
189 
190  MatrixDouble difference =
191  NN - outer_prod(row_data.getN(gg), col_data.getN(gg));
192  zero_entries(difference.data());
193 
194  my_split << difference << std::endl;
195  if (row_type != MBVERTEX) {
196  my_split << row_data.getDiffN(gg) << std::endl;
197  }
198 
199  if (row_type == MBVERTEX) {
200  my_split << row_data.getDiffN() << std::endl;
201  } else {
202  typedef ublas::array_adaptor<FieldData> storage_t;
203  storage_t st(nb_row_dofs * 3, &row_data.getDiffN()(gg, 0));
204  ublas::matrix<FieldData, ublas::row_major, storage_t>
205  digNatGaussPt(nb_row_dofs, 3, st);
206  my_split << std::endl << digNatGaussPt << std::endl;
207  }
208  }
209 
210  my_split << std::endl;
211 
213  }
214  };
215 
216  my_mult_H1_H1 op;
217 
218  ForcesAndSourcesCore_TestFE(MoFEM::Interface &_m_field)
219  : ForcesAndSourcesCore(_m_field), data_row(MBTET), data_col(MBTET){};
220 
224  }
225 
226  EntitiesFieldData data_row, data_col;
227 
230 
233 
234  CHKERR getEntitySense<MBEDGE>(data_row);
235  CHKERR getEntitySense<MBTRI>(data_row);
236  CHKERR getEntitySense<MBEDGE>(data_col);
237  CHKERR getEntitySense<MBTRI>(data_col);
238 
239  CHKERR getEntityDataOrder<MBEDGE>(data_row, H1);
240  CHKERR getEntityDataOrder<MBEDGE>(data_col, H1);
241  CHKERR getEntityDataOrder<MBTRI>(data_row, H1);
242  CHKERR getEntityDataOrder<MBTRI>(data_col, H1);
243  CHKERR getEntityDataOrder<MBTET>(data_row, H1);
244  CHKERR getEntityDataOrder<MBTET>(data_col, H1);
245  data_row.dataOnEntities[MBVERTEX][0].getBase() =
247 
248  const auto bn1 = mField.get_field_bit_number("FIELD1");
249  const auto bn2 = mField.get_field_bit_number("FIELD2");
250 
251  CHKERR getEntityFieldData(data_row, bn1, MBEDGE);
252  data_col.dataOnEntities[MBVERTEX][0].getBase() =
254  CHKERR getEntityFieldData(data_col, bn2, MBEDGE);
255  CHKERR getRowNodesIndices(data_row, bn1);
256  CHKERR getColNodesIndices(data_col, bn2);
257  CHKERR getEntityRowIndices(data_row, bn1, MBEDGE);
258  CHKERR getEntityColIndices(data_col, bn2, MBEDGE);
259  CHKERR getFaceNodes(data_row);
260  CHKERR getFaceNodes(data_col);
261 
262  MatrixDouble gauss_pts(4, 4);
263  for (int gg = 0; gg < 4; gg++) {
264  gauss_pts(0, gg) = G_TET_X4[gg];
265  gauss_pts(1, gg) = G_TET_Y4[gg];
266  gauss_pts(2, gg) = G_TET_Z4[gg];
267  gauss_pts(3, gg) = G_TET_W4[gg];
268  }
270  gauss_pts,
271  boost::shared_ptr<BaseFunctionCtx>(new EntPolynomialBaseCtx(
272  data_row, H1, AINSWORTH_LEGENDRE_BASE)));
274  gauss_pts,
275  boost::shared_ptr<BaseFunctionCtx>(new EntPolynomialBaseCtx(
276  data_col, H1, AINSWORTH_LEGENDRE_BASE)));
277 
278  CHKERR op.opLhs(data_row, data_col);
279 
281  }
282 
285 
287  }
288  };
289 
290  ForcesAndSourcesCore_TestFE fe1(m_field);
291  CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "TEST_FE", fe1);
292  }
293  CATCH_ERRORS;
294 
296 
297  return 0;
298 }

◆ zero_entries()

template<typename T >
void zero_entries ( T &  t)

Definition at line 19 of file forces_and_sources_getting_mult_H1_H1_atom_test.cpp.

19  {
20  for (auto &v : t)
21  if (std::abs(v) < eps)
22  v = 0;
23 }

Variable Documentation

◆ eps

constexpr double eps = 1e-6
constexpr

◆ help

char help[] = "...\n\n"
static
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:447
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.
Definition: UnknownInterface.hpp:93
MoFEM::ForcesAndSourcesCore::getEntityRowIndices
MoFEMErrorCode getEntityRowIndices(EntitiesFieldData &data, const int bit_number, const EntityType type_lo=MBVERTEX, const EntityType type_hi=MBPOLYHEDRON) const
Definition: ForcesAndSourcesCore.cpp:412
MoFEM::EntitiesFieldData::EntData
Data on single entity (This is passed as argument to DataOperator::doWork)
Definition: EntitiesFieldData.hpp:127
TeeStream
stream< TeeDevice > TeeStream
Definition: forces_and_sources_testing_contact_prism_element.cpp:10
MoFEM::CoreInterface::loop_finite_elements
virtual MoFEMErrorCode loop_finite_elements(const std::string problem_name, const std::string &fe_name, FEMethod &method, boost::shared_ptr< NumeredEntFiniteElement_multiIndex > fe_ptr=nullptr, MoFEMTypes bh=MF_EXIST, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr(), int verb=DEFAULT_VERBOSITY)=0
Make a loop over finite elements.
MoFEM::ProblemsManager::buildProblem
MoFEMErrorCode buildProblem(const std::string name, const bool square_matrix, int verb=VERBOSE)
build problem data structures
Definition: ProblemsManager.cpp:87
MoFEM::CoreTmp< 0 >
Core (interface) class.
Definition: Core.hpp:82
eps
constexpr double eps
Definition: forces_and_sources_getting_mult_H1_H1_atom_test.cpp:18
H1
@ H1
continuous field
Definition: definitions.h:85
MoFEM::DataOperator
base operator to do operations at Gauss Pt. level
Definition: DataOperators.hpp:24
EntityHandle
MoFEM::ProblemsManager
Problem manager is used to build and partition problems.
Definition: ProblemsManager.hpp:21
MoFEM::EntPolynomialBaseCtx
Class used to pass element data to calculate base functions on tet,triangle,edge.
Definition: EntPolynomialBaseCtx.hpp:22
G_TET_Z4
static const double G_TET_Z4[]
Definition: fem_tools.h:1120
MoFEM::ForcesAndSourcesCore::getFaceNodes
MoFEMErrorCode getFaceNodes(EntitiesFieldData &data) const
Get nodes on faces.
Definition: ForcesAndSourcesCore.cpp:921
MoFEM::CoreInterface::modify_finite_element_add_field_row
virtual MoFEMErrorCode modify_finite_element_add_field_row(const std::string &fe_name, const std::string name_row)=0
set field row which finite element use
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
MoFEM::Types::MatrixDouble
UBlasMatrix< double > MatrixDouble
Definition: Types.hpp:77
G_TET_Y4
static const double G_TET_Y4[]
Definition: fem_tools.h:1118
MoFEM::CoreTmp< 0 >::Finalize
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition: Core.cpp:112
MoFEM::EntitiesFieldData::EntData::getDiffN
MatrixDouble & getDiffN(const FieldApproximationBase base)
get derivatives of base functions
Definition: EntitiesFieldData.hpp:1316
MoFEM::ForcesAndSourcesCore::preProcess
virtual MoFEMErrorCode preProcess()
function is run at the beginning of loop
Definition: ForcesAndSourcesCore.cpp:1958
MoFEM::CoreInterface::add_ents_to_field_by_type
virtual MoFEMErrorCode add_ents_to_field_by_type(const Range &ents, const EntityType type, const std::string &name, int verb=DEFAULT_VERBOSITY)=0
Add entities to field meshset.
order
constexpr int order
Definition: dg_projection.cpp:18
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
TeeDevice
tee_device< std::ostream, std::ofstream > TeeDevice
Definition: forces_and_sources_testing_contact_prism_element.cpp:9
MoFEM::ForcesAndSourcesCore::getEntityFieldData
MoFEMErrorCode getEntityFieldData(EntitiesFieldData &data, const int bit_number, const EntityType type_lo=MBVERTEX, const EntityType type_hi=MBPOLYHEDRON) const
Definition: ForcesAndSourcesCore.cpp:770
MoFEM::Interface
DeprecatedCoreInterface Interface
Definition: Interface.hpp:1975
MoFEM::CoreInterface::add_ents_to_finite_element_by_type
virtual MoFEMErrorCode add_ents_to_finite_element_by_type(const EntityHandle entities, const EntityType type, const std::string &name, const bool recursive=true)=0
add entities to finite element
MoFEM::Exceptions::rval
static MoFEMErrorCodeGeneric< moab::ErrorCode > rval
Definition: Exceptions.hpp:74
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:535
MoFEM::CoreInterface::add_finite_element
virtual MoFEMErrorCode add_finite_element(const std::string &fe_name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
add finite element
MoFEM::CoreInterface::modify_finite_element_add_field_col
virtual MoFEMErrorCode modify_finite_element_add_field_col(const std::string &fe_name, const std::string name_row)=0
set field col which finite element use
MoFEM::CoreInterface::build_finite_elements
virtual MoFEMErrorCode build_finite_elements(int verb=DEFAULT_VERBOSITY)=0
Build finite elements.
MoFEM::ForcesAndSourcesCore::getColNodesIndices
MoFEMErrorCode getColNodesIndices(EntitiesFieldData &data, const int bit_number) const
get col node indices from FENumeredDofEntity_multiIndex
Definition: ForcesAndSourcesCore.cpp:328
MoFEM::CoreInterface::add_field
virtual MoFEMErrorCode add_field(const std::string &name, const FieldSpace space, const FieldApproximationBase base, const FieldCoefficientsNumber nb_of_coefficients, const TagType tag_type=MB_TAG_SPARSE, const enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
Add field.
double
MoFEM::EntitiesFieldData::EntData::getBase
FieldApproximationBase & getBase()
Get approximation base.
Definition: EntitiesFieldData.hpp:1300
MoFEM::ForcesAndSourcesCore::getEntityColIndices
MoFEMErrorCode getEntityColIndices(EntitiesFieldData &data, const int bit_number, const EntityType type_lo=MBVERTEX, const EntityType type_hi=MBPOLYHEDRON) const
Definition: ForcesAndSourcesCore.cpp:427
MoFEM::ForcesAndSourcesCore::operator()
virtual MoFEMErrorCode operator()()
function is run for every finite element
Definition: ForcesAndSourcesCore.cpp:1966
t
constexpr double t
plate stiffness
Definition: plate.cpp:59
MoFEM::ProblemsManager::partitionFiniteElements
MoFEMErrorCode partitionFiniteElements(const std::string name, bool part_from_moab=false, int low_proc=-1, int hi_proc=-1, int verb=VERBOSE)
partition finite elements
Definition: ProblemsManager.cpp:2167
MoFEM::CoreInterface::modify_problem_ref_level_add_bit
virtual MoFEMErrorCode modify_problem_ref_level_add_bit(const std::string &name_problem, const BitRefLevel &bit)=0
add ref level to problem
MoFEM::ForcesAndSourcesCore::getSpacesAndBaseOnEntities
MoFEMErrorCode getSpacesAndBaseOnEntities(EntitiesFieldData &data) const
Get field approximation space and base on entities.
Definition: ForcesAndSourcesCore.cpp:990
mesh_file_name
char mesh_file_name[255]
Definition: mesh_smoothing.cpp:23
MoFEM::ForcesAndSourcesCore::getRowNodesIndices
MoFEMErrorCode getRowNodesIndices(EntitiesFieldData &data, const int bit_number) const
get row node indices from FENumeredDofEntity_multiIndex
Definition: ForcesAndSourcesCore.cpp:311
MoFEM::ForcesAndSourcesCore
structure to get information form mofem into EntitiesFieldData
Definition: ForcesAndSourcesCore.hpp:22
G_TET_X4
static const double G_TET_X4[]
Definition: fem_tools.h:1116
G_TET_W4
static const double G_TET_W4[]
Definition: fem_tools.h:1122
v
const double v
phase velocity of light in medium (cm/ns)
Definition: initial_diffusion.cpp:40
MoFEM::CoreTmp< 0 >::Initialize
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
Definition: Core.cpp:72
CATCH_ERRORS
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:372
MoFEM::EntitiesFieldData::EntData::getN
MatrixDouble & getN(const FieldApproximationBase base)
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....
Definition: EntitiesFieldData.hpp:1305
MoFEM::Core
CoreTmp< 0 > Core
Definition: Core.hpp:1094
std
Definition: enable_if.hpp:5
MoFEM::ProblemsManager::partitionSimpleProblem
MoFEMErrorCode partitionSimpleProblem(const std::string name, int verb=VERBOSE)
partition problem dofs
Definition: ProblemsManager.cpp:1537
zero_entries
void zero_entries(T &t)
Definition: forces_and_sources_getting_mult_H1_H1_atom_test.cpp:19
AINSWORTH_LEGENDRE_BASE
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:60
MoFEM::PetscOptionsGetString
PetscErrorCode PetscOptionsGetString(PetscOptions *, const char pre[], const char name[], char str[], size_t size, PetscBool *set)
Definition: DeprecatedPetsc.hpp:172
MoFEM::CoreInterface::build_fields
virtual MoFEMErrorCode build_fields(int verb=DEFAULT_VERBOSITY)=0
MoFEM::CoreInterface::modify_finite_element_add_field_data
virtual MoFEMErrorCode modify_finite_element_add_field_data(const std::string &fe_name, const std::string name_filed)=0
set finite element field data
MoFEM::TetPolynomialBase::getValue
MoFEMErrorCode getValue(MatrixDouble &pts, boost::shared_ptr< BaseFunctionCtx > ctx_ptr)
Definition: TetPolynomialBase.cpp:1348
MoFEM::EntitiesFieldData::dataOnEntities
std::array< boost::ptr_vector< EntData >, MBMAXTYPE > dataOnEntities
Definition: EntitiesFieldData.hpp:56
MoFEM::BitRefManager
Managing BitRefLevels.
Definition: BitRefManager.hpp:21
MoFEM::CoreInterface::modify_problem_add_finite_element
virtual MoFEMErrorCode modify_problem_add_finite_element(const std::string name_problem, const std::string &fe_name)=0
add finite element to problem, this add entities assigned to finite element to a particular problem
MoFEM::Types::BitRefLevel
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition: Types.hpp:40
MoFEM::CoreInterface::build_adjacencies
virtual MoFEMErrorCode build_adjacencies(const Range &ents, int verb=DEFAULT_VERBOSITY)=0
build adjacencies
MoFEMFunctionBeginHot
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:440
MoFEM::CoreInterface::set_field_order
virtual MoFEMErrorCode set_field_order(const EntityHandle meshset, const EntityType type, const std::string &name, const ApproximationOrder order, int verb=DEFAULT_VERBOSITY)=0
Set order approximation of the entities in the field.
MoFEM::EntitiesFieldData
data structure for finite element entity
Definition: EntitiesFieldData.hpp:40
MoFEM::ForcesAndSourcesCore::postProcess
virtual MoFEMErrorCode postProcess()
function is run at the end of loop
Definition: ForcesAndSourcesCore.cpp:1981
MoFEM::ProblemsManager::partitionGhostDofs
MoFEMErrorCode partitionGhostDofs(const std::string name, int verb=VERBOSE)
determine ghost nodes
Definition: ProblemsManager.cpp:2339
MoFEM::CoreInterface::add_problem
virtual MoFEMErrorCode add_problem(const std::string &name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
Add problem.
help
static char help[]
Definition: forces_and_sources_getting_mult_H1_H1_atom_test.cpp:16
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
CHKERRG
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:483
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
MoFEM::TetPolynomialBase
Calculate base functions on tetrahedral.
Definition: TetPolynomialBase.hpp:19