v0.15.5
Loading...
Searching...
No Matches
FieldEvaluator.cpp
Go to the documentation of this file.
1/** \file FieldEvaluator.cpp
2 * \brief Field evaluator
3 *
4 */
5
6namespace MoFEM {
7
9 : cOre(const_cast<MoFEM::Core &>(core)) {
10
11 if (!LogManager::checkIfChannelExist("FieldEvaluatorWorld")) {
12 auto core_log = logging::core::get();
13
15 "FieldEvaluatorWorld"));
17 "FieldEvaluatorSync"));
19 "FieldEvaluatorSelf"));
20
21 LogManager::setLog("FieldEvaluatorWorld");
22 LogManager::setLog("FieldEvaluatorSync");
23 LogManager::setLog("FieldEvaluatorSelf");
24
25 MOFEM_LOG_TAG("FieldEvaluatorWorld", "FieldEvaluator");
26 MOFEM_LOG_TAG("FieldEvaluatorSync", "FieldEvaluator");
27 MOFEM_LOG_TAG("FieldEvaluatorSelf", "FieldEvaluator");
28 }
29
30 MOFEM_LOG("FieldEvaluatorWorld", Sev::noisy) << "Field evaluator intreface";
31}
32
34 boost::typeindex::type_index type_index, UnknownInterface **iface) const {
36 *iface = const_cast<FieldEvaluatorInterface *>(this);
38}
39
40template <int D>
42 boost::shared_ptr<SetIntegrationPtsMethodData> spd_ptr,
43 const std::string finite_element, BitRefLevel bit, BitRefLevel mask) {
44 CoreInterface &m_field = cOre;
46 EntityHandle fe_meshset;
47 fe_meshset = m_field.get_finite_element_meshset(finite_element);
48 Range entities;
49 CHKERR m_field.get_moab().get_entities_by_dimension(fe_meshset, D, entities,
50 true);
51 if (bit.any()) {
52 CHKERR m_field.getInterface<BitRefManager>()->filterEntitiesByRefLevel(
53 bit, mask, entities);
54 }
55 if (spd_ptr->treePtr) {
56 spd_ptr->treePtr->reset_tree();
57 } else {
58 spd_ptr->treePtr = boost::make_shared<AdaptiveKDTree>(&m_field.get_moab());
59 }
60 CHKERR spd_ptr->treePtr->build_tree(entities, &spd_ptr->rooTreeSet);
62}
63
64template <>
65MoFEMErrorCode FieldEvaluatorInterface::buildTree<2>(
66 boost::shared_ptr<SetIntegrationPtsMethodData> spd_ptr,
67 const std::string finite_element, BitRefLevel bit, BitRefLevel mask) {
68 return buildTreeImpl<2>(spd_ptr, finite_element, bit, mask);
69}
70
71template <>
72MoFEMErrorCode FieldEvaluatorInterface::buildTree<3>(
73 boost::shared_ptr<SetIntegrationPtsMethodData> spd_ptr,
74 const std::string finite_element, BitRefLevel bit, BitRefLevel mask) {
75 return buildTreeImpl<3>(spd_ptr, finite_element, bit, mask);
76}
77
79 boost::shared_ptr<SetIntegrationPtsMethodData> data_ptr)
80 : dataPtr(data_ptr) {}
81
83 ForcesAndSourcesCore *fe_raw_ptr, int order_row, int order_col,
84 int order_data) {
86
87 if (auto data_ptr = dataPtr.lock()) {
88
89 const auto nb_eval_points = data_ptr->nbEvalPoints;
90 const auto &shape_functions = data_ptr->shapeFunctions;
91 const auto &eval_pointentity_handle = data_ptr->evalPointEntityHandle;
92
93 if (auto fe_ptr = data_ptr->feMethodPtr) {
94
95 auto &fe = static_cast<ForcesAndSourcesCore &>(*fe_ptr);
96
97#ifndef NDEBUG
98 if (fe_ptr.get() != fe_raw_ptr)
99 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "Wrong FE ptr");
100#endif
101
102 const auto fe_ent = fe.numeredEntFiniteElementPtr->getEnt();
103 const auto fe_type = type_from_handle(fe_ent);
104 const auto fe_dim = moab::CN::Dimension(fe_type);
105
106 auto &gauss_pts = fe.gaussPts;
107 int nb_gauss_pts;
108
109 if (fe_dim == 3) {
110 gauss_pts.resize(4, nb_eval_points, false);
112 &shape_functions(0, 0), &shape_functions(0, 1),
113 &shape_functions(0, 2), &shape_functions(0, 3)};
115 &gauss_pts(0, 0), &gauss_pts(1, 0), &gauss_pts(2, 0)};
116
117 nb_gauss_pts = 0;
118 for (int nn = 0; nn != nb_eval_points; ++nn) {
119 if (eval_pointentity_handle[nn] == fe_ent) {
120 for (const int i : {0, 1, 2}) {
121 t_gauss_pts(i) = t_shape(i + 1);
122 }
123 gauss_pts(3, nb_gauss_pts) = nn;
124 ++t_gauss_pts;
125 ++nb_gauss_pts;
126 }
127 ++t_shape;
128 }
129 gauss_pts.resize(4, nb_gauss_pts, true);
130 } else if (fe_dim == 2) {
131 gauss_pts.resize(3, nb_eval_points, false);
133 &shape_functions(0, 0), &shape_functions(0, 1),
134 &shape_functions(0, 2)};
136 &gauss_pts(0, 0), &gauss_pts(1, 0)};
137 nb_gauss_pts = 0;
138 for (int nn = 0; nn != nb_eval_points; ++nn) {
139 if (eval_pointentity_handle[nn] == fe_ent) {
140 for (const int i : {0, 1}) {
141 t_gauss_pts(i) = t_shape(i + 1);
142 }
143 gauss_pts(2, nb_gauss_pts) = nn;
144 ++t_gauss_pts;
145 ++nb_gauss_pts;
146 }
147 ++t_shape;
148 }
149 gauss_pts.resize(3, nb_gauss_pts, true);
150 } else {
151 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
152 "Dimension not implemented");
153 }
154
155#ifndef NDEBUG
156 MOFEM_LOG("SELF", Sev::noisy)
157 << "nbEvalOPoints / nbGau_sspt_s: " << nb_eval_points << " / "
158 << nb_gauss_pts;
159 MOFEM_LOG("SELF", Sev::noisy) << "gauss pts: " << gauss_pts;
160#endif
161
162 } else
163 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
164 "Pointer to element does not exists");
165
166 } else
167 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
168 "Pointer to data does not exists");
169
171}
172
173template <int D>
175 const double *const point, const double distance, const std::string problem,
176 const std::string finite_element,
177 boost::shared_ptr<SetIntegrationPtsMethodData> data_ptr, int lower_rank,
178 int upper_rank, BitRefLevel bit, BitRefLevel mask,
179 boost::shared_ptr<CacheTuple> cache_ptr, MoFEMTypes bh,
180 VERBOSITY_LEVELS verb) {
181 CoreInterface &m_field = cOre;
183
184 std::vector<EntityHandle> leafs_out;
185 CHKERR data_ptr->treePtr->distance_search(point, distance, leafs_out);
186 Range tree_ents;
187 for (auto lit : leafs_out) //{
188 CHKERR m_field.get_moab().get_entities_by_dimension(lit, D, tree_ents,
189 true);
190
191 if (verb >= NOISY)
192 MOFEM_LOG("SELF", Sev::noisy) << "tree entities: " << tree_ents;
193
194 data_ptr->evalPointEntityHandle.resize(data_ptr->nbEvalPoints);
195 std::fill(data_ptr->evalPointEntityHandle.begin(),
196 data_ptr->evalPointEntityHandle.end(), 0);
197
198 std::vector<double> local_coords;
199 std::vector<double> shape;
200
201 auto nb_eval_points = data_ptr->nbEvalPoints;
202
203 for (auto tet : tree_ents) {
204 if ((type_from_handle(tet) != moab::MBTRI) &&
205 (type_from_handle(tet) != moab::MBTET)) {
206 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
207 "Wrong element type, Field Evaluator not implemented for Quads "
208 "and Hexes");
209 }
210 const EntityHandle *conn;
211 int num_nodes;
212 CHKERR m_field.get_moab().get_connectivity(tet, conn, num_nodes, true);
213
214 if constexpr (D == 3) {
215
216 local_coords.resize(3 * nb_eval_points);
217 shape.resize(4 * nb_eval_points);
218
219 std::array<double, 12> coords;
220 CHKERR m_field.get_moab().get_coords(conn, num_nodes, coords.data());
221
223 coords.data(), &data_ptr->evalPoints[0], nb_eval_points,
224 &local_coords[0]);
225 CHKERR Tools::shapeFunMBTET<3>(&shape[0], &local_coords[0],
226 &local_coords[1], &local_coords[2],
227 nb_eval_points);
228
229 FTensor::Index<'i', 4> i4;
231 &shape[0], &shape[1], &shape[2], &shape[3]};
233 &data_ptr->shapeFunctions(0, 0), &data_ptr->shapeFunctions(0, 1),
234 &data_ptr->shapeFunctions(0, 2), &data_ptr->shapeFunctions(0, 3)};
235
236 FTensor::Index<'j', 3> j3;
238 &local_coords[0], &local_coords[1], &local_coords[2]};
240 &(data_ptr->localCoords(0, 0)), &(data_ptr->localCoords(0, 1)),
241 &(data_ptr->localCoords(0, 2))};
242
243 for (int n = 0; n != nb_eval_points; ++n) {
244
245 const double eps = data_ptr->eps;
246 if (t_shape(0) >= 0 - eps && t_shape(0) <= 1 + eps &&
247
248 t_shape(1) >= 0 - eps && t_shape(1) <= 1 + eps &&
249
250 t_shape(2) >= 0 - eps && t_shape(2) <= 1 + eps &&
251
252 t_shape(3) >= 0 - eps && t_shape(3) <= 1 + eps) {
253
254 data_ptr->evalPointEntityHandle[n] = tet;
255 t_shape_data(i4) = t_shape(i4);
256 t_local_data(j3) = t_local(j3);
257 }
258
259 ++t_shape;
260 ++t_shape_data;
261 ++t_local;
262 ++t_local_data;
263 }
264 }
265
266 if constexpr (D == 2) {
267
268 local_coords.resize(2 * nb_eval_points);
269 shape.resize(3 * nb_eval_points);
270
271 std::array<double, 9> coords;
272 CHKERR m_field.get_moab().get_coords(conn, num_nodes, coords.data());
273
275 coords.data(), &data_ptr->evalPoints[0], nb_eval_points,
276 &local_coords[0]);
277 CHKERR Tools::shapeFunMBTRI<2>(&shape[0], &local_coords[0],
278 &local_coords[1], nb_eval_points);
279
280 FTensor::Index<'i', 3> i3;
282 &shape[0], &shape[1], &shape[2]};
284 &data_ptr->shapeFunctions(0, 0), &data_ptr->shapeFunctions(0, 1),
285 &data_ptr->shapeFunctions(0, 2)};
286
287 FTensor::Index<'j', 2> j2;
289 &local_coords[0], &local_coords[1]};
290 data_ptr->localCoords.resize(nb_eval_points, 2);
292 &(data_ptr->localCoords(0, 0)), &(data_ptr->localCoords(0, 1))};
293
294 for (int n = 0; n != nb_eval_points; ++n) {
295
296 const double eps = data_ptr->eps;
297 if (t_shape(0) >= 0 - eps && t_shape(0) <= 1 + eps &&
298
299 t_shape(1) >= 0 - eps && t_shape(1) <= 1 + eps &&
300
301 t_shape(2) >= 0 - eps && t_shape(2) <= 1 + eps) {
302
303 data_ptr->evalPointEntityHandle[n] = tet;
304 t_shape_data(i3) = t_shape(i3);
305 t_local_data(j2) = t_local(j2);
306 }
307
308 ++t_shape;
309 ++t_shape_data;
310 ++t_local;
311 ++t_local_data;
312 }
313 }
314 }
315
316 const Problem *prb_ptr;
317 CHKERR m_field.get_problem(problem, &prb_ptr);
318 boost::shared_ptr<NumeredEntFiniteElement_multiIndex> numered_fes(
320
321 Range in_tets;
322 in_tets.insert_list(data_ptr->evalPointEntityHandle.begin(),
323 data_ptr->evalPointEntityHandle.end());
324 in_tets = in_tets.subset_by_dimension(D);
325
326 if (verb >= NOISY)
327 MOFEM_LOG("SELF", Sev::noisy) << "in tets: " << in_tets << endl;
328
329 auto fe_miit =
330 m_field.get_finite_elements()->get<FiniteElement_name_mi_tag>().find(
331 finite_element);
332 if (fe_miit !=
333 m_field.get_finite_elements()->get<FiniteElement_name_mi_tag>().end()) {
334
335 for (auto peit = in_tets.pair_begin(); peit != in_tets.pair_end(); ++peit) {
336
337 auto lo =
338 prb_ptr->numeredFiniteElementsPtr->get<Unique_mi_tag>().lower_bound(
340 peit->first, (*fe_miit)->getFEUId()));
341 auto hi =
342 prb_ptr->numeredFiniteElementsPtr->get<Unique_mi_tag>().upper_bound(
344 peit->second, (*fe_miit)->getFEUId()));
345
346 if (bit.any()) {
347 if (verb >= VERY_NOISY)
348 MOFEM_LOG("FieldEvaluatorSelf", Sev::noisy);
349 for (auto it = lo; it != hi; ++it) {
350 auto fe_bit = (*it)->getBitRefLevel();
351 if ((fe_bit & mask) != fe_bit)
352 continue;
353 if ((fe_bit & bit).none())
354 continue;
355 numered_fes->insert(*it);
356 if (verb >= VERY_NOISY)
357 MOFEM_LOG("FieldEvaluatorSelf", Sev::noisy) << **it << endl;
358 }
359 } else {
360 numered_fes->insert(lo, hi);
361 if (verb >= VERY_NOISY) {
362 MOFEM_LOG("FieldEvaluatorSelf", Sev::noisy)
363 << "numbered elements:" << std::endl;
364 for (; lo != hi; ++lo)
365 MOFEM_LOG("FieldEvaluatorSelf", Sev::noisy) << **lo << endl;
366 }
367 }
368 }
369
370 if (auto fe_ptr = data_ptr->feMethodPtr) {
371
372 MOFEM_LOG_C("FieldEvaluatorSelf", Sev::noisy,
373 "Number elements %d to evaluate at proc %d",
374 numered_fes->size(), m_field.get_comm_rank());
375
376 if (!cache_ptr) {
377 cache_ptr = boost::make_shared<CacheTuple>();
378 CHKERR m_field.cache_problem_entities(prb_ptr->getName(), cache_ptr);
379
380 MOFEM_LOG("FieldEvaluatorSelf", Sev::noisy)
381 << "If you call that function many times in the loop consider to "
382 "set cache_ptr outside of the loop. Otherwise code can be slow.";
383 }
384
385 CHKERR m_field.loop_finite_elements(prb_ptr, finite_element, *fe_ptr,
386 lower_rank, upper_rank, numered_fes,
387 bh, cache_ptr, verb);
388
389 } else
390 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
391 "Pointer to element does not exists");
392 }
393
395}
396
397template <>
398MoFEMErrorCode FieldEvaluatorInterface::evalFEAtThePoint<2>(
399 const double *const point, const double distance, const std::string problem,
400 const std::string finite_element,
401 boost::shared_ptr<SetIntegrationPtsMethodData> data_ptr, int lower_rank,
402 int upper_rank, BitRefLevel bit, BitRefLevel mask,
403 boost::shared_ptr<CacheTuple> cache_ptr, MoFEMTypes bh,
404 VERBOSITY_LEVELS verb) {
405 return evalFEAtThePointImpl<2>(point, distance, problem, finite_element,
406 data_ptr, lower_rank, upper_rank, bit, mask,
407 cache_ptr, bh, verb);
408}
409
410template <>
411MoFEMErrorCode FieldEvaluatorInterface::evalFEAtThePoint<3>(
412 const double *const point, const double distance, const std::string problem,
413 const std::string finite_element,
414 boost::shared_ptr<SetIntegrationPtsMethodData> data_ptr, int lower_rank,
415 int upper_rank, BitRefLevel bit, BitRefLevel mask,
416 boost::shared_ptr<CacheTuple> cache_ptr, MoFEMTypes bh,
417 VERBOSITY_LEVELS verb) {
418 return evalFEAtThePointImpl<3>(point, distance, problem, finite_element,
419 data_ptr, lower_rank, upper_rank, bit, mask,
420 cache_ptr, bh, verb);
421}
422
423template <int D>
427 const std::string finite_element,
428 boost::shared_ptr<SetIntegrationPtsMethodData> data_ptr, int lower_rank,
429 int upper_rank, BitRefLevel bit, BitRefLevel mask,
430
431 MoFEMTypes bh, VERBOSITY_LEVELS verb) {
432
434
435 struct Op : public BaseOP {
436 Op(
437
439 const std::string finite_element,
440 boost::shared_ptr<SetIntegrationPtsMethodData> data_ptr, int lower_rank,
441 int upper_rank, BitRefLevel bit, BitRefLevel mask, double eps,
443
444 )
445 : BaseOP(NOSPACE, BaseOP::OPSPACE),
446 matDataPts(mat_data_pts), finiteElement(finite_element),
447 dataPtr(data_ptr), lowerRank(lower_rank), upperRank(upper_rank),
448 bitRefLevel(bit), maskRefLevel(mask), eps(eps), bhType(bh),
449 verbosity(verb) {}
450
451 MoFEMErrorCode doWork(int side, EntityType type,
454
455 MOFEM_LOG_CHANNEL("SELF");
456 MOFEM_LOG_TAG("SELF", "FieldEvaluator");
457
458 MoFEM::Interface &m_field = this->getPtrFE()->mField;
459 auto field_eval_ptr = m_field.getInterface<FieldEvaluatorInterface>();
460
461 auto nb_integration_pts = getGaussPts().size2();
462
463 auto prb_name = getPtrFE()->problemPtr->getName();
464 auto cache_ptr = getPtrFE()->getCacheWeakPtr();
465
466 auto coords_at_pts = getCoordsAtGaussPts();
467 for (int gg = 0; gg != nb_integration_pts; ++gg) {
468 double *coords_ptr = &coords_at_pts(gg, 0);
469
470 dataPtr->setEvalPoints(coords_ptr, 1);
471 CHKERR field_eval_ptr->evalFEAtThePoint<D>(
472 coords_ptr, eps, prb_name, finiteElement, dataPtr, lowerRank,
473 upperRank, bitRefLevel, maskRefLevel, cache_ptr.lock(), MF_EXIST,
474 QUIET);
475
476 // Get number of found finite elements, that consist gauss point.
477 // E.g. for edge flip expected is that only one element will be found.
478 auto loop_size = dataPtr->feMethodPtr->getLoopSize();
479 if (loop_size != 1) {
480 if (verbosity >= QUIET) {
482 LogManager::BitLineID | LogManager::BitScope);
483 MOFEM_LOG("SELF", Sev::warning)
484 << "Field evaluator wrong number of elements found: "
485 << loop_size;
486 }
487 }
488
489 for (auto &m : matDataPts) {
490 auto nb_rows = m.first->size1();
491 m.second->resize(nb_rows, nb_integration_pts, false);
492 if (gg == 0)
493 m.second->clear();
494 for (auto ll = 0; ll != loop_size; ++ll) {
495 for (int dd = 0; dd != nb_rows; ++dd) {
496 (*m.second)(dd, gg) += (*m.first)(dd, ll) / loop_size;
497 }
498 }
499 }
500 }
501
502 MOFEM_LOG_CHANNEL("SELF");
503
505 }
506
507 private:
509 std::string finiteElement;
510 boost::shared_ptr<SetIntegrationPtsMethodData> dataPtr;
511 int lowerRank;
512 int upperRank;
513 BitRefLevel bitRefLevel;
514 BitRefLevel maskRefLevel;
515 double eps;
516 MoFEMTypes bhType;
517 VERBOSITY_LEVELS verbosity;
518 };
519
520 double min[3];
521 double max[3];
522 unsigned int dep;
523 CHKERR data_ptr->treePtr->get_info(data_ptr->rooTreeSet, min, max, dep);
524 double dist = std::sqrt(pow(max[0] - min[0], 2) + pow(max[1] - min[1], 2) +
525 pow(max[2] - min[2], 2));
526
527 return new Op( mat_data_pts, finite_element, data_ptr,
528 lower_rank, upper_rank, bit, mask, 1e-6 * dist, bh, verb);
529}
530
533 boost::shared_ptr<MoFEM::ForcesAndSourcesCore> fe_method_ptr,
534 const double *eval_points, const int nb_eval_points, const double eps,
535 VERBOSITY_LEVELS verb)
536 : feMethodPtr(fe_method_ptr), evalPoints(eval_points),
537 nbEvalPoints(nb_eval_points), eps(eps), verb(verb) {
538 localCoords.resize(nbEvalPoints, 3);
539 shapeFunctions.resize(nbEvalPoints, 4);
540}
541
543 const double *ptr, const int nb_eval_points) {
544 evalPoints = ptr;
545 nbEvalPoints = nb_eval_points;
546 localCoords.resize(nbEvalPoints, 3, false);
547 shapeFunctions.resize(nbEvalPoints, 4, false);
548}
549
550template <>
552FieldEvaluatorInterface::getDataOperator<2>(
554 const std::string finite_element,
555 boost::shared_ptr<SetIntegrationPtsMethodData> data_ptr, int lower_rank,
556 int upper_rank, BitRefLevel bit, BitRefLevel mask, MoFEMTypes bh,
557 VERBOSITY_LEVELS verb) {
558 return getDataOperatorImpl<2>(mat_data_pts, finite_element,
559 data_ptr, lower_rank, upper_rank, bit, mask, bh,
560 verb);
561}
562
563template <>
565FieldEvaluatorInterface::getDataOperator<3>(
567 const std::string finite_element,
568 boost::shared_ptr<SetIntegrationPtsMethodData> data_ptr, int lower_rank,
569 int upper_rank, BitRefLevel bit, BitRefLevel mask, MoFEMTypes bh,
570 VERBOSITY_LEVELS verb) {
571 return getDataOperatorImpl<3>(mat_data_pts, finite_element,
572 data_ptr, lower_rank, upper_rank, bit, mask, bh,
573 verb);
574}
575
576} // namespace MoFEM
#define MOFEM_LOG_C(channel, severity, format,...)
static const double eps
VERBOSITY_LEVELS
Verbosity levels.
@ QUIET
@ NOISY
@ VERY_NOISY
MoFEMTypes
Those types control how functions respond on arguments, f.e. error handling.
@ MF_EXIST
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
@ NOSPACE
Definition definitions.h:83
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ MOFEM_DATA_INCONSISTENCY
Definition definitions.h:31
@ MOFEM_NOT_IMPLEMENTED
Definition definitions.h:32
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
multi_index_container< boost::shared_ptr< NumeredEntFiniteElement >, indexed_by< ordered_unique< tag< Unique_mi_tag >, const_mem_fun< NumeredEntFiniteElement::interface_type_EntFiniteElement, UId, &NumeredEntFiniteElement::getLocalUniqueId > >, ordered_non_unique< tag< Part_mi_tag >, member< NumeredEntFiniteElement, unsigned int, &NumeredEntFiniteElement::part > >, ordered_non_unique< tag< Ent_mi_tag >, const_mem_fun< NumeredEntFiniteElement::interface_type_RefEntity, EntityHandle, &NumeredEntFiniteElement::getEnt > >, ordered_non_unique< tag< Composite_Name_And_Part_mi_tag >, composite_key< NumeredEntFiniteElement, const_mem_fun< NumeredEntFiniteElement::interface_type_FiniteElement, boost::string_ref, &NumeredEntFiniteElement::getNameRef >, member< NumeredEntFiniteElement, unsigned int, &NumeredEntFiniteElement::part > > > > > NumeredEntFiniteElement_multiIndex
MultiIndex for entities for NumeredEntFiniteElement.
virtual const Problem * get_problem(const std::string problem_name) const =0
Get the problem object.
virtual const FiniteElement_multiIndex * get_finite_elements() const =0
Get the finite elements object.
virtual EntityHandle get_finite_element_meshset(const std::string name) const =0
static LoggerType & setLog(const std::string channel)
Set ans resset chanel logger.
#define MOFEM_LOG(channel, severity)
Log.
#define MOFEM_LOG_TAG(channel, tag)
Tag channel.
#define MOFEM_LOG_CHANNEL(channel)
Set and reset channel.
#define MOFEM_LOG_ATTRIBUTES(channel, bit)
Add attributes to channel.
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.
auto bit
set bit
FTensor::Index< 'i', SPACE_DIM > i
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
auto type_from_handle(const EntityHandle h)
get type from entity handle
FTensor::Index< 'm', 3 > m
Managing BitRefLevels.
virtual moab::Interface & get_moab()=0
virtual MoFEMErrorCode cache_problem_entities(const std::string prb_name, CacheTupleWeakPtr cache_ptr)=0
Cache variables.
virtual int get_comm_rank() const =0
Core (interface) class.
Definition Core.hpp:82
Deprecated interface functions.
UId getLocalUniqueIdCalculate() const
Generate UId for finite element entity.
Data on single entity (This is passed as argument to DataOperator::doWork)
void setEvalPoints(const double *ptr, const int nb_eval_points)
MoFEMErrorCode operator()(ForcesAndSourcesCore *fe_raw_ptr, int order_row, int order_col, int order_data)
Field evaluator interface.
MoFEMErrorCode query_interface(boost::typeindex::type_index type_index, UnknownInterface **iface) const
std::vector< std::pair< boost::shared_ptr< MatrixDouble >, boost::shared_ptr< MatrixDouble > > > MatDataPts
MoFEMErrorCode buildTreeImpl(boost::shared_ptr< SetIntegrationPtsMethodData > spd_ptr, const std::string finite_element, BitRefLevel bit, BitRefLevel mask)
FieldEvaluatorInterface(const MoFEM::Core &core)
MoFEMErrorCode evalFEAtThePointImpl(const double *const point, const double distance, const std::string problem, const std::string finite_element, boost::shared_ptr< SetIntegrationPtsMethodData > data_ptr, int lower_rank, int upper_rank, BitRefLevel bit, BitRefLevel mask, boost::shared_ptr< CacheTuple > cache_ptr, MoFEMTypes bh=MF_EXIST, VERBOSITY_LEVELS verb=QUIET)
ForcesAndSourcesCore::UserDataOperator * getDataOperatorImpl(MatDataPts mat_data_pts, const std::string finite_element, boost::shared_ptr< SetIntegrationPtsMethodData > data_ptr, int lower_rank, int upper_rank, BitRefLevel bit, BitRefLevel mask, MoFEMTypes bh=MF_EXIST, VERBOSITY_LEVELS verb=QUIET)
structure to get information from mofem into EntitiesFieldData
static boost::shared_ptr< SinkType > createSink(boost::shared_ptr< std::ostream > stream_ptr, std::string comm_filter)
Create a sink object.
static boost::shared_ptr< std::ostream > getStrmWorld()
Get the strm world object.
static boost::shared_ptr< std::ostream > getStrmSync()
Get the strm sync object.
static bool checkIfChannelExist(const std::string channel)
Check if channel exist.
static boost::shared_ptr< std::ostream > getStrmSelf()
Get the strm self object.
keeps basic data about problem
boost::shared_ptr< NumeredEntFiniteElement_multiIndex > numeredFiniteElementsPtr
store finite elements
static MoFEMErrorCode getLocalCoordinatesOnReferenceFourNodeTet(const double *elem_coords, const double *glob_coords, const int nb_nodes, double *local_coords)
Get the local coordinates on reference four node tet object.
Definition Tools.cpp:88
static MoFEMErrorCode getLocalCoordinatesOnReferenceThreeNodeTri(const double *elem_coords, const double *glob_coords, const int nb_nodes, double *local_coords)
Get the local coordinates on reference three node tri object.
Definition Tools.cpp:188
base class for all interface classes
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.