v0.14.0
HODataOperators.hpp
Go to the documentation of this file.
1 /** \file HODataOperators.hpp
2  * \brief Operators managing HO geometry
3 
4 */
5 
6 #ifndef __HO_DATA_OPERATORS_HPP__
7 #define __HO_DATA_OPERATORS_HPP__
8 
9 namespace MoFEM {
10 
11 /**
12  * @brief Calculate jacobian on Hex or other volume which is not simplex.
13  *
14  * Using not a field data, but geometrical coordinates of element.
15  *
16  * \note Use OpCalculateVectorFieldGradient to calculate Jacobian from field
17  * data.
18  *
19  */
22 
23  OpCalculateHOJacForVolume(boost::shared_ptr<MatrixDouble> jac_ptr);
24 
25  MoFEMErrorCode doWork(int side, EntityType type,
27 
28 private:
29  boost::shared_ptr<MatrixDouble> jacPtr;
30 };
31 
32 /** \deprecated use OpCalculateHOJacForVolume
33  */
35 
36 /**
37  * @brief Calculate HO coordinates at gauss points
38  *
39  */
40 template <int FIELD_DIM = 3>
42 
43  OpCalculateHOCoords(const std::string field_name)
45 
46  MoFEMErrorCode doWork(int side, EntityType type,
48 };
49 
50 /**
51  * @deprecated This class should be DIM = 3 specialization when default
52  * parameter is removed
53  *
54  */
56 
58 
59  MoFEMErrorCode doWork(int side, EntityType type,
61 };
62 
63 /**
64  * \brief Set inverse jacobian to base functions
65  *
66  * \deprecated Version with default variant DIM = 3 is deprecated, keep for
67  * back compatibility. Should be removed in future. Use of default variant make
68  * code implicit, what can be source of some hidden error. Explict interface is
69  * better, when user see and control outcome, and is aware of existing variants.
70  *
71  */
72 template <int DIM = 3>
74  using OpSetHOInvJacToScalarBasesImpl::OpSetHOInvJacToScalarBasesImpl;
75 };
76 
77 template <>
79  : public OpSetInvJacSpaceForFaceImpl<2, 1> {
80 
82 };
83 
84 /**
85  * \brief transform local reference derivatives of shape function to global
86  derivatives if higher order geometry is given
87  *
88 
89  * \ingroup mofem_forces_and_sources
90 */
92 
94  boost::shared_ptr<MatrixDouble> inv_jac_ptr)
95  : ForcesAndSourcesCore::UserDataOperator(space), invJacPtr(inv_jac_ptr) {
96 
97  if (!invJacPtr)
99  "Pointer for invJacPtr not allocated");
100 
101  doVertices = false;
102  if (space == HDIV)
103  doEdges = false;
104  }
105 
106 
107  MoFEMErrorCode doWork(int side, EntityType type,
109 
110 private:
111  boost::shared_ptr<MatrixDouble> invJacPtr;
114 };
115 
116 /**
117  * @brief Modify integration weights on face to take in account higher-order
118  * geometry
119  * @ingroup mofem_forces_and_sources_tri_element
120  *
121  */
126  MoFEMErrorCode doWork(int side, EntityType type,
128 };
129 
130 /**
131  * @brief Modify integration weights on face to take in account higher-order
132  * geometry
133  * @ingroup mofem_forces_and_sources_tri_element
134  *
135  */
140  MoFEMErrorCode doWork(int side, EntityType type,
142 };
143 
144 template <int SPACE_DIM>
146 
147 template <> struct OpSetHOWeightsOnSubDim<2> : public OpSetHOWeightsOnEdge {
149 };
150 
151 template <> struct OpSetHOWeightsOnSubDim<3> : public OpSetHOWeightsOnFace {
153 };
154 
155 /**
156  * \brief Set inverse jacobian to base functions
157  *
158  */
160 
161  OpSetHOWeights(boost::shared_ptr<VectorDouble> det_ptr)
163  detPtr(det_ptr) {
164  if (!detPtr)
166  "Pointer for detPtr not allocated");
167  }
168 
169  MoFEMErrorCode doWork(int side, EntityType type,
171 
172 private:
173  boost::shared_ptr<VectorDouble> detPtr;
174 };
175 
176 /** \brief Apply contravariant (Piola) transfer to Hdiv space for HO geometry
177 
178 * \ingroup mofem_forces_and_sources
179 */
182 
184  boost::shared_ptr<VectorDouble> det_ptr,
185  boost::shared_ptr<MatrixDouble> jac_ptr)
187  jacPtr(jac_ptr) {
188  doVertices = false;
189  if (space == HDIV)
190  doEdges = false;
191  }
192 
193  MoFEMErrorCode doWork(int side, EntityType type,
195 
196 private:
197  boost::shared_ptr<VectorDouble> detPtr;
198  boost::shared_ptr<MatrixDouble> jacPtr;
199 
202 };
203 
204 /** \brief Apply covariant (Piola) transfer to Hcurl space for HO geometry
205  */
208 
210  boost::shared_ptr<MatrixDouble> jac_inv_ptr)
212  jacInvPtr(jac_inv_ptr) {
213  doVertices = false;
214  if (space == HDIV)
215  doEdges = false;
216  if (!jacInvPtr)
218  "Pointer for jacPtr not allocated");
219  }
220 
221  MoFEMErrorCode doWork(int side, EntityType type,
223 
224 private:
225  boost::shared_ptr<MatrixDouble> jacInvPtr;
226 
229 };
230 
231 /** \brief Calculate jacobian for face element
232 
233  It is assumed that face element is XY plane. Applied
234  only for 2d problems and 2d problems embedded in 3d space
235 
236  \note If you push operators for HO normal befor this operator, HO geometry is
237  taken into account when you calculate jacobian.
238 
239 */
240 template <int DIM> struct OpCalculateHOJacForFaceImpl;
241 
242 template <>
245 
246  OpCalculateHOJacForFaceImpl(boost::shared_ptr<MatrixDouble> jac_ptr);
247 
248  MoFEMErrorCode doWork(int side, EntityType type,
250 
251 protected:
252  boost::shared_ptr<MatrixDouble> jacPtr;
253 };
254 
255 template <>
257 
259 
260  MoFEMErrorCode doWork(int side, EntityType type,
262 };
263 
266 
267 template <int DIM> struct OpCalculateHOJac;
268 
269 template <> struct OpCalculateHOJac<3> : public OpCalculateHOJacForVolume {
271 };
272 
273 template <> struct OpCalculateHOJac<2> : public OpCalculateHOJacForFaceImpl<2> {
275 };
276 
277 /** \brief Calculate normals at Gauss points of triangle element
278  * \ingroup mofem_forces_and_source
279  */
280 template <int FIELD_DIM = 3>
283 
286 
287  MoFEMErrorCode doWork(int side, EntityType type,
289 };
290 
291 /** \brief transform Hdiv base fluxes from reference element to physical
292  * triangle \ingroup mofem_forces_and_sources
293  *
294  * \note Matrix which keeps normal is assumed to have three columns, and number
295  * of rows should be equal to number of integration points.
296  *
297  */
300 
302  const FieldSpace space,
303  boost::shared_ptr<MatrixDouble> normals_at_gauss_pts = nullptr)
305  normalsAtGaussPts(normals_at_gauss_pts) {}
306 
307  MoFEMErrorCode doWork(int side, EntityType type,
309 
310 private:
311  boost::shared_ptr<MatrixDouble> normalsAtGaussPts;
312 };
313 
314 /** \brief transform Hcurl base fluxes from reference element to physical edge
315  * \ingroup mofem_forces_and_sources
316  */
319 
321  const FieldSpace space = HCURL,
322  boost::shared_ptr<MatrixDouble> tangent_at_pts = nullptr)
324  tangentAtGaussPts(tangent_at_pts) {}
325 
326  MoFEMErrorCode doWork(int side, EntityType type,
328 
329 private:
330  boost::shared_ptr<MatrixDouble> tangentAtGaussPts;
331 };
332 
333 /** \brief transform Hcurl base fluxes from reference element to physical
334  * triangle \ingroup mofem_forces_and_sources
335  */
338 
340  const FieldSpace space,
341  boost::shared_ptr<MatrixDouble> normals_at_pts = nullptr,
342  boost::shared_ptr<MatrixDouble> tangent1_at_pts = nullptr,
343  boost::shared_ptr<MatrixDouble> tangent2_at_pts = nullptr)
345  normalsAtPts(normals_at_pts), tangent1AtPts(tangent1_at_pts),
346  tangent2AtPts(tangent2_at_pts) {
347  if (normals_at_pts || tangent1_at_pts || tangent2_at_pts)
348  if (normals_at_pts && tangent1_at_pts && tangent2_at_pts)
350  "All elements in constructor have to set pointer");
351  }
352 
353  MoFEMErrorCode doWork(int side, EntityType type,
355 
356 private:
357  boost::shared_ptr<MatrixDouble> normalsAtPts;
358  boost::shared_ptr<MatrixDouble> tangent1AtPts;
359  boost::shared_ptr<MatrixDouble> tangent2AtPts;
360 
363 };
364 
365 /** \brief Calculate tangent vector on edge form HO geometry approximation
366  * \ingroup mofem_forces_and_sources
367  */
368 template <int FIELD_DIM = 3>
371 
373  std::string field_name,
374  boost::shared_ptr<MatrixDouble> tangents_at_pts = nullptr)
376  tangentsAtPts(tangents_at_pts) {}
377 
378  MoFEMErrorCode doWork(int side, EntityType type,
380 
381 private:
382  boost::shared_ptr<MatrixDouble> tangentsAtPts;
383 };
384 
385 /**
386  * @brief Scale base functions by inverses of measure of element
387  *
388  * @tparam OP
389  */
392 
394 
396  const FieldSpace space, boost::shared_ptr<VectorDouble> det_jac_ptr,
397  boost::shared_ptr<double> scale_det_ptr = nullptr);
398 
399  MoFEMErrorCode doWork(int side, EntityType type,
401 
402 private:
404  boost::shared_ptr<VectorDouble> detJacPtr;
405  boost::shared_ptr<double> scaleDetPtr = nullptr;
406 };
407 
408 /**
409  * @brief Add operators pushing bases from local to physical configuration
410  *
411  * @tparam FE_DIM dimension of element
412  * @tparam PROBLEM_DIM problem dimension
413  * @tparam SPACE_DIM space dimension
414  */
415 template <int FE_DIM, int PROBLEM_DIM, int SPACE_DIM> struct AddHOOps;
416 
417 template <> struct AddHOOps<2, 2, 2> {
418  AddHOOps() = delete;
419  static MoFEMErrorCode
420  add(boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pipeline,
421  std::vector<FieldSpace> spaces, std::string geom_field_name = "",
422  boost::shared_ptr<MatrixDouble> jac = nullptr,
423  boost::shared_ptr<VectorDouble> det = nullptr,
424  boost::shared_ptr<MatrixDouble> inv_jac = nullptr);
425 };
426 
427 template <> struct AddHOOps<2, 2, 3> {
428  AddHOOps() = delete;
429  static MoFEMErrorCode
430  add(boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pipeline,
431  std::vector<FieldSpace> spaces, std::string geom_field_name = "",
432  boost::shared_ptr<MatrixDouble> jac = nullptr,
433  boost::shared_ptr<VectorDouble> det = nullptr,
434  boost::shared_ptr<MatrixDouble> inv_jac = nullptr);
435 };
436 
437 template <> struct AddHOOps<1, 2, 2> {
438  AddHOOps() = delete;
439  static MoFEMErrorCode
440  add(boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pipeline,
441  std::vector<FieldSpace> spaces, std::string geom_field_name = "");
442 };
443 template <> struct AddHOOps<1, 3, 3> {
444  AddHOOps() = delete;
445  static MoFEMErrorCode
446  add(boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pipeline,
447  std::vector<FieldSpace> space, std::string geom_field_name = "");
448 };
449 
450 template <> struct AddHOOps<2, 3, 3> {
451  AddHOOps() = delete;
452  static MoFEMErrorCode
453  add(boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pipeline,
454  std::vector<FieldSpace> space, std::string geom_field_name = "");
455 };
456 
457 template <> struct AddHOOps<3, 3, 3> {
458  AddHOOps() = delete;
459  static MoFEMErrorCode
460  add(boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pipeline,
461  std::vector<FieldSpace> space, std::string geom_field_name = "",
462  boost::shared_ptr<MatrixDouble> jac = nullptr,
463  boost::shared_ptr<VectorDouble> det = nullptr,
464  boost::shared_ptr<MatrixDouble> inv_jac = nullptr);
465 };
466 
467 template <int FIELD_DIM>
473  const auto nb_dofs = data.getFieldData().size() / FIELD_DIM;
474  if (nb_dofs) {
475  if (type == MBVERTEX)
476  getCoordsAtGaussPts().clear();
477  auto t_base = data.getFTensor0N();
478  auto t_coords = getFTensor1CoordsAtGaussPts();
479  const auto nb_integration_pts = data.getN().size1();
480  const auto nb_base_functions = data.getN().size2();
481  for (auto gg = 0; gg != nb_integration_pts; ++gg) {
482  auto t_dof = data.getFTensor1FieldData<FIELD_DIM>();
483  size_t bb = 0;
484  for (; bb != nb_dofs; ++bb) {
485  t_coords(i) += t_base * t_dof(i);
486  ++t_dof;
487  ++t_base;
488  }
489  for (; bb < nb_base_functions; ++bb)
490  ++t_base;
491  ++t_coords;
492  }
493  }
495 }
496 
497 template <int FIELD_DIM>
502 
505 
506  auto get_ftensor1 = [](MatrixDouble &m) {
508  &m(0, 0), &m(0, 1), &m(0, 2));
509  };
510 
511  unsigned int nb_dofs = data.getFieldData().size();
512  if (nb_dofs != 0) {
513 
514  int nb_gauss_pts = data.getN().size1();
515  auto &tangent1_at_gauss_pts = getTangent1AtGaussPts();
516  auto &tangent2_at_gauss_pts = getTangent2AtGaussPts();
517  tangent1_at_gauss_pts.resize(nb_gauss_pts, 3, false);
518  tangent2_at_gauss_pts.resize(nb_gauss_pts, 3, false);
519 
520  switch (type) {
521  case MBVERTEX: {
522  tangent1_at_gauss_pts.clear();
523  tangent2_at_gauss_pts.clear();
524  }
525  case MBEDGE:
526  case MBTRI:
527  case MBQUAD: {
528 
529 #ifndef NDEBUG
530  if (2 * data.getN().size2() != data.getDiffN().size2()) {
531  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
532  "expected two direcatives in local element coordinates");
533  }
534  if (nb_dofs % FIELD_DIM != 0) {
535  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
536  "expected that number of dofs is multiplicative of field "
537  "dimension");
538  }
539 #endif
540 
541  if (nb_dofs > FIELD_DIM * data.getN().size2()) {
542  unsigned int nn = 0;
543  for (; nn != nb_dofs; nn++) {
544  if (!data.getFieldDofs()[nn]->getActive())
545  break;
546  }
547  if (nn > FIELD_DIM * data.getN().size2()) {
548  SETERRQ1(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
549  "Data inconsistency for base %s",
551  } else {
552  nb_dofs = nn;
553  if (!nb_dofs)
555  }
556  }
557  const int nb_base_functions = data.getN().size2();
558  auto t_base = data.getFTensor0N();
559  auto t_diff_base = data.getFTensor1DiffN<2>();
560  auto t_t1 = get_ftensor1(tangent1_at_gauss_pts);
561  auto t_t2 = get_ftensor1(tangent2_at_gauss_pts);
562  for (int gg = 0; gg != nb_gauss_pts; ++gg) {
563  auto t_data = data.getFTensor1FieldData<FIELD_DIM>();
564  int bb = 0;
565  for (; bb != nb_dofs / FIELD_DIM; ++bb) {
567  t_t1(i) += t_data(i) * t_diff_base(N0);
568  t_t2(i) += t_data(i) * t_diff_base(N1);
569  ++t_data;
570  ++t_base;
571  ++t_diff_base;
572  }
573  for (; bb != nb_base_functions; ++bb) {
574  ++t_base;
575  ++t_diff_base;
576  }
577  ++t_t1;
578  ++t_t2;
579  }
580  } break;
581  default:
582  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
583  }
584  }
585 
586  if (moab::CN::Dimension(type) == 2) {
587 
588  auto &normal_at_gauss_pts = getNormalsAtGaussPts();
589  auto &tangent1_at_gauss_pts = getTangent1AtGaussPts();
590  auto &tangent2_at_gauss_pts = getTangent2AtGaussPts();
591 
592  const auto nb_gauss_pts = tangent1_at_gauss_pts.size1();
593  normal_at_gauss_pts.resize(nb_gauss_pts, 3, false);
594 
595  auto t_normal = get_ftensor1(normal_at_gauss_pts);
596  auto t_t1 = get_ftensor1(tangent1_at_gauss_pts);
597  auto t_t2 = get_ftensor1(tangent2_at_gauss_pts);
598  for (unsigned int gg = 0; gg != nb_gauss_pts; ++gg) {
602  t_normal(j) = FTensor::levi_civita(i, j, k) * t_t1(k) * t_t2(i);
603  ++t_normal;
604  ++t_t1;
605  ++t_t2;
606  }
607  }
608 
610 }
611 
612 template <int FIELD_DIM>
617 
618  auto get_tangent = [&]() -> MatrixDouble & {
619  if (tangentsAtPts)
620  return *tangentsAtPts;
621  else
622  return getTangentAtGaussPts();
623  };
624 
625  auto &tangent = get_tangent();
626  int nb_gauss_pts = getGaussPts().size2();
627 
628  if (type == MBVERTEX) {
629  tangent.resize(nb_gauss_pts, 3, false);
630  tangent.clear();
631  }
632 
633  int nb_dofs = data.getFieldData().size();
634  if (nb_dofs != 0) {
635  const int nb_base_functions = data.getN().size2();
636  double *diff_base_function = &*data.getDiffN().data().begin();
637  auto tangent_at_gauss_pts =
638  getFTensor1FromPtr<FIELD_DIM, 3>(&*tangent.data().begin());
639 
641  int size = nb_dofs / FIELD_DIM;
642  if (nb_dofs % FIELD_DIM) {
643  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
644  "Inconsistent number of dofs and filed dimension");
645  }
646  for (int gg = 0; gg != nb_gauss_pts; ++gg) {
647  auto field_data = data.getFTensor1FieldData<FIELD_DIM>();
648  int bb = 0;
649  for (; bb < size; ++bb) {
650  tangent_at_gauss_pts(i) += field_data(i) * (*diff_base_function);
651  ++field_data;
652  ++diff_base_function;
653  }
654  for (; bb != nb_base_functions; ++bb) {
655  ++diff_base_function;
656  }
657  ++tangent_at_gauss_pts;
658  }
659  }
661 }
662 
663 /**
664  * @deprecated do not use this function, instead use AddHOOps.
665  *
666  * @tparam E
667  * @param field
668  * @param e
669  * @param h1
670  * @param hcurl
671  * @param hdiv
672  * @param l2
673  * @return MoFEMErrorCode
674  */
675 template <typename E>
676 MoFEMErrorCode addHOOpsVol(const std::string field, E &e, bool h1,
677  bool hcurl, bool hdiv, bool l2) {
678  std::vector<FieldSpace> spaces;
679  if (h1)
680  spaces.push_back(H1);
681  if (hcurl)
682  spaces.push_back(HCURL);
683  if (hdiv)
684  spaces.push_back(HDIV);
685  if (l2)
686  spaces.push_back(L2);
687  return AddHOOps<3, 3, 3>::add(e.getOpPtrVector(), spaces, field);
688 }
689 
690 /**
691  * @deprecated do not use this function, instead use AddHOOps.
692  *
693  * @tparam E
694  * @param field
695  * @param e
696  * @param hcurl
697  * @param hdiv
698  * @return MoFEMErrorCode
699  */
700 template <typename E>
701 MoFEMErrorCode addHOOpsFace3D(const std::string field, E &e,
702  bool hcurl, bool hdiv) {
703  std::vector<FieldSpace> spaces;
704  if (hcurl)
705  spaces.push_back(HCURL);
706  if (hdiv)
707  spaces.push_back(HDIV);
708  return AddHOOps<2, 3, 3>::add(e.getOpPtrVector(), spaces, field);
709 }
710 
711 }; // namespace MoFEM
712 
713 #endif
UBlasMatrix< double >
NOSPACE
@ NOSPACE
Definition: definitions.h:83
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:460
MoFEM::OpSetHOCovariantPiolaTransform::piolaN
MatrixDouble piolaN
Definition: HODataOperators.hpp:227
MoFEM::EdgeElementForcesAndSourcesCore
Edge finite element.
Definition: EdgeElementForcesAndSourcesCore.hpp:30
MoFEM::addHOOpsFace3D
MoFEMErrorCode addHOOpsFace3D(const std::string field, E &e, bool hcurl, bool hdiv)
Definition: HODataOperators.hpp:701
MoFEM::EntitiesFieldData::EntData
Data on single entity (This is passed as argument to DataOperator::doWork)
Definition: EntitiesFieldData.hpp:128
MoFEM::OpSetHOWeightsOnEdge::OpSetHOWeightsOnEdge
OpSetHOWeightsOnEdge()
Definition: HODataOperators.hpp:138
MoFEM::OpSetHOInvJacToScalarBasesImpl::doWork
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
Definition: HODataOperators.cpp:73
MoFEM::OpSetHOContravariantPiolaTransform::piolaDiffN
MatrixDouble piolaDiffN
Definition: HODataOperators.hpp:201
MoFEM::OpSetHOWeights::detPtr
boost::shared_ptr< VectorDouble > detPtr
Definition: HODataOperators.hpp:173
MoFEM::OpGetHOTangentsOnEdge
Calculate tangent vector on edge form HO geometry approximation.
Definition: HODataOperators.hpp:369
H1
@ H1
continuous field
Definition: definitions.h:85
MoFEM::OpSetHOContravariantPiolaTransform::jacPtr
boost::shared_ptr< MatrixDouble > jacPtr
Definition: HODataOperators.hpp:198
FTensor::Tensor1
Definition: Tensor1_value.hpp:8
MoFEM::ForcesAndSourcesCore::UserDataOperator::OPSPACE
@ OPSPACE
operator do Work is execute on space data
Definition: ForcesAndSourcesCore.hpp:570
MoFEM::OpSetHOWeights::OpSetHOWeights
OpSetHOWeights(boost::shared_ptr< VectorDouble > det_ptr)
Definition: HODataOperators.hpp:161
MoFEM::OpHOSetCovariantPiolaTransformOnFace3D::normalsAtPts
boost::shared_ptr< MatrixDouble > normalsAtPts
Definition: HODataOperators.hpp:357
MoFEM::addHOOpsVol
MoFEMErrorCode addHOOpsVol(const std::string field, E &e, bool h1, bool hcurl, bool hdiv, bool l2)
Definition: HODataOperators.hpp:676
CHK_THROW_MESSAGE
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
Definition: definitions.h:609
ApproximationBaseNames
const static char *const ApproximationBaseNames[]
Definition: definitions.h:72
L2
@ L2
field with C-1 continuity
Definition: definitions.h:88
MoFEM::OpHOSetContravariantPiolaTransformOnEdge3D::doWork
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
Definition: HODataOperators.cpp:609
MoFEM::OpSetHOContravariantPiolaTransform
Apply contravariant (Piola) transfer to Hdiv space for HO geometry.
Definition: HODataOperators.hpp:180
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
E
MoFEM::OpScaleBaseBySpaceInverseOfMeasure::doWork
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
Definition: HODataOperators.cpp:704
MoFEM::OpSetHOCovariantPiolaTransform::piolaDiffN
MatrixDouble piolaDiffN
Definition: HODataOperators.hpp:228
MoFEM::OpHOSetContravariantPiolaTransformOnFace3D
transform Hdiv base fluxes from reference element to physical triangle
Definition: HODataOperators.hpp:298
MoFEM::OpGetHOTangentsOnEdge::tangentsAtPts
boost::shared_ptr< MatrixDouble > tangentsAtPts
Definition: HODataOperators.hpp:382
MoFEM::OpSetInvJacToScalarBasesBasic
Definition: UserDataOperators.hpp:2835
MoFEM::EntitiesFieldData::EntData::getFieldData
const VectorDouble & getFieldData() const
get dofs values
Definition: EntitiesFieldData.hpp:1254
MoFEM::OpSetHOContravariantPiolaTransform::piolaN
MatrixDouble piolaN
Definition: HODataOperators.hpp:200
MoFEM::OpSetHOWeightsOnFace::OpSetHOWeightsOnFace
OpSetHOWeightsOnFace()
Definition: HODataOperators.hpp:124
MoFEM::OpSetHOInvJacVectorBase::diffNinvJac
MatrixDouble diffNinvJac
Definition: HODataOperators.hpp:113
MoFEM::OpCalculateHOJacForVolume::OpCalculateHOJacForVolume
OpCalculateHOJacForVolume(boost::shared_ptr< MatrixDouble > jac_ptr)
Definition: HODataOperators.cpp:9
MoFEM::EntitiesFieldData::EntData::getDiffN
MatrixDouble & getDiffN(const FieldApproximationBase base)
get derivatives of base functions
Definition: EntitiesFieldData.hpp:1329
MoFEM::OpScaleBaseBySpaceInverseOfMeasure::detJacPtr
boost::shared_ptr< VectorDouble > detJacPtr
Definition: HODataOperators.hpp:404
MoFEM::OpScaleBaseBySpaceInverseOfMeasure
Scale base functions by inverses of measure of element.
Definition: HODataOperators.hpp:390
FTensor::levi_civita
constexpr std::enable_if<(Dim0<=2 &&Dim1<=2), Tensor2_Expr< Levi_Civita< T >, T, Dim0, Dim1, i, j > >::type levi_civita(const Index< i, Dim0 > &, const Index< j, Dim1 > &)
levi_civita functions to make for easy adhoc use
Definition: Levi_Civita.hpp:617
MoFEM::DataOperator::doVertices
bool & doVertices
\deprectaed If false skip vertices
Definition: DataOperators.hpp:90
MoFEM::OpCalculateHOCoords
Calculate HO coordinates at gauss points.
Definition: HODataOperators.hpp:41
MoFEM::OpCalculateHOJacForFaceImpl
Calculate jacobian for face element.
Definition: HODataOperators.hpp:240
FIELD_DIM
constexpr int FIELD_DIM
Definition: child_and_parent.cpp:15
MoFEM::OpSetInvJacToScalarBasesBasic::OpSetInvJacToScalarBasesBasic
OpSetInvJacToScalarBasesBasic(FieldSpace space, boost::shared_ptr< MatrixDouble > inv_jac_ptr)
Definition: UserDataOperators.cpp:10
MoFEM::OpCalculateHOJac
Definition: HODataOperators.hpp:267
MoFEM::OpSetHOInvJacToScalarBasesImpl
Definition: HODataOperators.hpp:55
MoFEM::EntitiesFieldData::EntData::getFTensor0N
FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > getFTensor0N(const FieldApproximationBase base)
Get base function as Tensor0.
Definition: EntitiesFieldData.hpp:1502
MoFEM::OpSetHOCovariantPiolaTransform::OpSetHOCovariantPiolaTransform
OpSetHOCovariantPiolaTransform(const FieldSpace space, boost::shared_ptr< MatrixDouble > jac_inv_ptr)
Definition: HODataOperators.hpp:209
MoFEM::OpCalculateHOJacForVolume::doWork
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
Definition: HODataOperators.cpp:22
FieldSpace
FieldSpace
approximation spaces
Definition: definitions.h:82
MoFEM::OpSetHOContravariantPiolaTransform::detPtr
boost::shared_ptr< VectorDouble > detPtr
Definition: HODataOperators.hpp:197
FTensor::Number< 0 >
MoFEM::OpHOSetCovariantPiolaTransformOnFace3D::OpHOSetCovariantPiolaTransformOnFace3D
OpHOSetCovariantPiolaTransformOnFace3D(const FieldSpace space, boost::shared_ptr< MatrixDouble > normals_at_pts=nullptr, boost::shared_ptr< MatrixDouble > tangent1_at_pts=nullptr, boost::shared_ptr< MatrixDouble > tangent2_at_pts=nullptr)
Definition: HODataOperators.hpp:339
MoFEM::OpCalculateHOJacForFaceImpl< 2 >::jacPtr
boost::shared_ptr< MatrixDouble > jacPtr
Definition: HODataOperators.hpp:252
MoFEM::OpCalculateHOJacForFaceImpl< 3 >
Definition: HODataOperators.hpp:256
MoFEM::OpHOSetContravariantPiolaTransformOnEdge3D::OpHOSetContravariantPiolaTransformOnEdge3D
OpHOSetContravariantPiolaTransformOnEdge3D(const FieldSpace space=HCURL, boost::shared_ptr< MatrixDouble > tangent_at_pts=nullptr)
Definition: HODataOperators.hpp:320
MoFEM
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
MoFEM::OpCalculateHOJacForVolume
Calculate jacobian on Hex or other volume which is not simplex.
Definition: HODataOperators.hpp:20
MoFEM::OpSetHOCovariantPiolaTransform
Apply covariant (Piola) transfer to Hcurl space for HO geometry.
Definition: HODataOperators.hpp:206
MoFEM::ForcesAndSourcesCore::UserDataOperator
Definition: ForcesAndSourcesCore.hpp:549
MoFEM::FaceElementForcesAndSourcesCore::UserDataOperator
default operator for TRI element
Definition: FaceElementForcesAndSourcesCore.hpp:94
MoFEM::OpScaleBaseBySpaceInverseOfMeasure::fieldSpace
FieldSpace fieldSpace
Definition: HODataOperators.hpp:403
MoFEM::OpSetHOInvJacToScalarBases
Set inverse jacobian to base functions.
Definition: HODataOperators.hpp:73
convert.type
type
Definition: convert.py:64
MoFEM::OpSetHOWeights::doWork
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
Definition: HODataOperators.cpp:205
MoFEM::OpHOSetCovariantPiolaTransformOnFace3D::piolaN
MatrixDouble piolaN
Definition: HODataOperators.hpp:361
MoFEM::OpSetHOContravariantPiolaTransform::doWork
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
Definition: HODataOperators.cpp:229
MoFEM::EntitiesFieldData::EntData::getBase
FieldApproximationBase & getBase()
Get approximation base.
Definition: EntitiesFieldData.hpp:1313
MoFEM::OpSetHOCovariantPiolaTransform::jacInvPtr
boost::shared_ptr< MatrixDouble > jacInvPtr
Definition: HODataOperators.hpp:225
MoFEM::OpSetHOInvJacVectorBase::doWork
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
Definition: HODataOperators.cpp:109
MoFEM::OpSetHOWeightsOnFace
Modify integration weights on face to take in account higher-order geometry.
Definition: HODataOperators.hpp:122
MoFEM::FaceElementForcesAndSourcesCore
Face finite element.
Definition: FaceElementForcesAndSourcesCore.hpp:23
MoFEM::OpCalculateHOCoords::doWork
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
Definition: HODataOperators.hpp:469
MoFEM::OpCalculateHOJacForFaceImpl< 2 >
Definition: HODataOperators.hpp:243
MoFEM::OpSetHOInvJacVectorBase::OpSetHOInvJacVectorBase
OpSetHOInvJacVectorBase(const FieldSpace space, boost::shared_ptr< MatrixDouble > inv_jac_ptr)
Definition: HODataOperators.hpp:93
MoFEM::OpGetHONormalsOnFace
Calculate normals at Gauss points of triangle element.
Definition: HODataOperators.hpp:281
MoFEM::OpHOSetCovariantPiolaTransformOnFace3D::diffPiolaN
MatrixDouble diffPiolaN
Definition: HODataOperators.hpp:362
MoFEM::EntitiesFieldData::EntData::getFieldDofs
const VectorDofs & getFieldDofs() const
get dofs data stature FEDofEntity
Definition: EntitiesFieldData.hpp:1269
MoFEM::OpHOSetCovariantPiolaTransformOnFace3D::tangent1AtPts
boost::shared_ptr< MatrixDouble > tangent1AtPts
Definition: HODataOperators.hpp:358
MoFEM::OpHOSetContravariantPiolaTransformOnFace3D::OpHOSetContravariantPiolaTransformOnFace3D
OpHOSetContravariantPiolaTransformOnFace3D(const FieldSpace space, boost::shared_ptr< MatrixDouble > normals_at_gauss_pts=nullptr)
Definition: HODataOperators.hpp:301
MoFEM::ForcesAndSourcesCore::UserDataOperator
friend class UserDataOperator
Definition: ForcesAndSourcesCore.hpp:482
i
FTensor::Index< 'i', SPACE_DIM > i
Definition: hcurl_divergence_operator_2d.cpp:27
MoFEM::OpGetHOTangentsOnEdge::doWork
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
Definition: HODataOperators.hpp:614
MoFEM::OpSetHOWeightsOnSubDim
Definition: HODataOperators.hpp:145
MoFEM::OpSetHOInvJacVectorBase::invJacPtr
boost::shared_ptr< MatrixDouble > invJacPtr
Definition: HODataOperators.hpp:111
MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator
Definition: VolumeElementForcesAndSourcesCore.hpp:108
field_name
constexpr auto field_name
Definition: poisson_2d_homogeneous.cpp:13
FTensor::Index< 'i', FIELD_DIM >
MoFEM::AddHOOps
Add operators pushing bases from local to physical configuration.
Definition: HODataOperators.hpp:415
MoFEM::EdgeElementForcesAndSourcesCore::UserDataOperator
default operator for EDGE element
Definition: EdgeElementForcesAndSourcesCore.hpp:68
MoFEM::ForcesAndSourcesCore
structure to get information form mofem into EntitiesFieldData
Definition: ForcesAndSourcesCore.hpp:22
MoFEM::OpCalculateHOJacForVolume::jacPtr
boost::shared_ptr< MatrixDouble > jacPtr
Definition: HODataOperators.hpp:29
MoFEM::OpSetHOWeightsOnEdge::doWork
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
Definition: HODataOperators.cpp:180
MoFEM::DataOperator::doEdges
bool & doEdges
\deprectaed If false skip edges
Definition: DataOperators.hpp:91
MoFEM::OpHOSetCovariantPiolaTransformOnFace3D
transform Hcurl base fluxes from reference element to physical triangle
Definition: HODataOperators.hpp:336
MoFEM::EntitiesFieldData::EntData::getFTensor1DiffN
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1DiffN(const FieldApproximationBase base)
Get derivatives of base functions.
Definition: EntitiesFieldData.cpp:526
MoFEM::OpHOSetContravariantPiolaTransformOnEdge3D::tangentAtGaussPts
boost::shared_ptr< MatrixDouble > tangentAtGaussPts
Definition: HODataOperators.hpp:330
MoFEM::EntitiesFieldData::EntData::getFTensor1FieldData
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1FieldData()
Return FTensor of rank 1, i.e. vector from filed data coefficients.
Definition: EntitiesFieldData.hpp:1288
MoFEM::OpGetHOTangentsOnEdge::OpGetHOTangentsOnEdge
OpGetHOTangentsOnEdge(std::string field_name, boost::shared_ptr< MatrixDouble > tangents_at_pts=nullptr)
Definition: HODataOperators.hpp:372
MoFEM::OpHOSetCovariantPiolaTransformOnFace3D::tangent2AtPts
boost::shared_ptr< MatrixDouble > tangent2AtPts
Definition: HODataOperators.hpp:359
MoFEM::EntitiesFieldData::EntData::getN
MatrixDouble & getN(const FieldApproximationBase base)
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....
Definition: EntitiesFieldData.hpp:1318
MoFEM::OpSetHOCovariantPiolaTransform::doWork
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
Definition: HODataOperators.cpp:334
MoFEM::OpSetHOInvJacVectorBase::diffHdivInvJac
MatrixDouble diffHdivInvJac
Definition: HODataOperators.hpp:112
j
FTensor::Index< 'j', 3 > j
Definition: matrix_function.cpp:19
MoFEM::OpSetHOWeights
Set inverse jacobian to base functions.
Definition: HODataOperators.hpp:159
MoFEM::OpScaleBaseBySpaceInverseOfMeasure::OpScaleBaseBySpaceInverseOfMeasure
OpScaleBaseBySpaceInverseOfMeasure(const FieldSpace space, boost::shared_ptr< VectorDouble > det_jac_ptr, boost::shared_ptr< double > scale_det_ptr=nullptr)
Definition: HODataOperators.cpp:693
HCURL
@ HCURL
field with continuous tangents
Definition: definitions.h:86
MOFEM_DATA_INCONSISTENCY
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
MoFEM::OpGetHONormalsOnFace::OpGetHONormalsOnFace
OpGetHONormalsOnFace(std::string field_name)
Definition: HODataOperators.hpp:284
MoFEM::OpHOSetContravariantPiolaTransformOnFace3D::normalsAtGaussPts
boost::shared_ptr< MatrixDouble > normalsAtGaussPts
Definition: HODataOperators.hpp:311
m
FTensor::Index< 'm', 3 > m
Definition: shallow_wave.cpp:80
MoFEM::OpSetHOInvJacVectorBase
transform local reference derivatives of shape function to global derivatives if higher order geometr...
Definition: HODataOperators.hpp:91
MoFEM::OpSetHOWeightsOnEdge
Modify integration weights on face to take in account higher-order geometry.
Definition: HODataOperators.hpp:136
MoFEMFunctionBeginHot
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:453
MoFEM::OpCalculateHOCoords::OpCalculateHOCoords
OpCalculateHOCoords(const std::string field_name)
Definition: HODataOperators.hpp:43
k
FTensor::Index< 'k', 3 > k
Definition: matrix_function.cpp:20
MoFEM::OpHOSetContravariantPiolaTransformOnFace3D::doWork
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
Definition: HODataOperators.cpp:459
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
MoFEM::OpSetHOContravariantPiolaTransform::OpSetHOContravariantPiolaTransform
OpSetHOContravariantPiolaTransform(const FieldSpace space, boost::shared_ptr< VectorDouble > det_ptr, boost::shared_ptr< MatrixDouble > jac_ptr)
Definition: HODataOperators.hpp:183
HDIV
@ HDIV
field with continuous normal traction
Definition: definitions.h:87
MOFEM_NOT_IMPLEMENTED
@ MOFEM_NOT_IMPLEMENTED
Definition: definitions.h:32
MoFEM::OpSetHOWeightsOnFace::doWork
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
Definition: HODataOperators.cpp:153
MoFEM::OpSetInvJacSpaceForFaceImpl
Transform local reference derivatives of shape functions to global derivatives.
Definition: UserDataOperators.hpp:2833
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
MoFEM::OpHOSetCovariantPiolaTransformOnFace3D::doWork
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
Definition: HODataOperators.cpp:512
MoFEM::OpGetHONormalsOnFace::doWork
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
Definition: HODataOperators.hpp:499
MoFEM::OpScaleBaseBySpaceInverseOfMeasure::scaleDetPtr
boost::shared_ptr< double > scaleDetPtr
Definition: HODataOperators.hpp:405
MoFEM::ForcesAndSourcesCore::UserDataOperator::OPROW
@ OPROW
operator doWork function is executed on FE rows
Definition: ForcesAndSourcesCore.hpp:567
MoFEM::OpHOSetContravariantPiolaTransformOnEdge3D
transform Hcurl base fluxes from reference element to physical edge
Definition: HODataOperators.hpp:317