v0.14.0
Functions | Variables
check_base_functions_derivatives_on_tet.cpp File Reference
#include <MoFEM.hpp>

Go to the source code of this file.

Functions

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

Variables

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

Function Documentation

◆ main()

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

Definition at line 14 of file check_base_functions_derivatives_on_tet.cpp.

14  {
15 
16  MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
17 
18  try {
19 
20  // Select space
21  enum spaces { H1TET, HDIVTET, HCURLTET, LASTSPACEOP };
22 
23  const char *list_spaces[] = {"h1tet", "hdivtet", "hcurltet"};
24 
25  PetscBool flg;
26  PetscInt choise_space_value = H1TET;
27  CHKERR PetscOptionsGetEList(PETSC_NULL, NULL, "-space", list_spaces,
28  LASTSPACEOP, &choise_space_value, &flg);
29 
30  if (flg != PETSC_TRUE) {
31  SETERRQ(PETSC_COMM_SELF, MOFEM_IMPOSSIBLE_CASE, "space not set");
32  }
33 
34  FieldSpace space = LASTSPACE;
35  if (choise_space_value == H1TET) {
36  space = H1;
37  } else if (choise_space_value == HDIVTET) {
38  space = HDIV;
39  } else if (choise_space_value == HCURLTET) {
40  space = HCURL;
41  }
42 
43  // Select base
44  enum bases { AINSWORTH, DEMKOWICZ, LASBASETOP };
45 
46  const char *list_bases[] = {"ainsworth", "demkowicz"};
47 
48  PetscInt choice_base_value = AINSWORTH;
49  CHKERR PetscOptionsGetEList(PETSC_NULL, NULL, "-base", list_bases,
50  LASBASETOP, &choice_base_value, &flg);
51 
52  if (flg != PETSC_TRUE) {
53  SETERRQ(PETSC_COMM_SELF, MOFEM_IMPOSSIBLE_CASE, "base not set");
54  }
55 
57  if (choice_base_value == AINSWORTH) {
59  } else if (choice_base_value == DEMKOWICZ) {
60  base = DEMKOWICZ_JACOBI_BASE;
61  }
62 
63  moab::Core mb_instance;
64  moab::Interface &moab = mb_instance;
65  int rank;
66  MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
67 
68  // create one tet
69  double tet_coords[] = {0, 0, 0, 0.5, 0, 0, 0, 0.5, 0, 0, 0, 0.5};
70  EntityHandle nodes[4];
71  for (int nn = 0; nn < 4; nn++) {
72  CHKERR moab.create_vertex(&tet_coords[3 * nn], nodes[nn]);
73  }
74  EntityHandle tet;
75  CHKERR moab.create_element(MBTET, nodes, 4, tet);
76 
77  // Create MoFEM database
78  MoFEM::Core core(moab);
79  MoFEM::Interface &m_field = core;
80 
81  // set entities bit level
82  BitRefLevel bit_level0;
83  bit_level0.set(0);
84  CHKERR m_field.getInterface<BitRefManager>()->setBitRefLevelByDim(
85  0, 3, bit_level0);
86 
87 
88  // Fields
89  CHKERR m_field.add_field("FIELD", space, base, 1);
90 
91  // FE TET
92  CHKERR m_field.add_finite_element("TET_FE");
93 
94  // Define rows/cols and element data
95  CHKERR m_field.modify_finite_element_add_field_row("TET_FE", "FIELD");
96  CHKERR m_field.modify_finite_element_add_field_col("TET_FE", "FIELD");
97  CHKERR m_field.modify_finite_element_add_field_data("TET_FE", "FIELD");
98 
99 
100  // Problem
101  CHKERR m_field.add_problem("TEST_PROBLEM");
102 
103 
104  // set finite elements for problem
105  CHKERR m_field.modify_problem_add_finite_element("TEST_PROBLEM", "TET_FE");
106 
107  // set refinement level for problem
108  CHKERR m_field.modify_problem_ref_level_add_bit("TEST_PROBLEM", bit_level0);
109 
110 
111  // meshset consisting all entities in mesh
112  EntityHandle root_set = moab.get_root_set();
113  // add entities to field
114  CHKERR m_field.add_ents_to_field_by_type(root_set, MBTET, "FIELD");
115 
116  // add entities to finite element
117  CHKERR
118  m_field.add_ents_to_finite_element_by_type(root_set, MBTET, "TET_FE");
119 
120 
121  // set app. order
122  int order = 5;
123  if (space == H1) {
124  CHKERR m_field.set_field_order(root_set, MBTET, "FIELD", order);
125  CHKERR m_field.set_field_order(root_set, MBTRI, "FIELD", order);
126  CHKERR m_field.set_field_order(root_set, MBEDGE, "FIELD", order);
127  CHKERR m_field.set_field_order(root_set, MBVERTEX, "FIELD", 1);
128 
129  }
130  if (space == HCURL) {
131  CHKERR m_field.set_field_order(root_set, MBTET, "FIELD", order);
132  CHKERR m_field.set_field_order(root_set, MBTRI, "FIELD", order);
133  CHKERR m_field.set_field_order(root_set, MBEDGE, "FIELD", order);
134 
135  }
136  if (space == HDIV) {
137  CHKERR m_field.set_field_order(root_set, MBTET, "FIELD", order);
138  CHKERR m_field.set_field_order(root_set, MBTRI, "FIELD", order);
139  }
140 
141  /****/
142  // build database
143  // build field
144  CHKERR m_field.build_fields();
145 
146  // build finite elemnts
147  CHKERR m_field.build_finite_elements();
148 
149  // build adjacencies
150  CHKERR m_field.build_adjacencies(bit_level0);
151 
152 
153  /****/
154  ProblemsManager *prb_mng_ptr;
155  CHKERR m_field.getInterface(prb_mng_ptr);
156 
157  // build problem
158  CHKERR prb_mng_ptr->buildProblem("TEST_PROBLEM", true);
159 
160  // mesh partitioning
161  // partition
162  CHKERR prb_mng_ptr->partitionSimpleProblem("TEST_PROBLEM");
163  CHKERR prb_mng_ptr->partitionFiniteElements("TEST_PROBLEM");
164 
165  // what are ghost nodes, see Petsc Manual
166  CHKERR prb_mng_ptr->partitionGhostDofs("TEST_PROBLEM");
167 
168  typedef tee_device<std::ostream, std::ofstream> TeeDevice;
169  typedef stream<TeeDevice> TeeStream;
170 
171  std::ofstream ofs("forces_and_sources_checking_derivatives.txt");
172  TeeDevice my_tee(std::cout, ofs);
173  TeeStream my_split(my_tee);
174 
175  struct OpCheckingDerivatives
177 
178  TeeStream &mySplit;
179  OpCheckingDerivatives(TeeStream &my_split)
181  "FIELD", UserDataOperator::OPROW),
182  mySplit(my_split) {}
183 
184  MoFEMErrorCode doWork(int side, EntityType type,
187 
188  if (data.getFieldData().size() == 0)
190 
191  mySplit << "Type " << type << " side " << side << endl;
192 
193  if (data.getFieldDofs()[0]->getSpace() == H1) {
194 
195  // mySplit << std::fixed << data.getN() << std::endl;
196  // mySplit << std::fixed << data.getDiffN() << std::endl;
197 
198  const int nb_dofs = data.getN().size2();
199  for (int dd = 0; dd != nb_dofs; dd++) {
200  const double dksi = (data.getN()(1, dd) - data.getN()(0, dd)) / eps;
201  const double deta = (data.getN()(3, dd) - data.getN()(2, dd)) / eps;
202  const double dzeta =
203  (data.getN()(5, dd) - data.getN()(4, dd)) / eps;
204  mySplit << "DKsi " << dksi - data.getDiffN()(6, 3 * dd + 0)
205  << std::endl;
206  mySplit << "DEta " << deta - data.getDiffN()(6, 3 * dd + 1)
207  << std::endl;
208  mySplit << "DZeta " << dzeta - data.getDiffN()(6, 3 * dd + 2)
209  << std::endl;
210  if (fabs(dksi - data.getDiffN()(6, 3 * dd + 0)) > eps_diff) {
211  SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
212  "H1 inconsistent dKsi derivative");
213  }
214  if (fabs(deta - data.getDiffN()(6, 3 * dd + 1)) > eps_diff) {
215  SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
216  "H1 inconsistent dEta derivative");
217  }
218  if (fabs(dzeta - data.getDiffN()(6, 3 * dd + 2)) > eps_diff) {
219  SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
220  "H1 inconsistent dZeta derivative");
221  }
222  }
223  }
224 
225  if (data.getFieldDofs()[0]->getSpace() == HDIV ||
226  data.getFieldDofs()[0]->getSpace() == HCURL) {
227 
229  &data.getN()(0, 0), &data.getN()(0, 1), &data.getN()(0, 2), 3);
231  &data.getN()(1, 0), &data.getN()(1, 1), &data.getN()(1, 2), 3);
233  &data.getN()(2, 0), &data.getN()(2, 1), &data.getN()(2, 2), 3);
235  &data.getN()(3, 0), &data.getN()(3, 1), &data.getN()(3, 2), 3);
236  FTensor::Tensor1<double *, 3> base_zeta_m(
237  &data.getN()(4, 0), &data.getN()(4, 1), &data.getN()(4, 2), 3);
238  FTensor::Tensor1<double *, 3> base_zeta_p(
239  &data.getN()(5, 0), &data.getN()(5, 1), &data.getN()(5, 2), 3);
240 
242  &data.getDiffN()(6, 0), &data.getDiffN()(6, 3),
243  &data.getDiffN()(6, 6), &data.getDiffN()(6, 1),
244  &data.getDiffN()(6, 4), &data.getDiffN()(6, 7),
245  &data.getDiffN()(6, 2), &data.getDiffN()(6, 5),
246  &data.getDiffN()(6, 8), 9);
247 
252 
253  const int nb_dofs = data.getN().size2() / 3;
254  for (int dd = 0; dd != nb_dofs; dd++) {
256  dksi(i) = (base_ksi_p(i) - base_ksi_m(i)) / eps;
258  deta(i) = (base_eta_p(i) - base_eta_m(i)) / eps;
260  dzeta(i) = (base_zeta_p(i) - base_zeta_m(i)) / eps;
261 
262  dksi(i) -= diff_base(i, N0);
263  deta(i) -= diff_base(i, N1);
264  dzeta(i) -= diff_base(i, N2);
265 
266  mySplit << "dKsi " << dksi(0) << " " << dksi(1) << " " << dksi(2)
267  << " " << sqrt(dksi(i) * dksi(i)) << endl;
268  mySplit << "dEta " << deta(0) << " " << deta(1) << " " << deta(2)
269  << " " << sqrt(deta(i) * deta(i)) << endl;
270  mySplit << "dZeta " << dzeta(0) << " " << dzeta(1) << " "
271  << dzeta(2) << " " << sqrt(dzeta(i) * dzeta(i)) << endl;
272 
273  if (sqrt(dksi(i) * dksi(i)) > eps_diff) {
274  SETERRQ2(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
275  "%s inconsistent dKsi derivative for type %d",
276  FieldSpaceNames[data.getFieldDofs()[0]->getSpace()],
277  type);
278  }
279  if (sqrt(deta(i) * deta(i)) > eps_diff) {
280  SETERRQ2(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
281  "%s inconsistent dEta derivative for type %d",
282  FieldSpaceNames[data.getFieldDofs()[0]->getSpace()],
283  type);
284  }
285  if (sqrt(dzeta(i) * dzeta(i)) > eps_diff) {
286  SETERRQ2(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
287  "%s inconsistent dZeta derivative for type %d",
288  FieldSpaceNames[data.getFieldDofs()[0]->getSpace()],
289  type);
290  }
291 
292  ++base_ksi_m;
293  ++base_ksi_p;
294  ++base_eta_m;
295  ++base_eta_p;
296  ++base_zeta_m;
297  ++base_zeta_p;
298  ++diff_base;
299  }
300  }
301 
303  }
304  };
305 
306  struct MyFE : public VolumeElementForcesAndSourcesCore {
307 
308  MyFE(MoFEM::Interface &m_field)
310  int getRule(int order) { return -1; };
311 
314 
315  const double ksi = G_TET_X1[0];
316  const double eta = G_TET_Y1[0];
317  const double zeta = G_TET_Z1[0];
318 
319  gaussPts.resize(4, 7);
320  gaussPts.clear();
321 
322  gaussPts(0, 0) = ksi - eps;
323  gaussPts(1, 0) = eta;
324  gaussPts(2, 0) = zeta;
325  gaussPts(0, 1) = ksi + eps;
326  gaussPts(1, 1) = eta;
327  gaussPts(2, 1) = zeta;
328 
329  gaussPts(0, 2) = ksi;
330  gaussPts(1, 2) = eta - eps;
331  gaussPts(2, 2) = zeta;
332  gaussPts(0, 3) = ksi;
333  gaussPts(1, 3) = eta + eps;
334  gaussPts(2, 3) = zeta;
335 
336  gaussPts(0, 4) = ksi;
337  gaussPts(1, 4) = eta;
338  gaussPts(2, 4) = zeta - eps;
339  gaussPts(0, 5) = ksi;
340  gaussPts(1, 5) = eta;
341  gaussPts(2, 5) = zeta + eps;
342 
343  gaussPts(0, 6) = ksi;
344  gaussPts(1, 6) = eta;
345  gaussPts(2, 6) = zeta;
346 
347  for (unsigned int ii = 0; ii != gaussPts.size2(); ii++) {
348  gaussPts(3, ii) = 1;
349  }
350 
351  cerr << gaussPts << endl;
352 
354  }
355  };
356 
357  MyFE tet_fe(m_field);
358 
359  tet_fe.getOpPtrVector().push_back(new OpCheckingDerivatives(my_split));
360  CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "TET_FE", tet_fe);
361 
362  }
363  CATCH_ERRORS;
364 
366 
367 }

