v0.14.0
forces_and_sources_testing_users_base.cpp
Go to the documentation of this file.
1 /** \file forces_and_sources_testing_users_base.cpp
2  * \example forces_and_sources_testing_users_base.cpp
3  *
4  * Primarily this is used for testing if the code can handle user base. It is
5  * also, an example of how to build and use user approximation base. This is a
6  * test, so we used RT base by Demkowicz recipe.
7  *
8  * Note that triple defines approximation element; element entity type,
9  * approximation space and approximation base. Entity type determines the
10  * integration method; approximation space determines the adjacency of the
11  * matrix and approximation base determines together with space the regularity
12  * of approximation.
13  *
14  */
15 
16 
17 
18 #include <MoFEM.hpp>
19 #include <Hdiv.hpp>
20 
21 namespace bio = boost::iostreams;
22 using bio::stream;
23 using bio::tee_device;
24 
25 using namespace MoFEM;
26 
27 // used for profiling
28 static PetscBool quiet = PETSC_FALSE;
29 static PetscBool base_cache = PETSC_FALSE;
30 
31 static char help[] = "...\n\n";
32 
33 /**
34  * @brief Class used to calculate base functions at integration points
35  *
36  */
38 
39  SomeUserPolynomialBase() = default;
40  ~SomeUserPolynomialBase() = default;
41 
42  /**
43  * @brief Return interface to this class when one ask for for tetrahedron,
44  * otherisw return interface class for generic class.
45  *
46  * @param iface interface class
47  * @return MoFEMErrorCode
48  */
49  MoFEMErrorCode query_interface(boost::typeindex::type_index type_index,
50  UnknownInterface **iface) const {
52  *iface = const_cast<SomeUserPolynomialBase *>(this);
54  }
55 
56  /**
57  * @brief Calculate base functions at intergeneration points
58  *
59  * @param pts
60  * @param ctx_ptr
61  * @return MoFEMErrorCode
62  */
64  boost::shared_ptr<BaseFunctionCtx> ctx_ptr) {
66 
67  cTx = ctx_ptr->getInterface<EntPolynomialBaseCtx>();
68 
69  int nb_gauss_pts = pts.size2();
70  if (!nb_gauss_pts) {
72  }
73 
74  if (pts.size1() < 3) {
75  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
76  "Wrong dimension of pts, should be at least 3 rows with "
77  "coordinates");
78  }
79 
80  switch (cTx->sPace) {
81  case HDIV:
82  CHKERR getValueHdivForCGGBubble(pts);
83  break;
84  default:
85  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "Not yet implemented");
86  }
87 
89  }
90 
91 private:
93 
95 
98 
99  const FieldApproximationBase base = cTx->bAse;
100  // This should be used only in case USER_BASE is selected
101  if (cTx->bAse != USER_BASE) {
102  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
103  "Wrong base, should be USER_BASE");
104  }
105 
106  // This is example, simply use Demkowicz HDiv base to generate base
107  // functions
108 
109  EntitiesFieldData &data = cTx->dAta;
110  int nb_gauss_pts = pts.size2();
111 
112  // calculate shape functions, i.e. barycentric coordinates
113  shapeFun.resize(nb_gauss_pts, 4, false);
114  CHKERR ShapeMBTET(&*shapeFun.data().begin(), &pts(0, 0), &pts(1, 0),
115  &pts(2, 0), nb_gauss_pts);
116  // derivatives of shape functions
117  double diff_shape_fun[12];
118  CHKERR ShapeDiffMBTET(diff_shape_fun);
119 
120  int volume_order = data.dataOnEntities[MBTET][0].getOrder();
121 
122  int p_f[4];
123  double *phi_f[4];
124  double *diff_phi_f[4];
125 
126  // Calculate base function on tet faces
127  for (int ff = 0; ff != 4; ff++) {
128  int face_order = data.dataOnEntities[MBTRI][ff].getOrder();
129  int order = volume_order > face_order ? volume_order : face_order;
130  data.dataOnEntities[MBTRI][ff].getN(base).resize(
131  nb_gauss_pts, 3 * NBFACETRI_DEMKOWICZ_HDIV(order), false);
132  data.dataOnEntities[MBTRI][ff].getDiffN(base).resize(
133  nb_gauss_pts, 9 * NBFACETRI_DEMKOWICZ_HDIV(order), false);
134  p_f[ff] = order;
135  phi_f[ff] = &*data.dataOnEntities[MBTRI][ff].getN(base).data().begin();
136  diff_phi_f[ff] =
137  &*data.dataOnEntities[MBTRI][ff].getDiffN(base).data().begin();
139  continue;
141  &data.facesNodes(ff, 0), order, &*shapeFun.data().begin(),
142  diff_shape_fun, phi_f[ff], diff_phi_f[ff], nb_gauss_pts, 4);
143  }
144 
145  // Calculate base functions in tet interior
146  if (NBVOLUMETET_DEMKOWICZ_HDIV(volume_order) > 0) {
147  data.dataOnEntities[MBTET][0].getN(base).resize(
148  nb_gauss_pts, 3 * NBVOLUMETET_DEMKOWICZ_HDIV(volume_order), false);
149  data.dataOnEntities[MBTET][0].getDiffN(base).resize(
150  nb_gauss_pts, 9 * NBVOLUMETET_DEMKOWICZ_HDIV(volume_order), false);
151  double *phi_v = &*data.dataOnEntities[MBTET][0].getN(base).data().begin();
152  double *diff_phi_v =
153  &*data.dataOnEntities[MBTET][0].getDiffN(base).data().begin();
155  volume_order, &*shapeFun.data().begin(), diff_shape_fun, p_f, phi_f,
156  diff_phi_f, phi_v, diff_phi_v, nb_gauss_pts);
157  }
158 
159  // Set size of face base correctly
160  for (int ff = 0; ff != 4; ff++) {
161  int face_order = data.dataOnEntities[MBTRI][ff].getOrder();
162  data.dataOnEntities[MBTRI][ff].getN(base).resize(
163  nb_gauss_pts, 3 * NBFACETRI_DEMKOWICZ_HDIV(face_order), true);
164  data.dataOnEntities[MBTRI][ff].getDiffN(base).resize(
165  nb_gauss_pts, 9 * NBFACETRI_DEMKOWICZ_HDIV(face_order), true);
166  }
167 
169  }
170 };
171 
172 int main(int argc, char *argv[]) {
173 
174  // Initialise MoFEM, MPI and petsc
175  MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
176 
177  try {
178 
179  // create moab
180  moab::Core mb_instance;
181  // get interface to moab databse
182  moab::Interface &moab = mb_instance;
183 
184  // get file
185  PetscBool flg = PETSC_TRUE;
186  char mesh_file_name[255];
187 #if PETSC_VERSION_GE(3, 6, 4)
188  CHKERR PetscOptionsGetString(PETSC_NULL, "", "-my_file", mesh_file_name,
189  255, &flg);
190 #else
191  CHKERR PetscOptionsGetString(PETSC_NULL, PETSC_NULL, "-my_file",
192  mesh_file_name, 255, &flg);
193 #endif
194  if (flg != PETSC_TRUE) {
195  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
196  "*** ERROR -my_file (MESH FILE NEEDED)");
197  }
198 
199  // create MoFEM database
200  MoFEM::Core core(moab);
201  // get interface to moab database
202  MoFEM::Interface &m_field = core;
203 
204  // load mesh file
205  const char *option;
206  option = "";
207  CHKERR moab.load_file(mesh_file_name, 0, option);
208 
209  // set bit refinement level
210  CHKERR m_field.getInterface<BitRefManager>()->setBitRefLevelByDim(
211  0, 3, BitRefLevel().set(0));
212 
213  // Create fields, field "FIELD_CGG" has user base, it means that recipe how
214  // to construct approximation is provided by user. Is set that user provided
215  // base is in h-div space.
216  CHKERR m_field.add_field("FILED_CGG", HDIV, USER_BASE, 1);
217  CHKERR m_field.add_field("FILED_RT", HDIV, DEMKOWICZ_JACOBI_BASE, 1);
218 
219  // get access to "FIELD_CGG" data structure
220  auto field_ptr = m_field.get_field_structure("FILED_CGG");
221  // get table associating number of dofs to entities depending on
222  // approximation order set on those entities.
223  auto field_order_table =
224  const_cast<Field *>(field_ptr)->getFieldOrderTable();
225 
226  // function set zero number of dofs
227  auto get_cgg_bubble_order_zero = [](int p) { return 0; };
228  // function set non-zero number of dofs on tetrahedrons
229  auto get_cgg_bubble_order_face = [](int p) {
230  return NBFACETRI_DEMKOWICZ_HDIV(p);
231  };
232  auto get_cgg_bubble_order_tet = [](int p) {
233  return NBVOLUMETET_DEMKOWICZ_HDIV(p);
234  };
235  field_order_table[MBVERTEX] = get_cgg_bubble_order_zero;
236  field_order_table[MBEDGE] = get_cgg_bubble_order_zero;
237  field_order_table[MBTRI] = get_cgg_bubble_order_face;
238  field_order_table[MBTET] = get_cgg_bubble_order_tet;
239  const_cast<Field *>(field_ptr)->rebuildDofsOrderMap();
240 
241  auto &dof_order_map = field_ptr->getDofOrderMap(MBTET);
242  for(auto d = 0; d!=10; ++d) {
243  MOFEM_LOG("WORLD", Sev::noisy) << "dof " << dof_order_map[d];
244  }
245 
246  CHKERR m_field.add_finite_element("FE");
247 
248  // define rows/cols and element data
249  CHKERR m_field.modify_finite_element_add_field_row("FE", "FILED_CGG");
250  CHKERR m_field.modify_finite_element_add_field_col("FE", "FILED_CGG");
251  CHKERR m_field.modify_finite_element_add_field_data("FE", "FILED_CGG");
252  CHKERR m_field.modify_finite_element_add_field_row("FE", "FILED_RT");
253  CHKERR m_field.modify_finite_element_add_field_col("FE", "FILED_RT");
254  CHKERR m_field.modify_finite_element_add_field_data("FE", "FILED_RT");
255 
256  // add problem
257  CHKERR m_field.add_problem("PROBLEM");
258 
259  // set finite elements for problem
260  CHKERR m_field.modify_problem_add_finite_element("PROBLEM", "FE");
261  // set refinement level for problem
262  CHKERR m_field.modify_problem_ref_level_add_bit("PROBLEM",
263  BitRefLevel().set(0));
264 
265  // meshset consisting all entities in mesh
266  EntityHandle root_set = moab.get_root_set();
267  // add entities to field
268  CHKERR m_field.add_ents_to_field_by_type(root_set, MBTET, "FILED_CGG");
269  CHKERR m_field.add_ents_to_field_by_type(root_set, MBTET, "FILED_RT");
270  // add entities to finite element
271  CHKERR m_field.add_ents_to_finite_element_by_type(root_set, MBTET, "FE");
272 
273  // set app. order
274  int order = 3;
275  CHKERR PetscOptionsGetInt(PETSC_NULL, "", "-order", &order, PETSC_NULL);
276  CHKERR m_field.set_field_order(root_set, MBTRI, "FILED_CGG", order);
277  CHKERR m_field.set_field_order(root_set, MBTET, "FILED_CGG", order);
278  CHKERR m_field.set_field_order(root_set, MBTRI, "FILED_RT", order);
279  CHKERR m_field.set_field_order(root_set, MBTET, "FILED_RT", order);
280 
281  /****/
282  // build database
283  // build field
284  CHKERR m_field.build_fields();
285  // build finite elemnts
286  CHKERR m_field.build_finite_elements();
287  // build adjacencies
288  CHKERR m_field.build_adjacencies(BitRefLevel().set(0));
289 
290  // build problem
291  CHKERR m_field.getInterface<ProblemsManager>()->buildProblem("PROBLEM",
292  true);
293  // dofs partitioning
294  CHKERR m_field.getInterface<ProblemsManager>()->partitionSimpleProblem(
295  "PROBLEM");
296  CHKERR m_field.getInterface<ProblemsManager>()->partitionFiniteElements(
297  "PROBLEM");
298  // what are ghost nodes, see Petsc Manual
299  CHKERR m_field.getInterface<ProblemsManager>()->partitionGhostDofs(
300  "PROBLEM");
301 
302  typedef tee_device<std::ostream, std::ofstream> TeeDevice;
303  typedef stream<TeeDevice> TeeStream;
304 
305  std::ofstream ofs("forces_and_sources_testing_users_base.txt");
306  TeeDevice my_tee(std::cout, ofs);
307  TeeStream my_split(my_tee);
308 
309  /**
310  * Simple user data operator which main purpose is to print values
311  * of base functions at integration points.
312  *
313  */
315 
316  TeeStream &my_split;
317  MyOp1(const std::string &row_field, const std::string &col_field,
318  TeeStream &_my_split, char type)
320  row_field, col_field, type),
321  my_split(_my_split) {
322  sYmm = false;
323  }
324 
325  MoFEMErrorCode doWork(int side, EntityType type,
328  if (data.getIndices().empty()) {
330  }
331  if (!quiet) {
332  my_split << rowFieldName << endl;
333  my_split << "side: " << side << " type: " << type << std::endl;
334  my_split << data << endl;
335  my_split << data.getN() << endl;
336  my_split << endl;
337  }
339  }
340 
341  MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type,
342  EntityType col_type,
343  EntitiesFieldData::EntData &row_data,
344  EntitiesFieldData::EntData &col_data) {
346  if (row_data.getIndices().empty())
348  if (col_data.getIndices().empty())
350  if (!quiet) {
351  my_split << rowFieldName << " : " << colFieldName << endl;
352  my_split << "row side: " << row_side << " row_type: " << row_type
353  << std::endl;
354  my_split << "col side: " << col_side << " col_type: " << col_type
355  << std::endl;
356  my_split << row_data.getIndices().size() << " : "
357  << col_data.getIndices().size() << endl;
358  my_split << endl;
359  }
361  }
362  };
363 
364  // create finite element instance
366  // set class needed to construct user approximation base
367  fe1.getUserPolynomialBase() =
368  boost::shared_ptr<BaseFunction>(new SomeUserPolynomialBase());
369 
370  CHKERR PetscOptionsGetBool(PETSC_NULL, "", "-quiet", &quiet, PETSC_NULL);
371 
372  // push user data operators
373  fe1.getOpPtrVector().push_back(
374  new MyOp1("FILED_CGG", "FILED_CGG", my_split,
376  fe1.getOpPtrVector().push_back(
377  new MyOp1("FILED_CGG", "FILED_RT", my_split,
379 
380  CHKERR PetscOptionsGetBool(PETSC_NULL, "", "-base_cache", &base_cache,
381  PETSC_NULL);
382 
383  if (base_cache) {
386  }
389  }
390  }
391 
392  // iterate over finite elements, and execute user data operators on each
393  // of them
394  CHKERR m_field.loop_finite_elements("PROBLEM", "FE", fe1);
395 
396  if(base_cache) {
399  }
400 
401  }
402  CATCH_ERRORS;
403 
405 
406  return 0;
407 }
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:460
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference 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:128
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::CoreTmp< 0 >
Core (interface) class.
Definition: Core.hpp:82
MoFEM::Hdiv_Demkowicz_Face_MBTET_ON_FACE
MoFEMErrorCode Hdiv_Demkowicz_Face_MBTET_ON_FACE(int *faces_nodes, int p, double *N, double *diffN, double *phi_f, double *diff_phi_f, int gdim, int nb)
Definition: Hdiv.cpp:617
EntityHandle
MoFEM::Field::getDofOrderMap
const std::array< ApproximationOrder, MAX_DOFS_ON_ENTITY > & getDofOrderMap(const EntityType type) const
get hash-map relating dof index on entity with its order
Definition: FieldMultiIndices.hpp:268
MoFEM::ProblemsManager
Problem manager is used to build and partition problems.
Definition: ProblemsManager.hpp:21
quiet
static PetscBool quiet
Definition: forces_and_sources_testing_users_base.cpp:28
MoFEM::EntPolynomialBaseCtx
Class used to pass element data to calculate base functions on tet,triangle,edge.
Definition: EntPolynomialBaseCtx.hpp:22
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::CoreInterface::get_field_structure
virtual const Field * get_field_structure(const std::string &name, enum MoFEMTypes bh=MF_EXIST) const =0
get field structure
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
MoFEM::Types::MatrixDouble
UBlasMatrix< double > MatrixDouble
Definition: Types.hpp:77
MoFEM::EntPolynomialBaseCtx::dAta
EntitiesFieldData & dAta
Definition: EntPolynomialBaseCtx.hpp:36
MoFEM.hpp
help
static char help[]
Definition: forces_and_sources_testing_users_base.cpp:31
MoFEM::Hdiv_Demkowicz_Interior_MBTET
MoFEMErrorCode Hdiv_Demkowicz_Interior_MBTET(int p, double *N, double *diffN, int p_face[], double *phi_f[4], double *diff_phi_f[4], double *phi_v, double *diff_phi_v, int gdim)
Definition: Hdiv.cpp:762
MoFEM::BaseFunction
Base class if inherited used to calculate base functions.
Definition: BaseFunction.hpp:40
MoFEM::CoreTmp< 0 >::Finalize
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition: Core.cpp:112
MoFEM::ForcesAndSourcesCore::UserDataOperator::OPROWCOL
@ OPROWCOL
operator doWork is executed on FE rows &columns
Definition: ForcesAndSourcesCore.hpp:569
MoFEM::EntitiesFieldData::facesNodes
MatrixInt facesNodes
nodes on finite element faces
Definition: EntitiesFieldData.hpp:45
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.
MoFEM::Field
Provide data structure for (tensor) field approximation.
Definition: FieldMultiIndices.hpp:51
USER_BASE
@ USER_BASE
user implemented approximation base
Definition: definitions.h:68
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:2002
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:548
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::ForcesAndSourcesCore::getUserPolynomialBase
auto & getUserPolynomialBase()
Get the User Polynomial Base object.
Definition: ForcesAndSourcesCore.hpp:97
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::TetPolynomialBase::swichCacheHdivBaseInteriorDemkowicz
static bool swichCacheHdivBaseInteriorDemkowicz(const void *ptr)
Definition: TetPolynomialBase.cpp:2007
SomeUserPolynomialBase::getValue
MoFEMErrorCode getValue(MatrixDouble &pts, boost::shared_ptr< BaseFunctionCtx > ctx_ptr)
Calculate base functions at intergeneration points.
Definition: forces_and_sources_testing_users_base.cpp:63
MoFEM::CoreInterface::build_finite_elements
virtual MoFEMErrorCode build_finite_elements(int verb=DEFAULT_VERBOSITY)=0
Build finite elements.
NBFACETRI_DEMKOWICZ_HDIV
#define NBFACETRI_DEMKOWICZ_HDIV(P)
Definition: h1_hdiv_hcurl_l2.h:139
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.
base_cache
static PetscBool base_cache
Definition: forces_and_sources_testing_users_base.cpp:29
convert.type
type
Definition: convert.py:64
MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator
friend class UserDataOperator
Definition: VolumeElementForcesAndSourcesCore.hpp:105
MoFEM::EntitiesFieldData::EntData::getIndices
const VectorInt & getIndices() const
Get global indices of dofs on entity.
Definition: EntitiesFieldData.hpp:1204
SomeUserPolynomialBase::shapeFun
MatrixDouble shapeFun
Definition: forces_and_sources_testing_users_base.cpp:94
NBVOLUMETET_DEMKOWICZ_HDIV
#define NBVOLUMETET_DEMKOWICZ_HDIV(P)
Definition: h1_hdiv_hcurl_l2.h:140
MoFEM::VolumeElementForcesAndSourcesCore
Volume finite element base.
Definition: VolumeElementForcesAndSourcesCore.hpp:26
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::EntPolynomialBaseCtx::bAse
const FieldApproximationBase bAse
Definition: EntPolynomialBaseCtx.hpp:39
MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator
Definition: VolumeElementForcesAndSourcesCore.hpp:108
mesh_file_name
char mesh_file_name[255]
Definition: mesh_smoothing.cpp:23
MoFEM::UnknownInterface
base class for all interface classes
Definition: UnknownInterface.hpp:34
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
MOFEM_LOG
#define MOFEM_LOG(channel, severity)
Log.
Definition: LogManager.hpp:308
CATCH_ERRORS
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:385
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:1308
SomeUserPolynomialBase::query_interface
MoFEMErrorCode query_interface(boost::typeindex::type_index type_index, UnknownInterface **iface) const
Return interface to this class when one ask for for tetrahedron, otherisw return interface class for ...
Definition: forces_and_sources_testing_users_base.cpp:49
MoFEM::Core
CoreTmp< 0 > Core
Definition: Core.hpp:1140
Hdiv.hpp
Implementation of H-curl base function.
ShapeDiffMBTET
PetscErrorCode ShapeDiffMBTET(double *diffN)
calculate derivatives of shape functions
Definition: fem_tools.c:319
MoFEM::ForcesAndSourcesCore::getOpPtrVector
boost::ptr_deque< UserDataOperator > & getOpPtrVector()
Use to push back operator for row operator.
Definition: ForcesAndSourcesCore.hpp:83
MOFEM_DATA_INCONSISTENCY
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
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
sdf_hertz_2d_axisymm_plane.d
float d
Definition: sdf_hertz_2d_axisymm_plane.py:4
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::EntitiesFieldData::dataOnEntities
std::array< boost::ptr_vector< EntData >, MBMAXTYPE > dataOnEntities
Definition: EntitiesFieldData.hpp:57
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:453
SomeUserPolynomialBase::getValueHdivForCGGBubble
MoFEMErrorCode getValueHdivForCGGBubble(MatrixDouble &pts)
Definition: forces_and_sources_testing_users_base.cpp:96
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
SomeUserPolynomialBase
Class used to calculate base functions at integration points.
Definition: forces_and_sources_testing_users_base.cpp:37
main
int main(int argc, char *argv[])
Definition: forces_and_sources_testing_users_base.cpp:172
MoFEM::CoreInterface::add_problem
virtual MoFEMErrorCode add_problem(const std::string &name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
Add problem.
MoFEM::PetscOptionsGetInt
PetscErrorCode PetscOptionsGetInt(PetscOptions *, const char pre[], const char name[], PetscInt *ivalue, PetscBool *set)
Definition: DeprecatedPetsc.hpp:142
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
SomeUserPolynomialBase::cTx
EntPolynomialBaseCtx * cTx
Definition: forces_and_sources_testing_users_base.cpp:92
HDIV
@ HDIV
field with continuous normal traction
Definition: definitions.h:87
MOFEM_NOT_IMPLEMENTED
@ MOFEM_NOT_IMPLEMENTED
Definition: definitions.h:32
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
MoFEM::TetPolynomialBase::swichCacheHDivBaseFaceDemkowicz
static bool swichCacheHDivBaseFaceDemkowicz(const void *ptr)
Definition: TetPolynomialBase.cpp:2002
MoFEM::ForcesAndSourcesCore::UserDataOperator::OPROW
@ OPROW
operator doWork function is executed on FE rows
Definition: ForcesAndSourcesCore.hpp:567
ShapeMBTET
PetscErrorCode ShapeMBTET(double *N, const double *G_X, const double *G_Y, const double *G_Z, int DIM)
calculate shape functions
Definition: fem_tools.c:306
MoFEM::PetscOptionsGetBool
PetscErrorCode PetscOptionsGetBool(PetscOptions *, const char pre[], const char name[], PetscBool *bval, PetscBool *set)
Definition: DeprecatedPetsc.hpp:182