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

Testing projection child and parent.

/**
* \example child_and_parent.cpp
*
* Testing projection child and parent.
*
*/
#include <MoFEM.hpp>
using namespace MoFEM;
static char help[] = "...\n\n";
constexpr char FIELD_NAME[] = "U";
constexpr int FIELD_DIM = 1;
constexpr int SPACE_DIM = 2;
template <int DIM> struct ElementsAndOps {};
template <> struct ElementsAndOps<2> {
using DomainEleOp = DomainEle::UserDataOperator;
using BoundaryEleOp = BoundaryEle::UserDataOperator;
};
template <> struct ElementsAndOps<3> {
using DomainEle = VolumeElementForcesAndSourcesCore;
using DomainEleOp = DomainEle::UserDataOperator;
};
using DomainEleOp = DomainEle::UserDataOperator;
using BoundaryEleOp = BoundaryEle::UserDataOperator;
template <int FIELD_DIM> struct ApproxFieldFunction;
template <> struct ApproxFieldFunction<1> {
double operator()(const double x, const double y, const double z) {
return x * x + y * y + x * y + pow(x, 3) + pow(y, 3) + pow(x, 4) +
pow(y, 4);
}
};
struct AtomTest {
AtomTest(MoFEM::Interface &m_field) : mField(m_field) {}
private:
checkResults(boost::function<bool(FEMethod *fe_method_ptr)> test_bit);
struct CommonData {
boost::shared_ptr<VectorDouble> approxVals;
};
template <int FIELD_DIM> struct OpError;
};
template <> struct AtomTest::OpError<1> : public DomainEleOp {
boost::shared_ptr<CommonData> commonDataPtr;
OpError(boost::shared_ptr<CommonData> &common_data_ptr)
: DomainEleOp(FIELD_NAME, OPROW), commonDataPtr(common_data_ptr) {}
MoFEMErrorCode doWork(int side, EntityType type, EntData &data) {
if (const size_t nb_dofs = data.getIndices().size()) {
const int nb_integration_pts = getGaussPts().size2();
auto t_w = getFTensor0IntegrationWeight();
auto t_val = getFTensor0FromVec(*(commonDataPtr->approxVals));
auto t_coords = getFTensor1CoordsAtGaussPts();
VectorDouble nf(nb_dofs, false);
nf.clear();
FTensor::Index<'i', 3> i;
const double volume = getMeasure();
auto t_row_base = data.getFTensor0N();
double error = 0;
for (int gg = 0; gg != nb_integration_pts; ++gg) {
const double alpha = t_w * volume;
double diff = t_val - AtomTest::approxFunction(t_coords(0), t_coords(1),
t_coords(2));
error += alpha * pow(diff, 2);
for (size_t r = 0; r != nb_dofs; ++r) {
nf[r] += alpha * t_row_base * diff;
++t_row_base;
}
++t_w;
++t_val;
++t_coords;
}
const int index = 0;
CHKERR VecSetValue(commonDataPtr->L2Vec, index, error, ADD_VALUES);
CHKERR VecSetValues(commonDataPtr->resVec, data, &nf[0], ADD_VALUES);
}
}
};
template <typename ELE_OP, typename PARENT_ELE>
struct OpCheckGaussCoords : public ELE_OP {
MoFEMErrorCode doWork(int side, EntityType type,
MatrixDouble parent_coords;
PARENT_ELE parent_fe(this->getPtrFE()->mField);
auto op = new ELE_OP(NOSPACE, ELE_OP::OPSPACE);
op->doWorkRhsHook = [&](DataOperator *op_ptr, int side,
EntityType type,
MOFEM_LOG("SELF", Sev::noisy)
<< "parent_coords in op "
<< static_cast<ELE_OP *>(op_ptr)->getCoordsAtGaussPts();
parent_coords = static_cast<ELE_OP *>(op_ptr)->getCoordsAtGaussPts();
};
parent_fe.getOpPtrVector().push_back(op);
MOFEM_LOG("SELF", Sev::noisy) << "fe name " << this->getFEName();
CHKERR this->loopParent(this->getFEName(), &parent_fe);
MOFEM_LOG("SELF", Sev::noisy) << "parent_coords " << parent_coords;
MatrixDouble child_coords = this->getCoordsAtGaussPts();
MOFEM_LOG("SELF", Sev::noisy) << "child_coords " << child_coords;
child_coords -= parent_coords;
MOFEM_LOG("SELF", Sev::noisy) << "Corrds diffs" << child_coords;
double n = 0;
for (auto d : child_coords.data())
n += d * d;
if (sqrt(n) > 1e-12)
SETERRQ1(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
"Parent and child global coords at integration points are "
"diffrent norm = %3.2e",
sqrt(n));
}
};
//! [Run programme]
auto test_bit_child = [](FEMethod *fe_ptr) {
const auto &bit = fe_ptr->numeredEntFiniteElementPtr->getBitRefLevel();
MOFEM_LOG("SELF", Sev::noisy) << bit << " " << bit.test(0);
return bit.test(1);
};
}
//! [Run programme]
//! [Read mesh]
MOFEM_LOG("WORLD", Sev::verbose) << "Dim " << simpleInterface->getDim();
auto bit_level0 = simpleInterface->getBitRefLevel();
auto &moab = mField.get_moab();
auto refine_mesh = [&](auto bit_level1) {
auto meshset_level0_ptr = get_temp_meshset_ptr(moab);
CHKERR mField.getInterface<BitRefManager>()->getEntitiesByRefLevel(
bit_level0, BitRefLevel().set(), *meshset_level0_ptr);
// random mesh refinement
auto meshset_ref_edges_ptr = get_temp_meshset_ptr(moab);
Range edges_to_refine;
CHKERR moab.get_entities_by_type(*meshset_level0_ptr, MBEDGE,
edges_to_refine);
int ii = 0;
for (Range::iterator eit = edges_to_refine.begin();
eit != edges_to_refine.end(); eit++, ii++) {
int numb = ii % 2;
if (numb == 0) {
CHKERR moab.add_entities(*meshset_ref_edges_ptr, &*eit, 1);
}
}
CHKERR refine->addVerticesInTheMiddleOfEdges(*meshset_ref_edges_ptr,
bit_level1, false, VERBOSE);
if (simpleInterface->getDim() == 3) {
CHKERR refine->refineTets(*meshset_level0_ptr, bit_level1, VERBOSE);
} else if (simpleInterface->getDim() == 2) {
CHKERR refine->refineTris(*meshset_level0_ptr, bit_level1, VERBOSE);
} else {
SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
"Dimension not handled by test");
}
};
BitRefLevel bit_level1;
bit_level1.set(1);
CHKERR refine_mesh(bit_level1);
}
//! [Read mesh]
//! [Set up problem]
// Add field
constexpr int order = 4;
CHKERR mField.getInterface<ProblemsManager>()->removeDofsOnEntities(
simpleInterface->getProblemName(), FIELD_NAME, BitRefLevel().set(0),
BitRefLevel().set(0));
}
//! [Set up problem]
boost::shared_ptr<DomainEle> domainChildLhs, domainChildRhs;
//! [Push operators to pipeline]
auto rule = [](int, int, int p) -> int { return 2 * p; };
CHKERR pipeline_mng->setDomainLhsIntegrationRule(rule);
CHKERR pipeline_mng->setDomainRhsIntegrationRule(rule);
auto test_bit_parent = [](FEMethod *fe_ptr) {
const auto &bit = fe_ptr->numeredEntFiniteElementPtr->getBitRefLevel();
MOFEM_LOG("SELF", Sev::noisy) << bit << " " << bit.test(0);
return bit.test(0);
};
pipeline_mng->getDomainLhsFE()->exeTestHook = test_bit_parent;
pipeline_mng->getDomainRhsFE()->exeTestHook = test_bit_parent;
auto beta = [](const double, const double, const double) { return 1; };
// Make aliased shared pointer, and create child element
domainChildLhs = boost::make_shared<DomainEle>(mField);
domainChildLhs->getRuleHook = rule;
domainChildLhs->getOpPtrVector().push_back(
domainChildRhs = boost::make_shared<DomainEle>(mField);
domainChildLhs->getRuleHook = rule;
domainChildRhs->getOpPtrVector().push_back(
auto parent_op_lhs = new DomainEleOp(NOSPACE, DomainEleOp::OPSPACE);
parent_op_lhs->doWorkRhsHook = [&](DataOperator *op_ptr, int side,
EntityType type,
auto domain_op = static_cast<DomainEleOp *>(op_ptr);
MOFEM_LOG("SELF", Sev::noisy) << "LHS Pipeline FE";
SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID, "FE not allocated");
auto &bit =
domain_op->getFEMethod()->numeredEntFiniteElementPtr->getBitRefLevel();
if (bit == BitRefLevel().set(0)) {
CHKERR domain_op->loopChildren(domain_op->getFEName(),
domainChildLhs.get(), VERBOSE, Sev::noisy);
} else {
CHKERR domain_op->loopThis(domain_op->getFEName(), domainChildLhs.get(),
VERBOSE, Sev::noisy);
}
};
auto parent_op_rhs = new DomainEleOp(NOSPACE, DomainEleOp::OPSPACE);
parent_op_rhs->doWorkRhsHook = [&](DataOperator *op_ptr, int side,
EntityType type,
auto domain_op = static_cast<DomainEleOp *>(op_ptr);
MOFEM_LOG("SELF", Sev::noisy) << "RHS Pipeline FE";
SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID, "FE not allocated");
auto &bit =
domain_op->getFEMethod()->numeredEntFiniteElementPtr->getBitRefLevel();
if (bit == BitRefLevel().set(0)) {
CHKERR domain_op->loopChildren(domain_op->getFEName(),
domainChildRhs.get(), VERBOSE, Sev::noisy);
} else if ((bit & BitRefLevel().set(0)).any()) {
CHKERR domain_op->loopThis(domain_op->getFEName(), domainChildRhs.get(),
VERBOSE, Sev::noisy);
}
};
pipeline_mng->getOpDomainLhsPipeline().push_back(parent_op_lhs);
pipeline_mng->getOpDomainRhsPipeline().push_back(parent_op_rhs);
}
//! [Push operators to pipeline]
//! [Solve]
MOFEM_LOG("WORLD", Sev::inform) << "Solve problem";
auto solver = pipeline_mng->createKSP();
CHKERR KSPSetFromOptions(solver);
CHKERR KSPSetUp(solver);
auto dm = simpleInterface->getDM();
auto D = createDMVector(dm);
auto F = vectorDuplicate(D);
CHKERR KSPSolve(solver, F, D);
CHKERR VecGhostUpdateBegin(D, INSERT_VALUES, SCATTER_FORWARD);
CHKERR VecGhostUpdateEnd(D, INSERT_VALUES, SCATTER_FORWARD);
CHKERR DMoFEMMeshToLocalVector(dm, D, INSERT_VALUES, SCATTER_REVERSE);
}
//! [Solve]
auto &moab = mField.get_moab();
auto bit_level0 = BitRefLevel().set(0);
auto bit_level1 = BitRefLevel().set(1);
auto bit_level2 = BitRefLevel().set(2);
auto refine_mesh = [&]() {
auto meshset_level1_ptr = get_temp_meshset_ptr(moab);
CHKERR mField.getInterface<BitRefManager>()->getEntitiesByDimAndRefLevel(
bit_level1, BitRefLevel().set(), simpleInterface->getDim(),
*meshset_level1_ptr);
// random mesh refinement
auto meshset_ref_edges_ptr = get_temp_meshset_ptr(moab);
Range edges_to_refine;
CHKERR mField.getInterface<BitRefManager>()->getEntitiesByTypeAndRefLevel(
bit_level1, BitRefLevel().set(), MBEDGE, edges_to_refine);
CHKERR refine->addVerticesInTheMiddleOfEdges(edges_to_refine, bit_level2,
if (simpleInterface->getDim() == 3) {
CHKERR refine->refineTets(*meshset_level1_ptr, bit_level2, VERBOSE);
} else if (simpleInterface->getDim() == 2) {
CHKERR refine->refineTris(*meshset_level1_ptr, bit_level2, VERBOSE);
} else {
SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
"Dimension not handled by test");
}
Range meshsets;
CHKERR moab.get_entities_by_type(0, MBENTITYSET, meshsets, true);
for (auto m : meshsets) {
->updateMeshsetByEntitiesChildren(m, bit_level2, m, MBMAXTYPE, false);
}
};
CHKERR refine_mesh();
simpleInterface->getBitRefLevel() = bit_level1 | bit_level2;
CHKERR mField.getInterface<ProblemsManager>()->removeDofsOnEntities(
simpleInterface->getProblemName(), FIELD_NAME, BitRefLevel().set(),
bit_level0 | bit_level1);
auto project_data = [&]() {
pipeline_mng->getDomainLhsFE().reset();
pipeline_mng->getDomainRhsFE().reset();
pipeline_mng->getBoundaryRhsFE().reset();
auto rule = [](int, int, int p) -> int { return 2 * p; };
CHKERR pipeline_mng->setDomainLhsIntegrationRule(rule);
CHKERR pipeline_mng->setDomainRhsIntegrationRule(rule);
auto test_bit_ref = [](FEMethod *fe_ptr) {
const auto &bit = fe_ptr->numeredEntFiniteElementPtr->getBitRefLevel();
MOFEM_LOG("SELF", Sev::noisy) << "ref : " << bit << " " << bit.test(2);
return bit.test(2);
};
pipeline_mng->getDomainLhsFE()->exeTestHook = test_bit_ref;
pipeline_mng->getDomainRhsFE()->exeTestHook = test_bit_ref;
pipeline_mng->getBoundaryRhsFE()->exeTestHook = test_bit_ref;
auto beta = [](const double, const double, const double) { return 1; };
auto field_vals_ptr = boost::make_shared<VectorDouble>();
auto domainParentRhs = boost::make_shared<DomainParentEle>(mField);
domainParentRhs->getOpPtrVector().push_back(
pipeline_mng->getOpDomainLhsPipeline().push_back(
pipeline_mng->getOpDomainRhsPipeline().push_back(
new OpRunParent(domainParentRhs, bit_level2, bit_level2,
domainParentRhs, bit_level2, BitRefLevel().set()));
pipeline_mng->getOpDomainRhsPipeline().push_back(
new OpDomainTimesScalarField(FIELD_NAME, field_vals_ptr, beta));
pipeline_mng->getOpDomainRhsPipeline().push_back(
pipeline_mng->getOpBoundaryRhsPipeline().push_back(
CHKERR checkResults([](FEMethod *fe_ptr) { return true; });
};
CHKERR project_data();
}
//! [Postprocess results]
//! [Check results]
boost::function<bool(FEMethod *fe_method_ptr)> test_bit) {
pipeline_mng->getDomainLhsFE().reset();
pipeline_mng->getDomainRhsFE().reset();
pipeline_mng->getBoundaryRhsFE().reset();
auto rule = [](int, int, int p) -> int { return 2 * p + 1; };
CHKERR pipeline_mng->setDomainRhsIntegrationRule(rule);
pipeline_mng->getDomainRhsFE()->exeTestHook = test_bit;
auto common_data_ptr = boost::make_shared<CommonData>();
common_data_ptr->resVec = createDMVector(simpleInterface->getDM());
common_data_ptr->L2Vec = createVectorMPI(
mField.get_comm(), (!mField.get_comm_rank()) ? 1 : 0, 1);
common_data_ptr->approxVals = boost::make_shared<VectorDouble>();
pipeline_mng->getOpDomainRhsPipeline().push_back(
common_data_ptr->approxVals));
pipeline_mng->getOpDomainRhsPipeline().push_back(
new OpError<FIELD_DIM>(common_data_ptr));
CHKERR VecZeroEntries(common_data_ptr->L2Vec);
CHKERR VecZeroEntries(common_data_ptr->resVec);
CHKERR pipeline_mng->loopFiniteElements();
CHKERR VecAssemblyBegin(common_data_ptr->L2Vec);
CHKERR VecAssemblyEnd(common_data_ptr->L2Vec);
CHKERR VecAssemblyBegin(common_data_ptr->resVec);
CHKERR VecAssemblyEnd(common_data_ptr->resVec);
double nrm2;
CHKERR VecNorm(common_data_ptr->resVec, NORM_2, &nrm2);
const double *array;
CHKERR VecGetArrayRead(common_data_ptr->L2Vec, &array);
MOFEM_LOG_C("WORLD", Sev::inform, "Error %6.4e Vec norm %6.4e\n",
std::sqrt(array[0]), nrm2);
CHKERR VecRestoreArrayRead(common_data_ptr->L2Vec, &array);
constexpr double eps = 1e-8;
if (nrm2 > eps)
SETERRQ1(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
"Not converged solution err = %6.4e", nrm2);
}
//! [Check results]
int main(int argc, char *argv[]) {
// Initialisation of MoFEM/PETSc and MOAB data structures
MoFEM::Core::Initialize(&argc, &argv, NULL, help);
try {
//! [Register MoFEM discrete manager in PETSc]
DMType dm_name = "DMMOFEM";
//! [Register MoFEM discrete manager in PETSc
//! [Create MoAB]
moab::Core mb_instance; ///< mesh database
moab::Interface &moab = mb_instance; ///< mesh database interface
//! [Create MoAB]
//! [Create MoFEM]
MoFEM::Core core(moab); ///< finite element database
MoFEM::Interface &m_field = core; ///< finite element database insterface
//! [Create MoFEM]
//! [AtomTest]
AtomTest ex(m_field);
CHKERR ex.runProblem();
//! [AtomTest]
}
}
#define MOFEM_LOG_C(channel, severity, format,...)
static char help[]
int main()
static const double eps
ElementsAndOps< SPACE_DIM >::DomianParentEle DomainParentEle
ElementsAndOps< SPACE_DIM >::BoundaryParentEle BoundaryParentEle
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpSource< 1, FIELD_DIM > OpDomainSource
constexpr char FIELD_NAME[]
DomainEle::UserDataOperator DomainEleOp
Finire element operator type.
constexpr int SPACE_DIM
boost::shared_ptr< DomainEle > domainChildRhs
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< 1, FIELD_DIM > OpDomainMass
constexpr int FIELD_DIM
boost::shared_ptr< DomainEle > domainChildLhs
[Set up problem]
ElementsAndOps< SPACE_DIM >::BoundaryEle BoundaryEle
ElementsAndOps< SPACE_DIM >::DomainEle DomainEle
[Define dimension]
@ VERBOSE
#define CATCH_ERRORS
Catch errors.
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base nme:nme847.
Definition definitions.h:60
@ H1
continuous field
Definition definitions.h:85
@ NOSPACE
Definition definitions.h:83
#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
@ MOFEM_DATA_INCONSISTENCY
Definition definitions.h:31
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
constexpr int order
@ F
PetscErrorCode DMoFEMMeshToLocalVector(DM dm, Vec l, InsertMode mode, ScatterMode scatter_mode)
set local (or ghosted) vector values on mesh for partition only
Definition DMMoFEM.cpp:523
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
Definition DMMoFEM.cpp:43
auto createDMVector(DM dm)
Get smart vector from DM.
Definition DMMoFEM.hpp:1099
MoFEMErrorCode loopFiniteElements(SmartPetscObj< DM > dm=nullptr)
Iterate finite elements.
boost::ptr_deque< UserDataOperator > & getOpDomainLhsPipeline()
Get the Op Domain Lhs Pipeline object.
SmartPetscObj< KSP > createKSP(SmartPetscObj< DM > dm=nullptr)
Create KSP (linear) solver.
boost::ptr_deque< UserDataOperator > & getOpBoundaryRhsPipeline()
Get the Op Boundary Rhs Pipeline object.
boost::ptr_deque< UserDataOperator > & getOpDomainRhsPipeline()
Get the Op Domain Rhs Pipeline object.
#define MOFEM_LOG(channel, severity)
Log.
auto test_bit_child
lambda function used to select elements on which finite element pipelines are executed.
auto bit
set bit
FTensor::Index< 'i', SPACE_DIM > i
constexpr char FIELD_NAME[]
constexpr int FIELD_DIM
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< BASE_DIM, FIELD_DIM > OpDomainMass
OPerator to integrate mass matrix for least square approximation.
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpSource< BASE_DIM, FIELD_DIM > OpDomainSource
Operator to integrate the right hand side matrix for the problem.
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpBaseTimesScalar< 1 > OpDomainTimesScalarField
double D
const double n
refractive index of diffusive medium
PipelineManager::ElementsAndOpsByDim< FE_DIM >::DomianParentEle DomianParentEle
Definition level_set.cpp:38
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition Types.hpp:40
implementation of Data Operators for Forces and Sources
Definition Common.hpp:10
SmartPetscObj< Vec > vectorDuplicate(Vec vec)
Create duplicate vector of smart vector.
auto createVectorMPI(MPI_Comm comm, PetscInt n, PetscInt N)
Create MPI Vector.
static auto getFTensor0FromVec(ublas::vector< T, A > &data)
Get tensor rank 0 (scalar) form data vector.
auto get_temp_meshset_ptr(moab::Interface &moab)
Create smart pointer to temporary meshset.
MoFEMErrorCode VecSetValues(Vec V, const EntitiesFieldData::EntData &data, const double *ptr, InsertMode iora)
Assemble PETSc vector.
int r
Definition sdf.py:8
FTensor::Index< 'm', 3 > m
boost::shared_ptr< VectorDouble > approxVals
SmartPetscObj< Vec > L2Vec
SmartPetscObj< Vec > resVec
Operator to evaluate errors.
boost::shared_ptr< CommonData > commonDataPtr
OpError(boost::shared_ptr< MatrixDouble > data_ptr)
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
Simple * simpleInterface
MoFEMErrorCode checkResults()
[Check results]
MoFEMErrorCode solveSystem()
[Push operators to pipeline]
static ApproxFieldFunction< FIELD_DIM > approxFunction
AtomTest(MoFEM::Interface &m_field)
MoFEMErrorCode setupProblem()
[Read mesh]
MoFEMErrorCode refineResults()
[Solve]
MoFEMErrorCode readMesh()
[Run programme]
MoFEMErrorCode assembleSystem()
[Push operators to pipeline]
MoFEM::Interface & mField
MoFEMErrorCode runProblem()
[Run programme]
Managing BitRefLevels.
virtual moab::Interface & get_moab()=0
virtual MPI_Comm & get_comm() const =0
virtual int get_comm_rank() const =0
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
base operator to do operations at Gauss Pt. level
Deprecated interface functions.
Data on single entity (This is passed as argument to DataOperator::doWork)
FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > getFTensor0N(const FieldApproximationBase base)
Get base function as Tensor0.
const VectorInt & getIndices() const
Get global indices of dofs on entity.
Mesh refinement interface.
Get value at integration points for scalar field.
Operator to execute finite element instance on parent element. This operator is typically used to pro...
PipelineManager interface.
boost::shared_ptr< FEMethod > & getDomainRhsFE()
boost::shared_ptr< FEMethod > & getDomainLhsFE()
MoFEMErrorCode setDomainRhsIntegrationRule(RuleHookFun rule)
MoFEMErrorCode setBoundaryRhsIntegrationRule(RuleHookFun rule)
boost::shared_ptr< FEMethod > & getBoundaryRhsFE()
MoFEMErrorCode setDomainLhsIntegrationRule(RuleHookFun rule)
Problem manager is used to build and partition problems.
Simple interface for fast problem set-up.
Definition Simple.hpp:27
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
int getDim() const
Get the problem dimension.
Definition Simple.hpp:331
MoFEMErrorCode loadFile(const std::string options, const std::string mesh_file_name, LoadFileFunc loadFunc=defaultLoadFileFunc)
Load mesh file.
Definition Simple.cpp:194
MoFEMErrorCode reSetUp(bool only_dm=false)
Rebuild internal MoFEM data structures.
Definition Simple.cpp:791
MoFEMErrorCode addBoundaryField(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 boundary.
Definition Simple.cpp:358
MoFEMErrorCode getOptions()
get options
Definition Simple.cpp:180
MoFEMErrorCode getDM(DM *dm)
Get DM.
Definition Simple.cpp:829
MoFEMErrorCode setFieldOrder(const std::string field_name, const int order, const Range *ents=NULL)
Set field order.
Definition Simple.cpp:599
BitRefLevel & getBitRefLevelMask()
Get the BitRefLevel.
Definition Simple.hpp:380
MoFEMErrorCode setUp(const PetscBool is_partitioned=PETSC_TRUE)
Setup problem.
Definition Simple.cpp:765
const std::string getProblemName() const
Get the Problem Name.
Definition Simple.hpp:408
BitRefLevel & getBitRefLevel()
Get the BitRefLevel.
Definition Simple.hpp:373
intrusive_ptr for managing petsc objects
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)