v0.10.0
hcurl_divergence_operator_2d.cpp

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

/**
* \file hcurl_divergence_operator_2d.cpp
* \example hcurl_divergence_operator_2d.cpp
*
* Testing Hcurl base, transfromed to Hdiv base in 2d using Green theorem
*
*/
/* 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 |
FaceElementForcesAndSourcesCore::NO_CONTRAVARIANT_TRANSFORM_HDIV |
FaceElementForcesAndSourcesCore::NO_COVARIANT_TRANSFORM_HCURL>;
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_diff_base_fun = data.getFTensor2DiffN<3, 2>();
for (int gg = 0; gg != nb_gauss_pts; gg++) {
const double val = getArea() * getGaussPts()(2, gg);
for (int bb = 0; bb != nb_dofs; bb++) {
dIv += val * (t_diff_base_fun(0, 0) + t_diff_base_fun(1, 1));
++t_diff_base_fun;
}
}
}
};
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();
FTensor::Tensor1<double, 3> t_normal(-getDirection()[1], getDirection()[0],
0.);
auto t_base_fun = data.getFTensor1N<3>();
for (int gg = 0; gg != nb_gauss_pts; gg++) {
const double val = getGaussPts()(1, gg);
for (int bb = 0; bb != nb_dofs; bb++) {
fLux += val * t_normal(i) * t_base_fun(i);
++t_base_fun;
}
}
}
};
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, "");
ParallelComm *pcomm = ParallelComm::get_pcomm(&moab, MYPCOMM_INDEX);
if (pcomm == NULL)
pcomm = new ParallelComm(&moab, PETSC_COMM_WORLD);
// 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_IMPOSIBLE_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_type(0, MBTRI, "FIELD1");
CHKERR m_field.add_ents_to_field_by_type(0, MBQUAD, "FIELD1");
// Set order
CHKERR m_field.set_field_order(0, MBTRI, "FIELD1", order);
CHKERR m_field.set_field_order(0, MBQUAD, "FIELD1", order);
CHKERR m_field.set_field_order(0, MBEDGE, "FIELD1", order);
// Add entities to elements
CHKERR m_field.add_ents_to_finite_element_by_type(0, MBTRI, "FACE_FE");
CHKERR m_field.add_ents_to_finite_element_by_type(0, MBQUAD, "FACE_FE");
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);
CHKERR m_field.add_ents_to_finite_element_by_type(faces_skin, MBEDGE,
"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;
MatrixDouble inv_jac(2, 2), jac(2, 2);
fe_face.getOpPtrVector().push_back(new OpCalculateJacForFace(jac));
fe_face.getOpPtrVector().push_back(new OpCalculateInvJacForFace(inv_jac));
fe_face.getOpPtrVector().push_back(new OpMakeHdivFromHcurl());
fe_face.getOpPtrVector().push_back(
fe_face.getOpPtrVector().push_back(new OpSetInvJacHcurlFace(inv_jac));
fe_face.getOpPtrVector().push_back(
new OpMakeHighOrderGeometryWeightsOnFace());
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;
fe_edge.getOpPtrVector().push_back(
new OpSetContravariantPiolaTransformOnEdge());
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();
PetscPrintf(PETSC_COMM_WORLD, "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);
}
}
static Index< 'p', 3 > p
ForcesAndSourcesCore::UserDataOperator UserDataOperator
MoFEM::FaceElementForcesAndSourcesCore FaceEle
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:441
FieldApproximationBase
approximation base
Definition: definitions.h:150
@ 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
@ HCURL
field with continuous tangents
Definition: definitions.h:178
#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_IMPOSIBLE_CASE
Definition: definitions.h:127
@ MOFEM_ATOM_TEST_INVALID
Definition: definitions.h:132
@ MOFEM_INVALID_DATA
Definition: definitions.h:128
#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
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.
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
FaceEle::UserDataOperator FaceEleOp
int main(int argc, char *argv[])
static char help[]
FTensor::Index< 'i', 3 > i
double tol
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
OpSetInvJacHcurlFaceImpl< 2 > OpSetInvJacHcurlFace
OpCalculateJacForFaceImpl< 2 > OpCalculateJacForFace
PetscErrorCode PetscOptionsGetInt(PetscOptions *, const char pre[], const char name[], PetscInt *ivalue, PetscBool *set)
OpCalculateInvJacForFaceImpl< 2 > OpCalculateInvJacForFace
CoreTmp< 0 > Core
Definition: Core.hpp:1129
OpSetContravariantPiolaTransformFaceImpl< 2 > OpSetContravariantPiolaTransformFace
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.
virtual moab::Interface & get_moab()=0
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.