v0.14.0
Loading...
Searching...
No Matches
continuity_check_on_skeleton_3d.cpp

Testing integration on skeleton for 3D.

Testing integration on skeleton for 3DChecking 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.
*
*/
#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);
auto core_log = logging::core::get();
core_log->add_sink(
LogManager::createSink(LogManager::getStrmSelf(), "ATOM_TEST"));
LogManager::setLog("ATOM_TEST");
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;
// 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_dim(root_set, 3, "F2");
// set app. order
int order = 3;
CHKERR m_field.set_field_order(root_set, MBTET, "F2", order);
CHKERR m_field.set_field_order(root_set, MBHEX, "F2", order);
CHKERR m_field.set_field_order(root_set, MBTRI, "F2", order);
CHKERR m_field.set_field_order(root_set, MBQUAD, "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_dim(root_set, 3, "V1");
Range faces;
CHKERR m_field.getInterface<BitRefManager>()->getEntitiesByDimAndRefLevel(
bit_level0, BitRefLevel().set(), 2, faces);
CHKERR m_field.add_ents_to_finite_element_by_dim(faces, 2, "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 {
VectorDouble dotNormalFace;
VectorDouble dotNormalEleLeft;
VectorDouble dotNormalEleRight;
};
struct SkeletonFE
: public FaceElementForcesAndSourcesCore::UserDataOperator {
VolumeElementForcesAndSourcesCoreOnSide volSideFe;
struct OpVolSide
: public VolumeElementForcesAndSourcesCoreOnSide::UserDataOperator {
CommonData &elemData;
OpVolSide(CommonData &elem_data)
: VolumeElementForcesAndSourcesCoreOnSide::UserDataOperator(
"F2", UserDataOperator::OPROW),
elemData(elem_data) {}
MoFEMErrorCode doWork(int side, EntityType type,
if (CN::Dimension(type) == 3) {
MatrixDouble diff =
getCoordsAtGaussPts() - getFaceCoordsAtGaussPts();
MOFEM_LOG("ATOM_TEST", Sev::noisy)
<< "getCoordsAtGaussPts() " << getCoordsAtGaussPts();
MOFEM_LOG("ATOM_TEST", Sev::noisy)
<< "getFaceCoordsAtGaussPts() " << getFaceCoordsAtGaussPts();
constexpr double eps = 1e-12;
if (norm_inf(diff) > eps) {
MOFEM_LOG("ATOM_TEST", Sev::error) << "diff " << diff;
SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
"Coordinates at integration pts are different");
}
}
const size_t nb_dofs = data.getFieldData().size();
if (nb_dofs) {
const size_t nb_integration_pts = getGaussPts().size2();
FTensor::Index<'i', 3> i;
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);
}
VectorDouble *ptr_dot_elem_data = nullptr;
if (getFEMethod()->nInTheLoop == 0)
ptr_dot_elem_data = &elemData.dotNormalEleLeft;
else
ptr_dot_elem_data = &elemData.dotNormalEleRight;
auto &dot_elem_data = *ptr_dot_elem_data;
if (dot_elem_data.size() == 0) {
dot_elem_data.resize(nb_integration_pts, false);
dot_elem_data.clear();
}
auto t_hdiv_base = data.getFTensor1N<3>();
for (size_t gg = 0; gg != nb_integration_pts; ++gg) {
auto t_data = data.getFTensor0FieldData();
for (size_t bb = 0; bb != nb_dofs; ++bb) {
dot_elem_data(gg) += (t_to_do_dot(i) * t_hdiv_base(i)) * t_data;
++t_hdiv_base;
++t_data;
}
}
}
}
};
SkeletonFE(MoFEM::Interface &m_field, CommonData &elem_data)
: FaceElementForcesAndSourcesCore::UserDataOperator(
"F2", UserDataOperator::OPROW),
volSideFe(m_field), elemData(elem_data) {
auto jac_ptr = boost::make_shared<MatrixDouble>();
auto inv_jac_ptr = boost::make_shared<MatrixDouble>();
auto det_ptr = boost::make_shared<VectorDouble>();
volSideFe.getOpPtrVector().push_back(new OpCalculateHOJac<3>(jac_ptr));
volSideFe.getOpPtrVector().push_back(
new OpInvertMatrix<3>(jac_ptr, det_ptr, inv_jac_ptr));
volSideFe.getOpPtrVector().push_back(new OpSetHOWeights(det_ptr));
volSideFe.getOpPtrVector().push_back(
volSideFe.getOpPtrVector().push_back(
new OpSetHOContravariantPiolaTransform(HDIV, det_ptr, jac_ptr));
volSideFe.getOpPtrVector().push_back(
new SkeletonFE::OpVolSide(elemData));
}
const size_t nb_integration_pts = getGaussPts().size2();
if (side == 0 && type == MBEDGE) {
elemData.dotNormalFace.resize(nb_integration_pts, false);
elemData.dotNormalFace.clear();
}
const size_t nb_dofs = data.getFieldData().size();
if (nb_dofs) {
FTensor::Index<'i', 3> i;
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>();
for (size_t gg = 0; gg != nb_integration_pts; ++gg) {
auto t_data = data.getFTensor0FieldData();
for (size_t bb = 0; bb != nb_dofs; ++bb) {
elemData.dotNormalFace(gg) +=
(t_to_do_dot(i) * t_hdiv_base(i)) * t_data;
++t_hdiv_base;
++t_data;
}
}
}
if (CN::Dimension(type) == 2) {
elemData.dotNormalEleLeft.resize(0, false);
elemData.dotNormalEleRight.resize(0, false);
CHKERR loopSideVolumes("V1", volSideFe);
auto check_continuity_of_base = [&](auto &vol_dot_data) {
if (vol_dot_data.size() != elemData.dotNormalFace.size())
SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
"Inconsistent number of integration points");
const double eps = 1e-12;
for (size_t gg = 0; gg != nb_integration_pts; ++gg) {
const double error =
std::abs(vol_dot_data(gg) - elemData.dotNormalFace(gg));
MOFEM_LOG("ATOM_TEST", Sev::noisy) << "Error: " << error;
if (error > eps)
SETERRQ2(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
"Inconsistency %3.4e != %3.4e", vol_dot_data(gg),
elemData.dotNormalFace(gg));
}
};
if (elemData.dotNormalEleLeft.size() != 0)
CHKERR check_continuity_of_base(elemData.dotNormalEleLeft);
else if (elemData.dotNormalEleRight.size() != 0)
CHKERR check_continuity_of_base(elemData.dotNormalEleRight);
}
}
};
CommonData elem_data;
face_fe.getRuleHook = [](int, int, int) { return 1; };
face_fe.getOpPtrVector().push_back(
face_fe.getOpPtrVector().push_back(
face_fe.getOpPtrVector().push_back(new SkeletonFE(m_field, elem_data));
auto field_ents_ptr = m_field.get_field_ents();
auto cache_ptr = boost::make_shared<CacheTuple>();
CHKERR m_field.cache_problem_entities("P1", cache_ptr);
for (auto &ent_ptr : (*field_ents_ptr)) {
MOFEM_LOG("ATOM_TEST", Sev::verbose) << *ent_ptr;
for(auto &v : ent_ptr->getEntFieldData())
v = 1;
CHKERR m_field.loop_finite_elements("P1", "S2", face_fe, 0, 1, nullptr,
MF_ZERO, cache_ptr);
for (auto &v : ent_ptr->getEntFieldData())
v = 0;
}
}
return 0;
}
static char help[]
int main()
static const double eps
#define CATCH_ERRORS
Catch errors.
@ MF_ZERO
Definition definitions.h:98
FieldApproximationBase
approximation base
Definition definitions.h:58
@ LASTBASE
Definition definitions.h:69
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition definitions.h:60
@ DEMKOWICZ_JACOBI_BASE
Definition definitions.h:66
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
FieldSpace
approximation spaces
Definition definitions.h:82
@ LASTSPACE
FieldSpace in [ 0, LASTSPACE )
Definition definitions.h:89
@ HCURL
field with continuous tangents
Definition definitions.h:86
@ HDIV
field with continuous normal traction
Definition definitions.h:87
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ MOFEM_ATOM_TEST_INVALID
Definition definitions.h:40
#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 ...
constexpr int order
virtual const FieldEntity_multiIndex * get_field_ents() const =0
Get the field ents object.
virtual MoFEMErrorCode add_ents_to_finite_element_by_dim(const EntityHandle entities, const int dim, const std::string &name, const bool recursive=true)=0
add entities to 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
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 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 add_ents_to_field_by_dim(const Range &ents, const int dim, const std::string &name, int verb=DEFAULT_VERBOSITY)=0
Add entities to field meshset.
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.
#define MOFEM_LOG(channel, severity)
Log.
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.
MoFEMErrorCode partitionGhostDofs(const std::string name, int verb=VERBOSE)
determine ghost nodes
MoFEMErrorCode buildProblem(const std::string name, const bool square_matrix, int verb=VERBOSE)
build problem data structures
MoFEMErrorCode partitionProblem(const std::string name, int verb=VERBOSE)
partition problem dofs (collective)
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
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
FTensor::Index< 'i', SPACE_DIM > i
const double v
phase velocity of light in medium (cm/ns)
char mesh_file_name[255]
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
implementation of Data Operators for Forces and Sources
Definition Common.hpp:10
Managing BitRefLevels.
virtual MoFEMErrorCode cache_problem_entities(const std::string prb_name, CacheTupleWeakPtr cache_ptr)=0
Cache variables.
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:82
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
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition Core.cpp:112
Deprecated interface functions.
auto getFTensor1Normal()
get edge normal NOTE: it should be used only in 2D analysis
Data on single entity (This is passed as argument to DataOperator::doWork)
const VectorDouble & getFieldData() const
get dofs values
FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > getFTensor0FieldData()
Resturn scalar files as a FTensor of rank 0.
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1N(FieldApproximationBase base)
Get base functions for Hdiv/Hcurl spaces.
FieldSpace & getSpace()
Get field space.
MatrixDouble & getGaussPts()
matrix of integration (Gauss) points for Volume Element
transform Hdiv base fluxes from reference element to physical triangle
transform Hcurl base fluxes from reference element to physical triangle
Apply contravariant (Piola) transfer to Hdiv space for HO geometry.
Apply covariant (Piola) transfer to Hcurl space for HO geometry.
Set inverse jacobian to base functions.
Problem manager is used to build and partition problems.
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
SkeletonFE(MoFEM::Interface &m_field, CommonData &elem_data)