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

Go to the source code of this file.

Functions

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

Variables

static char help []
 

Function Documentation

◆ main()

int main ( int  argc,
char *  argv[] 
)
Examples
elasticity_mixed_formulation.cpp.

Definition at line 22 of file elasticity_mixed_formulation.cpp.

22  {
23 
24  const string default_options = "-ksp_type gmres \n"
25  "-pc_type lu \n"
26  "-pc_factor_mat_solver_type mumps \n"
27  "-mat_mumps_icntl_20 0 \n"
28  "-ksp_monitor\n";
29 
30  string param_file = "param_file.petsc";
31  if (!static_cast<bool>(ifstream(param_file))) {
32  std::ofstream file(param_file.c_str(), std::ios::ate);
33  if (file.is_open()) {
34  file << default_options;
35  file.close();
36  }
37  }
38 
39  MoFEM::Core::Initialize(&argc, &argv, param_file.c_str(), help);
40 
41  try {
42 
43  // Create mesh database
44  moab::Core mb_instance; // create database
45  moab::Interface &moab = mb_instance; // create interface to database
46 
47  // Create moab communicator
48  // Create separate MOAB communicator, it is duplicate of PETSc communicator.
49  // NOTE That this should eliminate potential communication problems between
50  // MOAB and PETSC functions.
51  ParallelComm *pcomm = ParallelComm::get_pcomm(&moab, MYPCOMM_INDEX);
52  auto moab_comm_wrap =
53  boost::make_shared<WrapMPIComm>(PETSC_COMM_WORLD, false);
54  if (pcomm == NULL)
55  pcomm = new ParallelComm(&moab, moab_comm_wrap->get_comm());
56 
57  // Get command line options
58  char mesh_file_name[255];
59  PetscBool flg_file;
60  int order_p = 2; // default approximation order_p
61  int order_u = 3; // default approximation order_u
62  PetscBool is_partitioned = PETSC_FALSE;
63  PetscBool flg_test = PETSC_FALSE; // true check if error is numerical error
64 
65  CHKERR PetscOptionsBegin(PETSC_COMM_WORLD, "", "Mix elastic problem",
66  "none");
67 
68  CHKERR PetscOptionsString("-my_file", "mesh file name", "", "mesh.h5m",
69  mesh_file_name, 255, &flg_file);
70  // Set approximation order
71  CHKERR PetscOptionsInt("-my_order_p", "approximation order_p", "", order_p,
72  &order_p, PETSC_NULL);
73  CHKERR PetscOptionsInt("-my_order_u", "approximation order_u", "", order_u,
74  &order_u, PETSC_NULL);
75 
76  CHKERR PetscOptionsBool("-is_partitioned", "is_partitioned?", "",
77  is_partitioned, &is_partitioned, PETSC_NULL);
78 
79  // Set testing (used by CTest)
80  CHKERR PetscOptionsBool("-test", "if true is ctest", "", flg_test,
81  &flg_test, PETSC_NULL);
82  ierr = PetscOptionsEnd();
83 
84  if (flg_file != PETSC_TRUE) {
85  SETERRQ(PETSC_COMM_SELF, 1, "*** ERROR -my_file (MESH FILE NEEDED)");
86  }
87 
88  // Read whole mesh or part of it if partitioned
89  if (is_partitioned == PETSC_TRUE) {
90  // This is a case of distributed mesh and algebra. In that case each
91  // processor keeps only part of the problem.
92  const char *option;
93  option = "PARALLEL=READ_PART;"
94  "PARALLEL_RESOLVE_SHARED_ENTS;"
95  "PARTITION=PARALLEL_PARTITION;";
96  CHKERR moab.load_file(mesh_file_name, 0, option);
97  } else {
98  // In this case, we have distributed algebra, i.e. assembly of vectors and
99  // matrices is in parallel, but whole mesh is stored on all processors.
100  // Solver and matrix scales well, however problem set-up of problem is
101  // not fully parallel.
102  const char *option;
103  option = "";
104  CHKERR moab.load_file(mesh_file_name, 0, option);
105  }
106 
107  // Create MoFEM database and link it to MoAB
108  MoFEM::Core core(moab);
109  MoFEM::Interface &m_field = core;
110 
111  // Print boundary conditions and material parameters
112  MeshsetsManager *meshsets_mng_ptr;
113  CHKERR m_field.getInterface(meshsets_mng_ptr);
114  CHKERR meshsets_mng_ptr->printDisplacementSet();
115  CHKERR meshsets_mng_ptr->printForceSet();
116  CHKERR meshsets_mng_ptr->printMaterialsSet();
117 
118  BitRefLevel bit_level0;
119  bit_level0.set(0);
120  CHKERR m_field.getInterface<BitRefManager>()->setBitRefLevelByDim(
121  0, 3, bit_level0);
122 
123  // **** ADD FIELDS **** //
124  CHKERR m_field.add_field("MESH_NODE_POSITIONS", H1, AINSWORTH_LEGENDRE_BASE,
125  3);
126  CHKERR m_field.add_ents_to_field_by_type(0, MBTET, "MESH_NODE_POSITIONS");
127  CHKERR m_field.set_field_order(0, MBVERTEX, "MESH_NODE_POSITIONS", 1);
128 
129  CHKERR m_field.add_field("U", H1, AINSWORTH_LEGENDRE_BASE, 3);
130  CHKERR m_field.add_ents_to_field_by_type(0, MBTET, "U");
131  CHKERR m_field.set_field_order(0, MBVERTEX, "U", 1);
132  CHKERR m_field.set_field_order(0, MBEDGE, "U", order_u);
133  CHKERR m_field.set_field_order(0, MBTRI, "U", order_u);
134  CHKERR m_field.set_field_order(0, MBTET, "U", order_u);
135 
136  CHKERR m_field.add_field("P", H1, AINSWORTH_LEGENDRE_BASE, 1);
137  CHKERR m_field.add_ents_to_field_by_type(0, MBTET, "P");
138  CHKERR m_field.set_field_order(0, MBVERTEX, "P", 1);
139  CHKERR m_field.set_field_order(0, MBEDGE, "P", order_p);
140  CHKERR m_field.set_field_order(0, MBTRI, "P", order_p);
141  CHKERR m_field.set_field_order(0, MBTET, "P", order_p);
142 
143  CHKERR m_field.build_fields();
144 
145  // CHKERR m_field.getInterface<FieldBlas>()->setField(
146  // 0, MBVERTEX, "P"); // initial p = 0 everywhere
147  {
148  Projection10NodeCoordsOnField ent_method_material(m_field,
149  "MESH_NODE_POSITIONS");
150  CHKERR m_field.loop_dofs("MESH_NODE_POSITIONS", ent_method_material);
151  }
152  // setup elements for loading
154  CHKERR MetaNodalForces::addElement(m_field, "U");
155  CHKERR MetaEdgeForces::addElement(m_field, "U");
156 
157  // **** ADD ELEMENTS **** //
158 
159  // Add finite element (this defines element, declaration comes later)
160  CHKERR m_field.add_finite_element("ELASTIC");
161  CHKERR m_field.modify_finite_element_add_field_row("ELASTIC", "U");
162  CHKERR m_field.modify_finite_element_add_field_col("ELASTIC", "U");
163  CHKERR m_field.modify_finite_element_add_field_data("ELASTIC", "U");
164 
165  CHKERR m_field.modify_finite_element_add_field_row("ELASTIC", "P");
166  CHKERR m_field.modify_finite_element_add_field_col("ELASTIC", "P");
167  CHKERR m_field.modify_finite_element_add_field_data("ELASTIC", "P");
169  "MESH_NODE_POSITIONS");
170 
171  // Add entities to that element
172  CHKERR m_field.add_ents_to_finite_element_by_type(0, MBTET, "ELASTIC");
173  // build finite elements
174  CHKERR m_field.build_finite_elements();
175  // build adjacencies between elements and degrees of freedom
176  CHKERR m_field.build_adjacencies(bit_level0);
177 
178  // **** BUILD DM **** //
179  DM dm;
180  DMType dm_name = "DM_ELASTIC_MIX";
181  // Register DM problem
182  CHKERR DMRegister_MoFEM(dm_name);
183  CHKERR DMCreate(PETSC_COMM_WORLD, &dm);
184  CHKERR DMSetType(dm, dm_name);
185  // Create DM instance
186  CHKERR DMMoFEMCreateMoFEM(dm, &m_field, dm_name, bit_level0);
187  // Configure DM form line command options (DM itself, solvers,
188  // pre-conditioners, ... )
189  CHKERR DMSetFromOptions(dm);
190  // Add elements to dm (only one here)
191  CHKERR DMMoFEMAddElement(dm, "ELASTIC");
192  if (m_field.check_finite_element("PRESSURE_FE"))
193  CHKERR DMMoFEMAddElement(dm, "PRESSURE_FE");
194  if (m_field.check_finite_element("FORCE_FE"))
195  CHKERR DMMoFEMAddElement(dm, "FORCE_FE");
196  CHKERR DMMoFEMSetIsPartitioned(dm, is_partitioned);
197  // setup the DM
198  CHKERR DMSetUp(dm);
199 
200  DataAtIntegrationPts commonData(m_field);
201  CHKERR commonData.getParameters();
202 
203  boost::shared_ptr<FEMethod> nullFE;
204  boost::shared_ptr<VolumeElementForcesAndSourcesCore> feLhs(
205  new VolumeElementForcesAndSourcesCore(m_field));
206  boost::shared_ptr<VolumeElementForcesAndSourcesCore> feRhs(
207  new VolumeElementForcesAndSourcesCore(m_field));
208 
210 
211  // loop over blocks
212  for (auto &sit : commonData.setOfBlocksData) {
213 
214  feLhs->getOpPtrVector().push_back(
215  new OpAssembleP(commonData, sit.second));
216  feLhs->getOpPtrVector().push_back(
217  new OpAssembleK(commonData, sit.second));
218  feLhs->getOpPtrVector().push_back(
219  new OpAssembleG(commonData, sit.second));
220 
221  auto u_ptr = boost::make_shared<MatrixDouble>();
222 
223  post_proc.getOpPtrVector().push_back(
224  new OpCalculateVectorFieldValues<3>("U", u_ptr));
225  post_proc.getOpPtrVector().push_back(
226  new OpCalculateScalarFieldValues("P", commonData.pPtr));
227  post_proc.getOpPtrVector().push_back(
229  commonData.gradDispPtr));
230 
232 
233  post_proc.getOpPtrVector().push_back(
234 
235  new OpPPMap(
236 
237  post_proc.getPostProcMesh(), post_proc.getMapGaussPts(),
238  {{"P", commonData.pPtr}},
239 
240  {{"U", u_ptr}},
241 
242  {},
243 
244  {}));
245 
246  post_proc.getOpPtrVector().push_back(new OpPostProcStress(
247  post_proc.getPostProcMesh(), post_proc.getMapGaussPts(), commonData,
248  sit.second));
249  }
250 
251  Mat Aij; // Stiffness matrix
252  Vec d, F_ext; // Vector of DOFs and the RHS
253 
254  {
256  CHKERR VecZeroEntries(d);
257  CHKERR VecDuplicate(d, &F_ext);
258  CHKERR DMCreateMatrix_MoFEM(dm, &Aij);
259  CHKERR VecGhostUpdateBegin(d, INSERT_VALUES, SCATTER_FORWARD);
260  CHKERR VecGhostUpdateEnd(d, INSERT_VALUES, SCATTER_FORWARD);
261  CHKERR DMoFEMMeshToLocalVector(dm, d, INSERT_VALUES, SCATTER_REVERSE);
262  CHKERR MatZeroEntries(Aij);
263  }
264 
265  // Assign global matrix/vector
266  feLhs->ksp_B = Aij;
267  feLhs->ksp_f = F_ext;
268 
269  boost::shared_ptr<DirichletDisplacementBc> dirichlet_bc_ptr(
270  new DirichletDisplacementBc(m_field, "U", Aij, d, F_ext));
271  dirichlet_bc_ptr->snes_ctx = FEMethod::CTX_SNESNONE;
272  dirichlet_bc_ptr->ts_ctx = FEMethod::CTX_TSNONE;
273  CHKERR DMoFEMPreProcessFiniteElements(dm, dirichlet_bc_ptr.get());
274  CHKERR DMoFEMLoopFiniteElements(dm, "ELASTIC", feLhs);
275 
276  // Assemble pressure and traction forces.
277  boost::ptr_map<std::string, NeumannForcesSurface> neumann_forces;
278  CHKERR MetaNeumannForces::setMomentumFluxOperators(m_field, neumann_forces,
279  F_ext, "U");
280 
281  {
282  boost::ptr_map<std::string, NeumannForcesSurface>::iterator mit =
283  neumann_forces.begin();
284  for (; mit != neumann_forces.end(); mit++) {
285  CHKERR DMoFEMLoopFiniteElements(dm, mit->first.c_str(),
286  &mit->second->getLoopFe());
287  }
288  }
289  // Assemble forces applied to nodes
290  boost::ptr_map<std::string, NodalForce> nodal_forces;
291  CHKERR MetaNodalForces::setOperators(m_field, nodal_forces, F_ext, "U");
292 
293  {
294  boost::ptr_map<std::string, NodalForce>::iterator fit =
295  nodal_forces.begin();
296  for (; fit != nodal_forces.end(); fit++) {
297  CHKERR DMoFEMLoopFiniteElements(dm, fit->first.c_str(),
298  &fit->second->getLoopFe());
299  }
300  }
301  // Assemble edge forces
302  boost::ptr_map<std::string, EdgeForce> edge_forces;
303  CHKERR MetaEdgeForces::setOperators(m_field, edge_forces, F_ext, "U");
304 
305  {
306  boost::ptr_map<std::string, EdgeForce>::iterator fit =
307  edge_forces.begin();
308  for (; fit != edge_forces.end(); fit++) {
309  CHKERR DMoFEMLoopFiniteElements(dm, fit->first.c_str(),
310  &fit->second->getLoopFe());
311  }
312  }
313 
314  CHKERR DMoFEMPostProcessFiniteElements(dm, dirichlet_bc_ptr.get());
315  CHKERR DMMoFEMKSPSetComputeOperators(dm, "ELASTIC", feLhs, nullFE, nullFE);
316  CHKERR VecAssemblyBegin(F_ext);
317  CHKERR VecAssemblyEnd(F_ext);
318 
319  // **** SOLVE **** //
320 
321  KSP solver;
322 
323  CHKERR KSPCreate(PETSC_COMM_WORLD, &solver);
324  CHKERR KSPSetFromOptions(solver);
325  CHKERR KSPSetOperators(solver, Aij, Aij);
326  CHKERR KSPSetUp(solver);
327  CHKERR KSPSolve(solver, F_ext, d);
328 
329  // VecView(F_ext, PETSC_VIEWER_STDOUT_WORLD);
330  // CHKERR VecView(d, PETSC_VIEWER_STDOUT_WORLD); // Print out the results
331  CHKERR DMoFEMMeshToGlobalVector(dm, d, INSERT_VALUES, SCATTER_REVERSE);
332 
333  // Save data on mesh
334  CHKERR DMoFEMLoopFiniteElements(dm, "ELASTIC", &post_proc);
335  PetscPrintf(PETSC_COMM_WORLD, "Output file: %s\n", "out.h5m");
336  CHKERR post_proc.writeFile("out.h5m");
337 
338  CHKERR MatDestroy(&Aij);
339  CHKERR VecDestroy(&d);
340  CHKERR VecDestroy(&F_ext);
341  CHKERR DMDestroy(&dm);
342  }
343  CATCH_ERRORS;
344 
345  // finish work cleaning memory, getting statistics, etc
347  return 0;
348 }

