v0.14.0
forces_and_sources_testing_flat_prism_element.cpp
/** \file forces_and_sources_testing_flat_prism_element.cpp
* \brief test for flat prism element
* \example forces_and_sources_testing_flat_prism_element.cpp
*
*/
#include <MoFEM.hpp>
namespace bio = boost::iostreams;
using bio::stream;
using bio::tee_device;
using namespace MoFEM;
static char help[] = "...\n\n";
int main(int argc, char *argv[]) {
MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
try {
moab::Core mb_instance;
moab::Interface &moab = mb_instance;
PetscBool flg = PETSC_TRUE;
char mesh_file_name[255];
#if PETSC_VERSION_GE(3, 6, 4)
CHKERR PetscOptionsGetString(PETSC_NULL, "", "-my_file", mesh_file_name,
255, &flg);
#else
CHKERR PetscOptionsGetString(PETSC_NULL, PETSC_NULL, "-my_file",
mesh_file_name, 255, &flg);
#endif
if (flg != PETSC_TRUE) {
SETERRQ(PETSC_COMM_SELF, 1, "*** ERROR -my_file (MESH FILE NEEDED)");
}
const char *option;
option = "";
CHKERR moab.load_file(mesh_file_name, 0, option);
// Create MoFEM (Joseph) database
MoFEM::Core core(moab);
MoFEM::Interface &m_field = core;
PrismInterface *interface;
CHKERR m_field.getInterface(interface);
// set entitities bit level
CHKERR m_field.getInterface<BitRefManager>()->setBitRefLevelByDim(
0, 3, BitRefLevel().set(0));
std::vector<BitRefLevel> bit_levels;
bit_levels.push_back(BitRefLevel().set(0));
int ll = 1;
// for(_IT_CUBITMESHSETS_BY_BCDATA_TYPE_FOR_LOOP_(m_field,SIDESET|INTERFACESET,cit))
// {
CHKERR PetscPrintf(PETSC_COMM_WORLD, "Insert Interface %d\n",
cit->getMeshsetId());
EntityHandle cubit_meshset = cit->getMeshset();
{
// get tet enties form back bit_level
EntityHandle ref_level_meshset = 0;
CHKERR moab.create_meshset(MESHSET_SET, ref_level_meshset);
;
->getEntitiesByTypeAndRefLevel(bit_levels.back(),
BitRefLevel().set(), MBTET,
ref_level_meshset);
->getEntitiesByTypeAndRefLevel(bit_levels.back(),
BitRefLevel().set(), MBPRISM,
ref_level_meshset);
Range ref_level_tets;
CHKERR moab.get_entities_by_handle(ref_level_meshset, ref_level_tets,
true);
;
// get faces and test to split
CHKERR interface->getSides(cubit_meshset, bit_levels.back(), true, 0);
// set new bit level
bit_levels.push_back(BitRefLevel().set(ll++));
// split faces and
CHKERR interface->splitSides(ref_level_meshset, bit_levels.back(),
cubit_meshset, true, true, 0);
// clean meshsets
CHKERR moab.delete_entities(&ref_level_meshset, 1);
;
}
// update cubit meshsets
for (_IT_CUBITMESHSETS_FOR_LOOP_(m_field, ciit)) {
EntityHandle cubit_meshset = ciit->meshset;
->updateMeshsetByEntitiesChildren(cubit_meshset, bit_levels.back(),
cubit_meshset, MBMAXTYPE, true);
}
}
// Fields
CHKERR m_field.add_field("FIELD1", H1, AINSWORTH_LEGENDRE_BASE, 3);
CHKERR m_field.add_field("MESH_NODE_POSITIONS", H1, AINSWORTH_LEGENDRE_BASE,
3);
CHKERR m_field.add_field("FIELD2", NOFIELD, NOBASE, 3);
{
// Creating and adding no field entities.
const double coords[] = {0, 0, 0};
EntityHandle no_field_vertex;
CHKERR m_field.get_moab().create_vertex(coords, no_field_vertex);
;
Range range_no_field_vertex;
range_no_field_vertex.insert(no_field_vertex);
CHKERR m_field.getInterface<BitRefManager>()->setBitRefLevel(
range_no_field_vertex, BitRefLevel().set());
EntityHandle meshset = m_field.get_field_meshset("FIELD2");
CHKERR m_field.get_moab().add_entities(meshset, range_no_field_vertex);
;
}
// FE
CHKERR m_field.add_finite_element("TEST_FE1");
CHKERR m_field.add_finite_element("TEST_FE2");
// Define rows/cols and element data
CHKERR m_field.modify_finite_element_add_field_row("TEST_FE1", "FIELD1");
CHKERR m_field.modify_finite_element_add_field_col("TEST_FE1", "FIELD1");
CHKERR m_field.modify_finite_element_add_field_data("TEST_FE1", "FIELD1");
"MESH_NODE_POSITIONS");
CHKERR m_field.modify_finite_element_add_field_row("TEST_FE2", "FIELD1");
// CHKERR m_field.modify_finite_element_add_field_row("TEST_FE2","FIELD2");
// CHKERR m_field.modify_finite_element_add_field_col("TEST_FE2","FIELD1");
CHKERR m_field.modify_finite_element_add_field_col("TEST_FE2", "FIELD2");
CHKERR m_field.modify_finite_element_add_field_data("TEST_FE2", "FIELD1");
CHKERR m_field.modify_finite_element_add_field_data("TEST_FE2", "FIELD2");
// Problem
CHKERR m_field.add_problem("TEST_PROBLEM");
// set finite elements for problem
CHKERR m_field.modify_problem_add_finite_element("TEST_PROBLEM",
"TEST_FE1");
CHKERR m_field.modify_problem_add_finite_element("TEST_PROBLEM",
"TEST_FE2");
// set refinement level for problem
CHKERR m_field.modify_problem_ref_level_add_bit("TEST_PROBLEM",
bit_levels.back());
// meshset consisting all entities in mesh
EntityHandle root_set = moab.get_root_set();
// add entities to field
CHKERR m_field.add_ents_to_field_by_type(root_set, MBTET, "FIELD1");
CHKERR m_field.add_ents_to_field_by_type(root_set, MBTET,
"MESH_NODE_POSITIONS");
// add entities to finite element
CHKERR m_field.add_ents_to_finite_element_by_type(root_set, MBPRISM,
"TEST_FE1", 10);
CHKERR m_field.add_ents_to_finite_element_by_type(root_set, MBPRISM,
"TEST_FE2", 10);
// set app. order
// see Hierarchic Finite Element Bases on Unstructured Tetrahedral Meshes
// (Mark Ainsworth & Joe Coyle)
int order = 3;
CHKERR m_field.set_field_order(root_set, MBTET, "FIELD1", order);
CHKERR m_field.set_field_order(root_set, MBTRI, "FIELD1", order);
CHKERR m_field.set_field_order(root_set, MBEDGE, "FIELD1", order);
CHKERR m_field.set_field_order(root_set, MBVERTEX, "FIELD1", 1);
CHKERR m_field.set_field_order(root_set, MBTET, "MESH_NODE_POSITIONS", 2);
CHKERR m_field.set_field_order(root_set, MBTRI, "MESH_NODE_POSITIONS", 2);
CHKERR m_field.set_field_order(root_set, MBEDGE, "MESH_NODE_POSITIONS", 2);
CHKERR m_field.set_field_order(root_set, MBVERTEX, "MESH_NODE_POSITIONS",
1);
/****/
// build database
// build field
CHKERR m_field.build_fields();
// set FIELD1 from positions of 10 node tets
Projection10NodeCoordsOnField ent_method_field1(m_field, "FIELD1");
CHKERR m_field.loop_dofs("FIELD1", ent_method_field1);
Projection10NodeCoordsOnField ent_method_mesh_positions(
m_field, "MESH_NODE_POSITIONS");
CHKERR m_field.loop_dofs("MESH_NODE_POSITIONS", ent_method_mesh_positions);
// build finite elemnts
// build adjacencies
CHKERR m_field.build_adjacencies(bit_levels.back());
// build problem
ProblemsManager *prb_mng_ptr;
CHKERR m_field.getInterface(prb_mng_ptr);
CHKERR prb_mng_ptr->buildProblem("TEST_PROBLEM", false);
/****/
// mesh partitioning
// partition
CHKERR prb_mng_ptr->partitionSimpleProblem("TEST_PROBLEM");
CHKERR prb_mng_ptr->partitionFiniteElements("TEST_PROBLEM");
// what are ghost nodes, see Petsc Manual
CHKERR prb_mng_ptr->partitionGhostDofs("TEST_PROBLEM");
typedef tee_device<std::ostream, std::ofstream> TeeDevice;
typedef stream<TeeDevice> TeeStream;
std::ofstream ofs("forces_and_sources_testing_flat_prism_element.txt");
TeeDevice my_tee(std::cout, ofs);
TeeStream my_split(my_tee);
struct MyOp
TeeStream &mySplit;
MyOp(TeeStream &mySplit, const char type)
"FIELD1", "FIELD1", type),
mySplit(mySplit) {}
MoFEMErrorCode doWork(int side, EntityType type,
if (data.getFieldData().empty())
const double eps = 1e-4;
for (DoubleAllocator::iterator it = getNormal().data().begin();
it != getNormal().data().end(); it++) {
*it = fabs(*it) < eps ? 0.0 : *it;
}
for (DoubleAllocator::iterator it =
getNormalsAtGaussPtsF3().data().begin();
it != getNormalsAtGaussPtsF3().data().end(); it++) {
*it = fabs(*it) < eps ? 0.0 : *it;
}
for (DoubleAllocator::iterator it =
getTangent1AtGaussPtF3().data().begin();
it != getTangent1AtGaussPtF3().data().end(); it++) {
*it = fabs(*it) < eps ? 0.0 : *it;
}
for (DoubleAllocator::iterator it =
getTangent2AtGaussPtF3().data().begin();
it != getTangent2AtGaussPtF3().data().end(); it++) {
*it = fabs(*it) < eps ? 0.0 : *it;
}
for (DoubleAllocator::iterator it =
getNormalsAtGaussPtsF4().data().begin();
it != getNormalsAtGaussPtsF4().data().end(); it++) {
*it = fabs(*it) < eps ? 0.0 : *it;
}
for (DoubleAllocator::iterator it =
getTangent1AtGaussPtF4().data().begin();
it != getTangent1AtGaussPtF4().data().end(); it++) {
*it = fabs(*it) < eps ? 0.0 : *it;
}
for (DoubleAllocator::iterator it =
getTangent2AtGaussPtF4().data().begin();
it != getTangent2AtGaussPtF4().data().end(); it++) {
*it = fabs(*it) < eps ? 0.0 : *it;
}
mySplit << "NH1" << std::endl;
mySplit << "side: " << side << " type: " << type << std::endl;
mySplit << data << std::endl;
mySplit << std::setprecision(3) << getCoords() << std::endl;
mySplit << std::setprecision(3) << getCoordsAtGaussPts() << std::endl;
mySplit << std::setprecision(3) << getArea(0) << std::endl;
mySplit << std::setprecision(3) << getArea(1) << std::endl;
mySplit << std::setprecision(3) << "normal F3 " << getNormalF3()
<< std::endl;
mySplit << std::setprecision(3) << "normal F4 " << getNormalF4()
<< std::endl;
mySplit << std::setprecision(3) << "normal at Gauss pt F3 "
<< getNormalsAtGaussPtsF3() << std::endl;
mySplit << std::setprecision(3) << getTangent1AtGaussPtF3()
<< std::endl;
mySplit << std::setprecision(3) << getTangent2AtGaussPtF3()
<< std::endl;
mySplit << std::setprecision(3) << "normal at Gauss pt F4 "
<< getNormalsAtGaussPtsF4() << std::endl;
mySplit << std::setprecision(3) << getTangent1AtGaussPtF4()
<< std::endl;
mySplit << std::setprecision(3) << getTangent2AtGaussPtF4()
<< std::endl;
}
MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type,
EntityType col_type,
if (row_data.getFieldData().empty())
mySplit << "NH1NH1" << std::endl;
mySplit << "row side: " << row_side << " row_type: " << row_type
<< std::endl;
mySplit << row_data << std::endl;
mySplit << "NH1NH1" << std::endl;
mySplit << "col side: " << col_side << " col_type: " << col_type
<< std::endl;
mySplit << row_data << std::endl;
}
};
struct MyOp2
TeeStream &mySplit;
MyOp2(TeeStream &my_split, const char type)
"FIELD1", "FIELD2", type),
mySplit(my_split) {}
MoFEMErrorCode doWork(int side, EntityType type,
if (type != MBENTITYSET)
mySplit << "NOFIELD" << std::endl;
mySplit << "side: " << side << " type: " << type << std::endl;
mySplit << data << std::endl;
}
MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type,
EntityType col_type,
unSetSymm();
if (col_type != MBENTITYSET)
mySplit << "NOFILEDH1" << std::endl;
mySplit << "row side: " << row_side << " row_type: " << row_type
<< std::endl;
mySplit << row_data << std::endl;
mySplit << "col side: " << col_side << " col_type: " << col_type
<< std::endl;
mySplit << col_data << std::endl;
}
};
fe1.getOpPtrVector().push_back(
fe1.getOpPtrVector().push_back(
CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "TEST_FE1", fe1);
fe2.getOpPtrVector().push_back(
fe2.getOpPtrVector().push_back(
CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "TEST_FE2", fe2);
}
return 0;
}
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::FlatPrismElementForcesAndSourcesCore::UserDataOperator
friend class UserDataOperator
Definition: FlatPrismElementForcesAndSourcesCore.hpp:167
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
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
MoFEM::CoreInterface::get_field_meshset
virtual EntityHandle get_field_meshset(const std::string name) const =0
get field meshset
MoFEM::ForcesAndSourcesCore::UserDataOperator::OPROWCOL
@ OPROWCOL
operator doWork is executed on FE rows &columns
Definition: ForcesAndSourcesCore.hpp:569
main
int main(int argc, char *argv[])
Definition: forces_and_sources_testing_flat_prism_element.cpp:19
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
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
help
static char help[]
Definition: forces_and_sources_testing_flat_prism_element.cpp:17
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::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::ForcesAndSourcesCore::UserDataOperator::OPCOL
@ OPCOL
operator doWork function is executed on FE columns
Definition: ForcesAndSourcesCore.hpp:568
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::FlatPrismElementForcesAndSourcesCore::UserDataOperator
default operator for Flat Prism element
Definition: FlatPrismElementForcesAndSourcesCore.hpp:34
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
mesh_file_name
char mesh_file_name[255]
Definition: mesh_smoothing.cpp:23
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::FlatPrismElementForcesAndSourcesCore
FlatPrism finite element.
Definition: FlatPrismElementForcesAndSourcesCore.hpp:27
CATCH_ERRORS
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:385
_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
eps
static const double eps
Definition: check_base_functions_derivatives_on_tet.cpp:11
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::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
_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
MoFEM::CoreInterface::set_field_order
virtual MoFEMErrorCode set_field_order(const EntityHandle meshset, const EntityType type, const std::string &name, const ApproximationOrder order, int verb=DEFAULT_VERBOSITY)=0
Set order approximation of the entities in the field.
MoFEM::ProblemsManager::partitionGhostDofs
MoFEMErrorCode partitionGhostDofs(const std::string name, int verb=VERBOSE)
determine ghost nodes
Definition: ProblemsManager.cpp:2339
MoFEM::CoreInterface::add_problem
virtual MoFEMErrorCode add_problem(const std::string &name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
Add problem.
NOFIELD
@ NOFIELD
scalar or vector of scalars describe (no true field)
Definition: definitions.h:84
MoFEM::ForcesAndSourcesCore::UserDataOperator::OPROW
@ OPROW
operator doWork function is executed on FE rows
Definition: ForcesAndSourcesCore.hpp:567