v0.10.0
continuity_check_on_skeleton_3d.cpp

Testing integration on skeleton for 3D Checking continuity of hdiv and hcurl spaces on faces, and testing methods for integration on the skeleton.

/** \file continuity_check_on_skeleton_3d.cpp
* \example continuity_check_on_skeleton_3d.cpp
*
* \brief Testing integration on skeleton for 3D
*
* Checking continuity of hdiv and hcurl spaces on faces, and testing methods
* for integration on the skeleton.
*
*/
/* This file is part of MoFEM.
* MoFEM is free software: you can redistribute it and/or modify it under
* the terms of the GNU Lesser General Public License as published by the
* Free Software Foundation, either version 3 of the License, or (at your
* option) any later version.
*
* MoFEM is distributed in the hope that it will be useful, but WITHOUT
* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
* FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public
* License for more details.
*
* You should have received a copy of the GNU Lesser General Public
* License along with MoFEM. If not, see <http://www.gnu.org/licenses/>. */
#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)");
}
ParallelComm *pcomm = ParallelComm::get_pcomm(&moab, MYPCOMM_INDEX);
if (pcomm == NULL)
pcomm = new ParallelComm(&moab, PETSC_COMM_WORLD);
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;
// set entireties bit level
BitRefLevel bit_level0 = BitRefLevel().set(0);
CHKERR m_field.getInterface<BitRefManager>()->setBitRefLevelByDim(
0, 3, bit_level0);
// 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 { HDIV, HCURL, LASTSPACEOP };
const char *list_spaces[] = {"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)
return space;
}
return LASTSPACE;
};
CHKERR m_field.add_field("F2", get_space(), get_base(), 1);
// 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");
// 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);
CHKERR m_field.build_fields();
// add elements
CHKERR m_field.add_finite_element("V1");
CHKERR m_field.add_finite_element("S2");
CHKERR m_field.add_ents_to_finite_element_by_type(root_set, MBTET, "V1");
Range faces;
CHKERR m_field.getInterface<BitRefManager>()->getEntitiesByTypeAndRefLevel(
bit_level0, BitRefLevel().set(), MBTRI, faces);
CHKERR m_field.add_ents_to_finite_element_by_type(faces, MBTRI, "S2");
CHKERR m_field.build_adjacencies(bit_level0);
// Problems
CHKERR m_field.add_problem("P1");
// set refinement level for problem
CHKERR m_field.modify_problem_ref_level_add_bit("P1", bit_level0);
// build problems
ProblemsManager *prb_mng_ptr;
CHKERR m_field.getInterface(prb_mng_ptr);
CHKERR prb_mng_ptr->buildProblem("P1", true);
CHKERR prb_mng_ptr->partitionProblem("P1");
CHKERR prb_mng_ptr->partitionFiniteElements("P1");
CHKERR prb_mng_ptr->partitionGhostDofs("P1");
struct CommonData {
MatrixDouble dotNormalFace;
MatrixDouble dotNormalEleLeft;
MatrixDouble dotNormalEleRight;
};
struct SkeletonFE
struct OpVolSide
CommonData &elemData;
OpVolSide(CommonData &elem_data)
"F2", UserDataOperator::OPROW),
elemData(elem_data) {}
MoFEMErrorCode doWork(int side, EntityType type,
if (type == MBTRI && side == getFaceSideNumber()) {
MatrixDouble diff =
getCoordsAtGaussPts() - getFaceCoordsAtGaussPts();
const 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();
auto t_to_do_dot = getFTensor1Normal();
if (data.getSpace() == HCURL) {
auto s1 = getFTensor1Tangent1();
auto s2 = getFTensor1Tangent2();
t_to_do_dot(i) = s1(i) + s2(i);
}
auto t_hdiv_base = data.getFTensor1N<3>();
MatrixDouble *ptr_dot_elem_data = nullptr;
if (getFEMethod()->nInTheLoop == 0)
ptr_dot_elem_data = &elemData.dotNormalEleLeft;
else
ptr_dot_elem_data = &elemData.dotNormalEleRight;
MatrixDouble &dot_elem_data = *ptr_dot_elem_data;
dot_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) {
dot_elem_data(gg, bb) = t_to_do_dot(i) * t_hdiv_base(i);
++t_hdiv_base;
}
}
}
}
};
CommonData &elemData;
SkeletonFE(MoFEM::Interface &m_field, CommonData &elem_data)
"F2", UserDataOperator::OPROW),
volSideFe(m_field), elemData(elem_data) {
volSideFe.getOpPtrVector().push_back(
new SkeletonFE::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();
auto t_to_do_dot = getFTensor1Normal();
if (data.getSpace() == HCURL) {
auto s1 = getFTensor1Tangent1();
auto s2 = getFTensor1Tangent2();
t_to_do_dot(i) = s1(i) + s2(i);
}
auto t_hdiv_base = data.getFTensor1N<3>();
elemData.dotNormalFace.resize(nb_integration_pts, nb_dofs, false);
elemData.dotNormalEleLeft.resize(nb_integration_pts, 0, false);
elemData.dotNormalEleRight.resize(nb_integration_pts, 0, false);
for (size_t gg = 0; gg != nb_integration_pts; ++gg) {
for (size_t bb = 0; bb != nb_dofs; ++bb) {
elemData.dotNormalFace(gg, bb) = t_to_do_dot(i) * t_hdiv_base(i);
++t_hdiv_base;
}
}
CHKERR loopSideVolumes("V1", volSideFe);
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");
const 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));
}
};
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);
}
}
};
CommonData elem_data;
face_fe.getOpPtrVector().push_back(new SkeletonFE(m_field, elem_data));
CHKERR m_field.loop_finite_elements("P1", "S2", face_fe);
}
return 0;
}
ForcesAndSourcesCore::UserDataOperator UserDataOperator
int main(int argc, char *argv[])
static char help[]
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:441
FieldApproximationBase
approximation base
Definition: definitions.h:150
@ LASTBASE
Definition: definitions.h:161
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:152
@ DEMKOWICZ_JACOBI_BASE
Definition: definitions.h:158
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:516
FieldSpace
approximation spaces
Definition: definitions.h:174
@ LASTSPACE
FieldSpace in [ 0, LASTSPACE )
Definition: definitions.h:181
@ HCURL
field with continuous tangents
Definition: definitions.h:178
@ HDIV
field with continuous normal traction
Definition: definitions.h:179
#define MYPCOMM_INDEX
default communicator number PCOMM
Definition: definitions.h:292
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:415
@ MOFEM_ATOM_TEST_INVALID
Definition: definitions.h:132
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:485
#define CHKERR
Inline error check.
Definition: definitions.h:604
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:509
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 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 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_col(const std::string &fe_name, const std::string &name_row)=0
set field col 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.
VolumeElementForcesAndSourcesCoreOnSideSwitch< 0 > VolumeElementForcesAndSourcesCoreOnSide
Volume element used to integrate on skeleton.
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
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
FTensor::Index< 'i', 3 > i
char mesh_file_name[255]
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:67
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition: Types.hpp:51
ublas::matrix< double, ublas::row_major, DoubleAllocator > MatrixDouble
Definition: Types.hpp:76
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:21
CoreTmp< 0 > Core
Definition: Core.hpp:1129
PetscErrorCode PetscOptionsGetEList(PetscOptions *, const char pre[], const char name[], const char *const *list, PetscInt next, PetscInt *value, PetscBool *set)
PetscErrorCode PetscOptionsGetString(PetscOptions *, const char pre[], const char name[], char str[], size_t size, PetscBool *set)
DeprecatedCoreInterface Interface
Definition: Interface.hpp:1943
DataForcesAndSourcesCore::EntData EntData
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.
Core (interface) class.
Definition: Core.hpp:77
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
Definition: Core.cpp:60
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition: Core.cpp:100
Deprecated interface functions.
DEPRECATED MoFEMErrorCode loop_finite_elements(const Problem *problem_ptr, const std::string &fe_name, FEMethod &method, int lower_rank, int upper_rank, MoFEMTypes bh, int verb=DEFAULT_VERBOSITY)
MoFEMErrorCode getInterface(const MOFEMuuid &uuid, IFACE *&iface) const
Get interface by uuid and return reference to pointer of interface.
VolumeElementForcesAndSourcesCoreOnSideBase::UserDataOperator UserDataOperator