13static char help[] =
"...\n\n";
47 double operator()(
const double x,
const double y,
const double z) {
48 return x * x + y * y + x * y + pow(x, 3) + pow(y, 3) + pow(x, 4) +
83 template <
int FIELD_DIM>
struct OpError;
91 OpError(boost::shared_ptr<CommonData> &common_data_ptr)
96 if (
const size_t nb_dofs = data.
getIndices().size()) {
98 const int nb_integration_pts = getGaussPts().size2();
99 auto t_w = getFTensor0IntegrationWeight();
101 auto t_coords = getFTensor1CoordsAtGaussPts();
107 const double volume = getMeasure();
111 for (
int gg = 0; gg != nb_integration_pts; ++gg) {
113 const double alpha = t_w * volume;
116 error += alpha * pow(diff, 2);
118 for (
size_t r = 0; r != nb_dofs; ++r) {
119 nf[r] += alpha * t_row_base * diff;
137template <
typename ELE_OP,
typename PARENT_ELE>
147 PARENT_ELE parent_fe(this->getPtrFE()->mField);
155 <<
"parent_coords in op "
156 <<
static_cast<ELE_OP *
>(op_ptr)->getCoordsAtGaussPts();
158 parent_coords =
static_cast<ELE_OP *
>(op_ptr)->getCoordsAtGaussPts();
161 parent_fe.getOpPtrVector().push_back(op);
163 MOFEM_LOG(
"SELF", Sev::noisy) <<
"fe name " << this->getFEName();
164 CHKERR this->loopParent(this->getFEName(), &parent_fe);
165 MOFEM_LOG(
"SELF", Sev::noisy) <<
"parent_coords " << parent_coords;
167 MatrixDouble child_coords = this->getCoordsAtGaussPts();
168 MOFEM_LOG(
"SELF", Sev::noisy) <<
"child_coords " << child_coords;
170 child_coords -= parent_coords;
172 MOFEM_LOG(
"SELF", Sev::noisy) <<
"Corrds diffs" << child_coords;
175 for (
auto d : child_coords.data())
180 "Parent and child global coords at integration points are "
181 "diffrent norm = %3.2e",
197 const auto &
bit = fe_ptr->numeredEntFiniteElementPtr->getBitRefLevel();
222 auto refine_mesh = [&](
auto bit_level1) {
229 bit_level0,
BitRefLevel().set(), *meshset_level0_ptr);
233 Range edges_to_refine;
234 CHKERR moab.get_entities_by_type(*meshset_level0_ptr, MBEDGE,
237 for (Range::iterator eit = edges_to_refine.begin();
238 eit != edges_to_refine.end(); eit++, ii++) {
241 CHKERR moab.add_entities(*meshset_ref_edges_ptr, &*eit, 1);
244 CHKERR refine->addVerticesInTheMiddleOfEdges(*meshset_ref_edges_ptr,
247 CHKERR refine->refineTets(*meshset_level0_ptr, bit_level1,
VERBOSE);
249 CHKERR refine->refineTris(*meshset_level0_ptr, bit_level1,
VERBOSE);
252 "Dimension not handled by test");
260 CHKERR refine_mesh(bit_level1);
276 constexpr int order = 4;
295 auto rule = [](int, int,
int p) ->
int {
return 2 * p; };
300 auto test_bit_parent = [](
FEMethod *fe_ptr) {
301 const auto &
bit = fe_ptr->numeredEntFiniteElementPtr->getBitRefLevel();
323 parent_op_lhs->doWorkRhsHook = [&](
DataOperator *op_ptr,
int side,
326 auto domain_op =
static_cast<DomainEleOp *
>(op_ptr);
329 MOFEM_LOG(
"SELF", Sev::noisy) <<
"LHS Pipeline FE";
335 domain_op->getFEMethod()->numeredEntFiniteElementPtr->getBitRefLevel();
337 CHKERR domain_op->loopChildren(domain_op->getFEName(),
347 parent_op_rhs->doWorkRhsHook = [&](
DataOperator *op_ptr,
int side,
350 auto domain_op =
static_cast<DomainEleOp *
>(op_ptr);
353 MOFEM_LOG(
"SELF", Sev::noisy) <<
"RHS Pipeline FE";
359 domain_op->getFEMethod()->numeredEntFiniteElementPtr->getBitRefLevel();
361 CHKERR domain_op->loopChildren(domain_op->getFEName(),
382 MOFEM_LOG(
"WORLD", Sev::inform) <<
"Solve problem";
385 CHKERR KSPSetFromOptions(solver);
393 CHKERR VecGhostUpdateBegin(
D, INSERT_VALUES, SCATTER_FORWARD);
394 CHKERR VecGhostUpdateEnd(
D, INSERT_VALUES, SCATTER_FORWARD);
409 auto refine_mesh = [&]() {
417 *meshset_level1_ptr);
421 Range edges_to_refine;
423 bit_level1,
BitRefLevel().set(), MBEDGE, edges_to_refine);
425 CHKERR refine->addVerticesInTheMiddleOfEdges(edges_to_refine, bit_level2,
428 CHKERR refine->refineTets(*meshset_level1_ptr, bit_level2,
VERBOSE);
430 CHKERR refine->refineTris(*meshset_level1_ptr, bit_level2,
VERBOSE);
433 "Dimension not handled by test");
437 CHKERR moab.get_entities_by_type(0, MBENTITYSET, meshsets,
true);
438 for (
auto m : meshsets) {
440 ->updateMeshsetByEntitiesChildren(
m, bit_level2,
m, MBMAXTYPE,
false);
455 bit_level0 | bit_level1);
457 auto project_data = [&]() {
466 auto rule = [](int, int,
int p) ->
int {
return 2 * p; };
472 auto test_bit_ref = [](
FEMethod *fe_ptr) {
473 const auto &
bit = fe_ptr->numeredEntFiniteElementPtr->getBitRefLevel();
474 MOFEM_LOG(
"SELF", Sev::noisy) <<
"ref : " <<
bit <<
" " <<
bit.test(2);
483 auto field_vals_ptr = boost::make_shared<VectorDouble>();
485 auto domainParentRhs = boost::make_shared<DomainParentEle>(
mField);
486 domainParentRhs->getOpPtrVector().push_back(
492 new OpRunParent(domainParentRhs, bit_level2, bit_level2,
493 domainParentRhs, bit_level2,
BitRefLevel().set()));
525 boost::function<
bool(
FEMethod *fe_method_ptr)> test_bit) {
532 auto rule = [](int, int,
int p) ->
int {
return 2 * p + 1; };
536 auto common_data_ptr = boost::make_shared<CommonData>();
540 common_data_ptr->approxVals = boost::make_shared<VectorDouble>();
544 common_data_ptr->approxVals));
548 CHKERR VecZeroEntries(common_data_ptr->L2Vec);
549 CHKERR VecZeroEntries(common_data_ptr->resVec);
553 CHKERR VecAssemblyBegin(common_data_ptr->L2Vec);
554 CHKERR VecAssemblyEnd(common_data_ptr->L2Vec);
555 CHKERR VecAssemblyBegin(common_data_ptr->resVec);
556 CHKERR VecAssemblyEnd(common_data_ptr->resVec);
558 CHKERR VecNorm(common_data_ptr->resVec, NORM_2, &nrm2);
560 CHKERR VecGetArrayRead(common_data_ptr->L2Vec, &array);
561 MOFEM_LOG_C(
"WORLD", Sev::inform,
"Error %6.4e Vec norm %6.4e\n",
562 std::sqrt(array[0]), nrm2);
563 CHKERR VecRestoreArrayRead(common_data_ptr->L2Vec, &array);
565 constexpr double eps = 1e-8;
568 "Not converged solution err = %6.4e", nrm2);
573int main(
int argc,
char *argv[]) {
581 DMType dm_name =
"DMMOFEM";
586 moab::Core mb_instance;
587 moab::Interface &moab = mb_instance;
#define MOFEM_LOG_C(channel, severity, format,...)
ElementsAndOps< SPACE_DIM >::DomianParentEle DomainParentEle
ElementsAndOps< SPACE_DIM >::BoundaryParentEle BoundaryParentEle
constexpr char FIELD_NAME[]
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpSource< 1, FIELD_DIM > OpDomainSource
boost::shared_ptr< DomainEle > domainChildRhs
boost::shared_ptr< DomainEle > domainChildLhs
[Set up problem]
ElementsAndOps< SPACE_DIM >::DomainEle DomainEle
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< 1, FIELD_DIM > OpDomainMass
ElementsAndOps< SPACE_DIM >::BoundaryEle BoundaryEle
#define CATCH_ERRORS
Catch errors.
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ MOFEM_ATOM_TEST_INVALID
@ MOFEM_DATA_INCONSISTENCY
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
PetscErrorCode DMoFEMMeshToLocalVector(DM dm, Vec l, InsertMode mode, ScatterMode scatter_mode)
set local (or ghosted) vector values on mesh for partition only
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
auto createDMVector(DM dm)
Get smart vector from DM.
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.
FTensor::Index< 'i', SPACE_DIM > i
constexpr char FIELD_NAME[]
DomainEle::UserDataOperator DomainEleOp
Finire element operator type.
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
const double n
refractive index of diffusive medium
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
implementation of Data Operators for Forces and Sources
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.
FTensor::Index< 'm', 3 > m
double operator()(const double x, const double y, const double z)
Collected data use d by operator to evaluate errors for the test.
boost::shared_ptr< VectorDouble > approxVals
SmartPetscObj< Vec > L2Vec
SmartPetscObj< Vec > resVec
OpError(boost::shared_ptr< CommonData > &common_data_ptr)
boost::shared_ptr< CommonData > commonDataPtr
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
Operator to evaluate errors.
boost::shared_ptr< CommonData > commonDataPtr
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]
virtual moab::Interface & get_moab()=0
virtual MPI_Comm & get_comm() const =0
virtual int get_comm_rank() const =0
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
base operator to do operations at Gauss Pt. level
Deprecated interface functions.
Base face element used to integrate on skeleton.
default operator for EDGE element
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.
structure for User Loop Methods on finite elements
Base face element used to integrate on skeleton.
default operator for TRI element
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()
MoFEM::FaceElementForcesAndSourcesCore FaceEle
MoFEMErrorCode setDomainRhsIntegrationRule(RuleHookFun rule)
MoFEMErrorCode setBoundaryRhsIntegrationRule(RuleHookFun rule)
MoFEM::EdgeElementForcesAndSourcesCore EdgeEle
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.
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.
int getDim() const
Get the problem dimension.
MoFEMErrorCode loadFile(const std::string options, const std::string mesh_file_name, LoadFileFunc loadFunc=defaultLoadFileFunc)
Load mesh file.
MoFEMErrorCode reSetUp(bool only_dm=false)
Rebuild internal MoFEM data structures.
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.
MoFEMErrorCode getOptions()
get options
MoFEMErrorCode getDM(DM *dm)
Get DM.
MoFEMErrorCode setFieldOrder(const std::string field_name, const int order, const Range *ents=NULL)
Set field order.
BitRefLevel & getBitRefLevelMask()
Get the BitRefLevel.
MoFEMErrorCode setUp(const PetscBool is_partitioned=PETSC_TRUE)
Setup problem.
const std::string getProblemName() const
Get the Problem Name.
BitRefLevel & getBitRefLevel()
Get the BitRefLevel.
intrusive_ptr for managing petsc objects
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Volume finite element base.
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)