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 static char help[] = "...\n\n";
28 
29 /**
30  * @brief Class used to calculate base functions at integration points
31  *
32  */
34 
35  SomeUserPolynomialBase() = default;
36  ~SomeUserPolynomialBase() = default;
37 
38  /**
39  * @brief Return interface to this class when one ask for for tetrahedron,
40  * otherisw return interface class for generic class.
41  *
42  * @param iface interface class
43  * @return MoFEMErrorCode
44  */
45  MoFEMErrorCode query_interface(boost::typeindex::type_index type_index,
46  UnknownInterface **iface) const {
48  *iface = const_cast<SomeUserPolynomialBase *>(this);
50  }
51 
52  /**
53  * @brief Calculate base functions at intergeneration points
54  *
55  * @param pts
56  * @param ctx_ptr
57  * @return MoFEMErrorCode
58  */
60  boost::shared_ptr<BaseFunctionCtx> ctx_ptr) {
62 
63  cTx = ctx_ptr->getInterface<EntPolynomialBaseCtx>();
64 
65  int nb_gauss_pts = pts.size2();
66  if (!nb_gauss_pts) {
68  }
69 
70  if (pts.size1() < 3) {
71  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
72  "Wrong dimension of pts, should be at least 3 rows with "
73  "coordinates");
74  }
75 
76  switch (cTx->sPace) {
77  case HDIV:
78  CHKERR getValueHdivForCGGBubble(pts);
79  break;
80  default:
81  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "Not yet implemented");
82  }
83 
85  }
86 
87 private:
89 
91 
94 
95  const FieldApproximationBase base = cTx->bAse;
96  // This should be used only in case USER_BASE is selected
97  if (cTx->bAse != USER_BASE) {
98  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
99  "Wrong base, should be USER_BASE");
100  }
101 
102  // This is example, simply use Demkowicz HDiv base to generate base
103  // functions
104 
105  EntitiesFieldData &data = cTx->dAta;
106  int nb_gauss_pts = pts.size2();
107 
108  // calculate shape functions, i.e. barycentric coordinates
109  shapeFun.resize(nb_gauss_pts, 4, false);
110  CHKERR ShapeMBTET(&*shapeFun.data().begin(), &pts(0, 0), &pts(1, 0),
111  &pts(2, 0), nb_gauss_pts);
112  // derivatives of shape functions
113  double diff_shape_fun[12];
114  CHKERR ShapeDiffMBTET(diff_shape_fun);
115 
116  int volume_order = data.dataOnEntities[MBTET][0].getOrder();
117 
118  int p_f[4];
119  double *phi_f[4];
120  double *diff_phi_f[4];
121 
122  // Calculate base function on tet faces
123  for (int ff = 0; ff != 4; ff++) {
124  int face_order = data.dataOnEntities[MBTRI][ff].getOrder();
125  int order = volume_order > face_order ? volume_order : face_order;
126  data.dataOnEntities[MBTRI][ff].getN(base).resize(
127  nb_gauss_pts, 3 * NBFACETRI_DEMKOWICZ_HDIV(order), false);
128  data.dataOnEntities[MBTRI][ff].getDiffN(base).resize(
129  nb_gauss_pts, 9 * NBFACETRI_DEMKOWICZ_HDIV(order), false);
130  p_f[ff] = order;
131  phi_f[ff] = &*data.dataOnEntities[MBTRI][ff].getN(base).data().begin();
132  diff_phi_f[ff] =
133  &*data.dataOnEntities[MBTRI][ff].getDiffN(base).data().begin();
135  continue;
137  &data.facesNodes(ff, 0), order, &*shapeFun.data().begin(),
138  diff_shape_fun, phi_f[ff], diff_phi_f[ff], nb_gauss_pts, 4);
139  }
140 
141  // Calculate base functions in tet interior
142  if (NBVOLUMETET_DEMKOWICZ_HDIV(volume_order) > 0) {
143  data.dataOnEntities[MBTET][0].getN(base).resize(
144  nb_gauss_pts, 3 * NBVOLUMETET_DEMKOWICZ_HDIV(volume_order), false);
145  data.dataOnEntities[MBTET][0].getDiffN(base).resize(
146  nb_gauss_pts, 9 * NBVOLUMETET_DEMKOWICZ_HDIV(volume_order), false);
147  double *phi_v = &*data.dataOnEntities[MBTET][0].getN(base).data().begin();
148  double *diff_phi_v =
149  &*data.dataOnEntities[MBTET][0].getDiffN(base).data().begin();
151  volume_order, &*shapeFun.data().begin(), diff_shape_fun, p_f, phi_f,
152  diff_phi_f, phi_v, diff_phi_v, nb_gauss_pts);
153  }
154 
155  // Set size of face base correctly
156  for (int ff = 0; ff != 4; ff++) {
157  int face_order = data.dataOnEntities[MBTRI][ff].getOrder();
158  data.dataOnEntities[MBTRI][ff].getN(base).resize(
159  nb_gauss_pts, 3 * NBFACETRI_DEMKOWICZ_HDIV(face_order), true);
160  data.dataOnEntities[MBTRI][ff].getDiffN(base).resize(
161  nb_gauss_pts, 9 * NBFACETRI_DEMKOWICZ_HDIV(face_order), true);
162  }
163 
165  }
166 };
167 
168 int main(int argc, char *argv[]) {
169 
170  // Initialise MoFEM, MPI and petsc
171  MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
172 
173  try {
174 
175  // create moab
176  moab::Core mb_instance;
177  // get interface to moab databse
178  moab::Interface &moab = mb_instance;
179 
180  // get file
181  PetscBool flg = PETSC_TRUE;
182  char mesh_file_name[255];
183 #if PETSC_VERSION_GE(3, 6, 4)
184  CHKERR PetscOptionsGetString(PETSC_NULL, "", "-my_file", mesh_file_name,
185  255, &flg);
186 #else
187  CHKERR PetscOptionsGetString(PETSC_NULL, PETSC_NULL, "-my_file",
188  mesh_file_name, 255, &flg);
189 #endif
190  if (flg != PETSC_TRUE) {
191  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
192  "*** ERROR -my_file (MESH FILE NEEDED)");
193  }
194 
195  // create MoFEM database
196  MoFEM::Core core(moab);
197  // get interface to moab database
198  MoFEM::Interface &m_field = core;
199 
200  // load mesh file
201  const char *option;
202  option = "";
203  CHKERR moab.load_file(mesh_file_name, 0, option);
204 
205  // set bit refinement level
206  CHKERR m_field.getInterface<BitRefManager>()->setBitRefLevelByDim(
207  0, 3, BitRefLevel().set(0));
208 
209  // Create fields, field "FIELD_CGG" has user base, it means that recipe how
210  // to construct approximation is provided by user. Is set that user provided
211  // base is in h-div space.
212  CHKERR m_field.add_field("FILED_CGG", HDIV, USER_BASE, 1);
213  CHKERR m_field.add_field("FILED_RT", HDIV, DEMKOWICZ_JACOBI_BASE, 1);
214 
215  // get access to "FIELD_CGG" data structure
216  auto field_ptr = m_field.get_field_structure("FILED_CGG");
217  // get table associating number of dofs to entities depending on
218  // approximation order set on those entities.
219  auto field_order_table =
220  const_cast<Field *>(field_ptr)->getFieldOrderTable();
221 
222  // function set zero number of dofs
223  auto get_cgg_bubble_order_zero = [](int p) { return 0; };
224  // function set non-zero number of dofs on tetrahedrons
225  auto get_cgg_bubble_order_face = [](int p) {
226  return NBFACETRI_DEMKOWICZ_HDIV(p);
227  };
228  auto get_cgg_bubble_order_tet = [](int p) {
229  return NBVOLUMETET_DEMKOWICZ_HDIV(p);
230  };
231  field_order_table[MBVERTEX] = get_cgg_bubble_order_zero;
232  field_order_table[MBEDGE] = get_cgg_bubble_order_zero;
233  field_order_table[MBTRI] = get_cgg_bubble_order_face;
234  field_order_table[MBTET] = get_cgg_bubble_order_tet;
235  const_cast<Field *>(field_ptr)->rebuildDofsOrderMap();
236 
237  auto &dof_order_map = field_ptr->getDofOrderMap(MBTET);
238  for(auto d = 0; d!=10; ++d) {
239  MOFEM_LOG("WORLD", Sev::noisy) << "dof " << dof_order_map[d];
240  }
241 
242  CHKERR m_field.add_finite_element("FE");
243 
244  // define rows/cols and element data
245  CHKERR m_field.modify_finite_element_add_field_row("FE", "FILED_CGG");
246  CHKERR m_field.modify_finite_element_add_field_col("FE", "FILED_CGG");
247  CHKERR m_field.modify_finite_element_add_field_data("FE", "FILED_CGG");
248  CHKERR m_field.modify_finite_element_add_field_row("FE", "FILED_RT");
249  CHKERR m_field.modify_finite_element_add_field_col("FE", "FILED_RT");
250  CHKERR m_field.modify_finite_element_add_field_data("FE", "FILED_RT");
251 
252  // add problem
253  CHKERR m_field.add_problem("PROBLEM");
254 
255  // set finite elements for problem
256  CHKERR m_field.modify_problem_add_finite_element("PROBLEM", "FE");
257  // set refinement level for problem
258  CHKERR m_field.modify_problem_ref_level_add_bit("PROBLEM",
259  BitRefLevel().set(0));
260 
261  // meshset consisting all entities in mesh
262  EntityHandle root_set = moab.get_root_set();
263  // add entities to field
264  CHKERR m_field.add_ents_to_field_by_type(root_set, MBTET, "FILED_CGG");
265  CHKERR m_field.add_ents_to_field_by_type(root_set, MBTET, "FILED_RT");
266  // add entities to finite element
267  CHKERR m_field.add_ents_to_finite_element_by_type(root_set, MBTET, "FE");
268 
269  // set app. order
270  int order = 3;
271  CHKERR m_field.set_field_order(root_set, MBTRI, "FILED_CGG", order);
272  CHKERR m_field.set_field_order(root_set, MBTET, "FILED_CGG", order);
273  CHKERR m_field.set_field_order(root_set, MBTRI, "FILED_RT", order);
274  CHKERR m_field.set_field_order(root_set, MBTET, "FILED_RT", order);
275 
276  /****/
277  // build database
278  // build field
279  CHKERR m_field.build_fields();
280  // build finite elemnts
281  CHKERR m_field.build_finite_elements();
282  // build adjacencies
283  CHKERR m_field.build_adjacencies(BitRefLevel().set(0));
284 
285  // build problem
286  CHKERR m_field.getInterface<ProblemsManager>()->buildProblem("PROBLEM",
287  true);
288  // dofs partitioning
289  CHKERR m_field.getInterface<ProblemsManager>()->partitionSimpleProblem(
290  "PROBLEM");
291  CHKERR m_field.getInterface<ProblemsManager>()->partitionFiniteElements(
292  "PROBLEM");
293  // what are ghost nodes, see Petsc Manual
294  CHKERR m_field.getInterface<ProblemsManager>()->partitionGhostDofs(
295  "PROBLEM");
296 
297  typedef tee_device<std::ostream, std::ofstream> TeeDevice;
298  typedef stream<TeeDevice> TeeStream;
299 
300  std::ofstream ofs("forces_and_sources_testing_users_base.txt");
301  TeeDevice my_tee(std::cout, ofs);
302  TeeStream my_split(my_tee);
303 
304  /**
305  * Simple user data operator which main purpose is to print values
306  * of base functions at intergation points.
307  *
308  */
310 
311  TeeStream &my_split;
312  MyOp1(const std::string &row_filed, const std::string &col_field,
313  TeeStream &_my_split, char type)
315  row_filed, col_field, type),
316  my_split(_my_split) {
317  sYmm = false;
318  }
319 
320  MoFEMErrorCode doWork(int side, EntityType type,
323  if (data.getIndices().empty()) {
325  }
326  my_split << rowFieldName << endl;
327  my_split << "side: " << side << " type: " << type << std::endl;
328  my_split << data << endl;
329  my_split << data.getN() << endl;
330  my_split << endl;
332  }
333 
334  MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type,
335  EntityType col_type,
336  EntitiesFieldData::EntData &row_data,
337  EntitiesFieldData::EntData &col_data) {
339  if (row_data.getIndices().empty())
341  if (col_data.getIndices().empty())
343  my_split << rowFieldName << " : " << colFieldName << endl;
344  my_split << "row side: " << row_side << " row_type: " << row_type
345  << std::endl;
346  my_split << "col side: " << col_side << " col_type: " << col_type
347  << std::endl;
348  my_split << row_data.getIndices().size() << " : "
349  << col_data.getIndices().size() << endl;
350  my_split << endl;
352  }
353  };
354 
355  // create finite element instance
357  // set class needed to construct user approximation base
358  fe1.getUserPolynomialBase() =
359  boost::shared_ptr<BaseFunction>(new SomeUserPolynomialBase());
360 
361  // push user data oprators
362  fe1.getOpPtrVector().push_back(
363  new MyOp1("FILED_CGG", "FILED_CGG", my_split,
365  fe1.getOpPtrVector().push_back(
366  new MyOp1("FILED_CGG", "FILED_RT", my_split,
368 
369  // iterate over finite elements, and execute user data operators on each
370  // of them
371  CHKERR m_field.loop_finite_elements("PROBLEM", "FE", fe1);
372  }
373  CATCH_ERRORS;
374 
376 
377  return 0;
378 }
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::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
sdf_hertz.d
float d
Definition: sdf_hertz.py:5
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:254
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
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:27
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:46
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: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::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
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:59
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.
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:1201
SomeUserPolynomialBase::shapeFun
MatrixDouble shapeFun
Definition: forces_and_sources_testing_users_base.cpp:90
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:38
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: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
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:45
MoFEM::Core
CoreTmp< 0 > Core
Definition: Core.hpp:1094
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
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: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
SomeUserPolynomialBase::getValueHdivForCGGBubble
MoFEMErrorCode getValueHdivForCGGBubble(MatrixDouble &pts)
Definition: forces_and_sources_testing_users_base.cpp:92
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:33
main
int main(int argc, char *argv[])
Definition: forces_and_sources_testing_users_base.cpp:168
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
SomeUserPolynomialBase::cTx
EntPolynomialBaseCtx * cTx
Definition: forces_and_sources_testing_users_base.cpp:88
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:346
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