v0.14.0
forces_and_sources_testing_contact_prism_element.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 typedef tee_device<std::ostream, std::ofstream> TeeDevice;
10 typedef stream<TeeDevice> TeeStream;
11 
12 using namespace MoFEM;
13 
14 static char help[] = "...\n\n";
15 
17 
18  const char faceType;
19  MyOp(const char type, const char face_type)
21  "FIELD1", "FIELD1", type, face_type),
22  faceType(face_type) {}
23 
24  MoFEMErrorCode doWork(int side, EntityType type,
27 
28  if (data.getFieldData().empty())
30 
31  MOFEM_LOG("ATOM_TEST", Sev::inform) << "NH1";
32  MOFEM_LOG("ATOM_TEST", Sev::inform)
33  << "side: " << side << " type: " << type;
34  MOFEM_LOG("ATOM_TEST", Sev::inform) << data;
35 
36  if (faceType == FACEMASTER) {
37  MOFEM_LOG("ATOM_TEST", Sev::inform)
38  << "coords Master " << std::fixed << std::setprecision(2)
39  << getCoordsMaster();
40  MOFEM_LOG("ATOM_TEST", Sev::inform)
41  << "area Master " << std::fixed << std::setprecision(2)
42  << getAreaMaster();
43  MOFEM_LOG("ATOM_TEST", Sev::inform)
44  << "normal Master " << std::fixed << std::setprecision(2)
45  << getNormalMaster();
46  MOFEM_LOG("ATOM_TEST", Sev::inform)
47  << "coords at Gauss Pts Master " << std::fixed
48  << std::setprecision(2) << getCoordsAtGaussPtsMaster();
49  } else {
50  MOFEM_LOG("ATOM_TEST", Sev::inform)
51  << "coords Slave " << std::fixed << std::setprecision(2)
52  << getCoordsSlave();
53  MOFEM_LOG("ATOM_TEST", Sev::inform)
54  << "area Slave " << std::fixed << std::setprecision(2)
55  << getAreaSlave();
56  MOFEM_LOG("ATOM_TEST", Sev::inform)
57  << "normal Slave " << std::fixed << std::setprecision(2)
58  << getNormalSlave();
59  MOFEM_LOG("ATOM_TEST", Sev::inform)
60  << "coords at Gauss Pts Slave " << std::fixed
61  << std::setprecision(2) << getCoordsAtGaussPtsSlave();
62  }
64  }
65 
66  MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type,
67  EntityType col_type,
69  EntitiesFieldData::EntData &col_data) {
71 
72  if (row_data.getFieldData().empty())
74 
75  MOFEM_LOG("ATOM_TEST", Sev::inform) << "NH1NH1";
76  MOFEM_LOG("ATOM_TEST", Sev::inform)
77  << "row side: " << row_side << " row_type: " << row_type;
78  MOFEM_LOG("ATOM_TEST", Sev::inform) << row_data;
79  MOFEM_LOG("ATOM_TEST", Sev::inform) << "NH1NH1";
80  MOFEM_LOG("ATOM_TEST", Sev::inform)
81  << "col side: " << col_side << " col_type: " << col_type;
82  MOFEM_LOG("ATOM_TEST", Sev::inform) << col_data;
83 
85  }
86 };
87 
89 
90  CallingOp(const char type)
91  : ForcesAndSourcesCore::UserDataOperator("FIELD1", "FIELD1", type) {}
92 
93  MoFEMErrorCode doWork(int side, EntityType type,
96 
97  if (data.getFieldData().empty())
99 
100  MOFEM_LOG("ATOM_TEST", Sev::inform) << "Calling Operator NH1";
101  MOFEM_LOG("ATOM_TEST", Sev::inform)
102  << "side: " << side << " type: " << type;
103  MOFEM_LOG("ATOM_TEST", Sev::inform) << data;
104 
106  }
107 
108  MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type,
109  EntityType col_type,
110  EntitiesFieldData::EntData &row_data,
111  EntitiesFieldData::EntData &col_data) {
113 
114  if (row_data.getFieldData().empty())
116 
117  MOFEM_LOG("ATOM_TEST", Sev::inform) << "Calling Operator NH1NH1";
118  MOFEM_LOG("ATOM_TEST", Sev::inform)
119  << "row side: " << row_side << " row_type: " << row_type;
120  MOFEM_LOG("ATOM_TEST", Sev::inform) << row_data;
121  MOFEM_LOG("ATOM_TEST", Sev::inform) << "NH1NH1";
122  MOFEM_LOG("ATOM_TEST", Sev::inform)
123  << "col side: " << col_side << " col_type: " << col_type;
124  MOFEM_LOG("ATOM_TEST", Sev::inform) << col_data;
125 
127  }
128 };
129 
130 struct MyOp2
132 
133  MyOp2(const char type, const char face_type)
135  "FIELD1", "FIELD2", type, face_type) {}
136 
137  MoFEMErrorCode doWork(int side, EntityType type,
140 
141  if (type != MBENTITYSET)
143 
144  MOFEM_LOG("ATOM_TEST", Sev::inform) << "NPFIELD";
145  MOFEM_LOG("ATOM_TEST", Sev::inform)
146  << "side: " << side << " type: " << type;
147  MOFEM_LOG("ATOM_TEST", Sev::inform) << data;
149  }
150 
151  MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type,
152  EntityType col_type,
153  EntitiesFieldData::EntData &row_data,
154  EntitiesFieldData::EntData &col_data) {
156 
157  unSetSymm();
158 
159  if (col_type != MBENTITYSET)
161 
162  MOFEM_LOG("ATOM_TEST", Sev::inform) << "NOFILEDH1";
163  MOFEM_LOG("ATOM_TEST", Sev::inform)
164  << "row side: " << row_side << " row_type: " << row_type;
165  MOFEM_LOG("ATOM_TEST", Sev::inform) << row_data;
166  MOFEM_LOG("ATOM_TEST", Sev::inform)
167  << "col side: " << col_side << " col_type: " << col_type;
168  MOFEM_LOG("ATOM_TEST", Sev::inform) << col_data;
169 
171  }
172 };
173 
174 int main(int argc, char *argv[]) {
175 
176  MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
177 
178  try {
179 
180  moab::Core mb_instance;
181  moab::Interface &moab = mb_instance;
182 
183  PetscBool flg = PETSC_TRUE;
184  char mesh_file_name[255];
185  CHKERR PetscOptionsGetString(PETSC_NULL, "", "-my_file", mesh_file_name,
186  255, &flg);
187  if (flg != PETSC_TRUE)
188  SETERRQ(PETSC_COMM_SELF, MOFEM_INVALID_DATA, "error -my_file (mesh file not given");
189 
190  enum spaces { MYH1, MYHDIV, MYLASTSPACEOP };
191  const char *list_spaces[] = {"h1", "hdiv"};
192  PetscInt choice_space_value = H1;
193  CHKERR PetscOptionsGetEList(PETSC_NULL, NULL, "-space", list_spaces,
194  MYLASTSPACEOP, &choice_space_value, &flg);
195  bool is_hdiv = (choice_space_value == MYHDIV) ? true : false;
196 
197  const char *option;
198  option = "";
199  CHKERR moab.load_file(mesh_file_name, 0, option);
200 
201  // Create channel "ATOM_TEST" and add sink to file for that channel
202  auto add_atom_logging = [is_hdiv]() {
204  auto get_log_file_name = [is_hdiv]() {
205  if (is_hdiv)
206  return "forces_and_sources_testing_contact_prism_element_HDIV.txt";
207  else
208  return "forces_and_sources_testing_contact_prism_element.txt";
209  };
210 
211  auto core_log = logging::core::get();
212  core_log->add_sink(
214  LogManager::setLog("ATOM_TEST");
215  MOFEM_LOG_TAG("ATOM_TEST", "atom test");
216 
217  // Add log to file
218  logging::add_file_log(keywords::file_name = get_log_file_name(),
219  keywords::filter =
220  MoFEM::LogKeywords::channel == "ATOM_TEST");
222  };
223 
224  CHKERR add_atom_logging();
225 
226  // Create MoFEM (Joseph) database
227  MoFEM::Core core(moab);
228  MoFEM::Interface &m_field = core;
229  PrismInterface *interface;
230  CHKERR m_field.getInterface(interface);
231 
232  // set entities bit level
233  CHKERR m_field.getInterface<BitRefManager>()->setBitRefLevelByDim(
234  0, 3, BitRefLevel().set(0));
235  std::vector<BitRefLevel> bit_levels;
236  bit_levels.push_back(BitRefLevel().set(0));
237 
238  int ll = 1;
239  for (_IT_CUBITMESHSETS_BY_SET_TYPE_FOR_LOOP_(m_field, SIDESET, cit)) {
240 
241  CHKERR PetscPrintf(PETSC_COMM_WORLD, "Insert Interface %d\n",
242  cit->getMeshsetId());
243  EntityHandle cubit_meshset = cit->getMeshset();
244  {
245  // get tet enties form back bit_level
246  EntityHandle ref_level_meshset = 0;
247  CHKERR moab.create_meshset(MESHSET_SET, ref_level_meshset);
249  ->getEntitiesByTypeAndRefLevel(bit_levels.back(),
250  BitRefLevel().set(), MBTET,
251  ref_level_meshset);
253  ->getEntitiesByTypeAndRefLevel(bit_levels.back(),
254  BitRefLevel().set(), MBPRISM,
255  ref_level_meshset);
256  Range ref_level_tets;
257  CHKERR moab.get_entities_by_handle(ref_level_meshset, ref_level_tets,
258  true);
259  // get faces and test to split
260  CHKERR interface->getSides(cubit_meshset, bit_levels.back(), true, 0);
261  // set new bit level
262  bit_levels.push_back(BitRefLevel().set(ll++));
263  // split faces and
264  CHKERR interface->splitSides(ref_level_meshset, bit_levels.back(),
265  cubit_meshset, true, true, 0);
266  // clean meshsets
267  CHKERR moab.delete_entities(&ref_level_meshset, 1);
268  }
269  // update cubit meshsets
270  for (_IT_CUBITMESHSETS_FOR_LOOP_(m_field, ciit)) {
271  EntityHandle cubit_meshset = ciit->meshset;
273  ->updateMeshsetByEntitiesChildren(cubit_meshset, bit_levels.back(),
274  cubit_meshset, MBMAXTYPE, true);
275  }
276  }
277 
278  // Fields
279  if (is_hdiv)
280  CHKERR m_field.add_field("FIELD1", HDIV, DEMKOWICZ_JACOBI_BASE, 3);
281  else
282  CHKERR m_field.add_field("FIELD1", H1, AINSWORTH_LEGENDRE_BASE, 3);
283 
284  CHKERR m_field.add_field("MESH_NODE_POSITIONS", H1, AINSWORTH_LEGENDRE_BASE,
285  3);
286  CHKERR m_field.add_field("FIELD2", NOFIELD, NOBASE, 3);
287 
288  auto set_no_field_vertex = [&]() {
290  // Creating and adding no field entities.
291  const double coords[] = {0, 0, 0};
292  EntityHandle no_field_vertex;
293  CHKERR m_field.get_moab().create_vertex(coords, no_field_vertex);
294  Range range_no_field_vertex;
295  range_no_field_vertex.insert(no_field_vertex);
296  CHKERR m_field.getInterface<BitRefManager>()->setBitRefLevel(
297  range_no_field_vertex, BitRefLevel().set());
298  EntityHandle meshset = m_field.get_field_meshset("FIELD2");
299  CHKERR m_field.get_moab().add_entities(meshset, range_no_field_vertex);
301  };
302 
303  CHKERR set_no_field_vertex();
304 
305  // FE
306  CHKERR m_field.add_finite_element("TEST_FE1");
307  CHKERR m_field.add_finite_element("TEST_FE2");
308 
309  // Define rows/cols and element data
310  CHKERR m_field.modify_finite_element_add_field_row("TEST_FE1", "FIELD1");
311  CHKERR m_field.modify_finite_element_add_field_col("TEST_FE1", "FIELD1");
312  CHKERR m_field.modify_finite_element_add_field_data("TEST_FE1", "FIELD1");
313  CHKERR m_field.modify_finite_element_add_field_data("TEST_FE1",
314  "MESH_NODE_POSITIONS");
315 
316  CHKERR m_field.modify_finite_element_add_field_row("TEST_FE2", "FIELD1");
317  CHKERR m_field.modify_finite_element_add_field_col("TEST_FE2", "FIELD2");
318  CHKERR m_field.modify_finite_element_add_field_data("TEST_FE2", "FIELD1");
319  CHKERR m_field.modify_finite_element_add_field_data("TEST_FE2", "FIELD2");
320 
321  // Problem
322  CHKERR m_field.add_problem("TEST_PROBLEM");
323 
324  // set finite elements for problem
325  CHKERR m_field.modify_problem_add_finite_element("TEST_PROBLEM",
326  "TEST_FE1");
327  CHKERR m_field.modify_problem_add_finite_element("TEST_PROBLEM",
328  "TEST_FE2");
329  // set refinement level for problem
330  CHKERR m_field.modify_problem_ref_level_add_bit("TEST_PROBLEM",
331  bit_levels.back());
332 
333  // meshset consisting all entities in mesh
334  EntityHandle root_set = moab.get_root_set();
335  // add entities to field
336  CHKERR m_field.add_ents_to_field_by_type(root_set, MBTET, "FIELD1");
337  CHKERR m_field.add_ents_to_field_by_type(root_set, MBTET,
338  "MESH_NODE_POSITIONS");
339  // add entities to finite element
340  CHKERR m_field.add_ents_to_finite_element_by_type(root_set, MBPRISM,
341  "TEST_FE1");
342  CHKERR m_field.add_ents_to_finite_element_by_type(root_set, MBPRISM,
343  "TEST_FE2");
344 
345  // set app. order
346  constexpr int order = 3;
347 
348  if (is_hdiv) {
349  CHKERR m_field.set_field_order(root_set, MBTET, "FIELD1", 0);
350  CHKERR m_field.set_field_order(root_set, MBTRI, "FIELD1", order);
351  } else {
352  CHKERR m_field.set_field_order(root_set, MBTET, "FIELD1", order);
353  CHKERR m_field.set_field_order(root_set, MBTRI, "FIELD1", order);
354  CHKERR m_field.set_field_order(root_set, MBEDGE, "FIELD1", order);
355  CHKERR m_field.set_field_order(root_set, MBVERTEX, "FIELD1", 1);
356  }
357 
358  CHKERR m_field.set_field_order(root_set, MBTET, "MESH_NODE_POSITIONS", 2);
359  CHKERR m_field.set_field_order(root_set, MBTRI, "MESH_NODE_POSITIONS", 2);
360  CHKERR m_field.set_field_order(root_set, MBEDGE, "MESH_NODE_POSITIONS", 2);
361  CHKERR m_field.set_field_order(root_set, MBVERTEX, "MESH_NODE_POSITIONS",
362  1);
363 
364  // build field
365  CHKERR m_field.build_fields();
366 
367  if (!is_hdiv) {
368  // set FIELD1 from positions of 10 node tets
369  Projection10NodeCoordsOnField ent_method_field1(m_field, "FIELD1");
370  CHKERR m_field.loop_dofs("FIELD1", ent_method_field1);
371  }
372  Projection10NodeCoordsOnField ent_method_mesh_positions(
373  m_field, "MESH_NODE_POSITIONS");
374  CHKERR m_field.loop_dofs("MESH_NODE_POSITIONS", ent_method_mesh_positions);
375 
376  // build finite elemnts
377  CHKERR m_field.build_finite_elements();
378  // build adjacencies
379  CHKERR m_field.build_adjacencies(bit_levels.back());
380  // build problem
381  ProblemsManager *prb_mng_ptr;
382  CHKERR m_field.getInterface(prb_mng_ptr);
383  CHKERR prb_mng_ptr->buildProblem("TEST_PROBLEM", false);
384 
385  // mesh partitioning
386  // partition
387  CHKERR prb_mng_ptr->partitionSimpleProblem("TEST_PROBLEM");
388  CHKERR prb_mng_ptr->partitionFiniteElements("TEST_PROBLEM");
389  // what are ghost nodes, see Petsc Manual
390  CHKERR prb_mng_ptr->partitionGhostDofs("TEST_PROBLEM");
391 
393  using ContactDataOp =
395 
397  fe1.getOpPtrVector().push_back(
398  new MyOp(UMDataOp::OPROW, ContactDataOp::FACEMASTER));
399  fe1.getOpPtrVector().push_back(
400  new MyOp(UMDataOp::OPROW, ContactDataOp::FACESLAVE));
401  fe1.getOpPtrVector().push_back(
402  new MyOp(UMDataOp::OPROWCOL, ContactDataOp::FACEMASTERMASTER));
403  fe1.getOpPtrVector().push_back(
404  new MyOp(UMDataOp::OPROWCOL, ContactDataOp::FACEMASTERSLAVE));
405  fe1.getOpPtrVector().push_back(
406  new MyOp(UMDataOp::OPROWCOL, ContactDataOp::FACESLAVEMASTER));
407  fe1.getOpPtrVector().push_back(
408  new MyOp(UMDataOp::OPROWCOL, ContactDataOp::FACESLAVESLAVE));
409  fe1.getOpPtrVector().push_back(new CallingOp(UMDataOp::OPCOL));
410  fe1.getOpPtrVector().push_back(new CallingOp(UMDataOp::OPROW));
411  fe1.getOpPtrVector().push_back(new CallingOp(UMDataOp::OPROWCOL));
412  CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "TEST_FE1", fe1);
413 
415  fe2.getOpPtrVector().push_back(
416  new MyOp2(UMDataOp::OPCOL, ContactDataOp::FACEMASTER));
417  fe2.getOpPtrVector().push_back(
418  new MyOp2(UMDataOp::OPCOL, ContactDataOp::FACESLAVE));
419  fe2.getOpPtrVector().push_back(
420  new MyOp2(UMDataOp::OPROWCOL, ContactDataOp::FACEMASTERMASTER));
421  fe2.getOpPtrVector().push_back(
422  new MyOp2(UMDataOp::OPROWCOL, ContactDataOp::FACEMASTERSLAVE));
423  fe2.getOpPtrVector().push_back(
424  new MyOp2(UMDataOp::OPROWCOL, ContactDataOp::FACESLAVEMASTER));
425  fe2.getOpPtrVector().push_back(
426  new MyOp2(UMDataOp::OPROWCOL, ContactDataOp::FACESLAVESLAVE));
427  CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "TEST_FE2", fe2);
428  }
429  CATCH_ERRORS;
430 
432 
433  return 0;
434 }
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
SIDESET
@ SIDESET
Definition: definitions.h:160
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::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
MyOp2::MyOp2
MyOp2(const char type, const char face_type)
Definition: forces_and_sources_testing_contact_prism_element.cpp:133
MoFEM::ProblemsManager
Problem manager is used to build and partition problems.
Definition: ProblemsManager.hpp:21
MoFEM::PrismInterface
Create interface from given surface and insert flat prisms in-between.
Definition: PrismInterface.hpp:23
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::EntitiesFieldData::EntData::getFieldData
const VectorDouble & getFieldData() const
get dofs values
Definition: EntitiesFieldData.hpp:1244
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
MyOp::faceType
const char faceType
Definition: forces_and_sources_testing_contact_prism_element.cpp:18
MoFEM::CoreInterface::get_field_meshset
virtual EntityHandle get_field_meshset(const std::string name) const =0
get field meshset
MoFEM::LogManager::createSink
static boost::shared_ptr< SinkType > createSink(boost::shared_ptr< std::ostream > stream_ptr, std::string comm_filter)
Create a sink object.
Definition: LogManager.cpp:298
MoFEM::CoreInterface::add_ents_to_field_by_type
virtual MoFEMErrorCode add_ents_to_field_by_type(const Range &ents, const EntityType type, const std::string &name, int verb=DEFAULT_VERBOSITY)=0
Add entities to field meshset.
order
constexpr int order
Definition: dg_projection.cpp:18
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
TeeDevice
tee_device< std::ostream, std::ofstream > TeeDevice
Definition: forces_and_sources_testing_contact_prism_element.cpp:9
MoFEM::Interface
DeprecatedCoreInterface Interface
Definition: Interface.hpp:2002
MyOp::doWork
MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
Definition: forces_and_sources_testing_contact_prism_element.cpp:66
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
CallingOp
Definition: forces_and_sources_testing_contact_prism_element.cpp:88
MoFEM::ContactPrismElementForcesAndSourcesCore
ContactPrism finite element.
Definition: ContactPrismElementForcesAndSourcesCore.hpp:27
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::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::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.
MyOp::MyOp
MyOp(const char type, const char face_type)
Definition: forces_and_sources_testing_contact_prism_element.cpp:19
convert.type
type
Definition: convert.py:64
MoFEM::PrismInterface::splitSides
MoFEMErrorCode splitSides(const EntityHandle meshset, const BitRefLevel &bit, const int msId, const CubitBCType cubit_bc_type, const bool add_interface_entities, const bool recursive=false, int verb=QUIET)
Split nodes and other entities of tetrahedra on both sides of the interface and insert flat prisms in...
Definition: PrismInterface.cpp:519
MoFEM::ContactPrismElementForcesAndSourcesCore::UserDataOperator
default operator for Contact Prism element
Definition: ContactPrismElementForcesAndSourcesCore.hpp:208
MoFEM::LogManager::getStrmSelf
static boost::shared_ptr< std::ostream > getStrmSelf()
Get the strm self object.
Definition: LogManager.cpp:340
MOFEM_LOG_TAG
#define MOFEM_LOG_TAG(channel, tag)
Tag channel.
Definition: LogManager.hpp:339
MoFEM::ForcesAndSourcesCore::UserDataOperator
friend class UserDataOperator
Definition: ForcesAndSourcesCore.hpp:482
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
MyOp::doWork
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Definition: forces_and_sources_testing_contact_prism_element.cpp:24
mesh_file_name
char mesh_file_name[255]
Definition: mesh_smoothing.cpp:23
MoFEM::ForcesAndSourcesCore
structure to get information form mofem into EntitiesFieldData
Definition: ForcesAndSourcesCore.hpp:22
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
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
_IT_CUBITMESHSETS_BY_SET_TYPE_FOR_LOOP_
#define _IT_CUBITMESHSETS_BY_SET_TYPE_FOR_LOOP_(MESHSET_MANAGER, CUBITBCTYPE, IT)
Iterator that loops over a specific Cubit MeshSet having a particular BC meshset in a moFEM field.
Definition: MeshsetsManager.hpp:71
MoFEM::Core
CoreTmp< 0 > Core
Definition: Core.hpp:1140
MoFEM::ProblemsManager::partitionSimpleProblem
MoFEMErrorCode partitionSimpleProblem(const std::string name, int verb=VERBOSE)
partition problem dofs
Definition: ProblemsManager.cpp:1537
MoFEM::PrismInterface::getSides
MoFEMErrorCode getSides(const int msId, const CubitBCType cubit_bc_type, const BitRefLevel mesh_bit_level, const bool recursive, int verb=QUIET)
Store tetrahedra from each side of the interface separately in two child meshsets of the parent meshs...
Definition: PrismInterface.cpp:56
MyOp2
Definition: forces_and_sources_testing_contact_prism_element.cpp:130
AINSWORTH_LEGENDRE_BASE
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:60
MyOp
Operator used to check consistency between local coordinates and global cooridnates for integrated po...
Definition: field_evaluator.cpp:21
MoFEM::ForcesAndSourcesCore::getOpPtrVector
boost::ptr_deque< UserDataOperator > & getOpPtrVector()
Use to push back operator for row operator.
Definition: ForcesAndSourcesCore.hpp:83
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::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
MyOp2::doWork
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
Definition: forces_and_sources_testing_contact_prism_element.cpp:137
_IT_CUBITMESHSETS_FOR_LOOP_
#define _IT_CUBITMESHSETS_FOR_LOOP_(MESHSET_MANAGER, IT)
Iterator that loops over all the Cubit MeshSets in a moFEM field.
Definition: MeshsetsManager.hpp:34
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
CallingOp::doWork
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
Definition: forces_and_sources_testing_contact_prism_element.cpp:93
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::LogManager::setLog
static LoggerType & setLog(const std::string channel)
Set ans resset chanel logger.
Definition: LogManager.cpp:389
MyOp2::doWork
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: forces_and_sources_testing_contact_prism_element.cpp:151
MoFEM::ProblemsManager::partitionGhostDofs
MoFEMErrorCode partitionGhostDofs(const std::string name, int verb=VERBOSE)
determine ghost nodes
Definition: ProblemsManager.cpp:2339
CallingOp::CallingOp
CallingOp(const char type)
Definition: forces_and_sources_testing_contact_prism_element.cpp:90
main
int main(int argc, char *argv[])
Definition: forces_and_sources_testing_contact_prism_element.cpp:174
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:429
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:359
MOFEM_INVALID_DATA
@ MOFEM_INVALID_DATA
Definition: definitions.h:36
CallingOp::doWork
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: forces_and_sources_testing_contact_prism_element.cpp:108
NOFIELD
@ NOFIELD
scalar or vector of scalars describe (no true field)
Definition: definitions.h:84
help
static char help[]
Definition: forces_and_sources_testing_contact_prism_element.cpp:14