v0.10.0
continuity_check_on_skeleton_with_simple_2d_for_h1.cpp

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

/**
* \file continuity_check_on_skeleton_with_simple_2d_for_h1.cpp
* \ingroup mofem_simple_interface
\example continuity_check_on_skeleton_with_simple_2d_for_h1.cpp
*
* \brief Integration on skeleton for 2d
*
* Teting integration on skeleton and checking of continuity of hcurl space on
* edges.
*/
/* 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";
FaceElementForcesAndSourcesCore::NO_HO_GEOMETRY>;
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) {
MOFEM_LOG("ATOM", Sev::error)
<< "Quad coords: " << getCoordsAtGaussPts();
MOFEM_LOG("ATOM", Sev::error)
<< "Edge coords: " << getEdgeCoordsAtGaussPts();
MOFEM_LOG("ATOM", Sev::error) << "Difference: " << diff;
SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
"Coordinates at integration pts are different");
}
const size_t nb_dofs = data.getN().size2();
const size_t nb_integration_pts = data.getN().size1();
auto t_base = data.getFTensor0N();
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_base;
++t_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();
const size_t nb_integration_pts = data.getN().size1();
MOFEM_LOG("ATOM", Sev::noisy)
<< "Cords at integration points" << getCoordsAtGaussPts();
auto t_base = data.getFTensor0N();
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_base;
++t_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", H1, base, 1);
// CHKERR simple_interface->addDomainField("TEST_FIELD", L2,
// AINSWORTH_LEGENDRE_BASE, 1);
CHKERR simple_interface->addSkeletonField("FIELD", H1, base, 1);
// set fields order
CHKERR simple_interface->setFieldOrder("FIELD", 2);
// CHKERR simple_interface->setFieldOrder("TEST_FIELD", 1);
// 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));
skeleton_fe->getOpPtrVector().push_back(
new SkeletonFE(m_field, elem_data));
// iterate skeleton finite elements
CHKERR DMoFEMLoopFiniteElements(dm, simple_interface->getSkeletonFEName(),
skeleton_fe);
}
}
// finish work cleaning memory, getting statistics, etc.
return 0;
}
int main(int argc, char *argv[])
MoFEM::EdgeElementForcesAndSourcesCoreSwitch< EdgeElementForcesAndSourcesCore::NO_HO_GEOMETRY > EdgeEle
#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
@ H1
continuous field
Definition: definitions.h:177
#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
PetscErrorCode DMoFEMLoopFiniteElements(DM dm, const char fe_name[], MoFEM::FEMethod *method)
Executes FEMethod for finite elements in DM.
Definition: DMMMoFEM.cpp:507
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
Definition: DMMMoFEM.cpp:48
static LoggerType & setLog(const std::string channel)
Set ans resset chanel logger.
Definition: LogManager.cpp:368
#define MOFEM_LOG(channel, severity)
Log.
Definition: LogManager.hpp:312
#define MOFEM_LOG_ATTRIBUTES(channel, bit)
Add attributes to channel.
Definition: LogManager.hpp:300
#define MOFEM_LOG_FUNCTION()
Set scope.
Definition: LogManager.hpp:329
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:67
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)
DeprecatedCoreInterface Interface
Definition: Interface.hpp:1943
DataForcesAndSourcesCore::EntData EntData
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.
boost::ptr_vector< UserDataOperator > & getOpPtrVector()
Use to push back operator for row operator.
static boost::shared_ptr< SinkType > createSink(boost::shared_ptr< std::ostream > stream_ptr, std::string comm_filter)
Create a sink object.
Definition: LogManager.cpp:283
static boost::shared_ptr< std::ostream > getStrmSelf()
Get the strm self object.
Definition: LogManager.cpp:323
MoFEMErrorCode getInterface(const MOFEMuuid &uuid, IFACE *&iface) const
Get interface by uuid and return reference to pointer of interface.