v0.15.5
Loading...
Searching...
No Matches
ForcesAndSourcesCore.hpp
Go to the documentation of this file.
1/** \file ForcesAndSourcesCore.hpp
2 * \brief Generic finite element implementation using operator pipelines
3 *
4 * This file provides the core framework for finite element computations in MoFEM
5 * based on a pipeline of user data operators. Instead of monolithic finite element
6 * implementations, ForcesAndSourcesCore treats finite elements as sequences of
7 * modular operators that can be combined and reused.
8 *
9 * ## Core Concept: Operator Pipeline Architecture
10 *
11 * **Finite Element = Pipeline of Operators**
12 * Rather than traditional finite element classes with hardcoded functionality,
13 * MoFEM implements finite elements as pipelines of small, focused operators.
14 * Each operator performs a specific task (calculate gradients, apply boundary
15 * conditions, assemble matrices, etc.) and can be combined with others to
16 * create complex finite element behaviors.
17 *
18 * **ForcesAndSourcesCore**: Generic base class that manages the operator pipeline
19 * and provides the framework for element assembly. Users inherit from this class
20 * and add specific operators to implement their particular problem.
21 *
22 * ## Operator Types and Usage:
23 *
24 * **OPROW**: Operators that work with test functions (rows of system matrix)
25 * - Used for: Loading vectors, boundary conditions, source terms
26 * - Example: Applying Neumann boundary conditions, body forces
27 * - Operates during: Right-hand side vector assembly
28 *
29 * **OPCOL**: Operators that work with trial functions (columns of system matrix)
30 * - Used for: Geometric operations, coordinate transformations
31 * - Example: Calculating coordinates, setting up geometry
32 * - Operates during: Matrix structure setup and geometric computations
33 *
34 * **OPROWCOL**: Operators that work with both test and trial functions
35 * - Used for: Stiffness matrices, mass matrices, coupling terms
36 * - Example: Linear elasticity stiffness, thermal conductivity
37 * - Operates during: System matrix assembly (most common type)
38 *
39 * **OPSPACE**: Operators that work on the finite element space itself
40 * - Used for: Basis function transformations, coordinate mappings
41 * - Example: Piola transforms for Hdiv/Hcurl spaces, Jacobian operations
42 * - Operates during: Pre-processing before matrix/vector assembly
43 *
44 *
45 */
46
47
48
49#ifndef __FORCES_AND_SOURCES_CORE__HPP__
50#define __FORCES_AND_SOURCES_CORE__HPP__
51
52using namespace boost::numeric;
53
54namespace MoFEM {
55
56/** \brief structure to get information from mofem into EntitiesFieldData
57 * \ingroup mofem_forces_and_sources
58 *
59 */
61
63
65 typedef boost::function<int(int order_row, int order_col, int order_data)>
67
68 typedef boost::function<MoFEMErrorCode(ForcesAndSourcesCore *fe_raw_ptr,
69 int order_row, int order_col,
70 int order_data)>
72
73 /**
74 * \brief Hook to get rule
75 *
76 * \todo check preferred format how works with gcc and clang,
77 * see
78 * <http://www.boost.org/doc/libs/1_64_0/doc/html/function/tutorial.html#idp247873024>
79 */
81
82 /**
83 * @brief Set function to calculate integration rule
84 *
85 */
87
88 /** \brief Data operator to do calculations at integration points.
89 * \ingroup mofem_forces_and_sources
90
91 Is inherited and implemented by user to do calculations. It can be used in
92 many different ways but typically is used to integrate matrices (f.e.
93 stiffness matrix) and the right hand vector (f.e. force vector).
94
95 Note: It is assumed that operator is executed for symmetric problem. That
96 means that is executed for not repeating entities on finite element. For
97 example on triangle we have nodes, 3 edges and one face. Because of symmetry
98 calculations are for: nodes-nodes, nodes-edge_0, nodes-edge_1, nodes-edge_2,
99 nodes-face,
100 edge_0-edge_0, edge_0-edge_1, edge_0-edge_2,
101 edge_1-edge_1, edge_1-edge_2,
102 edge_2-edge_2,
103 edge_0-face, edge_1-face, edge_2-face,
104 face - face
105
106 In case of non symmetric problem in addition calculations of lower off
107 diagonal terms. F.e. edge_1-edge_0, edge_2-edge_0, edge_2-edge_1,
108
109 In that case class variable UserDataOperator::sYmm = false;
110
111 Note: By default sYmm is set for symmetric problems
112
113 */
114 struct UserDataOperator;
115
116 /** \brief Use to push back operator for row operator
117
118 It can be used to calculate nodal forces or other quantities on the mesh.
119
120 */
121 boost::ptr_deque<UserDataOperator> &getOpPtrVector() { return opPtrVector; }
122
123 /**
124 * @brief Get the Entity Polynomial Base object
125 *
126 * @return boost::shared_ptr<BaseFunction>&&
127 */
129
130 /**
131 * @brief Get the User Polynomial Base object
132 *
133 * @return boost::shared_ptr<BaseFunction>&
134 */
136
137 /**
138 * @brief Matrix of integration points
139 *
140 * Columns is equal to number of integration points, numver of rows
141 * depends on dimension of finite element entity, for example for
142 * tetrahedron rows are x,y,z,weight. Last row is integration weight.
143 *
144 * FIXME: that should be moved to private class data and acessed only by
145 * member function
146 */
148
149 virtual MoFEMErrorCode preProcess();
150 virtual MoFEMErrorCode operator()();
151 virtual MoFEMErrorCode postProcess();
152
153public:
154 /** \brief Get max order of approximation for data fields
155
156 getMaxDataOrder() returns maximal order on entities, for
157 all data on the element. So for example if finite element is triangle, and
158 triangle base function have order 4 and on edges base function have order
159 2, this function returns 4.
160
161 If finite element has for example 2 or more approximated fields, for
162 example Pressure (order 3) and displacement field (order 5), this function
163 returns 5.
164
165 */
166 int getMaxDataOrder() const;
167
168 /// \brief Get max order of approximation for field in rows
169 int getMaxRowOrder() const;
170
171 /// \brief Get max order of approximation for field in columns
172 int getMaxColOrder() const;
173
174 /**
175 * @brief Get the entity data
176 *
177 * @param space
178 * @param type
179 * @param side
180 * @return EntitiesFieldData::EntData&
181 */
182 auto &getEntData(const FieldSpace space, const EntityType type,
183 const int side) {
184 return dataOnElement[space]->dataOnEntities[type][side];
185 }
186
187 /**
188 * @brief Get data on entities and space
189 *
190 * Entities data are stored by space, by entity type, and entity side.
191 *
192 * @return std::array<boost::shared_ptr<EntitiesFieldData>, LASTSPACE>
193 */
195
196 /**
197 * @brief Get derived data on entities and space
198 *
199 * Entities data are stored by space, by entity type, and entity side. Derived
200 * data is used to store data on columns, so it shares information about shape
201 * functions wih rows.
202 *
203 * @return std::array<boost::shared_ptr<EntitiesFieldData>, LASTSPACE>
204 */
206
207protected:
208 /**
209 * \brief get sense (orientation) of entity
210 * @param type type of entity
211 * @param data entity data
212 * @return error code
213 */
215 const EntityType type,
216 boost::ptr_vector<EntitiesFieldData::EntData> &data) const;
217
218 /**
219 * @brief Get the entity data order
220 *
221 * @param type
222 * @param space
223 * @param data
224 * @return MoFEMErrorCode
225 */
227 const EntityType type, const FieldSpace space,
228 boost::ptr_vector<EntitiesFieldData::EntData> &data) const;
229
230 /**
231 * @brief Get the entity sense (orientation)
232 *
233 * @tparam type
234 * @param data
235 * @return MoFEMErrorCode
236 */
237 template <EntityType type>
239 return getEntitySense(type, data.dataOnEntities[type]);
240 }
241
242 /**
243 * @brief Get the entity data order for given space
244 *
245 * @tparam type
246 * @param data
247 * @param space
248 * @return MoFEMErrorCode
249 */
250 template <EntityType type>
252 const FieldSpace space) const {
253 return getEntityDataOrder(type, space, data.dataOnEntities[type]);
254 }
255
256 /** \name Indices */
257
258 /**@{*/
259
260 /// \brief get node indices
261 template <typename EXTRACTOR>
263 getNodesIndices(const int bit_number,
264 FieldEntity_vector_view &ents_field, VectorInt &nodes_indices,
265 VectorInt &local_nodes_indices, EXTRACTOR &&extractor) const;
266
267 /// \brief get row node indices from FENumeredDofEntity_multiIndex
269 const int bit_number) const;
270
271 /// \brief get col node indices from FENumeredDofEntity_multiIndex
273 const int bit_number) const;
274
275 template <typename EXTRACTOR>
276 MoFEMErrorCode getEntityIndices(EntitiesFieldData &data, const int bit_number,
277 FieldEntity_vector_view &ents_field,
278 const EntityType type_lo,
279 const EntityType type_hi,
280 EXTRACTOR &&extractor) const;
281
283 getEntityRowIndices(EntitiesFieldData &data, const int bit_number,
284 const EntityType type_lo = MBVERTEX,
285 const EntityType type_hi = MBPOLYHEDRON) const;
286
288 getEntityColIndices(EntitiesFieldData &data, const int bit_number,
289 const EntityType type_lo = MBVERTEX,
290 const EntityType type_hi = MBPOLYHEDRON) const;
291
292
293 /// \brief get col NoField indices
295 const int bit_number) const;
296
297 /// \brief get col NoField indices
299 const int bit_number) const;
300
301 /**@}*/
302
303 /** \name Data */
304
305 /**@{*/
306
307 /** Get bit ref level in entities, and set it to data
308 */
310
311 /**
312 * \brief Get data on nodes
313 * @param data Data structure
314 * @param field_name Field name
315 * @return Error code
316 */
318 const int bit_number) const;
319
321 getEntityFieldData(EntitiesFieldData &data, const int bit_number,
322 const EntityType type_lo = MBVERTEX,
323 const EntityType type_hi = MBPOLYHEDRON) const;
324
325 /**
326 * @brief Get field data on entities where no field is defined
327 *
328 * This is not standards, since op_data, dataOnElement[NOFIELD], has to be set
329 * for rows and columns for each NOFIELD, also size of sides is not determined.
330 *
331 * @param data @param bit_number @return MoFEMErrorCode
332 */
333 template <typename EXTRACTOR>
335 const int bit_number,
336 EXTRACTOR &&extractor) const;
337
338 /**@}*/
339
340 /// \brief Get nodes on faces
342
343 /// \brief Get field approximation space and base on entities
346
347 /** \name Data form NumeredDofEntity_multiIndex */
348
349 /**@{*/
350
351 /// \brief get indices of nodal indices which are declared for problem but
352 /// not this particular element
354 const NumeredDofEntity_multiIndex &dofs,
355 VectorInt &nodes_indices) const;
356
357 /// \brief get indices by type (generic function) which are declared for
358 /// problem but not this particular element
360 const NumeredDofEntity_multiIndex &dofs,
361 EntityType type, int side_number,
362 VectorInt &indices) const;
363
365 VectorInt &nodes_indices) const;
367 EntityType type, int side_number,
368 VectorInt &indices) const;
370 VectorInt &nodes_indices) const;
372 EntityType type, int side_number,
373 VectorInt &indices) const;
374
375 /**@}*/
376
377 /**
378 * \brief another variant of getRule
379 * @param order_row order of base function on row
380 * @param order_col order of base function on columns
381 * @param order_data order of base function approximating data
382 * @return integration rule
383 *
384 * This function is overloaded by the user. The integration rule
385 * is set such that specific operator implemented by the user is
386 * integrated accurately. For example if user implement bilinear operator
387 * \f[
388 * b(u,v) =
389 * \int_\mathcal{T}
390 * \frac{\partial u_i}{\partial x_j}\frac{\partial v_i}{\partial x_j}
391 * \textrm{d}\mathcal{T}
392 * \f]
393 * then if \f$u\f$ and \f$v\f$ are polynomial of given \em order, then
394 * exact integral would be \code int getRule(int order) { return
395 * 2*(order-1); }; \endcode
396 *
397 * The integration points and weights are set appropriately for given
398 * entity type and integration rule from \ref quad.c
399 *
400 * Method \ref ForcesAndSourcesCore::getRule takes at argument takes
401 * maximal polynomial order set on the element on all fields defined on
402 * the element. If a user likes to have more control, another variant of
403 * this function can be called which distinguishing between field orders
404 * on rows, columns and data, the i.e. first argument of a bilinear form,
405 * the second argument of bilinear form and field coefficients on the
406 * element.
407 *
408 * \note If user set rule to -1 or any other negative integer, then method
409 * \ref ForcesAndSourcesCore::setGaussPts is called. In that method user
410 * can implement own (specific) integration method.
411 *
412 * \bug this function should be const
413 */
414 virtual int getRule(int order_row, int order_col, int order_data);
415
416 /** \brief set user specific integration rule
417
418 This function allows for user defined integration rule. The key is to
419 called matrix gaussPts, which is used by other MoFEM procedures. Matrix
420 has number of rows equal to problem dimension plus one, where last index
421 is used to store weight values. %Number of columns is equal to number of
422 integration points.
423
424 \note This function is called if method \ref
425 ForcesAndSourcesCore::getRule is returning integer -1 or any other
426 negative integer.
427
428 User sets
429 \code
430 MatrixDouble gaussPts;
431 \endcode
432 where
433 \code
434 gaussPts.resize(dim+1,nb_gauss_pts);
435 \endcode
436 number rows represents local coordinates of integration points
437 in reference element, where last index in row is for integration weight.
438
439 */
440 virtual MoFEMErrorCode setGaussPts(int order_row, int order_col,
441 int order_data);
442
443 /**
444 * \brief Calculate base functions
445 * @return Error code
446 */
449
450 /**
451 * \brief Calculate base functions
452 * @return Error code
453 */
455
456 /**
457 * @brief Calculate Bernstein-Bezier base
458 *
459 * @return MoFEMErrorCode
460 */
462
463 /**
464 * @brief Create a entity data on element object
465 *
466 * @return MoFEMErrorCode
467 */
468 MoFEMErrorCode createDataOnElement(EntityType type);
469
470 /**
471 * @brief Iterate user data operators
472 *
473 * @return MoFEMErrorCode
474 */
476
477 /**@{*/
478
479 /** \name Deprecated (do not use) */
480
481 /** \deprecated Use getRule(int row_order, int col_order, int data order)
482 */
483 virtual int getRule(int order);
484
485 /** \deprecated setGaussPts(int row_order, int col_order, int data order);
486 */
487 virtual MoFEMErrorCode setGaussPts(int order);
488
489 /**@}*/
490
491 /**
492 * @brief Entity data on element entity rows fields
493 *
494 */
495 const std::array<boost::shared_ptr<EntitiesFieldData>, LASTSPACE>
497
498 /**
499 * @brief Entity data on element entity columns fields
500 *
501 */
502 const std::array<boost::shared_ptr<EntitiesFieldData>, LASTSPACE>
504
510
511 /**
512 * @brief Vector of finite element users data operators
513 *
514 */
515 boost::ptr_deque<UserDataOperator> opPtrVector;
516
517 friend class UserDataOperator;
518
519protected:
520 /**
521 * @brief Last evaluated type of element entity
522 *
523 */
525
526private:
527 /**
528 * @brief Pointer to entity polynomial base
529 *
530 */
531 boost::shared_ptr<BaseFunction> elementPolynomialBasePtr;
532
533 /**
534 * @brief Pointer to user polynomial base
535 */
536 boost::shared_ptr<BaseFunction> userPolynomialBasePtr;
537
538 /**
539 * @brief Element to integrate on the sides
540 *
541 */
543
544 /**
545 * @brief Set the pointer to face element on the side
546 *
547 * \note Function is is used by face element, while it iterates over
548 * elements on the side
549 *
550 * @param side_fe_ptr
551 * @return MoFEMErrorCode
552 */
554
555 /**u
556 * @brief Element to integrate parent or child
557 *
558 */
560
561 /**
562 * @brief Set the pointer to face element refined
563 *
564 * @param refine_fe_ptr
565 * @return MoFEMErrorCode
566 */
568
574
575 template <int DIM> friend struct OpCopyGeomDataToE;
576 template <typename E> friend struct OpBrokenLoopSide;
577
578protected:
579 MatrixDouble coordsAtGaussPts; ///< coordinated at gauss points
580 double elementMeasure; ///< Depending on dimension of elements, stores length,
581 ///< area or volume of element.
582};
583
585
586 /**
587 * \brief Controls loop over entities on element
588 *
589 * - OPROW is used if row vector is assembled
590 * - OPCOL is usually used if column vector is assembled
591 * - OPROWCOL is usually used for assemble matrices.
592 * - OPSPACE no field is defined for such operator. Is usually used to modify
593 * base
594 *
595 *
596 * For typical problem like Bubnov-Galerkin OPROW and OPCOL are the same. In
597 * more general case for example for non-square matrices columns and rows
598 * could have different numeration and/or different set of fields.
599 *
600 */
601 enum OpType {
602 OPROW = 1 << 0, ///< operator doWork function is executed on FE rows
603 OPCOL = 1 << 1, ///< operator doWork function is executed on FE columns
604 OPROWCOL = 1 << 2, ///< operator doWork is executed on FE rows &columns
605 OPSPACE = 1 << 3, ///< operator do Work is execute on space data
606 OPLAST = 1 << 3 ///< @deprecated would be removed
607 };
608
609 static const char *const OpTypeNames[];
610
611 char opType;
612 std::string rowFieldName;
613 std::string colFieldName;
615
616 /**
617 * @brief Constructor for operators working on finite element spaces
618 *
619 * This constructor is used when operators modify basis functions on
620 * approximation spaces. The operator runs for all data on the specified
621 * space without access to specific field data. Commonly used for:
622 * - Coordinate transformations (Piola transforms)
623 * - Integration weight modifications
624 * - Basis function scaling or rotation
625 *
626 * @param space Field space to operate on (H1, HDIV, HCURL, L2, etc.)
627 * @param type Operator type (OPSPACE, OPROW, OPCOL, OPROWCOL)
628 * @param symm Whether the operator maintains matrix symmetry (default: true)
629 */
630 UserDataOperator(const FieldSpace space, const char type = OPSPACE,
631 const bool symm = true);
632
633 /**
634 * @brief Constructor for operators working on a single field
635 *
636 * This constructor creates operators that work with data from a specific
637 * field. Used for operations like:
638 * - Calculating field values at integration points
639 * - Computing field gradients or derivatives
640 * - Applying boundary conditions to specific fields
641 * - Post-processing field data
642 *
643 * @param field_name Name of the field to operate on
644 * @param type Operator type:
645 * - OPROW: Works with test functions (RHS assembly)
646 * - OPCOL: Works with trial functions (geometry setup)
647 * - OPROWCOL: Works with both (matrix assembly)
648 * @param symm Whether the operator maintains matrix symmetry (default: true)
649 */
650 UserDataOperator(const std::string field_name, const char type,
651 const bool symm = true);
652
653 /**
654 * @brief Constructor for operators working on two fields (bilinear forms)
655 *
656 * This constructor creates operators that couple two different fields,
657 * essential for:
658 * - Mixed formulations (pressure-velocity coupling)
659 * - Multi-physics problems (thermal-mechanical coupling)
660 * - Contact mechanics (displacement-Lagrange multiplier)
661 * - Interface problems between different field types
662 *
663 * @param row_field_name Name of the field for test functions (rows)
664 * @param col_field_name Name of the field for trial functions (columns)
665 * @param type Operator type (typically OPROWCOL for coupling)
666 * @param symm Whether the operator maintains matrix symmetry (default: true)
667 * Note: Coupling between different fields is often non-symmetric
668 */
669 UserDataOperator(const std::string row_field_name,
670 const std::string col_field_name, const char type,
671 const bool symm = true);
672
673 /** \brief Return raw pointer to NumeredEntFiniteElement
674 */
675 inline boost::shared_ptr<const NumeredEntFiniteElement>
677
678 /**
679 * \brief Return finite element entity handle
680 * @return Finite element entity handle
681 */
682 inline EntityHandle getFEEntityHandle() const;
683
684 /**
685 * @brief Get dimension of finite element
686 *
687 * @return int
688 */
689 inline int getFEDim() const;
690
691 /**
692 * @brief Get dimension of finite element
693 *
694 * @return int
695 */
696 inline EntityType getFEType() const;
697
698 /**
699 * @brief Get the side number pointer
700 *
701 * \note For vertex is expectation. Side basses in argument of function doWork
702 * is zero. For other entity types side can be used as argument of this
703 * function.
704 *
705 * @param side_number
706 * @param type
707 * @return boost::weak_ptr<SideNumber>
708 */
709 inline boost::weak_ptr<SideNumber> getSideNumberPtr(const int side_number,
710 const EntityType type);
711
712 /**
713 * @brief Get the side entity
714 *
715 * \note For vertex is expectation. Side basses in argument of function
716 * doWork is zero. For other entity types side can be used as argument of
717 * this function.
718 *
719 * \code
720 * MoFEMErrorCode doWork(int side, EntityType type,
721 * EntitiesFieldData::EntData &data) {
722 * MoFEMFunctionBegin;
723 *
724 * if (type == MBVERTEX) {
725 * for (int n = 0; n != number_of_nodes; ++n)
726 * EntityHandle ent = getSideEntity(n, type);
727 *
728 * // Do somthing
729 *
730 * } else {
731 * EntityHandle ent = getSideEntity(side, type);
732 *
733 * // Do somthing
734 *
735 * }
736 *
737 * MoFEMFunctionReturn(0);
738 * }
739 * \endcode
740 *
741 * @param side_number
742 * @param type
743 */
744 inline EntityHandle getSideEntity(const int side_number,
745 const EntityType type);
746
747 /**
748 * @brief Get the number of nodes on finite element
749 *
750 * @return int
751 */
752 inline int getNumberOfNodesOnElement() const;
753
754 /** \brief Get row indices
755
756 Field could be or not declared for this element but is declared for
757 problem
758
759 \param field_name
760 \param type entity type
761 \param side side number, any number if type is MBVERTEX
762 \return indices
763
764 NOTE: Using those indices to assemble matrix will result in error if new
765 non-zero values need to be created.
766
767 */
768 MoFEMErrorCode getProblemRowIndices(const std::string filed_name,
769 const EntityType type, const int side,
770 VectorInt &indices) const;
771
772 /** \brief Get col indices
773
774 Field could be or not declared for this element but is declared for
775 problem
776
777 \param field_name
778 \param type entity type
779 \param side side number, any number if type is MBVERTEX
780 \return indices
781
782 NOTE: Using those indices to assemble matrix will result in error if new
783 non-zero values need to be created.
784
785 */
786 MoFEMErrorCode getProblemColIndices(const std::string filed_name,
787 const EntityType type, const int side,
788 VectorInt &indices) const;
789
790 /** \brief Return raw pointer to Finite Element Method object
791 */
792 inline const FEMethod *getFEMethod() const;
793
794 /**
795 * \brief Get operator types
796 * @return Return operator type
797 */
798 inline int getOpType() const;
799
800 /**
801 * \brief Set operator type
802 * @param Operator type
803 */
804 inline void setOpType(const OpType type);
805
806 /**
807 * \brief Add operator type
808 */
809 inline void addOpType(const OpType type);
810
811 /**
812 * \brief get number of finite element in the loop
813 * @return number of finite element
814 */
815 inline int getNinTheLoop() const;
816
817 /**
818 * \brief get size of elements in the loop
819 * @return loop size
820 */
821 inline int getLoopSize() const;
822
823 /** \brief Get name of the element
824 */
825 inline std::string getFEName() const;
826
827 /** \name Accessing KSP */
828
829 /**@{*/
830
831 inline const PetscData::Switches &getDataCtx() const;
832
833 inline const KspMethod::KSPContext getKSPCtx() const;
834
835 inline const SnesMethod::SNESContext getSNESCtx() const;
836
837 inline const TSMethod::TSContext getTSCtx() const;
838
839 /**@}*/
840
841 /**@{*/
842
843 inline Vec getKSPf() const;
844
845 inline Mat getKSPA() const;
846
847 inline Mat getKSPB() const;
848
849 /**@}*/
850
851 /** \name Accessing SNES */
852
853 /**@{*/
854
855 inline Vec getSNESf() const;
856
857 inline Vec getSNESx() const;
858
859 inline Mat getSNESA() const;
860
861 inline Mat getSNESB() const;
862
863 /**@}*/
864
865 /** \name Accessing TS */
866
867 /**@{*/
868
869 inline Vec getTSu() const;
870
871 inline Vec getTSu_t() const;
872
873 inline Vec getTSu_tt() const;
874
875 inline Vec getTSf() const;
876
877 inline Mat getTSA() const;
878
879 inline Mat getTSB() const;
880
881 inline int getTSstep() const;
882
883 inline double getTStime() const;
884
885 inline double getTStimeStep() const;
886
887 inline double getTSa() const;
888
889 inline double getTSaa() const;
890
891 /**@}*/
892
893 /**@{*/
894
895 /**@}*/
896
897 /** \name Base functions and integration points */
898
899 /**@{*/
900
901 /** \brief matrix of integration (Gauss) points for Volume Element
902 *
903 * For triangle: columns 0,1 are x,y coordinates respectively and column
904 * 2 is a weight value for example getGaussPts()(1,13) returns y
905 * coordinate of 13th Gauss point on particular volume element
906 *
907 * For tetrahedron: columns 0,1,2 are x,y,z coordinates respectively and
908 * column 3 is a weight value for example getGaussPts()(1,13) returns y
909 * coordinate of 13th Gauss point on particular volume element
910 *
911 */
912 inline MatrixDouble &getGaussPts();
913
914 /**
915 * @brief Get integration weights
916 *
917 * \code
918 * auto t_w = getFTensor0IntegrationWeight();
919 * for(int gg = 0; gg!=getGaussPts.size2(); ++gg) {
920 * // integrate something
921 * ++t_w;
922 * }
923 * \endcode
924 *
925 * @return FTensor::Tensor0<FTensor::PackPtr<double *, 1>>
926 */
927 inline auto getFTensor0IntegrationWeight();
928
929 /**@}*/
930
931 /** \name Coordinates and access to internal data */
932
933 /**@{*/
934
935 /** \brief Gauss points and weight, matrix (nb. of points x 3)
936
937 Column 0-2 integration points coordinate x, y and z, respectively. At rows
938 are integration points.
939
940 */
942
943 /**
944 * \brief Get coordinates at integration points assuming linear geometry
945 *
946 * \code
947 * auto t_coords = getFTensor1CoordsAtGaussPts();
948 * for(int gg = 0;gg!=nb_int_ptrs;gg++) {
949 * // do something
950 * ++t_coords;
951 * }
952 * \endcode
953 *
954 */
955 inline auto getFTensor1CoordsAtGaussPts();
956
957 /**@}*/
958
959 /**@{*/
960
961 /** \name Measures (area, volume, length, etc.) */
962
963 /**
964 * \brief get measure of element
965 * @return volume
966 */
967 inline double getMeasure() const;
968
969 /**
970 * \brief get measure of element
971 * @return volume
972 */
973 inline double &getMeasure();
974
975 /**}*/
976
977 /**@{*/
978
979 /** \name Loops */
980
981 using AdjCache =
982 std::map<EntityHandle,
983 std::vector<boost::weak_ptr<NumeredEntFiniteElement>>>;
984
985 /**
986 * @brief User calls this function to loop over elements on the side of
987 * face. This function calls finite element with its operator to do
988 * calculations.
989 *
990 * @param fe_name name of the side element
991 * @param side_fe pointer to the side element instance
992 * @param dim dimension the of side element
993 * @param ent_for_side entity handle for which adjacent volume or face will
994 * be accessed
995 * @param verb
996 * @param sev
997 * @param adj_cache
998 * @return MoFEMErrorCode
999 */
1000 MoFEMErrorCode loopSide(const string &fe_name, ForcesAndSourcesCore *side_fe,
1001 const size_t dim, const EntityHandle ent_for_side = 0,
1002 boost::shared_ptr<Range> fe_range = nullptr,
1003 const int verb = QUIET,
1004 const LogManager::SeverityLevel sev = Sev::noisy,
1005 AdjCache *adj_cache = nullptr
1006
1007 );
1008
1009 /**
1010 * @brief User calls this function to loop over the same element using a
1011 * different set of integration points. This function calls finite element
1012 * with its operator to do calculations.
1013 *
1014 * @param fe_name
1015 * @param this_fe
1016 * @param verb
1017 * @param sev
1018 * @return MoFEMErrorCode
1019 */
1020 MoFEMErrorCode loopThis(const string &fe_name,
1021 ForcesAndSourcesCore *this_fe,
1022 const int verb = QUIET,
1023 const LogManager::SeverityLevel sev = Sev::noisy);
1024
1025 /**
1026 * @brief User calls this function to loop over parent elements. This function
1027 * calls finite element with its operator to do calculations.
1028 *
1029 * @param fe_name
1030 * @param parent_fe
1031 * @param verb
1032 * @param sev
1033 * @return MoFEMErrorCode
1034 */
1035 MoFEMErrorCode loopParent(const string &fe_name,
1036 ForcesAndSourcesCore *parent_fe,
1037 const int verb = QUIET,
1038 const LogManager::SeverityLevel sev = Sev::noisy);
1039
1040 /**
1041 * @brief User calls this function to loop over parent elements. This function
1042 * calls finite element with its operator to do calculations.
1043 *
1044 * @param fe_name
1045 * @param child_fe
1046 * @param verb
1047 * @param sev
1048 * @return MoFEMErrorCode
1049 */
1050 MoFEMErrorCode loopChildren(const string &fe_name,
1051 ForcesAndSourcesCore *child_fe,
1052 const int verb = QUIET,
1053 const LogManager::SeverityLevel sev = Sev::noisy);
1054
1055 /**
1056 * @brief Iterate over range of elements
1057 *
1058 * @param fe_name
1059 * @param range_fe
1060 * @param fe_range
1061 * @param verb
1062 * @param sev
1063 * @return * MoFEMErrorCode
1064 */
1065 MoFEMErrorCode loopRange(const string &fe_name,
1066 ForcesAndSourcesCore *range_fe,
1067 boost::shared_ptr<Range> fe_range,
1068 const int verb = QUIET,
1069 const LogManager::SeverityLevel sev = Sev::noisy);
1070
1071 /**@}*/
1072
1073 inline ForcesAndSourcesCore *getPtrFE() const;
1074
1075 inline ForcesAndSourcesCore *getSidePtrFE() const;
1076
1077 inline ForcesAndSourcesCore *getRefinePtrFE() const;
1078
1079 /**@{*/
1080
1081 /** \name Database */
1082
1084 inline moab::Interface &getMoab() { return getMField().get_moab(); }
1085
1086 /**@}*/
1087
1088 /**@{ */
1089
1090 /** \name Pipelines */
1091
1092 virtual boost::weak_ptr<ForcesAndSourcesCore> getSubPipelinePtr() const {
1093 return boost::weak_ptr<ForcesAndSourcesCore>();
1094 };
1095
1096 /**@} */
1097
1098protected:
1100
1102
1103private:
1108};
1109
1110/// \deprecated Used ForcesAndSourcesCore instead
1112
1113boost::shared_ptr<const NumeredEntFiniteElement>
1117
1119 return getNumeredEntFiniteElementPtr()->getEnt();
1120}
1121
1125
1129
1130boost::weak_ptr<SideNumber>
1132 const int side_number, const EntityType type) {
1133 auto &side_table_by_side_and_type =
1134 ptrFE->numeredEntFiniteElementPtr->getSideNumberTable().get<1>();
1135 auto side_it =
1136 side_table_by_side_and_type.find(boost::make_tuple(type, side_number));
1137 if (side_it != side_table_by_side_and_type.end())
1138 return *side_it;
1139 else
1140 return boost::weak_ptr<SideNumber>();
1141}
1142
1145 const EntityType type) {
1146 if (auto side_ptr = getSideNumberPtr(side_number, type).lock())
1147 return side_ptr->ent;
1148 else
1149 return 0;
1150}
1151
1153 return ptrFE->getNumberOfNodes();
1154}
1155
1157 return ptrFE;
1158}
1159
1161
1163 opType = type;
1164}
1165
1167 opType |= type;
1168}
1169
1171 return getFEMethod()->getNinTheLoop();
1172}
1173
1175 return getFEMethod()->getLoopSize();
1176}
1177
1179 return getFEMethod()->getFEName();
1180}
1181
1182const PetscData::Switches &
1184 return getFEMethod()->data_ctx;
1185}
1186
1189 return getFEMethod()->ksp_ctx;
1190}
1191
1194 return getFEMethod()->snes_ctx;
1195}
1196
1199 return getFEMethod()->ts_ctx;
1200}
1201
1203#ifndef NDEBUG
1204 if (getFEMethod()->ksp_f == PETSC_NULLPTR)
1205 THROW_MESSAGE("KSP not set F vector");
1206#endif
1207 return getFEMethod()->ksp_f;
1208}
1209
1211#ifndef NDEBUG
1212 if (getFEMethod()->ksp_A == PETSC_NULLPTR)
1213 THROW_MESSAGE("KSP not set A vector");
1214#endif
1215 return getFEMethod()->ksp_A;
1216}
1217
1219#ifndef NDEBUG
1220 if (getFEMethod()->ksp_B == PETSC_NULLPTR)
1221 THROW_MESSAGE("KSP not set B vector");
1222#endif
1223 return getFEMethod()->ksp_B;
1224}
1225
1227#ifndef NDEBUG
1228 if (getFEMethod()->snes_f == PETSC_NULLPTR)
1229 THROW_MESSAGE("SNES not set F vector");
1230#endif
1231 return getFEMethod()->snes_f;
1232}
1233
1235#ifndef NDEBUG
1236 if (getFEMethod()->snes_x == PETSC_NULLPTR)
1237 THROW_MESSAGE("SNESnot set X vector");
1238#endif
1239 return getFEMethod()->snes_x;
1240}
1241
1243#ifndef NDEBUG
1244 if (getFEMethod()->snes_A == PETSC_NULLPTR)
1245 THROW_MESSAGE("SNES not set A vector");
1246#endif
1247 return getFEMethod()->snes_A;
1248}
1249
1251#ifndef NDEBUG
1252 if (getFEMethod()->snes_B == PETSC_NULLPTR)
1253 THROW_MESSAGE("SNES not set A matrix");
1254#endif
1255 return getFEMethod()->snes_B;
1256}
1257
1259#ifndef NDEBUG
1260 if (getFEMethod()->ts_u == PETSC_NULLPTR)
1261 THROW_MESSAGE("TS not set U vector");
1262#endif
1263 return getFEMethod()->ts_u;
1264}
1265
1267#ifndef NDEBUG
1268 if (getFEMethod()->ts_u_t == PETSC_NULLPTR)
1269 THROW_MESSAGE("TS not set U_t vector");
1270#endif
1271 return getFEMethod()->ts_u_t;
1272}
1273
1275#ifndef NDEBUG
1276 if (getFEMethod()->ts_u_tt == PETSC_NULLPTR)
1277 THROW_MESSAGE("TS not set U_tt vector");
1278#endif
1279 return getFEMethod()->ts_u_tt;
1280}
1281
1283#ifndef NDEBUG
1284 if (getFEMethod()->ts_F == PETSC_NULLPTR)
1285 THROW_MESSAGE("TS not set F vector");
1286#endif
1287 return getFEMethod()->ts_F;
1288}
1289
1291#ifndef NDEBUG
1292 if (getFEMethod()->ts_A == PETSC_NULLPTR)
1293 THROW_MESSAGE("TS not set A matrix");
1294#endif
1295 return getFEMethod()->ts_A;
1296}
1297
1299#ifndef NDEBUG
1300 if (getFEMethod()->ts_B == PETSC_NULLPTR)
1301 THROW_MESSAGE("TS not set B matrix");
1302#endif
1303 return getFEMethod()->ts_B;
1304}
1305
1307#ifndef NDEBUG
1308 if ((getFEMethod()->data_ctx & PetscData::PetscData::CtxSetTime).none())
1309 THROW_MESSAGE("TS not set time");
1310#endif
1311 return getFEMethod()->ts_step;
1312}
1313
1315#ifndef NDEBUG
1316 if ((getFEMethod()->data_ctx & PetscData::PetscData::CtxSetTime).none())
1317 THROW_MESSAGE("TS not set time");
1318#endif
1319 return getFEMethod()->ts_t;
1320}
1321
1323#ifndef NDEBUG
1324 if ((getFEMethod()->data_ctx & PetscData::PetscData::CtxSetTime).none())
1325 THROW_MESSAGE("TS not set time");
1326#endif
1327 return getFEMethod()->ts_dt;
1328}
1329
1331#ifndef NDEBUG
1332 if ((getFEMethod()->data_ctx & (PetscData::CtxSetA | PetscData::CtxSetB))
1333 .none() ||
1334 (getFEMethod()->data_ctx & (PetscData::CtxSetX_T)).none())
1335 THROW_MESSAGE("TS not set B matrix");
1336#endif
1337 return getFEMethod()->ts_a;
1338}
1339
1341#ifndef NDEBUG
1342 if ((getFEMethod()->data_ctx & (PetscData::CtxSetA | PetscData::CtxSetB))
1343 .none() ||
1344 (getFEMethod()->data_ctx & (PetscData::CtxSetX_TT)).none())
1345 THROW_MESSAGE("TS not set B matrix");
1346#endif
1347 return getFEMethod()->ts_aa;
1348}
1349
1353
1356 &(getGaussPts()(getGaussPts().size1() - 1, 0)));
1357}
1358
1362
1367
1372
1376
1379 &getCoordsAtGaussPts()(0, 0), &getCoordsAtGaussPts()(0, 1),
1380 &getCoordsAtGaussPts()(0, 2));
1381}
1382
1384 return static_cast<ForcesAndSourcesCore *>(ptrFE)->elementMeasure;
1385}
1386
1390
1391/**
1392 * @brief Element used to execute operators on side of the element
1393 *
1394 * @tparam E template for side element type
1395 *
1396 */
1397template <typename E>
1399
1401
1402 OpLoopSide(MoFEM::Interface &m_field, const std::string fe_name,
1403 const int side_dim, boost::shared_ptr<Range> fe_range,
1404 const LogManager::SeverityLevel sev = Sev::noisy,
1405 boost::shared_ptr<AdjCache> adj_cache = nullptr)
1406 : UserDataOperator(NOSPACE, OPSPACE), sideFEPtr(new E(m_field)),
1407 sideFEName(fe_name), sideDim(side_dim), sevLevel(sev),
1408 adjCache(adj_cache), feRange(fe_range) {}
1409
1410 /**
1411 * @brief Construct a new Op Loop Side object
1412 *
1413 * @param m_field
1414 * @param fe_name name of side (domain element)
1415 * @param side_dim dimension
1416 * @param sev severity level
1417 * @param adj_cache if set to null cache of the element is used
1418 */
1419 OpLoopSide(MoFEM::Interface &m_field, const std::string fe_name,
1420 const int side_dim,
1421 const LogManager::SeverityLevel sev = Sev::noisy,
1422 boost::shared_ptr<AdjCache> adj_cache = nullptr)
1423 : OpLoopSide(m_field, fe_name, side_dim, nullptr, sev, adj_cache) {}
1424
1425 MoFEMErrorCode doWork(int side, EntityType type,
1427 return loopSide(sideFEName, sideFEPtr.get(), sideDim, 0, feRange, VERBOSE,
1428 sevLevel, adjCache.get());
1429 };
1430
1431 boost::ptr_deque<UserDataOperator> &getOpPtrVector() {
1432 return sideFEPtr->getOpPtrVector();
1433 }
1434
1435 boost::shared_ptr<E> &getSideFEPtr() { return sideFEPtr; }
1436
1437 boost::weak_ptr<ForcesAndSourcesCore> getSubPipelinePtr() const override {
1438 return sideFEPtr;
1439 }
1440
1441protected:
1442 const std::string sideFEName;
1443 const int sideDim;
1444 boost::shared_ptr<E> sideFEPtr;
1446 boost::shared_ptr<AdjCache> adjCache;
1447 boost::shared_ptr<Range> feRange;
1448};
1449
1450/**
1451 * @brief Iterate over range of (sub)elements
1452 *
1453 * @tparam E
1454 */
1455template <typename E>
1457
1459
1460 OpLoopRange(MoFEM::Interface &m_field, const std::string fe_name,
1461 boost::shared_ptr<Range> fe_range,
1462 const LogManager::SeverityLevel sev = Sev::noisy)
1463 : UserDataOperator(NOSPACE, OPSPACE), rangeFEPtr(new E(m_field)),
1464 rangeFEName(fe_name), sevLevel(sev), feRange(fe_range) {}
1465
1466 MoFEMErrorCode doWork(int side, EntityType type,
1468 return loopRange(rangeFEName, rangeFEPtr.get(), feRange, VERBOSE, sevLevel);
1469 };
1470
1471 boost::ptr_deque<UserDataOperator> &getOpPtrVector() {
1472 return rangeFEPtr->getOpPtrVector();
1473 }
1474
1475 boost::shared_ptr<E> &getRangeFEPtr() { return rangeFEPtr; }
1476
1477 boost::weak_ptr<ForcesAndSourcesCore> getSubPipelinePtr() const override {
1478 return rangeFEPtr;
1479 }
1480
1481protected:
1482 const std::string rangeFEName;
1483 boost::shared_ptr<E> rangeFEPtr;
1485 boost::shared_ptr<Range> feRange;
1486};
1487
1488/**
1489 * \brief Copy geometry-related data from one element to other
1490 *
1491 * That can be used to copy high order geometry data from coarse element to
1492 * children. That is often a case when higher order geometry is defined only on
1493 * coarse elements.
1494 *
1495 * \note Integration points have to be located at the same geometric positions
1496 *
1497 * FIXME: Write atom test
1498 */
1499template <int DIM>
1501
1502} // namespace MoFEM
1503
1504#endif //__FORCES_AND_SOURCES_CORE__HPP__
1505
1506/**
1507 * \defgroup mofem_forces_and_sources Forces and sources
1508 * \ingroup mofem
1509 *
1510 * \brief Manages complexities related to assembly of vector and matrices at
1511 * single finite element level.
1512 *
1513 **/
@ QUIET
@ VERBOSE
FieldApproximationBase
approximation base
Definition definitions.h:58
FieldSpace
approximation spaces
Definition definitions.h:82
@ LASTSPACE
FieldSpace in [ 0, LASTSPACE )
Definition definitions.h:89
@ NOSPACE
Definition definitions.h:83
#define DEPRECATED
Definition definitions.h:17
#define THROW_MESSAGE(msg)
Throw MoFEM exception.
constexpr int order
multi_index_container< boost::shared_ptr< NumeredDofEntity >, indexed_by< ordered_unique< tag< Unique_mi_tag >, const_mem_fun< NumeredDofEntity::interface_type_DofEntity, UId, &NumeredDofEntity::getLocalUniqueId > >, ordered_non_unique< tag< Part_mi_tag >, member< NumeredDofEntity, unsigned int, &NumeredDofEntity::pArt > >, ordered_non_unique< tag< Idx_mi_tag >, member< NumeredDofEntity, DofIdx, &NumeredDofEntity::dofIdx > >, ordered_non_unique< tag< PetscGlobalIdx_mi_tag >, member< NumeredDofEntity, DofIdx, &NumeredDofEntity::petscGloablDofIdx > >, ordered_non_unique< tag< PetscLocalIdx_mi_tag >, member< NumeredDofEntity, DofIdx, &NumeredDofEntity::petscLocalDofIdx > >, ordered_non_unique< tag< Ent_mi_tag >, const_mem_fun< NumeredDofEntity::interface_type_DofEntity, EntityHandle, &NumeredDofEntity::getEnt > > > > NumeredDofEntity_multiIndex
MultiIndex container keeps NumeredDofEntity.
SeverityLevel
Severity levels.
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
UBlasVector< int > VectorInt
Definition Types.hpp:67
implementation of Data Operators for Forces and Sources
Definition Common.hpp:10
auto type_from_handle(const EntityHandle h)
get type from entity handle
DEPRECATED typedef ForcesAndSourcesCore ForcesAndSurcesCore
auto dimension_from_handle(const EntityHandle h)
get entity dimension form handle
std::vector< boost::weak_ptr< FieldEntity > > FieldEntity_vector_view
constexpr auto field_name
virtual moab::Interface & get_moab()=0
base operator to do operations at Gauss Pt. level
Deprecated interface functions.
Data on single entity (This is passed as argument to DataOperator::doWork)
data structure for finite element entity
std::array< boost::ptr_vector< EntData >, MBMAXTYPE > dataOnEntities
Structure for user loop methods on finite elements.
EntityHandle getFEEntityHandle() const
Get the entity handle of the current finite element.
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
Shared pointer to finite element database structure.
Base face element used to integrate on skeleton.
void setOpType(const OpType type)
Set operator type.
EntityHandle getFEEntityHandle() const
Return finite element entity handle.
MoFEMErrorCode getProblemColIndices(const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
Get col indices.
int getLoopSize() const
get size of elements in the loop
EntityType getFEType() const
Get dimension of finite element.
virtual boost::weak_ptr< ForcesAndSourcesCore > getSubPipelinePtr() const
MoFEMErrorCode loopSide(const string &fe_name, ForcesAndSourcesCore *side_fe, const size_t dim, const EntityHandle ent_for_side=0, boost::shared_ptr< Range > fe_range=nullptr, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy, AdjCache *adj_cache=nullptr)
User calls this function to loop over elements on the side of face. This function calls finite elemen...
MoFEMErrorCode loopRange(const string &fe_name, ForcesAndSourcesCore *range_fe, boost::shared_ptr< Range > fe_range, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
Iterate over range of elements.
int getNinTheLoop() const
get number of finite element in the loop
auto getFTensor1CoordsAtGaussPts()
Get coordinates at integration points assuming linear geometry.
int getNumberOfNodesOnElement() const
Get the number of nodes on finite element.
boost::weak_ptr< SideNumber > getSideNumberPtr(const int side_number, const EntityType type)
Get the side number pointer.
const SnesMethod::SNESContext getSNESCtx() const
MoFEMErrorCode loopParent(const string &fe_name, ForcesAndSourcesCore *parent_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
User calls this function to loop over parent elements. This function calls finite element with its op...
EntityHandle getSideEntity(const int side_number, const EntityType type)
Get the side entity.
MatrixDouble & getCoordsAtGaussPts()
Gauss points and weight, matrix (nb. of points x 3)
auto getFTensor0IntegrationWeight()
Get integration weights.
double getMeasure() const
get measure of element
std::string getFEName() const
Get name of the element.
MoFEMErrorCode getProblemRowIndices(const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
Get row indices.
boost::shared_ptr< const NumeredEntFiniteElement > getNumeredEntFiniteElementPtr() const
Return raw pointer to NumeredEntFiniteElement.
std::map< EntityHandle, std::vector< boost::weak_ptr< NumeredEntFiniteElement > > > AdjCache
void addOpType(const OpType type)
Add operator type.
OpType
Controls loop over entities on element.
@ OPCOL
operator doWork function is executed on FE columns
@ OPROW
operator doWork function is executed on FE rows
@ OPROWCOL
operator doWork is executed on FE rows &columns
@ OPSPACE
operator do Work is execute on space data
MoFEMErrorCode loopChildren(const string &fe_name, ForcesAndSourcesCore *child_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
User calls this function to loop over parent elements. This function calls finite element with its op...
virtual MoFEMErrorCode setPtrFE(ForcesAndSourcesCore *ptr)
const FEMethod * getFEMethod() const
Return raw pointer to Finite Element Method object.
int getFEDim() const
Get dimension of finite element.
MatrixDouble & getGaussPts()
matrix of integration (Gauss) points for Volume Element
MoFEMErrorCode loopThis(const string &fe_name, ForcesAndSourcesCore *this_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
User calls this function to loop over the same element using a different set of integration points....
structure to get information from mofem into EntitiesFieldData
int getMaxRowOrder() const
Get max order of approximation for field in rows.
MoFEMErrorCode calHierarchicalBaseFunctionsOnElement()
Calculate base functions.
auto & getDataOnElementBySpaceArray()
Get data on entities and space.
MoFEMErrorCode getEntityFieldData(EntitiesFieldData &data, const int bit_number, const EntityType type_lo=MBVERTEX, const EntityType type_hi=MBPOLYHEDRON) const
MoFEMErrorCode loopOverOperators()
Iterate user data operators.
auto & getUserPolynomialBase()
Get the User Polynomial Base object.
MoFEMErrorCode getEntitySense(EntitiesFieldData &data) const
Get the entity sense (orientation)
MoFEMErrorCode getEntityDataOrder(const EntityType type, const FieldSpace space, boost::ptr_vector< EntitiesFieldData::EntData > &data) const
Get the entity data order.
MoFEMErrorCode getProblemNodesColIndices(const std::string &field_name, VectorInt &nodes_indices) const
MoFEMErrorCode getRowNodesIndices(EntitiesFieldData &data, const int bit_number) const
get row node indices from FENumeredDofEntity_multiIndex
const std::array< boost::shared_ptr< EntitiesFieldData >, LASTSPACE > derivedDataOnElement
Entity data on element entity columns fields.
MoFEMErrorCode getNoFieldRowIndices(EntitiesFieldData &data, const int bit_number) const
get col NoField indices
virtual MoFEMErrorCode setGaussPts(int order_row, int order_col, int order_data)
set user specific integration rule
MoFEMErrorCode getNoFieldEntityFieldData(EntitiesFieldData &data, const int bit_number, EXTRACTOR &&extractor) const
Get field data on entities where no field is defined.
MoFEMErrorCode getEntityIndices(EntitiesFieldData &data, const int bit_number, FieldEntity_vector_view &ents_field, const EntityType type_lo, const EntityType type_hi, EXTRACTOR &&extractor) const
int getMaxDataOrder() const
Get max order of approximation for data fields.
auto & getElementPolynomialBase()
Get the Entity Polynomial Base object.
ForcesAndSourcesCore * sidePtrFE
Element to integrate on the sides.
MoFEMErrorCode getProblemTypeRowIndices(const std::string &field_name, EntityType type, int side_number, VectorInt &indices) const
MoFEMErrorCode getColNodesIndices(EntitiesFieldData &data, const int bit_number) const
get col node indices from FENumeredDofEntity_multiIndex
MoFEMErrorCode getSpacesAndBaseOnEntities(EntitiesFieldData &data) const
Get field approximation space and base on entities.
MoFEMErrorCode getNodesIndices(const int bit_number, FieldEntity_vector_view &ents_field, VectorInt &nodes_indices, VectorInt &local_nodes_indices, EXTRACTOR &&extractor) const
get node indices
MatrixDouble coordsAtGaussPts
coordinated at gauss points
virtual MoFEMErrorCode operator()()
Main operator function executed for each loop iteration.
boost::function< MoFEMErrorCode(ForcesAndSourcesCore *fe_raw_ptr, int order_row, int order_col, int order_data)> GaussHookFun
MoFEMErrorCode setSideFEPtr(const ForcesAndSourcesCore *side_fe_ptr)
Set the pointer to face element on the side.
virtual MoFEMErrorCode preProcess()
Pre-processing function executed at loop initialization.
boost::ptr_deque< UserDataOperator > opPtrVector
Vector of finite element users data operators.
MoFEMErrorCode getNodesFieldData(EntitiesFieldData &data, const int bit_number) const
Get data on nodes.
int getMaxColOrder() const
Get max order of approximation for field in columns.
const std::array< boost::shared_ptr< EntitiesFieldData >, LASTSPACE > dataOnElement
Entity data on element entity rows fields.
auto & getDerivedDataOnElementBySpaceArray()
Get derived data on entities and space.
MoFEMErrorCode getEntityRowIndices(EntitiesFieldData &data, const int bit_number, const EntityType type_lo=MBVERTEX, const EntityType type_hi=MBPOLYHEDRON) const
MoFEMErrorCode getProblemTypeIndices(const std::string &field_name, const NumeredDofEntity_multiIndex &dofs, EntityType type, int side_number, VectorInt &indices) const
get indices by type (generic function) which are declared for problem but not this particular element
MatrixDouble gaussPts
Matrix of integration points.
auto & getEntData(const FieldSpace space, const EntityType type, const int side)
Get the entity data.
MoFEMErrorCode getEntityDataOrder(EntitiesFieldData &data, const FieldSpace space) const
Get the entity data order for given space.
MoFEMErrorCode getProblemNodesRowIndices(const std::string &field_name, VectorInt &nodes_indices) const
MoFEMErrorCode getEntityColIndices(EntitiesFieldData &data, const int bit_number, const EntityType type_lo=MBVERTEX, const EntityType type_hi=MBPOLYHEDRON) const
MoFEMErrorCode getFaceNodes(EntitiesFieldData &data) const
Get nodes on faces.
EntityType lastEvaluatedElementEntityType
Last evaluated type of element entity.
boost::function< int(int order_row, int order_col, int order_data)> RuleHookFun
virtual MoFEMErrorCode postProcess()
Post-processing function executed at loop completion.
boost::ptr_deque< UserDataOperator > & getOpPtrVector()
Use to push back operator for row operator.
GaussHookFun setRuleHook
Set function to calculate integration rule.
boost::shared_ptr< BaseFunction > elementPolynomialBasePtr
Pointer to entity polynomial base.
MoFEMErrorCode getProblemTypeColIndices(const std::string &field_name, EntityType type, int side_number, VectorInt &indices) const
ForcesAndSourcesCore * refinePtrFE
Element to integrate parent or child.
MoFEMErrorCode getNoFieldColIndices(EntitiesFieldData &data, const int bit_number) const
get col NoField indices
MoFEMErrorCode getProblemNodesIndices(const std::string &field_name, const NumeredDofEntity_multiIndex &dofs, VectorInt &nodes_indices) const
get indices of nodal indices which are declared for problem but not this particular element
virtual int getRule(int order_row, int order_col, int order_data)
another variant of getRule
MoFEMErrorCode createDataOnElement(EntityType type)
Create a entity data on element object.
MoFEMErrorCode getEntitySense(const EntityType type, boost::ptr_vector< EntitiesFieldData::EntData > &data) const
get sense (orientation) of entity
MoFEMErrorCode calBernsteinBezierBaseFunctionsOnElement()
Calculate Bernstein-Bezier base.
boost::shared_ptr< BaseFunction > userPolynomialBasePtr
Pointer to user polynomial base.
RuleHookFun getRuleHook
Hook to get rule.
MoFEMErrorCode setRefineFEPtr(const ForcesAndSourcesCore *refine_fe_ptr)
Set the pointer to face element refined.
Mat & ksp_B
Reference to preconditioner matrix in KSP context.
KSPContext
Context enumeration for KSP solver phases.
Mat & ksp_A
Reference to system matrix in KSP context.
Vec & ksp_f
Reference to residual vector in KSP context.
Operator for broken loop side.
Copy geometry-related data from one element to other.
Iterate over range of (sub)elements.
OpLoopRange(MoFEM::Interface &m_field, const std::string fe_name, boost::shared_ptr< Range > fe_range, const LogManager::SeverityLevel sev=Sev::noisy)
boost::shared_ptr< E > & getRangeFEPtr()
boost::weak_ptr< ForcesAndSourcesCore > getSubPipelinePtr() const override
boost::ptr_deque< UserDataOperator > & getOpPtrVector()
const std::string rangeFEName
const LogManager::SeverityLevel sevLevel
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
boost::shared_ptr< Range > feRange
boost::shared_ptr< E > rangeFEPtr
Element used to execute operators on side of the element.
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
boost::shared_ptr< Range > feRange
boost::shared_ptr< E > & getSideFEPtr()
boost::weak_ptr< ForcesAndSourcesCore > getSubPipelinePtr() const override
const std::string sideFEName
boost::shared_ptr< AdjCache > adjCache
boost::shared_ptr< E > sideFEPtr
OpLoopSide(MoFEM::Interface &m_field, const std::string fe_name, const int side_dim, boost::shared_ptr< Range > fe_range, const LogManager::SeverityLevel sev=Sev::noisy, boost::shared_ptr< AdjCache > adj_cache=nullptr)
boost::ptr_deque< UserDataOperator > & getOpPtrVector()
const LogManager::SeverityLevel sevLevel
OpLoopSide(MoFEM::Interface &m_field, const std::string fe_name, const int side_dim, const LogManager::SeverityLevel sev=Sev::noisy, boost::shared_ptr< AdjCache > adj_cache=nullptr)
Construct a new Op Loop Side object.
static constexpr Switches CtxSetA
Jacobian matrix switch.
static constexpr Switches CtxSetX_TT
Second time derivative switch.
std::bitset< 8 > Switches
Bitset type for context switches.
Switches data_ctx
Current data context switches.
static constexpr Switches CtxSetX_T
First time derivative switch.
static constexpr Switches CtxSetB
Preconditioner matrix switch.
SNESContext
Context enumeration for SNES solver phases.
Vec & snes_f
Reference to residual vector.
Vec & snes_x
Reference to current solution state vector.
Mat & snes_B
Reference to preconditioner of Jacobian matrix.
Mat & snes_A
Reference to Jacobian matrix.
TSContext
Context enumeration for TS solver phases.
Mat & ts_A
Reference to Jacobian matrix: dF/dU + a*dF/dU_t + aa*dF/dU_tt.
Vec & ts_F
Reference to residual vector F(t,U,U_t,U_tt)
Vec & ts_u_tt
Reference to second time derivative of state vector d²U/dt²
Vec & ts_u_t
Reference to first time derivative of state vector dU/dt.
Vec & ts_u
Reference to current state vector U(t)
Mat & ts_B
Reference to preconditioner matrix for ts_A.
Base volume element used to integrate on contact surface (could be extended to other volume elements)...