v0.14.0
continuity_check_on_contact_prism_side_ele.cpp

Testing integration on volume side on contact elementChecking continuity of H1 ( later hdiv and hcurl spaces on faces), and testing methods for integration on volume side on contact element.

/** \file continuity_check_on_contact_prism_side_ele.cpp
* \example continuity_check_on_contact_prism_side_ele.cpp
*
* \brief Testing integration on volume side on contact element
*
* Checking continuity of H1 ( later hdiv and hcurl spaces on faces), and
* testing methods for integration on volume side on contact element.
*
*/
#include <MoFEM.hpp>
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;
int rank;
MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
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 database
MoFEM::Core core(moab);
// Access to database through interface
MoFEM::Interface &m_field = core;
enum side_contact { MASTERSIDE, SLAVESIDE };
// Fields
auto get_base = []() -> FieldApproximationBase {
enum bases { AINSWORTH, DEMKOWICZ, LASTBASEOP };
const char *list_bases[] = {"ainsworth", "demkowicz"};
PetscBool flg;
PetscInt choice_base_value = AINSWORTH;
CHKERR PetscOptionsGetEList(PETSC_NULL, NULL, "-base", list_bases,
LASTBASEOP, &choice_base_value, &flg);
if (flg == PETSC_TRUE) {
if (choice_base_value == AINSWORTH)
else if (choice_base_value == DEMKOWICZ)
return base;
}
return LASTBASE;
};
auto get_space = []() -> FieldSpace {
enum spaces { H1, HDIV, HCURL, LASTSPACEOP };
const char *list_spaces[] = {"h1", "hdiv", "hcurl"};
PetscBool flg;
PetscInt choice_space_value = HDIV;
CHKERR PetscOptionsGetEList(PETSC_NULL, NULL, "-space", list_spaces,
LASTSPACEOP, &choice_space_value, &flg);
if (flg == PETSC_TRUE) {
if (choice_space_value == HDIV)
else if (choice_space_value == HCURL)
else if (choice_space_value == H1)
space = FieldSpace::H1;
return space;
}
return LASTSPACE;
};
CHKERR m_field.add_field("F2", get_space(), get_base(), 1);
auto add_prism_interface = [&](Range &tets, Range &prisms,
Range &master_tris, Range &slave_tris,
EntityHandle &meshset_tets,
EntityHandle &meshset_prisms,
std::vector<BitRefLevel> &bit_levels) {
PrismInterface *interface;
CHKERR m_field.getInterface(interface);
int ll = 1;
if (cit->getName().compare(0, 11, "INT_CONTACT") == 0) {
CHKERR PetscPrintf(PETSC_COMM_WORLD, "Insert %s (id: %d)\n",
cit->getName().c_str(), cit->getMeshsetId());
EntityHandle cubit_meshset = cit->getMeshset();
Range tris;
CHKERR moab.get_entities_by_type(cubit_meshset, MBTRI, tris, true);
master_tris.merge(tris);
{
// get tet entities from 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 tets 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 tets
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);
}
}
}
for (unsigned int ll = 0; ll != bit_levels.size() - 1; ll++) {
CHKERR m_field.delete_ents_by_bit_ref(bit_levels[ll], bit_levels[ll],
true);
}
CHKERR m_field.getInterface<BitRefManager>()->shiftRightBitRef(
bit_levels.size() - 1);
CHKERR moab.create_meshset(MESHSET_SET, meshset_tets);
CHKERR moab.create_meshset(MESHSET_SET, meshset_prisms);
->getEntitiesByTypeAndRefLevel(bit_levels[0], BitRefLevel().set(),
MBTET, meshset_tets);
CHKERR moab.get_entities_by_handle(meshset_tets, tets, true);
->getEntitiesByTypeAndRefLevel(bit_levels[0], BitRefLevel().set(),
MBPRISM, meshset_prisms);
CHKERR moab.get_entities_by_handle(meshset_prisms, prisms);
Range tris;
CHKERR moab.get_adjacencies(prisms, 2, false, tris,
moab::Interface::UNION);
tris = tris.subset_by_type(MBTRI);
slave_tris = subtract(tris, master_tris);
};
Range all_tets, contact_prisms, master_tris, slave_tris;
EntityHandle meshset_tets, meshset_prisms;
std::vector<BitRefLevel> bit_levels;
bit_levels.push_back(BitRefLevel().set(0));
CHKERR m_field.getInterface<BitRefManager>()->setBitRefLevelByDim(
0, 3, bit_levels[0]);
CHKERR add_prism_interface(all_tets, contact_prisms, master_tris,
slave_tris, meshset_tets, meshset_prisms,
bit_levels);
// 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, "F2");
CHKERR m_field.add_ents_to_field_by_type(root_set, MBPRISM, "F2");
// set app. order
int order = 2;
CHKERR m_field.set_field_order(root_set, MBTET, "F2", order);
CHKERR m_field.set_field_order(root_set, MBTRI, "F2", order);
CHKERR m_field.set_field_order(root_set, MBEDGE, "F2", order);
if (get_space() == FieldSpace::H1) {
CHKERR m_field.set_field_order(root_set, MBVERTEX, "F2", 1);
}
CHKERR m_field.build_fields();
// add elements
CHKERR m_field.add_finite_element("V1");
CHKERR m_field.add_finite_element("C2");
CHKERR m_field.add_ents_to_finite_element_by_type(root_set, MBTET, "V1");
Range prism;
CHKERR m_field.getInterface<BitRefManager>()->getEntitiesByTypeAndRefLevel(
bit_levels[0], BitRefLevel().set(), MBPRISM, prism);
CHKERR m_field.add_ents_to_finite_element_by_type(root_set, MBPRISM, "C2");
CHKERR m_field.build_adjacencies(bit_levels[0]);
// Problems
CHKERR m_field.add_problem("P1");
// set refinement level for problem
CHKERR m_field.modify_problem_ref_level_add_bit("P1", bit_levels[0]);
// build problems
ProblemsManager *prb_mng_ptr;
CHKERR m_field.getInterface(prb_mng_ptr);
CHKERR prb_mng_ptr->buildProblem("P1", true);
struct CommonData {
MatrixDouble dotNormalFace;
MatrixDouble dotNormalEleLeft;
MatrixDouble dotNormalEleRight;
MatrixDouble shapeFunH1Values;
MatrixDouble shapeFunH1VolSide;
};
struct OnContactSideMaster
CommonData &elemData;
OpVolSide(CommonData &elem_data)
elemData(elem_data) {}
MoFEMErrorCode doWork(int side, EntityType type,
if (type == MBTRI && side == getFaceSideNumber()) {
MatrixDouble diff =
getCoordsAtGaussPts() - getMasterCoordsAtGaussPts();
constexpr double eps = 1e-12;
if (norm_inf(diff) > eps)
SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
"coordinates at integration pts are different");
const size_t nb_dofs = data.getN().size2() / 3;
const size_t nb_integration_pts = data.getN().size1();
if (data.getSpace() == H1) {
MatrixDouble *ptr_elem_data = nullptr;
ptr_elem_data = &elemData.shapeFunH1VolSide;
MatrixDouble &elem_data = *ptr_elem_data;
elem_data.resize(nb_integration_pts, nb_dofs, false);
for (size_t gg = 0; gg != nb_integration_pts; ++gg) {
FTensor::Tensor0<double *> t_base(&data.getN()(gg, 0));
for (size_t bb = 0; bb != nb_dofs; ++bb) {
elem_data(gg, bb) = t_base;
++t_base;
}
}
}
}
}
};
CommonData &elemData;
OnContactSideMaster(MoFEM::Interface &m_field, CommonData &elem_data)
FACEMASTER),
volSideFe(m_field), elemData(elem_data) {
volSideFe.getOpPtrVector().push_back(
new OnContactSideMaster::OpVolSide(elemData));
}
MoFEMErrorCode doWork(int side, EntityType type,
if (type == MBTRI && side == 0) {
const size_t nb_dofs = data.getN().size2() / 3;
const size_t nb_integration_pts = data.getN().size1();
if (data.getSpace() == H1) {
elemData.shapeFunH1Values.resize(nb_integration_pts, nb_dofs,
false);
elemData.shapeFunH1VolSide.resize(nb_integration_pts, nb_dofs,
false);
for (size_t gg = 0; gg != nb_integration_pts; ++gg) {
FTensor::Tensor0<double *> t_base(&data.getN()(gg, 0));
for (size_t bb = 0; bb != nb_dofs; ++bb) {
elemData.shapeFunH1Values(gg, bb) = t_base;
++t_base;
}
}
}
std::string side_fe_name = "V1";
const EntityHandle tri_master = getSideEntity(3, MBTRI);
CHKERR loopSideVolumes(side_fe_name, volSideFe, 3, tri_master);
auto check_continuity_of_base = [&](auto &vol_dot_data) {
if (vol_dot_data.size1() != elemData.dotNormalFace.size1())
SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
"Inconsistent number of integration points");
if (vol_dot_data.size2() != elemData.dotNormalFace.size2())
SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
"Inconsistent number of base functions");
constexpr double eps = 1e-12;
for (size_t gg = 0; gg != vol_dot_data.size1(); ++gg)
for (size_t bb = 0; bb != vol_dot_data.size2(); ++bb) {
const double error = std::abs(vol_dot_data(gg, bb) -
elemData.dotNormalFace(gg, bb));
if (error > eps)
SETERRQ2(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
"Inconsistency %3.4e != %3.4e", vol_dot_data(gg, bb),
elemData.dotNormalFace(gg, bb));
}
};
auto check_continuity_of_h1_base = [&](auto &vol_data) {
if (vol_data.size1() != elemData.shapeFunH1Values.size1())
SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
"Inconsistent number of integration points");
if (vol_data.size2() != elemData.shapeFunH1Values.size2())
SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
"Inconsistent number of base functions");
constexpr double eps = 1e-12;
for (size_t gg = 0; gg != vol_data.size1(); ++gg)
for (size_t bb = 0; bb != vol_data.size2(); ++bb) {
const double error = std::abs(
vol_data(gg, bb) - elemData.shapeFunH1Values(gg, bb));
if (error > eps)
SETERRQ2(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
"Inconsistency %3.4e != %3.4e", vol_data(gg, bb),
elemData.shapeFunH1Values(gg, bb));
}
};
if (elemData.dotNormalEleLeft.size2() != 0)
CHKERR check_continuity_of_base(elemData.dotNormalEleLeft);
else if (elemData.dotNormalEleRight.size2() != 0)
CHKERR check_continuity_of_base(elemData.dotNormalEleRight);
else if (elemData.shapeFunH1VolSide.size2() != 0)
CHKERR check_continuity_of_h1_base(elemData.shapeFunH1VolSide);
}
}
};
struct OnContactSideSlave
CommonData &elemData;
OpVolSide(CommonData &elem_data)
elemData(elem_data) {}
MoFEMErrorCode doWork(int side, EntityType type,
if (type == MBTRI && side == getFaceSideNumber()) {
MatrixDouble diff =
getCoordsAtGaussPts() - getMasterCoordsAtGaussPts();
constexpr double eps = 1e-12;
if (norm_inf(diff) > eps)
SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
"coordinates at integration pts are different");
const size_t nb_dofs = data.getN().size2() / 3;
const size_t nb_integration_pts = data.getN().size1();
if (data.getSpace() == H1) {
MatrixDouble *ptr_elem_data = nullptr;
ptr_elem_data = &elemData.shapeFunH1VolSide;
MatrixDouble &elem_data = *ptr_elem_data;
elem_data.resize(nb_integration_pts, nb_dofs, false);
for (size_t gg = 0; gg != nb_integration_pts; ++gg) {
FTensor::Tensor0<double *> t_base(&data.getN()(gg, 0));
for (size_t bb = 0; bb != nb_dofs; ++bb) {
elem_data(gg, bb) = t_base;
++t_base;
}
}
}
}
}
};
CommonData &elemData;
OnContactSideSlave(MoFEM::Interface &m_field, CommonData &elem_data)
FACESLAVE),
volSideFe(m_field), elemData(elem_data) {
volSideFe.getOpPtrVector().push_back(
new OnContactSideSlave::OpVolSide(elemData));
}
MoFEMErrorCode doWork(int side, EntityType type,
if (type == MBTRI && side == 0) {
const size_t nb_dofs = data.getN().size2() / 3;
const size_t nb_integration_pts = data.getN().size1();
if (data.getSpace() == H1) {
elemData.shapeFunH1Values.resize(nb_integration_pts, nb_dofs,
false);
elemData.shapeFunH1VolSide.resize(nb_integration_pts, nb_dofs,
false);
for (size_t gg = 0; gg != nb_integration_pts; ++gg) {
FTensor::Tensor0<double *> t_base(&data.getN()(gg, 0));
for (size_t bb = 0; bb != nb_dofs; ++bb) {
elemData.shapeFunH1Values(gg, bb) = t_base;
++t_base;
}
}
}
std::string side_fe_name = "V1";
const EntityHandle tri_slave = getSideEntity(4, MBTRI);
CHKERR loopSideVolumes(side_fe_name, volSideFe, 3, tri_slave);
auto check_continuity_of_base = [&](auto &vol_dot_data) {
if (vol_dot_data.size1() != elemData.dotNormalFace.size1())
SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
"Inconsistent number of integration points");
if (vol_dot_data.size2() != elemData.dotNormalFace.size2())
SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
"Inconsistent number of base functions");
constexpr double eps = 1e-12;
for (size_t gg = 0; gg != vol_dot_data.size1(); ++gg)
for (size_t bb = 0; bb != vol_dot_data.size2(); ++bb) {
const double error = std::abs(vol_dot_data(gg, bb) -
elemData.dotNormalFace(gg, bb));
if (error > eps)
SETERRQ2(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
"Inconsistency %3.4e != %3.4e", vol_dot_data(gg, bb),
elemData.dotNormalFace(gg, bb));
}
};
auto check_continuity_of_h1_base = [&](auto &vol_data) {
if (vol_data.size1() != elemData.shapeFunH1Values.size1())
SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
"Inconsistent number of integration points");
if (vol_data.size2() != elemData.shapeFunH1Values.size2())
SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
"Inconsistent number of base functions");
constexpr double eps = 1e-12;
for (size_t gg = 0; gg != vol_data.size1(); ++gg)
for (size_t bb = 0; bb != vol_data.size2(); ++bb) {
const double error = std::abs(
vol_data(gg, bb) - elemData.shapeFunH1Values(gg, bb));
if (error > eps)
SETERRQ2(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
"Inconsistency %3.4e != %3.4e", vol_data(gg, bb),
elemData.shapeFunH1Values(gg, bb));
}
};
if (elemData.dotNormalEleLeft.size2() != 0)
CHKERR check_continuity_of_base(elemData.dotNormalEleLeft);
else if (elemData.dotNormalEleRight.size2() != 0)
CHKERR check_continuity_of_base(elemData.dotNormalEleRight);
else if (elemData.shapeFunH1VolSide.size2() != 0)
CHKERR check_continuity_of_h1_base(elemData.shapeFunH1VolSide);
}
}
};
CommonData elem_data;
ContactPrismElementForcesAndSourcesCore contact_prism_fe_master(m_field);
ContactPrismElementForcesAndSourcesCore contact_prism_fe_slave(m_field);
// OnContactSideMaster
contact_prism_fe_master.getOpPtrVector().push_back(
new OnContactSideMaster(m_field, elem_data));
// OnContactSideSlave
contact_prism_fe_slave.getOpPtrVector().push_back(
new OnContactSideSlave(m_field, elem_data));
CHKERR m_field.loop_finite_elements("P1", "C2", contact_prism_fe_master);
CHKERR m_field.loop_finite_elements("P1", "C2", contact_prism_fe_slave);
}
return 0;
}
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
SIDESET
@ SIDESET
Definition: definitions.h:147
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::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
LASTBASE
@ LASTBASE
Definition: definitions.h:69
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
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
MoFEM::Types::MatrixDouble
UBlasMatrix< double > MatrixDouble
Definition: Types.hpp:77
MoFEM.hpp
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::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
MoFEM::Interface
DeprecatedCoreInterface Interface
Definition: Interface.hpp:1975
FieldSpace
FieldSpace
approximation spaces
Definition: definitions.h:82
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::ContactPrismElementForcesAndSourcesCore
ContactPrism finite element.
Definition: ContactPrismElementForcesAndSourcesCore.hpp:27
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
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::VolumeElementForcesAndSourcesCoreOnContactPrismSide
Base volume element used to integrate on contact surface (could be extended to other volume elements)
Definition: VolumeElementForcesAndSourcesCoreOnContactPrismSide.hpp:23
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.
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::VolumeElementForcesAndSourcesCore::UserDataOperator
friend class UserDataOperator
Definition: VolumeElementForcesAndSourcesCore.hpp:105
MoFEM::ContactPrismElementForcesAndSourcesCore::UserDataOperator
default operator for Contact Prism element
Definition: ContactPrismElementForcesAndSourcesCore.hpp:208
LASTSPACE
@ LASTSPACE
FieldSpace in [ 0, LASTSPACE )
Definition: definitions.h:89
MoFEM::CoreInterface::delete_ents_by_bit_ref
virtual MoFEMErrorCode delete_ents_by_bit_ref(const BitRefLevel bit, const BitRefLevel mask, const bool remove_parent=false, int verb=DEFAULT_VERBOSITY, MoFEMTypes mf=MF_ZERO)=0
delete entities form mofem and moab database
main
int main(int argc, char *argv[])
Definition: continuity_check_on_contact_prism_side_ele.cpp:16
MoFEM::ForcesAndSourcesCore::UserDataOperator
friend class UserDataOperator
Definition: ForcesAndSourcesCore.hpp:482
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
CATCH_ERRORS
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:372
DEMKOWICZ_JACOBI_BASE
@ DEMKOWICZ_JACOBI_BASE
Definition: definitions.h:66
FTensor::Tensor0
Definition: Tensor0.hpp:16
_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::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
MoFEM::Core
CoreTmp< 0 > Core
Definition: Core.hpp:1094
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
CommonData
Definition: continuity_check_on_skeleton_with_simple_2d_for_h1.cpp:22
AINSWORTH_LEGENDRE_BASE
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:60
MoFEM::EntitiesFieldData::EntData::getSpace
FieldSpace & getSpace()
Get field space.
Definition: EntitiesFieldData.hpp:1302
MoFEM::ForcesAndSourcesCore::getOpPtrVector
boost::ptr_deque< UserDataOperator > & getOpPtrVector()
Use to push back operator for row operator.
Definition: ForcesAndSourcesCore.hpp:83
HCURL
@ HCURL
field with continuous tangents
Definition: definitions.h:86
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
FieldApproximationBase
FieldApproximationBase
approximation base
Definition: definitions.h:58
_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
help
static char help[]
Definition: continuity_check_on_contact_prism_side_ele.cpp:14
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_ATOM_TEST_INVALID
@ MOFEM_ATOM_TEST_INVALID
Definition: definitions.h: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
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::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
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:346
MoFEM::ForcesAndSourcesCore::UserDataOperator::OPROW
@ OPROW
operator doWork function is executed on FE rows
Definition: ForcesAndSourcesCore.hpp:567