Variable Documentation

◆ help

char help[]
static
Initial value:
= "-my_order_p approximation order_p \n"
"-my_order_u approximation order_u \n"
"-is_partitioned is_partitioned? \n"
Examples
elasticity_mixed_formulation.cpp.

Definition at line 18 of file elasticity_mixed_formulation.cpp.

MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Definition: UnknownInterface.hpp:93
MYPCOMM_INDEX
#define MYPCOMM_INDEX
default communicator number PCOMM
Definition: definitions.h:228
MetaNeumannForces::addNeumannBCElements
static MoFEMErrorCode addNeumannBCElements(MoFEM::Interface &m_field, const std::string field_name, const std::string mesh_nodals_positions="MESH_NODE_POSITIONS", Range *intersect_ptr=NULL)
Declare finite element.
Definition: SurfacePressure.cpp:1974
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::CoreTmp< 0 >
Core (interface) class.
Definition: Core.hpp:82
H1
@ H1
continuous field
Definition: definitions.h:85
MetaNodalForces::addElement
static MoFEMErrorCode addElement(MoFEM::Interface &m_field, const std::string field_name, Range *intersect_ptr=NULL)
Add element taking information from NODESET.
Definition: NodalForce.hpp:92
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::OpCalculateVectorFieldValues
Get values at integration pts for tensor filed rank 1, i.e. vector field.
Definition: UserDataOperators.hpp:466
MoFEM::DMoFEMMeshToLocalVector
PetscErrorCode DMoFEMMeshToLocalVector(DM dm, Vec l, InsertMode mode, ScatterMode scatter_mode)
set local (or ghosted) vector values on mesh for partition only
Definition: DMMoFEM.cpp:523
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::DMMoFEMKSPSetComputeOperators
PetscErrorCode DMMoFEMKSPSetComputeOperators(DM dm, const char fe_name[], MoFEM::FEMethod *method, MoFEM::BasicMethod *pre_only, MoFEM::BasicMethod *post_only)
Set KSP operators and push mofem finite element methods.
Definition: DMMoFEM.cpp:678
DirichletDisplacementBc
Set Dirichlet boundary conditions on displacements.
Definition: DirichletBC.hpp:57
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::DMMoFEMAddElement
PetscErrorCode DMMoFEMAddElement(DM dm, std::string fe_name)
add element to dm
Definition: DMMoFEM.cpp:497
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
OpPostProcStress
Definition: ElasticityMixedFormulation.hpp:429
MoFEM::Interface
DeprecatedCoreInterface Interface
Definition: Interface.hpp:2010
MetaNodalForces::setOperators
static MoFEMErrorCode setOperators(MoFEM::Interface &m_field, boost::ptr_map< std::string, NodalForce > &nodal_forces, Vec F, const std::string field_name)
Set integration point operators.
Definition: NodalForce.hpp:128
MoFEM::DMCreateGlobalVector_MoFEM
PetscErrorCode DMCreateGlobalVector_MoFEM(DM dm, Vec *g)
DMShellSetCreateGlobalVector.
Definition: DMMoFEM.cpp:1167
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
MoFEM::DMCreateMatrix_MoFEM
PetscErrorCode DMCreateMatrix_MoFEM(DM dm, Mat *M)
Definition: DMMoFEM.cpp:1197
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
MoFEM::MeshsetsManager::printForceSet
MoFEMErrorCode printForceSet() const
print meshsets with force boundary conditions data structure
Definition: MeshsetsManager.cpp:304
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
OpAssembleK
Assemble K *.
Definition: ElasticityMixedFormulation.hpp:311
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
MetaNeumannForces::setMomentumFluxOperators
static MoFEMErrorCode setMomentumFluxOperators(MoFEM::Interface &m_field, boost::ptr_map< std::string, NeumannForcesSurface > &neumann_forces, Vec F, const std::string field_name, const std::string mesh_nodals_positions="MESH_NODE_POSITIONS")
Set operators to finite elements calculating right hand side vector.
Definition: SurfacePressure.cpp:2069
MoFEM::CoreInterface::build_finite_elements
virtual MoFEMErrorCode build_finite_elements(int verb=DEFAULT_VERBOSITY)=0
Build finite elements.
OpPPMap
OpPostProcMapInMoab< SPACE_DIM, SPACE_DIM > OpPPMap
Definition: photon_diffusion.cpp:29
MoFEM::DMoFEMMeshToGlobalVector
PetscErrorCode DMoFEMMeshToGlobalVector(DM dm, Vec g, InsertMode mode, ScatterMode scatter_mode)
set ghosted vector values on all existing mesh entities
Definition: DMMoFEM.cpp:535
MoFEM::OpCalculateScalarFieldValues
Get value at integration points for scalar field.
Definition: UserDataOperators.hpp:82
MoFEM::DMRegister_MoFEM
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
Definition: DMMoFEM.cpp:43
MoFEM::PostProcBrokenMeshInMoab< VolumeElementForcesAndSourcesCore >
Definition: PostProcBrokenMeshInMoabBase.hpp:670
OpAssembleP
Assemble P *.
Definition: ElasticityMixedFormulation.hpp:114
MoFEM::DMoFEMPreProcessFiniteElements
PetscErrorCode DMoFEMPreProcessFiniteElements(DM dm, MoFEM::FEMethod *method)
execute finite element method for each element in dm (problem)
Definition: DMMoFEM.cpp:546
MoFEM::DMoFEMPostProcessFiniteElements
PetscErrorCode DMoFEMPostProcessFiniteElements(DM dm, MoFEM::FEMethod *method)
execute finite element method for each element in dm (problem)
Definition: DMMoFEM.cpp:556
MoFEM::DMMoFEMCreateMoFEM
PetscErrorCode DMMoFEMCreateMoFEM(DM dm, MoFEM::Interface *m_field_ptr, const char problem_name[], const MoFEM::BitRefLevel bit_level, const MoFEM::BitRefLevel bit_mask=MoFEM::BitRefLevel().set())
Must be called by user to set MoFEM data structures.
Definition: DMMoFEM.cpp:114
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_field)=0
set finite element field data
MoFEM::VolumeElementForcesAndSourcesCore
Volume finite element base.
Definition: VolumeElementForcesAndSourcesCore.hpp:26
MoFEM::MeshsetsManager::printMaterialsSet
MoFEMErrorCode printMaterialsSet() const
print meshsets with material data structure set on it
Definition: MeshsetsManager.cpp:325
mesh_file_name
char mesh_file_name[255]
Definition: mesh_smoothing.cpp:22
MoFEM::OpCalculateVectorFieldGradient
Get field gradients at integration pts for scalar filed rank 0, i.e. vector field.
Definition: UserDataOperators.hpp:1535
MetaEdgeForces::addElement
static MoFEMErrorCode addElement(MoFEM::Interface &m_field, const std::string field_name, Range *intersect_ptr=NULL)
Add element taking information from NODESET.
Definition: EdgeForce.hpp:62
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
CATCH_ERRORS
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:385
MoFEM::MeshsetsManager::printDisplacementSet
MoFEMErrorCode printDisplacementSet() const
print meshsets with displacement boundary conditions data structure
Definition: MeshsetsManager.cpp:290
MoFEM::Core
CoreTmp< 0 > Core
Definition: Core.hpp:1148
help
static char help[]
Definition: elasticity_mixed_formulation.cpp:18
MoFEM::Exceptions::ierr
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Exceptions.hpp:76
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::MeshsetsManager
Interface for managing meshsets containing materials and boundary conditions.
Definition: MeshsetsManager.hpp:104
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::BitRefManager
Managing BitRefLevels.
Definition: BitRefManager.hpp:21
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
MoFEM::CoreInterface::check_finite_element
virtual bool check_finite_element(const std::string &name) const =0
Check if finite element is in database.
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::DMoFEMLoopFiniteElements
PetscErrorCode DMoFEMLoopFiniteElements(DM dm, const char fe_name[], MoFEM::FEMethod *method, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr())
Executes FEMethod for finite elements in DM.
Definition: DMMoFEM.cpp:586
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.
MoFEM::DMMoFEMSetIsPartitioned
PetscErrorCode DMMoFEMSetIsPartitioned(DM dm, PetscBool is_partitioned)
Definition: DMMoFEM.cpp:1123
OpAssembleG
Assemble G *.
Definition: ElasticityMixedFormulation.hpp:215
MetaEdgeForces::setOperators
static MoFEMErrorCode setOperators(MoFEM::Interface &m_field, boost::ptr_map< std::string, EdgeForce > &edge_forces, Vec F, const std::string field_name, std::string mesh_node_positions="MESH_NODE_POSITIONS")
Set integration point operators.
Definition: EdgeForce.hpp:97
MoFEM::OpPostProcMapInMoab
Post post-proc data at points from hash maps.
Definition: PostProcBrokenMeshInMoabBase.hpp:698