v0.15.0
Loading...
Searching...
No Matches
ForcesAndSourcesCore.cpp
Go to the documentation of this file.
1/** \file ForcesAndSourcesCore.cpp
2
3\brief Implementation of Elements on Entities for Forces and Sources
4*/
5
6#ifdef __cplusplus
7extern "C" {
8#endif
9#include <cblas.h>
10#include <lapack_wrap.h>
11// #include <gm_rule.h>
12#ifdef __cplusplus
13}
14#endif
15
16namespace MoFEM {
17
18static auto cmp_uid_lo(const boost::weak_ptr<FieldEntity> &a, const UId &b) {
19 if (auto a_ptr = a.lock()) {
20 if (a_ptr->getLocalUniqueId() < b)
21 return true;
22 else
23 return false;
24 } else {
25 return false;
26 }
27}
28
29static auto cmp_uid_hi(const UId &b, const boost::weak_ptr<FieldEntity> &a) {
30 if (auto a_ptr = a.lock()) {
31 if (b < a_ptr->getLocalUniqueId())
32 return true;
33 else
34 return false;
35 } else {
36 return true;
37 }
38}
39
41 :
42
43 mField(m_field), getRuleHook(0), setRuleHook(0),
44 dataOnElement{
45
46 boost::make_shared<EntitiesFieldData>(MBENTITYSET), // NOSPACE,
47 boost::make_shared<EntitiesFieldData>(MBENTITYSET), // NOFIELD
48 boost::make_shared<EntitiesFieldData>(MBENTITYSET), // H1
49 boost::make_shared<EntitiesFieldData>(MBENTITYSET), // HCURL
50 boost::make_shared<EntitiesFieldData>(MBENTITYSET), // HDIV
51 boost::make_shared<EntitiesFieldData>(MBENTITYSET) // L2
52
53 },
54 derivedDataOnElement{
55
56 nullptr,
57 boost::make_shared<DerivedEntitiesFieldData>(
58 dataOnElement[NOFIELD]), // NOFIELD
59 boost::make_shared<DerivedEntitiesFieldData>(dataOnElement[H1]), // H1
60 boost::make_shared<DerivedEntitiesFieldData>(
61 dataOnElement[HCURL]), // HCURL
62 boost::make_shared<DerivedEntitiesFieldData>(
63 dataOnElement[HDIV]), // HDIV
64 boost::make_shared<DerivedEntitiesFieldData>(dataOnElement[L2]) // L2
65
66 },
67 dataNoField(*dataOnElement[NOFIELD].get()),
68 dataH1(*dataOnElement[H1].get()), dataHcurl(*dataOnElement[HCURL].get()),
69 dataHdiv(*dataOnElement[HDIV].get()), dataL2(*dataOnElement[L2].get()),
70 lastEvaluatedElementEntityType(MBMAXTYPE), sidePtrFE(nullptr),
71 refinePtrFE(nullptr) {
72
73 dataOnElement[NOSPACE]->dataOnEntities[MBENTITYSET].push_back(
75
76 dataOnElement[NOFIELD]->dataOnEntities[MBENTITYSET].push_back(
78 derivedDataOnElement[NOFIELD]->dataOnEntities[MBENTITYSET].push_back(
80}
81
82// ** Sense **
83
85 const EntityType type,
86 boost::ptr_vector<EntitiesFieldData::EntData> &data) const {
88
89 auto &side_table = numeredEntFiniteElementPtr->getSideNumberTable().get<0>();
90 auto sit = side_table.lower_bound(get_id_for_min_type(type));
91 if (sit != side_table.end()) {
92 auto hi_sit = side_table.upper_bound(get_id_for_max_type(type));
93 for (; sit != hi_sit; ++sit) {
94 if (const auto side_number = (*sit)->side_number; side_number >= 0) {
95 const int brother_side_number = (*sit)->brother_side_number;
96 const int sense = (*sit)->sense;
97
98 data[side_number].getSense() = sense;
99 if (brother_side_number != -1)
100 data[brother_side_number].getSense() = sense;
101 }
102 }
103 }
105}
106
107// ** Order **
108
109template <typename ENTMULTIINDEX>
110static inline int getMaxOrder(const ENTMULTIINDEX &multi_index) {
111 int max_order = 0;
112 for (auto ent_field_weak_ptr : multi_index)
113 if (auto e = ent_field_weak_ptr.lock()) {
114 const int order = e->getMaxOrder();
115 max_order = (max_order < order) ? order : max_order;
116 }
117 return max_order;
118}
119
121 int max_order = 0;
122 for (auto e : getDataFieldEnts()) {
123 if (auto ptr = e.lock()) {
124 const int order = ptr->getMaxOrder();
125 max_order = (max_order < order) ? order : max_order;
126 }
127 }
128 return max_order;
129}
130
134
138
140 const EntityType type, const FieldSpace space,
141 boost::ptr_vector<EntitiesFieldData::EntData> &data) const {
143
144 auto set_order = [&]() {
146 auto &side_table = numeredEntFiniteElementPtr->getSideNumberTable();
147
148 for (unsigned int s = 0; s != data.size(); ++s)
149 data[s].getOrder() = 0;
150
151 const FieldEntity_vector_view &data_field_ent = getDataFieldEnts();
152
153 for (auto r = fieldsPtr->get<BitFieldId_space_mi_tag>().equal_range(space);
154 r.first != r.second; ++r.first) {
155
156 const auto field_bit_number = (*r.first)->getBitNumber();
157 const auto lo_uid = FieldEntity::getLocalUniqueIdCalculate(
158 field_bit_number, get_id_for_min_type(type));
159 auto lo = std::lower_bound(data_field_ent.begin(), data_field_ent.end(),
160 lo_uid, cmp_uid_lo);
161 if (lo != data_field_ent.end()) {
162 const auto hi_uid = FieldEntity::getLocalUniqueIdCalculate(
163 field_bit_number, get_id_for_max_type(type));
164 auto hi =
165 std::upper_bound(lo, data_field_ent.end(), hi_uid, cmp_uid_hi);
166 for (; lo != hi; ++lo) {
167
168 if (auto ptr = lo->lock()) {
169
170 auto &e = *ptr;
171 auto sit = side_table.find(e.getEnt());
172 if (sit != side_table.end()) {
173 auto &side = *sit;
174 if (const auto side_number = side->side_number;
175 side_number >= 0) {
176 ApproximationOrder ent_order = e.getMaxOrder();
177 auto &dat = data[side_number];
178 dat.getOrder() =
179 dat.getOrder() > ent_order ? dat.getOrder() : ent_order;
180 }
181 } else
182 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
183 "Entity on side of the element not found");
184 }
185 }
186 }
187 }
188
190 };
191
192 auto set_order_on_brother = [&]() {
194 auto &side_table =
195 numeredEntFiniteElementPtr->getSideNumberTable().get<0>();
196 auto sit = side_table.lower_bound(get_id_for_min_type(type));
197 if (sit != side_table.end()) {
198 auto hi_sit = side_table.upper_bound(get_id_for_max_type(type));
199 for (; sit != hi_sit; ++sit) {
200 const int brother_side_number = (*sit)->brother_side_number;
201 if (brother_side_number != -1) {
202 const int side_number = (*sit)->side_number;
203 data[brother_side_number].getOrder() = data[side_number].getOrder();
204 }
205 }
206 }
208 };
209
210 CHKERR set_order();
211 CHKERR set_order_on_brother();
212
214}
215
216// ** Indices **
217
218template <typename EXTRACTOR>
220 const int bit_number, FieldEntity_vector_view &ents_field,
221 VectorInt &nodes_indices, VectorInt &local_nodes_indices,
222 EXTRACTOR &&extractor) const {
224
225 auto field_it = fieldsPtr->get<BitFieldId_mi_tag>().find(
226 BitFieldId().set(bit_number - 1));
227 if (field_it != fieldsPtr->get<BitFieldId_mi_tag>().end()) {
228
229#ifndef NDEBUG
230 if ((*field_it)->getBitNumber() != bit_number)
231 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "Wrong bit number");
232#endif
233 const auto lo_uid = FieldEntity::getLocalUniqueIdCalculate(
234 bit_number, get_id_for_min_type<MBVERTEX>());
235 auto lo = std::lower_bound(ents_field.begin(), ents_field.end(), lo_uid,
236 cmp_uid_lo);
237 if (lo != ents_field.end()) {
238 const auto hi_uid = FieldEntity::getLocalUniqueIdCalculate(
239 bit_number, get_id_for_max_type<MBVERTEX>());
240 auto hi = std::upper_bound(lo, ents_field.end(), hi_uid, cmp_uid_hi);
241
242 const int num_nodes = getNumberOfNodes();
243 const int nb_dofs_on_vert = (*field_it)->getNbOfCoeffs();
244 const int max_nb_dofs = nb_dofs_on_vert * num_nodes;
245
246 int nb_dofs = 0;
247 for (auto it = lo; it != hi; ++it) {
248 if (auto e = it->lock()) {
249 if (auto cache = extractor(e).lock()) {
250 if (cache->loHi[0] != cache->loHi[1]) {
251 nb_dofs += std::distance(cache->loHi[0], cache->loHi[1]);
252 break;
253 }
254 }
255 }
256 }
257
258 if (nb_dofs) {
259 nodes_indices.resize(max_nb_dofs, false);
260 local_nodes_indices.resize(max_nb_dofs, false);
261 } else {
262 nodes_indices.resize(0, false);
263 local_nodes_indices.resize(0, false);
264 }
265
266 if (nb_dofs != max_nb_dofs) {
267 std::fill(nodes_indices.begin(), nodes_indices.end(), -1);
268 std::fill(local_nodes_indices.begin(), local_nodes_indices.end(), -1);
269 }
270
271 for (auto it = lo; it != hi; ++it) {
272 if (auto e = it->lock()) {
273 auto side_ptr = e->getSideNumberPtr();
274 if (const auto side_number = side_ptr->side_number;
275 side_number >= 0) {
276 const auto brother_side_number = side_ptr->brother_side_number;
277 if (auto cache = extractor(e).lock()) {
278 for (auto dit = cache->loHi[0]; dit != cache->loHi[1]; ++dit) {
279 auto &dof = **dit;
280 const int idx = dof.getPetscGlobalDofIdx();
281 const int local_idx = dof.getPetscLocalDofIdx();
282 const int pos =
283 side_number * nb_dofs_on_vert + dof.getDofCoeffIdx();
284 nodes_indices[pos] = idx;
285 local_nodes_indices[pos] = local_idx;
286 if (brother_side_number != -1) {
287 const int elem_idx = brother_side_number * nb_dofs_on_vert +
288 (*dit)->getDofCoeffIdx();
289 nodes_indices[elem_idx] = idx;
290 local_nodes_indices[elem_idx] = local_idx;
291 }
292 }
293 }
294 }
295 }
296 }
297 } else {
298 nodes_indices.resize(0, false);
299 local_nodes_indices.resize(0, false);
300 }
301
302 } else {
303 nodes_indices.resize(0, false);
304 local_nodes_indices.resize(0, false);
305 }
306
308}
309
312 const int bit_number) const {
313
314 struct Extractor {
315 boost::weak_ptr<EntityCacheNumeredDofs>
316 operator()(boost::shared_ptr<FieldEntity> &e) {
317 return e->entityCacheRowDofs;
318 }
319 };
320
321 return getNodesIndices(bit_number, getRowFieldEnts(),
322 data.dataOnEntities[MBVERTEX][0].getIndices(),
323 data.dataOnEntities[MBVERTEX][0].getLocalIndices(),
324 Extractor());
325}
326
329 const int bit_number) const {
330
331 struct Extractor {
332 boost::weak_ptr<EntityCacheNumeredDofs>
333 operator()(boost::shared_ptr<FieldEntity> &e) {
334 return e->entityCacheColDofs;
335 }
336 };
337
338 return getNodesIndices(bit_number, getColFieldEnts(),
339 data.dataOnEntities[MBVERTEX][0].getIndices(),
340 data.dataOnEntities[MBVERTEX][0].getLocalIndices(),
341 Extractor());
342}
343
344template <typename EXTRACTOR>
346 EntitiesFieldData &data, const int bit_number,
347 FieldEntity_vector_view &ents_field, const EntityType type_lo,
348 const EntityType type_hi, EXTRACTOR &&extractor) const {
350
351 for (EntityType t = type_lo; t != type_hi; ++t) {
352 for (auto &dat : data.dataOnEntities[t]) {
353 dat.getIndices().resize(0, false);
354 dat.getLocalIndices().resize(0, false);
355 }
356 }
357
358 const auto lo_uid = FieldEntity::getLocalUniqueIdCalculate(
359 bit_number, get_id_for_min_type(type_lo));
360 auto lo = std::lower_bound(ents_field.begin(), ents_field.end(), lo_uid,
361 cmp_uid_lo);
362 if (lo != ents_field.end()) {
363 const auto hi_uid = FieldEntity::getLocalUniqueIdCalculate(
364 bit_number, get_id_for_max_type(type_hi));
365 auto hi = std::upper_bound(lo, ents_field.end(), hi_uid, cmp_uid_hi);
366 if (lo != hi) {
367
368 std::vector<boost::weak_ptr<FieldEntity>> brother_ents_vec;
369
370 for (auto it = lo; it != hi; ++it)
371 if (auto e = it->lock()) {
372
373 const EntityType type = e->getEntType();
374 auto side_ptr = e->getSideNumberPtr();
375 if (const auto side = side_ptr->side_number; side >= 0) {
376 const auto nb_dofs_on_ent = e->getNbDofsOnEnt();
377 const auto brother_side = side_ptr->brother_side_number;
378 auto &dat = data.dataOnEntities[type][side];
379 auto &ent_field_indices = dat.getIndices();
380 auto &ent_field_local_indices = dat.getLocalIndices();
381
382 ent_field_indices.resize(nb_dofs_on_ent, false);
383 ent_field_local_indices.resize(nb_dofs_on_ent, false);
384 std::fill(ent_field_indices.data().begin(),
385 ent_field_indices.data().end(), -1);
386 std::fill(ent_field_local_indices.data().begin(),
387 ent_field_local_indices.data().end(), -1);
388
389 if (auto cache = extractor(e).lock()) {
390 for (auto dit = cache->loHi[0]; dit != cache->loHi[1]; ++dit) {
391 const int idx = (*dit)->getEntDofIdx();
392 ent_field_indices[idx] = (*dit)->getPetscGlobalDofIdx();
393 ent_field_local_indices[idx] = (*dit)->getPetscLocalDofIdx();
394 }
395 }
396
397 if (brother_side != -1) {
398 auto &dat_brother = data.dataOnEntities[type][brother_side];
399 dat_brother.getIndices().resize(nb_dofs_on_ent, false);
400 dat_brother.getLocalIndices().resize(nb_dofs_on_ent, false);
401 noalias(dat_brother.getIndices()) = dat.getIndices();
402 noalias(dat_brother.getLocalIndices()) = dat.getLocalIndices();
403 }
404 }
405 }
406 }
407 }
408
410}
411
413 EntitiesFieldData &data, const int bit_number, const EntityType type_lo,
414 const EntityType type_hi) const {
415
416 struct Extractor {
417 boost::weak_ptr<EntityCacheNumeredDofs>
418 operator()(boost::shared_ptr<FieldEntity> &e) {
419 return e->entityCacheRowDofs;
420 }
421 };
422
423 return getEntityIndices(data, bit_number, getRowFieldEnts(), type_lo, type_hi,
424 Extractor());
425}
426
428 EntitiesFieldData &data, const int bit_number, const EntityType type_lo,
429 const EntityType type_hi) const {
430
431 struct Extractor {
432 boost::weak_ptr<EntityCacheNumeredDofs>
433 operator()(boost::shared_ptr<FieldEntity> &e) {
434 return e->entityCacheColDofs;
435 }
436 };
437
438 return getEntityIndices(data, bit_number, getColFieldEnts(), type_lo, type_hi,
439 Extractor());
440}
441
443 const int bit_number, boost::shared_ptr<FENumeredDofEntity_multiIndex> dofs,
444 VectorInt &indices) const {
446 auto dit = dofs->get<Unique_mi_tag>().lower_bound(
448 auto hi_dit = dofs->get<Unique_mi_tag>().upper_bound(
450 indices.resize(std::distance(dit, hi_dit));
451 for (; dit != hi_dit; dit++) {
452 int idx = (*dit)->getPetscGlobalDofIdx();
453 indices[(*dit)->getDofCoeffIdx()] = idx;
454 }
456}
457
460 const int bit_number) const {
462#ifndef NDEBUG
463 if (data.dataOnEntities[MBENTITYSET].size() == 0) {
464 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
465 "data.dataOnEntities[MBENTITYSET] is empty");
466 }
467#endif
469 data.dataOnEntities[MBENTITYSET][0].getIndices());
471}
472
475 const int bit_number) const {
477#ifndef NDEBUG
478 if (data.dataOnEntities[MBENTITYSET].size() == 0) {
479 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
480 "data.dataOnEntities[MBENTITYSET] is empty");
481 }
482#endif
484 data.dataOnEntities[MBENTITYSET][0].getIndices());
486}
487
488// ** Indices from problem **
489
491 const std::string &field_name, const NumeredDofEntity_multiIndex &dofs,
492 VectorInt &nodes_indices) const {
494
495 const Field *field_struture = mField.get_field_structure(field_name);
496 if (field_struture->getSpace() == H1) {
497
498 const int num_nodes = getNumberOfNodes();
499 nodes_indices.resize(field_struture->getNbOfCoeffs() * num_nodes, false);
500 std::fill(nodes_indices.begin(), nodes_indices.end(), -1);
501
502 auto &side_table = const_cast<SideNumber_multiIndex &>(
503 numeredEntFiniteElementPtr->getSideNumberTable());
504 auto siit = side_table.get<1>().lower_bound(boost::make_tuple(MBVERTEX, 0));
505 auto hi_siit =
506 side_table.get<1>().lower_bound(boost::make_tuple(MBVERTEX, 10000));
507
508 int nn = 0;
509 for (; siit != hi_siit; siit++, nn++) {
510 if (siit->get()->side_number >= 0) {
511 auto bit_number = mField.get_field_bit_number(field_name);
512 const EntityHandle ent = siit->get()->ent;
513 auto dit = dofs.get<Unique_mi_tag>().lower_bound(
515 auto hi_dit = dofs.get<Unique_mi_tag>().upper_bound(
517 for (; dit != hi_dit; dit++) {
518 nodes_indices[siit->get()->side_number * (*dit)->getNbOfCoeffs() +
519 (*dit)->getDofCoeffIdx()] =
520 (*dit)->getPetscGlobalDofIdx();
521 }
522 }
523 }
524 } else {
525 nodes_indices.resize(0, false);
526 }
527
529}
530
532 const std::string &field_name, const NumeredDofEntity_multiIndex &dofs,
533 EntityType type, int side_number, VectorInt &indices) const {
535
536 indices.resize(0);
537
538 auto &side_table = const_cast<SideNumber_multiIndex &>(
539 numeredEntFiniteElementPtr->getSideNumberTable());
540 auto siit =
541 side_table.get<1>().lower_bound(boost::make_tuple(type, side_number));
542 auto hi_siit =
543 side_table.get<1>().upper_bound(boost::make_tuple(type, side_number));
544
545 for (; siit != hi_siit; siit++) {
546 if (siit->get()->side_number >= 0) {
547
548 const EntityHandle ent = siit->get()->ent;
549 auto bit_number = mField.get_field_bit_number(field_name);
550 auto dit = dofs.get<Unique_mi_tag>().lower_bound(
552 auto hi_dit = dofs.get<Unique_mi_tag>().upper_bound(
554 indices.resize(std::distance(dit, hi_dit));
555 for (; dit != hi_dit; dit++) {
556
557 indices[(*dit)->getEntDofIdx()] = (*dit)->getPetscGlobalDofIdx();
558 }
559 }
560 }
561
563}
564
566 const std::string &field_name, VectorInt &nodes_indices) const {
568 nodes_indices);
569}
570
573 EntityType type, int side_number,
574 VectorInt &indices) const {
576 type, side_number, indices);
577}
578
580 const std::string &field_name, VectorInt &nodes_indices) const {
582 nodes_indices);
583}
584
587 EntityType type, int side_number,
588 VectorInt &indices) const {
590 type, side_number, indices);
591}
592
593// ** Data **
594
597
598 for (auto &data : dataOnElement) {
599 if (data) {
600 for (auto &dat : data->dataOnEntities) {
601 for (auto &ent_dat : dat) {
602 ent_dat.getEntDataBitRefLevel().clear();
603 }
604 }
605 }
606 }
607
608 auto &field_ents = getDataFieldEnts();
609 for (auto it : field_ents) {
610 if (auto e = it.lock()) {
611 const FieldSpace space = e->getSpace();
612 if (space > NOFIELD) {
613 const EntityType type = e->getEntType();
614 const signed char side = e->getSideNumberPtr()->side_number;
615 if (side >= 0) {
616 if (auto &data = dataOnElement[space]) {
617 if (type == MBVERTEX) {
618 auto &dat = data->dataOnEntities[type][0];
619 dat.getEntDataBitRefLevel().resize(getNumberOfNodes(), false);
620 dat.getEntDataBitRefLevel()[side] = e->getBitRefLevel();
621 } else {
622 auto &dat = data->dataOnEntities[type][side];
623 dat.getEntDataBitRefLevel().resize(1, false);
624 dat.getEntDataBitRefLevel()[0] = e->getBitRefLevel();
625 }
626 }
627 }
628 } else {
629 if (auto &data = dataOnElement[NOFIELD]) {
630 auto &dat = data->dataOnEntities[MBENTITYSET][0];
631 dat.getEntDataBitRefLevel().resize(1, false);
632 dat.getEntDataBitRefLevel()[0] = e->getBitRefLevel();
633 }
634 }
635 }
636 }
637
639};
640
643 const int bit_number) const {
644
645 auto get_nodes_field_data = [&](VectorDouble &nodes_data,
646 VectorFieldEntities &field_entities,
647 VectorDofs &nodes_dofs, FieldSpace &space,
649 VectorInt &bb_node_order) {
651
652 nodes_data.resize(0, false);
653 nodes_dofs.resize(0, false);
654 field_entities.resize(0, false);
655
656 auto field_it =
657 fieldsPtr->get<BitFieldId_mi_tag>().find(BitFEId().set(bit_number - 1));
658 if (field_it != fieldsPtr->get<BitFieldId_mi_tag>().end()) {
659
660#ifndef NDEBUG
661 if ((*field_it)->getBitNumber() != bit_number)
662 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "Wrong bit number");
663#endif
664 const int nb_dofs_on_vert = (*field_it)->getNbOfCoeffs();
665 space = (*field_it)->getSpace();
666 base = (*field_it)->getApproxBase();
667
668 auto &field_ents = getDataFieldEnts();
669 const auto lo_uid = FieldEntity::getLocalUniqueIdCalculate(
670 bit_number, get_id_for_min_type<MBVERTEX>());
671 auto lo = std::lower_bound(field_ents.begin(), field_ents.end(), lo_uid,
672 cmp_uid_lo);
673 if (lo != field_ents.end()) {
674 const auto hi_uid = FieldEntity::getLocalUniqueIdCalculate(
675 bit_number, get_id_for_max_type<MBVERTEX>());
676 auto hi = std::upper_bound(lo, field_ents.end(), hi_uid, cmp_uid_hi);
677 if (lo != hi) {
678
679 int nb_dofs = 0;
680 for (auto it = lo; it != hi; ++it) {
681 if (auto e = it->lock()) {
682 if (auto cache = e->entityCacheDataDofs.lock()) {
683 if (cache->loHi[0] != cache->loHi[1]) {
684 nb_dofs += std::distance(cache->loHi[0], cache->loHi[1]);
685 break;
686 }
687 }
688 }
689 }
690
691 if (nb_dofs) {
692
693 const int num_nodes = getNumberOfNodes();
694 bb_node_order.resize(num_nodes, false);
695 bb_node_order.clear();
696 const int max_nb_dofs = nb_dofs_on_vert * num_nodes;
697 nodes_data.resize(max_nb_dofs, false);
698 nodes_dofs.resize(max_nb_dofs, false);
699 field_entities.resize(num_nodes, false);
700 std::fill(nodes_data.begin(), nodes_data.end(), 0);
701 std::fill(nodes_dofs.begin(), nodes_dofs.end(), nullptr);
702 std::fill(field_entities.begin(), field_entities.end(), nullptr);
703
704 std::vector<boost::weak_ptr<FieldEntity>> brother_ents_vec;
705
706 for (auto it = lo; it != hi; ++it) {
707 if (auto e = it->lock()) {
708 const auto &sn = e->getSideNumberPtr();
709 // Some field entities on skeleton can have negative side
710 // number
711 if (const auto side_number = sn->side_number;
712 side_number >= 0) {
713 const int brother_side_number = sn->brother_side_number;
714
715 field_entities[side_number] = e.get();
716 if (brother_side_number != -1) {
717 brother_ents_vec.emplace_back(e);
718 field_entities[side_number] = field_entities[side_number];
719 }
720
721 bb_node_order[side_number] = e->getMaxOrder();
722 int pos = side_number * nb_dofs_on_vert;
723 auto ent_filed_data_vec = e->getEntFieldData();
724 if (auto cache = e->entityCacheDataDofs.lock()) {
725 for (auto dit = cache->loHi[0]; dit != cache->loHi[1];
726 ++dit) {
727 const auto dof_idx = (*dit)->getEntDofIdx();
728 nodes_data[pos + dof_idx] = ent_filed_data_vec[dof_idx];
729 nodes_dofs[pos + dof_idx] =
730 reinterpret_cast<FEDofEntity *>((*dit).get());
731 }
732 }
733 }
734 }
735 }
736
737 for (auto &it : brother_ents_vec) {
738 if (const auto e = it.lock()) {
739 const auto &sn = e->getSideNumberPtr();
740 const int side_number = sn->side_number;
741 const int brother_side_number = sn->brother_side_number;
742 bb_node_order[brother_side_number] = bb_node_order[side_number];
743 int pos = side_number * nb_dofs_on_vert;
744 int brother_pos = brother_side_number * nb_dofs_on_vert;
745 for (int ii = 0; ii != nb_dofs_on_vert; ++ii) {
746 nodes_data[brother_pos] = nodes_data[pos];
747 nodes_dofs[brother_pos] = nodes_dofs[pos];
748 ++pos;
749 ++brother_pos;
750 }
751 }
752 }
753 }
754 }
755 }
756 }
757
759 };
760
761 return get_nodes_field_data(
762 data.dataOnEntities[MBVERTEX][0].getFieldData(),
763 data.dataOnEntities[MBVERTEX][0].getFieldEntities(),
764 data.dataOnEntities[MBVERTEX][0].getFieldDofs(),
765 data.dataOnEntities[MBVERTEX][0].getSpace(),
766 data.dataOnEntities[MBVERTEX][0].getBase(),
767 data.dataOnEntities[MBVERTEX][0].getBBNodeOrder());
768}
769
771 EntitiesFieldData &data, const int bit_number, const EntityType type_lo,
772 const EntityType type_hi) const {
774 for (EntityType t = type_lo; t != type_hi; ++t) {
775 for (auto &dat : data.dataOnEntities[t]) {
776 dat.getOrder() = 0;
777 dat.getBase() = NOBASE;
778 dat.getSpace() = NOSPACE;
779 dat.getFieldData().resize(0, false);
780 dat.getFieldDofs().resize(0, false);
781 dat.getFieldEntities().resize(0, false);
782 }
783 }
784
785 auto &field_ents = getDataFieldEnts();
786 const auto lo_uid = FieldEntity::getLocalUniqueIdCalculate(
787 bit_number, get_id_for_min_type(type_lo));
788 auto lo = std::lower_bound(field_ents.begin(), field_ents.end(), lo_uid,
789 cmp_uid_lo);
790 if (lo != field_ents.end()) {
791 const auto hi_uid = FieldEntity::getLocalUniqueIdCalculate(
792 bit_number, get_id_for_max_type(type_hi));
793 auto hi = std::upper_bound(lo, field_ents.end(), hi_uid, cmp_uid_hi);
794 if (lo != hi) {
795
796 std::vector<boost::weak_ptr<FieldEntity>> brother_ents_vec;
797
798 for (auto it = lo; it != hi; ++it)
799 if (auto e = it->lock()) {
800 auto side_ptr = e->getSideNumberPtr();
801 if (const auto side = side_ptr->side_number; side >= 0) {
802 const EntityType type = e->getEntType();
803 auto &dat = data.dataOnEntities[type][side];
804 auto &ent_field = dat.getFieldEntities();
805 auto &ent_field_dofs = dat.getFieldDofs();
806 auto &ent_field_data = dat.getFieldData();
807
808 const int brother_side = side_ptr->brother_side_number;
809 if (brother_side != -1)
810 brother_ents_vec.emplace_back(e);
811
812 dat.getBase() = e->getApproxBase();
813 dat.getSpace() = e->getSpace();
814 const int ent_order = e->getMaxOrder();
815 dat.getOrder() =
816 dat.getOrder() > ent_order ? dat.getOrder() : ent_order;
817
818 auto ent_data = e->getEntFieldData();
819 ent_field_data.resize(ent_data.size(), false);
820 noalias(ent_field_data) = ent_data;
821 ent_field_dofs.resize(ent_data.size(), false);
822 std::fill(ent_field_dofs.begin(), ent_field_dofs.end(), nullptr);
823 ent_field.resize(1, false);
824 ent_field[0] = e.get();
825 if (auto cache = e->entityCacheDataDofs.lock()) {
826 for (auto dit = cache->loHi[0]; dit != cache->loHi[1]; ++dit) {
827 ent_field_dofs[(*dit)->getEntDofIdx()] =
828 reinterpret_cast<FEDofEntity *>((*dit).get());
829 }
830 }
831 }
832 }
833
834 for (auto &it : brother_ents_vec) {
835 if (const auto e = it.lock()) {
836 const EntityType type = e->getEntType();
837 const int side = e->getSideNumberPtr()->side_number;
838 const int brother_side = e->getSideNumberPtr()->brother_side_number;
839 auto &dat = data.dataOnEntities[type][side];
840 auto &dat_brother = data.dataOnEntities[type][brother_side];
841 dat_brother.getBase() = dat.getBase();
842 dat_brother.getSpace() = dat.getSpace();
843 dat_brother.getOrder() = dat.getOrder();
844 dat_brother.getFieldData() = dat.getFieldData();
845 dat_brother.getFieldDofs() = dat.getFieldDofs();
846 dat_brother.getFieldEntities() = dat.getFieldEntities();
847 }
848 }
849 }
850 }
851
853}
854
856 const int bit_number, VectorDouble &ent_field_data,
857 VectorDofs &ent_field_dofs, VectorFieldEntities &ent_field) const {
858
860
861 ent_field_data.resize(0, false);
862 ent_field_dofs.resize(0, false);
863 ent_field.resize(0, false);
864
865 auto &field_ents = getDataFieldEnts();
866 const auto lo_uid = FieldEntity::getLocalUniqueIdCalculate(
867 bit_number, get_id_for_min_type<MBVERTEX>());
868 auto lo = std::lower_bound(field_ents.begin(), field_ents.end(), lo_uid,
869 cmp_uid_lo);
870 if (lo != field_ents.end()) {
871
872 ent_field.resize(field_ents.size(), false);
873 std::fill(ent_field.begin(), ent_field.end(), nullptr);
874
875 const auto hi_uid = FieldEntity::getLocalUniqueIdCalculate(
876 bit_number, get_id_for_max_type<MBENTITYSET>());
877 auto hi = std::upper_bound(lo, field_ents.end(), hi_uid, cmp_uid_hi);
878 if (lo != hi) {
879
880 int side = 0;
881 for (auto it = lo; it != hi; ++it, ++side)
882 if (auto e = it->lock()) {
883
884 const auto size = e->getNbDofsOnEnt();
885 ent_field_data.resize(size, false);
886 ent_field_dofs.resize(size, false);
887 ent_field[side] = e.get();
888 noalias(ent_field_data) = e->getEntFieldData();
889
890 if (auto cache = e->entityCacheDataDofs.lock()) {
891 for (auto dit = cache->loHi[0]; dit != cache->loHi[1]; ++dit) {
892 ent_field_dofs[(*dit)->getEntDofIdx()] =
893 reinterpret_cast<FEDofEntity *>((*dit).get());
894 }
895 }
896 }
897 }
898 }
899
901}
902
905 const int bit_number) const {
907 if (data.dataOnEntities[MBENTITYSET].size() == 0)
908 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
909 "No space to insert data");
910
912 bit_number, data.dataOnEntities[MBENTITYSET][0].getFieldData(),
913 data.dataOnEntities[MBENTITYSET][0].getFieldDofs(),
914 data.dataOnEntities[MBENTITYSET][0].getFieldEntities());
916}
917
918// ** Face **
919
923 auto &face_nodes = data.facesNodes;
924 auto &face_nodes_order = data.facesNodesOrder;
925 auto &side_table = const_cast<SideNumber_multiIndex &>(
926 numeredEntFiniteElementPtr->getSideNumberTable());
927 const auto ent = numeredEntFiniteElementPtr->getEnt();
928 const auto type = numeredEntFiniteElementPtr->getEntType();
929 const auto nb_faces = CN::NumSubEntities(type, 2);
930 const EntityHandle *conn_ele;
931 int num_nodes_ele;
932 CHKERR mField.get_moab().get_connectivity(ent, conn_ele, num_nodes_ele, true);
933 auto side_ptr_it = side_table.get<1>().lower_bound(
934 boost::make_tuple(CN::TypeDimensionMap[2].first, 0));
935 auto hi_side_ptr_it = side_table.get<1>().upper_bound(
936 boost::make_tuple(CN::TypeDimensionMap[2].second, 100));
937
938 for (; side_ptr_it != hi_side_ptr_it; ++side_ptr_it) {
939 const auto side = (*side_ptr_it)->side_number;
940 const auto sense = (*side_ptr_it)->sense;
941 const auto offset = (*side_ptr_it)->offset;
942
943 EntityType face_type;
944 int nb_nodes_face;
945 auto face_indices =
946 CN::SubEntityVertexIndices(type, 2, side, face_type, nb_nodes_face);
947 face_nodes.resize(nb_faces, nb_nodes_face);
948 face_nodes_order.resize(nb_faces, nb_nodes_face);
949
950 if (sense == 1)
951 for (int n = 0; n != nb_nodes_face; ++n)
952 face_nodes_order(side, n) = (n + offset) % nb_nodes_face;
953 else
954 for (int n = 0; n != nb_nodes_face; ++n)
955 face_nodes_order(side, n) =
956 (nb_nodes_face - (n - offset) % nb_nodes_face) % nb_nodes_face;
957
958 for (int n = 0; n != nb_nodes_face; ++n)
959 face_nodes(side, n) = face_indices[face_nodes_order(side, n)];
960
961#ifndef NDEBUG
962 const auto face_ent = (*side_ptr_it)->ent;
963 auto check = [&]() {
965 const EntityHandle *conn_face;
966 // int nb_nodes_face;
967 CHKERR mField.get_moab().get_connectivity(face_ent, conn_face,
968 nb_nodes_face, true);
969 face_nodes.resize(nb_faces, nb_nodes_face);
970 for (int nn = 0; nn != nb_nodes_face; ++nn) {
971 if (face_nodes(side, nn) !=
972 std::distance(
973 conn_ele,
974 std::find(conn_ele, &conn_ele[num_nodes_ele], conn_face[nn]))) {
975 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
976 "Wrong face numeration");
977 }
978 }
980 };
981 CHKERR check();
982#endif
983 }
984
986}
987
988// ** Space and Base **
989
991 EntitiesFieldData &data) const {
993
994 if (nInTheLoop == 0) {
995 data.bAse.reset();
996 for (EntityType t = MBVERTEX; t != MBMAXTYPE; ++t) {
997 data.spacesOnEntities[t].reset();
998 data.basesOnEntities[t].reset();
999 }
1000 for (int s = 0; s != LASTSPACE; ++s) {
1001 data.basesOnSpaces[s].reset();
1002 data.brokenBasesOnSpaces[s].reset();
1003 }
1004 }
1005
1006 auto fe_type = numeredEntFiniteElementPtr->getEntType();
1007
1008 if (getDataFieldEntsPtr())
1009 for (auto e : getDataFieldEnts()) {
1010 if (auto ptr = e.lock()) {
1011 // get data from entity
1012 const EntityType type = ptr->getEntType();
1013 if (DefaultElementAdjacency::getDefTypeMap(fe_type, type)) {
1014 const FieldSpace space = ptr->getSpace();
1015 const FieldContinuity continuity = ptr->getContinuity();
1016 const FieldApproximationBase approx = ptr->getApproxBase();
1017
1018 // set data
1019 switch (continuity) {
1020 case CONTINUOUS:
1021 data.basesOnSpaces[space].set(approx);
1022 break;
1023 case DISCONTINUOUS:
1024 data.brokenBasesOnSpaces[space].set(approx);
1025 break;
1026 default:
1027 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1028 "Continuity not defined");
1029 }
1030
1031 data.bAse.set(approx);
1032 data.spacesOnEntities[type].set(space);
1033 data.basesOnEntities[type].set(approx);
1034 }
1035 }
1036 }
1037 else
1038 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1039 "data fields ents not allocated on element");
1040
1041 for (auto space = 0; space != LASTSPACE; ++space)
1042 if ((data.brokenBasesOnSpaces[space] & data.basesOnSpaces[space]).any())
1043 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1044 "Discontinuous and continuous bases on the same space");
1045
1047}
1048
1050 const FieldApproximationBase b) {
1052 if (dataOnElement[H1]->bAse.test(b)) {
1053
1054 switch (static_cast<FieldApproximationBase>(b)) {
1055 case NOBASE:
1057 // BERNSTEIN_BEZIER_BASE is not hierarchical base
1059 }
1060 default:
1061 break;
1062 }
1063
1064 auto get_elem_base = [&](auto base) {
1065 if (b != USER_BASE) {
1069 "Functions generating approximation base not defined");
1070 return getElementPolynomialBase();
1071 } else {
1072 if (!getUserPolynomialBase())
1075 "Functions generating user approximation base not defined");
1076 return getUserPolynomialBase();
1077 }
1078 };
1079
1080 auto get_ctx = [&](auto space, auto base, auto continuity) {
1081 return boost::make_shared<EntPolynomialBaseCtx>(
1082 *dataOnElement[space], static_cast<FieldSpace>(space), continuity,
1083 static_cast<FieldApproximationBase>(base), NOBASE);
1084 };
1085
1086 auto get_value = [&](auto elem_base, auto &gaussPts, auto &&ctx) {
1087 return elem_base->getValue(gaussPts, ctx);
1088 };
1089
1090 for (int space = H1; space != LASTSPACE; ++space) {
1091
1092#ifndef NDEBUG
1093 if (dataOnElement[H1]->basesOnSpaces[space].test(b) &&
1094 dataOnElement[H1]->brokenBasesOnSpaces[space].test(b)) {
1095 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1096 "Discontinuous and continuous bases on the same space");
1097 }
1098#endif // NDEBUG
1099
1100 if (dataOnElement[H1]->basesOnSpaces[space].test(b))
1101 CHKERR get_value(get_elem_base(b), gaussPts,
1102 get_ctx(space, b, CONTINUOUS));
1103 else if (dataOnElement[H1]->brokenBasesOnSpaces[space].test(b))
1104 CHKERR get_value(get_elem_base(b), gaussPts,
1105 get_ctx(space, b, DISCONTINUOUS));
1106 }
1107
1108 }
1110}
1111
1114 /// Use the some node base. Node base is usually used for construction other
1115 /// bases.
1116 for (int space = HCURL; space != LASTSPACE; ++space) {
1117 dataOnElement[space]->dataOnEntities[MBVERTEX][0].getNSharedPtr(NOBASE) =
1118 dataOnElement[H1]->dataOnEntities[MBVERTEX][0].getNSharedPtr(NOBASE);
1119 dataOnElement[space]->dataOnEntities[MBVERTEX][0].getDiffNSharedPtr(
1120 NOBASE) =
1121 dataOnElement[H1]->dataOnEntities[MBVERTEX][0].getDiffNSharedPtr(
1122 NOBASE);
1123 }
1124 for (int b = AINSWORTH_LEGENDRE_BASE; b != LASTBASE; b++) {
1126 static_cast<FieldApproximationBase>(b));
1127 }
1129}
1130
1134
1135 const auto ele_type = numeredEntFiniteElementPtr->getEntType();
1136
1137 auto get_nodal_base_data = [&](EntitiesFieldData &data, auto field_ptr) {
1139 auto &space = data.dataOnEntities[MBVERTEX][0].getSpace();
1140 auto &base = data.dataOnEntities[MBVERTEX][0].getBase();
1141 auto &bb_node_order = data.dataOnEntities[MBVERTEX][0].getBBNodeOrder();
1142
1143 auto &field_ents = getDataFieldEnts();
1144 auto bit_number = field_ptr->getBitNumber();
1145 const auto lo_uid = FieldEntity::getLocalUniqueIdCalculate(
1146 bit_number, get_id_for_min_type<MBVERTEX>());
1147 auto lo = std::lower_bound(field_ents.begin(), field_ents.end(), lo_uid,
1148 cmp_uid_lo);
1149 if (lo != field_ents.end()) {
1150 const auto hi_uid = FieldEntity::getLocalUniqueIdCalculate(
1151 bit_number, get_id_for_max_type<MBVERTEX>());
1152 auto hi = std::upper_bound(lo, field_ents.end(), hi_uid, cmp_uid_hi);
1153 if (lo != hi) {
1154
1155 for (auto it = lo; it != hi; ++it)
1156 if (auto first_e = it->lock()) {
1157 space = first_e->getSpace();
1158 base = first_e->getApproxBase();
1159 const int num_nodes = getNumberOfNodes();
1160 bb_node_order.resize(num_nodes, false);
1161 bb_node_order.clear();
1162
1163 std::vector<boost::weak_ptr<FieldEntity>> brother_ents_vec;
1164
1165 for (; it != hi; ++it) {
1166 if (auto e = it->lock()) {
1167 const auto &sn = e->getSideNumberPtr();
1168 const int side_number = sn->side_number;
1169 const int brother_side_number = sn->brother_side_number;
1170 if (brother_side_number != -1)
1171 brother_ents_vec.emplace_back(e);
1172 bb_node_order[side_number] = e->getMaxOrder();
1173 }
1174 }
1175
1176 for (auto &it : brother_ents_vec) {
1177 if (const auto e = it.lock()) {
1178 const auto &sn = e->getSideNumberPtr();
1179 const int side_number = sn->side_number;
1180 const int brother_side_number = sn->brother_side_number;
1181 bb_node_order[brother_side_number] = bb_node_order[side_number];
1182 }
1183 }
1184
1185 break;
1186 }
1187 }
1188 }
1189
1191 };
1192
1193 auto get_entity_base_data = [&](EntitiesFieldData &data, auto field_ptr,
1194 const EntityType type_lo,
1195 const EntityType type_hi) {
1197 for (EntityType t = MBEDGE; t != MBPOLYHEDRON; ++t) {
1198 for (auto &dat : data.dataOnEntities[t]) {
1199 dat.getOrder() = 0;
1200 dat.getBase() = NOBASE;
1201 dat.getSpace() = NOSPACE;
1202 dat.getFieldData().resize(0, false);
1203 dat.getFieldDofs().resize(0, false);
1204 }
1205 }
1206
1207 auto &field_ents = getDataFieldEnts();
1208 auto bit_number = field_ptr->getBitNumber();
1209 const auto lo_uid = FieldEntity::getLocalUniqueIdCalculate(
1210 bit_number, get_id_for_min_type(type_lo));
1211 auto lo = std::lower_bound(field_ents.begin(), field_ents.end(), lo_uid,
1212 cmp_uid_lo);
1213 if (lo != field_ents.end()) {
1214 const auto hi_uid = FieldEntity::getLocalUniqueIdCalculate(
1215 bit_number, get_id_for_max_type(type_hi));
1216 auto hi = std::upper_bound(lo, field_ents.end(), hi_uid, cmp_uid_hi);
1217
1218 std::vector<boost::weak_ptr<FieldEntity>> brother_ents_vec;
1219 for (; lo != hi; ++lo) {
1220 if (auto e = lo->lock()) {
1221 if (auto cache = e->entityCacheDataDofs.lock()) {
1222 if (cache->loHi[0] != cache->loHi[1]) {
1223 if (const auto side = e->getSideNumberPtr()->side_number;
1224 side >= 0) {
1225 const EntityType type = e->getEntType();
1226 auto &dat = data.dataOnEntities[type][side];
1227 const int brother_side =
1228 e->getSideNumberPtr()->brother_side_number;
1229 if (brother_side != -1)
1230 brother_ents_vec.emplace_back(e);
1231 dat.getBase() = e->getApproxBase();
1232 dat.getSpace() = e->getSpace();
1233 const auto ent_order = e->getMaxOrder();
1234 dat.getOrder() =
1235 dat.getOrder() > ent_order ? dat.getOrder() : ent_order;
1236 }
1237 }
1238 }
1239 }
1240 }
1241
1242 for (auto &ent_ptr : brother_ents_vec) {
1243 if (auto e = ent_ptr.lock()) {
1244 const EntityType type = e->getEntType();
1245 const int side = e->getSideNumberPtr()->side_number;
1246 const int brother_side = e->getSideNumberPtr()->brother_side_number;
1247 auto &dat = data.dataOnEntities[type][side];
1248 auto &dat_brother = data.dataOnEntities[type][brother_side];
1249 dat_brother.getBase() = dat.getBase();
1250 dat_brother.getSpace() = dat.getSpace();
1251 dat_brother.getOrder() = dat.getOrder();
1252 }
1253 }
1254 }
1256 };
1257
1258 for (auto &e : getDataFieldEnts()) {
1259 if (auto ent_data_ptr = e.lock()) {
1260 if (ent_data_ptr->getApproxBase() == AINSWORTH_BERNSTEIN_BEZIER_BASE) {
1261 auto space = ent_data_ptr->getSpace();
1262 for (EntityType t = MBVERTEX; t != MBPOLYHEDRON; ++t) {
1263 for (auto &dat : (*dataOnElement[space]).dataOnEntities[t]) {
1264 for (auto &ptr : dat.getBBAlphaIndicesByOrderArray())
1265 ptr.reset();
1266 for (auto &ptr : dat.getBBNByOrderArray())
1267 ptr.reset();
1268 for (auto &ptr : dat.getBBDiffNByOrderArray())
1269 ptr.reset();
1270 }
1271 }
1272 }
1273 }
1274 }
1275
1276 auto check_space = [&](const auto space) {
1277 switch (space) {
1278 case H1:
1279 for (auto t = MBVERTEX; t <= ele_type; ++t) {
1280 if (dataOnElement[H1]->spacesOnEntities[t].test(H1))
1281 return true;
1282 }
1283 return false;
1284 case HCURL:
1285 for (auto t = MBEDGE; t <= ele_type; ++t) {
1286 if (dataOnElement[HCURL]->spacesOnEntities[t].test(HCURL))
1287 return true;
1288 }
1289 return false;
1290 case HDIV:
1291 for (auto t = MBTRI; t <= ele_type; ++t) {
1292 if (dataOnElement[HDIV]->spacesOnEntities[t].test(HDIV))
1293 return true;
1294 }
1295 return false;
1296 case L2:
1297 return dataOnElement[L2]->spacesOnEntities[ele_type].test(L2);
1298 break;
1299 default:
1300 THROW_MESSAGE("Not implemented");
1301 }
1302 };
1303
1304 std::set<string> fields_list;
1305 for (auto &e : getDataFieldEnts()) {
1306 if (auto ent_data_ptr = e.lock()) {
1307 if (ent_data_ptr->getApproxBase() == AINSWORTH_BERNSTEIN_BEZIER_BASE) {
1308 auto field_name = ent_data_ptr->getName();
1309 if (fields_list.find(field_name) == fields_list.end()) {
1310 auto field_ptr = ent_data_ptr->getFieldRawPtr();
1311 auto space = ent_data_ptr->getSpace();
1312 CHKERR get_nodal_base_data(*dataOnElement[space], field_ptr);
1313 CHKERR get_entity_base_data(*dataOnElement[space], field_ptr, MBEDGE,
1314 ele_type);
1315 if (check_space(space)) {
1316#ifndef NDEBUG
1317 if (ent_data_ptr->getContinuity() != CONTINUOUS)
1318 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1319 "Broken space not implemented");
1320#endif // NDEBUG
1322 -> getValue(gaussPts,
1323 boost::make_shared<EntPolynomialBaseCtx>(
1324 *dataOnElement[space], field_name,
1325 static_cast<FieldSpace>(space), CONTINUOUS,
1327 fields_list.insert(field_name);
1328 }
1329 }
1330 }
1331 }
1332 }
1333
1335};
1336
1339
1340 // Data on elements for proper spaces
1341 for (int space = H1; space != LASTSPACE; ++space) {
1342 dataOnElement[space]->setElementType(type);
1343 derivedDataOnElement[space]->setElementType(type);
1344 }
1345
1347}
1348
1349#define LOG_FUNCTION_NAME_WITH_OP_NAME(OP, CHANNEL, SEV) \
1350 MOFEM_LOG(CHANNEL, SEV) \
1351 << "(Calling user data operator " \
1352 << boost::typeindex::type_id_runtime(OP).pretty_name() << " rowField " \
1353 << (OP).rowFieldName << " colField " << (OP).colFieldName << ") "
1354
1355#define CATCH_OP_ERRORS(OP) \
1356 catch (MoFEMExceptionInitial const &ex) { \
1357 LOG_FUNCTION_NAME_WITH_OP_NAME(OP, "SELF", Sev::error) \
1358 << "in " << PETSC_FUNCTION_NAME; \
1359 return PetscError(PETSC_COMM_SELF, ex.lINE, PETSC_FUNCTION_NAME, __FILE__, \
1360 ex.errorCode, PETSC_ERROR_INITIAL, "%s", ex.what()); \
1361 } \
1362 catch (MoFEMExceptionRepeat const &ex) { \
1363 LOG_FUNCTION_NAME_WITH_OP_NAME(OP, "SELF", Sev::error) \
1364 << "in " << PETSC_FUNCTION_NAME; \
1365 return PetscError(PETSC_COMM_SELF, ex.lINE, PETSC_FUNCTION_NAME, __FILE__, \
1366 ex.errorCode, PETSC_ERROR_REPEAT, " "); \
1367 } \
1368 catch (MoFEMException const &ex) { \
1369 LOG_FUNCTION_NAME_WITH_OP_NAME(OP, "SELF", Sev::error) \
1370 << "in " << PETSC_FUNCTION_NAME; \
1371 SETERRQ(PETSC_COMM_SELF, ex.errorCode, "%s", ex.errorMessage); \
1372 } \
1373 catch (std::exception const &ex) { \
1374 LOG_FUNCTION_NAME_WITH_OP_NAME(OP, "SELF", Sev::error) \
1375 << "\nError: " << ex.what() << " at " << __LINE__ << " : " __FILE__ \
1376 << " in " << PETSC_FUNCTION_NAME; \
1377 SETERRQ(PETSC_COMM_SELF, MOFEM_STD_EXCEPTION_THROW, "%s", \
1378 "std::exception"); \
1379 }
1380
1381#ifndef NDEBUG
1382 #define LOG_OP(OP) \
1383 LOG_FUNCTION_NAME_WITH_OP_NAME(OP, "WORLD", Sev::noisy) \
1384 << "\nLocation: at " << __LINE__ << " : " __FILE__ << " in " \
1385 << PETSC_FUNCTION_NAME;
1386#else
1387 #define LOG_OP(OP)
1388#endif
1389
1392
1393 using UDO = UserDataOperator;
1394
1395 std::array<std::string, 2> field_name;
1396 std::array<const Field *, 3> field_struture;
1397 std::array<int, 2>
1398 field_id; // note the this is field bit number (nor field bit)
1399 std::array<FieldSpace, 2> space;
1400 std::array<FieldApproximationBase, 2> base;
1401
1402 constexpr std::array<UDO::OpType, 2> types = {UDO::OPROW, UDO::OPCOL};
1403 std::array<int, 2> last_eval_field_id = {0, 0};
1404
1405 std::array<boost::shared_ptr<EntitiesFieldData>, 2> op_data;
1406
1407 auto swap_bases = [&](auto &op) {
1409 for (size_t ss = 0; ss != 2; ++ss) {
1410 if (op.getOpType() & types[ss] || op.getOpType() & UDO::OPROWCOL) {
1411 switch (base[ss]) {
1413 CHKERR op_data[ss]->baseSwap(field_name[ss],
1415 default:
1416 break;
1417 }
1418 }
1419 }
1420
1422 };
1423
1424 const EntityType type = numeredEntFiniteElementPtr->getEntType();
1425
1426 // evaluate entity data only, no field specific data provided or known
1427 auto evaluate_op_space = [&](auto &op) {
1429
1430 // reseat all data which all field dependent
1431 dataOnElement[op.sPace]->resetFieldDependentData();
1432 std::fill(last_eval_field_id.begin(), last_eval_field_id.end(), 0);
1433
1434 switch (op.sPace) {
1435 case NOSPACE:
1436 LOG_OP(op);
1437 try {
1438 CHKERR op.doWork(
1439 0, MBENTITYSET,
1440 dataOnElement[NOSPACE]->dataOnEntities[MBENTITYSET][0]);
1441 }
1442 CATCH_OP_ERRORS(op);
1443 break;
1444 case NOFIELD:
1445 case H1:
1446 case HCURL:
1447 case HDIV:
1448 case L2:
1449 LOG_OP(op);
1450 try {
1451
1452 for (EntityType t = MBVERTEX; t != MBMAXTYPE; ++t) {
1453 for (auto &e : dataOnElement[op.sPace]->dataOnEntities[t]) {
1454 e.getSpace() = op.sPace;
1455 }
1456 }
1457
1458 CHKERR op.opRhs(*dataOnElement[op.sPace], false);
1459 }
1460 CATCH_OP_ERRORS(op);
1461 break;
1462 default:
1463 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
1464 "Not implemented for this space <%s>", FieldSpaceNames[op.sPace]);
1465 }
1466
1468 };
1469
1470 // set entity data
1471 auto set_op_entities_data = [&](auto ss, auto &op) {
1473
1474#ifndef NDEBUG
1475 if ((op.getNumeredEntFiniteElementPtr()->getBitFieldIdData() &
1477 .none()) {
1478 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1479 "no data field < %s > on finite element < %s >",
1480 field_name[ss].c_str(), getFEName().c_str());
1481 }
1482#endif
1483
1484 op_data[ss] =
1485 !ss ? dataOnElement[space[ss]] : derivedDataOnElement[space[ss]];
1486
1487 for (auto &data : op_data[ss]->dataOnEntities[MBENTITYSET]) {
1488 CHKERR data.resetFieldDependentData();
1489 }
1490
1491 auto get_data_for_nodes = [&]() {
1493 if (!ss)
1494 CHKERR getRowNodesIndices(*op_data[ss], field_id[ss]);
1495 else
1496 CHKERR getColNodesIndices(*op_data[ss], field_id[ss]);
1497 CHKERR getNodesFieldData(*op_data[ss], field_id[ss]);
1499 };
1500
1501 // get data on entities but not nodes
1502 auto get_data_for_entities = [&]() {
1504 CHKERR getEntityFieldData(*op_data[ss], field_id[ss], MBEDGE);
1505 if (!ss)
1506 CHKERR getEntityRowIndices(*op_data[ss], field_id[ss], MBEDGE);
1507 else
1508 CHKERR getEntityColIndices(*op_data[ss], field_id[ss], MBEDGE);
1510 };
1511
1512 auto get_data_for_meshset = [&]() {
1514 if (!ss) {
1515 CHKERR getNoFieldRowIndices(*op_data[ss], field_id[ss]);
1516 } else {
1517 CHKERR getNoFieldColIndices(*op_data[ss], field_id[ss]);
1518 }
1519 CHKERR getNoFieldFieldData(*op_data[ss], field_id[ss]);
1521 };
1522
1523 switch (space[ss]) {
1524 case H1:
1525 CHKERR get_data_for_nodes();
1526 case HCURL:
1527 case HDIV:
1528 CHKERR get_data_for_entities();
1529 break;
1530 case L2:
1531 switch (type) {
1532 case MBVERTEX:
1533 CHKERR get_data_for_nodes();
1534 break;
1535 default:
1536 CHKERR get_data_for_entities();
1537 }
1538 break;
1539 case NOFIELD:
1540 // if (!getNinTheLoop()) {
1541 // NOFIELD data are the same for each element, can be
1542 // retrieved only once
1543 CHKERR get_data_for_meshset();
1544 // }
1545 break;
1546 default:
1547 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
1548 "not implemented for this space < %s >",
1549 FieldSpaceNames[space[ss]]);
1550 }
1552 };
1553
1554 // evaluate operators with field data
1555 auto evaluate_op_for_fields = [&](auto &op) {
1557
1558 if (op.getOpType() & UDO::OPROW) {
1559 LOG_OP(op);
1560 try {
1561 CHKERR op.opRhs(*op_data[0], false);
1562 }
1563 CATCH_OP_ERRORS(op);
1564 }
1565
1566 if (op.getOpType() & UDO::OPCOL) {
1567 LOG_OP(op);
1568 try {
1569 CHKERR op.opRhs(*op_data[1], false);
1570 }
1571 CATCH_OP_ERRORS(op);
1572 }
1573
1574 if (op.getOpType() & UDO::OPROWCOL) {
1575 LOG_OP(op);
1576 try {
1577 CHKERR op.opLhs(*op_data[0], *op_data[1]);
1578 }
1579 CATCH_OP_ERRORS(op);
1580 }
1582 };
1583
1584 // collect bit ref level on all entities, fields and spaces
1586
1587 auto oit = opPtrVector.begin();
1588 auto hi_oit = opPtrVector.end();
1589
1590 // iterate over all operators on element
1591 for (; oit != hi_oit; oit++) {
1592
1593 LOG_OP(*oit);
1594 try {
1595
1596 CHKERR oit->setPtrFE(this);
1597
1598 if ((oit->opType & UDO::OPSPACE) == UDO::OPSPACE) {
1599
1600 // run operator for space but specific field
1601 CHKERR evaluate_op_space(*oit);
1602
1603 } else if (
1604
1605 (oit->opType & (UDO::OPROW | UDO::OPCOL | UDO::OPROWCOL)) ==
1606 oit->opType
1607
1608 ) {
1609
1610 field_name[0] = oit->rowFieldName;
1611 field_name[1] = oit->colFieldName;
1612
1613 for (size_t ss = 0; ss != 2; ss++) {
1614
1615#ifndef NDEBUG
1616 if (field_name[ss].empty()) {
1617 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1618 "Not set Field name in operator %ld (0-row, 1-column) in "
1619 "operator %s",
1620 ss,
1621 (boost::typeindex::type_id_runtime(*oit).pretty_name())
1622 .c_str());
1623 }
1624#endif
1625
1626 field_struture[ss] = mField.get_field_structure(field_name[ss]);
1627 field_id[ss] = field_struture[ss]->getBitNumber();
1628 space[ss] = field_struture[ss]->getSpace();
1629 base[ss] = field_struture[ss]->getApproxBase();
1630 }
1631
1632 // not that if field name do not change between operators, entity field
1633 // data are nor rebuild
1634 for (size_t ss = 0; ss != 2; ss++) {
1635
1636 if (oit->getOpType() & types[ss] ||
1637 oit->getOpType() & UDO::OPROWCOL) {
1638 if (last_eval_field_id[ss] != field_id[ss]) {
1639 CHKERR set_op_entities_data(ss, *oit);
1640 last_eval_field_id[ss] = field_id[ss];
1641 }
1642 }
1643 }
1644
1645 CHKERR swap_bases(*oit);
1646
1647 // run operator for given field or row, column or both
1648 CHKERR evaluate_op_for_fields(*oit);
1649
1650 CHKERR swap_bases(*oit);
1651
1652 } else {
1653
1654 for (int i = 0; i != 5; ++i)
1655 if (oit->opType & (1 << i))
1656 MOFEM_LOG("SELF", Sev::error) << UDO::OpTypeNames[i];
1657 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1658 "Impossible operator type");
1659 }
1660 }
1661 CATCH_OP_ERRORS(*oit);
1662 }
1664}
1665
1666const char *const ForcesAndSourcesCore::UserDataOperator::OpTypeNames[] = {
1667 "OPROW", " OPCOL", "OPROWCOL", "OPSPACE", "OPLAST"};
1668
1669MoFEMErrorCode ForcesAndSourcesCore::UserDataOperator::getProblemRowIndices(
1670 const std::string field_name, const EntityType type, const int side,
1671 VectorInt &indices) const {
1673 if (ptrFE == NULL)
1674 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
1675
1676 switch (type) {
1677 case MBVERTEX:
1679 break;
1680 default:
1681 CHKERR ptrFE->getProblemTypeRowIndices(field_name, type, side, indices);
1682 }
1684}
1685
1686MoFEMErrorCode ForcesAndSourcesCore::UserDataOperator::getProblemColIndices(
1687 const std::string field_name, const EntityType type, const int side,
1688 VectorInt &indices) const {
1690 if (ptrFE == NULL)
1691 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
1692
1693 switch (type) {
1694 case MBVERTEX:
1695 CHKERR ptrFE->getProblemNodesColIndices(field_name, indices);
1696 break;
1697 default:
1698 CHKERR ptrFE->getProblemTypeColIndices(field_name, type, side, indices);
1699 }
1701}
1702
1709
1716
1717MoFEMErrorCode ForcesAndSourcesCore::UserDataOperator::loopSide(
1718 const string &fe_name, ForcesAndSourcesCore *side_fe, const size_t side_dim,
1719 const EntityHandle ent_for_side, boost::shared_ptr<Range> fe_range,
1720 const int verb, const LogManager::SeverityLevel sev, AdjCache *adj_cache) {
1722
1723 const auto *problem_ptr = getFEMethod()->problemPtr;
1724 auto &numered_fe =
1725 problem_ptr->numeredFiniteElementsPtr->get<Unique_mi_tag>();
1726
1727 auto fe_miit = ptrFE->mField.get_finite_elements()
1729 .find(fe_name);
1730 if (fe_miit != ptrFE->mField.get_finite_elements()
1732 .end()) {
1733
1734 const auto ent = ent_for_side ? ent_for_side : getFEEntityHandle();
1735
1736 side_fe->feName = fe_name;
1737
1738 CHKERR side_fe->setSideFEPtr(ptrFE);
1739 CHKERR side_fe->copyBasicMethod(*getFEMethod());
1740 CHKERR side_fe->copyPetscData(*getFEMethod());
1741 CHKERR side_fe->copyKsp(*getFEMethod());
1742 CHKERR side_fe->copySnes(*getFEMethod());
1743 CHKERR side_fe->copyTs(*getFEMethod());
1744
1745 side_fe->cacheWeakPtr = getFEMethod()->cacheWeakPtr;
1746
1747 CHKERR side_fe->preProcess();
1748
1749 std::vector<boost::weak_ptr<NumeredEntFiniteElement>> fe_vec;
1750 auto get_numered_fe_ptr = [&](auto &fe_uid, Range &&adjacent_ents)
1751 -> std::vector<boost::weak_ptr<NumeredEntFiniteElement>> & {
1752 fe_vec.reserve(adjacent_ents.size());
1753 for (auto fe_ent : adjacent_ents) {
1754 auto miit = numered_fe.find(
1756 if (miit != numered_fe.end()) {
1757 fe_vec.emplace_back(*miit);
1758 }
1759 }
1760 return fe_vec;
1761 };
1762
1763 auto get_bit_entity_adjacency = [&]() {
1764 Range adjacent_ents;
1766 ptrFE->mField.getInterface<BitRefManager>()->getAdjacenciesAny(
1767 ent, side_dim, adjacent_ents),
1768 "getAdjacenciesAny failed");
1769 return adjacent_ents;
1770 };
1771
1772 auto get_bit_meshset_entities = [&]() {
1773 auto &bit = getFEMethod()->numeredEntFiniteElementPtr->getBitRefLevel();
1774 Range ents = *fe_range;
1776 ptrFE->mField.getInterface<BitRefManager>()->filterEntitiesByRefLevel(
1777 bit, BitRefLevel().set(), ents),
1778 "filterEntitiesByRefLevel failed");
1779 return ents;
1780 };
1781
1782 auto get_adj = [&](auto &fe_uid, auto get_adj_fun)
1783 -> std::vector<boost::weak_ptr<NumeredEntFiniteElement>> & {
1784 if (adj_cache) {
1785 try {
1786 return (*adj_cache).at(ent);
1787 } catch (const std::out_of_range &) {
1788 return (*adj_cache)[ent] = get_numered_fe_ptr(fe_uid, get_adj_fun());
1789 }
1790 } else {
1791 return get_numered_fe_ptr(fe_uid, get_adj_fun());
1792 }
1793 };
1794
1795 if (type_from_handle(ent) == MBENTITYSET) {
1796 if (!fe_range)
1797 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1798 "No range of finite elements");
1799 }
1800 auto adj = (!fe_range)
1801 ? get_adj((*fe_miit)->getFEUId(), get_bit_entity_adjacency)
1802 : get_adj((*fe_miit)->getFEUId(), get_bit_meshset_entities);
1803
1804 if (verb >= VERBOSE && !adj.empty())
1805 MOFEM_LOG("SELF", sev) << "Number of side finite elements " << adj.size();
1806
1807 int nn = 0;
1808 side_fe->loopSize = adj.size();
1809 for (auto fe_weak_ptr : adj) {
1810 if (auto fe_ptr = fe_weak_ptr.lock()) {
1811 if (verb >= VERBOSE)
1812 MOFEM_LOG("SELF", sev)
1813 << "Side finite element " << "(" << nn << "): " << *fe_ptr;
1814 side_fe->nInTheLoop = nn;
1815 side_fe->numeredEntFiniteElementPtr = fe_ptr;
1816 CHKERR (*side_fe)();
1817 }
1818 ++nn;
1819 }
1820
1821 CHKERR side_fe->postProcess();
1822 }
1823
1825}
1826
1827MoFEMErrorCode ForcesAndSourcesCore::UserDataOperator::loopThis(
1828 const string &fe_name, ForcesAndSourcesCore *this_fe, const int verb,
1829 const LogManager::SeverityLevel sev) {
1831
1832
1833
1834 auto &fes =
1835 ptrFE->mField.get_finite_elements()->get<FiniteElement_name_mi_tag>();
1836 auto fe_miit = fes.find(fe_name);
1837 if (fe_miit != fes.end()) {
1838
1839 this_fe->feName = fe_name;
1840
1841 CHKERR this_fe->setRefineFEPtr(ptrFE);
1842 CHKERR this_fe->copyBasicMethod(*getFEMethod());
1843 CHKERR this_fe->copyPetscData(*getFEMethod());
1844 CHKERR this_fe->copyKsp(*getFEMethod());
1845 CHKERR this_fe->copySnes(*getFEMethod());
1846 CHKERR this_fe->copyTs(*getFEMethod());
1847
1848 this_fe->cacheWeakPtr = getFEMethod()->cacheWeakPtr;
1849
1850 this_fe->nInTheLoop = getNinTheLoop();
1851 this_fe->loopSize = getLoopSize();
1852
1853 CHKERR this_fe->preProcess();
1854
1855 if (fe_name == getFEName()) {
1856
1857 if (verb >= VERBOSE)
1858 MOFEM_LOG("SELF", sev)
1859 << "This finite element: " << *getNumeredEntFiniteElementPtr();
1860
1861 this_fe->numeredEntFiniteElementPtr = getNumeredEntFiniteElementPtr();
1862 CHKERR (*this_fe)();
1863 } else {
1864
1865 auto get_numered_fe_ptr = [&](auto &fe_uid, auto fe_ent) {
1866 auto &numered_fe =
1867 getFEMethod()
1868 ->problemPtr->numeredFiniteElementsPtr->get<Unique_mi_tag>();
1869 auto it = numered_fe.find(
1871 boost::shared_ptr<const NumeredEntFiniteElement> this_numered_fe_ptr;
1872 if (it != numered_fe.end()) {
1873 this_numered_fe_ptr = *it;
1874 }
1875 return this_numered_fe_ptr;
1876 };
1877
1878 auto this_numered_fe_ptr = get_numered_fe_ptr(
1879 (*fe_miit)->getFEUId(), getNumeredEntFiniteElementPtr()->getEnt());
1880 if (this_numered_fe_ptr) {
1881
1882 if (verb >= VERBOSE)
1883 MOFEM_LOG("SELF", sev)
1884 << "This finite element: " << *this_numered_fe_ptr;
1885
1886 this_fe->numeredEntFiniteElementPtr = this_numered_fe_ptr;
1887 CHKERR (*this_fe)();
1888 }
1889 }
1890 CHKERR this_fe->postProcess();
1891 }
1892
1894}
1895
1896MoFEMErrorCode ForcesAndSourcesCore::UserDataOperator::loopParent(
1897 const string &fe_name, ForcesAndSourcesCore *parent_fe, const int verb,
1898 const LogManager::SeverityLevel sev) {
1900
1901 auto &fes =
1902 ptrFE->mField.get_finite_elements()->get<FiniteElement_name_mi_tag>();
1903 auto fe_miit = fes.find(fe_name);
1904 if (fe_miit != fes.end()) {
1905
1906 const auto *problem_ptr = getFEMethod()->problemPtr;
1907 auto &numered_fe =
1908 problem_ptr->numeredFiniteElementsPtr->get<Unique_mi_tag>();
1909
1910 parent_fe->feName = fe_name;
1911 CHKERR parent_fe->setRefineFEPtr(ptrFE);
1912 CHKERR parent_fe->copyBasicMethod(*getFEMethod());
1913 CHKERR parent_fe->copyPetscData(*getFEMethod());
1914 CHKERR parent_fe->copyKsp(*getFEMethod());
1915 CHKERR parent_fe->copySnes(*getFEMethod());
1916 CHKERR parent_fe->copyTs(*getFEMethod());
1917
1918 parent_fe->cacheWeakPtr = getFEMethod()->cacheWeakPtr;
1919
1920 const auto parent_ent = getNumeredEntFiniteElementPtr()->getParentEnt();
1921 auto miit = numered_fe.find(EntFiniteElement::getLocalUniqueIdCalculate(
1922 parent_ent, (*fe_miit)->getFEUId()));
1923 if (miit != numered_fe.end()) {
1924 if (verb >= VERBOSE)
1925 MOFEM_LOG("SELF", sev) << "Parent finite element: " << **miit;
1926 parent_fe->loopSize = 1;
1927 parent_fe->nInTheLoop = 0;
1928 parent_fe->numeredEntFiniteElementPtr = *miit;
1929 CHKERR parent_fe->preProcess();
1930 CHKERR (*parent_fe)();
1931 CHKERR parent_fe->postProcess();
1932 } else {
1933 if (verb >= VERBOSE)
1934 MOFEM_LOG("SELF", sev) << "Parent finite element: no parent";
1935 parent_fe->loopSize = 0;
1936 parent_fe->nInTheLoop = 0;
1937 CHKERR parent_fe->preProcess();
1938 CHKERR parent_fe->postProcess();
1939 }
1940 }
1941
1943}
1944
1945MoFEMErrorCode ForcesAndSourcesCore::UserDataOperator::loopChildren(
1946 const string &fe_name, ForcesAndSourcesCore *child_fe, const int verb,
1947 const LogManager::SeverityLevel sev) {
1949
1950 auto fe_miit = ptrFE->mField.get_finite_elements()
1952 .find(fe_name);
1953 if (fe_miit != ptrFE->mField.get_finite_elements()
1955 .end()) {
1956
1957 const auto *problem_ptr = getFEMethod()->problemPtr;
1958 auto &ref_ents = *getPtrFE()->mField.get_ref_ents();
1959 auto &numered_fe =
1960 problem_ptr->numeredFiniteElementsPtr->get<Unique_mi_tag>();
1961
1962 const auto parent_ent = getNumeredEntFiniteElementPtr()->getEnt();
1963 const auto parent_type = getNumeredEntFiniteElementPtr()->getEntType();
1964 auto range =
1965 ref_ents.get<Composite_ParentEnt_And_EntType_mi_tag>().equal_range(
1966 boost::make_tuple(parent_type, parent_ent));
1967
1968 if (auto size = std::distance(range.first, range.second)) {
1969
1970 std::vector<EntityHandle> childs_vec;
1971 childs_vec.reserve(size);
1972 for (; range.first != range.second; ++range.first)
1973 childs_vec.emplace_back((*range.first)->getEnt());
1974
1975 Range childs;
1976
1977 if ((childs_vec.back() - childs_vec.front() + 1) == size)
1978 childs = Range(childs_vec.front(), childs_vec.back());
1979 else
1980 childs.insert_list(childs_vec.begin(), childs_vec.end());
1981
1982 child_fe->feName = fe_name;
1983
1984 CHKERR child_fe->setRefineFEPtr(ptrFE);
1985 CHKERR child_fe->copyBasicMethod(*getFEMethod());
1986 CHKERR child_fe->copyPetscData(*getFEMethod());
1987 CHKERR child_fe->copyKsp(*getFEMethod());
1988 CHKERR child_fe->copySnes(*getFEMethod());
1989 CHKERR child_fe->copyTs(*getFEMethod());
1990
1991 child_fe->cacheWeakPtr = getFEMethod()->cacheWeakPtr;
1992
1993 CHKERR child_fe->preProcess();
1994
1995 int nn = 0;
1996 child_fe->loopSize = size;
1997
1998 for (auto p = childs.pair_begin(); p != childs.pair_end(); ++p) {
1999
2000 auto miit =
2001 numered_fe.lower_bound(EntFiniteElement::getLocalUniqueIdCalculate(
2002 p->first, (*fe_miit)->getFEUId()));
2003 auto hi_miit =
2004 numered_fe.upper_bound(EntFiniteElement::getLocalUniqueIdCalculate(
2005 p->second, (*fe_miit)->getFEUId()));
2006
2007 for (; miit != hi_miit; ++miit) {
2008
2009 if (verb >= VERBOSE)
2010 MOFEM_LOG("SELF", sev)
2011 << "Child finite element " << "(" << nn << "): " << **miit;
2012
2013 child_fe->nInTheLoop = nn++;
2014 child_fe->numeredEntFiniteElementPtr = *miit;
2015 CHKERR (*child_fe)();
2016 }
2017 }
2018 }
2019
2020 CHKERR child_fe->postProcess();
2021 }
2022
2024}
2025
2026MoFEMErrorCode ForcesAndSourcesCore::UserDataOperator::loopRange(
2027 const string &fe_name, ForcesAndSourcesCore *range_fe,
2028 boost::shared_ptr<Range> fe_range, const int verb,
2029 const LogManager::SeverityLevel sev) {
2031
2032 auto &fes =
2033 ptrFE->mField.get_finite_elements()->get<FiniteElement_name_mi_tag>();
2034 auto fe_miit = fes.find(fe_name);
2035 if (fe_miit != fes.end()) {
2036
2037 const auto *problem_ptr = getFEMethod()->problemPtr;
2038 auto &numered_fe =
2039 problem_ptr->numeredFiniteElementsPtr->get<Unique_mi_tag>();
2040
2041 range_fe->feName = fe_name;
2042 CHKERR range_fe->setRefineFEPtr(ptrFE);
2043 CHKERR range_fe->copyBasicMethod(*getFEMethod());
2044 CHKERR range_fe->copyPetscData(*getFEMethod());
2045 CHKERR range_fe->copyKsp(*getFEMethod());
2046 CHKERR range_fe->copySnes(*getFEMethod());
2047 CHKERR range_fe->copyTs(*getFEMethod());
2048
2049 range_fe->cacheWeakPtr = getFEMethod()->cacheWeakPtr;
2050
2051 auto get_numered_fe_ptr = [&](auto &fe_uid, auto fe_range, auto execute) {
2053 if (fe_range) {
2054 int nn = 0;
2055 for (auto p = fe_range->pair_begin(); p != fe_range->pair_end(); ++p) {
2056 auto first = p->first;
2057 auto second = p->second;
2058 auto lo = numered_fe.lower_bound(
2060 auto hi = numered_fe.upper_bound(
2062 for (; lo != hi; ++lo) {
2063 CHKERR execute(lo, nn++);
2064 }
2065 }
2066 }
2068 };
2069
2070 auto execute = [&](auto lo, auto nn) {
2072 if (verb >= VERBOSE)
2073 MOFEM_LOG("SELF", sev) << "Range finite element: " << **lo;
2074 range_fe->nInTheLoop = nn;
2075 range_fe->numeredEntFiniteElementPtr = *lo;
2076 CHKERR range_fe->preProcess();
2077 CHKERR (*range_fe)();
2078 CHKERR range_fe->postProcess();
2080 };
2081
2082 CHKERR get_numered_fe_ptr((*fe_miit)->getFEUId(), fe_range, execute);
2083 }
2084
2086}
2087
2088int ForcesAndSourcesCore::getRule(int order_row, int order_col,
2089 int order_data) {
2090 return getRuleHook ? getRuleHook(order_row, order_col, order_data)
2091 : getRule(order_data);
2092}
2093
2095 int order_data) {
2096 return setRuleHook ? setRuleHook(this, order_row, order_col, order_data)
2097 : setGaussPts(order_data);
2098}
2099
2101
2102/** \deprecated setGaussPts(int row_order, int col_order, int data order);
2103 */
2106 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "Sorry, not implemented");
2108}
2109
2110ForcesAndSourcesCore::UserDataOperator::UserDataOperator(const FieldSpace space,
2111 const char type,
2112 const bool symm)
2113 : DataOperator(symm), opType(type), sPace(space), ptrFE(nullptr) {}
2114
2115ForcesAndSourcesCore::UserDataOperator::UserDataOperator(
2116 const std::string field_name, const char type, const bool symm)
2117 : DataOperator(symm), opType(type), rowFieldName(field_name),
2118 colFieldName(field_name), sPace(LASTSPACE), ptrFE(nullptr) {}
2119
2120ForcesAndSourcesCore::UserDataOperator::UserDataOperator(
2121 const std::string row_field_name, const std::string col_field_name,
2122 const char type, const bool symm)
2123 : DataOperator(symm), opType(type), rowFieldName(row_field_name),
2124 colFieldName(col_field_name), sPace(LASTSPACE), ptrFE(nullptr) {}
2125
2136 if (operatorHook) {
2137 ierr = operatorHook();
2138 CHKERRG(ierr);
2139 } else {
2140#ifndef NDEBUG
2141 MOFEM_LOG("SELF", Sev::warning)
2142 << "No method operator() overloaded on element entity on finite "
2143 "element <"
2144 << boost::typeindex::type_id_runtime(*this).pretty_name() << ">";
2145#endif
2146 }
2148}
2157
2159ForcesAndSourcesCore::UserDataOperator::setPtrFE(ForcesAndSourcesCore *ptr) {
2161 if (!(ptrFE = dynamic_cast<ForcesAndSourcesCore *>(ptr)))
2162 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
2163 "User operator and finite element do not work together");
2165}
2166
2167} // namespace MoFEM
#define CATCH_OP_ERRORS(OP)
#define LOG_OP(OP)
constexpr double a
@ VERBOSE
FieldApproximationBase
approximation base
Definition definitions.h:58
@ LASTBASE
Definition definitions.h:69
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition definitions.h:60
@ USER_BASE
user implemented approximation base
Definition definitions.h:68
@ NOBASE
Definition definitions.h:59
@ AINSWORTH_BERNSTEIN_BEZIER_BASE
Definition definitions.h:64
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
FieldSpace
approximation spaces
Definition definitions.h:82
@ L2
field with C-1 continuity
Definition definitions.h:88
@ LASTSPACE
FieldSpace in [ 0, LASTSPACE )
Definition definitions.h:89
@ NOFIELD
scalar or vector of scalars describe (no true field)
Definition definitions.h:84
@ H1
continuous field
Definition definitions.h:85
@ NOSPACE
Definition definitions.h:83
@ HCURL
field with continuous tangents
Definition definitions.h:86
@ HDIV
field with continuous normal traction
Definition definitions.h:87
FieldContinuity
Field continuity.
Definition definitions.h:99
@ CONTINUOUS
Regular field.
@ DISCONTINUOUS
Broken continuity (No effect on L2 space)
static const char *const FieldSpaceNames[]
Definition definitions.h:92
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
@ 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 ...
#define THROW_MESSAGE(msg)
Throw MoFEM exception.
constexpr int order
multi_index_container< boost::shared_ptr< NumeredDofEntity >, indexed_by< ordered_unique< tag< Unique_mi_tag >, const_mem_fun< NumeredDofEntity::interface_type_DofEntity, UId, &NumeredDofEntity::getLocalUniqueId > >, ordered_non_unique< tag< Part_mi_tag >, member< NumeredDofEntity, unsigned int, &NumeredDofEntity::pArt > >, ordered_non_unique< tag< Idx_mi_tag >, member< NumeredDofEntity, DofIdx, &NumeredDofEntity::dofIdx > >, ordered_non_unique< tag< PetscGlobalIdx_mi_tag >, member< NumeredDofEntity, DofIdx, &NumeredDofEntity::petscGloablDofIdx > >, ordered_non_unique< tag< PetscLocalIdx_mi_tag >, member< NumeredDofEntity, DofIdx, &NumeredDofEntity::petscLocalDofIdx > >, ordered_non_unique< tag< Ent_mi_tag >, const_mem_fun< NumeredDofEntity::interface_type_DofEntity, EntityHandle, &NumeredDofEntity::getEnt > > > > NumeredDofEntity_multiIndex
MultiIndex container keeps NumeredDofEntity.
multi_index_container< boost::shared_ptr< SideNumber >, indexed_by< ordered_unique< member< SideNumber, EntityHandle, &SideNumber::ent > >, ordered_non_unique< composite_key< SideNumber, const_mem_fun< SideNumber, EntityType, &SideNumber::getEntType >, member< SideNumber, signed char, &SideNumber::side_number > > > > > SideNumber_multiIndex
SideNumber_multiIndex for SideNumber.
virtual MoFEMErrorCode getAdjacenciesAny(const EntityHandle from_entity, const int to_dimension, Range &adj_entities) const
Get the adjacencies associated with a entity to entities of a specified dimension.
virtual const Field * get_field_structure(const std::string &name, enum MoFEMTypes bh=MF_EXIST) const =0
get field structure
#define MOFEM_LOG(channel, severity)
Log.
SeverityLevel
Severity levels.
auto bit
set bit
FTensor::Index< 'i', SPACE_DIM > i
const double n
refractive index of diffusive medium
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
std::bitset< BITFEID_SIZE > BitFEId
Finite element Id.
Definition Types.hpp:43
int ApproximationOrder
Approximation on the entity.
Definition Types.hpp:26
std::bitset< BITFIELDID_SIZE > BitFieldId
Field Id.
Definition Types.hpp:42
uint128_t UId
Unique Id.
Definition Types.hpp:31
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition Types.hpp:40
UBlasVector< int > VectorInt
Definition Types.hpp:67
implementation of Data Operators for Forces and Sources
Definition Common.hpp:10
auto type_from_handle(const EntityHandle h)
get type from entity handle
static MoFEMErrorCode get_value(MatrixDouble &pts_x, MatrixDouble &pts_t, TYPE *ctx)
static auto cmp_uid_lo(const boost::weak_ptr< FieldEntity > &a, const UId &b)
EntityHandle get_id_for_max_type()
static auto cmp_uid_hi(const UId &b, const boost::weak_ptr< FieldEntity > &a)
EntityHandle get_id_for_min_type()
ublas::vector< FieldEntity *, FieldEntAllocator > VectorFieldEntities
std::vector< boost::weak_ptr< FieldEntity > > FieldEntity_vector_view
Field_multiIndex::index< FieldName_mi_tag >::type::iterator field_it
static int getMaxOrder(const ENTMULTIINDEX &multi_index)
ublas::vector< FEDofEntity *, DofsAllocator > VectorDofs
constexpr double t
plate stiffness
Definition plate.cpp:58
constexpr auto field_name
int loopSize
local number oe methods to process
MoFEMErrorCode copyBasicMethod(const BasicMethod &basic)
Copy data from other base method to this base method.
int getLoopSize() const
get loop size
boost::function< MoFEMErrorCode()> preProcessHook
Hook function for pre-processing.
const Field_multiIndex * fieldsPtr
raw pointer to fields container
int nInTheLoop
number currently of processed method
const Problem * problemPtr
raw pointer to problem
boost::function< MoFEMErrorCode()> operatorHook
Hook function for operator.
boost::weak_ptr< CacheTuple > cacheWeakPtr
int getNinTheLoop() const
get number of evaluated element in the loop
boost::function< MoFEMErrorCode()> postProcessHook
Hook function for post-processing.
Managing BitRefLevels.
MoFEMErrorCode filterEntitiesByRefLevel(const BitRefLevel bit, const BitRefLevel mask, Range &ents, int verb=QUIET) const
filter entities by bit ref level
virtual FieldBitNumber get_field_bit_number(const std::string name) const =0
get field bit number
virtual moab::Interface & get_moab()=0
virtual BitFieldId get_field_id(const std::string &name) const =0
Get field Id.
base operator to do operations at Gauss Pt. level
static bool getDefTypeMap(const EntityType fe_type, const EntityType ent_type)
Deprecated interface functions.
this class derive data form other data structure
UId getLocalUniqueIdCalculate() const
Generate UId for finite element entity.
Data on single entity (This is passed as argument to DataOperator::doWork)
data structure for finite element entity
MatrixInt facesNodes
nodes on finite element faces
std::array< std::bitset< LASTSPACE >, MBMAXTYPE > spacesOnEntities
spaces on entity types
std::bitset< LASTBASE > bAse
bases on element
std::array< boost::ptr_vector< EntData >, MBMAXTYPE > dataOnEntities
std::array< std::bitset< LASTBASE >, LASTSPACE > basesOnSpaces
base on spaces
std::array< std::bitset< LASTBASE >, LASTSPACE > brokenBasesOnSpaces
base on spaces
MatrixInt facesNodesOrder
order of face nodes on element
std::array< std::bitset< LASTBASE >, MBMAXTYPE > basesOnEntities
bases on entity types
keeps information about indexed dofs for the finite element
std::string feName
Name of finite element.
auto & getRowFieldEnts() const
auto getFEName() const
get finite element name
const FieldEntity_vector_view & getDataFieldEnts() const
boost::shared_ptr< FieldEntity_vector_view > & getDataFieldEntsPtr() const
auto getColDofsPtr() const
EntityHandle getFEEntityHandle() const
auto getNumberOfNodes() const
auto getRowDofsPtr() const
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
auto & getColFieldEnts() const
static UId getHiBitNumberUId(const FieldBitNumber bit_number)
UId getLocalUniqueIdCalculate()
Get the Local Unique Id Calculate object.
static UId getLoBitNumberUId(const FieldBitNumber bit_number)
static UId getLoLocalEntityBitNumber(const char bit_number, const EntityHandle ent)
static UId getHiLocalEntityBitNumber(const char bit_number, const EntityHandle ent)
Provide data structure for (tensor) field approximation.
FieldCoefficientsNumber getNbOfCoeffs() const
Get number of field coefficients.
FieldSpace getSpace() const
Get field approximation space.
FieldBitNumber getBitNumber() const
Get number of set bit in Field ID. Each field has uid, get getBitNumber get number of bit set for giv...
std::map< EntityHandle, std::vector< boost::weak_ptr< NumeredEntFiniteElement > > > AdjCache
structure to get information form mofem into EntitiesFieldData
MoFEMErrorCode loopOverOperators()
Iterate user data operators.
MoFEMErrorCode getProblemNodesRowIndices(const std::string &field_name, VectorInt &nodes_indices) const
MoFEMErrorCode getEntityIndices(EntitiesFieldData &data, const int bit_number, FieldEntity_vector_view &ents_field, const EntityType type_lo, const EntityType type_hi, EXTRACTOR &&extractor) const
MoFEMErrorCode getEntityRowIndices(EntitiesFieldData &data, const int bit_number, const EntityType type_lo=MBVERTEX, const EntityType type_hi=MBPOLYHEDRON) const
MoFEMErrorCode getRowNodesIndices(EntitiesFieldData &data, const int bit_number) const
get row node indices from FENumeredDofEntity_multiIndex
int getMaxRowOrder() const
Get max order of approximation for field in rows.
auto & getUserPolynomialBase()
Get the User Polynomial Base object.
MoFEMErrorCode getSpacesAndBaseOnEntities(EntitiesFieldData &data) const
Get field approximation space and base on entities.
MoFEMErrorCode getProblemTypeColIndices(const std::string &field_name, EntityType type, int side_number, VectorInt &indices) const
MoFEMErrorCode getNoFieldIndices(const int bit_number, boost::shared_ptr< FENumeredDofEntity_multiIndex > dofs, VectorInt &nodes_indices) const
get NoField indices
virtual MoFEMErrorCode operator()()
function is run for every finite element
MoFEMErrorCode getFaceNodes(EntitiesFieldData &data) const
Get nodes on faces.
const std::array< boost::shared_ptr< EntitiesFieldData >, LASTSPACE > derivedDataOnElement
Entity data on element entity columns fields.
MoFEMErrorCode getNoFieldRowIndices(EntitiesFieldData &data, const int bit_number) const
get col NoField indices
virtual MoFEMErrorCode setGaussPts(int order_row, int order_col, int order_data)
set user specific integration rule
MoFEMErrorCode getNodesIndices(const int bit_number, FieldEntity_vector_view &ents_field, VectorInt &nodes_indices, VectorInt &local_nodes_indices, EXTRACTOR &&extractor) const
get node indices
MoFEMErrorCode getEntityFieldData(EntitiesFieldData &data, const int bit_number, const EntityType type_lo=MBVERTEX, const EntityType type_hi=MBPOLYHEDRON) const
MoFEMErrorCode setRefineFEPtr(const ForcesAndSourcesCore *refine_fe_ptr)
Set the pointer to face element refined.
auto & getElementPolynomialBase()
Get the Entity Polynomial Base object.
ForcesAndSourcesCore * sidePtrFE
Element to integrate on the sides.
ForcesAndSourcesCore(Interface &m_field)
virtual int getRule(int order_row, int order_col, int order_data)
another variant of getRule
boost::ptr_deque< UserDataOperator > opPtrVector
Vector of finite element users data operators.
MoFEMErrorCode getEntitySense(const EntityType type, boost::ptr_vector< EntitiesFieldData::EntData > &data) const
get sense (orientation) of entity
const std::array< boost::shared_ptr< EntitiesFieldData >, LASTSPACE > dataOnElement
Entity data on element entity rows fields.
MoFEMErrorCode getProblemNodesColIndices(const std::string &field_name, VectorInt &nodes_indices) const
MoFEMErrorCode calHierarchicalBaseFunctionsOnElement()
Calculate base functions.
MatrixDouble gaussPts
Matrix of integration points.
MoFEMErrorCode calBernsteinBezierBaseFunctionsOnElement()
Calculate Bernstein-Bezier base.
virtual MoFEMErrorCode preProcess()
function is run at the beginning of loop
MoFEMErrorCode getColNodesIndices(EntitiesFieldData &data, const int bit_number) const
get col node indices from FENumeredDofEntity_multiIndex
MoFEMErrorCode setSideFEPtr(const ForcesAndSourcesCore *side_fe_ptr)
Set the pointer to face element on the side.
int getMaxColOrder() const
Get max order of approximation for field in columns.
MoFEMErrorCode getProblemNodesIndices(const std::string &field_name, const NumeredDofEntity_multiIndex &dofs, VectorInt &nodes_indices) const
get indices of nodal indices which are declared for problem but not this particular element
GaussHookFun setRuleHook
Set function to calculate integration rule.
MoFEMErrorCode getNodesFieldData(EntitiesFieldData &data, const int bit_number) const
Get data on nodes.
MoFEMErrorCode getEntityColIndices(EntitiesFieldData &data, const int bit_number, const EntityType type_lo=MBVERTEX, const EntityType type_hi=MBPOLYHEDRON) const
MoFEMErrorCode getEntityDataOrder(const EntityType type, const FieldSpace space, boost::ptr_vector< EntitiesFieldData::EntData > &data) const
Get the entity data order.
ForcesAndSourcesCore * refinePtrFE
Element to integrate parent or child.
MoFEMErrorCode getNoFieldColIndices(EntitiesFieldData &data, const int bit_number) const
get col NoField indices
int getMaxDataOrder() const
Get max order of approximation for data fields.
MoFEMErrorCode getNoFieldFieldData(const int bit_number, VectorDouble &ent_field_data, VectorDofs &ent_field_dofs, VectorFieldEntities &ent_field) const
Get field data on nodes.
MoFEMErrorCode getProblemTypeRowIndices(const std::string &field_name, EntityType type, int side_number, VectorInt &indices) const
MoFEMErrorCode getProblemTypeIndices(const std::string &field_name, const NumeredDofEntity_multiIndex &dofs, EntityType type, int side_number, VectorInt &indices) const
get indices by type (generic function) which are declared for problem but not this particular element
RuleHookFun getRuleHook
Hook to get rule.
MoFEMErrorCode createDataOnElement(EntityType type)
Create a entity data on element object.
virtual MoFEMErrorCode postProcess()
function is run at the end of loop
MoFEMErrorCode copyKsp(const KspMethod &ksp)
copy data form another method
MoFEMErrorCode copyPetscData(const PetscData &petsc_data)
boost::shared_ptr< NumeredDofEntity_multiIndex > numeredRowDofsPtr
store DOFs on rows for this problem
boost::shared_ptr< NumeredDofEntity_multiIndex > numeredColDofsPtr
store DOFs on columns for this problem
MoFEMErrorCode copySnes(const SnesMethod &snes)
Copy snes data.
MoFEMErrorCode copyTs(const TSMethod &ts)
Copy TS solver data.
boost::shared_ptr< SideNumber > getSideNumberPtr() const