v0.14.0
Loading...
Searching...
No Matches
EntitiesFieldData.cpp
Go to the documentation of this file.
1/** \file EntitiesFieldData.cpp
2\brief Implementation for Data Structures in Forces and Sources
3
4*/
5
6
7
8namespace MoFEM {
9
10EntitiesFieldData::EntData::EntData(const bool allocate_base_matrices)
11 : sEnse(0), oRder(0), bAse(NOBASE), entDataBitRefLevel(),
12 N(baseFunctionsAndBaseDerivatives[ZeroDerivative]),
13 diffN(baseFunctionsAndBaseDerivatives[FirstDerivative]) {
14 if (allocate_base_matrices) {
15
16 for (auto d = 0; d != LastDerivative; ++d) {
17 for (int b = 0; b != LASTBASE; ++b) {
19 }
20 }
21 }
22}
23
24int EntitiesFieldData::EntData::getSense() const { return sEnse; }
25
26boost::shared_ptr<MatrixDouble> &
28 const BaseDerivatives direvatie) {
29 return baseFunctionsAndBaseDerivatives[direvatie][base];
30}
31
32boost::shared_ptr<MatrixDouble> &
34 return N[base];
35}
36
37boost::shared_ptr<MatrixDouble> &EntitiesFieldData::EntData::getDiffNSharedPtr(
38 const FieldApproximationBase base) {
39 return diffN[base];
40}
41
42static void constructor_data(EntitiesFieldData *data, const EntityType type) {
43
45
46 auto set_default = [&]() {
47 std::array<size_t, MBMAXTYPE> count;
48 std::fill(count.begin(), count.end(), 0);
49 const int dim_type = moab::CN::Dimension(type);
50 data->dataOnEntities[MBVERTEX].resize(1);
51 if (type != MBVERTEX) {
52 for (auto dd = dim_type; dd > 0; --dd) {
53 int nb_ents = moab::CN::NumSubEntities(type, dd);
54 for (int ii = 0; ii != nb_ents; ++ii) {
55 auto sub_ent_type = moab::CN::SubEntityType(type, dd, ii);
56 count[sub_ent_type] = nb_ents;
57 }
58 for (auto tt = moab::CN::TypeDimensionMap[dd].first;
59 tt <= moab::CN::TypeDimensionMap[dd].second; ++tt) {
60 data->dataOnEntities[tt].resize(count[tt]);
61 }
62 }
63 }
64 };
65
66 switch (type) {
67 case MBENTITYSET:
68 break;
69
70 default:
71 set_default();
72 }
73}
74
76 constructor_data(this, type);
77}
78
81 constructor_data(this, type);
83}
84
85static void
87 const boost::shared_ptr<EntitiesFieldData> &data_ptr) {
88
90 using DerivedEntData = DerivedEntitiesFieldData::DerivedEntData;
91
92 for (int tt = MBVERTEX; tt != MBMAXTYPE; ++tt) {
93 auto &ent_data = data_ptr->dataOnEntities[tt];
94 auto &derived_ent_data = derived_data->dataOnEntities[tt];
95 for (auto c = derived_ent_data.size(); c < ent_data.size(); ++c) {
96 boost::shared_ptr<EntData> ent_data_ptr(data_ptr, &ent_data[c]);
97 derived_ent_data.push_back(new DerivedEntData(ent_data_ptr));
98 }
99 derived_ent_data.resize(ent_data.size());
100 }
101}
102
104 const boost::shared_ptr<EntitiesFieldData> &data_ptr)
105 : EntitiesFieldData(), dataPtr(data_ptr) {
107}
108
113}
114
117 sPace = NOSPACE;
118 bAse = NOBASE;
119 fieldEntities.resize(0, false);
120 iNdices.resize(0, false);
121 localIndices.resize(0, false);
122 dOfs.resize(0, false);
123 fieldData.resize(0, false);
125}
126
129 for (EntityType t = MBVERTEX; t != MBMAXTYPE; t++)
130 for (auto &e : dataOnEntities[t])
131 CHKERR e.resetFieldDependentData();
133}
134
137 const FieldApproximationBase base) {
139 auto make_swap = [](boost::shared_ptr<MatrixDouble> &ptr,
140 boost::shared_ptr<MatrixDouble> &ptrBB,
141 boost::shared_ptr<MatrixDouble> &swap_ptr) {
142 if (swap_ptr) {
143 ptr = swap_ptr;
144 swap_ptr.reset();
145 } else {
146 swap_ptr = ptr;
147 ptr = ptrBB;
148 }
149 };
150 make_swap(getNSharedPtr(base), getBBNSharedPtr(field_name), swapBaseNPtr);
151 make_swap(getDiffNSharedPtr(base), getBBDiffNSharedPtr(field_name),
152 swapBaseDiffNPtr);
154}
155
157 const FieldApproximationBase base) {
159 // Note: Do not swap bases on entities sets
160 for (int tt = MBVERTEX; tt != MBENTITYSET; ++tt) {
161 auto &ent_data = dataOnEntities[tt];
162 for (auto &side_data : ent_data)
163 CHKERR side_data.baseSwap(field_name, base);
164 }
166}
167
169 const std::string &field_name, const FieldApproximationBase base) {
171 auto make_swap = [](boost::shared_ptr<MatrixDouble> &ptr,
172 boost::shared_ptr<MatrixDouble> &ptrBB,
173 boost::shared_ptr<MatrixDouble> &swap_ptr) {
174 if (swap_ptr) {
175 ptr = swap_ptr;
176 swap_ptr.reset();
177 } else {
178 swap_ptr = ptr;
179 ptr = ptrBB;
180 }
181 };
182 make_swap(getDerivedNSharedPtr(base), getBBNSharedPtr(field_name),
183 swapBaseNPtr);
184 make_swap(getDerivedDiffNSharedPtr(base), getBBDiffNSharedPtr(field_name),
185 swapBaseDiffNPtr);
187}
188
190 const boost::shared_ptr<EntitiesFieldData::EntData> &ent_data_ptr)
191 : EntitiesFieldData::EntData(false), entDataPtr(ent_data_ptr) {}
192
194 return entDataPtr->getSense();
195}
196
197boost::shared_ptr<MatrixDouble> &
199 const FieldApproximationBase base, const BaseDerivatives derivative) {
200 if (baseFunctionsAndBaseDerivatives[derivative][base])
201 return baseFunctionsAndBaseDerivatives[derivative][base];
202 else
203 return entDataPtr->getNSharedPtr(base, derivative);
204}
205
206boost::shared_ptr<MatrixDouble> &
208 const FieldApproximationBase base) {
209 if (N[base])
210 return N[base];
211 else
212 return entDataPtr->getNSharedPtr(base);
213}
214
215boost::shared_ptr<MatrixDouble> &
217 const FieldApproximationBase base) {
218 if (diffN[base])
219 return diffN[base];
220 else
221 return entDataPtr->getDiffNSharedPtr(base);
222}
223const boost::shared_ptr<MatrixDouble> &
225 const FieldApproximationBase base) const {
226 if (N[base])
227 return N[base];
228 else
229 return entDataPtr->getNSharedPtr(base);
230}
231const boost::shared_ptr<MatrixDouble> &
233 const FieldApproximationBase base) const {
234 if (diffN[base])
235 return diffN[base];
236 else
237 return entDataPtr->getDiffNSharedPtr(base);
238}
239
240std::ostream &operator<<(std::ostream &os,
242 os << "sEnse: " << e.getSense() << std::endl
243 << "oRder: " << e.getOrder() << std::endl
244 << "global indices: " << e.getIndices() << std::endl
245 << "local indices: " << e.getLocalIndices() << std::endl;
246 // FIXME: precision should not be set here
247 os << "fieldData: " << std::fixed << std::setprecision(2) << e.getFieldData()
248 << std::endl;
249 MatrixDouble base = const_cast<EntitiesFieldData::EntData &>(e).getN();
250 MatrixDouble diff_base =
251 const_cast<EntitiesFieldData::EntData &>(e).getDiffN();
252 const double eps = 1e-6;
253 for (unsigned int ii = 0; ii != base.size1(); ii++) {
254 for (unsigned int jj = 0; jj != base.size2(); jj++) {
255 if (fabs(base(ii, jj)) < eps)
256 base(ii, jj) = 0;
257 }
258 }
259 for (unsigned int ii = 0; ii != diff_base.size1(); ii++) {
260 for (unsigned int jj = 0; jj != diff_base.size2(); jj++) {
261 if (fabs(diff_base(ii, jj)) < eps)
262 diff_base(ii, jj) = 0;
263 }
264 }
265 os << "N: " << std::fixed << base << std::endl
266 << "diffN: " << std::fixed << diff_base;
267 return os;
268}
269
270std::ostream &operator<<(std::ostream &os, const EntitiesFieldData &e) {
271 for (EntityType t = MBVERTEX; t != MBMAXTYPE; ++t) {
272 for (unsigned int nn = 0; nn < e.dataOnEntities[t].size(); nn++) {
273 os << "dataOnEntities[" << moab::CN::EntityTypeName(t) << "][" << nn
274 << "]" << std::endl
275 << e.dataOnEntities[t][nn] << std::endl;
276 }
277 }
278 return os;
279}
280
281/** \name Specializations for H1/L2 */
282
283/**@{*/
284
285template <>
287EntitiesFieldData::EntData::getFTensor1FieldData<3>() {
288 if (dOfs[0]->getNbOfCoeffs() != 3) {
289 std::stringstream s;
290 s << "Wrong number of coefficients is " << dOfs[0]->getNbOfCoeffs();
291 s << " but you ask for tensor rank 1 dimension 3";
292 THROW_MESSAGE(s.str());
293 }
294 double *ptr = &*fieldData.data().begin();
296 &ptr[2]);
297}
298
299template <>
301EntitiesFieldData::EntData::getFTensor1FieldData<2>() {
302 if (dOfs[0]->getNbOfCoeffs() != 2) {
303 std::stringstream s;
304 s << "Wrong number of coefficients is " << dOfs[0]->getNbOfCoeffs();
305 s << " but you ask for tensor rank 1 dimension 2";
306 THROW_MESSAGE(s.str());
307 }
308 double *ptr = &*fieldData.data().begin();
309 return FTensor::Tensor1<FTensor::PackPtr<double *, 2>, 2>(ptr, &ptr[1]);
310}
311
312template <>
314EntitiesFieldData::EntData::getFTensor1FieldData<1>() {
315 if (dOfs[0]->getNbOfCoeffs() != 1) {
316 std::stringstream s;
317 s << "Wrong number of coefficients is " << dOfs[0]->getNbOfCoeffs();
318 s << " but you ask for tensor rank 1 dimension 1";
319 THROW_MESSAGE(s.str());
320 }
321 double *ptr = &*fieldData.data().begin();
323}
324
325template <>
327EntitiesFieldData::EntData::getFTensor2FieldData<3, 3>() {
328 if (dOfs[0]->getNbOfCoeffs() != 9) {
329 std::stringstream s;
330 s << "Wrong number of coefficients is " << dOfs[0]->getNbOfCoeffs();
331 s << " but you ask for tensor rank 2 dimensions 3 by 3 so 9 coefficients "
332 "is expected";
333 THROW_MESSAGE(s.str());
334 }
335 double *ptr = &*fieldData.data().begin();
337 ptr, &ptr[1], &ptr[2], &ptr[3], &ptr[4], &ptr[5], &ptr[6], &ptr[7],
338 &ptr[8]);
339}
340
341template <>
343EntitiesFieldData::EntData::getFTensor2SymmetricFieldData<3>() {
344 if (dOfs[0]->getNbOfCoeffs() != 6) {
345 std::stringstream s;
346 s << "Wrong number of coefficients is " << dOfs[0]->getNbOfCoeffs();
347 s << " but you ask for symmetric tensor rank 2 dimensions 3 by 3 so 6 "
348 "coefficients "
349 "is expected";
350 THROW_MESSAGE(s.str());
351 }
352 double *ptr = &*fieldData.data().begin();
354 ptr, &ptr[1], &ptr[2], &ptr[3], &ptr[4], &ptr[5]);
355}
356
357template <>
359EntitiesFieldData::EntData::getFTensor2SymmetricFieldData<2>() {
360 if (dOfs[0]->getNbOfCoeffs() != 3) {
361 std::stringstream s;
362 s << "Wrong number of coefficients is " << dOfs[0]->getNbOfCoeffs();
363 s << " but you ask for symmetric tensor rank 2 dimensions 2 by 2 so 3 "
364 "coefficients "
365 "is expected";
366 THROW_MESSAGE(s.str());
367 }
368 double *ptr = &*fieldData.data().begin();
370 ptr, &ptr[1], &ptr[2]);
371}
372
375 if (dOfs[0]->getNbOfCoeffs() != 1) {
376 std::stringstream s;
377 s << "Wrong number of coefficients is " << dOfs[0]->getNbOfCoeffs();
378 s << " but expected scalar field, tensor of rank 0";
379 THROW_MESSAGE(s.str());
380 }
382 &*fieldData.data().begin());
383}
384
385template <int Tensor_Dim>
388 const FieldApproximationBase base) {
389 std::stringstream s;
390 s << "Template for tensor dimension " << Tensor_Dim << " not implemented";
391 THROW_MESSAGE(s.str());
393}
394
395template <int Tensor_Dim>
398 return getFTensor1DiffN<Tensor_Dim>(bAse);
399}
400
401/**
402 * \brief Get spatial derivative of base function tensor for dimension 3d
403 */
404template <>
406EntitiesFieldData::EntData::getFTensor1DiffN<3>(
407 const FieldApproximationBase base) {
408 double *ptr = &*getDiffN(base).data().begin();
410 &ptr[2]);
411}
412
413/**
414 * \brief Get spatial derivative of base function tensor for dimension 3d
415 */
416template <>
418EntitiesFieldData::EntData::getFTensor1DiffN<3>() {
419 return getFTensor1DiffN<3>(bAse);
420}
421
422/**
423 * \brief Get spatial derivative of base function tensor for dimension 2d
424 */
425template <>
427EntitiesFieldData::EntData::getFTensor1DiffN<2>(
428 const FieldApproximationBase base) {
429 double *ptr = &*getDiffN(base).data().begin();
430 return FTensor::Tensor1<FTensor::PackPtr<double *, 2>, 2>(ptr, &ptr[1]);
431}
432
433/**
434 * \brief Get spatial derivative of base function tensor for dimension 2d
435 */
436template <>
438EntitiesFieldData::EntData::getFTensor1DiffN<2>() {
439 return getFTensor1DiffN<2>(bAse);
440}
441
442template <int Tensor_Dim>
445 const int gg, const int bb) {
446 std::stringstream s;
447 s << "Template for tensor dimension " << Tensor_Dim << " not implemented";
448 THROW_MESSAGE(s.str());
450}
451
452/**
453 * \brief Get spatial derivative of base function tensor for dimension 3d
454 */
455template <>
457EntitiesFieldData::EntData::getFTensor1DiffN<3>(
458 const FieldApproximationBase base, const int gg, const int bb) {
459 double *ptr = &getDiffN(base)(gg, 3 * bb);
461 &ptr[2]);
462}
463
464/**
465 * \brief Get spatial derivative of base function tensor for dimension 3d
466 */
467template <>
469EntitiesFieldData::EntData::getFTensor1DiffN<3>(const int gg, const int bb) {
470 return getFTensor1DiffN<3>(bAse, gg, bb);
471}
472
473/**
474 * \brief Get spatial derivative of base function tensor for dimension 2d
475 */
476template <>
478EntitiesFieldData::EntData::getFTensor1DiffN<2>(
479 const FieldApproximationBase base, const int gg, const int bb) {
480 double *ptr = &getDiffN(base)(gg, 2 * bb);
481 return FTensor::Tensor1<FTensor::PackPtr<double *, 2>, 2>(ptr, &ptr[1]);
482}
483
484/**
485 * \brief Get spatial derivative of base function tensor for dimension 2d
486 */
487template <>
489EntitiesFieldData::EntData::getFTensor1DiffN<2>(const int gg, const int bb) {
490 return getFTensor1DiffN<2>(bAse, gg, bb);
491}
492
493/**@}*/
494
495/** \name Specializations for HDiv/HCrul */
496
497/**@{*/
498
499template <int Tensor_Dim>
502 std::stringstream s;
503 s << "Template for tensor dimension " << Tensor_Dim << " not implemented";
504 THROW_MESSAGE(s.str());
506}
507
508template <int Tensor_Dim>
511 const int gg, const int bb) {
512 std::stringstream s;
513 s << "Template for tensor dimension " << Tensor_Dim << " not implemented";
514 THROW_MESSAGE(s.str());
516}
517
518template <int Tensor_Dim0, int Tensor_Dim1>
520 Tensor_Dim0, Tensor_Dim1>
522 std::stringstream s;
523 s << "Template for tensor dimension " << Tensor_Dim0 << "x" << Tensor_Dim1
524 << " not implemented";
525 THROW_MESSAGE(s.str());
527}
528
529template <int Tensor_Dim0, int Tensor_Dim1>
531 Tensor_Dim0, Tensor_Dim1>
533 const int gg, const int bb) {
534 std::stringstream s;
535 s << "Template for tensor dimension " << Tensor_Dim0 << "x" << Tensor_Dim1
536 << " not implemented";
537 THROW_MESSAGE(s.str());
539}
540
541template <>
543EntitiesFieldData::EntData::getFTensor1N<3>(FieldApproximationBase base) {
544 double *t_n_ptr = &*getN(base).data().begin();
545 return FTensor::Tensor1<FTensor::PackPtr<double *, 3>, 3>(t_n_ptr, // HVEC0
546 &t_n_ptr[HVEC1],
547 &t_n_ptr[HVEC2]);
548}
549
550template <>
552EntitiesFieldData::EntData::getFTensor1N<3>(FieldApproximationBase base,
553 const int gg, const int bb) {
554 double *t_n_ptr = &getN(base)(gg, 3 * bb);
555 return FTensor::Tensor1<FTensor::PackPtr<double *, 3>, 3>(t_n_ptr, // HVEC0
556 &t_n_ptr[HVEC1],
557 &t_n_ptr[HVEC2]);
558}
559
560template <>
562EntitiesFieldData::EntData::getFTensor2DiffN<3, 3>(
564 double *t_diff_n_ptr = &*getDiffN(base).data().begin();
566 t_diff_n_ptr, &t_diff_n_ptr[HVEC0_1], &t_diff_n_ptr[HVEC0_2],
567 &t_diff_n_ptr[HVEC1_0], &t_diff_n_ptr[HVEC1_1], &t_diff_n_ptr[HVEC1_2],
568 &t_diff_n_ptr[HVEC2_0], &t_diff_n_ptr[HVEC2_1], &t_diff_n_ptr[HVEC2_2]);
569}
570
571template <>
573EntitiesFieldData::EntData::getFTensor2DiffN<3, 3>(FieldApproximationBase base,
574 const int gg, const int bb) {
575 double *t_diff_n_ptr = &getDiffN(base)(gg, 9 * bb);
577 t_diff_n_ptr, &t_diff_n_ptr[HVEC0_1], &t_diff_n_ptr[HVEC0_2],
578 &t_diff_n_ptr[HVEC1_0], &t_diff_n_ptr[HVEC1_1], &t_diff_n_ptr[HVEC1_2],
579 &t_diff_n_ptr[HVEC2_0], &t_diff_n_ptr[HVEC2_1], &t_diff_n_ptr[HVEC2_2]);
580}
581
582template <>
584EntitiesFieldData::EntData::getFTensor2DiffN<3, 2>(
586 double *t_diff_n_ptr = &*getDiffN(base).data().begin();
588 t_diff_n_ptr, &t_diff_n_ptr[HVEC0_1], &t_diff_n_ptr[HVEC1_0],
589 &t_diff_n_ptr[HVEC1_1], &t_diff_n_ptr[HVEC2_0], &t_diff_n_ptr[HVEC2_1]);
590}
591
592template <>
594EntitiesFieldData::EntData::getFTensor2DiffN<3, 2>(FieldApproximationBase base,
595 const int gg, const int bb) {
596 double *t_diff_n_ptr = &getDiffN(base)(gg, 6 * bb);
598 t_diff_n_ptr, &t_diff_n_ptr[HVEC0_1], &t_diff_n_ptr[HVEC1_0],
599 &t_diff_n_ptr[HVEC1_1], &t_diff_n_ptr[HVEC2_0], &t_diff_n_ptr[HVEC2_1]);
600}
601
602template <>
605 double *ptr = &(getN(base, BaseDerivatives::SecondDerivative))(0, 0);
607
608 &ptr[2 * HVEC0_0 + 0], &ptr[2 * HVEC0_0 + 1], &ptr[2 * HVEC0_1 + 0],
609 &ptr[2 * HVEC0_1 + 1],
610
611 &ptr[2 * HVEC1_0 + 0], &ptr[2 * HVEC1_0 + 1], &ptr[2 * HVEC1_1 + 0],
612 &ptr[2 * HVEC1_1 + 1],
613
614 &ptr[2 * HVEC2_0 + 0], &ptr[2 * HVEC2_0 + 1], &ptr[2 * HVEC2_1 + 0],
615 &ptr[2 * HVEC2_1 + 1]
616
617 };
618}
619
620template <int Tensor_Dim0, int Tensor_Dim1>
622 Tensor_Dim0, Tensor_Dim1>
624 std::stringstream s;
625 s << "Template for tensor dimension " << Tensor_Dim0 << ", " << Tensor_Dim1
626 << " not implemented";
627 THROW_MESSAGE(s.str());
629 Tensor_Dim0, Tensor_Dim1>();
630}
631
632template <>
634EntitiesFieldData::EntData::getFTensor2N<3, 3>(FieldApproximationBase base) {
635 double *t_n_ptr = &*(getN(base).data().begin());
637
638 &t_n_ptr[HVEC0], &t_n_ptr[HVEC1], &t_n_ptr[HVEC2],
639
640 &t_n_ptr[3 + HVEC0], &t_n_ptr[3 + HVEC1], &t_n_ptr[3 + HVEC2],
641
642 &t_n_ptr[6 + HVEC0], &t_n_ptr[6 + HVEC1], &t_n_ptr[6 + HVEC2]
643
644 );
645}
646
647template <int Tensor_Dim0, int Tensor_Dim1>
649 Tensor_Dim0, Tensor_Dim1>
651 const int gg, const int bb) {
652 std::stringstream s;
653 s << "Template for tensor dimension " << Tensor_Dim0 << ", " << Tensor_Dim1
654 << " not implemented";
655 THROW_MESSAGE(s.str());
657 Tensor_Dim0, Tensor_Dim1>();
658}
659
660template <>
662EntitiesFieldData::EntData::getFTensor2N<3, 3>(FieldApproximationBase base,
663 const int gg, const int bb) {
664 double *t_n_ptr = &getN(base)(gg, 9 * bb);
666
667 &t_n_ptr[HVEC0], &t_n_ptr[HVEC1], &t_n_ptr[HVEC2],
668
669 &t_n_ptr[3 + HVEC0], &t_n_ptr[3 + HVEC1], &t_n_ptr[3 + HVEC2],
670
671 &t_n_ptr[6 + HVEC0], &t_n_ptr[6 + HVEC1], &t_n_ptr[6 + HVEC2]
672
673 );
674}
675
676/**@}*/
677
678/** \name Bernstein-Bezier base only functions */
679
680/**@{*/
681
682boost::shared_ptr<MatrixInt> &
684 const std::string &field_name) {
685 return bbAlphaIndices[field_name];
686}
687
688boost::shared_ptr<MatrixDouble> &
690 return bbN[field_name];
691}
692
693/**
694 * Get shared pointer to BB base base functions
695 */
696const boost::shared_ptr<MatrixDouble> &
698 const std::string &field_name) const {
699 return bbN.at(field_name);
700}
701
702/**
703 * Get shared pointer to BB derivatives of base base functions
704 */
705boost::shared_ptr<MatrixDouble> &
707 return bbDiffN[field_name];
708}
709
710/**
711 * Get shared pointer to derivatives of BB base base functions
712 */
713const boost::shared_ptr<MatrixDouble> &
715 const std::string &field_name) const {
716 return bbDiffN.at(field_name);
717}
718
719std::map<std::string, boost::shared_ptr<MatrixInt>> &
721 return bbAlphaIndices;
722}
723
724std::map<std::string, boost::shared_ptr<MatrixDouble>> &
726 return bbN;
727}
728
729std::map<std::string, boost::shared_ptr<MatrixDouble>> &
731 return bbDiffN;
732}
733
734boost::shared_ptr<MatrixInt> &
736 return bbAlphaIndicesByOrder[o];
737}
738
739boost::shared_ptr<MatrixDouble> &
741 return bbNByOrder[o];
742}
743
744boost::shared_ptr<MatrixDouble> &
746 return bbDiffNByOrder[o];
747}
748
749std::array<boost::shared_ptr<MatrixInt>,
752 return bbAlphaIndicesByOrder;
753}
754
755std::array<boost::shared_ptr<MatrixDouble>,
758 return bbNByOrder;
759}
760
761std::array<boost::shared_ptr<MatrixDouble>,
764 return bbDiffNByOrder;
765}
766
767boost::shared_ptr<MatrixInt> &
769 const std::string &field_name) {
770 return entDataPtr->getBBAlphaIndicesSharedPtr(field_name);
771}
772
773boost::shared_ptr<MatrixDouble> &
775 const std::string &field_name) {
776 return entDataPtr->getBBNSharedPtr(field_name);
777}
778
779const boost::shared_ptr<MatrixDouble> &
781 const std::string &field_name) const {
782 return entDataPtr->getBBNSharedPtr(field_name);
783}
784
785/**
786 * Get shared pointer to BB derivatives of base base functions
787 */
788boost::shared_ptr<MatrixDouble> &
790 const std::string &field_name) {
791 return entDataPtr->getBBDiffNSharedPtr(field_name);
792}
793
794/**
795 * Get shared pointer to derivatives of BB base base functions
796 */
797const boost::shared_ptr<MatrixDouble> &
799 const std::string &field_name) const {
800 return entDataPtr->getBBDiffNSharedPtr(field_name);
801}
802
803/**@}*/
804
806 return entDataBitRefLevel;
807}
808
809std::vector<BitRefLevel> &
811 return entDataPtr->getEntDataBitRefLevel();
812}
813
814} // namespace MoFEM
static Index< 'o', 3 > o
static const double eps
T data[Tensor_Dim]
T data[Tensor_Dim0][Tensor_Dim1]
FieldApproximationBase
approximation base
Definition: definitions.h:58
@ LASTBASE
Definition: definitions.h:69
@ NOBASE
Definition: definitions.h:59
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:447
@ NOSPACE
Definition: definitions.h:83
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
@ HVEC0
Definition: definitions.h:186
@ HVEC1
Definition: definitions.h:186
@ HVEC2
Definition: definitions.h:186
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
@ HVEC1_1
Definition: definitions.h:196
@ HVEC0_1
Definition: definitions.h:195
@ HVEC1_0
Definition: definitions.h:193
@ HVEC2_1
Definition: definitions.h:197
@ HVEC1_2
Definition: definitions.h:199
@ HVEC2_2
Definition: definitions.h:200
@ HVEC2_0
Definition: definitions.h:194
@ HVEC0_2
Definition: definitions.h:198
@ HVEC0_0
Definition: definitions.h:192
#define CHKERR
Inline error check.
Definition: definitions.h:535
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:440
#define THROW_MESSAGE(msg)
Throw MoFEM exception.
Definition: definitions.h:561
const double c
speed of light (cm/ns)
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
UBlasMatrix< double > MatrixDouble
Definition: Types.hpp:77
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
static void constructor_derived_data(DerivedEntitiesFieldData *derived_data, const boost::shared_ptr< EntitiesFieldData > &data_ptr)
static void constructor_data(EntitiesFieldData *data, const EntityType type)
constexpr double t
plate stiffness
Definition: plate.cpp:59
constexpr auto field_name
const int N
Definition: speed_test.cpp:3
Derived ata on single entity (This is passed as argument to DataOperator::doWork)
boost::shared_ptr< MatrixDouble > & getBBNSharedPtr(const std::string &field_name)
MoFEMErrorCode baseSwap(const std::string &field_name, const FieldApproximationBase base)
Used by Bernstein base to keep temporally pointer.
std::vector< BitRefLevel > & getEntDataBitRefLevel()
int getSense() const
get entity sense, need to calculate base functions with conforming approximation fields
boost::shared_ptr< MatrixInt > & getBBAlphaIndicesSharedPtr(const std::string &field_name)
boost::shared_ptr< MatrixDouble > & getBBDiffNSharedPtr(const std::string &field_name)
boost::shared_ptr< MatrixDouble > & getDiffNSharedPtr(const FieldApproximationBase base)
DerivedEntData(const boost::shared_ptr< EntitiesFieldData::EntData > &ent_data_ptr)
boost::shared_ptr< MatrixDouble > & getNSharedPtr(const FieldApproximationBase base, const BaseDerivatives derivative)
this class derive data form other data structure
MoFEMErrorCode setElementType(const EntityType type)
const boost::shared_ptr< EntitiesFieldData > dataPtr
DerivedEntitiesFieldData(const boost::shared_ptr< EntitiesFieldData > &data_ptr)
Data on single entity (This is passed as argument to DataOperator::doWork)
virtual std::array< boost::shared_ptr< MatrixDouble >, MaxBernsteinBezierOrder > & getBBDiffNByOrderArray()
virtual std::array< boost::shared_ptr< MatrixDouble >, MaxBernsteinBezierOrder > & getBBNByOrderArray()
auto getFTensor2N()
Get base functions for Hdiv space.
VectorDouble fieldData
Field data on entity.
ApproximationOrder getOrder() const
get approximation order
VectorInt localIndices
Local indices on entity.
virtual std::array< boost::shared_ptr< MatrixInt >, MaxBernsteinBezierOrder > & getBBAlphaIndicesByOrderArray()
virtual std::map< std::string, boost::shared_ptr< MatrixInt > > & getBBAlphaIndicesMap()
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1DiffN()
Get derivatives of base functions.
static constexpr size_t MaxBernsteinBezierOrder
VectorInt iNdices
Global indices on entity.
virtual boost::shared_ptr< MatrixInt > & getBBAlphaIndicesSharedPtr(const std::string &field_name)
const VectorDouble & getFieldData() const
get dofs values
virtual boost::shared_ptr< MatrixDouble > & getBBDiffNByOrderSharedPtr(const size_t o)
get BB base derivative by order
virtual std::map< std::string, boost::shared_ptr< MatrixDouble > > & getBBNMap()
get hash map of base function for BB base, key is a field name
virtual MoFEMErrorCode baseSwap(const std::string &field_name, const FieldApproximationBase base)
Swap bases functions.
const VectorInt & getLocalIndices() const
get local indices of dofs on entity
virtual std::map< std::string, boost::shared_ptr< MatrixDouble > > & getBBDiffNMap()
get hash map of direvarives base function for BB base, key is a field name
virtual std::vector< BitRefLevel > & getEntDataBitRefLevel()
virtual int getSense() const
get entity sense, need to calculate base functions with conforming approximation fields
VectorFieldEntities fieldEntities
Field entities.
FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > getFTensor0FieldData()
Resturn scalar files as a FTensor of rank 0.
virtual boost::shared_ptr< MatrixInt > & getBBAlphaIndicesByOrderSharedPtr(const size_t o)
get ALpha indices for BB base by order
FTensor::Tensor3< FTensor::PackPtr< double *, Tensor_Dim0 *Tensor_Dim1 *Tensor_Dim2 >, Tensor_Dim0, Tensor_Dim1, Tensor_Dim2 > getFTensor3Diff2N()
Get second derivatives of base functions for Hvec space.
auto getFTensor1N()
Get base functions for Hdiv space.
virtual boost::shared_ptr< MatrixDouble > & getBBNSharedPtr(const std::string &field_name)
std::array< std::array< boost::shared_ptr< MatrixDouble >, LASTBASE >, LastDerivative > baseFunctionsAndBaseDerivatives
FTensor::Tensor2< FTensor::PackPtr< double *, Tensor_Dim0 *Tensor_Dim1 >, Tensor_Dim0, Tensor_Dim1 > getFTensor2DiffN()
Get derivatives of base functions for Hdiv space.
virtual boost::shared_ptr< MatrixDouble > & getDiffNSharedPtr(const FieldApproximationBase base)
virtual boost::shared_ptr< MatrixDouble > & getBBNByOrderSharedPtr(const size_t o)
get BB base by order
virtual boost::shared_ptr< MatrixDouble > & getNSharedPtr(const FieldApproximationBase base, const BaseDerivatives direvatie)
const VectorInt & getIndices() const
Get global indices of dofs on entity.
EntData(const bool allocate_base_matrices=true)
FieldApproximationBase bAse
Field approximation base.
virtual boost::shared_ptr< MatrixDouble > & getBBDiffNSharedPtr(const std::string &field_name)
data structure for finite element entity
friend std::ostream & operator<<(std::ostream &os, const EntitiesFieldData &e)
MoFEMErrorCode resetFieldDependentData()
std::bitset< LASTBASE > bAse
bases on element
std::array< boost::ptr_vector< EntData >, MBMAXTYPE > dataOnEntities
virtual MoFEMErrorCode baseSwap(const std::string &field_name, const FieldApproximationBase base)
Swap approximation base.
virtual MoFEMErrorCode setElementType(const EntityType type)