v0.15.0
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 spd_ptr->treePtr.reset(new AdaptiveKDTree(&m_field.get_moab()));
56 CHKERR spd_ptr->treePtr->build_tree(entities, &spd_ptr->rooTreeSet);
58}
59
60template <>
61MoFEMErrorCode FieldEvaluatorInterface::buildTree<2>(
62 boost::shared_ptr<SetIntegrationPtsMethodData> spd_ptr,
63 const std::string finite_element, BitRefLevel bit, BitRefLevel mask) {
64 return buildTreeImpl<2>(spd_ptr, finite_element, bit, mask);
65}
66
67template <>
68MoFEMErrorCode FieldEvaluatorInterface::buildTree<3>(
69 boost::shared_ptr<SetIntegrationPtsMethodData> spd_ptr,
70 const std::string finite_element, BitRefLevel bit, BitRefLevel mask) {
71 return buildTreeImpl<3>(spd_ptr, finite_element, bit, mask);
72}
73
75 boost::shared_ptr<SetIntegrationPtsMethodData> data_ptr)
76 : dataPtr(data_ptr) {}
77
79 ForcesAndSourcesCore *fe_raw_ptr, int order_row, int order_col,
80 int order_data) {
82
83 if (auto data_ptr = dataPtr.lock()) {
84
85 const auto nb_eval_points = data_ptr->nbEvalPoints;
86 const auto &shape_functions = data_ptr->shapeFunctions;
87 const auto &eval_pointentity_handle = data_ptr->evalPointEntityHandle;
88
89 if (auto fe_ptr = data_ptr->feMethodPtr.lock()) {
90
91 auto &fe = static_cast<ForcesAndSourcesCore &>(*fe_ptr);
92
93#ifndef NDEBUG
94 if (fe_ptr.get() != fe_raw_ptr)
95 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "Wrong FE ptr");
96#endif
97
98 const auto fe_ent = fe.numeredEntFiniteElementPtr->getEnt();
99 const auto fe_type = type_from_handle(fe_ent);
100 const auto fe_dim = moab::CN::Dimension(fe_type);
101
102 auto &gauss_pts = fe.gaussPts;
103 int nb_gauss_pts;
104
105 if (fe_dim == 3) {
106 gauss_pts.resize(4, nb_eval_points, false);
108 &shape_functions(0, 0), &shape_functions(0, 1),
109 &shape_functions(0, 2), &shape_functions(0, 3)};
111 &gauss_pts(0, 0), &gauss_pts(1, 0), &gauss_pts(2, 0)};
112
113 nb_gauss_pts = 0;
114 for (int nn = 0; nn != nb_eval_points; ++nn) {
115 if (eval_pointentity_handle[nn] == fe_ent) {
116 for (const int i : {0, 1, 2}) {
117 t_gauss_pts(i) = t_shape(i + 1);
118 }
119 gauss_pts(3, nb_gauss_pts) = nn;
120 ++t_gauss_pts;
121 ++nb_gauss_pts;
122 }
123 ++t_shape;
124 }
125 gauss_pts.resize(4, nb_gauss_pts, true);
126 } else if (fe_dim == 2) {
127 gauss_pts.resize(3, nb_eval_points, false);
129 &shape_functions(0, 0), &shape_functions(0, 1),
130 &shape_functions(0, 2)};
132 &gauss_pts(0, 0), &gauss_pts(1, 0)};
133 nb_gauss_pts = 0;
134 for (int nn = 0; nn != nb_eval_points; ++nn) {
135 if (eval_pointentity_handle[nn] == fe_ent) {
136 for (const int i : {0, 1}) {
137 t_gauss_pts(i) = t_shape(i + 1);
138 }
139 gauss_pts(2, nb_gauss_pts) = nn;
140 ++t_gauss_pts;
141 ++nb_gauss_pts;
142 }
143 ++t_shape;
144 }
145 gauss_pts.resize(3, nb_gauss_pts, true);
146 } else {
147 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
148 "Dimension not implemented");
149 }
150
151#ifndef NDEBUG
152 MOFEM_LOG("SELF", Sev::noisy)
153 << "nbEvalOPoints / nbGau_sspt_s: " << nb_eval_points << " / "
154 << nb_gauss_pts;
155 MOFEM_LOG("SELF", Sev::noisy) << "gauss pts: " << gauss_pts;
156#endif
157
158 } else
159 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
160 "Pointer to element does not exists");
161
162 } else
163 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
164 "Pointer to data does not exists");
165
167}
168
169template <int D>
171 const double *const point, const double distance, const std::string problem,
172 const std::string finite_element,
173 boost::shared_ptr<SetIntegrationPtsMethodData> data_ptr, int lower_rank,
174 int upper_rank, BitRefLevel bit, BitRefLevel mask,
175 boost::shared_ptr<CacheTuple> cache_ptr, MoFEMTypes bh,
176 VERBOSITY_LEVELS verb) {
177 CoreInterface &m_field = cOre;
179
180 std::vector<EntityHandle> leafs_out;
181 CHKERR data_ptr->treePtr->distance_search(point, distance, leafs_out);
182 Range tree_ents;
183 for (auto lit : leafs_out) //{
184 CHKERR m_field.get_moab().get_entities_by_dimension(lit, D, tree_ents,
185 true);
186
187 if (verb >= NOISY)
188 MOFEM_LOG("SELF", Sev::noisy) << "tree entities: " << tree_ents;
189
190 data_ptr->evalPointEntityHandle.resize(data_ptr->nbEvalPoints);
191 std::fill(data_ptr->evalPointEntityHandle.begin(),
192 data_ptr->evalPointEntityHandle.end(), 0);
193
194 std::vector<double> local_coords;
195 std::vector<double> shape;
196
197 auto nb_eval_points = data_ptr->nbEvalPoints;
198
199 for (auto tet : tree_ents) {
200 if ((type_from_handle(tet) != moab::MBTRI) &&
201 (type_from_handle(tet) != moab::MBTET)) {
202 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
203 "Wrong element type, Field Evaluator not implemented for Quads "
204 "and Hexes");
205 }
206 const EntityHandle *conn;
207 int num_nodes;
208 CHKERR m_field.get_moab().get_connectivity(tet, conn, num_nodes, true);
209
210 if constexpr (D == 3) {
211
212 local_coords.resize(3 * nb_eval_points);
213 shape.resize(4 * nb_eval_points);
214
215 std::array<double, 12> coords;
216 CHKERR m_field.get_moab().get_coords(conn, num_nodes, coords.data());
217
219 coords.data(), &data_ptr->evalPoints[0], nb_eval_points,
220 &local_coords[0]);
221 CHKERR Tools::shapeFunMBTET<3>(&shape[0], &local_coords[0],
222 &local_coords[1], &local_coords[2],
223 nb_eval_points);
224
225 FTensor::Index<'i', 4> i4;
227 &shape[0], &shape[1], &shape[2], &shape[3]};
229 &data_ptr->shapeFunctions(0, 0), &data_ptr->shapeFunctions(0, 1),
230 &data_ptr->shapeFunctions(0, 2), &data_ptr->shapeFunctions(0, 3)};
231
232 FTensor::Index<'j', 3> j3;
234 &local_coords[0], &local_coords[1], &local_coords[2]};
236 &(data_ptr->localCoords(0, 0)), &(data_ptr->localCoords(0, 1)),
237 &(data_ptr->localCoords(0, 2))};
238
239 for (int n = 0; n != nb_eval_points; ++n) {
240
241 const double eps = data_ptr->eps;
242 if (t_shape(0) >= 0 - eps && t_shape(0) <= 1 + eps &&
243
244 t_shape(1) >= 0 - eps && t_shape(1) <= 1 + eps &&
245
246 t_shape(2) >= 0 - eps && t_shape(2) <= 1 + eps &&
247
248 t_shape(3) >= 0 - eps && t_shape(3) <= 1 + eps) {
249
250 data_ptr->evalPointEntityHandle[n] = tet;
251 t_shape_data(i4) = t_shape(i4);
252 t_local_data(j3) = t_local(j3);
253 }
254
255 ++t_shape;
256 ++t_shape_data;
257 ++t_local;
258 ++t_local_data;
259 }
260 }
261
262 if constexpr (D == 2) {
263
264 local_coords.resize(2 * nb_eval_points);
265 shape.resize(3 * nb_eval_points);
266
267 std::array<double, 9> coords;
268 CHKERR m_field.get_moab().get_coords(conn, num_nodes, coords.data());
269
271 coords.data(), &data_ptr->evalPoints[0], nb_eval_points,
272 &local_coords[0]);
273 CHKERR Tools::shapeFunMBTRI<2>(&shape[0], &local_coords[0],
274 &local_coords[1], nb_eval_points);
275
276 FTensor::Index<'i', 3> i3;
278 &shape[0], &shape[1], &shape[2]};
280 &data_ptr->shapeFunctions(0, 0), &data_ptr->shapeFunctions(0, 1),
281 &data_ptr->shapeFunctions(0, 2)};
282
283 FTensor::Index<'j', 2> j2;
285 &local_coords[0], &local_coords[1]};
286 data_ptr->localCoords.resize(nb_eval_points, 2);
288 &(data_ptr->localCoords(0, 0)), &(data_ptr->localCoords(0, 1))};
289
290 for (int n = 0; n != nb_eval_points; ++n) {
291
292 const double eps = data_ptr->eps;
293 if (t_shape(0) >= 0 - eps && t_shape(0) <= 1 + eps &&
294
295 t_shape(1) >= 0 - eps && t_shape(1) <= 1 + eps &&
296
297 t_shape(2) >= 0 - eps && t_shape(2) <= 1 + eps) {
298
299 data_ptr->evalPointEntityHandle[n] = tet;
300 t_shape_data(i3) = t_shape(i3);
301 t_local_data(j2) = t_local(j2);
302 }
303
304 ++t_shape;
305 ++t_shape_data;
306 ++t_local;
307 ++t_local_data;
308 }
309 }
310 }
311
312 const Problem *prb_ptr;
313 CHKERR m_field.get_problem(problem, &prb_ptr);
314 boost::shared_ptr<NumeredEntFiniteElement_multiIndex> numered_fes(
316
317 Range in_tets;
318 in_tets.insert_list(data_ptr->evalPointEntityHandle.begin(),
319 data_ptr->evalPointEntityHandle.end());
320 in_tets = in_tets.subset_by_dimension(D);
321
322 if (verb >= NOISY)
323 MOFEM_LOG("SELF", Sev::noisy) << "in tets: " << in_tets << endl;
324
325 auto fe_miit =
326 m_field.get_finite_elements()->get<FiniteElement_name_mi_tag>().find(
327 finite_element);
328 if (fe_miit !=
329 m_field.get_finite_elements()->get<FiniteElement_name_mi_tag>().end()) {
330
331 for (auto peit = in_tets.pair_begin(); peit != in_tets.pair_end(); ++peit) {
332
333 auto lo =
334 prb_ptr->numeredFiniteElementsPtr->get<Unique_mi_tag>().lower_bound(
336 peit->first, (*fe_miit)->getFEUId()));
337 auto hi =
338 prb_ptr->numeredFiniteElementsPtr->get<Unique_mi_tag>().upper_bound(
340 peit->second, (*fe_miit)->getFEUId()));
341
342 if (bit.any()) {
343 for (auto it = lo; it != hi; ++it) {
344 auto fe_bit = (*it)->getBitRefLevel();
345 if ((fe_bit & mask) != fe_bit)
346 continue;
347 if ((fe_bit & bit).none())
348 continue;
349 numered_fes->insert(*it);
350 }
351 } else {
352 numered_fes->insert(lo, hi);
353 }
354
355 if (verb >= VERY_NOISY)
356 std::cerr << "numered elements:" << std::endl;
357 for (; lo != hi; ++lo)
358 if (verb >= VERY_NOISY)
359 std::cerr << **lo << endl;
360 }
361 if (verb >= VERY_NOISY)
362 std::cerr << std::endl;
363
364 if (auto fe_ptr = data_ptr->feMethodPtr.lock()) {
365
366 MOFEM_LOG_C("FieldEvaluatorSync", Sev::verbose,
367 "Number elements %d to evaluate at proc %d",
368 numered_fes->size(), m_field.get_comm_rank());
370
371 if (!cache_ptr) {
372 cache_ptr = boost::make_shared<CacheTuple>();
373 CHKERR m_field.cache_problem_entities(prb_ptr->getName(), cache_ptr);
374
375 MOFEM_LOG("FieldEvaluatorSync", Sev::noisy)
376 << "If you call that function many times in the loop consider to "
377 "set "
378 "cache_ptr outside of the loop. Otherwise code can be slow.";
379 }
380
381 CHKERR m_field.loop_finite_elements(prb_ptr, finite_element, *fe_ptr,
382 lower_rank, upper_rank, numered_fes,
383 bh, cache_ptr, verb);
384
385 } else
386 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
387 "Pointer to element does not exists");
388 }
389
391}
392
393template <>
394MoFEMErrorCode FieldEvaluatorInterface::evalFEAtThePoint<2>(
395 const double *const point, const double distance, const std::string problem,
396 const std::string finite_element,
397 boost::shared_ptr<SetIntegrationPtsMethodData> data_ptr, int lower_rank,
398 int upper_rank, BitRefLevel bit, BitRefLevel mask,
399 boost::shared_ptr<CacheTuple> cache_ptr, MoFEMTypes bh,
400 VERBOSITY_LEVELS verb) {
401 return evalFEAtThePointImpl<2>(point, distance, problem, finite_element,
402 data_ptr, lower_rank, upper_rank, bit, mask,
403 cache_ptr, bh, verb);
404}
405
406template <>
407MoFEMErrorCode FieldEvaluatorInterface::evalFEAtThePoint<3>(
408 const double *const point, const double distance, const std::string problem,
409 const std::string finite_element,
410 boost::shared_ptr<SetIntegrationPtsMethodData> data_ptr, int lower_rank,
411 int upper_rank, BitRefLevel bit, BitRefLevel mask,
412 boost::shared_ptr<CacheTuple> cache_ptr, MoFEMTypes bh,
413 VERBOSITY_LEVELS verb) {
414 return evalFEAtThePointImpl<3>(point, distance, problem, finite_element,
415 data_ptr, lower_rank, upper_rank, bit, mask,
416 cache_ptr, bh, verb);
417}
418
419template <int D>
424 const std::string finite_element,
425 boost::shared_ptr<SetIntegrationPtsMethodData> data_ptr, int lower_rank,
426 int upper_rank, BitRefLevel bit, BitRefLevel mask,
427
428 MoFEMTypes bh, VERBOSITY_LEVELS verb) {
429
431
432 struct Op : public BaseOP {
433 Op(
434
437 const std::string finite_element,
438 boost::shared_ptr<SetIntegrationPtsMethodData> data_ptr, int lower_rank,
439 int upper_rank, BitRefLevel bit, BitRefLevel mask, double eps,
441
442 )
443 : BaseOP(NOSPACE, BaseOP::OPSPACE), vecDataPts(vec_data_pts),
444 matDataPts(mat_data_pts), finiteElement(finite_element),
445 dataPtr(data_ptr), lowerRank(lower_rank), upperRank(upper_rank),
446 bitRefLevel(bit), maskRefLevel(mask), eps(eps), bhType(bh),
447 verbosity(verb) {}
448
449 MoFEMErrorCode doWork(int side, EntityType type,
452
453 MOFEM_LOG_CHANNEL("SELF");
454 MOFEM_LOG_TAG("SELF", "FieldEvaluator");
455
456 auto &m_field = this->getPtrFE()->mField;
457 auto field_eval_ptr = m_field.getInterface<FieldEvaluatorInterface>();
458
459 auto nb_integration_pts = getGaussPts().size2();
460 auto &spatial_coords = getCoordsAtGaussPts();
461
462 auto prb_name = getPtrFE()->problemPtr->getName();
463 auto cache_ptr = getPtrFE()->getCacheWeakPtr();
464
465 auto coords_at_pts = getCoordsAtGaussPts();
466 for (int gg = 0; gg != nb_integration_pts; ++gg) {
467 double *coords_ptr = &coords_at_pts(gg, 0);
468
469 dataPtr->setEvalPoints(coords_ptr, 1);
470 CHKERR field_eval_ptr->evalFEAtThePoint<D>(
471 coords_ptr, eps, prb_name, finiteElement, dataPtr, lowerRank,
472 upperRank, bitRefLevel, maskRefLevel, cache_ptr.lock(), MF_EXIST,
473 QUIET);
474
475 // Get number of found finite elements, that consist gauss point.
476 // For edge flip expected is that only one element will be found.
477 auto loop_size = dataPtr->feMethodPtr.lock()->getLoopSize();
478 if (loop_size != 1) {
479 if (verbosity >= QUIET) {
481 LogManager::BitLineID | LogManager::BitScope);
482 MOFEM_LOG("SELF", Sev::warning)
483 << "Field evaluator wrong number of elements found %d",
484 loop_size;
485 }
486 }
487
488 for (auto &v : vecDataPts) {
489 v.second->resize(nb_integration_pts, false);
490 if (gg == 0)
491 v.second->clear();
492 for (auto ll = 0; ll != loop_size; ++ll)
493 (*v.second)(gg) += (*v.first)(ll) / loop_size;
494 }
495
496 for (auto &m : matDataPts) {
497 auto nb_rows = m.first->size1();
498 m.second->resize(nb_rows, nb_integration_pts, false);
499 if (gg == 0)
500 m.second->clear();
501 for (auto ll = 0; ll != loop_size; ++ll) {
502 for (int dd = 0; dd != nb_rows; ++dd) {
503 (*m.second)(dd, gg) += (*m.first)(dd, ll) / loop_size;
504 }
505 }
506 }
507 }
508
509 MOFEM_LOG_CHANNEL("SELF");
510
512 }
513
514 private:
517 std::string finiteElement;
518 boost::shared_ptr<SetIntegrationPtsMethodData> dataPtr;
519 int lowerRank;
520 int upperRank;
521 BitRefLevel bitRefLevel;
522 BitRefLevel maskRefLevel;
523 double eps;
524 MoFEMTypes bhType;
525 VERBOSITY_LEVELS verbosity;
526 };
527
528 double min[3];
529 double max[3];
530 unsigned int dep;
531 CHKERR data_ptr->treePtr->get_info(data_ptr->rooTreeSet, min, max, dep);
532 double dist = std::sqrt(pow(max[0] - min[0], 2) + pow(max[1] - min[1], 2) +
533 pow(max[2] - min[2], 2));
534
535 return new Op(vec_data_pts, mat_data_pts, finite_element, data_ptr,
536 lower_rank, upper_rank, bit, mask, 1e-6 * dist, bh, verb);
537}
538
541 boost::shared_ptr<MoFEM::ForcesAndSourcesCore> fe_method_ptr,
542 const double *eval_points, const int nb_eval_points, const double eps,
543 VERBOSITY_LEVELS verb)
544 : feMethodPtr(fe_method_ptr), evalPoints(eval_points),
545 nbEvalPoints(nb_eval_points), eps(eps), verb(verb) {
546 localCoords.resize(nbEvalPoints, 3);
547 shapeFunctions.resize(nbEvalPoints, 4);
548}
549
551 const double *ptr, const int nb_eval_points) {
552 evalPoints = ptr;
553 nbEvalPoints = nb_eval_points;
554 localCoords.resize(nbEvalPoints, 3, false);
555 shapeFunctions.resize(nbEvalPoints, 4, false);
556}
557
558template <>
560FieldEvaluatorInterface::getDataOperator<2>(
563 const std::string finite_element,
564 boost::shared_ptr<SetIntegrationPtsMethodData> data_ptr, int lower_rank,
565 int upper_rank, BitRefLevel bit, BitRefLevel mask, MoFEMTypes bh,
566 VERBOSITY_LEVELS verb) {
567 return getDataOperatorImpl<2>(vec_data_pts, mat_data_pts, finite_element,
568 data_ptr, lower_rank, upper_rank, bit, mask, bh,
569 verb);
570}
571
572template <>
574FieldEvaluatorInterface::getDataOperator<3>(
577 const std::string finite_element,
578 boost::shared_ptr<SetIntegrationPtsMethodData> data_ptr, int lower_rank,
579 int upper_rank, BitRefLevel bit, BitRefLevel mask, MoFEMTypes bh,
580 VERBOSITY_LEVELS verb) {
581 return getDataOperatorImpl<3>(vec_data_pts, mat_data_pts, finite_element,
582 data_ptr, lower_rank, upper_rank, bit, mask, bh,
583 verb);
584}
585
586} // namespace MoFEM
#define MOFEM_LOG_SYNCHRONISE(comm)
Synchronise "SYNC" channel.
#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 v
phase velocity of light in medium (cm/ns)
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 MPI_Comm & get_comm() const =0
virtual int get_comm_rank() const =0
Core (interface) class.
Definition Core.hpp:82
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
ForcesAndSourcesCore::UserDataOperator * getDataOperatorImpl(VecDataPts vec_data_pts, 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)
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)
std::vector< std::pair< boost::shared_ptr< VectorDouble >, boost::shared_ptr< VectorDouble > > > VecDataPts
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.