v0.15.5
Loading...
Searching...
No Matches
EntitiesFieldData.hpp
Go to the documentation of this file.
1/** \file EntitiesFieldData.hpp
2
3\brief Data structures for accessing information about finite element and its
4degrees of freedom.
5
6*/
7
8
9
10#ifndef __ENTITIES_FIELD_DATA_HPP__
11#define __ENTITIES_FIELD_DATA_HPP__
12
13using namespace boost::numeric;
14
15namespace MoFEM {
16
17using DofsAllocator = ublas::unbounded_array<
18
19 FEDofEntity *, std::allocator<FEDofEntity *>
20
21 >;
22
23using VectorDofs = ublas::vector<FEDofEntity *, DofsAllocator>;
24
25using FieldEntAllocator = ublas::unbounded_array<
26
27 FieldEntity *, std::allocator<FieldEntity *>
28
29 >;
30
31using VectorFieldEntities = ublas::vector<FieldEntity *, FieldEntAllocator>;
32
33/** \brief data structure for finite element entity
34 * \ingroup mofem_forces_and_sources_user_data_operators
35 *
36 * It keeps that about indices of degrees of freedom, dofs data, base functions
37 * functions, entity side number, type of entities, approximation order, etc.
38 *
39 */
41
42 struct EntData;
43
44 std::bitset<LASTBASE> bAse; ///< bases on element
45 MatrixInt facesNodes; ///< nodes on finite element faces
46 MatrixInt facesNodesOrder; ///< order of face nodes on element
47
48 std::array<std::bitset<LASTSPACE>, MBMAXTYPE>
49 spacesOnEntities; ///< spaces on entity types
50 std::array<std::bitset<LASTBASE>, MBMAXTYPE>
51 basesOnEntities; ///< bases on entity types
52 std::array<std::bitset<LASTBASE>, LASTSPACE>
53 basesOnSpaces; ///< base on spaces
54 std::array<std::bitset<LASTBASE>, LASTSPACE>
55 brokenBasesOnSpaces; ///< base on spaces
56 std::array<boost::ptr_vector<EntData>, MBMAXTYPE>
57 dataOnEntities; ///< data on nodes, base
58 ///< function, dofs
59 ///< values, etc.
60
61 /**
62 * @brief Construct EntitiesFieldData for given entity type
63 *
64 * @param type Entity type (MBVERTEX, MBEDGE, MBTRI, MBTET, etc.)
65 */
66 EntitiesFieldData(const EntityType type);
67 virtual ~EntitiesFieldData() = default;
68
69 /**
70 * @brief Set element type and initialize data structures
71 *
72 * @param type Entity type to set
73 * @return MoFEMErrorCode Error code (0 on success)
74 */
75 virtual MoFEMErrorCode setElementType(const EntityType type);
76
77 /**
78 * @brief Reset data associated with particular field name
79 *
80 * @return MoFEMErrorCode Error code (0 on success)
81 */
83
84 /**
85 * @brief Swap approximation base
86 *
87 * Bernstein-Bezier (BB) base is not hierarchical, and is calculated for
88 * particular field, since it all shape functions change with the order. BB
89 * base is precalculated for every field, and when user push operator with
90 * particular field using BB base, pointers to shape functions and
91 * BaseDerivatives of shape functions are set to particular location, once
92 * operator is executed, pointers are switch back to its original position.
93 *
94 * getNSharedPtr(base) <=== getBBNSharedPtr(field_name);
95 * // DO OPERATOR WORK
96 * getNSharedPtr(base) ==> getBBNSharedPtr(field_name);
97 *
98 * @param field_name Name of the field for which to swap base
99 * @param base Approximation base type to swap
100 * @return MoFEMErrorCode Error code (0 on success)
101 */
102 virtual MoFEMErrorCode baseSwap(const std::string &field_name,
103 const FieldApproximationBase base);
104
105 friend std::ostream &operator<<(std::ostream &os, const EntitiesFieldData &e);
106
107protected:
109};
110
111/** \brief this class derives data from other data structure
112 * \ingroup mofem_forces_and_sources_user_data_operators
113 *
114 *
115 * It behaves like normal data structure it is used to share base functions with
116 * other data structures. Dofs values, approx. order and
117 * indices are not shared.
118 *
119 * Shape functions, senses are shared with other data structure.
120 *
121 */
123
124 struct DerivedEntData;
125
126 /**
127 * @brief Construct DerivedEntitiesFieldData from existing EntitiesFieldData
128 *
129 * @param data_ptr Shared pointer to source EntitiesFieldData
130 */
132 const boost::shared_ptr<EntitiesFieldData> &data_ptr);
133
134 /**
135 * @brief Set element type for derived data
136 *
137 * @param type Entity type to set
138 * @return MoFEMErrorCode Error code (0 on success)
139 */
140 MoFEMErrorCode setElementType(const EntityType type);
141
142private:
143 const boost::shared_ptr<EntitiesFieldData> dataPtr;
144};
145
146/** \brief Data on single entity (This is passed as argument to
147 * DataOperator::doWork) \ingroup mofem_forces_and_sources_user_data_operators
148 * \nosubgrouping
149 *
150 * \todo Hdiv and Hcurl functions should be accessed through common interface.
151 */
153
162
163 /** \name Constructor and destructor */
164
165 /**@{*/
166
167 /**
168 * @brief Construct EntData object
169 *
170 * @param allocate_base_matrices If true, allocate base function matrices (default: true)
171 */
172 EntData(const bool allocate_base_matrices = true);
173 virtual ~EntData() = default;
174
175 /**@}*/
176
177 /** \name Sense, order and indices */
178
179 /**@{*/
180
181 /**
182 * @brief Get entity sense for conforming approximation fields
183 *
184 * Entity sense is needed to calculate base functions with conforming
185 * approximation fields.
186 *
187 * @return int Entity sense (orientation)
188 */
189 virtual int getSense() const;
190
191 /**
192 * @brief Get approximation order
193 *
194 * @return ApproximationOrder Approximation order for this entity
195 */
196 inline ApproximationOrder getOrder() const;
197
198 /**
199 * @brief Get global indices of degrees of freedom on entity
200 *
201 * @return const VectorInt& Reference to vector of global DOF indices
202 */
203 inline const VectorInt &getIndices() const;
204
205 /**
206 * @brief Get global indices of DOFs on entity up to given order
207 *
208 * @param order Maximum approximation order to include
209 * @return const VectorIntAdaptor Adaptor to subset of global indices
210 */
212
213 /**
214 * @brief Get local indices of degrees of freedom on entity
215 *
216 * @return const VectorInt& Reference to vector of local DOF indices
217 */
218 inline const VectorInt &getLocalIndices() const;
219
220 /**
221 * @brief Get local indices of DOFs on entity up to given order
222 *
223 * @param order Maximum approximation order to include
224 * @return const VectorIntAdaptor Adaptor to subset of local indices
225 */
227
228 /**
229 * @brief Get reference to entity sense for modification
230 *
231 * @return int& Reference to entity sense
232 */
233 inline int &getSense();
234
235 /**
236 * @brief Get reference to approximation order for modification
237 *
238 * @return ApproximationOrder& Reference to approximation order
239 */
241
242 /**
243 * @brief Get reference to global indices for modification
244 *
245 * @return VectorInt& Reference to global indices vector
246 */
247 inline VectorInt &getIndices();
248
249 /**
250 * @brief Get reference to local indices for modification
251 *
252 * @return VectorInt& Reference to local indices vector
253 */
254 inline VectorInt &getLocalIndices();
255
256 /**@}*/
257
258 /** \name Data on entity */
259
260 /**@{*/
261
262 /**
263 * @brief Get DOF values on entity
264 *
265 * @return const VectorDouble& Reference to vector of DOF values
266 */
267 inline const VectorDouble &getFieldData() const;
268
269 /**
270 * @brief Get DOF values on entity up to given order
271 *
272 * @param order Maximum approximation order to include
273 * @return const VectorAdaptor Adaptor to subset of DOF values
274 */
275 inline const VectorAdaptor getFieldDataUpToOrder(int order);
276
277 /**
278 * @brief Get DOF data structures (const version)
279 *
280 * @return const VectorDofs& Reference to vector of FEDofEntity pointers
281 */
282 inline const VectorDofs &getFieldDofs() const;
283
284 /**
285 * @brief Get DOF values on entity for modification
286 *
287 * @return VectorDouble& Reference to vector of DOF values
288 */
289 inline VectorDouble &getFieldData();
290
291 /**
292 * @brief Get field entities (const version)
293 *
294 * @return const VectorFieldEntities& Reference to vector of FieldEntity pointers
295 */
296 inline const VectorFieldEntities &getFieldEntities() const;
297
298 /**
299 * @brief Get field entities for modification
300 *
301 * @return VectorFieldEntities& Reference to vector of FieldEntity pointers
302 */
304
305 /**
306 * @brief Get entity bit reference level
307 *
308 * @return std::vector<BitRefLevel>& Reference to bit reference level vector
309 */
310 virtual std::vector<BitRefLevel> &getEntDataBitRefLevel();
311
312 /**
313 * @brief Return FTensor of rank 1, i.e. vector from field data coefficients
314 *
315 * \code
316 * auto t_vec = data.getFTensor1FieldData<3>();
317 * \endcode
318 *
319 * @tparam Tensor_Dim size of vector
320 * @return FTensor::Tensor1<FTensor::PackPtr<double *, Tensor_Dim>,
321 * Tensor_Dim>
322 */
323 template <int Tensor_Dim>
326
327 /**
328 * @brief Return FTensor rank 2, i.e. matrix from field data coefficients
329 *
330 * \code
331 * auto t_mat = data.getFTensor2FieldData<3,3>();
332 * \endcode
333 *
334 * @tparam Tensor_Dim0 First dimension of tensor
335 * @tparam Tensor_Dim1 Second dimension of tensor
336 * @return FTensor::Tensor2<FTensor::PackPtr<double *, Tensor_Dim0 *
337 * Tensor_Dim1>, Tensor_Dim0, Tensor_Dim1> Rank 2 tensor
338 */
339 template <int Tensor_Dim0, int Tensor_Dim1>
341 Tensor_Dim0, Tensor_Dim1>
343
344 /**
345 * @brief Return symmetric FTensor rank 2, i.e. matrix from field data
346 * coefficients
347 *
348 * \code
349 * auto t_mat = data.getFTensor2SymmetricFieldData<3>();
350 * \endcode
351 *
352 * @tparam Tensor_Dim Dimension of the symmetric tensor
353 * @return FTensor::Tensor2_symmetric<FTensor::PackPtr<double *, (Tensor_Dim
354 * * (Tensor_Dim + 1)) / 2>, Tensor_Dim> Symmetric rank 2 tensor
355 */
356 template <int Tensor_Dim>
358 FTensor::PackPtr<double *, (Tensor_Dim * (Tensor_Dim + 1)) / 2>,
359 Tensor_Dim>
361
362 /**
363 * @brief Return scalar files as a FTensor of rank 0
364 *
365 * @return FTensor::Tensor0<FTensor::PackPtr<double *,1> >
366 */
368
369 inline VectorDofs &getFieldDofs();
370
371 /**@}*/
372
373 /** \name Base and space */
374
375 /**@{*/
376
377 /**
378 * \brief Get approximation base
379 * @return Approximation base
380 */
382
383 /**
384 * \brief Get field space
385 * @return Field space
386 */
387 inline FieldSpace &getSpace();
388
389 /**
390 * @brief Get shared pointer to base functions with derivatives
391 *
392 * @param base Approximation base type
393 * @param derivative Derivative order (ZeroDerivative, FirstDerivative, etc.)
394 * @return boost::shared_ptr<MatrixDouble>& Shared pointer to base function matrix
395 */
396 virtual boost::shared_ptr<MatrixDouble> &
398 const BaseDerivatives derivative);
399
400 /**
401 * @brief Get shared pointer to base functions (zero derivative)
402 *
403 * @param base Approximation base type
404 * @return boost::shared_ptr<MatrixDouble>& Shared pointer to base function matrix
405 */
406 virtual boost::shared_ptr<MatrixDouble> &
408
409 /**
410 * @brief Get shared pointer to derivatives of base functions
411 *
412 * @param base Approximation base type
413 * @return boost::shared_ptr<MatrixDouble>& Shared pointer to derivative matrix
414 */
415 virtual boost::shared_ptr<MatrixDouble> &
417
418 /**@}*/
419
420 /** \name Get base functions for H1/L2 */
421
422 /**@{*/
423
424 /** \brief get base functions
425 * this return matrix (nb. of rows is equal to nb. of Gauss pts, nb. of
426 * columns is equal to number of base functions on this entity.
427 *
428 * \note Note that for vectorial base, like Hdiv or Hcurl, in columns are
429 * vectorial base functions. For tonsorial would be base
430 * functions. Interpretation depends on type of base, scalar, vectorial or
431 * tonsorial and dimension of problem.
432 *
433 */
434 inline MatrixDouble &getN(const FieldApproximationBase base);
435
436 /**
437 * @brief Get base functions for specific field (Bernstein-Bezier)
438 *
439 * @param field_name Name of the field
440 * @return MatrixDouble& Reference to base function matrix
441 */
442 inline MatrixDouble &getN(const std::string &field_name);
443
444 /**
445 * @brief Get base functions using default base
446 *
447 * @return MatrixDouble& Reference to base function matrix
448 */
449 inline MatrixDouble &getN();
450
451 /** \brief get derivatives of base functions
452 *
453 * Matrix at rows has nb. of Gauss pts, at columns it has derivative of
454 * base functions. Columns are structured as follows, [ dN1/dx, dN1/dy,
455 * dN1/dz, dN2/dx, dN2/dy, dN2/dz, ... ]
456 *
457 * Scalar base functions:
458 * Note that base functions are calculated in file H1.c
459 * Above description not apply for derivatives of nodal functions, since
460 * derivative of nodal functions in case of simplexes, EDGES, TRIANGLES and
461 * TETS are constant. So that matrix rows represents nb. of base
462 * functions, columns are derivatives. Nb. of columns depend on element
463 * dimension, for EDGES is one, for TRIS is 2 and TETS is 3.
464 *
465 * \note Note that for node element this function make no sense.
466 *
467 * Tensorial base functions:
468 * \note Note: In rows are integration pts, columns are formatted that that
469 * components of vectors and then derivatives, for example row for given
470 * integration points is formatted in array
471 * \f[
472 * t_{0,0}, t_{1,0}, t_{1,0}, t_{0,1}, t_{1,1}, t_{1,1}, t_{0,2}, t_{1,2},
473 * t_{1,2} \f] where comma express derivative, i.e. \f$t_{2,1} =
474 * \frac{\partial t_2}{\partial \xi_1}\f$
475 *
476 */
477 inline MatrixDouble &getDiffN(const FieldApproximationBase base);
478
479 /**
480 * @brief Get base function derivatives with specific derivative order
481 *
482 * @param base Approximation base type
483 * @param derivative Derivative order
484 * @return MatrixDouble& Reference to derivative matrix
485 */
486 inline MatrixDouble &getN(const FieldApproximationBase base,
487 const BaseDerivatives derivative);
488
489 /**
490 * @brief Get derivatives of base functions for specific field
491 *
492 * @param field_name Name of the field
493 * @return MatrixDouble& Reference to derivative matrix
494 */
495 inline MatrixDouble &getDiffN(const std::string &field_name);
496
497 /**
498 * @brief Get derivatives of base functions using default base
499 *
500 * @return MatrixDouble& Reference to derivative matrix
501 */
502 inline MatrixDouble &getDiffN();
503
504 /**
505 * @brief Get base function derivatives using default base
506 *
507 * @param derivative Derivative order
508 * @return MatrixDouble& Reference to derivative matrix
509 */
510 inline MatrixDouble &getN(const BaseDerivatives derivative);
511
512 /**
513 * @brief Get base functions at specific Gauss point
514 *
515 * @param base Approximation base type
516 * @param gg Gauss point number
517 * @return const VectorAdaptor Adaptor to base functions at Gauss point
518 */
519 inline const VectorAdaptor getN(const FieldApproximationBase base,
520 const int gg);
521
522 /**
523 * @brief Get base functions at specific Gauss point using default base
524 *
525 * @param gg Gauss point number
526 * @return const VectorAdaptor Adaptor to base functions at Gauss point
527 */
528 inline const VectorAdaptor getN(const int gg);
529
530 /**
531 * @brief Get derivative of base functions at specific Gauss point
532 *
533 * Returned matrix has base functions in rows and derivatives in columns.
534 *
535 * @param base Approximation base type
536 * @param gg Gauss point number
537 * @return const MatrixAdaptor Adaptor to derivative matrix at Gauss point
538 */
539 inline const MatrixAdaptor getDiffN(const FieldApproximationBase base,
540 const int gg);
541
542 /**
543 * @brief Get derivative of base functions at specific Gauss point using default base
544 *
545 * @param gg Gauss point number
546 * @return const MatrixAdaptor Adaptor to derivative matrix at Gauss point
547 */
548 inline const MatrixAdaptor getDiffN(const int gg);
549
550 /**
551 * @brief Get base functions at Gauss point with limited number of functions
552 *
553 * Note that for multi-field elements, different fields can have different
554 * approximation orders. Since hierarchical approximation basis is used,
555 * base functions are calculated once for element using maximal order.
556 *
557 * @param base Approximation base type
558 * @param gg Gauss point number
559 * @param nb_base_functions Number of base functions to return
560 * @return const VectorAdaptor Adaptor to limited base functions
561 */
562 inline const VectorAdaptor getN(const FieldApproximationBase base,
563 const int gg, const int nb_base_functions);
564
565 /**
566 * @brief Get base functions at Gauss point with limited number (default base)
567 *
568 * @param gg Gauss point number
569 * @param nb_base_functions Number of base functions to return
570 * @return const VectorAdaptor Adaptor to limited base functions
571 */
572 inline const VectorAdaptor getN(const int gg, const int nb_base_functions);
573
574 /**
575 * @brief Get derivatives of base functions at Gauss point with limited number
576 *
577 * Note that for multi-field elements, different fields can have different
578 * approximation orders. Since hierarchical approximation basis is used,
579 * base functions are calculated once for element using maximal order.
580 *
581 * @param base Approximation base type
582 * @param gg Gauss point number
583 * @param nb_base_functions Number of base functions
584 * @return const MatrixAdaptor Adaptor to limited derivative matrix
585 */
586 inline const MatrixAdaptor getDiffN(const FieldApproximationBase base,
587 const int gg,
588 const int nb_base_functions);
589
590 /**
591 * @brief Get derivatives of base functions at Gauss point with limited number (default base)
592 *
593 * @param gg Gauss point number
594 * @param nb_base_functions Number of base functions
595 * @return const MatrixAdaptor Adaptor to limited derivative matrix
596 */
597 inline const MatrixAdaptor getDiffN(const int gg,
598 const int nb_base_functions);
599
600 /**@}*/
601
602 /** \name Get base functions for vectorial approximation bases, i.e.
603 * Hdiv/Hcurl */
604
605 /**@{*/
606
607 /** \brief get Hdiv of base functions at Gauss pts
608 *
609 * \param base Approximation base
610 * \param gg nb. of Gauss point
611 *
612 */
613 template <int DIM>
614 inline const MatrixAdaptor getVectorN(const FieldApproximationBase base,
615 const int gg);
616
617 /** \brief get Hdiv of base functions at Gauss pts
618 *
619 * \param gg nb. of Gauss point
620 * \param number of of base functions
621 *
622 */
623 template <int DIM> inline const MatrixAdaptor getVectorN(const int gg);
624
625 /** \brief get DiffHdiv of base functions at Gauss pts
626 *
627 * \param base Approximation base
628 * \param gg nb. of Gauss point
629 * \param number of of base functions
630 *
631 */
632 template <int DIM0, int DIM1>
634 const int gg);
635
636 /** \brief get DiffHdiv of base functions at Gauss pts
637 *
638 * \param gg nb. of Gauss point
639 * \param number of of base functions
640 *
641 */
642 template <int DIM0, int DIM1>
643 inline const MatrixAdaptor getVectorDiffN(const int gg);
644
645 /** \brief get DiffHdiv of base functions at Gauss pts
646 *
647 * \param base Approximation base
648 * \param gg nb. of Gauss point
649 * \param number of of base functions
650 *
651 */
652 template <int DIM0, int DIM1>
654 const int dof, const int gg);
655
656 /** \brief get DiffHdiv of base functions at Gauss pts
657 *
658 * \param gg nb. of Gauss point
659 * \param number of of base functions
660 *
661 */
662 template <int DIM0, int DIM1>
663 inline const MatrixAdaptor getVectorDiffN(const int dof, const int gg);
664
665 /**@}*/
666
667 /** \name Get base functions with FTensor */
668
669 /**@{*/
670
671 /**
672 * \brief Get base function as Tensor0
673 *
674 * \param base
675 * \return Tensor0
676 *
677 */
680
681 /**
682 * \brief Get base function as Tensor0
683 *
684 * Return base functions for field base
685 *
686 * \return Tensor0
687 *
688 */
690
691 /**
692 * \brief Get base function as Tensor0 (Loop by integration points)
693 *
694 * \param base
695 * \param gg integration points
696 * \param bb base function
697 * \return Tensor0
698
699 Note that:
700 \code
701 t0 = data.getFTensor0N(base,bb);
702 ++t0
703 \endcode
704 Increment in above code will move pointer to base function in next
705 integration point.
706
707 *
708 */
710 getFTensor0N(const FieldApproximationBase base, const int gg, const int bb);
711
712 /**
713 * \brief Get base function as Tensor0 (Loop by integration points)
714 *
715 * Return base functions for field base
716 *
717 * \param bb base function
718 * \return Tensor0
719 *
720 */
722 getFTensor0N(const int gg, const int bb);
723
724 /**
725 * \brief Get derivatives of base functions
726 *
727 * For volume element like tetrahedral or prism,
728 * \code
729 * Tensor1<double*,3> diff_base = data.getFTensor1DiffN<3>();
730 * \endcode
731 *
732 * For face element like triangle or quad
733 * \code
734 * Tensor1<double*,2> diff_base = data.getFTensor1DiffN<2>();
735 * \endcode
736 *
737 * \param base functions
738 * \return Tensor rank 1 (vector)
739 *
740 */
741 template <int Tensor_Dim>
744
745 /**
746 * \brief Get derivatives of base functions
747 *
748 * For volume element like tetrahedral or prism,
749 * \code
750 * Tensor1<double*,3> diff_base = data.getFTensor1DiffN<3>();
751 * \endcode
752 *
753 * For face element like triangle or quad
754 * \code
755 * Tensor1<double*,2> diff_base = data.getFTensor1DiffN<2>();
756 * \endcode
757 *
758 * \return Tensor rank 1 (vector)
759 *
760 */
761 template <int Tensor_Dim>
764
765 /**
766 * \brief Get derivatives of base functions (Loop by integration points)
767 *
768 * For volume element like tetrahedral or prism,
769 * \code
770 * Tensor1<double*,3> diff_base = data.getFTensor1DiffN<3>(base,gg,bb);
771 * \endcode
772 * where bb is base function and gg is integration pt. Operator ++diff_base
773 * will move tensor pointer to next integration point.
774 *
775 * For face element like triangle or quad
776 * \code
777 * Tensor1<double*,2> diff_base = data.getFTensor1DiffN<2>(base,gg,bb);
778 * \endcode
779 *
780 * \return Tensor rank 1 (vector)
781 *
782 */
783 template <int Tensor_Dim>
785 getFTensor1DiffN(const FieldApproximationBase base, const int gg,
786 const int bb);
787
788 /**
789 * \brief Get derivatives of base functions (Loop by integration points)
790 *
791 * For volume element like tetrahedral or prism,
792 * \code
793 * Tensor1<double*,3> diff_base = data.getFTensor1DiffN<3>(gg,bb);
794 * \endcode
795 * where bb is base function and gg is integration pt. Operator ++diff_base
796 * will move tensor pointer to next base function.
797 *
798 * For face element like triangle or quad
799 * \code
800 * Tensor1<double*,2> diff_base = data.getFTensor1DiffN<2>(gg,bb);
801 * \endcode
802 *
803 * \return Tensor rank 1 (vector)
804 *
805 */
806 template <int Tensor_Dim>
808 getFTensor1DiffN(const int gg, const int bb);
809
810 /** \brief Get second derivatives of scalar base functions
811 */
812 template <int Tensor_Dim>
814 Tensor_Dim, Tensor_Dim>
816
817 /** \brief Get second derivatives of scalar base functions at integration
818 * pts
819 */
820 template <int Tensor_Dim>
822 Tensor_Dim, Tensor_Dim>
824 const int bb);
825
826 /** \brief Get second derivatives of scalar base functions
827 */
828 template <int Tensor_Dim>
830 Tensor_Dim, Tensor_Dim>
832
833 /** \brief Get second derivatives of scalar base functions at integration
834 * pts
835 */
836 template <int Tensor_Dim>
838 Tensor_Dim, Tensor_Dim>
839 getFTensor2DiffN2(const int gg, const int bb);
840
841
842 /** \brief Get base functions for Hdiv/Hcurl spaces
843
844 \note You probably like to use getFTensor1N(), in typical use base is
845 set automatically based on base set to field.
846
847 * @param base Approximation base
848
849 Example:
850 \code
851 FTensor::Index<'i',3> i;
852 int nb_dofs = data.getFieldData().size();
853 auto t_n_hdiv = data.getFTensor1N<3>();
854 for(int gg = 0;gg!=nb_gauss_pts;gg++) {
855 int ll = 0;
856 for(;ll!=nb_dofs;ll++) {
857 double dot_product = t_n_hdiv(i)*t_n_hdiv(i);
858 ++t_n_hdiv;
859 }
860 for(;ll!=data.getVectorN().size2()/3;ll++) {
861 ++t_n_hdiv;
862 }
863 }
864 \endcode
865
866 */
867 template <int Tensor_Dim>
870
871 /** \brief Get base functions for Hdiv space
872
873 Example:
874 \code
875 FTensor::Index<'i',3> i;
876 int nb_dofs = data.getFieldData().size();
877 auto t_n_hdiv = data.getFTensor1N<3>();
878 for(int gg = 0;gg!=nb_gauss_pts;gg++) {
879 int ll = 0;
880 for(;ll!=nb_dofs;ll++) {
881 double dot_product = t_n_hdiv(i)*t_n_hdiv(i);
882 ++t_n_hdiv;
883 }
884 for(;ll!=data.getVectorN().size2()/3;ll++) {
885 ++t_n_hdiv;
886 }
887 }
888 \endcode
889
890 */
891 template <int Tensor_Dim> auto getFTensor1N();
892
893 /** \brief Get derivatives of base functions for Hdiv space
894 */
895 template <int Tensor_Dim0, int Tensor_Dim1>
897 Tensor_Dim0, Tensor_Dim1>
899
900 /** \brief Get derivatives of base functions for Hdiv space at integration
901 * pts
902 */
903 template <int Tensor_Dim0, int Tensor_Dim1>
905 Tensor_Dim0, Tensor_Dim1>
906 getFTensor2DiffN(FieldApproximationBase base, const int gg, const int bb);
907
908 /** \brief Get derivatives of base functions for Hdiv space
909 */
910 template <int Tensor_Dim0, int Tensor_Dim1>
912 Tensor_Dim0, Tensor_Dim1>
914 return getFTensor2DiffN<Tensor_Dim0, Tensor_Dim1>(bAse);
915 }
916
917 /** \brief Get derivatives of base functions for Hdiv space at integration
918 * pts
919 */
920 template <int Tensor_Dim0, int Tensor_Dim1>
922 Tensor_Dim0, Tensor_Dim1>
923 getFTensor2DiffN(const int gg, const int bb) {
924 return getFTensor2DiffN<Tensor_Dim0, Tensor_Dim1>(bAse, gg, bb);
925 }
926
927 /** \brief Get second derivatives of base functions for Hvec space
928 */
929 template <int Tensor_Dim0, int Tensor_Dim1, int Tensor_Dim2>
932 Tensor_Dim0, Tensor_Dim1, Tensor_Dim2>
934
935 /** \brief Get second derivatives of base functions for Hvec space
936 */
937 template <int Tensor_Dim0, int Tensor_Dim1, int Tensor_Dim2>
938 inline FTensor::Tensor3<
940 Tensor_Dim0, Tensor_Dim1, Tensor_Dim2>
942 return getFTensor3Diff2N<Tensor_Dim0, Tensor_Dim1, Tensor_Dim2>(bAse);
943 }
944
945 /** \brief Get derivatives of base functions for tonsorial Hdiv space
946 */
947 template <int Tensor_Dim0, int Tensor_Dim1, int Tensor_Dim2>
950 Tensor_Dim0, Tensor_Dim1, Tensor_Dim2>
952
953 /** \brief Get derivatives of base functions for tonsorial Hdiv space at integration
954 * pts
955 */
956 template <int Tensor_Dim0, int Tensor_Dim1, int Tensor_Dim2>
959 Tensor_Dim0, Tensor_Dim1, Tensor_Dim2>
960 getFTensor3DiffN(FieldApproximationBase base, const int gg, const int bb);
961
962 /** \brief Get derivatives of base functions for tonsorial Hdiv space
963 */
964 template <int Tensor_Dim0, int Tensor_Dim1, int Tensor_Dim2>
965 inline FTensor::Tensor3<
967 Tensor_Dim0, Tensor_Dim1, Tensor_Dim2>
969 return getFTensor3DiffN<Tensor_Dim0, Tensor_Dim1, Tensor_Dim2>(bAse);
970 }
971
972 /** \brief Get derivatives of base functions for tonsorial Hdiv space at integration
973 * pts
974 */
975 template <int Tensor_Dim0, int Tensor_Dim1, int Tensor_Dim2>
976 inline FTensor::Tensor3<
978 Tensor_Dim0, Tensor_Dim1, Tensor_Dim2>
979 getFTensor3DiffN(const int gg, const int bb) {
980 return getFTensor3DiffN<Tensor_Dim0, Tensor_Dim1, Tensor_Dim2>(bAse, gg,
981 bb);
982 }
983
984 /**
985 * \brief Get Hdiv base functions at integration point
986
987 \code
988 FTensor::Index<'i',3> i;
989 for(int gg = 0;gg!=nb_gauss_pts;gg++) {
990 auto t_base = data.getFTensor1N(base,gg,bb);
991 for(int bb = 0;bb!=nb_base_functions;bb++) {
992 auto dot = t_base(i)*t_base(i);
993 }
994 }
995 \endcode
996
997 */
998 template <int Tensor_Dim>
1000 getFTensor1N(FieldApproximationBase base, const int gg, const int bb);
1001
1002 /**
1003 * \brief Get Hdiv base functions at integration point
1004
1005 \code
1006 FTensor::Index<'i',3> i;
1007 for(int gg = 0;gg!=nb_gauss_pts;gg++) {
1008 auto t_base = data.getFTensor1N(gg,0);
1009 for(int bb = 0;bb!=nb_base_functions;bb++) {
1010 double dot = t_base(i)*t_base(i);
1011 }
1012 }
1013 \endcode
1014
1015 */
1016 template <int Tensor_Dim>
1017 inline auto getFTensor1N(const int gg, const int bb);
1018
1019 /** \brief Get base functions for Hdiv/Hcurl spaces
1020
1021 \note You probably like to use getFTensor1N(), in typical use base is
1022 set automatically based on base set to field.
1023
1024 * @param base Approximation base
1025
1026 Example:
1027 \code
1028 FTensor::Index<'i',3> i;
1029 FTensor::Index<'i',3> j;
1030 int nb_dofs = data.getFieldData().size();
1031 auto t_n_hdiv = data.getFTensor2N<3,3>();
1032 for(int gg = 0;gg!=nb_gauss_pts;gg++) {
1033 int ll = 0;
1034 for(;ll!=nb_dofs;ll++) {
1035 double dot_product = t_n_hdiv(i,j)*t_n_hdiv(i,j);
1036 ++t_n_hdiv;
1037 }
1038 for(;ll!=data.getVectorN().size2()/3;ll++) {
1039 ++t_n_hdiv;
1040 }
1041 }
1042 \endcode
1043
1044 */
1045 template <int Tensor_Dim0, int Tensor_Dim1>
1047 Tensor_Dim0, Tensor_Dim1>
1049
1050 /** \brief Get base functions for Hdiv space
1051
1052 Example:
1053 \code
1054 FTensor::Index<'i',3> i;
1055 FTensor::Index<'j',3> j;
1056
1057 int nb_dofs = data.getFieldData().size();
1058 auto t_n_hdiv = data.getFTensor2N<3,3>();
1059 for(int gg = 0;gg!=nb_gauss_pts;gg++) {
1060 int ll = 0;
1061 for(;ll!=nb_dofs;ll++) {
1062 double dot_product = t_n_hdiv(i,j)*t_n_hdiv(i,j);
1063 ++t_n_hdiv;
1064 }
1065 for(;ll!=data.getVectorN().size2()/3;ll++) {
1066 ++t_n_hdiv;
1067 }
1068 }
1069 \endcode
1070
1071 */
1072 template <int Tensor_Dim0, int Tensor_Dim1> auto getFTensor2N();
1073
1074 /** \brief Get base functions for tensor Hdiv/Hcurl spaces
1075
1076 \note You probably like to use getFTensor2N(), in typical use base is
1077 set automatically based on base set to field.
1078
1079 @param base Approximation base
1080
1081 Example:
1082 \code
1083 FTensor::Index<'i',3> i;
1084 FTensor::Index<'j',3> i;
1085 int nb_dofs = data.getFieldData().size();
1086 for(int gg = 0;gg!=nb_gauss_pts;gg++) {
1087 auto t_n_hdiv = data.getFTensor2N<3>(base,gg,bb);
1088 int ll = 0;
1089 for(;ll!=nb_dofs;ll++) {
1090 double dot_product = t_n_hdiv(i,j)*t_n_hdiv(i,j);
1091 ++t_n_hdiv;
1092 }
1093 for(;ll!=data.getVectorN().size2()/3;ll++) {
1094 ++t_n_hdiv;
1095 }
1096 }
1097 \endcode
1098
1099 */
1100 template <int Tensor_Dim0, int Tensor_Dim1>
1102 Tensor_Dim0, Tensor_Dim1>
1103 getFTensor2N(FieldApproximationBase base, const int gg, const int bb);
1104
1105 /** \brief Get base functions for Hdiv space
1106
1107 Example:
1108 \code
1109 FTensor::Index<'i',3> i;
1110 FTensor::Index<'j',3> j;
1111 int nb_dofs = data.getFieldData().size();
1112 for(int gg = 0;gg!=nb_gauss_pts;++gg) {
1113 int ll = 0;
1114 auto t_n_hdiv = data.getFTensor2N<3,3>(gg,0);
1115 for(;ll!=nb_dofs;ll++) {
1116 double dot_product = t_n_hdiv(i)*t_n_hdiv(i);
1117 ++t_n_hdiv;
1118 }
1119 for(;ll!=data.getVectorN().size2()/3;ll++) {
1120 ++t_n_hdiv;
1121 }
1122 }
1123 \endcode
1124
1125 */
1126 template <int Tensor_Dim0, int Tensor_Dim1>
1127 auto getFTensor2N(const int gg, const int bb);
1128
1129 /**@}*/
1130
1131 /** \name Auxiliary functions */
1132
1133 /**@{*/
1134
1135 friend std::ostream &operator<<(std::ostream &os,
1137
1138 /**
1139 * Reset data associated with particular field name
1140 * @return error code
1141 */
1143
1144 /**@}*/
1145
1146 /** \name Bernstein-Bezier base only functions */
1147
1148 /**@{*/
1149
1150 /**
1151 * @brief Get orders at the nodes
1152 *
1153 * @return VectorInt&
1154 */
1155 inline VectorInt &getBBNodeOrder();
1156
1157 /**
1158 * @brief Get file BB indices
1159 *
1160 * @return MatrixInt&
1161 */
1162 inline MatrixInt &getBBAlphaIndices();
1163
1164 virtual boost::shared_ptr<MatrixInt> &
1165 getBBAlphaIndicesSharedPtr(const std::string &field_name);
1166
1167 /**
1168 * Get shared pointer to BB base base functions
1169 */
1170 virtual boost::shared_ptr<MatrixDouble> &
1171 getBBNSharedPtr(const std::string &field_name);
1172
1173 /**
1174 * Get shared pointer to BB base base functions
1175 */
1176 virtual const boost::shared_ptr<MatrixDouble> &
1177 getBBNSharedPtr(const std::string &field_name) const;
1178
1179 /**
1180 * Get shared pointer to BB derivatives of base base functions
1181 */
1182 virtual boost::shared_ptr<MatrixDouble> &
1183 getBBDiffNSharedPtr(const std::string &field_name);
1184
1185 /**
1186 * Get shared pointer to derivatives of BB base base functions
1187 */
1188 virtual const boost::shared_ptr<MatrixDouble> &
1189 getBBDiffNSharedPtr(const std::string &field_name) const;
1190
1191 virtual std::map<std::string, boost::shared_ptr<MatrixInt>> &
1193
1194 /**
1195 * @brief get hash map of base function for BB base, key is a field name
1196 *
1197 * @return std::map<std::string, boost::shared_ptr<MatrixDouble>>&
1198 */
1199 virtual std::map<std::string, boost::shared_ptr<MatrixDouble>> &getBBNMap();
1200
1201 /**
1202 * @brief get hash map of derivatives base function for BB base, key is a
1203 * field name
1204 *
1205 * @return std::map<std::string, boost::shared_ptr<MatrixDouble>>&
1206 */
1207 virtual std::map<std::string, boost::shared_ptr<MatrixDouble>> &
1208 getBBDiffNMap();
1209
1210 /**
1211 * @brief get ALpha indices for BB base by order
1212 *
1213 * @param o approximation order
1214 * @return boost::shared_ptr<MatrixInt>&
1215 */
1216 virtual boost::shared_ptr<MatrixInt> &
1217 getBBAlphaIndicesByOrderSharedPtr(const size_t o);
1218
1219 /**
1220 * @brief get BB base by order
1221 *
1222 * @param o
1223 * @return boost::shared_ptr<MatrixDouble>&
1224 */
1225 virtual boost::shared_ptr<MatrixDouble> &
1226 getBBNByOrderSharedPtr(const size_t o);
1227
1228 /**
1229 * @brief get BB base derivative by order
1230 *
1231 * @param o
1232 * @return boost::shared_ptr<MatrixDouble>&
1233 */
1234 virtual boost::shared_ptr<MatrixDouble> &
1235 getBBDiffNByOrderSharedPtr(const size_t o);
1236
1237 static constexpr size_t MaxBernsteinBezierOrder = BITFEID_SIZE;
1238
1239 virtual std::array<boost::shared_ptr<MatrixInt>, MaxBernsteinBezierOrder> &
1241
1242 virtual std::array<boost::shared_ptr<MatrixDouble>, MaxBernsteinBezierOrder> &
1244
1245 virtual std::array<boost::shared_ptr<MatrixDouble>, MaxBernsteinBezierOrder> &
1247
1248 /**
1249 * @brief Swap bases functions
1250 *
1251 * Some base are not hierarchical and depend on approximation order. Such case
1252 * demand special handling, that appropriate base order is set depending on
1253 * field, such that is accessible in operator.
1254 *
1255 * @note Base is not swap on meshsets
1256 *
1257 * @param field_name
1258 * @param base
1259 * @return MoFEMErrorCode
1260 */
1261 virtual MoFEMErrorCode baseSwap(const std::string &field_name,
1262 const FieldApproximationBase base);
1263
1264 /**@}*/
1265
1266 /** \name Broken spaces functions */
1267
1268protected:
1269 int sEnse; ///< Entity sense (orientation)
1270 ApproximationOrder oRder; ///< Entity order
1271 FieldSpace sPace; ///< Entity space
1272 FieldApproximationBase bAse; ///< Field approximation base
1273 VectorInt iNdices; ///< Global indices on entity
1274 VectorInt localIndices; ///< Local indices on entity
1275 VectorDofs dOfs; ///< DoFs on entity
1277 VectorDouble fieldData; ///< Field data on entity
1278 std::vector<BitRefLevel> entDataBitRefLevel; ///< Bit ref level in entity
1279
1280 std::vector<int> dofBrokenSideVec; ///< Map side to dofs number
1281 std::vector<EntityType>
1282 dofBrokenTypeVec; ///< Map type of entity to dof number
1283
1284 std::array<std::array<boost::shared_ptr<MatrixDouble>, LASTBASE>,
1287
1288 std::array<boost::shared_ptr<MatrixDouble>, LASTBASE> &N; ///< Base functions
1289 std::array<boost::shared_ptr<MatrixDouble>, LASTBASE>
1290 &diffN; ///< Derivatives of base functions
1291
1292 std::string bbFieldName; ///< field name
1293 VectorInt bbNodeOrder; ///< order of nodes
1294 std::map<std::string, boost::shared_ptr<MatrixDouble>> bbN;
1295 std::map<std::string, boost::shared_ptr<MatrixDouble>> bbDiffN;
1296 std::map<std::string, boost::shared_ptr<MatrixInt>>
1297 bbAlphaIndices; ///< Indices for Bernstein-Bezier (BB) base
1298
1299 std::array<boost::shared_ptr<MatrixDouble>, MaxBernsteinBezierOrder>
1300 bbNByOrder; ///< BB base functions by order
1301 std::array<boost::shared_ptr<MatrixDouble>, MaxBernsteinBezierOrder>
1302 bbDiffNByOrder; ///< BB base functions derivatives by order
1303 std::array<boost::shared_ptr<MatrixInt>, MaxBernsteinBezierOrder>
1304 bbAlphaIndicesByOrder; ///< BB alpha indices by order
1305
1306protected:
1307 /**
1308 * @brief Used by Bernstein base to keep temporary pointer
1309 *
1310 * @copydoc MoFEM::EntitiesFieldData::baseSwap
1311 */
1312 boost::shared_ptr<MatrixDouble> swapBaseNPtr;
1313
1314 /**
1315 * @brief Used by Bernstein base to keep temporary pointer
1316 *
1317 * @copydoc MoFEM::EntitiesFieldData::baseSwap
1318 */
1319 boost::shared_ptr<MatrixDouble> swapBaseDiffNPtr;
1320
1321 friend struct OpAddParentEntData;
1322
1323 template <typename OpBase> friend struct OpGetBrokenBaseSideData;
1324};
1325
1327
1328/** \brief Derived data on single entity (This is passed as argument to
1329 * DataOperator::doWork) \ingroup mofem_forces_and_sources_user_data_operators
1330 * \nosubgrouping
1331 *
1332 * DerivedEntData share part information with EntData except information about
1333 * base functions.
1334 *
1335 */
1338
1340 const boost::shared_ptr<EntitiesFieldData::EntData> &ent_data_ptr);
1341
1342 int getSense() const;
1343
1344 //// \brief get entity bit ref level
1345 std::vector<BitRefLevel> &getEntDataBitRefLevel();
1346
1347 boost::shared_ptr<MatrixDouble> &
1349 const BaseDerivatives derivative);
1350
1351 boost::shared_ptr<MatrixDouble> &
1353
1354 boost::shared_ptr<MatrixDouble> &
1356
1357 const boost::shared_ptr<MatrixDouble> &
1358 getNSharedPtr(const FieldApproximationBase base) const;
1359
1360 const boost::shared_ptr<MatrixDouble> &
1361 getDiffNSharedPtr(const FieldApproximationBase base) const;
1362
1363 inline boost::shared_ptr<MatrixDouble> &
1365
1366 inline boost::shared_ptr<MatrixDouble> &
1368
1369 boost::shared_ptr<MatrixInt> &
1370 getBBAlphaIndicesSharedPtr(const std::string &field_name);
1371
1372 /**
1373 * Get shared pointer to BB base base functions
1374 */
1375 boost::shared_ptr<MatrixDouble> &
1376 getBBNSharedPtr(const std::string &field_name);
1377
1378 /**
1379 * Get shared pointer to BB base base functions
1380 */
1381 const boost::shared_ptr<MatrixDouble> &
1382 getBBNSharedPtr(const std::string &field_name) const;
1383
1384 /**
1385 * Get shared pointer to BB derivatives of base base functions
1386 */
1387 boost::shared_ptr<MatrixDouble> &
1388 getBBDiffNSharedPtr(const std::string &field_name);
1389
1390 /**
1391 * Get shared pointer to derivatives of BB base base functions
1392 */
1393 const boost::shared_ptr<MatrixDouble> &
1394 getBBDiffNSharedPtr(const std::string &field_name) const;
1395
1396 /**
1397 * @copydoc MoFEM::EntitiesFieldData::EntData::swapBaseNPtr
1398 */
1399 MoFEMErrorCode baseSwap(const std::string &field_name,
1400 const FieldApproximationBase base);
1401
1402 /** \name Broken spaces functions */
1403
1404 /**@{*/
1405
1406protected:
1407 const boost::shared_ptr<EntitiesFieldData::EntData> entDataPtr;
1408};
1409
1413
1415 return iNdices;
1416}
1417
1418const VectorIntAdaptor
1420 unsigned int size = 0;
1421 if (auto dof = dOfs[0]) {
1422 size = dof->getOrderNbDofs(order) * dof->getNbOfCoeffs();
1423 size = size < iNdices.size() ? size : iNdices.size();
1424 }
1425 int *data = &*iNdices.data().begin();
1426 return VectorIntAdaptor(size, ublas::shallow_array_adaptor<int>(size, data));
1427}
1428
1430 return localIndices;
1431}
1432
1433const VectorIntAdaptor
1435 unsigned int size = 0;
1436 if (auto dof = dOfs[0]) {
1437 size = dof->getOrderNbDofs(order) * dof->getNbOfCoeffs();
1438 size = size < localIndices.size() ? size : localIndices.size();
1439 }
1440 int *data = &*localIndices.data().begin();
1441 return VectorIntAdaptor(size, ublas::shallow_array_adaptor<int>(size, data));
1442}
1443
1445
1447
1449
1451 return localIndices;
1452}
1453
1455 return fieldData;
1456}
1457
1458const VectorAdaptor
1460 unsigned int size = 0;
1461 if (auto dof = dOfs[0]) {
1462 size = dof->getOrderNbDofs(order) * dof->getNbOfCoeffs();
1463 size = size < fieldData.size() ? size : fieldData.size();
1464 }
1465 double *data = &*fieldData.data().begin();
1466 return getVectorAdaptor(data, size);
1467}
1468
1470 return dOfs;
1471}
1472
1474
1476
1480
1481const VectorFieldEntities &
1483 return fieldEntities;
1484}
1485
1486template <int Tensor_Dim>
1489 std::stringstream s;
1490 s << "Not implemented for this dimension dim = " << Tensor_Dim;
1491 THROW_MESSAGE(s.str());
1492}
1493
1494template <int Tensor_Dim0, int Tensor_Dim1>
1496 Tensor_Dim0, Tensor_Dim1>
1498 std::stringstream s;
1499 s << "Not implemented for this dimension dim0 = " << Tensor_Dim0;
1500 s << " and dim1 " << Tensor_Dim1;
1501 THROW_MESSAGE(s.str());
1502}
1503
1504template <int Tensor_Dim>
1506 FTensor::PackPtr<double *, (Tensor_Dim * (Tensor_Dim + 1)) / 2>, Tensor_Dim>
1508 std::stringstream s;
1509 s << "Not implemented for this dimension dim = " << Tensor_Dim;
1510 THROW_MESSAGE(s.str());
1511}
1512
1514
1516
1519#ifndef NDEBUG
1520 if (!getNSharedPtr(base)) {
1521 MOFEM_LOG_C("SELF", Sev::error, "Ptr to base %s functions is null",
1524 "Null pointer to base functions");
1525 }
1526#endif // NDEBUG
1527 return *(getNSharedPtr(base));
1528}
1529
1531#ifndef NDEBUG
1532 if (!getBBNSharedPtr(field_name)) {
1533 MOFEM_LOG_C("SELF", Sev::error, "Ptr to field %s base functions is null ",
1534 field_name.c_str());
1536 "Null pointer to base functions");
1537 }
1538#endif // NDEBUG
1539 return *(getBBNSharedPtr(field_name));
1540}
1541
1543
1546 return *(getDiffNSharedPtr(base));
1547}
1548
1551 const BaseDerivatives derivative) {
1552#ifndef NDEBUG
1553 if (!getNSharedPtr(base, derivative)) {
1554 MOFEM_LOG_C("SELF", Sev::error,
1555 "Ptr to base %s functions derivative %d is null",
1556 ApproximationBaseNames[base], derivative);
1557 THROW_MESSAGE("Null pointer");
1558 }
1559#endif
1560 return *(getNSharedPtr(base, derivative));
1561}
1562
1565 #ifndef NDEBUG
1566 if (!getBBDiffNSharedPtr(field_name)) {
1568 "Null pointer to base functions derivative");
1569 }
1570 #endif // NDEBUG
1571 return *(getBBDiffNSharedPtr(field_name));
1572}
1573
1575
1578 return getN(bAse, derivative);
1579}
1580
1581const VectorAdaptor
1583 const int gg) {
1584 int size = getN(base).size2();
1585 double *data = &getN(base)(gg, 0);
1586 return VectorAdaptor(size, ublas::shallow_array_adaptor<double>(size, data));
1587}
1588
1590 return getN(bAse, gg);
1591}
1592
1593const MatrixAdaptor
1595 const int gg) {
1596 // FIXME: That is bug, it will not work if number of integration pts is
1597 // equal to number of nodes on entity. User who not implementing low
1598 // level DataOperator will not experience this.
1599 if (getN(base).size1() == getDiffN(base).size1()) {
1600 int size = getN(base).size2();
1601 int dim = getDiffN(base).size2() / size;
1602 double *data = &getDiffN(base)(gg, 0);
1603 return MatrixAdaptor(
1604 getN(base).size2(), dim,
1605 ublas::shallow_array_adaptor<double>(getDiffN(base).size2(), data));
1606 } else {
1607 // in some cases, f.e. for derivatives of nodal base functions at only
1608 // one gauss point is needed
1609 return MatrixAdaptor(
1610 getN(base).size1(), getN(base).size2(),
1611 ublas::shallow_array_adaptor<double>(getDiffN(base).data().size(),
1612 &getDiffN(base).data()[0]));
1613 }
1614}
1615
1617 return getDiffN(bAse, gg);
1618}
1619
1620const VectorAdaptor
1622 const int gg, const int nb_base_functions) {
1623 (void)getN()(gg, nb_base_functions -
1624 1); // throw error if nb_base_functions is to big
1625 double *data = &getN(base)(gg, 0);
1626 return VectorAdaptor(nb_base_functions, ublas::shallow_array_adaptor<double>(
1627 nb_base_functions, data));
1628}
1629
1630const VectorAdaptor
1631EntitiesFieldData::EntData::getN(const int gg, const int nb_base_functions) {
1632 return getN(bAse, gg, nb_base_functions);
1633}
1634
1635const MatrixAdaptor
1637 const int gg,
1638 const int nb_base_functions) {
1639 // FIXME: That is bug, it will not work if number of integration pts is
1640 // equal to number of nodes on entity. User who not implementing low
1641 // level DataOperator will not experience this.
1642 if (getN(base).size1() == getDiffN(base).size1()) {
1643 (void)getN(base)(gg,
1644 nb_base_functions -
1645 1); // throw error if nb_base_functions is to big
1646 int dim = getDiffN(base).size2() / getN(base).size2();
1647 double *data = &getDiffN(base)(gg, 0);
1648 return MatrixAdaptor(
1649 nb_base_functions, dim,
1650 ublas::shallow_array_adaptor<double>(dim * nb_base_functions, data));
1651 } else {
1652 // in some cases, f.e. for derivatives of nodal base functions only one
1653 // gauss point is needed
1654 return MatrixAdaptor(
1655 getN(base).size1(), getN(base).size2(),
1656 ublas::shallow_array_adaptor<double>(getDiffN(base).data().size(),
1657 &getDiffN(base).data()[0]));
1658 }
1659}
1660
1661const MatrixAdaptor
1663 const int nb_base_functions) {
1664 return getDiffN(bAse, gg, nb_base_functions);
1665}
1666
1667template <int DIM>
1668const MatrixAdaptor
1670 const int gg) {
1671 if (PetscUnlikely(getN(base).size2() % DIM)) {
1672 THROW_MESSAGE("Wrong dimension");
1673 }
1674
1675 const int nb_base_functions = getN(base).size2() / DIM;
1676 double *data = &getN(base)(gg, 0);
1677 return MatrixAdaptor(
1678 nb_base_functions, DIM,
1679 ublas::shallow_array_adaptor<double>(DIM * nb_base_functions, data));
1680}
1681
1682template <int DIM>
1684 return getVectorN<DIM>(bAse, gg);
1685}
1686
1687template <int DIM0, int DIM1>
1688const MatrixAdaptor
1690 const int gg) {
1691 if (PetscUnlikely(getDiffN(base).size2() % (DIM0 * DIM1))) {
1692 THROW_MESSAGE("Wrong dimension");
1693 }
1694
1695 const int nb_base_functions = getN(base).size2() / (DIM0 * DIM1);
1696 double *data = &getN(base)(gg, 0);
1697 return MatrixAdaptor(nb_base_functions, DIM0 * DIM1,
1698 ublas::shallow_array_adaptor<double>(
1699 DIM0 * DIM1 * nb_base_functions, data));
1700}
1701
1702template <int DIM0, int DIM1>
1704 return getVectorDiffN<DIM0, DIM1>(bAse, gg);
1705}
1706
1707template <int DIM0, int DIM1>
1708const MatrixAdaptor
1710 const int dof, const int gg) {
1711 double *data =
1712 &EntitiesFieldData::EntData::getDiffN(base)(gg, DIM0 * DIM1 * dof);
1713 return MatrixAdaptor(DIM0, DIM1,
1714 ublas::shallow_array_adaptor<double>(DIM0 * DIM1, data));
1715}
1716
1717template <int DIM0, int DIM1>
1719 const int gg) {
1720 return getVectorDiffN<DIM0, DIM1>(bAse, dof, gg);
1721}
1722
1725 double *ptr = &*getN(base).data().begin();
1727};
1728
1731 return getFTensor0N(bAse);
1732};
1733
1736 const int gg, const int bb) {
1737 double *ptr = &getN(base)(gg, bb);
1739};
1740
1742EntitiesFieldData::EntData::getFTensor0N(const int gg, const int bb) {
1743 return getFTensor0N(bAse, gg, bb);
1744};
1745
1746template <int Tensor_Dim> auto EntitiesFieldData::EntData::getFTensor1N() {
1747 return getFTensor1N<Tensor_Dim>(bAse);
1748}
1749
1750template <int Tensor_Dim>
1751auto EntitiesFieldData::EntData::getFTensor1N(const int gg, const int bb) {
1752 return getFTensor1N<Tensor_Dim>(bAse, gg, bb);
1753}
1754
1755template <int Tensor_Dim0, int Tensor_Dim1>
1757 return getFTensor2N<Tensor_Dim0, Tensor_Dim1>(bAse);
1758}
1759
1760template <int Tensor_Dim0, int Tensor_Dim1>
1761auto EntitiesFieldData::EntData::getFTensor2N(const int gg, const int bb) {
1762 return getFTensor2N<Tensor_Dim0, Tensor_Dim1>(bAse, gg, bb);
1763}
1764
1765/** \name Bernstein-Bezier base only functions */
1766
1767/**@{*/
1768
1770
1772 return *getBBAlphaIndicesSharedPtr(bbFieldName);
1773}
1774
1775/**@}*/
1776
1777/** \name DerivedEntData */
1778
1779/**@{*/
1780
1781boost::shared_ptr<MatrixDouble> &
1786
1787boost::shared_ptr<MatrixDouble> &
1792
1793/**@}*/
1794
1795/**
1796 * @brief Assemble PETSc vector
1797 *
1798 * Function extract indices from entity data and assemble vector
1799 *
1800 * <a
1801 * href=https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecSetValues.html>See
1802 * PETSc documentation</a>
1803 *
1804 * @param V
1805 * @param data
1806 * @param ptr
1807 * @param iora
1808 * @return MoFEMErrorCode
1809 */
1810template <typename T = EntityStorage>
1812 const EntitiesFieldData::EntData &data,
1813 const double *ptr, InsertMode iora) {
1814 static_assert(!std::is_same<T, T>::value,
1815 "VecSetValues value for this data storage is not implemented");
1816 return MOFEM_NOT_IMPLEMENTED;
1817}
1818
1819template <>
1822 const double *ptr, InsertMode iora) {
1823 return VecSetValues(V, data.getIndices().size(), &*data.getIndices().begin(),
1824 ptr, iora);
1825}
1826
1827/**
1828 * @brief Assemble PETSc vector
1829 *
1830 * Function extract indices from entity data and assemble vector
1831 *
1832 * <a
1833 * href=https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Vec/VecSetValues.html>See
1834 * PETSc documentation</a>
1835 *
1836 * @param V
1837 * @param data
1838 * @param vec
1839 * @param iora
1840 * @return MoFEMErrorCode
1841 */
1842template <typename T = EntityStorage>
1844 const EntitiesFieldData::EntData &data,
1845 const VectorDouble &vec, InsertMode iora) {
1846 return VecSetValues<T>(V, data, &*vec.data().begin(), iora);
1847}
1848
1849/**
1850 * @brief Assemble PETSc matrix
1851 *
1852 * Function extract indices from entity data and assemble vector
1853 *
1854 * <a
1855 * href=https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatSetValues.html>See
1856 * PETSc documentation</a>
1857 *
1858 * @param M
1859 * @param row_data
1860 * @param col_data
1861 * @param ptr
1862 * @param iora
1863 * @return MoFEMErrorCode
1864 */
1865template <typename T = EntityStorage>
1867 const EntitiesFieldData::EntData &row_data,
1868 const EntitiesFieldData::EntData &col_data,
1869 const double *ptr, InsertMode iora) {
1870 static_assert(!std::is_same<T, T>::value,
1871 "MatSetValues value for this data storage is not implemented");
1872 return MOFEM_NOT_IMPLEMENTED;
1873}
1874
1875/**
1876 * @brief Assemble PETSc matrix
1877 *
1878 * Function extract indices from entity data and assemble vector
1879 *
1880 * <a
1881 * href=https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatSetValues.html>See
1882 * PETSc documentation</a>
1883 *
1884 * @param M
1885 * @param row_data
1886 * @param col_data
1887 * @param mat
1888 * @param iora
1889 * @return MoFEMErrorCode
1890 */
1891template <typename T = EntityStorage>
1893 const EntitiesFieldData::EntData &row_data,
1894 const EntitiesFieldData::EntData &col_data,
1895 const MatrixDouble &mat, InsertMode iora) {
1896 return MatSetValues<T>(M, row_data, col_data, &*mat.data().begin(), iora);
1897}
1898
1899template <>
1902 const EntitiesFieldData::EntData &col_data,
1903 const double *ptr, InsertMode iora) {
1904 return MatSetValues(
1905 M, row_data.getIndices().size(), &*row_data.getIndices().begin(),
1906 col_data.getIndices().size(), &*col_data.getIndices().begin(), ptr, iora);
1907}
1908
1909/** \name Specializations for tensor base function */
1910
1911/**@{*/
1912
1913template <>
1915EntitiesFieldData::EntData::getFTensor1N<3>(FieldApproximationBase base);
1916
1917template <>
1919EntitiesFieldData::EntData::getFTensor1N<3>(FieldApproximationBase base,
1920 const int gg, const int bb);
1921
1922template <>
1924EntitiesFieldData::EntData::getFTensor2N<3, 3>(FieldApproximationBase base);
1925
1926/**@}*/
1927
1928/** \name Specializations for derivatives of base functions */
1929
1930/**@{*/
1931
1932template <>
1934EntitiesFieldData::EntData::getFTensor1DiffN<3>(
1935 const FieldApproximationBase base);
1936template <>
1938EntitiesFieldData::EntData::getFTensor1DiffN<3>();
1939
1940template <>
1942EntitiesFieldData::EntData::getFTensor1DiffN<2>(
1943 const FieldApproximationBase base);
1944template <>
1946EntitiesFieldData::EntData::getFTensor1DiffN<2>();
1947
1948template <>
1950EntitiesFieldData::EntData::getFTensor2DiffN<3, 3>(FieldApproximationBase base);
1951template <>
1953EntitiesFieldData::EntData::getFTensor2DiffN<3, 3>(FieldApproximationBase base,
1954 const int gg, const int bb);
1955template <>
1957EntitiesFieldData::EntData::getFTensor2DiffN<3, 2>(FieldApproximationBase base);
1958template <>
1960EntitiesFieldData::EntData::getFTensor2DiffN<3, 2>(FieldApproximationBase base,
1961 const int gg, const int bb);
1962
1963template <>
1966
1967template <>
1970
1971template <>
1974 const int gg, const int bb);
1975
1976template <>
1978EntitiesFieldData::EntData::getFTensor2DiffN2<2>(
1979 const FieldApproximationBase base);
1980
1981template <>
1983EntitiesFieldData::EntData::getFTensor2DiffN2<2>(
1984 const FieldApproximationBase base, const int gg, const int bb);
1985
1986template <>
1988EntitiesFieldData::EntData::getFTensor2DiffN2<2>();
1989
1990template <>
1992EntitiesFieldData::EntData::getFTensor2DiffN2<2>(const int gg, const int bb);
1993
1994template <>
1996EntitiesFieldData::EntData::getFTensor2DiffN2<3>(
1997 const FieldApproximationBase base);
1998
1999template <>
2001EntitiesFieldData::EntData::getFTensor2DiffN2<3>(
2002 const FieldApproximationBase base, const int gg, const int bb);
2003
2004template <>
2006EntitiesFieldData::EntData::getFTensor2DiffN2<3>();
2007
2008template <>
2010EntitiesFieldData::EntData::getFTensor2DiffN2<3>(const int gg, const int bb);
2011
2012/**@}*/
2013
2014/** \name Specializations for field data */
2015
2016/**@{*/
2017
2018template <>
2020EntitiesFieldData::EntData::getFTensor1FieldData<3>();
2021
2022template <>
2024EntitiesFieldData::EntData::getFTensor1FieldData<2>();
2025
2026template <>
2028EntitiesFieldData::EntData::getFTensor1FieldData<1>();
2029
2030template <>
2032EntitiesFieldData::EntData::getFTensor2FieldData<1, 1>();
2033
2034template <>
2036EntitiesFieldData::EntData::getFTensor2FieldData<1, 2>();
2037
2038template <>
2040EntitiesFieldData::EntData::getFTensor2FieldData<1, 3>();
2041
2042template <>
2044EntitiesFieldData::EntData::getFTensor2FieldData<2, 2>();
2045
2046template <>
2048EntitiesFieldData::EntData::getFTensor2FieldData<3, 3>();
2049
2050template <>
2052EntitiesFieldData::EntData::getFTensor2SymmetricFieldData<3>();
2053
2054template <>
2056EntitiesFieldData::EntData::getFTensor2SymmetricFieldData<2>();
2057
2058/**@}*/
2059
2060/**
2061 * @deprecated Use EntitiesFieldData
2062 */
2064
2065/**
2066 * @deprecated Use DerivedEntitiesFieldData
2067 */
2069
2070} // namespace MoFEM
2071
2072#endif //__ENTITIES_FIELD_DATA_HPP__
2073
2074/**
2075 * \defgroup mofem_forces_and_sources_user_data_operators User data operator
2076 * data structures \ingroup
2077 *
2078 * \brief Users data structures and operator
2079 *
2080 * Data structures passed by argument to MoFEM::DataOperator::doWork and generic
2081 * user operators operating on those structures.
2082 *
2083 */
#define MOFEM_LOG_C(channel, severity, format,...)
FieldApproximationBase
approximation base
Definition definitions.h:58
@ LASTBASE
Definition definitions.h:69
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
FieldSpace
approximation spaces
Definition definitions.h:82
@ LASTSPACE
FieldSpace in [ 0, LASTSPACE )
Definition definitions.h:89
@ MOFEM_DATA_INCONSISTENCY
Definition definitions.h:31
@ MOFEM_NOT_IMPLEMENTED
Definition definitions.h:32
static const char *const ApproximationBaseNames[]
Definition definitions.h:72
#define DEPRECATED
Definition definitions.h:17
#define BITFEID_SIZE
max number of finite elements
#define THROW_MESSAGE(msg)
Throw MoFEM exception.
constexpr int order
constexpr int DIM1
Definition level_set.cpp:21
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
VectorShallowArrayAdaptor< int > VectorIntAdaptor
Definition Types.hpp:116
int ApproximationOrder
Approximation on the entity.
Definition Types.hpp:26
VectorShallowArrayAdaptor< double > VectorAdaptor
Definition Types.hpp:115
UBlasMatrix< double > MatrixDouble
Definition Types.hpp:77
UBlasVector< double > VectorDouble
Definition Types.hpp:68
MatrixShallowArrayAdaptor< double > MatrixAdaptor
Matrix adaptor.
Definition Types.hpp:132
UBlasMatrix< int > MatrixInt
Definition Types.hpp:76
UBlasVector< int > VectorInt
Definition Types.hpp:67
implementation of Data Operators for Forces and Sources
Definition Common.hpp:10
auto getVectorAdaptor(T1 ptr, const size_t n)
Get Vector adaptor.
Definition Templates.hpp:31
DEPRECATED typedef DerivedEntitiesFieldData DerivedDataForcesAndSourcesCore
ublas::vector< FieldEntity *, FieldEntAllocator > VectorFieldEntities
DEPRECATED typedef EntitiesFieldData DataForcesAndSourcesCore
MoFEMErrorCode MatSetValues(Mat M, const EntitiesFieldData::EntData &row_data, const EntitiesFieldData::EntData &col_data, const double *ptr, InsertMode iora)
Assemble PETSc matrix.
ublas::unbounded_array< FieldEntity *, std::allocator< FieldEntity * > > FieldEntAllocator
MoFEMErrorCode VecSetValues< EntityStorage >(Vec V, const EntitiesFieldData::EntData &data, const double *ptr, InsertMode iora)
MoFEMErrorCode MatSetValues< EntityStorage >(Mat M, const EntitiesFieldData::EntData &row_data, const EntitiesFieldData::EntData &col_data, const double *ptr, InsertMode iora)
ublas::unbounded_array< FEDofEntity *, std::allocator< FEDofEntity * > > DofsAllocator
MoFEMErrorCode VecSetValues(Vec V, const EntitiesFieldData::EntData &data, const double *ptr, InsertMode iora)
Assemble PETSc vector.
ublas::vector< FEDofEntity *, DofsAllocator > VectorDofs
constexpr auto field_name
const int N
Definition speed_test.cpp:3
Derived data on single entity (This is passed as argument to DataOperator::doWork)
boost::shared_ptr< MatrixDouble > & getBBNSharedPtr(const std::string &field_name)
boost::shared_ptr< MatrixDouble > & getDerivedDiffNSharedPtr(const FieldApproximationBase base)
const boost::shared_ptr< EntitiesFieldData::EntData > entDataPtr
MoFEMErrorCode baseSwap(const std::string &field_name, const FieldApproximationBase base)
Used by Bernstein base to keep temporary pointer.
std::vector< BitRefLevel > & getEntDataBitRefLevel()
Get entity bit reference level.
int getSense() const
Get entity sense for conforming approximation fields.
boost::shared_ptr< MatrixInt > & getBBAlphaIndicesSharedPtr(const std::string &field_name)
boost::shared_ptr< MatrixDouble > & getBBDiffNSharedPtr(const std::string &field_name)
boost::shared_ptr< MatrixDouble > & getDiffNSharedPtr(const FieldApproximationBase base)
Get shared pointer to derivatives of base functions.
boost::shared_ptr< MatrixDouble > & getNSharedPtr(const FieldApproximationBase base, const BaseDerivatives derivative)
Get shared pointer to base functions with derivatives.
boost::shared_ptr< MatrixDouble > & getDerivedNSharedPtr(const FieldApproximationBase base)
this class derives data from other data structure
MoFEMErrorCode setElementType(const EntityType type)
Set element type for derived data.
const boost::shared_ptr< EntitiesFieldData > dataPtr
Data on single entity (This is passed as argument to DataOperator::doWork)
std::vector< int > dofBrokenSideVec
Map side to dofs number.
FTensor::Tensor2< FTensor::PackPtr< double *, Tensor_Dim *Tensor_Dim >, Tensor_Dim, Tensor_Dim > getFTensor2DiffN2(const FieldApproximationBase base)
Get second derivatives of scalar base functions.
std::array< boost::shared_ptr< MatrixDouble >, MaxBernsteinBezierOrder > bbDiffNByOrder
BB base functions derivatives by order.
boost::shared_ptr< MatrixDouble > swapBaseDiffNPtr
Used by Bernstein base to keep temporary pointer.
virtual std::array< boost::shared_ptr< MatrixDouble >, MaxBernsteinBezierOrder > & getBBDiffNByOrderArray()
const VectorIntAdaptor getIndicesUpToOrder(int order)
Get global indices of DOFs on entity up to given order.
virtual std::array< boost::shared_ptr< MatrixDouble >, MaxBernsteinBezierOrder > & getBBNByOrderArray()
int sEnse
Entity sense (orientation)
auto getFTensor2N()
Get base functions for Hdiv space.
VectorDouble fieldData
Field data on entity.
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1DiffN(const int gg, const int bb)
Get derivatives of base functions (Loop by integration points)
MatrixDouble & getN()
Get base functions using default base.
FieldApproximationBase & getBase()
Get approximation base.
boost::shared_ptr< MatrixDouble > swapBaseNPtr
Used by Bernstein base to keep temporary pointer.
ApproximationOrder getOrder() const
Get approximation order.
std::map< std::string, boost::shared_ptr< MatrixInt > > bbAlphaIndices
Indices for Bernstein-Bezier (BB) base.
FTensor::Tensor3< FTensor::PackPtr< double *, Tensor_Dim0 *Tensor_Dim1 *Tensor_Dim2 >, Tensor_Dim0, Tensor_Dim1, Tensor_Dim2 > getFTensor3DiffN(const int gg, const int bb)
Get derivatives of base functions for tonsorial Hdiv space at integration pts.
friend std::ostream & operator<<(std::ostream &os, const EntitiesFieldData::EntData &e)
VectorInt localIndices
Local indices on entity.
const MatrixAdaptor getVectorDiffN(FieldApproximationBase base, const int gg)
get DiffHdiv of base functions at Gauss pts
const VectorFieldEntities & getFieldEntities() const
Get field entities (const version)
virtual std::array< boost::shared_ptr< MatrixInt >, MaxBernsteinBezierOrder > & getBBAlphaIndicesByOrderArray()
virtual std::map< std::string, boost::shared_ptr< MatrixInt > > & getBBAlphaIndicesMap()
const VectorIntAdaptor getLocalIndicesUpToOrder(int order)
Get local indices of DOFs on entity up to given order.
const VectorAdaptor getFieldDataUpToOrder(int order)
Get DOF values on entity up to given order.
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1DiffN()
Get derivatives of base functions.
std::vector< BitRefLevel > entDataBitRefLevel
Bit ref level in entity.
FTensor::Tensor2< FTensor::PackPtr< double *, Tensor_Dim0 *Tensor_Dim1 >, Tensor_Dim0, Tensor_Dim1 > getFTensor2DiffN(const int gg, const int bb)
Get derivatives of base functions for Hdiv space at integration pts.
static constexpr size_t MaxBernsteinBezierOrder
VectorInt iNdices
Global indices on entity.
virtual boost::shared_ptr< MatrixDouble > & getNSharedPtr(const FieldApproximationBase base, const BaseDerivatives derivative)
Get shared pointer to base functions with derivatives.
virtual boost::shared_ptr< MatrixInt > & getBBAlphaIndicesSharedPtr(const std::string &field_name)
const VectorDouble & getFieldData() const
Get DOF values on entity.
virtual boost::shared_ptr< MatrixDouble > & getBBDiffNByOrderSharedPtr(const size_t o)
get BB base derivative by order
FTensor::Tensor3< FTensor::PackPtr< double *, Tensor_Dim0 *Tensor_Dim1 *Tensor_Dim2 >, Tensor_Dim0, Tensor_Dim1, Tensor_Dim2 > getFTensor3Diff2N(FieldApproximationBase base)
Get second derivatives of base functions for Hvec space.
virtual std::map< std::string, boost::shared_ptr< MatrixDouble > > & getBBNMap()
get hash map of base function for BB base, key is a field name
const MatrixAdaptor getVectorN(const FieldApproximationBase base, const int gg)
get Hdiv of base functions at Gauss pts
virtual MoFEMErrorCode baseSwap(const std::string &field_name, const FieldApproximationBase base)
Swap bases functions.
const VectorInt & getLocalIndices() const
Get local indices of degrees of freedom on entity.
ApproximationOrder oRder
Entity order.
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1FieldData()
Return FTensor of rank 1, i.e. vector from field data coefficients.
std::vector< EntityType > dofBrokenTypeVec
Map type of entity to dof number.
FTensor::Tensor3< FTensor::PackPtr< double *, Tensor_Dim0 *Tensor_Dim1 *Tensor_Dim2 >, Tensor_Dim0, Tensor_Dim1, Tensor_Dim2 > getFTensor3DiffN()
Get derivatives of base functions for tonsorial Hdiv space.
virtual std::map< std::string, boost::shared_ptr< MatrixDouble > > & getBBDiffNMap()
get hash map of derivatives base function for BB base, key is a field name
virtual std::vector< BitRefLevel > & getEntDataBitRefLevel()
Get entity bit reference level.
virtual int getSense() const
Get entity sense for conforming approximation fields.
VectorFieldEntities fieldEntities
Field entities.
FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > getFTensor0N()
Get base function as Tensor0.
FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > getFTensor0FieldData()
Return scalar files as a FTensor of rank 0.
MatrixInt & getBBAlphaIndices()
Get file BB indices.
std::array< boost::shared_ptr< MatrixInt >, MaxBernsteinBezierOrder > bbAlphaIndicesByOrder
BB alpha indices by order.
MatrixDouble & getDiffN()
Get derivatives of base functions using default base.
VectorInt & getBBNodeOrder()
Get orders at the nodes.
virtual boost::shared_ptr< MatrixInt > & getBBAlphaIndicesByOrderSharedPtr(const size_t o)
get ALpha indices for BB base by order
std::array< boost::shared_ptr< MatrixDouble >, LASTBASE > & diffN
Derivatives of base functions.
std::map< std::string, boost::shared_ptr< MatrixDouble > > bbDiffN
const VectorDofs & getFieldDofs() const
Get DOF data structures (const version)
FTensor::Tensor2_symmetric< FTensor::PackPtr< double *,(Tensor_Dim *(Tensor_Dim+1))/2 >, Tensor_Dim > getFTensor2SymmetricFieldData()
Return symmetric FTensor rank 2, i.e. matrix from field data coefficients.
FTensor::Tensor3< FTensor::PackPtr< double *, Tensor_Dim0 *Tensor_Dim1 *Tensor_Dim2 >, Tensor_Dim0, Tensor_Dim1, Tensor_Dim2 > getFTensor3Diff2N()
Get second derivatives of base functions for Hvec space.
FieldSpace & getSpace()
Get field space.
auto getFTensor1N()
Get base functions for Hdiv space.
virtual boost::shared_ptr< MatrixDouble > & getBBNSharedPtr(const std::string &field_name)
std::array< boost::shared_ptr< MatrixDouble >, MaxBernsteinBezierOrder > bbNByOrder
BB base functions by order.
std::array< std::array< boost::shared_ptr< MatrixDouble >, LASTBASE >, LastDerivative > baseFunctionsAndBaseDerivatives
FTensor::Tensor2< FTensor::PackPtr< double *, Tensor_Dim0 *Tensor_Dim1 >, Tensor_Dim0, Tensor_Dim1 > getFTensor2DiffN()
Get derivatives of base functions for Hdiv space.
std::map< std::string, boost::shared_ptr< MatrixDouble > > bbN
virtual boost::shared_ptr< MatrixDouble > & getDiffNSharedPtr(const FieldApproximationBase base)
Get shared pointer to derivatives of base functions.
FTensor::Tensor2< FTensor::PackPtr< double *, Tensor_Dim *Tensor_Dim >, Tensor_Dim, Tensor_Dim > getFTensor2DiffN2()
Get second derivatives of scalar base functions.
FTensor::Tensor2< FTensor::PackPtr< double *, Tensor_Dim0 *Tensor_Dim1 >, Tensor_Dim0, Tensor_Dim1 > getFTensor2FieldData()
Return FTensor rank 2, i.e. matrix from field data coefficients.
virtual boost::shared_ptr< MatrixDouble > & getBBNByOrderSharedPtr(const size_t o)
get BB base by order
FTensor::Tensor2< FTensor::PackPtr< double *, Tensor_Dim *Tensor_Dim >, Tensor_Dim, Tensor_Dim > getFTensor2DiffN2(const int gg, const int bb)
Get second derivatives of scalar base functions at integration pts.
const VectorInt & getIndices() const
Get global indices of degrees of freedom on entity.
FieldApproximationBase bAse
Field approximation base.
FTensor::Tensor2< FTensor::PackPtr< double *, Tensor_Dim *Tensor_Dim >, Tensor_Dim, Tensor_Dim > getFTensor2DiffN2(const FieldApproximationBase base, const int gg, const int bb)
Get second derivatives of scalar base functions at integration pts.
virtual boost::shared_ptr< MatrixDouble > & getBBDiffNSharedPtr(const std::string &field_name)
std::array< boost::shared_ptr< MatrixDouble >, LASTBASE > & N
Base functions.
data structure for finite element entity
friend std::ostream & operator<<(std::ostream &os, const EntitiesFieldData &e)
MoFEMErrorCode resetFieldDependentData()
Reset data associated with particular field name.
MatrixInt facesNodes
nodes on finite element faces
std::array< std::bitset< LASTSPACE >, MBMAXTYPE > spacesOnEntities
spaces on entity types
std::bitset< LASTBASE > bAse
bases on element
std::array< boost::ptr_vector< EntData >, MBMAXTYPE > dataOnEntities
std::array< std::bitset< LASTBASE >, LASTSPACE > basesOnSpaces
base on spaces
virtual ~EntitiesFieldData()=default
virtual MoFEMErrorCode baseSwap(const std::string &field_name, const FieldApproximationBase base)
Swap approximation base.
virtual MoFEMErrorCode setElementType(const EntityType type)
Set element type and initialize data structures.
std::array< std::bitset< LASTBASE >, LASTSPACE > brokenBasesOnSpaces
base on spaces
MatrixInt facesNodesOrder
order of face nodes on element
std::array< std::bitset< LASTBASE >, MBMAXTYPE > basesOnEntities
bases on entity types
keeps information about indexed dofs for the finite element
Struct keeps handle to entity in the field.
Operator to project base functions from parent entity to child.