v0.15.0
Loading...
Searching...
No Matches
child_and_parent.cpp
Go to the documentation of this file.
1/**
2 * @file child_and_parent.cpp
3 * @brief projection from parent to child element
4 * @example mofem/atom_tests/child_and_parent.cpp
5 * @copyright Copyright (c) 2025
6 *
7 */
8
9#include <MoFEM.hpp>
10
11using namespace MoFEM;
12
13static char help[] = "...\n\n";
14
15constexpr char FIELD_NAME[] = "U";
16constexpr int FIELD_DIM = 1;
17constexpr int SPACE_DIM = 2;
18
19template <int DIM> struct ElementsAndOps {};
20
29
34
37using DomainEleOp = DomainEle::UserDataOperator;
39
41using BoundaryEleOp = BoundaryEle::UserDataOperator;
43
44template <int FIELD_DIM> struct ApproxFieldFunction;
45
46template <> struct ApproxFieldFunction<1> {
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) +
49 pow(y, 4);
50 }
51};
52
57
58struct AtomTest {
59
60 AtomTest(MoFEM::Interface &m_field) : mField(m_field) {}
61
63
64private:
67
69
75 checkResults(boost::function<bool(FEMethod *fe_method_ptr)> test_bit);
77 struct CommonData {
78 boost::shared_ptr<VectorDouble> approxVals;
81 };
82
83 template <int FIELD_DIM> struct OpError;
84};
85
88
89template <> struct AtomTest::OpError<1> : public DomainEleOp {
90 boost::shared_ptr<CommonData> commonDataPtr;
91 OpError(boost::shared_ptr<CommonData> &common_data_ptr)
92 : DomainEleOp(FIELD_NAME, OPROW), commonDataPtr(common_data_ptr) {}
93 MoFEMErrorCode doWork(int side, EntityType type, EntData &data) {
95
96 if (const size_t nb_dofs = data.getIndices().size()) {
97
98 const int nb_integration_pts = getGaussPts().size2();
99 auto t_w = getFTensor0IntegrationWeight();
100 auto t_val = getFTensor0FromVec(*(commonDataPtr->approxVals));
101 auto t_coords = getFTensor1CoordsAtGaussPts();
102
103 VectorDouble nf(nb_dofs, false);
104 nf.clear();
105
106 FTensor::Index<'i', 3> i;
107 const double volume = getMeasure();
108
109 auto t_row_base = data.getFTensor0N();
110 double error = 0;
111 for (int gg = 0; gg != nb_integration_pts; ++gg) {
112
113 const double alpha = t_w * volume;
114 double diff = t_val - AtomTest::approxFunction(t_coords(0), t_coords(1),
115 t_coords(2));
116 error += alpha * pow(diff, 2);
117
118 for (size_t r = 0; r != nb_dofs; ++r) {
119 nf[r] += alpha * t_row_base * diff;
120 ++t_row_base;
121 }
122
123 ++t_w;
124 ++t_val;
125 ++t_coords;
126 }
127
128 const int index = 0;
129 CHKERR VecSetValue(commonDataPtr->L2Vec, index, error, ADD_VALUES);
130 CHKERR VecSetValues(commonDataPtr->resVec, data, &nf[0], ADD_VALUES);
131 }
132
134 }
135};
136
137template <typename ELE_OP, typename PARENT_ELE>
138struct OpCheckGaussCoords : public ELE_OP {
140
141 MoFEMErrorCode doWork(int side, EntityType type,
144
145 MatrixDouble parent_coords;
146
147 PARENT_ELE parent_fe(this->getPtrFE()->mField);
148 auto op = new ELE_OP(NOSPACE, ELE_OP::OPSPACE);
149 op->doWorkRhsHook = [&](DataOperator *op_ptr, int side,
150 EntityType type,
153
154 MOFEM_LOG("SELF", Sev::noisy)
155 << "parent_coords in op "
156 << static_cast<ELE_OP *>(op_ptr)->getCoordsAtGaussPts();
157
158 parent_coords = static_cast<ELE_OP *>(op_ptr)->getCoordsAtGaussPts();
160 };
161 parent_fe.getOpPtrVector().push_back(op);
162
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;
166
167 MatrixDouble child_coords = this->getCoordsAtGaussPts();
168 MOFEM_LOG("SELF", Sev::noisy) << "child_coords " << child_coords;
169
170 child_coords -= parent_coords;
171
172 MOFEM_LOG("SELF", Sev::noisy) << "Corrds diffs" << child_coords;
173
174 double n = 0;
175 for (auto d : child_coords.data())
176 n += d * d;
177
178 if (sqrt(n) > 1e-12)
179 SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
180 "Parent and child global coords at integration points are "
181 "diffrent norm = %3.2e",
182 sqrt(n));
183
185 }
186};
187
188//! [Run programme]
195
196 auto test_bit_child = [](FEMethod *fe_ptr) {
197 const auto &bit = fe_ptr->numeredEntFiniteElementPtr->getBitRefLevel();
198 MOFEM_LOG("SELF", Sev::noisy) << bit << " " << bit.test(0);
199 return bit.test(1);
200 };
201
205}
206//! [Run programme]
207
208//! [Read mesh]
211
215
216 MOFEM_LOG("WORLD", Sev::verbose) << "Dim " << simpleInterface->getDim();
217
218 auto bit_level0 = simpleInterface->getBitRefLevel();
219
220 auto &moab = mField.get_moab();
221
222 auto refine_mesh = [&](auto bit_level1) {
224
225 auto refine = mField.getInterface<MeshRefinement>();
226
227 auto meshset_level0_ptr = get_temp_meshset_ptr(moab);
228 CHKERR mField.getInterface<BitRefManager>()->getEntitiesByRefLevel(
229 bit_level0, BitRefLevel().set(), *meshset_level0_ptr);
230
231 // random mesh refinement
232 auto meshset_ref_edges_ptr = get_temp_meshset_ptr(moab);
233 Range edges_to_refine;
234 CHKERR moab.get_entities_by_type(*meshset_level0_ptr, MBEDGE,
235 edges_to_refine);
236 int ii = 0;
237 for (Range::iterator eit = edges_to_refine.begin();
238 eit != edges_to_refine.end(); eit++, ii++) {
239 int numb = ii % 2;
240 if (numb == 0) {
241 CHKERR moab.add_entities(*meshset_ref_edges_ptr, &*eit, 1);
242 }
243 }
244 CHKERR refine->addVerticesInTheMiddleOfEdges(*meshset_ref_edges_ptr,
245 bit_level1, false, VERBOSE);
246 if (simpleInterface->getDim() == 3) {
247 CHKERR refine->refineTets(*meshset_level0_ptr, bit_level1, VERBOSE);
248 } else if (simpleInterface->getDim() == 2) {
249 CHKERR refine->refineTris(*meshset_level0_ptr, bit_level1, VERBOSE);
250 } else {
251 SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
252 "Dimension not handled by test");
253 }
254
256 };
257
258 BitRefLevel bit_level1;
259 bit_level1.set(1);
260 CHKERR refine_mesh(bit_level1);
263
265}
266//! [Read mesh]
267
268//! [Set up problem]
286//! [Set up problem]
287
288boost::shared_ptr<DomainEle> domainChildLhs, domainChildRhs;
289
290//! [Push operators to pipeline]
294
295 auto rule = [](int, int, int p) -> int { return 2 * p; };
296
297 CHKERR pipeline_mng->setDomainLhsIntegrationRule(rule);
298 CHKERR pipeline_mng->setDomainRhsIntegrationRule(rule);
299
300 auto test_bit_parent = [](FEMethod *fe_ptr) {
301 const auto &bit = fe_ptr->numeredEntFiniteElementPtr->getBitRefLevel();
302 MOFEM_LOG("SELF", Sev::noisy) << bit << " " << bit.test(0);
303 return bit.test(0);
304 };
305
306 pipeline_mng->getDomainLhsFE()->exeTestHook = test_bit_parent;
307 pipeline_mng->getDomainRhsFE()->exeTestHook = test_bit_parent;
308
309 auto beta = [](const double, const double, const double) { return 1; };
310
311 // Make aliased shared pointer, and create child element
312 domainChildLhs = boost::make_shared<DomainEle>(mField);
313 domainChildLhs->getRuleHook = rule;
314 domainChildLhs->getOpPtrVector().push_back(
316
317 domainChildRhs = boost::make_shared<DomainEle>(mField);
318 domainChildLhs->getRuleHook = rule;
319 domainChildRhs->getOpPtrVector().push_back(
321
322 auto parent_op_lhs = new DomainEleOp(NOSPACE, DomainEleOp::OPSPACE);
323 parent_op_lhs->doWorkRhsHook = [&](DataOperator *op_ptr, int side,
324 EntityType type,
326 auto domain_op = static_cast<DomainEleOp *>(op_ptr);
328
329 MOFEM_LOG("SELF", Sev::noisy) << "LHS Pipeline FE";
330
331 if (!domainChildLhs)
332 SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID, "FE not allocated");
333
334 auto &bit =
335 domain_op->getFEMethod()->numeredEntFiniteElementPtr->getBitRefLevel();
336 if (bit == BitRefLevel().set(0)) {
337 CHKERR domain_op->loopChildren(domain_op->getFEName(),
338 domainChildLhs.get(), VERBOSE, Sev::noisy);
339 } else {
340 CHKERR domain_op->loopThis(domain_op->getFEName(), domainChildLhs.get(),
341 VERBOSE, Sev::noisy);
342 }
344 };
345
346 auto parent_op_rhs = new DomainEleOp(NOSPACE, DomainEleOp::OPSPACE);
347 parent_op_rhs->doWorkRhsHook = [&](DataOperator *op_ptr, int side,
348 EntityType type,
350 auto domain_op = static_cast<DomainEleOp *>(op_ptr);
352
353 MOFEM_LOG("SELF", Sev::noisy) << "RHS Pipeline FE";
354
355 if (!domainChildRhs)
356 SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID, "FE not allocated");
357
358 auto &bit =
359 domain_op->getFEMethod()->numeredEntFiniteElementPtr->getBitRefLevel();
360 if (bit == BitRefLevel().set(0)) {
361 CHKERR domain_op->loopChildren(domain_op->getFEName(),
362 domainChildRhs.get(), VERBOSE, Sev::noisy);
363 } else if ((bit & BitRefLevel().set(0)).any()) {
364 CHKERR domain_op->loopThis(domain_op->getFEName(), domainChildRhs.get(),
365 VERBOSE, Sev::noisy);
366 }
368 };
369
370 pipeline_mng->getOpDomainLhsPipeline().push_back(parent_op_lhs);
371 pipeline_mng->getOpDomainRhsPipeline().push_back(parent_op_rhs);
372
374}
375//! [Push operators to pipeline]
376
377//! [Solve]
381
382 MOFEM_LOG("WORLD", Sev::inform) << "Solve problem";
383
384 auto solver = pipeline_mng->createKSP();
385 CHKERR KSPSetFromOptions(solver);
386 CHKERR KSPSetUp(solver);
387
388 auto dm = simpleInterface->getDM();
389 auto D = createDMVector(dm);
390 auto F = vectorDuplicate(D);
391
392 CHKERR KSPSolve(solver, F, D);
393 CHKERR VecGhostUpdateBegin(D, INSERT_VALUES, SCATTER_FORWARD);
394 CHKERR VecGhostUpdateEnd(D, INSERT_VALUES, SCATTER_FORWARD);
395 CHKERR DMoFEMMeshToLocalVector(dm, D, INSERT_VALUES, SCATTER_REVERSE);
397}
398
399//! [Solve]
402
403 auto &moab = mField.get_moab();
404
405 auto bit_level0 = BitRefLevel().set(0);
406 auto bit_level1 = BitRefLevel().set(1);
407 auto bit_level2 = BitRefLevel().set(2);
408
409 auto refine_mesh = [&]() {
411
412 auto refine = mField.getInterface<MeshRefinement>();
413
414 auto meshset_level1_ptr = get_temp_meshset_ptr(moab);
415 CHKERR mField.getInterface<BitRefManager>()->getEntitiesByDimAndRefLevel(
416 bit_level1, BitRefLevel().set(), simpleInterface->getDim(),
417 *meshset_level1_ptr);
418
419 // random mesh refinement
420 auto meshset_ref_edges_ptr = get_temp_meshset_ptr(moab);
421 Range edges_to_refine;
422 CHKERR mField.getInterface<BitRefManager>()->getEntitiesByTypeAndRefLevel(
423 bit_level1, BitRefLevel().set(), MBEDGE, edges_to_refine);
424
425 CHKERR refine->addVerticesInTheMiddleOfEdges(edges_to_refine, bit_level2,
426 VERBOSE);
427 if (simpleInterface->getDim() == 3) {
428 CHKERR refine->refineTets(*meshset_level1_ptr, bit_level2, VERBOSE);
429 } else if (simpleInterface->getDim() == 2) {
430 CHKERR refine->refineTris(*meshset_level1_ptr, bit_level2, VERBOSE);
431 } else {
432 SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
433 "Dimension not handled by test");
434 }
435
436 Range meshsets;
437 CHKERR moab.get_entities_by_type(0, MBENTITYSET, meshsets, true);
438 for (auto m : meshsets) {
440 ->updateMeshsetByEntitiesChildren(m, bit_level2, m, MBMAXTYPE, false);
441 }
442
444 };
445
446 CHKERR refine_mesh();
447
448 simpleInterface->getBitRefLevel() = bit_level1 | bit_level2;
450
452
453 CHKERR mField.getInterface<ProblemsManager>()->removeDofsOnEntities(
455 bit_level0 | bit_level1);
456
457 auto project_data = [&]() {
459
461
462 pipeline_mng->getDomainLhsFE().reset();
463 pipeline_mng->getDomainRhsFE().reset();
464 pipeline_mng->getBoundaryRhsFE().reset();
465
466 auto rule = [](int, int, int p) -> int { return 2 * p; };
467
468 CHKERR pipeline_mng->setDomainLhsIntegrationRule(rule);
469 CHKERR pipeline_mng->setDomainRhsIntegrationRule(rule);
470 CHKERR pipeline_mng->setBoundaryRhsIntegrationRule(rule);
471
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);
475 return bit.test(2);
476 };
477
478 pipeline_mng->getDomainLhsFE()->exeTestHook = test_bit_ref;
479 pipeline_mng->getDomainRhsFE()->exeTestHook = test_bit_ref;
480 pipeline_mng->getBoundaryRhsFE()->exeTestHook = test_bit_ref;
481
482 auto beta = [](const double, const double, const double) { return 1; };
483 auto field_vals_ptr = boost::make_shared<VectorDouble>();
484
485 auto domainParentRhs = boost::make_shared<DomainParentEle>(mField);
486 domainParentRhs->getOpPtrVector().push_back(
487 new OpCalculateScalarFieldValues(FIELD_NAME, field_vals_ptr));
488
489 pipeline_mng->getOpDomainLhsPipeline().push_back(
491 pipeline_mng->getOpDomainRhsPipeline().push_back(
492 new OpRunParent(domainParentRhs, bit_level2, bit_level2,
493 domainParentRhs, bit_level2, BitRefLevel().set()));
494
497 pipeline_mng->getOpDomainRhsPipeline().push_back(
498 new OpDomainTimesScalarField(FIELD_NAME, field_vals_ptr, beta));
499
500 pipeline_mng->getOpDomainRhsPipeline().push_back(
502
503 pipeline_mng->getOpBoundaryRhsPipeline().push_back(
505
507
508 simpleInterface->getBitRefLevel() = bit_level2;
511
512 CHKERR checkResults([](FEMethod *fe_ptr) { return true; });
513
515 };
516
517 CHKERR project_data();
518
520}
521//! [Postprocess results]
522
523//! [Check results]
525 boost::function<bool(FEMethod *fe_method_ptr)> test_bit) {
528 pipeline_mng->getDomainLhsFE().reset();
529 pipeline_mng->getDomainRhsFE().reset();
530 pipeline_mng->getBoundaryRhsFE().reset();
531
532 auto rule = [](int, int, int p) -> int { return 2 * p + 1; };
533 CHKERR pipeline_mng->setDomainRhsIntegrationRule(rule);
534 pipeline_mng->getDomainRhsFE()->exeTestHook = test_bit;
535
536 auto common_data_ptr = boost::make_shared<CommonData>();
537 common_data_ptr->resVec = createDMVector(simpleInterface->getDM());
538 common_data_ptr->L2Vec = createVectorMPI(
539 mField.get_comm(), (!mField.get_comm_rank()) ? 1 : 0, 1);
540 common_data_ptr->approxVals = boost::make_shared<VectorDouble>();
541
542 pipeline_mng->getOpDomainRhsPipeline().push_back(
544 common_data_ptr->approxVals));
545 pipeline_mng->getOpDomainRhsPipeline().push_back(
546 new OpError<FIELD_DIM>(common_data_ptr));
547
548 CHKERR VecZeroEntries(common_data_ptr->L2Vec);
549 CHKERR VecZeroEntries(common_data_ptr->resVec);
550
551 CHKERR pipeline_mng->loopFiniteElements();
552
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);
557 double nrm2;
558 CHKERR VecNorm(common_data_ptr->resVec, NORM_2, &nrm2);
559 const double *array;
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);
564
565 constexpr double eps = 1e-8;
566 if (nrm2 > eps)
567 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
568 "Not converged solution err = %6.4e", nrm2);
570}
571//! [Check results]
572
573int main(int argc, char *argv[]) {
574
575 // Initialisation of MoFEM/PETSc and MOAB data structures
576 MoFEM::Core::Initialize(&argc, &argv, NULL, help);
577
578 try {
579
580 //! [Register MoFEM discrete manager in PETSc]
581 DMType dm_name = "DMMOFEM";
582 CHKERR DMRegister_MoFEM(dm_name);
583 //! [Register MoFEM discrete manager in PETSc
584
585 //! [Create MoAB]
586 moab::Core mb_instance; ///< mesh database
587 moab::Interface &moab = mb_instance; ///< mesh database interface
588 //! [Create MoAB]
589
590 //! [Create MoFEM]
591 MoFEM::Core core(moab); ///< finite element database
592 MoFEM::Interface &m_field = core; ///< finite element database insterface
593 //! [Create MoFEM]
594
595 //! [AtomTest]
596 AtomTest ex(m_field);
597 CHKERR ex.runProblem();
598 //! [AtomTest]
599 }
601
603}
#define MOFEM_LOG_C(channel, severity, format,...)
int main()
static const double eps
ElementsAndOps< SPACE_DIM >::DomianParentEle DomainParentEle
ElementsAndOps< SPACE_DIM >::BoundaryParentEle BoundaryParentEle
static char help[]
constexpr char FIELD_NAME[]
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpSource< 1, FIELD_DIM > OpDomainSource
constexpr int SPACE_DIM
boost::shared_ptr< DomainEle > domainChildRhs
constexpr int FIELD_DIM
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
@ VERBOSE
#define CATCH_ERRORS
Catch errors.
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
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:514
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:1119
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
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
double D
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.
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.
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
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:118
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.
structure for User Loop Methods on finite elements
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.
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:261
int getDim() const
Get the problem dimension.
Definition Simple.hpp:333
MoFEMErrorCode loadFile(const std::string options, const std::string mesh_file_name, LoadFileFunc loadFunc=defaultLoadFileFunc)
Load mesh file.
Definition Simple.cpp:191
MoFEMErrorCode reSetUp(bool only_dm=false)
Rebuild internal MoFEM data structures.
Definition Simple.cpp:762
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:355
MoFEMErrorCode getOptions()
get options
Definition Simple.cpp:180
MoFEMErrorCode getDM(DM *dm)
Get DM.
Definition Simple.cpp:800
MoFEMErrorCode setFieldOrder(const std::string field_name, const int order, const Range *ents=NULL)
Set field order.
Definition Simple.cpp:575
BitRefLevel & getBitRefLevelMask()
Get the BitRefLevel.
Definition Simple.hpp:382
MoFEMErrorCode setUp(const PetscBool is_partitioned=PETSC_TRUE)
Setup problem.
Definition Simple.cpp:736
const std::string getProblemName() const
Get the Problem Name.
Definition Simple.hpp:410
BitRefLevel & getBitRefLevel()
Get the BitRefLevel.
Definition Simple.hpp:375
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)