v0.15.5
Loading...
Searching...
No Matches
HODataOperators.hpp
Go to the documentation of this file.
1/** \file HODataOperators.hpp
2 * \brief High-order geometry operators for finite element computations
3 *
4 * This file provides operators for managing high-order (curved) geometry in
5 * finite element analysis. These operators handle geometric transformations,
6 * coordinate mappings, and basis function modifications required when elements
7 * have curved boundaries or use high-order shape functions.
8 *
9 * ## Main Functionality:
10 *
11 * **Jacobian Calculations**: Compute transformation matrices between reference
12 * and physical coordinates for curved elements (OpCalculateHOJacForVolume,
13 * OpCalculateHOJacForFace, OpGetHOTangentsOnEdge, OpGetHONormalsOnFace)
14 *
15 * **Coordinate Transformations**: Calculate physical coordinates at integration
16 * points using high-order approximations (OpCalculateHOCoords)
17 *
18 * **Basis Function Transformations**: Apply Piola transforms and inverse Jacobian
19 * operations to basis functions for H1, Hdiv, and Hcurl spaces
20 * (OpSetHOContravariantPiolaTransform, OpSetHOCovariantPiolaTransform)
21 *
22 * **Integration Weight Corrections**: Adjust integration weights to account for
23 * geometric distortions in curved elements (OpSetHOWeights, OpSetHOWeightsOnFace)
24 *
25 * **Pipeline Management**: AddHOOps template specializations automatically add
26 * appropriate HO operators to finite element pipelines based on element
27 * dimension and problem type
28 *
29 * ## Usage Pattern:
30 * 1. Use AddHOOps<FE_DIM, PROBLEM_DIM, SPACE_DIM>::add() to automatically
31 * add required HO operators to your finite element pipeline
32 * 2. Individual operators can be used directly for specialized applications
33 * 3. Operators handle both internal computations and external data access
34 * via shared pointers to Jacobian matrices and determinants
35 *
36 * ## Key Benefits:
37 * - Accurate integration on curved domains
38 * - Proper handling of geometric nonlinearity
39 * - Seamless integration with MoFEM pipeline system
40 * - Support for mixed finite element spaces (H1, Hdiv, Hcurl)
41 *
42 * @note Essential for problems involving curved boundaries, contact mechanics,
43 * fluid-structure interaction, and any application requiring geometric
44 * accuracy beyond linear approximation.
45 */
46
47#ifndef __HO_DATA_OPERATORS_HPP__
48#define __HO_DATA_OPERATORS_HPP__
49
50namespace MoFEM {
51
52/**
53 * @brief Calculate Jacobian on Hex or other volume elements which are not simplex
54 *
55 * This structure calculates the Jacobian matrix for high-order volume elements
56 * using geometrical coordinates of the element rather than field data. It's
57 * specifically designed for hexahedral and other non-simplex volume elements.
58 *
59 * @note Use OpCalculateVectorFieldGradient to calculate Jacobian from field data.
60 */
63
64 /**
65 * @brief Constructor for HO Jacobian calculation on volumes
66 *
67 * @param jac_ptr Shared pointer to matrix for storing calculated Jacobian values
68 */
69 OpCalculateHOJacForVolume(boost::shared_ptr<MatrixDouble> jac_ptr);
70
71 MoFEMErrorCode doWork(int side, EntityType type,
73
74private:
75 boost::shared_ptr<MatrixDouble> jacPtr;
76};
77
78/** \deprecated use OpCalculateHOJacForVolume
79 */
81
82/**
83 * @brief Calculate high-order coordinates at Gauss points
84 *
85 * This template structure computes high-order coordinate values at integration
86 * (Gauss) points using field approximations. It supports different field
87 * dimensions and is commonly used in high-order finite element methods.
88 *
89 * @tparam FIELD_DIM Dimension of the coordinate field (default: 3)
90 */
91template <int FIELD_DIM = 3>
93
94 /**
95 * @brief Constructor for HO coordinates calculation
96 *
97 * @param field_name Name of the field containing coordinate information
98 */
101
102 MoFEMErrorCode doWork(int side, EntityType type,
104};
105
106/**
107 * @deprecated This class should be DIM = 3 specialization when default
108 * parameter is removed
109 *
110 */
118
119/**
120 * \brief Set inverse jacobian to base functions
121 *
122 * \deprecated Version with default variant DIM = 3 is deprecated, keep for
123 * back compatibility. Should be removed in future. Use of default variant make
124 * code implicit, what can be source of some hidden error. Explict interface is
125 * better, when user see and control outcome, and is aware of existing variants.
126 *
127 */
128template <int DIM = 3, int DERIVATIVE = 1>
130 using OpSetHOInvJacToScalarBasesImpl::OpSetHOInvJacToScalarBasesImpl;
131};
132
133template <>
135 : public OpSetHOInvJacToScalarBases<3, 1> {
136 using OpSetHOInvJacToScalarBases<3, 1>::OpSetHOInvJacToScalarBases;
137 MoFEMErrorCode doWork(int side, EntityType type,
139};
140
141template <>
143 : public OpSetInvJacSpaceForFaceImpl<2, 1> {
144
145 using OpSetInvJacSpaceForFaceImpl<2, 1>::OpSetInvJacSpaceForFaceImpl;
146};
147
148template <>
150 : public OpSetInvJacSpaceForFaceImpl<2, 2> {
151 using OpSetInvJacSpaceForFaceImpl<2, 2>::OpSetInvJacSpaceForFaceImpl;
152};
153
154
155
156/**
157 * @brief Transform local reference derivatives of shape functions to global derivatives
158 *
159 * This structure transforms derivatives of shape functions from local reference
160 * coordinates to global coordinates when higher-order geometry is used. It's
161 * essential for vector-based finite element spaces requiring coordinate transformations.
162 *
163 * @ingroup mofem_forces_and_sources
164 */
166
167 /**
168 * @brief Constructor for HO inverse Jacobian transformation for vector bases
169 *
170 * @param space Field space specification for the transformation
171 * @param inv_jac_ptr Shared pointer to matrix containing inverse Jacobian data
172 */
174 boost::shared_ptr<MatrixDouble> inv_jac_ptr)
175 : ForcesAndSourcesCore::UserDataOperator(space), invJacPtr(inv_jac_ptr) { if (!invJacPtr)
177 "Pointer for invJacPtr not allocated");
178
179 doVertices = false;
180 if (space == HDIV)
181 doEdges = false;
182 }
183
184
185 MoFEMErrorCode doWork(int side, EntityType type,
187
188private:
189 boost::shared_ptr<MatrixDouble> invJacPtr;
192};
193
194/**
195 * @brief Modify integration weights on face to take into account higher-order geometry
196 *
197 * This structure adjusts integration weights for face elements when higher-order
198 * geometry is present. It accounts for geometric distortions introduced by
199 * curved boundaries and high-order shape functions.
200 *
201 * @ingroup mofem_forces_and_sources_tri_element
202 */
205
206 /**
207 * @brief Default constructor for HO weights on face
208 */
211
212 MoFEMErrorCode doWork(int side, EntityType type,
214};
215
216/**
217 * @brief Modify integration weights on edge to take into account higher-order geometry
218 *
219 * This structure adjusts integration weights for edge elements when higher-order
220 * geometry is present. It accounts for geometric distortions introduced by
221 * curved edges and high-order shape functions along the edge.
222 *
223 * @ingroup mofem_forces_and_sources_tri_element
224 */
227
228 /**
229 * @brief Default constructor for HO weights on edge
230 */
233
234 MoFEMErrorCode doWork(int side, EntityType type,
236};
237
238template <int SPACE_DIM>
240
244
248
249/**
250 * @brief Set integration weights for higher-order elements
251 *
252 * This structure sets integration weights by computing the determinant of the
253 * Jacobian matrix for higher-order elements. It ensures proper integration
254 * when dealing with curved geometries and high-order shape functions.
255 */
257
258 /**
259 * @brief Constructor for HO weights calculation
260 *
261 * @param det_ptr Shared pointer to vector for storing determinant values
262 */
263 OpSetHOWeights(boost::shared_ptr<VectorDouble> det_ptr)
265 detPtr(det_ptr) {
266 if (!detPtr)
268 "Pointer for detPtr not allocated");
269 }
270
271 MoFEMErrorCode doWork(int side, EntityType type,
273
274private:
275 boost::shared_ptr<VectorDouble> detPtr;
276};
277
278/**
279 * @brief Apply contravariant (Piola) transformation to Hdiv space for HO geometry
280 *
281 * This structure applies the contravariant Piola transformation to basis functions
282 * in Hdiv spaces when higher-order geometry is present. This transformation is
283 * essential for preserving the divergence properties of vector fields across
284 * element boundaries in curved geometries.
285 *
286 * @ingroup mofem_forces_and_sources
287 */
290
291 /**
292 * @brief Constructor for contravariant Piola transformation
293 *
294 * @param space Field space specification (typically HDIV)
295 * @param det_ptr Shared pointer to determinant values for transformation
296 * @param jac_ptr Shared pointer to Jacobian matrix for transformation
297 */
299 boost::shared_ptr<VectorDouble> det_ptr,
300 boost::shared_ptr<MatrixDouble> jac_ptr)
302 jacPtr(jac_ptr) {
303 doVertices = false;
304 if (space == HDIV)
305 doEdges = false;
306 }
307
308 MoFEMErrorCode doWork(int side, EntityType type,
310
311private:
312 boost::shared_ptr<VectorDouble> detPtr;
313 boost::shared_ptr<MatrixDouble> jacPtr;
314
317};
318
319/**
320 * @brief Apply covariant (Piola) transformation to Hcurl space for HO geometry
321 *
322 * This structure applies the covariant Piola transformation to basis functions
323 * in Hcurl spaces when higher-order geometry is present. This transformation
324 * preserves the curl properties of vector fields and maintains tangential
325 * continuity across element boundaries in curved geometries.
326 */
329
330 /**
331 * @brief Constructor for covariant Piola transformation
332 *
333 * @param space Field space specification (typically HCURL)
334 * @param jac_inv_ptr Shared pointer to inverse Jacobian matrix for transformation
335 */
337 boost::shared_ptr<MatrixDouble> jac_inv_ptr)
339 jacInvPtr(jac_inv_ptr) {
340 doVertices = false;
341 if (space == HDIV)
342 doEdges = false;
343 if (!jacInvPtr)
345 "Pointer for jacPtr not allocated");
346 }
347
348 MoFEMErrorCode doWork(int side, EntityType type,
350
351private:
352 boost::shared_ptr<MatrixDouble> jacInvPtr;
353
356};
357
358/** \brief Calculate jacobian for face element
359
360 It is assumed that face element is XY plane. Applied
361 only for 2d problems and 2d problems embedded in 3d space
362
363 \note If you push operators for HO normal befor this operator, HO geometry is
364 taken into account when you calculate jacobian.
365
366*/
367template <int DIM> struct OpCalculateHOJacForFaceImpl;
368
369template <>
372
373 OpCalculateHOJacForFaceImpl(boost::shared_ptr<MatrixDouble> jac_ptr);
374
375 MoFEMErrorCode doWork(int side, EntityType type,
377
378protected:
379 boost::shared_ptr<MatrixDouble> jacPtr;
380};
381
382template <>
384
385 using OpCalculateHOJacForFaceImpl<2>::OpCalculateHOJacForFaceImpl;
386
387 MoFEMErrorCode doWork(int side, EntityType type,
389};
390
393
394template <int DIM> struct OpCalculateHOJac;
395
399
400template <> struct OpCalculateHOJac<2> : public OpCalculateHOJacForFaceImpl<2> {
401 using OpCalculateHOJacForFaceImpl<2>::OpCalculateHOJacForFaceImpl;
402};
403
404/**
405 * @brief Calculate normal vectors at Gauss points of face elements
406 *
407 * This template structure computes normal vectors at integration points for
408 * face elements using high-order geometry approximation. It's essential for
409 * boundary conditions, surface integrals, and flux calculations in finite
410 * element methods.
411 *
412 * @tparam FIELD_DIM Dimension of the field for normal calculation (default: 3)
413 *
414 * @ingroup mofem_forces_and_sources
415 */
416template <int FIELD_DIM = 3>
419
421
422 /**
423 * @brief Constructor for HO normal vectors calculation on faces
424 *
425 * @param field_name Field name where normal vectors will be stored
426 * @param add_field If true, add to existing field values; if false, overwrite (default: false)
427 */
428 OpGetHONormalsOnFace(std::string field_name, bool add_field = false)
429 : OP(field_name, OPROW), addField(add_field) {}
430
432 boost::shared_ptr<MatrixDouble> tangent1_at_pts,
433 boost::shared_ptr<MatrixDouble> tangent2_at_pts,
434 boost::shared_ptr<MatrixDouble> normals_at_pts,
435 bool add_field = false)
436 : OP(field_name, OPROW), tangent1AtPts(tangent1_at_pts),
437 tangent2AtPts(tangent2_at_pts), normalsAtPts(normals_at_pts),
438 addField(add_field) {}
439
440 MoFEMErrorCode doWork(int side, EntityType type,
442
443private:
447
451
455
456 boost::shared_ptr<MatrixDouble> tangent1AtPts;
457 boost::shared_ptr<MatrixDouble> tangent2AtPts;
458 boost::shared_ptr<MatrixDouble> normalsAtPts;
459
460 bool addField = false; ///< if true add field values, do not overwrite
461};
462
463/**
464 * @brief Transform Hdiv base fluxes from reference element to physical face in 3D
465 *
466 * This structure applies contravariant Piola transformation to Hdiv basis
467 * functions on 3D face elements. It transforms flux values from reference
468 * coordinates to physical coordinates while preserving divergence properties.
469 *
470 * @note Matrix which keeps normal vectors is assumed to have three columns,
471 * and number of rows should be equal to number of integration points.
472 *
473 * @ingroup mofem_forces_and_sources
474 */
477
478 /**
479 * @brief Constructor for contravariant Piola transformation on 3D faces
480 *
481 * @param space Field space specification (typically HDIV)
482 * @param normals_at_gauss_pts Shared pointer to matrix containing normal vectors at integration points (optional)
483 */
485 const FieldSpace space,
486 boost::shared_ptr<MatrixDouble> normals_at_gauss_pts = nullptr)
488 normalsAtGaussPts(normals_at_gauss_pts) {}
489
490 MoFEMErrorCode doWork(int side, EntityType type,
492
493private:
494 boost::shared_ptr<MatrixDouble> normalsAtGaussPts;
495};
496
497/**
498 * @brief Transform Hcurl base functions from reference element to physical edge in 3D
499 *
500 * This structure applies contravariant Piola transformation to Hcurl basis
501 * functions on 3D edge elements. It transforms curl-conforming basis functions
502 * from reference coordinates to physical coordinates while preserving tangential
503 * continuity across element boundaries.
504 *
505 * @ingroup mofem_forces_and_sources
506 */
509
510 /**
511 * @brief Constructor for contravariant Piola transformation on 3D edges
512 *
513 * @param space Field space specification (default: HCURL)
514 * @param tangent_at_pts Shared pointer to matrix containing tangent vectors at integration points (optional)
515 */
517 const FieldSpace space = HCURL,
518 boost::shared_ptr<MatrixDouble> tangent_at_pts = nullptr)
520 tangentAtGaussPts(tangent_at_pts) {}
521
522 MoFEMErrorCode doWork(int side, EntityType type,
524
525private:
526 boost::shared_ptr<MatrixDouble> tangentAtGaussPts;
527};
528
529/**
530 * @brief Transform Hcurl base functions from reference element to physical face in 3D
531 *
532 * This structure applies covariant Piola transformation to Hcurl basis functions
533 * on 3D face elements. It transforms curl-conforming basis functions from
534 * reference coordinates to physical coordinates using normal and tangent vectors,
535 * preserving the curl properties of the field.
536 *
537 * @ingroup mofem_forces_and_sources
538 */
541
542 /**
543 * @brief Constructor for covariant Piola transformation on 3D faces
544 *
545 * @param space Field space specification (typically HCURL)
546 * @param normals_at_pts Shared pointer to matrix containing normal vectors at integration points (optional)
547 * @param tangent1_at_pts Shared pointer to matrix containing first tangent vectors at integration points (optional)
548 * @param tangent2_at_pts Shared pointer to matrix containing second tangent vectors at integration points (optional)
549 */
551 const FieldSpace space,
552 boost::shared_ptr<MatrixDouble> normals_at_pts = nullptr,
553 boost::shared_ptr<MatrixDouble> tangent1_at_pts = nullptr,
554 boost::shared_ptr<MatrixDouble> tangent2_at_pts = nullptr)
556 normalsAtPts(normals_at_pts), tangent1AtPts(tangent1_at_pts),
557 tangent2AtPts(tangent2_at_pts) {
558 if (normals_at_pts || tangent1_at_pts || tangent2_at_pts)
559 if (normals_at_pts && tangent1_at_pts && tangent2_at_pts)
561 "All elements in constructor have to set pointer");
562 }
563
564 MoFEMErrorCode doWork(int side, EntityType type,
566
567private:
568 boost::shared_ptr<MatrixDouble> normalsAtPts;
569 boost::shared_ptr<MatrixDouble> tangent1AtPts;
570 boost::shared_ptr<MatrixDouble> tangent2AtPts;
571
574};
575
576/**
577 * @brief Calculate tangent vector on edge from HO geometry approximation
578 *
579 * This template structure computes tangent vectors along edges using high-order
580 * geometry approximation. It extracts tangent information at integration points
581 * and stores the results for further processing in finite element calculations.
582 *
583 * @tparam FIELD_DIM Dimension of the field (default: 3)
584 *
585 * @ingroup mofem_forces_and_sources
586 */
587template <int FIELD_DIM = 3>
590
591 /**
592 * @brief Constructor for HO tangents calculation on edges
593 *
594 * @param field_name Name of the field used for tangent calculation
595 * @param tangents_at_pts Shared pointer to matrix for storing calculated tangents (optional)
596 */
598 std::string field_name,
599 boost::shared_ptr<MatrixDouble> tangents_at_pts = nullptr)
601 tangentsAtPts(tangents_at_pts) {}
602
603 MoFEMErrorCode doWork(int side, EntityType type,
605
606private:
607 boost::shared_ptr<MatrixDouble> tangentsAtPts;
608};
609
610/**
611 * @brief Scale base functions by inverses of measure of element
612 *
613 * @tparam OP
614 */
617
619
621 const FieldSpace space, const FieldApproximationBase base,
622 boost::shared_ptr<VectorDouble> det_jac_ptr,
623 boost::shared_ptr<double> scale_det_ptr = nullptr);
624
625 MoFEMErrorCode doWork(int side, EntityType type,
627
628private:
631 boost::shared_ptr<VectorDouble> detJacPtr;
632 boost::shared_ptr<double> scaleDetPtr = nullptr;
633};
634
635/**
636 * @brief Add operators pushing bases from local to physical configuration
637 *
638 * @tparam FE_DIM dimension of element
639 * @tparam PROBLEM_DIM problem dimension
640 * @tparam SPACE_DIM space dimension
641 */
642template <int FE_DIM, int PROBLEM_DIM, int SPACE_DIM> struct AddHOOps;
643
644/**
645 * @brief Specialization for 2D face elements in 2D problems with 2D space
646 *
647 * This specialization adds high-order operators for 2D face elements operating
648 * in 2D problems within 2D coordinate space. Commonly used for planar face
649 * elements in 2D finite element analysis.
650 */
651template <> struct AddHOOps<2, 2, 2> {
652 AddHOOps() = delete;
653
654 /**
655 * @brief Add high-order operators to the pipeline for 2D face elements
656 *
657 * @param pipeline Pipeline to add operators to
658 * @param spaces Vector of field spaces to be handled
659 * @param geom_field_name Name of the geometry field (default: empty string)
660 * @param jac Shared pointer to Jacobian matrix. If provided, values will be set by this function.
661 * If null, Jacobian will be handled internally.
662 * @param det Shared pointer to determinant vector. If provided, values will be set by this function.
663 * If null, determinant will be handled internally.
664 * @param inv_jac Shared pointer to inverse Jacobian matrix. If provided, values will be set by this function.
665 * If null, inverse Jacobian will be handled internally.
666 * @return MoFEMErrorCode Error code indicating success or failure
667 */
668 static MoFEMErrorCode
669 add(boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pipeline,
670 std::vector<FieldSpace> spaces, std::string geom_field_name = "",
671 boost::shared_ptr<MatrixDouble> jac = nullptr,
672 boost::shared_ptr<VectorDouble> det = nullptr,
673 boost::shared_ptr<MatrixDouble> inv_jac = nullptr);
674};
675
676/**
677 * @brief Specialization for 2D face elements in 2D problems embedded in 3D space
678 *
679 * This specialization adds high-order operators for 2D face elements operating
680 * in 2D problems but embedded within 3D coordinate space. Typically used for
681 * surface elements or membranes in 3D environments.
682 */
683template <> struct AddHOOps<2, 2, 3> {
684 AddHOOps() = delete;
685
686 /**
687 * @brief Add high-order operators to the pipeline for 2D face elements embedded in 3D
688 *
689 * @param pipeline Pipeline to add operators to
690 * @param spaces Vector of field spaces to be handled
691 * @param geom_field_name Name of the geometry field (default: empty string)
692 * @param jac Shared pointer to Jacobian matrix. If provided, values will be set by this function.
693 * If null, Jacobian will be handled internally.
694 * @param det Shared pointer to determinant vector. If provided, values will be set by this function.
695 * If null, determinant will be handled internally.
696 * @param inv_jac Shared pointer to inverse Jacobian matrix. If provided, values will be set by this function.
697 * If null, inverse Jacobian will be handled internally.
698 * @return MoFEMErrorCode Error code indicating success or failure
699 */
700 static MoFEMErrorCode
701 add(boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pipeline,
702 std::vector<FieldSpace> spaces, std::string geom_field_name = "",
703 boost::shared_ptr<MatrixDouble> jac = nullptr,
704 boost::shared_ptr<VectorDouble> det = nullptr,
705 boost::shared_ptr<MatrixDouble> inv_jac = nullptr);
706};
707
708/**
709 * @brief Specialization for 1D edge elements in 2D problems with 2D space
710 *
711 * This specialization adds high-order operators for 1D edge elements operating
712 * in 2D problems within 2D coordinate space. Used for edge-based computations
713 * such as boundary conditions or line integrals in 2D finite element analysis.
714 */
715template <> struct AddHOOps<1, 2, 2> {
716 AddHOOps() = delete;
717
718 /**
719 * @brief Add high-order operators to the pipeline for 1D edge elements in 2D
720 *
721 * @param pipeline Pipeline to add operators to
722 * @param spaces Vector of field spaces to be handled
723 * @param geom_field_name Name of the geometry field (default: empty string)
724 * @return MoFEMErrorCode Error code indicating success or failure
725 *
726 * @note This specialization handles Jacobian, determinant, and inverse Jacobian internally.
727 */
728 static MoFEMErrorCode
729 add(boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pipeline,
730 std::vector<FieldSpace> spaces, std::string geom_field_name = "");
731};
732
733/**
734 * @brief Specialization for 1D edge elements in 3D problems with 3D space
735 *
736 * This specialization adds high-order operators for 1D edge elements operating
737 * in 3D problems within 3D coordinate space. Used for edge-based computations
738 * such as line integrals, boundary conditions, or wire-like structures in 3D
739 * finite element analysis.
740 */
741template <> struct AddHOOps<1, 3, 3> {
742 AddHOOps() = delete;
743
744 /**
745 * @brief Add high-order operators to the pipeline for 1D edge elements in 3D
746 *
747 * @param pipeline Pipeline to add operators to
748 * @param space Vector of field spaces to be handled
749 * @param geom_field_name Name of the geometry field (default: empty string)
750 * @return MoFEMErrorCode Error code indicating success or failure
751 *
752 * @note This specialization handles Jacobian, determinant, and inverse Jacobian internally.
753 */
754 static MoFEMErrorCode
755 add(boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pipeline,
756 std::vector<FieldSpace> space, std::string geom_field_name = "");
757};
758
759/**
760 * @brief Specialization for 2D face elements in 3D problems with 3D space
761 *
762 * This specialization adds high-order operators for 2D face elements operating
763 * in 3D problems within 3D coordinate space. Commonly used for surface elements,
764 * boundary conditions, membranes, or shell structures in 3D finite element
765 * analysis.
766 */
767template <> struct AddHOOps<2, 3, 3> {
768 AddHOOps() = delete;
769
770 /**
771 * @brief Add high-order operators to the pipeline for 2D face elements in 3D
772 *
773 * @param pipeline Pipeline to add operators to
774 * @param space Vector of field spaces to be handled
775 * @param geom_field_name Name of the geometry field (default: empty string)
776 * @return MoFEMErrorCode Error code indicating success or failure
777 *
778 * @note This specialization handles Jacobian, determinant, and inverse Jacobian internally.
779 */
780 static MoFEMErrorCode
781 add(boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pipeline,
782 std::vector<FieldSpace> space, std::string geom_field_name = "");
783};
784
785/**
786 * @brief Specialization for 3D volume elements in 3D problems with 3D space
787 *
788 * This specialization adds high-order operators for 3D volume elements operating
789 * in 3D problems within 3D coordinate space. This is the most common case for
790 * solid mechanics, heat transfer, and other 3D finite element applications
791 * involving volumetric elements like tetrahedra, hexahedra, and prisms.
792 */
793template <> struct AddHOOps<3, 3, 3> {
794 AddHOOps() = delete;
795
796 /**
797 * @brief Add high-order operators to the pipeline for 3D volume elements
798 *
799 * @param pipeline Pipeline to add operators to
800 * @param space Vector of field spaces to be handled
801 * @param geom_field_name Name of the geometry field (default: empty string)
802 * @param jac Shared pointer to Jacobian matrix. If provided, values will be set by this function.
803 * If null, Jacobian will be handled internally.
804 * @param det Shared pointer to determinant vector. If provided, values will be set by this function.
805 * If null, determinant will be handled internally.
806 * @param inv_jac Shared pointer to inverse Jacobian matrix. If provided, values will be set by this function.
807 * If null, inverse Jacobian will be handled internally.
808 * @return MoFEMErrorCode Error code indicating success or failure
809 */
810 static MoFEMErrorCode
811 add(boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pipeline,
812 std::vector<FieldSpace> space, std::string geom_field_name = "",
813 boost::shared_ptr<MatrixDouble> jac = nullptr,
814 boost::shared_ptr<VectorDouble> det = nullptr,
815 boost::shared_ptr<MatrixDouble> inv_jac = nullptr);
816};
817
818template <int FIELD_DIM>
824 const auto nb_dofs = data.getFieldData().size() / FIELD_DIM;
825 if (nb_dofs) {
826 if (type == MBVERTEX)
827 getCoordsAtGaussPts().clear();
828 auto t_base = data.getFTensor0N();
829 auto t_coords = getFTensor1CoordsAtGaussPts();
830 const auto nb_integration_pts = data.getN().size1();
831 const auto nb_base_functions = data.getN().size2();
832 for (auto gg = 0; gg != nb_integration_pts; ++gg) {
833 auto t_dof = data.getFTensor1FieldData<FIELD_DIM>();
834 size_t bb = 0;
835 for (; bb != nb_dofs; ++bb) {
836 t_coords(i) += t_base * t_dof(i);
837 ++t_dof;
838 ++t_base;
839 }
840 for (; bb < nb_base_functions; ++bb)
841 ++t_base;
842 ++t_coords;
843 }
844 }
846}
847
848template <int FIELD_DIM>
853
856
857 auto get_ftensor1 = [](MatrixDouble &m) {
859 &m(0, 0), &m(0, 1), &m(0, 2));
860 };
861
862 unsigned int nb_dofs = data.getFieldData().size();
863 if (nb_dofs != 0) {
864
865 int nb_gauss_pts = data.getN().size1();
866 auto &tangent1_at_gauss_pts = getTangent1AtGaussPts();
867 auto &tangent2_at_gauss_pts = getTangent2AtGaussPts();
868 tangent1_at_gauss_pts.resize(nb_gauss_pts, 3, false);
869 tangent2_at_gauss_pts.resize(nb_gauss_pts, 3, false);
870
871 switch (type) {
872 case MBVERTEX: {
873 if (!addField) {
874 tangent1_at_gauss_pts.clear();
875 tangent2_at_gauss_pts.clear();
876 }
877 }
878 case MBEDGE:
879 case MBTRI:
880 case MBQUAD: {
881
882#ifndef NDEBUG
883 if (2 * data.getN().size2() != data.getDiffN().size2()) {
884 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
885 "expected two derivative in local element coordinates");
886 }
887 if (nb_dofs % FIELD_DIM != 0) {
888 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
889 "expected that number of dofs is multiplicative of field "
890 "dimension");
891 }
892#endif
893
894 if (nb_dofs > FIELD_DIM * data.getN().size2()) {
895 unsigned int nn = 0;
896 for (; nn != nb_dofs; nn++) {
897 if (!data.getFieldDofs()[nn]->getActive())
898 break;
899 }
900 if (nn > FIELD_DIM * data.getN().size2()) {
901 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
902 "Data inconsistency for base %s",
904 } else {
905 nb_dofs = nn;
906 if (!nb_dofs)
908 }
909 }
910 const int nb_base_functions = data.getN().size2();
911 auto t_base = data.getFTensor0N();
912 auto t_diff_base = data.getFTensor1DiffN<2>();
913 auto t_t1 = get_ftensor1(tangent1_at_gauss_pts);
914 auto t_t2 = get_ftensor1(tangent2_at_gauss_pts);
915 for (int gg = 0; gg != nb_gauss_pts; ++gg) {
916 auto t_data = data.getFTensor1FieldData<FIELD_DIM>();
917 int bb = 0;
918 for (; bb != nb_dofs / FIELD_DIM; ++bb) {
920 t_t1(i) += t_data(i) * t_diff_base(N0);
921 t_t2(i) += t_data(i) * t_diff_base(N1);
922 ++t_data;
923 ++t_base;
924 ++t_diff_base;
925 }
926 for (; bb != nb_base_functions; ++bb) {
927 ++t_base;
928 ++t_diff_base;
929 }
930 ++t_t1;
931 ++t_t2;
932 }
933 } break;
934 default:
935 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
936 }
937 }
938
939 if (moab::CN::Dimension(type) == 2) {
940
941 auto &normal_at_gauss_pts = getNormalsAtGaussPts();
942 auto &tangent1_at_gauss_pts = getTangent1AtGaussPts();
943 auto &tangent2_at_gauss_pts = getTangent2AtGaussPts();
944
945 const auto nb_gauss_pts = tangent1_at_gauss_pts.size1();
946 normal_at_gauss_pts.resize(nb_gauss_pts, 3, false);
947
948 auto t_normal = get_ftensor1(normal_at_gauss_pts);
949 auto t_t1 = get_ftensor1(tangent1_at_gauss_pts);
950 auto t_t2 = get_ftensor1(tangent2_at_gauss_pts);
951 for (unsigned int gg = 0; gg != nb_gauss_pts; ++gg) {
952 FTensor::Index<'i', 3> i;
953 FTensor::Index<'j', 3> j;
954 FTensor::Index<'k', 3> k;
955 t_normal(j) = FTensor::levi_civita(i, j, k) * t_t1(k) * t_t2(i);
956 ++t_normal;
957 ++t_t1;
958 ++t_t2;
959 }
960 }
961
963}
964
965template <int FIELD_DIM>
970
971 auto get_tangent = [&]() -> MatrixDouble & {
972 if (tangentsAtPts)
973 return *tangentsAtPts;
974 else
975 return getTangentAtGaussPts();
976 };
977
978 auto &tangent = get_tangent();
979 int nb_gauss_pts = getGaussPts().size2();
980
981 if (type == MBVERTEX) {
982 tangent.resize(nb_gauss_pts, 3, false);
983 tangent.clear();
984 }
985
986 int nb_dofs = data.getFieldData().size();
987 if (nb_dofs != 0) {
988 const int nb_base_functions = data.getN().size2();
989 double *diff_base_function = &*data.getDiffN().data().begin();
990 auto tangent_at_gauss_pts =
991 getFTensor1FromPtr<FIELD_DIM, 3>(&*tangent.data().begin());
992
994 int size = nb_dofs / FIELD_DIM;
995 if (nb_dofs % FIELD_DIM) {
996 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
997 "Inconsistent number of dofs and field dimension");
998 }
999 for (int gg = 0; gg != nb_gauss_pts; ++gg) {
1000 auto field_data = data.getFTensor1FieldData<FIELD_DIM>();
1001 int bb = 0;
1002 for (; bb < size; ++bb) {
1003 tangent_at_gauss_pts(i) += field_data(i) * (*diff_base_function);
1004 ++field_data;
1005 ++diff_base_function;
1006 }
1007 for (; bb != nb_base_functions; ++bb) {
1008 ++diff_base_function;
1009 }
1010 ++tangent_at_gauss_pts;
1011 }
1012 }
1014}
1015
1016/**
1017 * @deprecated do not use this function, instead use AddHOOps.
1018 *
1019 * @tparam E
1020 * @param field
1021 * @param e
1022 * @param h1
1023 * @param hcurl
1024 * @param hdiv
1025 * @param l2
1026 * @return MoFEMErrorCode
1027 */
1028template <typename E>
1029MoFEMErrorCode addHOOpsVol(const std::string field, E &e, bool h1,
1030 bool hcurl, bool hdiv, bool l2) {
1031 std::vector<FieldSpace> spaces;
1032 if (h1)
1033 spaces.push_back(H1);
1034 if (hcurl)
1035 spaces.push_back(HCURL);
1036 if (hdiv)
1037 spaces.push_back(HDIV);
1038 if (l2)
1039 spaces.push_back(L2);
1040 return AddHOOps<3, 3, 3>::add(e.getOpPtrVector(), spaces, field);
1041}
1042
1043/**
1044 * @deprecated do not use this function, instead use AddHOOps.
1045 *
1046 * @tparam E
1047 * @param field
1048 * @param e
1049 * @param hcurl
1050 * @param hdiv
1051 * @return MoFEMErrorCode
1052 */
1053template <typename E>
1054MoFEMErrorCode addHOOpsFace3D(const std::string field, E &e,
1055 bool hcurl, bool hdiv) {
1056 std::vector<FieldSpace> spaces;
1057 if (hcurl)
1058 spaces.push_back(HCURL);
1059 if (hdiv)
1060 spaces.push_back(HDIV);
1061 return AddHOOps<2, 3, 3>::add(e.getOpPtrVector(), spaces, field);
1062}
1063
1064}; // namespace MoFEM
1065
1066#endif
constexpr int FIELD_DIM
FieldApproximationBase
approximation base
Definition definitions.h:58
#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
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ MOFEM_DATA_INCONSISTENCY
Definition definitions.h:31
@ MOFEM_NOT_IMPLEMENTED
Definition definitions.h:32
static const char *const ApproximationBaseNames[]
Definition definitions.h:72
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
FTensor::Index< 'i', SPACE_DIM > i
FTensor::Index< 'j', 3 > j
FTensor::Index< 'k', 3 > k
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
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
implementation of Data Operators for Forces and Sources
Definition Common.hpp:10
MoFEMErrorCode addHOOpsFace3D(const std::string field, E &e, bool hcurl, bool hdiv)
MoFEMErrorCode addHOOpsVol(const std::string field, E &e, bool h1, bool hcurl, bool hdiv, bool l2)
constexpr auto field_name
FTensor::Index< 'm', 3 > m
Add operators pushing bases from local to physical configuration.
bool & doVertices
\deprectaed If false skip vertices
bool & doEdges
\deprectaed If false skip edges
Data on single entity (This is passed as argument to DataOperator::doWork)
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1DiffN(const FieldApproximationBase base)
Get derivatives of base functions.
FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > getFTensor0N(const FieldApproximationBase base)
Get base function as Tensor0.
FieldApproximationBase & getBase()
Get approximation base.
MatrixDouble & getDiffN(const FieldApproximationBase base)
get derivatives of base functions
MatrixDouble & getN(const FieldApproximationBase base)
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....
const VectorDouble & getFieldData() const
Get DOF values on entity.
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1FieldData()
Return FTensor of rank 1, i.e. vector from field data coefficients.
const VectorDofs & getFieldDofs() const
Get DOF data structures (const version)
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.
@ OPROW
operator doWork function is executed on FE rows
@ OPSPACE
operator do Work is execute on space data
structure to get information from mofem into EntitiesFieldData
Calculate high-order coordinates at Gauss points.
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
OpCalculateHOCoords(const std::string field_name)
Constructor for HO coordinates calculation.
boost::shared_ptr< MatrixDouble > jacPtr
Calculate jacobian for face element.
Calculate Jacobian on Hex or other volume elements which are not simplex.
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.
Calculate normal vectors at Gauss points of face elements.
boost::shared_ptr< MatrixDouble > normalsAtPts
MatrixDouble & getTangent1AtGaussPts()
boost::shared_ptr< MatrixDouble > tangent1AtPts
boost::shared_ptr< MatrixDouble > tangent2AtPts
OpGetHONormalsOnFace(std::string field_name, boost::shared_ptr< MatrixDouble > tangent1_at_pts, boost::shared_ptr< MatrixDouble > tangent2_at_pts, boost::shared_ptr< MatrixDouble > normals_at_pts, bool add_field=false)
OpGetHONormalsOnFace(std::string field_name, bool add_field=false)
Constructor for HO normal vectors calculation on faces.
MatrixDouble & getTangent2AtGaussPts()
MatrixDouble & getNormalsAtGaussPts()
bool addField
if true add field values, do not overwrite
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
Calculate tangent vector on edge from HO geometry approximation.
boost::shared_ptr< MatrixDouble > tangentsAtPts
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
OpGetHOTangentsOnEdge(std::string field_name, boost::shared_ptr< MatrixDouble > tangents_at_pts=nullptr)
Constructor for HO tangents calculation on edges.
Transform Hcurl base functions from reference element to physical edge in 3D.
boost::shared_ptr< MatrixDouble > tangentAtGaussPts
OpHOSetContravariantPiolaTransformOnEdge3D(const FieldSpace space=HCURL, boost::shared_ptr< MatrixDouble > tangent_at_pts=nullptr)
Constructor for contravariant Piola transformation on 3D edges.
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.
OpHOSetContravariantPiolaTransformOnFace3D(const FieldSpace space, boost::shared_ptr< MatrixDouble > normals_at_gauss_pts=nullptr)
Constructor for contravariant Piola transformation on 3D faces.
Transform Hcurl base functions from reference element to physical face in 3D.
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)
Constructor for covariant Piola transformation on 3D faces.
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
Scale base functions by inverses of measure of element.
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.
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
OpSetHOContravariantPiolaTransform(const FieldSpace space, boost::shared_ptr< VectorDouble > det_ptr, boost::shared_ptr< MatrixDouble > jac_ptr)
Constructor for contravariant Piola transformation.
Apply covariant (Piola) transformation to Hcurl space for HO geometry.
OpSetHOCovariantPiolaTransform(const FieldSpace space, boost::shared_ptr< MatrixDouble > jac_inv_ptr)
Constructor for covariant Piola transformation.
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.
OpSetHOInvJacVectorBase(const FieldSpace space, boost::shared_ptr< MatrixDouble > inv_jac_ptr)
Constructor for HO inverse Jacobian transformation for vector bases.
Modify integration weights on edge 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.
OpSetHOWeightsOnEdge()
Default constructor for HO weights on edge.
Modify integration weights on face to take into account higher-order geometry.
OpSetHOWeightsOnFace()
Default constructor for HO weights on face.
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.
OpSetHOWeights(boost::shared_ptr< VectorDouble > det_ptr)
Constructor for HO weights calculation.
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
Transform local reference derivatives of shape functions to global derivatives.
OpSetInvJacToScalarBasesBasic(FieldSpace space, boost::shared_ptr< MatrixDouble > inv_jac_ptr)