v0.14.0
hcurl_divergence_operator_2d.cpp

Testing Hcurl base, transfromed to Hdiv base in 2d using Green theorem.

Note 0: This is low-level implementation. Note 1: Generic implementation for Quad/Tri mesh of arbitrary order.

/**
* \file hcurl_divergence_operator_2d.cpp
* \example hcurl_divergence_operator_2d.cpp
*
* Testing Hcurl base, transfromed to Hdiv base in 2d using Green theorem.
*
* Note 0: This is low-level implementation.
* Note 1: Generic implementation for Quad/Tri mesh of arbitrary order.
*
*/
#include <MoFEM.hpp>
using namespace MoFEM;
static char help[] = "...\n\n";
constexpr int SPACE_DIM = 2;
struct OpDivergence : public FaceEleOp {
double &dIv;
OpDivergence(double &div) : FaceEleOp("FIELD1", OPROW), dIv(div) {}
MoFEMErrorCode doWork(int side, EntityType type,
const int nb_dofs = data.getIndices().size();
if (nb_dofs == 0)
const int nb_gauss_pts = data.getN().size1();
auto t_w = getFTensor0IntegrationWeight();
auto t_diff_base_fun = data.getFTensor2DiffN<3, 2>();
const auto area = getMeasure();
for (int gg = 0; gg != nb_gauss_pts; gg++) {
const double val = area * t_w;
for (int bb = 0; bb != nb_dofs; bb++) {
dIv += val * t_diff_base_fun(i, i);
++t_diff_base_fun;
}
++t_w;
}
}
};
struct OpFlux : public EdgeEleOp {
double &fLux;
OpFlux(double &flux) : EdgeEleOp("FIELD1", OPROW), fLux(flux) {}
MoFEMErrorCode doWork(int side, EntityType type,
const int nb_dofs = data.getIndices().size();
if (nb_dofs == 0)
const int nb_gauss_pts = data.getN().size1();
auto t_normal = getFTensor1Normal();
auto t_base_fun = data.getFTensor1N<3>();
auto t_w = getFTensor0IntegrationWeight();
for (int gg = 0; gg != nb_gauss_pts; gg++) {
for (int bb = 0; bb != nb_dofs; bb++) {
fLux += t_w * t_normal(i) * t_base_fun(i);
++t_base_fun;
}
++t_w;
}
}
};
int main(int argc, char *argv[]) {
MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
try {
moab::Core mb_instance;
moab::Interface &moab = mb_instance;
PetscBool flg_file = PETSC_TRUE;
char mesh_file_name[255];
CHKERR PetscOptionsGetString(PETSC_NULL, "", "-my_file", mesh_file_name,
255, &flg_file);
if (flg_file != PETSC_TRUE)
SETERRQ(PETSC_COMM_SELF, MOFEM_INVALID_DATA,
"*** ERROR -my_file (MESH FILE NEEDED)");
// Read mesh to MOAB
CHKERR moab.load_file(mesh_file_name, 0, "");
// Create MoFEM instance
MoFEM::Core core(moab);
MoFEM::Interface &m_field = core;
// set entities bit level
BitRefLevel bit_level0 = BitRefLevel().set(0);
CHKERR m_field.getInterface<BitRefManager>()->setBitRefLevelByDim(
0, 2, bit_level0);
// Declare elements
enum bases { AINSWORTH, DEMKOWICZ, LASBASETOP };
const char *list_bases[] = {"ainsworth", "demkowicz"};
PetscBool flg;
PetscInt choice_base_value = AINSWORTH;
CHKERR PetscOptionsGetEList(PETSC_NULL, NULL, "-base", list_bases,
LASBASETOP, &choice_base_value, &flg);
if (flg != PETSC_TRUE)
SETERRQ(PETSC_COMM_SELF, MOFEM_IMPOSSIBLE_CASE, "base not set");
if (choice_base_value == AINSWORTH)
else if (choice_base_value == DEMKOWICZ)
int order = 5;
CHKERR PetscOptionsGetInt(PETSC_NULL, "", "-order", &order, PETSC_NULL);
CHKERR m_field.add_field("FIELD1", HCURL, base, 1);
CHKERR m_field.add_finite_element("FACE_FE");
// Define rows/cols and element data
CHKERR m_field.modify_finite_element_add_field_row("FACE_FE", "FIELD1");
CHKERR m_field.modify_finite_element_add_field_col("FACE_FE", "FIELD1");
CHKERR m_field.modify_finite_element_add_field_data("FACE_FE", "FIELD1");
CHKERR m_field.add_finite_element("EDGE_FE");
// Define rows/cols and element data
CHKERR m_field.modify_finite_element_add_field_row("EDGE_FE", "FIELD1");
CHKERR m_field.modify_finite_element_add_field_col("EDGE_FE", "FIELD1");
CHKERR m_field.modify_finite_element_add_field_data("EDGE_FE", "FIELD1");
// Problem
CHKERR m_field.add_problem("TEST_PROBLEM");
// set finite elements for problem
CHKERR m_field.modify_problem_add_finite_element("TEST_PROBLEM", "FACE_FE");
CHKERR m_field.modify_problem_add_finite_element("TEST_PROBLEM", "EDGE_FE");
// set refinement level for problem
CHKERR m_field.modify_problem_ref_level_add_bit("TEST_PROBLEM", bit_level0);
// Add entities
CHKERR m_field.add_ents_to_field_by_dim(0, SPACE_DIM, "FIELD1");
// Set order
for (auto t : {MBEDGE, MBTRI, MBQUAD})
CHKERR m_field.set_field_order(0, t, "FIELD1", order);
// Add entities to elements
auto set_edge_elements_entities_on_mesh_skin = [&]() {
Range faces;
CHKERR moab.get_entities_by_dimension(0, 2, faces, false);
Skinner skin(&m_field.get_moab());
Range faces_skin;
CHKERR skin.find_skin(0, faces, false, faces_skin);
faces_skin, SPACE_DIM - 1, "EDGE_FE");
};
CHKERR set_edge_elements_entities_on_mesh_skin();
// Build database
CHKERR m_field.build_fields();
// build finite elemnts
// build adjacencies
CHKERR m_field.build_adjacencies(bit_level0);
// build problem
ProblemsManager *prb_mng_ptr;
CHKERR m_field.getInterface(prb_mng_ptr);
CHKERR prb_mng_ptr->buildProblem("TEST_PROBLEM", true);
// Partition
CHKERR prb_mng_ptr->partitionSimpleProblem("TEST_PROBLEM");
CHKERR prb_mng_ptr->partitionFiniteElements("TEST_PROBLEM");
CHKERR prb_mng_ptr->partitionGhostDofs("TEST_PROBLEM");
// integration rule
auto rule = [&](int, int, int p) { return 2 * p; };
auto calculate_divergence = [&]() {
double div = 0;
FaceEle fe_face(m_field);
fe_face.getRuleHook = rule;
CHKERR AddHOOps<2, 2, 2>::add(fe_face.getOpPtrVector(), {HDIV});
fe_face.getOpPtrVector().push_back(new OpDivergence(div));
CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "FACE_FE", fe_face);
return div;
};
auto calculate_flux = [&]() {
double flux = 0;
EdgeEle fe_edge(m_field);
fe_edge.getRuleHook = rule;
CHKERR AddHOOps<1, 2, 2>::add(fe_edge.getOpPtrVector(), {HDIV});
fe_edge.getOpPtrVector().push_back(new OpFlux(flux));
CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "EDGE_FE", fe_edge);
return flux;
};
const double div = calculate_divergence();
const double flux = calculate_flux();
MOFEM_LOG_CHANNEL("WOLD"); // reset channel
MOFEM_LOG_C("WORLD", Sev::inform,
"Div = %4.3e Flux = %3.4e Error = %4.3e\n", div, flux,
div - flux);
constexpr double tol = 1e-8;
if (std::abs(div - flux) > tol)
SETERRQ2(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
"Test failed (div != flux) %3.4e != %3.4e", div, flux);
}
}
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:460
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference 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:128
MoFEM::CoreInterface::loop_finite_elements
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.
MoFEM::ProblemsManager::buildProblem
MoFEMErrorCode buildProblem(const std::string name, const bool square_matrix, int verb=VERBOSE)
build problem data structures
Definition: ProblemsManager.cpp:87
MoFEM::CoreTmp< 0 >
Core (interface) class.
Definition: Core.hpp:82
MOFEM_LOG_CHANNEL
#define MOFEM_LOG_CHANNEL(channel)
Set and reset channel.
Definition: LogManager.hpp:284
MoFEM::ProblemsManager
Problem manager is used to build and partition problems.
Definition: ProblemsManager.hpp:21
MoFEM::CoreInterface::modify_finite_element_add_field_row
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
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
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_IMPOSSIBLE_CASE
@ MOFEM_IMPOSSIBLE_CASE
Definition: definitions.h:35
MoFEM::CoreInterface::add_ents_to_field_by_dim
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.
order
constexpr int order
Definition: dg_projection.cpp:18
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
MoFEM::Interface
DeprecatedCoreInterface Interface
Definition: Interface.hpp:2002
MoFEM::EntitiesFieldData::EntData::getFTensor2DiffN
FTensor::Tensor2< FTensor::PackPtr< double *, Tensor_Dim0 *Tensor_Dim1 >, Tensor_Dim0, Tensor_Dim1 > getFTensor2DiffN(FieldApproximationBase base)
Get derivatives of base functions for Hdiv space.
Definition: EntitiesFieldData.cpp:660
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
MoFEM::FaceElementForcesAndSourcesCore::UserDataOperator
friend class UserDataOperator
Definition: FaceElementForcesAndSourcesCore.hpp:86
MoFEM::CoreInterface::add_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
help
static char help[]
Definition: hcurl_divergence_operator_2d.cpp:18
MoFEM::CoreInterface::get_moab
virtual moab::Interface & get_moab()=0
MoFEM
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
MoFEM::CoreInterface::modify_finite_element_add_field_col
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
MoFEM::CoreInterface::build_finite_elements
virtual MoFEMErrorCode build_finite_elements(int verb=DEFAULT_VERBOSITY)=0
Build finite elements.
MoFEM::FaceElementForcesAndSourcesCore::UserDataOperator
default operator for TRI element
Definition: FaceElementForcesAndSourcesCore.hpp:94
MOFEM_LOG_C
#define MOFEM_LOG_C(channel, severity, format,...)
Definition: LogManager.hpp:311
MoFEM::CoreInterface::add_field
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.
convert.type
type
Definition: convert.py:64
OpFlux
Definition: hcurl_divergence_operator_2d.cpp:56
MoFEM::EntitiesFieldData::EntData::getIndices
const VectorInt & getIndices() const
Get global indices of dofs on entity.
Definition: EntitiesFieldData.hpp:1204
SPACE_DIM
constexpr int SPACE_DIM
Definition: hcurl_divergence_operator_2d.cpp:26
MoFEM::FaceElementForcesAndSourcesCore
Face finite element.
Definition: FaceElementForcesAndSourcesCore.hpp:23
MoFEM::CoreInterface::add_ents_to_finite_element_by_dim
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
t
constexpr double t
plate stiffness
Definition: plate.cpp:59
MoFEM::ForcesAndSourcesCore::UserDataOperator
friend class UserDataOperator
Definition: ForcesAndSourcesCore.hpp:482
i
FTensor::Index< 'i', SPACE_DIM > i
Definition: hcurl_divergence_operator_2d.cpp:27
MoFEM::ProblemsManager::partitionFiniteElements
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
Definition: ProblemsManager.cpp:2167
FaceEleOp
FaceEle::UserDataOperator FaceEleOp
Definition: hcurl_divergence_operator_2d.cpp:23
MoFEM::CoreInterface::modify_problem_ref_level_add_bit
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', SPACE_DIM >
mesh_file_name
char mesh_file_name[255]
Definition: mesh_smoothing.cpp:23
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: hcurl_divergence_operator_2d.cpp:83
Range
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
CATCH_ERRORS
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:385
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:1308
MoFEM::Core
CoreTmp< 0 > Core
Definition: Core.hpp:1140
MoFEM::ProblemsManager::partitionSimpleProblem
MoFEMErrorCode partitionSimpleProblem(const std::string name, int verb=VERBOSE)
partition problem dofs
Definition: ProblemsManager.cpp:1537
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
AINSWORTH_LEGENDRE_BASE
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:60
HCURL
@ HCURL
field with continuous tangents
Definition: definitions.h:86
MoFEM::PetscOptionsGetString
PetscErrorCode PetscOptionsGetString(PetscOptions *, const char pre[], const char name[], char str[], size_t size, PetscBool *set)
Definition: DeprecatedPetsc.hpp:172
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::CoreInterface::build_fields
virtual MoFEMErrorCode build_fields(int verb=DEFAULT_VERBOSITY)=0
MoFEM::CoreInterface::modify_finite_element_add_field_data
virtual MoFEMErrorCode modify_finite_element_add_field_data(const std::string &fe_name, const std::string name_filed)=0
set finite element field data
FieldApproximationBase
FieldApproximationBase
approximation base
Definition: definitions.h:58
MoFEM::BitRefManager
Managing BitRefLevels.
Definition: BitRefManager.hpp:21
MoFEM::CoreInterface::modify_problem_add_finite_element
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
MoFEM::Types::BitRefLevel
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition: Types.hpp:40
MOFEM_ATOM_TEST_INVALID
@ MOFEM_ATOM_TEST_INVALID
Definition: definitions.h:40
MoFEM::CoreInterface::build_adjacencies
virtual MoFEMErrorCode build_adjacencies(const Range &ents, int verb=DEFAULT_VERBOSITY)=0
build adjacencies
MoFEM::CoreInterface::set_field_order
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.
MoFEM::ProblemsManager::partitionGhostDofs
MoFEMErrorCode partitionGhostDofs(const std::string name, int verb=VERBOSE)
determine ghost nodes
Definition: ProblemsManager.cpp:2339
MoFEM::CoreInterface::add_problem
virtual MoFEMErrorCode add_problem(const std::string &name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
Add problem.
convert.int
int
Definition: convert.py:64
MoFEM::PetscOptionsGetInt
PetscErrorCode PetscOptionsGetInt(PetscOptions *, const char pre[], const char name[], PetscInt *ivalue, PetscBool *set)
Definition: DeprecatedPetsc.hpp:142
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
MOFEM_INVALID_DATA
@ MOFEM_INVALID_DATA
Definition: definitions.h:36
OpDivergence
Definition: hcurl_divergence_operator_2d.cpp:29
tol
double tol
Definition: mesh_smoothing.cpp:27