v0.15.5
Loading...
Searching...
No Matches
HODataOperators.cpp
Go to the documentation of this file.
1/** \file HODataOperators.cpp
2
3\brief Set of operators for high-order geometry approximation.
4
5*/
6
7namespace MoFEM {
8
10 boost::shared_ptr<MatrixDouble> jac_ptr)
12 jacPtr(jac_ptr) {
13
14 for (auto t = MBEDGE; t != MBMAXTYPE; ++t)
15 doEntities[t] = false;
16
17 if (!jacPtr)
18 THROW_MESSAGE("Jac pointer not allocated");
19}
20
22OpCalculateHOJacForVolume::doWork(int side, EntityType type,
25
26 const auto nb_base_functions = data.getN(NOBASE).size2();
27 if (nb_base_functions) {
28
29 const auto nb_gauss_pts = getGaussPts().size2();
30 auto t_diff_base = data.getFTensor1DiffN<3>(NOBASE);
31 auto &coords = getCoords();
32
33#ifndef NDEBUG
34 if (nb_gauss_pts != data.getDiffN().size1())
35 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
36 "Inconsistent number base functions and gauss points");
37 if (nb_base_functions != data.getDiffN().size2() / 3)
38 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
39 "Inconsistent number of base functions");
40 if (coords.size() != 3 * nb_base_functions)
41 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
42 "Number of vertex coordinates and base functions is inconsistent "
43 "%ld != %ld",
44 coords.size() / 3, nb_base_functions);
45#endif
46
47 jacPtr->resize(9, nb_gauss_pts, false);
48 jacPtr->clear();
49 auto t_jac = getFTensor2FromMat<3, 3>(*jacPtr);
50 auto t_vol_inv_jac = getInvJac();
51
52 FTensor::Index<'i', 3> i;
53 FTensor::Index<'j', 3> j;
54 FTensor::Index<'k', 3> k;
55
56 for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
58 &coords[0], &coords[1], &coords[2]);
59 for (size_t bb = 0; bb != nb_base_functions; ++bb) {
60 t_jac(i, j) += t_coords(i) * (t_vol_inv_jac(k, j) * t_diff_base(k));
61 ++t_diff_base;
62 ++t_coords;
63 }
64 ++t_jac;
65 }
66 }
67
70
75
76 if (getFEDim() == 3) {
77
78 auto transform_base = [&](MatrixDouble &diff_n) {
79 return applyTransform<3, 3, 3, 3>(diff_n);
80 };
81
82 for (int b = AINSWORTH_LEGENDRE_BASE; b != LASTBASE; b++) {
84 CHKERR transform_base(data.getDiffN(base));
85 }
86
87 switch (type) {
88 case MBVERTEX:
89 for (auto &m : data.getBBDiffNMap())
90 CHKERR transform_base(*(m.second));
91 break;
92 default:
93 for (auto &ptr : data.getBBDiffNByOrderArray())
94 if (ptr)
95 CHKERR transform_base(*ptr);
96 }
97
98 } else {
99 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "Use different operator");
100 }
101
103}
104
109
110 if (getFEDim() != 3)
111 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
112 "This operator can be used only with element with volume dimension");
113
114 auto apply_transform = [&](MatrixDouble &diff2_n) {
116#ifndef NDEBUG
117 if (diff2_n.size2() % 9)
118 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
119 "Inconsistent number of base functions for second derivatives");
120#endif
121 size_t nb_functions = diff2_n.size2() / 9;
122 if (nb_functions) {
123#ifndef NDEBUG
124 if (diff2_n.size1() != getGaussPts().size2())
125 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
126 "Wrong number of Gauss Pts %ld != %ld", diff2_n.size1(),
127 getGaussPts().size2());
128#endif
129
130 size_t nb_gauss_pts = diff2_n.size1();
131
132 diffNinvJac.resize(diff2_n.size1(), diff2_n.size2(), false);
133
134 FTensor::Index<'i', 3> i;
135 FTensor::Index<'j', 3> j;
136 FTensor::Index<'k', 3> k;
137 FTensor::Index<'l', 3> l;
138 auto t_diff2_n = getFTensor2FromPtr<3, 3>(&*diffNinvJac.data().begin());
139 auto t_diff2_n_ref = getFTensor2FromPtr<3, 3>(&*diff2_n.data().begin());
140 auto t_inv_jac = getFTensor2FromMat<3, 3>(*invJacPtr);
141 for (size_t gg = 0; gg != nb_gauss_pts; ++gg, ++t_inv_jac) {
142 for (size_t dd = 0; dd != nb_functions; ++dd) {
143 t_diff2_n(i, j) =
144 t_inv_jac(k, i) * (t_inv_jac(l, j) * t_diff2_n_ref(k, l));
145 ++t_diff2_n;
146 ++t_diff2_n_ref;
147 }
148 }
149
150 diff2_n.swap(diffNinvJac);
151 }
153 };
154
155 if (!(type == MBVERTEX && sPace == L2)) {
156
157 for (int b = AINSWORTH_LEGENDRE_BASE; b != USER_BASE; b++) {
158 FieldApproximationBase base = static_cast<FieldApproximationBase>(b);
159 CHKERR apply_transform(
160 data.getN(base, BaseDerivatives::SecondDerivative));
161 }
162
163 auto error = [&](auto &m) {
165 if (m.size2())
166 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
167 "Operator do not work for Bernstein-Bezier. This "
168 "functionality is "
169 "not implemented.");
171 };
172
173 switch (type) {
174 case MBVERTEX:
175 for (auto &m : data.getBBDiffNMap()) {
176 CHKERR error(*(m.second));
177 CHKERR apply_transform(*(m.second));
178 }
179 break;
180 default:
181 for (auto &ptr : data.getBBDiffNByOrderArray())
182 if (ptr) {
183 CHKERR error(*ptr);
184 CHKERR apply_transform(*ptr);
185 }
186 }
187 }
188
190}
191
193OpSetHOInvJacVectorBase::doWork(int side, EntityType type,
196
197 FTensor::Index<'i', 3> i;
198 FTensor::Index<'j', 3> j;
199 FTensor::Index<'k', 3> k;
200
201 for (int b = AINSWORTH_LEGENDRE_BASE; b != LASTBASE; b++) {
202
203 FieldApproximationBase base = static_cast<FieldApproximationBase>(b);
204
205 diffHdivInvJac.resize(data.getDiffN(base).size1(),
206 data.getDiffN(base).size2(), false);
207
208 unsigned int nb_gauss_pts = data.getDiffN(base).size1();
209 unsigned int nb_base_functions = data.getDiffN(base).size2() / 9;
210 if (nb_base_functions == 0)
211 continue;
212
213 auto t_diff_n = data.getFTensor2DiffN<3, 3>(base);
214 double *t_inv_diff_n_ptr = &*diffHdivInvJac.data().begin();
216 t_inv_diff_n_ptr, &t_inv_diff_n_ptr[HVEC0_1],
217 &t_inv_diff_n_ptr[HVEC0_2], &t_inv_diff_n_ptr[HVEC1_0],
218 &t_inv_diff_n_ptr[HVEC1_1], &t_inv_diff_n_ptr[HVEC1_2],
219 &t_inv_diff_n_ptr[HVEC2_0], &t_inv_diff_n_ptr[HVEC2_1],
220 &t_inv_diff_n_ptr[HVEC2_2]);
221 auto t_inv_jac = getFTensor2FromMat<3, 3>(*invJacPtr);
222 for (unsigned int gg = 0; gg != nb_gauss_pts; ++gg) {
223 for (unsigned int bb = 0; bb != nb_base_functions; ++bb) {
224 t_inv_diff_n(i, j) = t_inv_jac(k, j) * t_diff_n(i, k);
225 ++t_diff_n;
226 ++t_inv_diff_n;
227 }
228 ++t_inv_jac;
229 }
230
231 data.getDiffN(base).swap(diffHdivInvJac);
232 }
233
235}
236
240
241 /* Note that is different treatment to OpSetHOWeight, that is because volume
242 * element by default calculate HO jacobian on physical element, whereas face
243 * element uses reference element. */
244
245 const size_t nb_int_pts = getGaussPts().size2();
246 if (getNormalsAtGaussPts().size1()) {
247 if (getNormalsAtGaussPts().size1() == nb_int_pts) {
248 double a = getMeasure();
249 if (getNumeredEntFiniteElementPtr()->getEntType() == MBTRI)
250 a *= 2;
251 auto t_w = getFTensor0IntegrationWeight();
252 auto t_normal = getFTensor1NormalsAtGaussPts();
253 FTensor::Index<'i', 3> i;
254 for (size_t gg = 0; gg != nb_int_pts; ++gg) {
255 t_w *= sqrt(t_normal(i) * t_normal(i)) / a;
256 ++t_w;
257 ++t_normal;
258 }
259 } else {
260 SETERRQ(PETSC_COMM_SELF, MOFEM_IMPOSSIBLE_CASE,
261 "Number of rows in getNormalsAtGaussPts should be equal to "
262 "number of integration points, but is not, i.e. %ld != %ld",
263 getNormalsAtGaussPts().size1(), nb_int_pts);
264 }
265 }
267}
268
272 const size_t nb_int_pts = getGaussPts().size2();
273 if (getTangentAtGaussPts().size1()) {
274 if (getTangentAtGaussPts().size1() == nb_int_pts) {
275 double a = getMeasure();
276 auto t_w = getFTensor0IntegrationWeight();
277 auto t_tangent = getFTensor1TangentAtGaussPts();
278 FTensor::Index<'i', 3> i;
279 for (size_t gg = 0; gg != nb_int_pts; ++gg) {
280 t_w *= sqrt(t_tangent(i) * t_tangent(i)) / a;
281 ++t_w;
282 ++t_tangent;
283 }
284 } else {
285 SETERRQ(PETSC_COMM_SELF, MOFEM_IMPOSSIBLE_CASE,
286 "Number of rows in getTangentAtGaussPts should be equal to "
287 "number of integration points, but is not, i.e. %ld != %ld",
288 getTangentAtGaussPts().size1(), nb_int_pts);
289 }
290 }
292}
293
294MoFEMErrorCode OpSetHOWeights::doWork(int side, EntityType type,
297
298 const auto nb_integration_pts = detPtr->size();
299
300#ifndef NDEBUG
301 if (nb_integration_pts != getGaussPts().size2())
302 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
303 "Inconsistent number of data points");
304#endif
305
306 auto t_w = getFTensor0IntegrationWeight();
307 auto t_det = getFTensor0FromVec(*detPtr);
308 for (size_t gg = 0; gg != nb_integration_pts; ++gg) {
309 t_w *= t_det;
310 ++t_w;
311 ++t_det;
312 }
313
315}
316
321
322 FTensor::Index<'i', 3> i;
323 FTensor::Index<'j', 3> j;
324 FTensor::Index<'k', 3> k;
325
326#ifndef NDEBUG
327 if (!detPtr)
328 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
329 "Pointer for detPtr not allocated");
330 if (!jacPtr)
331 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
332 "Pointer for jacPtr not allocated");
333#endif
334
335 for (int b = AINSWORTH_LEGENDRE_BASE; b != LASTBASE; ++b) {
336
337 FieldApproximationBase base = static_cast<FieldApproximationBase>(b);
338
339 auto nb_gauss_pts = data.getNSharedPtr(base) ? data.getN(base).size1() : 0;
340 auto nb_base_functions =
341 data.getNSharedPtr(base) ? data.getN(base).size2() / 3 : 0;
342 auto nb_diff_base_functions =
343 data.getDiffNSharedPtr(base) ? data.getDiffN(base).size2() / 9 : 0;
344
345#ifndef NDEBUG
346 if (nb_diff_base_functions) {
347 if (nb_diff_base_functions != nb_base_functions)
348 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
349 "Wrong number base functions %ld != %ld",
350 nb_diff_base_functions, nb_base_functions);
351 if (data.getDiffN(base).size1() != nb_gauss_pts)
352 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
353 "Wrong number integration points");
354 }
355#endif
356
357 if (nb_gauss_pts && nb_base_functions) {
358
359 piolaN.resize(nb_gauss_pts, data.getN(base).size2(), false);
360
361 auto t_n = data.getFTensor1N<3>(base);
362 double *t_transformed_n_ptr = &*piolaN.data().begin();
364 t_transformed_n_ptr, // HVEC0
365 &t_transformed_n_ptr[HVEC1], &t_transformed_n_ptr[HVEC2]);
366
367 auto t_det = getFTensor0FromVec(*detPtr);
368 auto t_jac = getFTensor2FromMat<3, 3>(*jacPtr);
369
370 for (unsigned int gg = 0; gg != nb_gauss_pts; ++gg) {
371 for (unsigned int bb = 0; bb != nb_base_functions; ++bb) {
372 const double a = 1. / t_det;
373 t_transformed_n(i) = a * (t_jac(i, k) * t_n(k));
374 ++t_n;
375 ++t_transformed_n;
376 }
377 ++t_det;
378 ++t_jac;
379 }
380
381 data.getN(base).swap(piolaN);
382 }
383
384 if (nb_gauss_pts && nb_diff_base_functions) {
385
386 piolaDiffN.resize(nb_gauss_pts, data.getDiffN(base).size2(), false);
387
388 auto t_diff_n = data.getFTensor2DiffN<3, 3>(base);
389 double *t_transformed_diff_n_ptr = &*piolaDiffN.data().begin();
391 t_transformed_diff_n(t_transformed_diff_n_ptr,
392 &t_transformed_diff_n_ptr[HVEC0_1],
393 &t_transformed_diff_n_ptr[HVEC0_2],
394 &t_transformed_diff_n_ptr[HVEC1_0],
395 &t_transformed_diff_n_ptr[HVEC1_1],
396 &t_transformed_diff_n_ptr[HVEC1_2],
397 &t_transformed_diff_n_ptr[HVEC2_0],
398 &t_transformed_diff_n_ptr[HVEC2_1],
399 &t_transformed_diff_n_ptr[HVEC2_2]);
400
401 auto t_det = getFTensor0FromVec(*detPtr);
402 for (unsigned int gg = 0; gg != nb_gauss_pts; ++gg) {
403 for (unsigned int bb = 0; bb != nb_diff_base_functions; ++bb) {
404 const double a = 1. / t_det;
405 t_transformed_diff_n(i, k) = a * t_diff_n(i, k);
406 ++t_diff_n;
407 ++t_transformed_diff_n;
408 }
409 ++t_det;
410 }
411
412 data.getDiffN(base).swap(piolaDiffN);
413 }
414 }
415
417}
418
423
424 FTensor::Index<'i', 3> i;
425 FTensor::Index<'j', 3> j;
426 FTensor::Index<'k', 3> k;
427
428 for (int b = AINSWORTH_LEGENDRE_BASE; b != LASTBASE; b++) {
429
430 FieldApproximationBase base = static_cast<FieldApproximationBase>(b);
431
432 unsigned int nb_gauss_pts = data.getN(base).size1();
433 unsigned int nb_base_functions = data.getN(base).size2() / 3;
434 piolaN.resize(nb_gauss_pts, data.getN(base).size2(), false);
435 piolaDiffN.resize(nb_gauss_pts, data.getDiffN(base).size2(), false);
436
437 auto t_n = data.getFTensor1N<3>(base);
438 double *t_transformed_n_ptr = &*piolaN.data().begin();
440 t_transformed_n_ptr, // HVEC0
441 &t_transformed_n_ptr[HVEC1], &t_transformed_n_ptr[HVEC2]);
442 auto t_diff_n = data.getFTensor2DiffN<3, 3>(base);
443 double *t_transformed_diff_n_ptr = &*piolaDiffN.data().begin();
444 FTensor::Tensor2<FTensor::PackPtr<double *, 9>, 3, 3> t_transformed_diff_n(
445 t_transformed_diff_n_ptr, &t_transformed_diff_n_ptr[HVEC0_1],
446 &t_transformed_diff_n_ptr[HVEC0_2], &t_transformed_diff_n_ptr[HVEC1_0],
447 &t_transformed_diff_n_ptr[HVEC1_1], &t_transformed_diff_n_ptr[HVEC1_2],
448 &t_transformed_diff_n_ptr[HVEC2_0], &t_transformed_diff_n_ptr[HVEC2_1],
449 &t_transformed_diff_n_ptr[HVEC2_2]);
450
451 auto t_inv_jac = getFTensor2FromMat<3, 3>(*jacInvPtr);
452
453 for (unsigned int gg = 0; gg != nb_gauss_pts; ++gg) {
454 for (unsigned int bb = 0; bb != nb_base_functions; ++bb) {
455 t_transformed_n(i) = t_inv_jac(k, i) * t_n(k);
456 t_transformed_diff_n(i, k) = t_inv_jac(j, i) * t_diff_n(j, k);
457 ++t_n;
458 ++t_transformed_n;
459 ++t_diff_n;
460 ++t_transformed_diff_n;
461 }
462 ++t_inv_jac;
463 }
464
465 data.getN(base).swap(piolaN);
466 data.getDiffN(base).swap(piolaDiffN);
467 }
468
470}
471
473 boost::shared_ptr<MatrixDouble> jac_ptr)
475 jacPtr(jac_ptr) {}
476
480
482
483 const auto nb_gauss_pts = getGaussPts().size2();
484
485 jacPtr->resize(4, nb_gauss_pts, false);
486 auto t_jac = getFTensor2FromMat<2, 2>(*jacPtr);
487
488 FTensor::Index<'i', 2> i;
489 FTensor::Index<'j', 2> j;
490
491 auto t_t1 = getFTensor1Tangent1AtGaussPts();
492 auto t_t2 = getFTensor1Tangent2AtGaussPts();
495
496 for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
497 t_jac(i, N0) = t_t1(i);
498 t_jac(i, N1) = t_t2(i);
499 ++t_t1;
500 ++t_t2;
501 ++t_jac;
502 }
503
505};
506
511
512 FTensor::Index<'i', 3> i;
513 FTensor::Index<'j', 3> j;
514 FTensor::Index<'k', 3> k;
515
516 size_t nb_gauss_pts = getGaussPts().size2();
517 jacPtr->resize(9, nb_gauss_pts, false);
518
519 auto t_t1 = getFTensor1Tangent1AtGaussPts();
520 auto t_t2 = getFTensor1Tangent2AtGaussPts();
521 auto t_normal = getFTensor1NormalsAtGaussPts();
522
526
527 auto t_jac = getFTensor2FromMat<3, 3>(*jacPtr);
528 for (size_t gg = 0; gg != nb_gauss_pts; ++gg, ++t_jac) {
529
530 t_jac(i, N0) = t_t1(i);
531 t_jac(i, N1) = t_t2(i);
532
533 const double l = sqrt(t_normal(j) * t_normal(j));
534
535 t_jac(i, N2) = t_normal(i) / l;
536
537 ++t_t1;
538 ++t_t2;
539 ++t_normal;
540 }
541
543}
544
546 int side, EntityType type, EntitiesFieldData::EntData &data) {
547 FTensor::Index<'i', 3> i;
549
550 if (moab::CN::Dimension(type) != 2)
552
553 auto get_normals_ptr = [&]() {
555 return &*normalsAtGaussPts->data().begin();
556 else
557 return &*getNormalsAtGaussPts().data().begin();
558 };
559
560 auto apply_transform_nonlinear_geometry = [&](auto base, auto nb_gauss_pts,
561 auto nb_base_functions) {
563
564 auto ptr = get_normals_ptr();
566 &ptr[0], &ptr[1], &ptr[2]);
567
568 auto t_base = data.getFTensor1N<3>(base);
569 for (int gg = 0; gg != nb_gauss_pts; ++gg) {
570 const auto l2 = t_normal(i) * t_normal(i);
571 for (int bb = 0; bb != nb_base_functions; ++bb) {
572 const auto v = t_base(0);
573 t_base(i) = (v / l2) * t_normal(i);
574 ++t_base;
575 }
576 ++t_normal;
577 }
579 };
580
581 for (int b = AINSWORTH_LEGENDRE_BASE; b != LASTBASE; b++) {
582
583 FieldApproximationBase base = static_cast<FieldApproximationBase>(b);
584 const auto &base_functions = data.getN(base);
585 const auto nb_gauss_pts = base_functions.size1();
586
587 if (nb_gauss_pts) {
588
589 const auto nb_base_functions = base_functions.size2() / 3;
590 CHKERR apply_transform_nonlinear_geometry(base, nb_gauss_pts,
591 nb_base_functions);
592 }
593 }
594
596}
597
599 int side, EntityType type, EntitiesFieldData::EntData &data) {
601
602 const auto type_dim = moab::CN::Dimension(type);
603 if (type_dim != 1 && type_dim != 2)
605
606 FTensor::Index<'i', 3> i;
607 FTensor::Index<'j', 3> j;
608 FTensor::Index<'k', 2> k;
609
610 auto get_jac = [&]() {
612 double *ptr_n = &*normalsAtPts->data().begin();
613 double *ptr_t1 = &*tangent1AtPts->data().begin();
614 double *ptr_t2 = &*tangent2AtPts->data().begin();
616 &ptr_t1[0], &ptr_t2[0], &ptr_n[0],
617
618 &ptr_t1[1], &ptr_t2[1], &ptr_n[1],
619
620 &ptr_t1[2], &ptr_t2[2], &ptr_n[2]);
621 } else {
622 double *ptr_n = &*getNormalsAtGaussPts().data().begin();
623 double *ptr_t1 = &*getTangent1AtGaussPts().data().begin();
624 double *ptr_t2 = &*getTangent2AtGaussPts().data().begin();
626 &ptr_t1[0], &ptr_t2[0], &ptr_n[0],
627
628 &ptr_t1[1], &ptr_t2[1], &ptr_n[1],
629
630 &ptr_t1[2], &ptr_t2[2], &ptr_n[2]);
631 }
632 };
633
634 for (int b = AINSWORTH_LEGENDRE_BASE; b != LASTBASE; ++b) {
635
636 FieldApproximationBase base = static_cast<FieldApproximationBase>(b);
637
638 auto &baseN = data.getN(base);
639 auto &diffBaseN = data.getDiffN(base);
640
641 int nb_dofs = baseN.size2() / 3;
642 int nb_gauss_pts = baseN.size1();
643
644 piolaN.resize(baseN.size1(), baseN.size2());
645 diffPiolaN.resize(diffBaseN.size1(), diffBaseN.size2());
646
647 if (nb_dofs > 0 && nb_gauss_pts > 0) {
648
649 auto t_h_curl = data.getFTensor1N<3>(base);
650 auto t_diff_h_curl = data.getFTensor2DiffN<3, 2>(base);
651
652 FTensor::Tensor1<FTensor::PackPtr<double *, 3>, 3> t_transformed_h_curl(
653 &piolaN(0, HVEC0), &piolaN(0, HVEC1), &piolaN(0, HVEC2));
654
656 t_transformed_diff_h_curl(
660
661 int cc = 0;
662 double det;
664
665 // HO geometry is set, so jacobian is different at each gauss point
666 auto t_m_at_pts = get_jac();
667 for (int gg = 0; gg != nb_gauss_pts; ++gg) {
668 CHKERR determinantTensor3by3(t_m_at_pts, det);
669 CHKERR invertTensor3by3(t_m_at_pts, det, t_inv_m);
670 for (int ll = 0; ll != nb_dofs; ll++) {
671 t_transformed_h_curl(i) = t_inv_m(j, i) * t_h_curl(j);
672 t_transformed_diff_h_curl(i, k) = t_inv_m(j, i) * t_diff_h_curl(j, k);
673 ++t_h_curl;
674 ++t_transformed_h_curl;
675 ++t_diff_h_curl;
676 ++t_transformed_diff_h_curl;
677 ++cc;
678 }
679 ++t_m_at_pts;
680 }
681
682#ifndef NDEBUG
683 if (cc != nb_gauss_pts * nb_dofs)
684 SETERRQ(PETSC_COMM_SELF, MOFEM_IMPOSSIBLE_CASE, "Data inconsistency");
685#endif
686
687 baseN.swap(piolaN);
688 diffBaseN.swap(diffPiolaN);
689 }
690 }
691
693}
694
696 int side, EntityType type, EntitiesFieldData::EntData &data) {
697 FTensor::Index<'i', 3> i;
699
700 if (type != MBEDGE)
702
703 const auto nb_gauss_pts = getGaussPts().size2();
704
706 if (tangentAtGaussPts->size1() != nb_gauss_pts)
707 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
708 "Wrong number of integration points %ld!=%ld",
709 tangentAtGaussPts->size1(), nb_gauss_pts);
710
711 auto get_base_at_pts = [&](auto base) {
713 &data.getN(base)(0, HVEC0), &data.getN(base)(0, HVEC1),
714 &data.getN(base)(0, HVEC2));
715 return t_h_curl;
716 };
717
718 auto get_tangent_ptr = [&]() {
719 if (tangentAtGaussPts) {
720 return &*tangentAtGaussPts->data().begin();
721 } else {
722 return &*getTangentAtGaussPts().data().begin();
723 }
724 };
725
726 auto get_tangent_at_pts = [&]() {
727 auto ptr = get_tangent_ptr();
729 &ptr[0], &ptr[1], &ptr[2]);
730 return t_m_at_pts;
731 };
732
733 auto calculate_squared_edge_length = [&]() {
734 std::vector<double> l1;
735 if (nb_gauss_pts) {
736 l1.resize(nb_gauss_pts);
737 auto t_m_at_pts = get_tangent_at_pts();
738 for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
739 l1[gg] = t_m_at_pts(i) * t_m_at_pts(i);
740 ++t_m_at_pts;
741 }
742 }
743 return l1;
744 };
745
746 auto l1 = calculate_squared_edge_length();
747
748 for (int b = AINSWORTH_LEGENDRE_BASE; b != LASTBASE; b++) {
749
750 FieldApproximationBase base = static_cast<FieldApproximationBase>(b);
751 const auto nb_dofs = data.getN(base).size2() / 3;
752
753 if (nb_gauss_pts && nb_dofs) {
754 auto t_h_curl = get_base_at_pts(base);
755 int cc = 0;
756 auto t_m_at_pts = get_tangent_at_pts();
757 for (int gg = 0; gg != nb_gauss_pts; ++gg) {
758 const double l0 = l1[gg];
759 for (int ll = 0; ll != nb_dofs; ++ll) {
760 const double val = t_h_curl(0);
761 const double a = val / l0;
762 t_h_curl(i) = t_m_at_pts(i) * a;
763 ++t_h_curl;
764 ++cc;
765 }
766 ++t_m_at_pts;
767 }
768
769#ifndef NDEBUG
770 if (cc != nb_gauss_pts * nb_dofs)
771 SETERRQ(PETSC_COMM_SELF, MOFEM_IMPOSSIBLE_CASE, "Data inconsistency");
772#endif
773 }
774 }
775
777}
778
780 const FieldSpace space, const FieldApproximationBase base,
781 boost::shared_ptr<VectorDouble> det_jac_ptr,
782 boost::shared_ptr<double> scale_det_ptr)
783 : OP(space), fieldSpace(space), fieldBase(base), detJacPtr(det_jac_ptr),
784 scaleDetPtr(scale_det_ptr) {
785 if (!detJacPtr) {
786 CHK_THROW_MESSAGE(MOFEM_DATA_INCONSISTENCY, "detJacPtr not set");
787 }
788}
789
794
795 double scale_det = 1.0;
796 if (scaleDetPtr)
797 scale_det = *scaleDetPtr;
798
799 auto scale = [&](auto base) {
801 auto &base_fun = data.getN(base);
802 if (base_fun.size2()) {
803
804 auto &det_vec = *detJacPtr;
805 const auto nb_int_pts = base_fun.size1();
806
807 if (nb_int_pts != det_vec.size())
808 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
809 "Number of integration pts in detJacPtr does not mush "
810 "number of integration pts in base function");
811
812 for (auto gg = 0; gg != nb_int_pts; ++gg) {
813 auto row = ublas::matrix_row<MatrixDouble>(base_fun, gg);
814 row *= scale_det / det_vec[gg];
815 }
816
817 auto &diff_base_fun = data.getDiffN(base);
818 if (diff_base_fun.size2()) {
819 for (auto gg = 0; gg != nb_int_pts; ++gg) {
820 auto row = ublas::matrix_row<MatrixDouble>(diff_base_fun, gg);
821 row *= scale_det / det_vec[gg];
822 }
823 }
824 }
826 };
827
828 if (this->getFEDim() == 3) {
829 switch (fieldSpace) {
830 case L2:
831 if (type >= MBTET)
833 break;
834 default:
835 SETERRQ(PETSC_COMM_SELF, MOFEM_IMPOSSIBLE_CASE, "space not handled");
836 }
837 } else if (this->getFEDim() == 2) {
838 switch (fieldSpace) {
839 case L2:
840 if (type >= MBTRI)
842 break;
843 default:
844 SETERRQ(PETSC_COMM_SELF, MOFEM_IMPOSSIBLE_CASE, "space not handled");
845 }
846 } else if (this->getFEDim() == 1) {
847 switch (fieldSpace) {
848 case L2:
849 if (type >= MBEDGE)
851 break;
852 default:
853 SETERRQ(PETSC_COMM_SELF, MOFEM_IMPOSSIBLE_CASE, "space not handled");
854 }
855 } else {
856 SETERRQ(PETSC_COMM_SELF, MOFEM_IMPOSSIBLE_CASE, "dimension not handled");
857 }
858
860}
861
863 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pipeline,
864 std::vector<FieldSpace> spaces, std::string geom_field_name,
865 boost::shared_ptr<MatrixDouble> jac_ptr,
866 boost::shared_ptr<VectorDouble> det_ptr,
867 boost::shared_ptr<MatrixDouble> inv_jac_ptr) {
869
870 if (!jac_ptr)
871 jac_ptr = boost::make_shared<MatrixDouble>();
872 if (!det_ptr)
873 det_ptr = boost::make_shared<VectorDouble>();
874 if (!inv_jac_ptr)
875 inv_jac_ptr = boost::make_shared<MatrixDouble>();
876
877 if (geom_field_name.empty()) {
878
879 pipeline.push_back(new OpCalculateHOJac<2>(jac_ptr));
880
881 } else {
882
883 pipeline.push_back(new OpCalculateHOCoords<2>(geom_field_name));
884 pipeline.push_back(
885 new OpCalculateVectorFieldGradient<2, 2>(geom_field_name, jac_ptr));
886 pipeline.push_back(new OpGetHONormalsOnFace<2>(geom_field_name));
887 }
888
889 pipeline.push_back(new OpInvertMatrix<2>(jac_ptr, det_ptr, inv_jac_ptr));
890 pipeline.push_back(new OpSetHOWeightsOnFace());
891
892 for (auto s : spaces) {
893 switch (s) {
894 case NOSPACE:
895 break;
896 case H1:
897 case L2:
898 pipeline.push_back(new OpSetHOInvJacToScalarBases<2>(s, inv_jac_ptr));
899 break;
900 case HCURL:
901 pipeline.push_back(new OpSetCovariantPiolaTransformOnFace2D(inv_jac_ptr));
902 pipeline.push_back(new OpSetInvJacHcurlFace(inv_jac_ptr));
903 break;
904 case HDIV:
905 pipeline.push_back(new OpMakeHdivFromHcurl());
906 pipeline.push_back(new OpSetContravariantPiolaTransformOnFace2D(jac_ptr));
907 pipeline.push_back(new OpSetInvJacHcurlFace(inv_jac_ptr));
908 break;
909 default:
910 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
911 "Space %s not yet implemented", FieldSpaceNames[s]);
912 }
913 }
914
916}
917
919 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pipeline,
920 std::vector<FieldSpace> spaces, std::string geom_field_name,
921 boost::shared_ptr<MatrixDouble> jac_ptr,
922 boost::shared_ptr<VectorDouble> det_ptr,
923 boost::shared_ptr<MatrixDouble> inv_jac_ptr) {
925
926 if (!jac_ptr)
927 jac_ptr = boost::make_shared<MatrixDouble>();
928 if (!det_ptr)
929 det_ptr = boost::make_shared<VectorDouble>();
930 if (!inv_jac_ptr)
931 inv_jac_ptr = boost::make_shared<MatrixDouble>();
932
933 if (geom_field_name.empty()) {
934
935 } else {
936
937 pipeline.push_back(new OpCalculateHOCoords<3>(geom_field_name));
938 pipeline.push_back(new OpGetHONormalsOnFace<3>(geom_field_name));
939 }
940
941 pipeline.push_back(new OpCalculateHOJacForFaceEmbeddedIn3DSpace(jac_ptr));
942 pipeline.push_back(new OpInvertMatrix<3>(jac_ptr, det_ptr, inv_jac_ptr));
943 pipeline.push_back(new OpSetHOWeightsOnFace());
944
945 for (auto s : spaces) {
946 switch (s) {
947 case NOSPACE:
948 break;
949 case H1:
950 pipeline.push_back(
952 break;
953 case HDIV:
954 pipeline.push_back(new OpMakeHdivFromHcurl());
955 pipeline.push_back(
957 jac_ptr));
958 pipeline.push_back(
960 break;
961 case L2:
962 pipeline.push_back(
964 break;
965 default:
966 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
967 "Space %s not yet implemented", FieldSpaceNames[s]);
968 }
969 }
970
972}
973
975 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pipeline,
976 std::vector<FieldSpace> spaces, std::string geom_field_name) {
978
979 if (geom_field_name.empty()) {
980
981 } else {
982
983 pipeline.push_back(new OpCalculateHOCoords<2>(geom_field_name));
984 pipeline.push_back(new OpGetHOTangentsOnEdge<2>(geom_field_name));
985 }
986
987 for (auto s : spaces) {
988 switch (s) {
989 case NOSPACE:
990 break;
991 case HCURL:
992 pipeline.push_back(new OpHOSetContravariantPiolaTransformOnEdge3D(HCURL));
993 break;
994 case HDIV:
995 pipeline.push_back(new OpSetContravariantPiolaTransformOnEdge2D());
996 break;
997 default:
998 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
999 "Space %s not yet implemented", FieldSpaceNames[s]);
1000 }
1001 }
1002
1004}
1005
1007 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pipeline,
1008 std::vector<FieldSpace> spaces, std::string geom_field_name) {
1010
1011 if (geom_field_name.empty()) {
1012
1013 } else {
1014
1015 pipeline.push_back(new OpCalculateHOCoords<3>(geom_field_name));
1016 pipeline.push_back(new OpGetHOTangentsOnEdge<3>(geom_field_name));
1017 }
1018
1019 for (auto s : spaces) {
1020 switch (s) {
1021 case HCURL:
1022 pipeline.push_back(new OpHOSetContravariantPiolaTransformOnEdge3D(HCURL));
1023 break;
1024 default:
1025 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
1026 "Space %s not yet implemented", FieldSpaceNames[s]);
1027 }
1028 }
1029
1031}
1032
1034 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pipeline,
1035 std::vector<FieldSpace> spaces, std::string geom_field_name) {
1037
1038 if (geom_field_name.empty()) {
1039 } else {
1040
1041 pipeline.push_back(new OpCalculateHOCoords<3>(geom_field_name));
1042 pipeline.push_back(new OpGetHONormalsOnFace<3>(geom_field_name));
1043 }
1044
1045 for (auto s : spaces) {
1046 switch (s) {
1047 case NOSPACE:
1048 break;
1049 case HCURL:
1050 pipeline.push_back(new OpHOSetCovariantPiolaTransformOnFace3D(HCURL));
1051 break;
1052 case HDIV:
1053 pipeline.push_back(new OpHOSetContravariantPiolaTransformOnFace3D(HDIV));
1054 break;
1055 case L2:
1056 break;
1057 default:
1058 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
1059 "Space %s not yet implemented", FieldSpaceNames[s]);
1060 break;
1061 }
1062 }
1063
1065}
1066
1068 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pipeline,
1069 std::vector<FieldSpace> spaces, std::string geom_field_name,
1070 boost::shared_ptr<MatrixDouble> jac, boost::shared_ptr<VectorDouble> det,
1071 boost::shared_ptr<MatrixDouble> inv_jac) {
1073
1074 if (!jac)
1075 jac = boost::make_shared<MatrixDouble>();
1076 if (!det)
1077 det = boost::make_shared<VectorDouble>();
1078 if (!inv_jac)
1079 inv_jac = boost::make_shared<MatrixDouble>();
1080
1081 if (geom_field_name.empty()) {
1082
1083 pipeline.push_back(new OpCalculateHOJac<3>(jac));
1084
1085 } else {
1086
1087 pipeline.push_back(new OpCalculateHOCoords<3>(geom_field_name));
1088 pipeline.push_back(
1089 new OpCalculateVectorFieldGradient<3, 3>(geom_field_name, jac));
1090 }
1091
1092 pipeline.push_back(new OpInvertMatrix<3>(jac, det, inv_jac));
1093 pipeline.push_back(new OpSetHOWeights(det));
1094
1095 for (auto s : spaces) {
1096 switch (s) {
1097 case NOSPACE:
1098 break;
1099 case H1:
1100 pipeline.push_back(new OpSetHOInvJacToScalarBases<3>(H1, inv_jac));
1101 break;
1102 case HCURL:
1103 pipeline.push_back(new OpSetHOCovariantPiolaTransform(HCURL, inv_jac));
1104 pipeline.push_back(new OpSetHOInvJacVectorBase(HCURL, inv_jac));
1105 break;
1106 case HDIV:
1107 pipeline.push_back(
1109 break;
1110 case L2:
1111 pipeline.push_back(new OpSetHOInvJacToScalarBases<3>(L2, inv_jac));
1112 break;
1113 default:
1114 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
1115 "Space %s not yet implemented", FieldSpaceNames[s]);
1116 break;
1117 }
1118 }
1119
1121}
1122
1123} // namespace MoFEM
static MoFEMErrorCode get_jac(EntitiesFieldData::EntData &col_data, int gg, MatrixDouble &jac_stress, MatrixDouble &jac)
constexpr double a
T data[Tensor_Dim]
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
#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
@ 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
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 ...
@ HVEC0
@ HVEC1
@ HVEC2
@ MOFEM_IMPOSSIBLE_CASE
Definition definitions.h:35
@ 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()
@ HVEC1_1
@ HVEC0_1
@ HVEC1_0
@ HVEC2_1
@ HVEC1_2
@ HVEC2_2
@ HVEC2_0
@ HVEC0_2
@ HVEC0_0
#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.
FTensor::Index< 'i', SPACE_DIM > i
const double v
phase velocity of light in medium (cm/ns)
FTensor::Index< 'l', 3 > l
FTensor::Index< 'j', 3 > j
FTensor::Index< 'k', 3 > k
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
implementation of Data Operators for Forces and Sources
Definition Common.hpp:10
OpSetInvJacHcurlFaceImpl< 2 > OpSetInvJacHcurlFace
OpSetContravariantPiolaTransformOnFace2DImpl< 3 > OpSetContravariantPiolaTransformOnFace2DEmbeddedIn3DSpace
OpSetInvJacHcurlFaceImpl< 3 > OpSetInvJacHcurlFaceEmbeddedIn3DSpace
MoFEMErrorCode invertTensor3by3(ublas::matrix< T, L, A > &jac_data, ublas::vector< T, A > &det_data, ublas::matrix< T, L, A > &inv_jac_data)
Calculate inverse of tensor rank 2 at integration points.
OpSetCovariantPiolaTransformOnFace2DImpl< 2 > OpSetCovariantPiolaTransformOnFace2D
OpSetContravariantPiolaTransformOnFace2DImpl< 2 > OpSetContravariantPiolaTransformOnFace2D
static auto getFTensor0FromVec(ublas::vector< T, A > &data)
Get tensor rank 0 (scalar) form data vector.
OpCalculateHOJacForFaceImpl< 3 > OpCalculateHOJacForFaceEmbeddedIn3DSpace
static auto determinantTensor3by3(T &t)
Calculate the determinant of a 3x3 matrix or a tensor of rank 2.
constexpr double t
plate stiffness
Definition plate.cpp:58
FTensor::Index< 'm', 3 > m
Add operators pushing bases from local to physical configuration.
std::array< bool, MBMAXTYPE > doEntities
If true operator is executed for entity.
virtual MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
Operator for bi-linear form, usually to calculate values on left hand side.
MatrixDouble & getTangentAtGaussPts()
get tangent vector to edge curve at integration points
FTensor::Tensor1< FTensor::PackPtr< double *, 3 >, DIM > getFTensor1TangentAtGaussPts()
Get tangent vector at integration points.
Data on single entity (This is passed as argument to DataOperator::doWork)
FTensor::Tensor2< FTensor::PackPtr< double *, Tensor_Dim0 *Tensor_Dim1 >, Tensor_Dim0, Tensor_Dim1 > getFTensor2DiffN(FieldApproximationBase base)
Get derivatives of base functions for Hdiv space.
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1DiffN(const FieldApproximationBase base)
Get derivatives of base functions.
virtual std::array< boost::shared_ptr< MatrixDouble >, MaxBernsteinBezierOrder > & getBBDiffNByOrderArray()
MatrixDouble & getDiffN(const FieldApproximationBase base)
get derivatives of base functions
virtual boost::shared_ptr< MatrixDouble > & getNSharedPtr(const FieldApproximationBase base, const BaseDerivatives derivative)
Get shared pointer to base functions with derivatives.
MatrixDouble & getN(const FieldApproximationBase base)
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....
virtual std::map< std::string, boost::shared_ptr< MatrixDouble > > & getBBDiffNMap()
get hash map of derivatives base function for BB base, key is a field name
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1N(FieldApproximationBase base)
Get base functions for Hdiv/Hcurl spaces.
virtual boost::shared_ptr< MatrixDouble > & getDiffNSharedPtr(const FieldApproximationBase base)
Get shared pointer to derivatives of base functions.
MatrixDouble & getNormalsAtGaussPts()
if higher order geometry return normals at Gauss pts.
MatrixDouble & getTangent2AtGaussPts()
if higher order geometry return tangent vector to triangle at Gauss pts.
MatrixDouble & getTangent1AtGaussPts()
if higher order geometry return tangent vector to triangle at Gauss pts.
auto getFTensor0IntegrationWeight()
Get integration weights.
double getMeasure() const
get measure of element
boost::shared_ptr< const NumeredEntFiniteElement > getNumeredEntFiniteElementPtr() const
Return raw pointer to NumeredEntFiniteElement.
int getFEDim() const
Get dimension of finite element.
MatrixDouble & getGaussPts()
matrix of integration (Gauss) points for Volume Element
Calculate high-order coordinates at Gauss points.
Calculate jacobian for face element.
boost::shared_ptr< MatrixDouble > jacPtr
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
OpCalculateHOJacForVolume(boost::shared_ptr< MatrixDouble > jac_ptr)
Constructor for HO Jacobian calculation on volumes.
Get field gradients at integration pts for scalar field rank 0, i.e. vector field.
Calculate normal vectors at Gauss points of face elements.
Calculate tangent vector on edge from HO geometry approximation.
Transform Hcurl base functions from reference element to physical edge in 3D.
boost::shared_ptr< MatrixDouble > tangentAtGaussPts
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
Transform Hdiv base fluxes from reference element to physical face in 3D.
boost::shared_ptr< MatrixDouble > normalsAtGaussPts
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
Transform Hcurl base functions from reference element to physical face in 3D.
boost::shared_ptr< MatrixDouble > normalsAtPts
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
boost::shared_ptr< MatrixDouble > tangent2AtPts
boost::shared_ptr< MatrixDouble > tangent1AtPts
Operator for inverting matrices at integration points.
Make Hdiv space from Hcurl space in 2d.
boost::shared_ptr< VectorDouble > detJacPtr
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
OpScaleBaseBySpaceInverseOfMeasure(const FieldSpace space, const FieldApproximationBase base, boost::shared_ptr< VectorDouble > det_jac_ptr, boost::shared_ptr< double > scale_det_ptr=nullptr)
Apply contravariant (Piola) transformation to Hdiv space for HO geometry.
boost::shared_ptr< VectorDouble > detPtr
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
boost::shared_ptr< MatrixDouble > jacPtr
Apply covariant (Piola) transformation to Hcurl space for HO geometry.
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
boost::shared_ptr< MatrixDouble > jacInvPtr
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
Set inverse jacobian to base functions.
Transform local reference derivatives of shape functions to global derivatives.
boost::shared_ptr< MatrixDouble > invJacPtr
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
Modify integration weights on face to take into account higher-order geometry.
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
Set integration weights for higher-order elements.
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
boost::shared_ptr< VectorDouble > detPtr
FTensor::Tensor2< double *, 3, 3 > & getInvJac()
get element inverse Jacobian
double scale
Definition plastic.cpp:123