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.
static char help[] =
"...\n\n";
int main(
int argc,
char *argv[]) {
try {
moab::Core mb_instance;
moab::Interface &moab = mb_instance;
int rank;
MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
PetscBool flg = PETSC_TRUE;
#if PETSC_VERSION_GE(3, 6, 4)
255, &flg);
#else
CHKERR PetscOptionsGetString(PETSC_NULL, PETSC_NULL,
"-my_file",
#endif
if (flg != PETSC_TRUE) {
SETERRQ(PETSC_COMM_SELF, 1, "*** ERROR -my_file (MESH FILE NEEDED)");
}
const char *option;
option = "";
enum side_contact { MASTERSIDE, SLAVESIDE };
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;
}
};
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)
return space;
}
};
auto add_prism_interface = [&](
Range &tets,
Range &prisms,
std::vector<BitRefLevel> &bit_levels) {
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());
CHKERR moab.get_entities_by_type(cubit_meshset, MBTRI, tris,
true);
master_tris.merge(tris);
{
CHKERR moab.create_meshset(MESHSET_SET, ref_level_meshset);
->getEntitiesByTypeAndRefLevel(bit_levels.back(),
ref_level_meshset);
->getEntitiesByTypeAndRefLevel(bit_levels.back(),
ref_level_meshset);
CHKERR moab.get_entities_by_handle(ref_level_meshset,
ref_level_tets, true);
0);
cubit_meshset, true, true, 0);
CHKERR moab.delete_entities(&ref_level_meshset, 1);
}
->updateMeshsetByEntitiesChildren(
cubit_meshset, bit_levels.back(), cubit_meshset, MBMAXTYPE,
true);
}
}
}
for (unsigned int ll = 0; ll != bit_levels.size() - 1; ll++) {
true);
}
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);
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;
std::vector<BitRefLevel> bit_levels;
0, 3, bit_levels[0]);
CHKERR add_prism_interface(all_tets, contact_prisms, master_tris,
slave_tris, meshset_tets, meshset_prisms,
bit_levels);
}
};
struct OnContactSideMaster
"F2", UserDataOperator::OPROW),
elemData(elem_data) {}
if (type == MBTRI && side == getFaceSideNumber()) {
getCoordsAtGaussPts() - getMasterCoordsAtGaussPts();
constexpr double eps = 1e-12;
if (norm_inf(diff) >
eps)
"coordinates at integration pts are different");
const size_t nb_dofs = data.
getN().size2() / 3;
const size_t nb_integration_pts = data.
getN().size1();
ptr_elem_data = &elemData.shapeFunH1VolSide;
elem_data.resize(nb_integration_pts, nb_dofs, false);
for (size_t gg = 0; gg != nb_integration_pts; ++gg) {
for (size_t bb = 0; bb != nb_dofs; ++bb) {
elem_data(gg, bb) = t_base;
++t_base;
}
}
}
}
}
};
"F2", UserDataOperator::OPROW,
FACEMASTER),
volSideFe(m_field), elemData(elem_data) {
new OnContactSideMaster::OpVolSide(elemData));
}
if (type == MBTRI && side == 0) {
const size_t nb_dofs = data.
getN().size2() / 3;
const size_t nb_integration_pts = data.
getN().size1();
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) {
for (size_t bb = 0; bb != nb_dofs; ++bb) {
elemData.shapeFunH1Values(gg, bb) = t_base;
++t_base;
}
}
}
std::string side_fe_name = "V1";
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())
"Inconsistent number of integration points");
if (vol_dot_data.size2() != elemData.dotNormalFace.size2())
"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));
"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())
"Inconsistent number of integration points");
if (vol_data.size2() != elemData.shapeFunH1Values.size2())
"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));
"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
"F2", UserDataOperator::OPROW),
elemData(elem_data) {}
if (type == MBTRI && side == getFaceSideNumber()) {
getCoordsAtGaussPts() - getMasterCoordsAtGaussPts();
constexpr double eps = 1e-12;
if (norm_inf(diff) >
eps)
"coordinates at integration pts are different");
const size_t nb_dofs = data.
getN().size2() / 3;
const size_t nb_integration_pts = data.
getN().size1();
ptr_elem_data = &elemData.shapeFunH1VolSide;
elem_data.resize(nb_integration_pts, nb_dofs, false);
for (size_t gg = 0; gg != nb_integration_pts; ++gg) {
for (size_t bb = 0; bb != nb_dofs; ++bb) {
elem_data(gg, bb) = t_base;
++t_base;
}
}
}
}
}
};
"F2", UserDataOperator::OPROW,
FACESLAVE),
volSideFe(m_field), elemData(elem_data) {
new OnContactSideSlave::OpVolSide(elemData));
}
if (type == MBTRI && side == 0) {
const size_t nb_dofs = data.
getN().size2() / 3;
const size_t nb_integration_pts = data.
getN().size1();
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) {
for (size_t bb = 0; bb != nb_dofs; ++bb) {
elemData.shapeFunH1Values(gg, bb) = t_base;
++t_base;
}
}
}
std::string side_fe_name = "V1";
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())
"Inconsistent number of integration points");
if (vol_dot_data.size2() != elemData.dotNormalFace.size2())
"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));
"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())
"Inconsistent number of integration points");
if (vol_data.size2() != elemData.shapeFunH1Values.size2())
"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));
"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);
}
}
};
new OnContactSideMaster(m_field, elem_data));
new OnContactSideSlave(m_field, elem_data));
}
return 0;
}
ForcesAndSourcesCore::UserDataOperator UserDataOperator
#define CATCH_ERRORS
Catch errors.
FieldApproximationBase
approximation base
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
FieldSpace
approximation spaces
@ LASTSPACE
FieldSpace in [ 0, LASTSPACE )
@ HCURL
field with continuous tangents
@ HDIV
field with continuous normal traction
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ MOFEM_ATOM_TEST_INVALID
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
virtual MoFEMErrorCode add_finite_element(const std::string &fe_name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
add finite element
virtual MoFEMErrorCode build_finite_elements(int verb=DEFAULT_VERBOSITY)=0
Build finite elements.
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
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
virtual MoFEMErrorCode modify_finite_element_add_field_data(const std::string &fe_name, const std::string name_filed)=0
set finite element field data
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
virtual MoFEMErrorCode build_fields(int verb=DEFAULT_VERBOSITY)=0
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.
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.
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.
#define _IT_CUBITMESHSETS_FOR_LOOP_(MESHSET_MANAGER, IT)
Iterator that loops over all the Cubit MeshSets in a moFEM field.
#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.
MoFEMErrorCode buildProblem(const std::string name, const bool square_matrix, int verb=VERBOSE)
build problem data structures
virtual MoFEMErrorCode add_problem(const std::string &name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
Add problem.
virtual MoFEMErrorCode modify_problem_ref_level_add_bit(const std::string &name_problem, const BitRefLevel &bit)=0
add ref level to problem
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
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
UBlasMatrix< double > MatrixDouble
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
implementation of Data Operators for Forces and Sources
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
virtual MoFEMErrorCode build_adjacencies(const Range &ents, int verb=DEFAULT_VERBOSITY)=0
build adjacencies
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.
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Deprecated interface functions.
Data on single entity (This is passed as argument to DataOperator::doWork)
MatrixDouble & getN(const FieldApproximationBase base)
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....
FieldSpace & getSpace()
Get field space.
boost::ptr_deque< UserDataOperator > & getOpPtrVector()
Use to push back operator for row operator.
Create interface from given surface and insert flat prisms in-between.
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...
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...
Problem manager is used to build and partition problems.
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.