Variable Documentation

◆ eps

const double eps = 1e-6
static

◆ eps_diff

const double eps_diff = 1e-8
static

Definition at line 12 of file check_base_functions_derivatives_on_tet.cpp.

◆ help

char help[] = "...\n\n"
static

Definition at line 9 of file check_base_functions_derivatives_on_tet.cpp.

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::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
H1
@ H1
continuous field
Definition: definitions.h:85
FTensor::Tensor1
Definition: Tensor1_value.hpp:8
EntityHandle
MoFEM::ProblemsManager
Problem manager is used to build and partition problems.
Definition: ProblemsManager.hpp:21
eps_diff
static const double eps_diff
Definition: check_base_functions_derivatives_on_tet.cpp:12
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
NOBASE
@ NOBASE
Definition: definitions.h:59
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
MoFEM::ForcesAndSourcesCore::setGaussPts
virtual MoFEMErrorCode setGaussPts(int order_row, int order_col, int order_data)
set user specific integration rule
Definition: ForcesAndSourcesCore.cpp:1926
MoFEM::EntitiesFieldData::EntData::getFieldData
const VectorDouble & getFieldData() const
get dofs values
Definition: EntitiesFieldData.hpp:1241
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::DataOperator::doWork
virtual MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
Operator for bi-linear form, usually to calculate values on left hand side.
Definition: DataOperators.hpp:40
MOFEM_IMPOSSIBLE_CASE
@ MOFEM_IMPOSSIBLE_CASE
Definition: definitions.h:35
zeta
double zeta
Viscous hardening.
Definition: plastic.cpp:177
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.
FTensor::Tensor2< double *, 3, 3 >
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::Interface
DeprecatedCoreInterface Interface
Definition: Interface.hpp:1975
FieldSpace
FieldSpace
approximation spaces
Definition: definitions.h:82
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
G_TET_Z1
static const double G_TET_Z1[]
Definition: fem_tools.h:1114
FTensor::Number< 0 >
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::UserDataOperator
Definition: ForcesAndSourcesCore.hpp:549
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.
convert.type
type
Definition: convert.py:64
G_TET_Y1
static const double G_TET_Y1[]
Definition: fem_tools.h:1113
eta
double eta
Definition: free_surface.cpp:170
G_TET_X1
static const double G_TET_X1[]
Definition: fem_tools.h:1112
LASTSPACE
@ LASTSPACE
FieldSpace in [ 0, LASTSPACE )
Definition: definitions.h:89
MoFEM::ForcesAndSourcesCore::getRule
virtual int getRule(int order_row, int order_col, int order_data)
another variant of getRule
Definition: ForcesAndSourcesCore.cpp:1920
MoFEM::EntitiesFieldData::EntData::getFieldDofs
const VectorDofs & getFieldDofs() const
get dofs data stature FEDofEntity
Definition: EntitiesFieldData.hpp:1256
MoFEM::VolumeElementForcesAndSourcesCore
Volume finite element base.
Definition: VolumeElementForcesAndSourcesCore.hpp:26
i
FTensor::Index< 'i', SPACE_DIM > i
Definition: hcurl_divergence_operator_2d.cpp:27
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::VolumeElementForcesAndSourcesCore::UserDataOperator
Definition: VolumeElementForcesAndSourcesCore.hpp:108
FTensor::Index< 'i', 3 >
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
FTensor::dd
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
FieldSpaceNames
const static char *const FieldSpaceNames[]
Definition: definitions.h:92
CATCH_ERRORS
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:372
DEMKOWICZ_JACOBI_BASE
@ DEMKOWICZ_JACOBI_BASE
Definition: definitions.h:66
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
MoFEM::ProblemsManager::partitionSimpleProblem
MoFEMErrorCode partitionSimpleProblem(const std::string name, int verb=VERBOSE)
partition problem dofs
Definition: ProblemsManager.cpp:1537
eps
static const double eps
Definition: check_base_functions_derivatives_on_tet.cpp:11
help
static char help[]
Definition: check_base_functions_derivatives_on_tet.cpp:9
AINSWORTH_LEGENDRE_BASE
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:60
HCURL
@ HCURL
field with continuous tangents
Definition: definitions.h:86
MoFEM::ForcesAndSourcesCore::gaussPts
MatrixDouble gaussPts
Matrix of integration points.
Definition: ForcesAndSourcesCore.hpp:109
MoFEM::PetscOptionsGetEList
PetscErrorCode PetscOptionsGetEList(PetscOptions *, const char pre[], const char name[], const char *const *list, PetscInt next, PetscInt *value, PetscBool *set)
Definition: DeprecatedPetsc.hpp:203
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
FieldApproximationBase
FieldApproximationBase
approximation base
Definition: definitions.h:58
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_ATOM_TEST_INVALID
@ MOFEM_ATOM_TEST_INVALID
Definition: definitions.h: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::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.
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
HDIV
@ HDIV
field with continuous normal traction
Definition: definitions.h:87
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346