v0.14.0
continuity_check_on_skeleton_with_simple_2d_for_hcurl.cpp

Integration on skeleton for 2dTesting integration on skeleton and checking of continuity of hcurl space on edges.

/**
* \file continuity_check_on_skeleton_with_simple_2d_for_hcurl.cpp
* \ingroup mofem_simple_interface
* \example continuity_check_on_skeleton_with_simple_2d_for_hcurl.cpp
*
* \brief Integration on skeleton for 2d
*
* Testing integration on skeleton and checking of continuity of hcurl space on
* edges.
*/
#include <MoFEM.hpp>
using namespace MoFEM;
static char help[] = "...\n\n";
struct CommonData {
MatrixDouble dotEdge;
MatrixDouble dotEleLeft;
MatrixDouble dotEleRight;
};
struct SkeletonFE : public EdgeEleOp {
struct OpFaceSide : public FaceEleOnSideOp {
CommonData &elemData;
OpFaceSide(CommonData &elem_data)
: FaceEleOnSideOp("FIELD", UserDataOperator::OPROW), elemData(elem_data) {}
MoFEMErrorCode doWork(int side, EntityType type,
if (type == MBEDGE && side == getEdgeSideNumber()) {
MatrixDouble diff = getCoordsAtGaussPts() - getEdgeCoordsAtGaussPts();
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_tangent = getFTensor1Direction();
auto t_hcurl_base = data.getFTensor1N<3>();
MatrixDouble *ptr_dot_elem_data = nullptr;
if (getFEMethod()->nInTheLoop == 0)
ptr_dot_elem_data = &elemData.dotEleLeft;
else
ptr_dot_elem_data = &elemData.dotEleRight;
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_tangent(i) * t_hcurl_base(i);
++t_hcurl_base;
}
}
}
}
};
CommonData &elemData;
FaceEleOnSide faceSideFe;
SkeletonFE(MoFEM::Interface &m_field, CommonData &elem_data)
faceSideFe(m_field), elemData(elem_data) {
faceSideFe.getOpPtrVector().push_back(new SkeletonFE::OpFaceSide(elemData));
}
MoFEMErrorCode doWork(int side, EntityType type,
if (type == MBEDGE) {
const size_t nb_dofs = data.getN().size2() / 3;
const size_t nb_integration_pts = data.getN().size1();
auto t_tangent = getFTensor1Direction();
auto t_hcurl_base = data.getFTensor1N<3>();
elemData.dotEdge.resize(nb_integration_pts, nb_dofs, false);
elemData.dotEleLeft.resize(nb_integration_pts, 0, false);
elemData.dotEleRight.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.dotEdge(gg, bb) = t_tangent(i) * t_hcurl_base(i);
++t_hcurl_base;
}
}
CHKERR loopSideFaces("dFE", faceSideFe);
auto check_continuity_of_base = [&](auto &vol_dot_data) {
if (vol_dot_data.size1() != elemData.dotEdge.size1())
SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
"Inconsistent number of integration points");
if (vol_dot_data.size2() != elemData.dotEdge.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.dotEdge(gg, bb));
if (error > eps)
SETERRQ4(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
"Inconsistency (%d, %d) %3.4e != %3.4e", gg, bb,
vol_dot_data(gg, bb), elemData.dotEdge(gg, bb));
else
MOFEM_LOG("ATOM", Sev::noisy) << "Ok";
}
};
if (elemData.dotEleLeft.size2() != 0)
CHKERR check_continuity_of_base(elemData.dotEleLeft);
else if (elemData.dotEleRight.size2() != 0)
CHKERR check_continuity_of_base(elemData.dotEleRight);
}
}
};
int main(int argc, char *argv[]) {
// initialize petsc
MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
try {
// Create MoAB database
moab::Core moab_core;
moab::Interface &moab = moab_core;
// Create MoFEM database and link it to MoAB
MoFEM::Core mofem_core(moab);
MoFEM::Interface &m_field = mofem_core;
auto core_log = logging::core::get();
core_log->add_sink(
// Register DM Manager
DMType dm_name = "DMMOFEM";
// Simple interface
Simple *simple_interface;
CHKERR m_field.getInterface(simple_interface);
{
// get options from command line
CHKERR simple_interface->getOptions();
// load mesh file
CHKERR simple_interface->loadFile("");
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;
};
// add fields
auto base = get_base();
CHKERR simple_interface->addDomainField("FIELD", HCURL, base, 1);
CHKERR simple_interface->addSkeletonField("FIELD", HCURL, base, 1);
// set fields order
CHKERR simple_interface->setFieldOrder("FIELD", 3);
// setup problem
CHKERR simple_interface->setUp();
// get dm
auto dm = simple_interface->getDM();
// create elements
CommonData elem_data;
boost::shared_ptr<EdgeEle> skeleton_fe =
boost::shared_ptr<EdgeEle>(new EdgeEle(m_field));
CHKERR AddHOOps<1, 2, 2>::add(skeleton_fe->getOpPtrVector(), {HCURL});
skeleton_fe->getOpPtrVector().push_back(
new SkeletonFE(m_field, elem_data));
// iterate skeleton finite elements
skeleton_fe);
}
}
// finish work cleaning memory, getting statistics, etc.
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::EdgeElementForcesAndSourcesCore
Edge finite element.
Definition: EdgeElementForcesAndSourcesCore.hpp:30
MoFEM::EntitiesFieldData::EntData
Data on single entity (This is passed as argument to DataOperator::doWork)
Definition: EntitiesFieldData.hpp:127
CommonData::dotEdge
MatrixDouble dotEdge
Definition: continuity_check_on_skeleton_with_simple_2d_for_h1.cpp:23
MoFEM::CoreTmp< 0 >
Core (interface) class.
Definition: Core.hpp:82
LASTBASE
@ LASTBASE
Definition: definitions.h:69
SkeletonFE::OpFaceSide
Definition: continuity_check_on_skeleton_with_simple_2d_for_h1.cpp:30
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
CommonData::dotEleLeft
MatrixDouble dotEleLeft
Definition: continuity_check_on_skeleton_with_simple_2d_for_h1.cpp:24
MoFEM::Types::MatrixDouble
UBlasMatrix< double > MatrixDouble
Definition: Types.hpp:77
MoFEM::Simple::getSkeletonFEName
const std::string getSkeletonFEName() const
Get the Skeleton FE Name.
Definition: Simple.hpp:355
MoFEM::Simple::loadFile
MoFEMErrorCode loadFile(const std::string options, const std::string mesh_file_name, LoadFileFunc loadFunc=defaultLoadFileFunc)
Load mesh file.
Definition: Simple.cpp:194
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::LogManager::createSink
static boost::shared_ptr< SinkType > createSink(boost::shared_ptr< std::ostream > stream_ptr, std::string comm_filter)
Create a sink object.
Definition: LogManager.cpp:298
MoFEM::Simple
Simple interface for fast problem set-up.
Definition: Simple.hpp:27
SkeletonFE
Definition: continuity_check_on_skeleton_with_simple_2d_for_h1.cpp:28
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
MoFEM::Simple::addSkeletonField
MoFEMErrorCode addSkeletonField(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_ZERO, int verb=-1)
Add field on skeleton.
Definition: Simple.cpp:300
MoFEM::Interface
DeprecatedCoreInterface Interface
Definition: Interface.hpp:1975
MoFEM::Simple::getOptions
MoFEMErrorCode getOptions()
get options
Definition: Simple.cpp:180
MoFEM::Simple::getDM
MoFEMErrorCode getDM(DM *dm)
Get DM.
Definition: Simple.cpp:667
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:535
MoFEM::FaceElementForcesAndSourcesCore::UserDataOperator
friend class UserDataOperator
Definition: FaceElementForcesAndSourcesCore.hpp:86
MoFEM
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
MoFEM::ForcesAndSourcesCore::UserDataOperator
Definition: ForcesAndSourcesCore.hpp:549
convert.type
type
Definition: convert.py:64
CommonData::dotEleRight
MatrixDouble dotEleRight
Definition: continuity_check_on_skeleton_with_simple_2d_for_h1.cpp:25
help
static char help[]
Definition: continuity_check_on_skeleton_with_simple_2d_for_hcurl.cpp:15
MoFEM::Simple::addDomainField
MoFEMErrorCode addDomainField(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_ZERO, int verb=-1)
Add field on domain.
Definition: Simple.cpp:264
MoFEM::DMRegister_MoFEM
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
Definition: DMMoFEM.cpp:47
MoFEM::LogManager::getStrmSelf
static boost::shared_ptr< std::ostream > getStrmSelf()
Get the strm self object.
Definition: LogManager.cpp:340
EdgeEle
MoFEM::EdgeElementForcesAndSourcesCore EdgeEle
Definition: continuity_check_on_skeleton_with_simple_2d_for_hcurl.cpp:18
MoFEM::ForcesAndSourcesCore::UserDataOperator
friend class UserDataOperator
Definition: ForcesAndSourcesCore.hpp:482
i
FTensor::Index< 'i', SPACE_DIM > i
Definition: hcurl_divergence_operator_2d.cpp:27
FTensor::Index< 'i', 3 >
MoFEM::Simple::setFieldOrder
MoFEMErrorCode setFieldOrder(const std::string field_name, const int order, const Range *ents=NULL)
Set field order.
Definition: Simple.cpp:473
MoFEM::AddHOOps
Add operators pushing bases from local to physical configuration.
Definition: HODataOperators.hpp:503
MoFEM::EdgeElementForcesAndSourcesCore::UserDataOperator
default operator for EDGE element
Definition: EdgeElementForcesAndSourcesCore.hpp:68
main
int main(int argc, char *argv[])
Definition: continuity_check_on_skeleton_with_simple_2d_for_hcurl.cpp:145
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_LOG
#define MOFEM_LOG(channel, severity)
Log.
Definition: LogManager.hpp:308
CATCH_ERRORS
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:372
DEMKOWICZ_JACOBI_BASE
@ DEMKOWICZ_JACOBI_BASE
Definition: definitions.h:66
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
OpFaceSide
Definition: hello_world.cpp:94
eps
static const double eps
Definition: check_base_functions_derivatives_on_tet.cpp:11
MoFEM::EntitiesFieldData::EntData::getFTensor1N
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1N(FieldApproximationBase base)
Get base functions for Hdiv/Hcurl spaces.
Definition: EntitiesFieldData.cpp:640
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::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::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::FaceElementForcesAndSourcesCoreOnSide::UserDataOperator
default operator for Face element
Definition: FaceElementForcesAndSourcesCoreOnSide.hpp:92
FieldApproximationBase
FieldApproximationBase
approximation base
Definition: definitions.h:58
MOFEM_ATOM_TEST_INVALID
@ MOFEM_ATOM_TEST_INVALID
Definition: definitions.h:40
MoFEMFunctionBeginHot
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:440
MoFEM::LogManager::setLog
static LoggerType & setLog(const std::string channel)
Set ans resset chanel logger.
Definition: LogManager.cpp:389
MoFEM::Simple::setUp
MoFEMErrorCode setUp(const PetscBool is_partitioned=PETSC_TRUE)
Setup problem.
Definition: Simple.cpp:611
MoFEM::DMoFEMLoopFiniteElements
PetscErrorCode DMoFEMLoopFiniteElements(DM dm, const char fe_name[], MoFEM::FEMethod *method, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr())
Executes FEMethod for finite elements in DM.
Definition: DMMoFEM.cpp:590
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
MoFEM::FaceElementForcesAndSourcesCoreOnSide
Base face element used to integrate on skeleton.
Definition: FaceElementForcesAndSourcesCoreOnSide.hpp:19
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346