v0.14.0
hdiv_continuity_check.cpp
Go to the documentation of this file.
1 
2 
3 #include <MoFEM.hpp>
4 
5 namespace bio = boost::iostreams;
6 using bio::stream;
7 using bio::tee_device;
8 
9 using namespace MoFEM;
10 
11 static char help[] = "...\n\n";
12 
13 static const double face_coords[4][9] = {{0, 0, 0, 1, 0, 0, 0, 0, 1},
14 
15  {1, 0, 0, 0, 1, 0, 0, 0, 1},
16 
17  {0, 0, 0, 0, 1, 0, 0, 0, 1},
18 
19  {0, 0, 0, 1, 0, 0, 0, 1, 0}};
20 
21 int main(int argc, char *argv[]) {
22 
23  MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
24 
25  try {
26 
27  enum bases { HDIV_AINSWORTH, HDIV_DEMKOWICZ, LASTOP };
28 
29  const char *list[] = {"hdiv_ainsworth", "hdiv_demkowicz"};
30 
31  PetscBool flg;
32  PetscInt choise_value = HDIV_AINSWORTH;
33  CHKERR PetscOptionsGetEList(PETSC_NULL, NULL, "-base", list, LASTOP,
34  &choise_value, &flg);
35  if (flg != PETSC_TRUE) {
36  SETERRQ(PETSC_COMM_SELF, MOFEM_IMPOSSIBLE_CASE, "base not set");
37  }
38 
39  moab::Core mb_instance;
40  moab::Interface &moab = mb_instance;
41  int rank;
42  MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
43 
44  char mesh_file_name[255];
45 #if PETSC_VERSION_GE(3, 6, 4)
46  CHKERR PetscOptionsGetString(PETSC_NULL, "", "-my_file", mesh_file_name,
47  255, &flg);
48 #else
49  CHKERR PetscOptionsGetString(PETSC_NULL, PETSC_NULL, "-my_file",
50  mesh_file_name, 255, &flg);
51 #endif
52  if (flg != PETSC_TRUE) {
53  SETERRQ(PETSC_COMM_SELF, 1, "*** ERROR -my_file (MESH FILE NEEDED)");
54  }
55  const char *option;
56  option = "";
57  CHKERR moab.load_file(mesh_file_name, 0, option);
58 
59  MoFEM::Core core(moab);
60  MoFEM::Interface &m_field = core;
61 
62  // set entities bit level
63  BitRefLevel bit_level0;
64  bit_level0.set(0);
65  EntityHandle meshset_level0;
66  CHKERR moab.create_meshset(MESHSET_SET, meshset_level0);
67  CHKERR m_field.getInterface<BitRefManager>()->setBitRefLevelByDim(
68  0, 3, bit_level0);
69 
70  // Fields
71  CHKERR m_field.add_field("MESH_NODE_POSITIONS", H1, AINSWORTH_LEGENDRE_BASE,
72  3);
73  switch (choise_value) {
74  case HDIV_AINSWORTH:
75  CHKERR m_field.add_field("HDIV", HDIV, AINSWORTH_LEGENDRE_BASE, 1);
76  break;
77  case HDIV_DEMKOWICZ:
78  CHKERR m_field.add_field("HDIV", HDIV, DEMKOWICZ_JACOBI_BASE, 1);
79  break;
80  }
81 
82  // FE
83  CHKERR m_field.add_finite_element("TET_FE");
84  CHKERR m_field.add_finite_element("TRI_FE");
85  CHKERR m_field.add_finite_element("SKIN_FE");
86 
87  // Define rows/cols and element data
88  CHKERR m_field.modify_finite_element_add_field_row("TET_FE", "HDIV");
89  CHKERR m_field.modify_finite_element_add_field_col("TET_FE", "HDIV");
90  CHKERR m_field.modify_finite_element_add_field_data("TET_FE", "HDIV");
92  "MESH_NODE_POSITIONS");
93 
94  CHKERR m_field.modify_finite_element_add_field_row("SKIN_FE", "HDIV");
95  CHKERR m_field.modify_finite_element_add_field_col("SKIN_FE", "HDIV");
96  CHKERR m_field.modify_finite_element_add_field_data("SKIN_FE", "HDIV");
98  "MESH_NODE_POSITIONS");
99 
100  CHKERR m_field.modify_finite_element_add_field_row("TRI_FE", "HDIV");
101  CHKERR m_field.modify_finite_element_add_field_col("TRI_FE", "HDIV");
102  CHKERR m_field.modify_finite_element_add_field_data("TRI_FE", "HDIV");
104  "MESH_NODE_POSITIONS");
105 
106  // Problem
107  CHKERR m_field.add_problem("TEST_PROBLEM");
108 
109  // set finite elements for problem
110  CHKERR m_field.modify_problem_add_finite_element("TEST_PROBLEM", "TET_FE");
111  CHKERR m_field.modify_problem_add_finite_element("TEST_PROBLEM", "SKIN_FE");
112  CHKERR m_field.modify_problem_add_finite_element("TEST_PROBLEM", "TRI_FE");
113 
114  // set refinement level for problem
115  CHKERR m_field.modify_problem_ref_level_add_bit("TEST_PROBLEM", bit_level0);
116 
117  // meshset consisting all entities in mesh
118  EntityHandle root_set = moab.get_root_set();
119  // add entities to field
120  CHKERR m_field.add_ents_to_field_by_type(root_set, MBTET, "HDIV");
121 
122  // add entities to finite element
123  CHKERR m_field.add_ents_to_finite_element_by_type(root_set, MBTET,
124  "TET_FE");
125 
126  Range tets;
127  CHKERR m_field.getInterface<BitRefManager>()->getEntitiesByTypeAndRefLevel(
128  BitRefLevel().set(0), BitRefLevel().set(), MBTET, tets);
129  Skinner skin(&moab);
130  Range skin_faces; // skin faces from 3d ents
131  CHKERR skin.find_skin(0, tets, false, skin_faces);
132  CHKERR m_field.add_ents_to_finite_element_by_type(skin_faces, MBTRI,
133  "SKIN_FE");
134 
135  Range faces;
136  CHKERR m_field.getInterface<BitRefManager>()->getEntitiesByTypeAndRefLevel(
137  BitRefLevel().set(0), BitRefLevel().set(), MBTRI, faces);
138  faces = subtract(faces, skin_faces);
139  CHKERR m_field.add_ents_to_finite_element_by_type(faces, MBTRI, "TRI_FE");
140 
141  // set app. order
142  int order = 4;
143  CHKERR m_field.set_field_order(root_set, MBTET, "HDIV", order);
144  CHKERR m_field.set_field_order(root_set, MBTRI, "HDIV", order);
145 
146  CHKERR m_field.add_ents_to_field_by_type(0, MBTET, "MESH_NODE_POSITIONS");
147  CHKERR m_field.set_field_order(0, MBVERTEX, "MESH_NODE_POSITIONS", 1);
148  CHKERR m_field.set_field_order(0, MBEDGE, "MESH_NODE_POSITIONS", 2);
149  CHKERR m_field.set_field_order(0, MBTRI, "MESH_NODE_POSITIONS", 2);
150  CHKERR m_field.set_field_order(0, MBTET, "MESH_NODE_POSITIONS", 2);
151 
152  // build field
153  CHKERR m_field.build_fields();
154  // build finite elemnts
155  CHKERR m_field.build_finite_elements();
156  // build adjacencies
157  CHKERR m_field.build_adjacencies(bit_level0);
158  // build problem
159  ProblemsManager *prb_mng_ptr;
160  CHKERR m_field.getInterface(prb_mng_ptr);
161  CHKERR prb_mng_ptr->buildProblem("TEST_PROBLEM", true);
162 
163  // project geometry form 10 node tets on higher order approx. functions
164  Projection10NodeCoordsOnField ent_method(m_field, "MESH_NODE_POSITIONS");
165  CHKERR m_field.loop_dofs("MESH_NODE_POSITIONS", ent_method);
166 
167  // mesh partitioning
168  // partition
169  CHKERR prb_mng_ptr->partitionSimpleProblem("TEST_PROBLEM");
170  CHKERR prb_mng_ptr->partitionFiniteElements("TEST_PROBLEM");
171  // what are ghost nodes, see Petsc Manual
172  CHKERR prb_mng_ptr->partitionGhostDofs("TEST_PROBLEM");
173 
174  Vec v;
175  CHKERR m_field.getInterface<VecManager>()->vecCreateGhost("TEST_PROBLEM",
176  ROW, &v);
177  CHKERR VecSetRandom(v, PETSC_NULL);
178  CHKERR m_field.getInterface<VecManager>()->setLocalGhostVector(
179  "TEST_PROBLEM", ROW, v, INSERT_VALUES, SCATTER_REVERSE);
180  CHKERR VecDestroy(&v);
181 
182  typedef tee_device<std::ostream, std::ofstream> TeeDevice;
183  typedef stream<TeeDevice> TeeStream;
184  std::ofstream ofs("forces_and_sources_hdiv_continuity_check.txt");
185  TeeDevice my_tee(std::cout, ofs);
186  TeeStream my_split(my_tee);
187 
188  struct OpTetFluxes
190 
191  MoFEM::Interface &m_field;
192  Tag tH;
193 
194  OpTetFluxes(MoFEM::Interface &m_field, Tag th)
196  "HDIV", UserDataOperator::OPROW),
197  m_field(m_field), tH(th) {}
198 
199  MoFEMErrorCode doWork(int side, EntityType type,
202 
203  if (data.getFieldData().size() == 0)
205 
206  if (type == MBTRI) {
207 
208  boost::shared_ptr<const NumeredEntFiniteElement> mofem_fe =
209  getNumeredEntFiniteElementPtr();
210  SideNumber_multiIndex &side_table = mofem_fe->getSideNumberTable();
211  EntityHandle face = side_table.get<1>()
212  .find(boost::make_tuple(type, side))
213  ->get()
214  ->ent;
215  int sense = side_table.get<1>()
216  .find(boost::make_tuple(type, side))
217  ->get()
218  ->sense;
219 
220  VectorDouble t(3, 0);
221  int nb_dofs = data.getN().size2() / 3;
222  for (int dd = 0; dd != nb_dofs; dd++) {
223  for (int ddd = 0; ddd != 3; ddd++)
224  t(ddd) +=
225  data.getVectorN<3>(side)(dd, ddd) * data.getFieldData()[dd];
226  }
227 
228  double *t_ptr;
229  CHKERR m_field.get_moab().tag_get_by_ptr(tH, &face, 1,
230  (const void **)&t_ptr);
231  for (int dd = 0; dd < 3; dd++)
232  t_ptr[dd] += sense * t[dd];
233  }
234 
236  }
237  };
238 
239  struct MyTetFE : public VolumeElementForcesAndSourcesCore {
240 
241  MyTetFE(MoFEM::Interface &m_field)
243  int getRule(int order) { return -1; };
244 
245  MatrixDouble N;
246  MoFEMErrorCode setGaussPts(int order) {
248 
249  N.resize(1, 3);
250  CHKERR ShapeMBTRI(&N(0, 0), G_TRI_X1, G_TRI_Y1, 1);
251 
252  gaussPts.resize(4, 4);
253  for (int ff = 0; ff < 4; ff++) {
254  for (int dd = 0; dd < 3; dd++) {
255  gaussPts(dd, ff) =
256  cblas_ddot(3, &N(0, 0), 1, &face_coords[ff][dd], 3);
257  }
258  gaussPts(3, ff) = G_TRI_W1[0];
259  }
260 
262  }
263  };
264 
265  struct OpFacesSkinFluxes
267 
268  MoFEM::Interface &m_field;
269  Tag tH1, tH2;
270  TeeStream &mySplit;
271 
272  OpFacesSkinFluxes(MoFEM::Interface &m_field, Tag th1, Tag th2,
273  TeeStream &my_split)
275  "HDIV", UserDataOperator::OPROW),
276  m_field(m_field), tH1(th1), tH2(th2), mySplit(my_split) {}
277 
278  MoFEMErrorCode doWork(int side, EntityType type,
281 
282  if (type != MBTRI)
284 
285  EntityHandle face = getNumeredEntFiniteElementPtr()->getEnt();
286 
287  double *t_ptr;
288  CHKERR m_field.get_moab().tag_get_by_ptr(tH1, &face, 1,
289  (const void **)&t_ptr);
290  double *tn_ptr;
291  CHKERR m_field.get_moab().tag_get_by_ptr(tH2, &face, 1,
292  (const void **)&tn_ptr);
293 
294  *tn_ptr = getNormalsAtGaussPts()(0, 0) * t_ptr[0] +
295  getNormalsAtGaussPts()(0, 1) * t_ptr[1] +
296  getNormalsAtGaussPts()(0, 2) * t_ptr[2];
297 
298  int nb_dofs = data.getN().size2() / 3;
299  int dd = 0;
300  for (; dd < nb_dofs; dd++) {
301  *tn_ptr += -getNormalsAtGaussPts()(0, 0) *
302  data.getN()(0, 3 * dd + 0) * data.getFieldData()[dd] -
303  getNormalsAtGaussPts()(0, 1) * data.getN()(0, 3 * dd + 1) *
304  data.getFieldData()[dd] -
305  getNormalsAtGaussPts()(0, 2) * data.getN()(0, 3 * dd + 2) *
306  data.getFieldData()[dd];
307  }
308 
309  const double eps = 1e-8;
310  if (fabs(*tn_ptr) > eps) {
311  SETERRQ1(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
312  "HDiv continuity failed %6.4e", *tn_ptr);
313  }
314 
315  mySplit.precision(5);
316 
317  mySplit << face << " " /*<< std::fixed*/ << fabs(*tn_ptr) << std::endl;
318 
320  }
321  };
322 
323  struct OpFacesFluxes
325 
326  MoFEM::Interface &m_field;
327  Tag tH1, tH2;
328  TeeStream &mySplit;
329 
330  OpFacesFluxes(MoFEM::Interface &m_field, Tag th1, Tag th2,
331  TeeStream &my_split)
333  "HDIV", UserDataOperator::OPROW),
334  m_field(m_field), tH1(th1), tH2(th2), mySplit(my_split) {}
335 
336  MoFEMErrorCode doWork(int side, EntityType type,
339 
340  if (type != MBTRI)
342 
343  EntityHandle face = getNumeredEntFiniteElementPtr()->getEnt();
344 
345  double *t_ptr;
346  CHKERR m_field.get_moab().tag_get_by_ptr(tH1, &face, 1,
347  (const void **)&t_ptr);
348  double *tn_ptr;
349  CHKERR m_field.get_moab().tag_get_by_ptr(tH2, &face, 1,
350  (const void **)&tn_ptr);
351 
352  *tn_ptr = getNormalsAtGaussPts()(0, 0) * t_ptr[0] +
353  getNormalsAtGaussPts()(0, 1) * t_ptr[1] +
354  getNormalsAtGaussPts()(0, 2) * t_ptr[2];
355 
356  const double eps = 1e-8;
357  if (fabs(*tn_ptr) > eps) {
358  SETERRQ1(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
359  "HDiv continuity failed %6.4e", *tn_ptr);
360  }
361 
362  mySplit.precision(5);
363 
364  mySplit << face << " " /*<< std::fixed*/ << fabs(*tn_ptr) << std::endl;
365 
367  }
368  };
369 
370  struct MyTriFE : public FaceElementForcesAndSourcesCore {
371 
372  MyTriFE(MoFEM::Interface &m_field)
373  : FaceElementForcesAndSourcesCore(m_field) {}
374  int getRule(int order) { return -1; };
375 
376  MoFEMErrorCode setGaussPts(int order) {
378 
379  gaussPts.resize(3, 1);
380  gaussPts(0, 0) = G_TRI_X1[0];
381  gaussPts(1, 0) = G_TRI_Y1[0];
382  gaussPts(2, 0) = G_TRI_W1[0];
383 
385  }
386  };
387 
388  MyTetFE tet_fe(m_field);
389  MyTriFE tri_fe(m_field);
390  MyTriFE skin_fe(m_field);
391 
392  Tag th1;
393  double def_val[] = {0, 0, 0};
394  CHKERR moab.tag_get_handle("T", 3, MB_TYPE_DOUBLE, th1,
395  MB_TAG_CREAT | MB_TAG_SPARSE, &def_val);
396 
397  auto material_grad_mat = boost::make_shared<MatrixDouble>();
398  auto material_det_vec = boost::make_shared<VectorDouble>();
399  auto material_inv_grad_mat = boost::make_shared<MatrixDouble>();
400 
401  tet_fe.getOpPtrVector().push_back(new OpCalculateVectorFieldGradient<3, 3>(
402  "MESH_NODE_POSITIONS", material_grad_mat));
403  tet_fe.getOpPtrVector().push_back(new OpInvertMatrix<3>(
404  material_grad_mat, material_det_vec, material_inv_grad_mat));
405  tet_fe.getOpPtrVector().push_back(new OpSetHOWeights(material_det_vec));
406  tet_fe.getOpPtrVector().push_back(new OpSetHOContravariantPiolaTransform(
407  HDIV, material_det_vec, material_grad_mat));
408  tet_fe.getOpPtrVector().push_back(
409  new OpSetHOInvJacVectorBase(HDIV, material_inv_grad_mat));
410  tet_fe.getOpPtrVector().push_back(new OpTetFluxes(m_field, th1));
411 
412  Tag th2;
413  CHKERR moab.tag_get_handle("TN", 1, MB_TYPE_DOUBLE, th2,
414  MB_TAG_CREAT | MB_TAG_SPARSE, &def_val);
415  tri_fe.getOpPtrVector().push_back(
416  new OpGetHONormalsOnFace("MESH_NODE_POSITIONS"));
417  tri_fe.getOpPtrVector().push_back(
419  tri_fe.getOpPtrVector().push_back(
420  new OpFacesFluxes(m_field, th1, th2, my_split));
421 
422  skin_fe.getOpPtrVector().push_back(
423  new OpGetHONormalsOnFace("MESH_NODE_POSITIONS"));
424  skin_fe.getOpPtrVector().push_back(
426  skin_fe.getOpPtrVector().push_back(
427  new OpFacesSkinFluxes(m_field, th1, th2, my_split));
428 
429  for (Range::iterator fit = faces.begin(); fit != faces.end(); fit++) {
430  CHKERR moab.tag_set_data(th1, &*fit, 1, &def_val);
431  CHKERR moab.tag_set_data(th2, &*fit, 1, &def_val);
432  }
433 
434  CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "TET_FE", tet_fe);
435  my_split << "internal\n";
436  CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "TRI_FE", tri_fe);
437  my_split << "skin\n";
438  CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "SKIN_FE", skin_fe);
439 
440  EntityHandle meshset;
441  CHKERR moab.create_meshset(MESHSET_SET, meshset);
442  CHKERR m_field.getInterface<BitRefManager>()->getEntitiesByTypeAndRefLevel(
443  BitRefLevel().set(0), BitRefLevel().set(), MBTRI, meshset);
444  CHKERR moab.write_file("out.vtk", "VTK", "", &meshset, 1);
445  }
446  CATCH_ERRORS;
447 
449 }
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
help
static char help[]
Definition: hdiv_continuity_check.cpp:11
G_TRI_X1
static const double G_TRI_X1[]
Definition: fem_tools.h:312
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::CoreInterface::loop_dofs
virtual MoFEMErrorCode loop_dofs(const Problem *problem_ptr, const std::string &field_name, RowColData rc, DofMethod &method, int lower_rank, int upper_rank, int verb=DEFAULT_VERBOSITY)=0
Make a loop over dofs.
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
EntityHandle
MoFEM::ProblemsManager
Problem manager is used to build and partition problems.
Definition: ProblemsManager.hpp:21
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::OpSetHOContravariantPiolaTransform
Apply contravariant (Piola) transfer to Hdiv space for HO geometry.
Definition: HODataOperators.hpp:180
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
MoFEM::Types::MatrixDouble
UBlasMatrix< double > MatrixDouble
Definition: Types.hpp:77
MoFEM::OpHOSetContravariantPiolaTransformOnFace3D
transform Hdiv base fluxes from reference element to physical triangle
Definition: HODataOperators.hpp:298
MoFEM::th
Tag th
Definition: Projection10NodeCoordsOnField.cpp:122
MoFEM::EntitiesFieldData::EntData::getFieldData
const VectorDouble & getFieldData() const
get dofs values
Definition: EntitiesFieldData.hpp:1241
OpFacesFluxes
Definition: hdiv_divergence_operator.cpp:31
MoFEM.hpp
MoFEM::Projection10NodeCoordsOnField
Projection of edge entities with one mid-node on hierarchical basis.
Definition: Projection10NodeCoordsOnField.hpp:24
MoFEM::CoreTmp< 0 >::Finalize
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition: Core.cpp:112
MOFEM_IMPOSSIBLE_CASE
@ MOFEM_IMPOSSIBLE_CASE
Definition: definitions.h:35
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
ROW
@ ROW
Definition: definitions.h:123
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
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:535
MoFEM::FaceElementForcesAndSourcesCore::UserDataOperator
friend class UserDataOperator
Definition: FaceElementForcesAndSourcesCore.hpp:86
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::get_moab
virtual moab::Interface & get_moab()=0
MoFEM
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
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::FaceElementForcesAndSourcesCore::UserDataOperator
default operator for TRI element
Definition: FaceElementForcesAndSourcesCore.hpp:94
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
SideNumber_multiIndex
multi_index_container< boost::shared_ptr< SideNumber >, indexed_by< ordered_unique< member< SideNumber, EntityHandle, &SideNumber::ent > >, ordered_non_unique< composite_key< SideNumber, const_mem_fun< SideNumber, EntityType, &SideNumber::getEntType >, member< SideNumber, signed char, &SideNumber::side_number > > > > > SideNumber_multiIndex
SideNumber_multiIndex for SideNumber.
Definition: RefEntsMultiIndices.hpp:101
MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator
friend class UserDataOperator
Definition: VolumeElementForcesAndSourcesCore.hpp:105
G_TRI_Y1
static const double G_TRI_Y1[]
Definition: fem_tools.h:313
MoFEM::FaceElementForcesAndSourcesCore
Face finite element.
Definition: FaceElementForcesAndSourcesCore.hpp:23
MoFEM::OpGetHONormalsOnFace
Calculate normals at Gauss points of triangle element.
Definition: HODataOperators.hpp:281
t
constexpr double t
plate stiffness
Definition: plate.cpp:59
MoFEM::VolumeElementForcesAndSourcesCore
Volume finite element base.
Definition: VolumeElementForcesAndSourcesCore.hpp:26
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::VecManager
Vector manager is used to create vectors \mofem_vectors.
Definition: VecManager.hpp:23
MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator
Definition: VolumeElementForcesAndSourcesCore.hpp:108
mesh_file_name
char mesh_file_name[255]
Definition: mesh_smoothing.cpp:23
MoFEM::OpCalculateVectorFieldGradient
Get field gradients at integration pts for scalar filed rank 0, i.e. vector field.
Definition: UserDataOperators.hpp:1536
main
int main(int argc, char *argv[])
Definition: hdiv_continuity_check.cpp:21
MoFEM::EntitiesFieldData::EntData::getVectorN
const MatrixAdaptor getVectorN(const FieldApproximationBase base, const int gg)
get Hdiv of base functions at Gauss pts
Definition: EntitiesFieldData.hpp:1434
N
const int N
Definition: speed_test.cpp:3
face_coords
static const double face_coords[4][9]
Definition: hdiv_continuity_check.cpp:13
v
const double v
phase velocity of light in medium (cm/ns)
Definition: initial_diffusion.cpp:40
Range
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
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
G_TRI_W1
static const double G_TRI_W1[]
Definition: fem_tools.h:314
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
MoFEM::OpSetHOWeights
Set inverse jacobian to base functions.
Definition: HODataOperators.hpp:159
AINSWORTH_LEGENDRE_BASE
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:60
EigenMatrix::Vec
const FTensor::Tensor2< T, Dim, Dim > Vec
Definition: MatrixFunction.hpp:66
MoFEM::PetscOptionsGetString
PetscErrorCode PetscOptionsGetString(PetscOptions *, const char pre[], const char name[], char str[], size_t size, PetscBool *set)
Definition: DeprecatedPetsc.hpp:172
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::OpInvertMatrix
Definition: UserDataOperators.hpp:3254
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::Types::VectorDouble
UBlasVector< double > VectorDouble
Definition: Types.hpp:68
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::OpSetHOInvJacVectorBase
transform local reference derivatives of shape function to global derivatives if higher order geometr...
Definition: HODataOperators.hpp:91
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.
FractureMechanics::LASTOP
@ LASTOP
Definition: CrackFrontElement.hpp:23
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
MoFEM::ForcesAndSourcesCore::UserDataOperator::OPROW
@ OPROW
operator doWork function is executed on FE rows
Definition: ForcesAndSourcesCore.hpp:567
ShapeMBTRI
PetscErrorCode ShapeMBTRI(double *N, const double *X, const double *Y, const int G_DIM)
calculate shape functions on triangle
Definition: fem_tools.c:182