v0.15.0
Loading...
Searching...
No Matches
UserDataOperators.hpp
Go to the documentation of this file.
1/** \file UserDataOperators.hpp
2 * \brief User data Operators
3
4*/
5
6#ifndef __USER_DATA_OPERATORS_HPP__
7 #define __USER_DATA_OPERATORS_HPP__
8
9namespace MoFEM {
10
11/** \name Get values at Gauss pts */
12
13/**@{*/
14
15/** \name Scalar values */
16
17/**@{*/
18
19/** \brief Scalar field values at integration points
20 *
21 */
22template <class T, class A>
25
27 const std::string field_name,
28 boost::shared_ptr<ublas::vector<T, A>> data_ptr,
29 const EntityType zero_type = MBVERTEX)
31 dataPtr(data_ptr), zeroType(zero_type) {
32 if (!dataPtr)
33 THROW_MESSAGE("Pointer is not set");
34 }
35
37 const std::string field_name,
38 boost::shared_ptr<ublas::vector<T, A>> data_ptr,
39 SmartPetscObj<Vec> data_vec, const EntityType zero_type = MBVERTEX)
41 dataPtr(data_ptr), zeroType(zero_type), dataVec(data_vec) {
42 if (!dataPtr)
43 THROW_MESSAGE("Pointer is not set");
44 }
45
46 /**
47 * \brief calculate values of scalar field at integration points
48 * @param side side entity number
49 * @param type side entity type
50 * @param data entity data
51 * @return error code
52 */
53 MoFEMErrorCode doWork(int side, EntityType type,
55
56protected:
57 boost::shared_ptr<ublas::vector<T, A>> dataPtr;
61};
62
63/**
64 * \brief Specialization of member function
65 *
66 */
67template <class T, class A>
69 int side, EntityType type, EntitiesFieldData::EntData &data) {
71 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "Not implemented for T = %s",
72 typeid(T).name() // boost::core::demangle(typeid(T).name()).c_str()
73 );
75}
76
77/**
78 * \brief Get value at integration points for scalar field
79 *
80 * \ingroup mofem_forces_and_sources_user_data_operators
81 */
83 : public OpCalculateScalarFieldValues_General<double, DoubleAllocator> {
84
86 double, DoubleAllocator>::OpCalculateScalarFieldValues_General;
87
88 /**
89 * \brief calculate values of scalar field at integration points
90 * @param side side entity number
91 * @param type side entity type
92 * @param data entity data
93 * @return error code
94 */
95 MoFEMErrorCode doWork(int side, EntityType type,
98 VectorDouble &vec = *dataPtr;
99 const size_t nb_gauss_pts = getGaussPts().size2();
100 if (type == zeroType || vec.size() != nb_gauss_pts) {
101 vec.resize(nb_gauss_pts, false);
102 vec.clear();
103 }
104
105 const size_t nb_dofs = data.getFieldData().size();
106
107 if (nb_dofs) {
108
109 if (dataVec.use_count()) {
110 dotVector.resize(nb_dofs, false);
111 const double *array;
112 CHKERR VecGetArrayRead(dataVec, &array);
113 const auto &local_indices = data.getLocalIndices();
114 for (int i = 0; i != local_indices.size(); ++i)
115 if (local_indices[i] != -1)
116 dotVector[i] = array[local_indices[i]];
117 else
118 dotVector[i] = 0;
119 CHKERR VecRestoreArrayRead(dataVec, &array);
120 data.getFieldData().swap(dotVector);
121 }
122
123 const size_t nb_base_functions = data.getN().size2();
124 auto base_function = data.getFTensor0N();
125 auto values_at_gauss_pts = getFTensor0FromVec(vec);
126 for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
127 auto field_data = data.getFTensor0FieldData();
128 size_t bb = 0;
129 for (; bb != nb_dofs; ++bb) {
130 values_at_gauss_pts += field_data * base_function;
131 ++field_data;
132 ++base_function;
133 }
134 // It is possible to have more base functions than dofs
135 for (; bb < nb_base_functions; ++bb)
136 ++base_function;
137 ++values_at_gauss_pts;
138 }
139
140 if (dataVec.use_count()) {
141 data.getFieldData().swap(dotVector);
142 }
143 }
144
146 }
147};
148
149/**
150 * @brief Get rate of scalar field at integration points
151 *
152 * \ingroup mofem_forces_and_sources_user_data_operators
153 */
154
155template <PetscData::DataContext CTX>
158
160 const std::string field_name, boost::shared_ptr<VectorDouble> data_ptr,
161 const EntityType zero_at_type = MBVERTEX)
164 dataPtr(data_ptr), zeroAtType(zero_at_type) {
165 if (!dataPtr)
166 THROW_MESSAGE("Pointer is not set");
167 }
168
169 MoFEMErrorCode doWork(int side, EntityType type,
172
173 const size_t nb_gauss_pts = getGaussPts().size2();
174
175 VectorDouble &vec = *dataPtr;
176 if (type == zeroAtType || vec.size() != nb_gauss_pts) {
177 vec.resize(nb_gauss_pts, false);
178 vec.clear();
179 }
180
181 auto &local_indices = data.getLocalIndices();
182 const size_t nb_dofs = local_indices.size();
183 if (nb_dofs) {
184
185 const double *array;
186
187 auto get_array = [&](const auto ctx, auto vec) {
189 #ifndef NDEBUG
190 if ((getFEMethod()->data_ctx & ctx).none()) {
191 MOFEM_LOG_CHANNEL("SELF");
192 MOFEM_LOG("SELF", Sev::error)
193 << "In this case field degrees of freedom are read from vector. "
194 "That usually happens when time solver is used, and acces to "
195 "first or second rates is needed. You probably not set ts_u, "
196 "ts_u_t, or ts_u_tt and associated data structure, i.e. "
197 "data_ctx to CTX_SET_X, CTX_SET_X_T, or CTX_SET_X_TT "
198 "respectively";
199 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "Vector not set!");
200 }
201 #endif
202 CHKERR VecGetArrayRead(vec, &array);
204 };
205
206 auto restore_array = [&](auto vec) {
207 return VecRestoreArrayRead(vec, &array);
208 };
209
210 switch (CTX) {
212 CHKERR get_array(PetscData::CtxSetX, getFEMethod()->ts_u);
213 break;
215 CHKERR get_array(PetscData::CtxSetX_T, getFEMethod()->ts_u_t);
216 break;
218 CHKERR get_array(PetscData::CtxSetX_TT, getFEMethod()->ts_u_tt);
219 break;
220 default:
221 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
222 "That case is not implemented");
223 }
224
225 std::array<double, MAX_DOFS_ON_ENTITY> dot_dofs_vector;
226 for (int i = 0; i != local_indices.size(); ++i)
227 if (local_indices[i] != -1)
228 dot_dofs_vector[i] = array[local_indices[i]];
229 else
230 dot_dofs_vector[i] = 0;
231
232 switch (CTX) {
234 CHKERR restore_array(getFEMethod()->ts_u);
235 break;
237 CHKERR restore_array(getFEMethod()->ts_u_t);
238 break;
240 CHKERR restore_array(getFEMethod()->ts_u_tt);
241 break;
242 default:
243 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
244 "That case is not implemented");
245 }
246
247 const size_t nb_base_functions = data.getN().size2();
248 auto base_function = data.getFTensor0N();
249 auto values_at_gauss_pts = getFTensor0FromVec(vec);
250
251 for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
252 size_t bb = 0;
253 for (; bb != nb_dofs; ++bb) {
254 values_at_gauss_pts += dot_dofs_vector[bb] * base_function;
255 ++base_function;
256 }
257 // Number of dofs can be smaller than number of Tensor_Dim x base
258 // functions
259 for (; bb < nb_base_functions; ++bb)
260 ++base_function;
261 ++values_at_gauss_pts;
262 }
263 }
265 }
266
267private:
268 boost::shared_ptr<VectorDouble> dataPtr;
270};
271
276
277/**
278 * \deprecated Name inconsistent with other operators
279 *
280 */
282
283/**@}*/
284
285/** \name Vector field values at integration points */
286
287/**@{*/
288
289/** \brief Calculate field values for tenor field rank 1, i.e. vector field
290 *
291 */
292template <int Tensor_Dim, class T, class L, class A>
295
297 const std::string field_name,
298 boost::shared_ptr<ublas::matrix<T, L, A>> data_ptr,
299 const EntityType zero_type = MBVERTEX)
302 dataPtr(data_ptr), zeroType(zero_type) {
303 if (!dataPtr)
304 THROW_MESSAGE("Pointer is not set");
305 }
306
307 /**
308 * \brief calculate values of vector field at integration points
309 * @param side side entity number
310 * @param type side entity type
311 * @param data entity data
312 * @return error code
313 */
314 MoFEMErrorCode doWork(int side, EntityType type,
316
317protected:
318 boost::shared_ptr<ublas::matrix<T, L, A>> dataPtr;
320};
321
322template <int Tensor_Dim, class T, class L, class A>
325 int side, EntityType type, EntitiesFieldData::EntData &data) {
327 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
328 "Not implemented for T = %s and dim = %d",
329 typeid(T).name(), // boost::core::demangle(typeid(T).name()),
330 Tensor_Dim);
332}
333
334/** \brief Calculate field values (template specialization) for tensor field
335 * rank 1, i.e. vector field
336 *
337 */
338template <int Tensor_Dim>
340 ublas::row_major, DoubleAllocator>
342
344 boost::shared_ptr<MatrixDouble> data_ptr,
345 const EntityType zero_type = MBVERTEX)
348 dataPtr(data_ptr), zeroType(zero_type) {
349 if (!dataPtr)
350 THROW_MESSAGE("Pointer is not set");
351 }
352
354 boost::shared_ptr<MatrixDouble> data_ptr,
355 SmartPetscObj<Vec> data_vec,
356 const EntityType zero_type = MBVERTEX)
359 dataPtr(data_ptr), dataVec(data_vec), zeroType(zero_type) {
360 if (!dataPtr)
361 THROW_MESSAGE("Pointer is not set");
362 }
363
364 MoFEMErrorCode doWork(int side, EntityType type,
366
367protected:
368 boost::shared_ptr<MatrixDouble> dataPtr;
372};
373
374/**
375 * \brief Member function specialization calculating values for tenor field rank
376 *
377 */
378template <int Tensor_Dim>
380 Tensor_Dim, double, ublas::row_major,
381 DoubleAllocator>::doWork(int side, EntityType type,
384
385 const size_t nb_gauss_pts = getGaussPts().size2();
386 auto &mat = *dataPtr;
387 if (type == zeroType || mat.size2() != nb_gauss_pts) {
388 mat.resize(Tensor_Dim, nb_gauss_pts, false);
389 mat.clear();
390 }
391
392 const size_t nb_dofs = data.getFieldData().size();
393 if (nb_dofs) {
394
395 if (dataVec.use_count()) {
396 dotVector.resize(nb_dofs, false);
397 const double *array;
398 CHKERR VecGetArrayRead(dataVec, &array);
399 const auto &local_indices = data.getLocalIndices();
400 for (int i = 0; i != local_indices.size(); ++i)
401 if (local_indices[i] != -1)
402 dotVector[i] = array[local_indices[i]];
403 else
404 dotVector[i] = 0;
405 CHKERR VecRestoreArrayRead(dataVec, &array);
406 data.getFieldData().swap(dotVector);
407 }
408
409 if (nb_gauss_pts) {
410 const size_t nb_base_functions = data.getN().size2();
411 auto base_function = data.getFTensor0N();
412 auto values_at_gauss_pts = getFTensor1FromMat<Tensor_Dim>(mat);
413 FTensor::Index<'I', Tensor_Dim> I;
414 const size_t size = nb_dofs / Tensor_Dim;
415 if (nb_dofs % Tensor_Dim) {
416 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
417 "Nb. of DOFs is inconsistent with Tensor_Dim");
418 }
419 for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
420 auto field_data = data.getFTensor1FieldData<Tensor_Dim>();
421
422 #ifndef NDEBUG
423 if (field_data.l2() != field_data.l2()) {
424 MOFEM_LOG("SELF", Sev::error) << "field data: " << field_data;
425 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
426 "Wrong number in coefficients");
427 }
428 #endif
429
430 size_t bb = 0;
431 for (; bb != size; ++bb) {
432
433 #ifndef SINGULARITY
434 #ifndef NDEBUG
435 if (base_function != base_function) {
436 MOFEM_LOG("SELF", Sev::error) << "base function: " << base_function;
437 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
438 "Wrong number number in base functions");
439 }
440 #endif
441 #endif
442
443 values_at_gauss_pts(I) += field_data(I) * base_function;
444 ++field_data;
445 ++base_function;
446 }
447 // Number of dofs can be smaller than number of Tensor_Dim x base
448 // functions
449 for (; bb < nb_base_functions; ++bb)
450 ++base_function;
451 ++values_at_gauss_pts;
452 }
453 }
454
455 if (dataVec.use_count()) {
456 data.getFieldData().swap(dotVector);
457 }
458 }
460}
461
462/** \brief Get values at integration pts for tensor field rank 1, i.e. vector
463 * field
464 *
465 * \ingroup mofem_forces_and_sources_user_data_operators
466 */
467template <int Tensor_Dim>
470 Tensor_Dim, double, ublas::row_major, DoubleAllocator> {
471
473 Tensor_Dim, double, ublas::row_major,
474 DoubleAllocator>::OpCalculateVectorFieldValues_General;
475};
476
477/**@}*/
478
479/** \name Vector field values at integration points */
480
481/**@{*/
482
483/** \brief Calculate field values (template specialization) for tensor field
484 * rank 1, i.e. vector field
485 *
486 */
487template <int Tensor_Dim, CoordinateTypes COORDINATE_SYSTEM = CARTESIAN>
490
492 const std::string field_name, boost::shared_ptr<VectorDouble> data_ptr,
493 const EntityType zero_type = MBVERTEX)
496 dataPtr(data_ptr), zeroType(zero_type) {
497 if (!dataPtr)
498 THROW_MESSAGE("Pointer is not set");
499 }
500
501 MoFEMErrorCode doWork(int side, EntityType type,
504
505 // When we move to C++17 add if constexpr()
506 if constexpr (COORDINATE_SYSTEM == POLAR || COORDINATE_SYSTEM == SPHERICAL)
507 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
508 "%s coordiante not implemented",
509 CoordinateTypesNames[COORDINATE_SYSTEM]);
510
511 const size_t nb_gauss_pts = getGaussPts().size2();
512 auto &vec = *dataPtr;
513 if (type == zeroType) {
514 vec.resize(nb_gauss_pts, false);
515 vec.clear();
516 }
517
518 const size_t nb_dofs = data.getFieldData().size();
519 if (nb_dofs) {
520
521 if (nb_gauss_pts) {
522 const size_t nb_base_functions = data.getN().size2();
523 auto values_at_gauss_pts = getFTensor0FromVec(vec);
524 FTensor::Index<'I', Tensor_Dim> I;
525 const size_t size = nb_dofs / Tensor_Dim;
526 #ifndef NDEBUG
527 if (nb_dofs % Tensor_Dim) {
528 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
529 "Number of dofs should multiple of dimensions");
530 }
531 #endif
532
533 // When we move to C++17 add if constexpr()
534 if constexpr (COORDINATE_SYSTEM == CARTESIAN) {
535 auto diff_base_function = data.getFTensor1DiffN<Tensor_Dim>();
536 for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
537 auto field_data = data.getFTensor1FieldData<Tensor_Dim>();
538 size_t bb = 0;
539 for (; bb != size; ++bb) {
540 values_at_gauss_pts += field_data(I) * diff_base_function(I);
541 ++field_data;
542 ++diff_base_function;
543 }
544 // Number of dofs can be smaller than number of Tensor_Dim x base
545 // functions
546 for (; bb < nb_base_functions; ++bb)
547 ++diff_base_function;
548 ++values_at_gauss_pts;
549 }
550 }
551
552 // When we move to C++17 add if constexpr()
553 if constexpr (COORDINATE_SYSTEM == CYLINDRICAL) {
554 auto t_coords = getFTensor1CoordsAtGaussPts();
555 auto values_at_gauss_pts = getFTensor0FromVec(vec);
556 auto base_function = data.getFTensor0N();
557 auto diff_base_function = data.getFTensor1DiffN<Tensor_Dim>();
558 for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
559 auto field_data = data.getFTensor1FieldData<Tensor_Dim>();
560 size_t bb = 0;
561 for (; bb != size; ++bb) {
562 values_at_gauss_pts += field_data(I) * diff_base_function(I);
563 values_at_gauss_pts +=
564 base_function * (field_data(0) / t_coords(0));
565 ++field_data;
566 ++base_function;
567 ++diff_base_function;
568 }
569 // Number of dofs can be smaller than number of Tensor_Dim x base
570 // functions
571 for (; bb < nb_base_functions; ++bb) {
572 ++base_function;
573 ++diff_base_function;
574 }
575 ++values_at_gauss_pts;
576 ++t_coords;
577 }
578 }
579 }
580 }
582 }
583
584protected:
585 boost::shared_ptr<VectorDouble> dataPtr;
587};
588
589/** \brief Approximate field values for given petsc vector
590 *
591 * \note Look at PetscData to see what vectors could be extracted with that user
592 * data operator.
593 *
594 * \ingroup mofem_forces_and_sources_user_data_operators
595 */
596template <int Tensor_Dim, PetscData::DataContext CTX>
599
601 const std::string field_name, boost::shared_ptr<MatrixDouble> data_ptr,
602 const EntityType zero_at_type = MBVERTEX, bool throw_error = true)
605 dataPtr(data_ptr), zeroAtType(zero_at_type), throwError(throw_error) {
606 if (!dataPtr)
607 THROW_MESSAGE("Pointer is not set");
608 }
609
610 MoFEMErrorCode doWork(int side, EntityType type,
613
614 auto &local_indices = data.getLocalIndices();
615 const size_t nb_dofs = local_indices.size();
616 const size_t nb_gauss_pts = getGaussPts().size2();
617
618 MatrixDouble &mat = *dataPtr;
619 if (type == zeroAtType) {
620 mat.resize(Tensor_Dim, nb_gauss_pts, false);
621 mat.clear();
622 }
623 if (!nb_dofs)
625
626 if (!throwError) {
627 if ((getFEMethod()->data_ctx & PetscData::Switches(CTX)).none()) {
629 }
630 }
631
632 const double *array;
633
634 auto get_array = [&](const auto ctx, auto vec) {
636 #ifndef NDEBUG
637 if ((getFEMethod()->data_ctx & ctx).none()) {
638 MOFEM_LOG_CHANNEL("SELF");
639 MOFEM_LOG("SELF", Sev::error)
640 << "In this case field degrees of freedom are read from vector. "
641 "That usually happens when time solver is used, and access to "
642 "first or second rates is needed. You probably not set ts_u, "
643 "ts_u_t, or ts_u_tt and associated data structure, i.e. "
644 "data_ctx to CTX_SET_X, CTX_SET_DX, CTX_SET_X_T, or "
645 "CTX_SET_X_TT respectively";
646 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "Vector not set");
647 }
648 #endif
649 CHKERR VecGetArrayRead(vec, &array);
651 };
652
653 auto restore_array = [&](auto vec) {
654 return VecRestoreArrayRead(vec, &array);
655 };
656
657 switch (CTX) {
659 CHKERR get_array(PetscData::CtxSetX, getFEMethod()->ts_u);
660 break;
662 CHKERR get_array(PetscData::CtxSetX, getFEMethod()->dx);
663 break;
665 CHKERR get_array(PetscData::CtxSetX_T, getFEMethod()->ts_u_t);
666 break;
668 CHKERR get_array(PetscData::CtxSetX_TT, getFEMethod()->ts_u_tt);
669 break;
670 default:
671 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
672 "That case is not implemented");
673 }
674
675 dotVector.resize(local_indices.size());
676 for (int i = 0; i != local_indices.size(); ++i)
677 if (local_indices[i] != -1)
678 dotVector[i] = array[local_indices[i]];
679 else
680 dotVector[i] = 0;
681
682 switch (CTX) {
684 CHKERR restore_array(getFEMethod()->ts_u);
685 break;
687 CHKERR restore_array(getFEMethod()->dx);
688 break;
690 CHKERR restore_array(getFEMethod()->ts_u_t);
691 break;
693 CHKERR restore_array(getFEMethod()->ts_u_tt);
694 break;
695 default:
696 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
697 "That case is not implemented");
698 }
699
700 const size_t nb_base_functions = data.getN().size2();
701 auto base_function = data.getFTensor0N();
702 auto values_at_gauss_pts = getFTensor1FromMat<Tensor_Dim>(mat);
703
704 FTensor::Index<'I', Tensor_Dim> I;
705 const size_t size = nb_dofs / Tensor_Dim;
706 if (nb_dofs % Tensor_Dim) {
707 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "Data inconsistency");
708 }
709 for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
710 auto field_data = getFTensor1FromArray<Tensor_Dim, Tensor_Dim>(dotVector);
711 size_t bb = 0;
712 for (; bb != size; ++bb) {
713 values_at_gauss_pts(I) += field_data(I) * base_function;
714 ++field_data;
715 ++base_function;
716 }
717 // Number of dofs can be smaller than number of Tensor_Dim x base
718 // functions
719 for (; bb < nb_base_functions; ++bb)
720 ++base_function;
721 ++values_at_gauss_pts;
722 }
724 }
725
726protected:
727 boost::shared_ptr<MatrixDouble> dataPtr;
731};
732
733/** \brief Get rate of values at integration pts for tensor field
734 * rank 1, i.e. vector field
735 *
736 * \ingroup mofem_forces_and_sources_user_data_operators
737 */
738template <int Tensor_Dim>
742
743/** \brief Get second rate of values at integration pts for tensor
744 * field rank 1, i.e. vector field
745 *
746 * \ingroup mofem_forces_and_sources_user_data_operators
747 */
748template <int Tensor_Dim>
752
753/** \brief Get second time second update vector at integration pts for tensor
754 * field rank 1, i.e. vector field
755 *
756 * \ingroup mofem_forces_and_sources_user_data_operators
757 */
758template <int Tensor_Dim>
762
763/**@}*/
764
765/** \name Tensor field values at integration points */
766
767/**@{*/
768
769/** \brief Calculate field values for tenor field rank 2.
770 *
771 */
772template <int Tensor_Dim0, int Tensor_Dim1, class T, class L, class A>
775
777 const std::string field_name,
778 boost::shared_ptr<ublas::matrix<T, L, A>> data_ptr,
779 const EntityType zero_type = MBVERTEX)
782 dataPtr(data_ptr), zeroType(zero_type) {
783 if (!dataPtr)
784 THROW_MESSAGE("Pointer is not set");
785 }
786
787 MoFEMErrorCode doWork(int side, EntityType type,
789
790protected:
791 boost::shared_ptr<ublas::matrix<T, L, A>> dataPtr;
793};
794
795template <int Tensor_Dim0, int Tensor_Dim1, class T, class L, class A>
798 doWork(int side, EntityType type, EntitiesFieldData::EntData &data) {
800 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
801 "Not implemented for T = %s, dim0 = %d and dim1 = %d",
802 typeid(T).name(), // boost::core::demangle(typeid(T).name()),
803 Tensor_Dim0, Tensor_Dim1);
805}
806
807template <int Tensor_Dim0, int Tensor_Dim1>
808struct OpCalculateTensor2FieldValues_General<Tensor_Dim0, Tensor_Dim1, double,
809 ublas::row_major, DoubleAllocator>
811
813 const std::string field_name,
814 boost::shared_ptr<
815 ublas::matrix<double, ublas::row_major, DoubleAllocator>>
816 data_ptr,
817 const EntityType zero_type = MBVERTEX)
820 dataPtr(data_ptr), zeroType(zero_type) {
821 if (!dataPtr)
822 THROW_MESSAGE("Pointer is not set");
823 }
824
826 const std::string field_name,
827 boost::shared_ptr<
828 ublas::matrix<double, ublas::row_major, DoubleAllocator>>
829 data_ptr,
830 SmartPetscObj<Vec> data_vec, const EntityType zero_type = MBVERTEX)
833 dataPtr(data_ptr), zeroType(zero_type), dataVec(data_vec) {
834 if (!dataPtr)
835 THROW_MESSAGE("Pointer is not set");
836 }
837
838 MoFEMErrorCode doWork(int side, EntityType type,
840
841protected:
842 boost::shared_ptr<ublas::matrix<double, ublas::row_major, DoubleAllocator>>
847};
848
849template <int Tensor_Dim0, int Tensor_Dim1>
851 Tensor_Dim0, Tensor_Dim1, double, ublas::row_major,
852 DoubleAllocator>::doWork(int side, EntityType type,
855
856 MatrixDouble &mat = *dataPtr;
857 const size_t nb_gauss_pts = data.getN().size1();
858 if (type == zeroType) {
859 mat.resize(Tensor_Dim0 * Tensor_Dim1, nb_gauss_pts, false);
860 mat.clear();
861 }
862
863 const size_t nb_dofs = data.getFieldData().size();
864
865 if (dataVec.use_count()) {
866 dotVector.resize(nb_dofs, false);
867 const double *array;
868 CHKERR VecGetArrayRead(dataVec, &array);
869 const auto &local_indices = data.getLocalIndices();
870 for (int i = 0; i != local_indices.size(); ++i)
871 if (local_indices[i] != -1)
872 dotVector[i] = array[local_indices[i]];
873 else
874 dotVector[i] = 0;
875 CHKERR VecRestoreArrayRead(dataVec, &array);
876 data.getFieldData().swap(dotVector);
877 }
878
879 if (nb_dofs) {
880 const size_t nb_base_functions = data.getN().size2();
881 auto base_function = data.getFTensor0N();
882 auto values_at_gauss_pts =
883 getFTensor2FromMat<Tensor_Dim0, Tensor_Dim1>(mat);
884 FTensor::Index<'i', Tensor_Dim0> i;
885 FTensor::Index<'j', Tensor_Dim1> j;
886 const size_t size = nb_dofs / (Tensor_Dim0 * Tensor_Dim1);
887 for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
888 auto field_data = data.getFTensor2FieldData<Tensor_Dim0, Tensor_Dim1>();
889 size_t bb = 0;
890 for (; bb != size; ++bb) {
891 values_at_gauss_pts(i, j) += field_data(i, j) * base_function;
892 ++field_data;
893 ++base_function;
894 }
895 for (; bb < nb_base_functions; ++bb)
896 ++base_function;
897 ++values_at_gauss_pts;
898 }
899
900 if (dataVec.use_count()) {
901 data.getFieldData().swap(dotVector);
902 }
903 }
905}
906
907/** \brief Get values at integration pts for tensor field rank 2, i.e. matrix
908 * field
909 *
910 * \ingroup mofem_forces_and_sources_user_data_operators
911 */
912template <int Tensor_Dim0, int Tensor_Dim1>
915 Tensor_Dim0, Tensor_Dim1, double, ublas::row_major, DoubleAllocator> {
916
918 Tensor_Dim0, Tensor_Dim1, double, ublas::row_major,
919 DoubleAllocator>::OpCalculateTensor2FieldValues_General;
920};
921
922/** \brief Get time direvarive values at integration pts for tensor field rank
923 * 2, i.e. matrix field
924 *
925 * \ingroup mofem_forces_and_sources_user_data_operators
926 */
927template <int Tensor_Dim0, int Tensor_Dim1>
930
932 boost::shared_ptr<MatrixDouble> data_ptr,
933 const EntityType zero_at_type = MBVERTEX)
936 dataPtr(data_ptr), zeroAtType(zero_at_type) {
937 if (!dataPtr)
938 THROW_MESSAGE("Pointer is not set");
939 }
940
941 MoFEMErrorCode doWork(int side, EntityType type,
944
945 const size_t nb_gauss_pts = getGaussPts().size2();
946 MatrixDouble &mat = *dataPtr;
947 if (type == zeroAtType) {
948 mat.resize(Tensor_Dim0 * Tensor_Dim1, nb_gauss_pts, false);
949 mat.clear();
950 }
951 const auto &local_indices = data.getLocalIndices();
952 const size_t nb_dofs = local_indices.size();
953 if (nb_dofs) {
954 dotVector.resize(nb_dofs, false);
955 const double *array;
956 CHKERR VecGetArrayRead(getFEMethod()->ts_u_t, &array);
957 for (size_t i = 0; i != local_indices.size(); ++i)
958 if (local_indices[i] != -1)
959 dotVector[i] = array[local_indices[i]];
960 else
961 dotVector[i] = 0;
962 CHKERR VecRestoreArrayRead(getFEMethod()->ts_u_t, &array);
963
964 const size_t nb_base_functions = data.getN().size2();
965
966 auto base_function = data.getFTensor0N();
967 auto values_at_gauss_pts =
968 getFTensor2FromMat<Tensor_Dim0, Tensor_Dim1>(mat);
969 FTensor::Index<'i', Tensor_Dim0> i;
970 FTensor::Index<'j', Tensor_Dim1> j;
971 const size_t size = nb_dofs / (Tensor_Dim0 * Tensor_Dim1);
972 for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
973 auto field_data =
974 getFTensor2FromPtr<Tensor_Dim0, Tensor_Dim1>(&*dotVector.begin());
975 size_t bb = 0;
976 for (; bb != size; ++bb) {
977 values_at_gauss_pts(i, j) += field_data(i, j) * base_function;
978 ++field_data;
979 ++base_function;
980 }
981 for (; bb < nb_base_functions; ++bb)
982 ++base_function;
983 ++values_at_gauss_pts;
984 }
985 }
987 }
988
989protected:
990 boost::shared_ptr<MatrixDouble> dataPtr; ///< Data computed into this matrix
991 EntityType zeroAtType; ///< Zero values at Gauss point at this type
992 VectorDouble dotVector; ///< Keeps temporary values of time derivatives
993};
994
995/**
996 * @brief Calculate symmetric tensor field values at integration pts.
997 *
998 * @tparam Tensor_Dim
999
1000 * \ingroup mofem_forces_and_sources_user_data_operators
1001 */
1002template <int Tensor_Dim>
1005
1007 const std::string field_name, boost::shared_ptr<MatrixDouble> data_ptr,
1008 const EntityType zero_type = MBEDGE, const int zero_side = 0)
1011 dataPtr(data_ptr), zeroType(zero_type), zeroSide(zero_side) {
1012 if (!dataPtr)
1013 THROW_MESSAGE("Pointer is not set");
1014 }
1015
1017 const std::string field_name, boost::shared_ptr<MatrixDouble> data_ptr,
1018 SmartPetscObj<Vec> data_vec, const EntityType zero_type = MBEDGE,
1019 const int zero_side = 0)
1022 dataPtr(data_ptr), zeroType(zero_type), zeroSide(zero_side),
1023 dataVec(data_vec) {
1024 if (!dataPtr)
1025 THROW_MESSAGE("Pointer is not set");
1026 }
1027
1028 MoFEMErrorCode doWork(int side, EntityType type,
1031 MatrixDouble &mat = *dataPtr;
1032 const int nb_gauss_pts = getGaussPts().size2();
1033 if (type == this->zeroType && side == zeroSide) {
1034 mat.resize((Tensor_Dim * (Tensor_Dim + 1)) / 2, nb_gauss_pts, false);
1035 mat.clear();
1036 }
1037 const int nb_dofs = data.getFieldData().size();
1038 if (!nb_dofs)
1040
1041 if (dataVec.use_count()) {
1042 dotVector.resize(nb_dofs, false);
1043 const double *array;
1044 CHKERR VecGetArrayRead(dataVec, &array);
1045 const auto &local_indices = data.getLocalIndices();
1046 for (int i = 0; i != local_indices.size(); ++i)
1047 if (local_indices[i] != -1)
1048 dotVector[i] = array[local_indices[i]];
1049 else
1050 dotVector[i] = 0;
1051 CHKERR VecRestoreArrayRead(dataVec, &array);
1052 data.getFieldData().swap(dotVector);
1053 }
1054
1055 const int nb_base_functions = data.getN().size2();
1056 auto base_function = data.getFTensor0N();
1057 auto values_at_gauss_pts = getFTensor2SymmetricFromMat<Tensor_Dim>(mat);
1058 FTensor::Index<'i', Tensor_Dim> i;
1059 FTensor::Index<'j', Tensor_Dim> j;
1060 const int size = nb_dofs / ((Tensor_Dim * (Tensor_Dim + 1)) / 2);
1061 for (int gg = 0; gg != nb_gauss_pts; ++gg) {
1062 auto field_data = data.getFTensor2SymmetricFieldData<Tensor_Dim>();
1063 int bb = 0;
1064 for (; bb != size; ++bb) {
1065 values_at_gauss_pts(i, j) += field_data(i, j) * base_function;
1066 ++field_data;
1067 ++base_function;
1068 }
1069 for (; bb < nb_base_functions; ++bb)
1070 ++base_function;
1071 ++values_at_gauss_pts;
1072 }
1073
1074 if (dataVec.use_count()) {
1075 data.getFieldData().swap(dotVector);
1076 }
1077
1079 }
1080
1081protected:
1082 boost::shared_ptr<MatrixDouble> dataPtr;
1084 const int zeroSide;
1087};
1088
1089/**
1090 * @brief Calculate symmetric tensor field rates ant integratio pts.
1091 *
1092 * @tparam Tensor_Dim
1093 *
1094 * \ingroup mofem_forces_and_sources_user_data_operators
1095 */
1096template <int Tensor_Dim>
1099
1101 const std::string field_name, boost::shared_ptr<MatrixDouble> data_ptr,
1102 const EntityType zero_type = MBEDGE, const int zero_side = 0)
1105 dataPtr(data_ptr), zeroType(zero_type), zeroSide(zero_side) {
1106 if (!dataPtr)
1107 THROW_MESSAGE("Pointer is not set");
1108 }
1109
1110 MoFEMErrorCode doWork(int side, EntityType type,
1113 const int nb_gauss_pts = getGaussPts().size2();
1114 MatrixDouble &mat = *dataPtr;
1115 constexpr auto symm_size = (Tensor_Dim * (Tensor_Dim + 1)) / 2;
1116 if (type == zeroType && side == zeroSide) {
1117 mat.resize(symm_size, nb_gauss_pts, false);
1118 mat.clear();
1119 }
1120 auto &local_indices = data.getLocalIndices();
1121 const int nb_dofs = local_indices.size();
1122 if (!nb_dofs)
1124
1125 #ifndef NDEBUG
1126 if ((getFEMethod()->data_ctx & PetscData::CtxSetX_T).none()) {
1127 MOFEM_LOG_CHANNEL("SELF");
1128 MOFEM_LOG("SELF", Sev::error)
1129 << "In this case field degrees of freedom are read from vector. "
1130 "That usually happens when time solver is used, and acces to "
1131 "first rates is needed. You probably not set "
1132 "ts_u_t and associated data structure data_ctx to CTX_SET_X_T "
1133 "respectively";
1134 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "Vector not set!");
1135 }
1136 #endif
1137
1138 dotVector.resize(nb_dofs, false);
1139 const double *array;
1140 CHKERR VecGetArrayRead(getFEMethod()->ts_u_t, &array);
1141 for (int i = 0; i != local_indices.size(); ++i)
1142 if (local_indices[i] != -1)
1143 dotVector[i] = array[local_indices[i]];
1144 else
1145 dotVector[i] = 0;
1146 CHKERR VecRestoreArrayRead(getFEMethod()->ts_u_t, &array);
1147
1148 const int nb_base_functions = data.getN().size2();
1149
1150 auto base_function = data.getFTensor0N();
1151 auto values_at_gauss_pts = getFTensor2SymmetricFromMat<Tensor_Dim>(mat);
1152 FTensor::Index<'i', Tensor_Dim> i;
1153 FTensor::Index<'j', Tensor_Dim> j;
1154 const int size = nb_dofs / symm_size;
1155 for (int gg = 0; gg != nb_gauss_pts; ++gg) {
1156 auto field_data = getFTensorDotData<Tensor_Dim>();
1157 int bb = 0;
1158 for (; bb != size; ++bb) {
1159 values_at_gauss_pts(i, j) += field_data(i, j) * base_function;
1160 ++field_data;
1161 ++base_function;
1162 }
1163 for (; bb < nb_base_functions; ++bb)
1164 ++base_function;
1165 ++values_at_gauss_pts;
1166 }
1167
1169 }
1170
1171protected:
1172 boost::shared_ptr<MatrixDouble> dataPtr;
1174 const int zeroSide;
1176
1177 template <int Dim> inline auto getFTensorDotData() {
1178 static_assert(Dim || !Dim, "not implemented");
1179 }
1180};
1181
1182template <>
1183template <>
1184inline auto
1187 &dotVector[0], &dotVector[1], &dotVector[2], &dotVector[3], &dotVector[4],
1188 &dotVector[5]);
1189}
1190
1191template <>
1192template <>
1193inline auto
1196 &dotVector[0], &dotVector[1], &dotVector[2]);
1197}
1198
1199/**@}*/
1200
1201/** \name Gradients and Hessian of scalar fields at integration points */
1202
1203/**@{*/
1204
1205/**
1206 * \brief Evaluate field gradient values for scalar field, i.e. gradient is
1207 * tensor rank 1 (vector)
1208 *
1209 */
1210template <int Tensor_Dim, class T, class L, class A>
1212 : public OpCalculateVectorFieldValues_General<Tensor_Dim, T, L, A> {
1213
1215 Tensor_Dim, T, L, A>::OpCalculateVectorFieldValues_General;
1216};
1217
1218/** \brief Evaluate field gradient values for scalar field, i.e. gradient is
1219 * tensor rank 1 (vector), specialization
1220 *
1221 */
1222template <int Tensor_Dim>
1224 ublas::row_major, DoubleAllocator>
1226 Tensor_Dim, double, ublas::row_major, DoubleAllocator> {
1227
1229 Tensor_Dim, double, ublas::row_major,
1230 DoubleAllocator>::OpCalculateVectorFieldValues_General;
1231
1232 /**
1233 * \brief calculate gradient values of scalar field at integration points
1234 * @param side side entity number
1235 * @param type side entity type
1236 * @param data entity data
1237 * @return error code
1238 */
1239 MoFEMErrorCode doWork(int side, EntityType type,
1241};
1242
1243/**
1244 * \brief Member function specialization calculating scalar field gradients for
1245 * tenor field rank 1
1246 *
1247 */
1248template <int Tensor_Dim>
1250 Tensor_Dim, double, ublas::row_major,
1251 DoubleAllocator>::doWork(int side, EntityType type,
1254
1255 const size_t nb_gauss_pts = this->getGaussPts().size2();
1256 auto &mat = *this->dataPtr;
1257 if (type == this->zeroType) {
1258 mat.resize(Tensor_Dim, nb_gauss_pts, false);
1259 mat.clear();
1260 }
1261
1262 const int nb_dofs = data.getFieldData().size();
1263 if (nb_dofs) {
1264
1265 if (nb_gauss_pts) {
1266 const int nb_base_functions = data.getN().size2();
1267 auto diff_base_function = data.getFTensor1DiffN<Tensor_Dim>();
1268 auto gradients_at_gauss_pts = getFTensor1FromMat<Tensor_Dim>(mat);
1269
1270 #ifndef NDEBUG
1271 if (nb_dofs > nb_base_functions)
1272 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1273 "Number of base functions inconsistent with number of DOFs "
1274 "(%d > %d)",
1275 nb_dofs, nb_base_functions);
1276
1277 if (data.getDiffN().size2() != nb_base_functions * Tensor_Dim)
1278 SETERRQ(
1279 PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1280 "Number of base functions inconsistent with number of derivatives "
1281 "(%lu != %d)",
1282 data.getDiffN().size2(), nb_base_functions);
1283
1284 if (data.getDiffN().size1() != nb_gauss_pts)
1285 SETERRQ(
1286 PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1287 "Number of base functions inconsistent with number of integration "
1288 "pts (%lu != %lu)",
1289 data.getDiffN().size2(), nb_gauss_pts);
1290
1291 #endif
1292
1293 FTensor::Index<'I', Tensor_Dim> I;
1294 for (int gg = 0; gg < nb_gauss_pts; ++gg) {
1295 auto field_data = data.getFTensor0FieldData();
1296 int bb = 0;
1297 for (; bb != nb_dofs; ++bb) {
1298 gradients_at_gauss_pts(I) += field_data * diff_base_function(I);
1299 ++field_data;
1300 ++diff_base_function;
1301 }
1302 // Number of dofs can be smaller than number of base functions
1303 for (; bb < nb_base_functions; ++bb)
1304 ++diff_base_function;
1305 ++gradients_at_gauss_pts;
1306 }
1307 }
1308 }
1309
1311}
1312
1313/** \brief Get field gradients at integration pts for scalar field rank 0, i.e.
1314 * vector field
1315 *
1316 * \ingroup mofem_forces_and_sources_user_data_operators
1317 */
1318template <int Tensor_Dim>
1321 Tensor_Dim, double, ublas::row_major, DoubleAllocator> {
1323 Tensor_Dim, double, ublas::row_major,
1324 DoubleAllocator>::OpCalculateScalarFieldGradient_General;
1325};
1326
1327/** \brief Evaluate field gradient values for scalar field, i.e. gradient is
1328 * tensor rank 1 (vector), specialization
1329 *
1330 */
1331template <int Tensor_Dim>
1334 Tensor_Dim, double, ublas::row_major, DoubleAllocator> {
1335
1337 Tensor_Dim, double, ublas::row_major,
1338 DoubleAllocator>::OpCalculateVectorFieldValues_General;
1339
1340 /**
1341 * \brief calculate gradient values of scalar field at integration points
1342 * @param side side entity number
1343 * @param type side entity type
1344 * @param data entity data
1345 * @return error code
1346 */
1347 MoFEMErrorCode doWork(int side, EntityType type,
1349};
1350
1351template <int Tensor_Dim>
1353 int side, EntityType type, EntitiesFieldData::EntData &data) {
1355
1356 const size_t nb_gauss_pts = this->getGaussPts().size2();
1357 auto &mat = *this->dataPtr;
1358 if (type == this->zeroType) {
1359 mat.resize(Tensor_Dim * Tensor_Dim, nb_gauss_pts, false);
1360 mat.clear();
1361 }
1362
1363 const int nb_dofs = data.getFieldData().size();
1364 if (nb_dofs) {
1365
1366 if (nb_gauss_pts) {
1367 const int nb_base_functions = data.getN().size2();
1368
1369 auto &hessian_base = data.getN(BaseDerivatives::SecondDerivative);
1370 #ifndef NDEBUG
1371 if (hessian_base.size1() != nb_gauss_pts) {
1372 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1373 "Wrong number of integration pts (%ld != %d)",
1374 hessian_base.size1(), nb_gauss_pts);
1375 }
1376 if (hessian_base.size2() != nb_base_functions * Tensor_Dim * Tensor_Dim) {
1377 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1378 "Wrong number of base functions (%ld != %d)",
1379 hessian_base.size2() / (Tensor_Dim * Tensor_Dim),
1380 nb_base_functions);
1381 }
1382 if (hessian_base.size2() < nb_dofs * Tensor_Dim * Tensor_Dim) {
1383 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1384 "Wrong number of base functions (%ld < %d)",
1385 hessian_base.size2(), nb_dofs * Tensor_Dim * Tensor_Dim);
1386 }
1387 #endif
1388
1389 auto t_diff2_base_function = getFTensor2FromPtr<Tensor_Dim, Tensor_Dim>(
1390 &*hessian_base.data().begin());
1391
1392 auto hessian_at_gauss_pts =
1393 getFTensor2FromMat<Tensor_Dim, Tensor_Dim>(mat);
1394
1395 FTensor::Index<'I', Tensor_Dim> I;
1396 FTensor::Index<'J', Tensor_Dim> J;
1397 for (int gg = 0; gg != nb_gauss_pts; ++gg) {
1398 auto field_data = data.getFTensor0FieldData();
1399 int bb = 0;
1400 for (; bb != nb_dofs; ++bb) {
1401 hessian_at_gauss_pts(I, J) +=
1402 field_data * t_diff2_base_function(I, J);
1403 ++field_data;
1404 ++t_diff2_base_function;
1405 }
1406 // Number of dofs can be smaller than number of base functions
1407 for (; bb < nb_base_functions; ++bb) {
1408 ++t_diff2_base_function;
1409 }
1410
1411 ++hessian_at_gauss_pts;
1412 }
1413 }
1414 }
1415
1417}
1418
1419/**}*/
1420
1421/** \name Gradients and hessian of tensor fields at integration points */
1422
1423/**@{*/
1424
1425/**
1426 * \brief Evaluate field gradient values for vector field, i.e. gradient is
1427 * tensor rank 2
1428 *
1429 */
1430template <int Tensor_Dim0, int Tensor_Dim1, class T, class L, class A>
1432 : public OpCalculateTensor2FieldValues_General<Tensor_Dim0, Tensor_Dim1, T,
1433 L, A> {
1434
1436 Tensor_Dim0, Tensor_Dim1, T, L, A>::OpCalculateTensor2FieldValues_General;
1437};
1438
1439template <int Tensor_Dim0, int Tensor_Dim1>
1440struct OpCalculateVectorFieldGradient_General<Tensor_Dim0, Tensor_Dim1, double,
1441 ublas::row_major, DoubleAllocator>
1443 Tensor_Dim0, Tensor_Dim1, double, ublas::row_major, DoubleAllocator> {
1444
1446 Tensor_Dim0, Tensor_Dim1, double, ublas::row_major,
1447 DoubleAllocator>::OpCalculateTensor2FieldValues_General;
1448
1449 /**
1450 * \brief calculate values of vector field at integration points
1451 * @param side side entity number
1452 * @param type side entity type
1453 * @param data entity data
1454 * @return error code
1455 */
1456 MoFEMErrorCode doWork(int side, EntityType type,
1458};
1459
1460/**
1461 * \brief Member function specialization calculating vector field gradients for
1462 * tenor field rank 2
1463 *
1464 */
1465template <int Tensor_Dim0, int Tensor_Dim1>
1467 Tensor_Dim0, Tensor_Dim1, double, ublas::row_major,
1468 DoubleAllocator>::doWork(int side, EntityType type,
1471 if (!this->dataPtr)
1472 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1473 "Data pointer not allocated");
1474
1475 const size_t nb_gauss_pts = this->getGaussPts().size2();
1476 auto &mat = *this->dataPtr;
1477 if (type == this->zeroType) {
1478 mat.resize(Tensor_Dim0 * Tensor_Dim1, nb_gauss_pts, false);
1479 mat.clear();
1480 }
1481
1482 if (nb_gauss_pts) {
1483 const size_t nb_dofs = data.getFieldData().size();
1484
1485 if (nb_dofs) {
1486
1487 if (this->dataVec.use_count()) {
1488 this->dotVector.resize(nb_dofs, false);
1489 const double *array;
1490 CHKERR VecGetArrayRead(this->dataVec, &array);
1491 const auto &local_indices = data.getLocalIndices();
1492 for (int i = 0; i != local_indices.size(); ++i)
1493 if (local_indices[i] != -1)
1494 this->dotVector[i] = array[local_indices[i]];
1495 else
1496 this->dotVector[i] = 0;
1497 CHKERR VecRestoreArrayRead(this->dataVec, &array);
1498 data.getFieldData().swap(this->dotVector);
1499 }
1500
1501 const int nb_base_functions = data.getN().size2();
1502 auto diff_base_function = data.getFTensor1DiffN<Tensor_Dim1>();
1503 auto gradients_at_gauss_pts =
1504 getFTensor2FromMat<Tensor_Dim0, Tensor_Dim1>(mat);
1505 FTensor::Index<'I', Tensor_Dim0> I;
1506 FTensor::Index<'J', Tensor_Dim1> J;
1507 int size = nb_dofs / Tensor_Dim0;
1508 if (nb_dofs % Tensor_Dim0) {
1509 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1510 "Data inconsistency");
1511 }
1512 for (int gg = 0; gg < nb_gauss_pts; ++gg) {
1513 auto field_data = data.getFTensor1FieldData<Tensor_Dim0>();
1514
1515 #ifndef NDEBUG
1516 if (field_data.l2() != field_data.l2()) {
1517 MOFEM_LOG("SELF", Sev::error) << "field data " << field_data;
1518 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1519 "Wrong number in coefficients");
1520 }
1521 #endif
1522
1523 int bb = 0;
1524 for (; bb < size; ++bb) {
1525 #ifndef SINGULARITY
1526 #ifndef NDEBUG
1527 if (diff_base_function.l2() != diff_base_function.l2()) {
1528 MOFEM_LOG("SELF", Sev::error)
1529 << "diff_base_function: " << diff_base_function;
1530 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1531 "Wrong number number in base functions");
1532 }
1533 #endif
1534 #endif
1535
1536 gradients_at_gauss_pts(I, J) += field_data(I) * diff_base_function(J);
1537 ++field_data;
1538 ++diff_base_function;
1539 }
1540 // Number of dofs can be smaller than number of Tensor_Dim0 x base
1541 // functions
1542 for (; bb != nb_base_functions; ++bb)
1543 ++diff_base_function;
1544 ++gradients_at_gauss_pts;
1545 }
1546
1547 if (this->dataVec.use_count()) {
1548 data.getFieldData().swap(this->dotVector);
1549 }
1550 }
1551 }
1553}
1554
1555/** \brief Get field gradients at integration pts for scalar field rank 0, i.e.
1556 * vector field
1557 *
1558 * \ingroup mofem_forces_and_sources_user_data_operators
1559 */
1560template <int Tensor_Dim0, int Tensor_Dim1>
1563 Tensor_Dim0, Tensor_Dim1, double, ublas::row_major, DoubleAllocator> {
1564
1566 Tensor_Dim0, Tensor_Dim1, double, ublas::row_major,
1567 DoubleAllocator>::OpCalculateVectorFieldGradient_General;
1568};
1569
1570/** \brief Get field gradients time derivative at integration pts for scalar
1571 * field rank 0, i.e. vector field
1572 *
1573 * \ingroup mofem_forces_and_sources_user_data_operators
1574 */
1575template <int Tensor_Dim0, int Tensor_Dim1>
1578
1580 boost::shared_ptr<MatrixDouble> data_ptr,
1581 const EntityType zero_at_type = MBVERTEX)
1584 dataPtr(data_ptr), zeroAtType(zero_at_type) {
1585 if (!dataPtr)
1586 THROW_MESSAGE("Pointer is not set");
1587 }
1588
1589 MoFEMErrorCode doWork(int side, EntityType type,
1592
1593 const auto &local_indices = data.getLocalIndices();
1594 const int nb_dofs = local_indices.size();
1595 const int nb_gauss_pts = this->getGaussPts().size2();
1596
1597 auto &mat = *this->dataPtr;
1598 if (type == this->zeroAtType) {
1599 mat.resize(Tensor_Dim0 * Tensor_Dim1, nb_gauss_pts, false);
1600 mat.clear();
1601 }
1602 if (!nb_dofs)
1604
1605 dotVector.resize(nb_dofs, false);
1606 const double *array;
1607 CHKERR VecGetArrayRead(getFEMethod()->ts_u_t, &array);
1608 for (int i = 0; i != local_indices.size(); ++i)
1609 if (local_indices[i] != -1)
1610 dotVector[i] = array[local_indices[i]];
1611 else
1612 dotVector[i] = 0;
1613 CHKERR VecRestoreArrayRead(getFEMethod()->ts_u_t, &array);
1614
1615 const int nb_base_functions = data.getN().size2();
1616 auto diff_base_function = data.getFTensor1DiffN<Tensor_Dim1>();
1617 auto gradients_at_gauss_pts =
1618 getFTensor2FromMat<Tensor_Dim0, Tensor_Dim1>(mat);
1619 FTensor::Index<'I', Tensor_Dim0> I;
1620 FTensor::Index<'J', Tensor_Dim1> J;
1621 int size = nb_dofs / Tensor_Dim0;
1622 if (nb_dofs % Tensor_Dim0) {
1623 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "Data inconsistency");
1624 }
1625
1626 for (int gg = 0; gg < nb_gauss_pts; ++gg) {
1627 auto field_data = getFTensor1FromPtr<Tensor_Dim0>(&*dotVector.begin());
1628 int bb = 0;
1629 for (; bb < size; ++bb) {
1630 gradients_at_gauss_pts(I, J) += field_data(I) * diff_base_function(J);
1631 ++field_data;
1632 ++diff_base_function;
1633 }
1634 // Number of dofs can be smaller than number of Tensor_Dim0 x base
1635 // functions
1636 for (; bb != nb_base_functions; ++bb)
1637 ++diff_base_function;
1638 ++gradients_at_gauss_pts;
1639 }
1641 }
1642
1643private:
1644 boost::shared_ptr<MatrixDouble> dataPtr; ///< Data computed into this matrix
1645 EntityType zeroAtType; ///< Zero values at Gauss point at this type
1646 VectorDouble dotVector; ///< Keeps temporary values of time derivatives
1647};
1648
1649/**
1650 * \brief Evaluate field gradient values for symmetric 2nd order tensor field,
1651 * i.e. gradient is tensor rank 3
1652 *
1653 */
1654template <int Tensor_Dim0, int Tensor_Dim1, class T, class L, class A>
1656 : public OpCalculateTensor2FieldValues_General<Tensor_Dim0, Tensor_Dim1, T,
1657 L, A> {
1658
1660 const std::string field_name, boost::shared_ptr<MatrixDouble> data_ptr,
1661 const EntityType zero_type = MBVERTEX)
1662 : OpCalculateTensor2FieldValues_General<Tensor_Dim0, Tensor_Dim1, T, L,
1663 A>(field_name, data_ptr,
1664 zero_type) {}
1665};
1666
1667template <int Tensor_Dim0, int Tensor_Dim1>
1669 Tensor_Dim0, Tensor_Dim1, double, ublas::row_major, DoubleAllocator>
1671 Tensor_Dim0, Tensor_Dim1, double, ublas::row_major, DoubleAllocator> {
1672
1674 const std::string field_name, boost::shared_ptr<MatrixDouble> data_ptr,
1675 const EntityType zero_type = MBVERTEX)
1676 : OpCalculateTensor2FieldValues_General<Tensor_Dim0, Tensor_Dim1, double,
1677 ublas::row_major,
1679 field_name, data_ptr, zero_type) {}
1680
1681 /**
1682 * \brief calculate values of vector field at integration points
1683 * @param side side entity number
1684 * @param type side entity type
1685 * @param data entity data
1686 * @return error code
1687 */
1688 MoFEMErrorCode doWork(int side, EntityType type,
1690};
1691
1692/**
1693 * \brief Member function specialization calculating tensor field gradients for
1694 * symmetric tensor field rank 2
1695 *
1696 */
1697template <int Tensor_Dim0, int Tensor_Dim1>
1698MoFEMErrorCode OpCalculateTensor2SymmetricFieldGradient_General<
1699 Tensor_Dim0, Tensor_Dim1, double, ublas::row_major,
1700 DoubleAllocator>::doWork(int side, EntityType type,
1703 if (!this->dataPtr)
1704 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1705 "Data pointer not allocated");
1706
1707 const size_t nb_gauss_pts = this->getGaussPts().size2();
1708 constexpr size_t msize = (Tensor_Dim0 * (Tensor_Dim0 + 1)) / 2;
1709 auto &mat = *this->dataPtr;
1710 if (type == this->zeroType) {
1711 mat.resize(msize * Tensor_Dim1, nb_gauss_pts, false);
1712 mat.clear();
1713 }
1714
1715 if (nb_gauss_pts) {
1716 const size_t nb_dofs = data.getFieldData().size();
1717
1718 if (nb_dofs) {
1719
1720 const int nb_base_functions = data.getN().size2();
1721 auto diff_base_function = data.getFTensor1DiffN<Tensor_Dim1>();
1722 auto gradients_at_gauss_pts =
1723 getFTensor3DgFromMat<Tensor_Dim0, Tensor_Dim1>(mat);
1724 FTensor::Index<'I', Tensor_Dim0> I;
1725 FTensor::Index<'J', Tensor_Dim0> J;
1726 FTensor::Index<'K', Tensor_Dim1> K;
1727 int size = nb_dofs / msize;
1728 if (nb_dofs % msize) {
1729 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1730 "Data inconsistency");
1731 }
1732 for (int gg = 0; gg < nb_gauss_pts; ++gg) {
1733 auto field_data = data.getFTensor2SymmetricFieldData<Tensor_Dim0>();
1734 int bb = 0;
1735 for (; bb < size; ++bb) {
1736 gradients_at_gauss_pts(I, J, K) +=
1737 field_data(I, J) * diff_base_function(K);
1738 ++field_data;
1739 ++diff_base_function;
1740 }
1741 // Number of dofs can be smaller than number of Tensor_Dim0 x base
1742 // functions
1743 for (; bb != nb_base_functions; ++bb)
1744 ++diff_base_function;
1745 ++gradients_at_gauss_pts;
1746 }
1747 }
1748 }
1750}
1751
1752/** \brief Get field gradients at integration pts for symmetric tensorial field
1753 * rank 2
1754 *
1755 * \ingroup mofem_forces_and_sources_user_data_operators
1756 */
1757template <int Tensor_Dim0, int Tensor_Dim1>
1760 Tensor_Dim0, Tensor_Dim1, double, ublas::row_major, DoubleAllocator> {
1761
1763 const std::string field_name, boost::shared_ptr<MatrixDouble> data_ptr,
1764 const EntityType zero_type = MBVERTEX)
1766 Tensor_Dim0, Tensor_Dim1, double, ublas::row_major,
1767 DoubleAllocator>(field_name, data_ptr, zero_type) {}
1768};
1769
1770template <int Tensor_Dim0, int Tensor_Dim1>
1773 Tensor_Dim0, Tensor_Dim1, double, ublas::row_major, DoubleAllocator> {
1774
1776 Tensor_Dim0, Tensor_Dim1, double, ublas::row_major,
1777 DoubleAllocator>::OpCalculateTensor2FieldValues_General;
1778
1779 /**
1780 * \brief calculate values of vector field at integration points
1781 * @param side side entity number
1782 * @param type side entity type
1783 * @param data entity data
1784 * @return error code
1785 */
1786 MoFEMErrorCode doWork(int side, EntityType type,
1788};
1789
1790/**
1791 * \brief Member function specialization calculating vector field gradients for
1792 * tenor field rank 2
1793 *
1794 */
1795template <int Tensor_Dim0, int Tensor_Dim1>
1797 int side, EntityType type, EntitiesFieldData::EntData &data) {
1799 if (!this->dataPtr)
1800 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1801 "Data pointer not allocated");
1802
1803 const size_t nb_gauss_pts = this->getGaussPts().size2();
1804 auto &mat = *this->dataPtr;
1805 if (type == this->zeroType) {
1806 mat.resize(Tensor_Dim0 * Tensor_Dim1 * Tensor_Dim1, nb_gauss_pts, false);
1807 mat.clear();
1808 }
1809
1810 if (nb_gauss_pts) {
1811 const size_t nb_dofs = data.getFieldData().size();
1812
1813 if (nb_dofs) {
1814
1815 if (this->dataVec.use_count()) {
1816 this->dotVector.resize(nb_dofs, false);
1817 const double *array;
1818 CHKERR VecGetArrayRead(this->dataVec, &array);
1819 const auto &local_indices = data.getLocalIndices();
1820 for (int i = 0; i != local_indices.size(); ++i)
1821 if (local_indices[i] != -1)
1822 this->dotVector[i] = array[local_indices[i]];
1823 else
1824 this->dotVector[i] = 0;
1825 CHKERR VecRestoreArrayRead(this->dataVec, &array);
1826 data.getFieldData().swap(this->dotVector);
1827 }
1828
1829 const int nb_base_functions = data.getN().size2();
1830
1831 auto &hessian_base = data.getN(BaseDerivatives::SecondDerivative);
1832 #ifndef NDEBUG
1833 if (hessian_base.size1() != nb_gauss_pts) {
1834 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1835 "Wrong number of integration pts (%d != %d)",
1836 hessian_base.size1(), nb_gauss_pts);
1837 }
1838 if (hessian_base.size2() !=
1839 nb_base_functions * Tensor_Dim1 * Tensor_Dim1) {
1840 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1841 "Wrong number of base functions (%d != %d)",
1842 hessian_base.size2() / (Tensor_Dim1 * Tensor_Dim1),
1843 nb_base_functions);
1844 }
1845 if (hessian_base.size2() < nb_dofs * Tensor_Dim1 * Tensor_Dim1) {
1846 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1847 "Wrong number of base functions (%d < %d)",
1848 hessian_base.size2(), nb_dofs * Tensor_Dim1 * Tensor_Dim1);
1849 }
1850 #endif
1851
1852 auto t_diff2_base_function = getFTensor2FromPtr<Tensor_Dim1, Tensor_Dim1>(
1853 &*hessian_base.data().begin());
1854
1855 auto t_hessian_at_gauss_pts =
1856 getFTensor3FromMat<Tensor_Dim0, Tensor_Dim1, Tensor_Dim1>(mat);
1857
1858 FTensor::Index<'I', Tensor_Dim0> I;
1859 FTensor::Index<'J', Tensor_Dim1> J;
1860 FTensor::Index<'K', Tensor_Dim1> K;
1861
1862 int size = nb_dofs / Tensor_Dim0;
1863 #ifndef NDEBUG
1864 if (nb_dofs % Tensor_Dim0) {
1865 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1866 "Data inconsistency");
1867 }
1868 #endif
1869
1870 for (int gg = 0; gg < nb_gauss_pts; ++gg) {
1871 auto field_data = data.getFTensor1FieldData<Tensor_Dim0>();
1872 int bb = 0;
1873 for (; bb < size; ++bb) {
1874 t_hessian_at_gauss_pts(I, J, K) +=
1875 field_data(I) * t_diff2_base_function(J, K);
1876 ++field_data;
1877 ++t_diff2_base_function;
1878 }
1879 // Number of dofs can be smaller than number of Tensor_Dim0 x base
1880 // functions
1881 for (; bb != nb_base_functions; ++bb)
1882 ++t_diff2_base_function;
1883 ++t_hessian_at_gauss_pts;
1884 }
1885
1886 if (this->dataVec.use_count()) {
1887 data.getFieldData().swap(this->dotVector);
1888 }
1889 }
1890 }
1892}
1893
1894/**@}*/
1895
1896/** \name Transform tensors and vectors */
1897
1898/**@{*/
1899
1900/**
1901 * @brief Calculate \f$ \pmb\sigma_{ij} = \mathbf{D}_{ijkl} \pmb\varepsilon_{kl}
1902 * \f$
1903 *
1904 * @tparam DIM
1905 *
1906 * \ingroup mofem_forces_and_sources_user_data_operators
1907 */
1908template <int DIM_01, int DIM_23, int S = 0>
1911
1914
1915 /**
1916 * @deprecated Do not use this constructor
1917 */
1920 boost::shared_ptr<MatrixDouble> in_mat,
1921 boost::shared_ptr<MatrixDouble> out_mat,
1922 boost::shared_ptr<MatrixDouble> d_mat)
1923 : UserOp(field_name, OPROW), inMat(in_mat), outMat(out_mat), dMat(d_mat) {
1924 // Only is run for vertices
1925 std::fill(&doEntities[MBEDGE], &doEntities[MBMAXTYPE], false);
1926 if (!inMat)
1927 THROW_MESSAGE("Pointer for in mat is null");
1928 if (!outMat)
1929 THROW_MESSAGE("Pointer for out mat is null");
1930 if (!dMat)
1931 THROW_MESSAGE("Pointer for tensor mat is null");
1932 }
1933
1934 OpTensorTimesSymmetricTensor(boost::shared_ptr<MatrixDouble> in_mat,
1935 boost::shared_ptr<MatrixDouble> out_mat,
1936 boost::shared_ptr<MatrixDouble> d_mat)
1937 : UserOp(NOSPACE, OPSPACE), inMat(in_mat), outMat(out_mat), dMat(d_mat) {
1938 // Only is run for vertices
1939 if (!inMat)
1940 THROW_MESSAGE("Pointer for in mat is null");
1941 if (!outMat)
1942 THROW_MESSAGE("Pointer for out mat is null");
1943 if (!dMat)
1944 THROW_MESSAGE("Pointer for tensor mat is null");
1945 }
1946
1947 MoFEMErrorCode doWork(int side, EntityType type, EntData &data) {
1949 const size_t nb_gauss_pts = getGaussPts().size2();
1950 auto t_D = getFTensor4DdgFromMat<DIM_01, DIM_23, S>(*(dMat));
1951 auto t_in = getFTensor2SymmetricFromMat<DIM_01>(*(inMat));
1952 outMat->resize((DIM_23 * (DIM_23 + 1)) / 2, nb_gauss_pts, false);
1953 auto t_out = getFTensor2SymmetricFromMat<DIM_23>(*(outMat));
1954 for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
1955 t_out(i, j) = t_D(i, j, k, l) * t_in(k, l);
1956 ++t_in;
1957 ++t_out;
1958 }
1960 }
1961
1962private:
1963 FTensor::Index<'i', DIM_01> i;
1964 FTensor::Index<'j', DIM_01> j;
1965 FTensor::Index<'k', DIM_23> k;
1966 FTensor::Index<'l', DIM_23> l;
1967
1968 boost::shared_ptr<MatrixDouble> inMat;
1969 boost::shared_ptr<MatrixDouble> outMat;
1970 boost::shared_ptr<MatrixDouble> dMat;
1971};
1972
1973template <int DIM>
1975
1978
1979 /**
1980 * @deprecated Do not use this constructor
1981 */
1983 boost::shared_ptr<MatrixDouble> in_mat,
1984 boost::shared_ptr<MatrixDouble> out_mat)
1985 : UserOp(field_name, OPROW), inMat(in_mat), outMat(out_mat) {
1986 // Only is run for vertices
1987 std::fill(&doEntities[MBEDGE], &doEntities[MBMAXTYPE], false);
1988 if (!inMat)
1989 THROW_MESSAGE("Pointer not set for in matrix");
1990 if (!outMat)
1991 THROW_MESSAGE("Pointer not set for in matrix");
1992 }
1993
1994 OpSymmetrizeTensor(boost::shared_ptr<MatrixDouble> in_mat,
1995 boost::shared_ptr<MatrixDouble> out_mat)
1996 : UserOp(NOSPACE, OPSPACE), inMat(in_mat), outMat(out_mat) {
1997 // Only is run for vertices
1998 if (!inMat)
1999 THROW_MESSAGE("Pointer not set for in matrix");
2000 if (!outMat)
2001 THROW_MESSAGE("Pointer not set for in matrix");
2002 }
2003
2004 MoFEMErrorCode doWork(int side, EntityType type, EntData &data) {
2006 const size_t nb_gauss_pts = getGaussPts().size2();
2007 auto t_in = getFTensor2FromMat<DIM, DIM>(*(inMat));
2008 outMat->resize((DIM * (DIM + 1)) / 2, nb_gauss_pts, false);
2009 auto t_out = getFTensor2SymmetricFromMat<DIM>(*(outMat));
2010 for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
2011 t_out(i, j) = (t_in(i, j) || t_in(j, i)) / 2;
2012 ++t_in;
2013 ++t_out;
2014 }
2016 }
2017
2018private:
2021 boost::shared_ptr<MatrixDouble> inMat;
2022 boost::shared_ptr<MatrixDouble> outMat;
2023};
2024
2026
2029
2030 /**
2031 * @deprecated Do not use this constructor
2032 */
2033 DEPRECATED OpScaleMatrix(const std::string field_name, const double scale,
2034 boost::shared_ptr<MatrixDouble> in_mat,
2035 boost::shared_ptr<MatrixDouble> out_mat)
2036 : UserOp(field_name, OPROW), inMat(in_mat), outMat(out_mat) {
2037 scalePtr = boost::make_shared<double>(scale);
2038 // Only is run for vertices
2039 std::fill(&doEntities[MBEDGE], &doEntities[MBMAXTYPE], false);
2040 if (!inMat)
2041 THROW_MESSAGE("Pointer not set for in matrix");
2042 if (!outMat)
2043 THROW_MESSAGE("Pointer not set for in matrix");
2044 }
2045
2046 OpScaleMatrix(boost::shared_ptr<double> scale_ptr,
2047 boost::shared_ptr<MatrixDouble> in_mat,
2048 boost::shared_ptr<MatrixDouble> out_mat)
2049 : UserOp(NOSPACE, OPSPACE), scalePtr(scale_ptr), inMat(in_mat),
2050 outMat(out_mat) {
2051 if (!scalePtr)
2052 THROW_MESSAGE("Pointer not set for scale");
2053 if (!inMat)
2054 THROW_MESSAGE("Pointer not set for in matrix");
2055 if (!outMat)
2056 THROW_MESSAGE("Pointer not set for in matrix");
2057 }
2058
2059 MoFEMErrorCode doWork(int side, EntityType type, EntData &data) {
2061 outMat->resize(inMat->size1(), inMat->size2(), false);
2062 noalias(*outMat) = (*scalePtr) * (*inMat);
2064 }
2065
2066private:
2067 boost::shared_ptr<double> scalePtr;
2068 boost::shared_ptr<MatrixDouble> inMat;
2069 boost::shared_ptr<MatrixDouble> outMat;
2070};
2071
2072/**@}*/
2073
2074/** \name H-div/H-curls (Vectorial bases) values at integration points */
2075
2076/**@{*/
2077
2078/** \brief Get vector field for H-div approximation
2079 * \ingroup mofem_forces_and_sources_user_data_operators
2080 */
2081template <int Base_Dim, int Field_Dim, class T, class L, class A>
2083
2084/** \brief Get vector field for H-div approximation
2085 * \ingroup mofem_forces_and_sources_user_data_operators
2086 */
2087template <int Field_Dim>
2089 ublas::row_major, DoubleAllocator>
2091
2093 boost::shared_ptr<MatrixDouble> data_ptr,
2094 SmartPetscObj<Vec> data_vec,
2095 const EntityType zero_type = MBEDGE,
2096 const int zero_side = 0)
2099 dataPtr(data_ptr), dataVec(data_vec), zeroType(zero_type),
2100 zeroSide(zero_side) {
2101 if (!dataPtr)
2102 THROW_MESSAGE("Pointer is not set");
2103 }
2104
2106 boost::shared_ptr<MatrixDouble> data_ptr,
2107 const EntityType zero_type = MBEDGE,
2108 const int zero_side = 0)
2110 field_name, data_ptr, SmartPetscObj<Vec>(), zero_type, zero_side) {}
2111
2112 /**
2113 * \brief Calculate values of vector field at integration points
2114 * @param side side entity number
2115 * @param type side entity type
2116 * @param data entity data
2117 * @return error code
2118 */
2119 MoFEMErrorCode doWork(int side, EntityType type,
2121
2122private:
2123 boost::shared_ptr<MatrixDouble> dataPtr;
2126 const int zeroSide;
2128};
2129
2130template <int Field_Dim>
2132 3, Field_Dim, double, ublas::row_major,
2133 DoubleAllocator>::doWork(int side, EntityType type,
2136 const size_t nb_integration_points = this->getGaussPts().size2();
2137 if (type == zeroType && side == zeroSide) {
2138 dataPtr->resize(Field_Dim, nb_integration_points, false);
2139 dataPtr->clear();
2140 }
2141 const size_t nb_dofs = data.getFieldData().size();
2142 if (!nb_dofs)
2144
2145 if (dataVec.use_count()) {
2146 dotVector.resize(nb_dofs, false);
2147 const double *array;
2148 CHKERR VecGetArrayRead(dataVec, &array);
2149 const auto &local_indices = data.getLocalIndices();
2150 for (int i = 0; i != local_indices.size(); ++i)
2151 if (local_indices[i] != -1)
2152 dotVector[i] = array[local_indices[i]];
2153 else
2154 dotVector[i] = 0;
2155 CHKERR VecRestoreArrayRead(dataVec, &array);
2156 data.getFieldData().swap(dotVector);
2157 }
2158
2159 const size_t nb_base_functions = data.getN().size2() / 3;
2160 FTensor::Index<'i', Field_Dim> i;
2161 auto t_n_hdiv = data.getFTensor1N<3>();
2162 auto t_data = getFTensor1FromMat<Field_Dim>(*dataPtr);
2163 for (size_t gg = 0; gg != nb_integration_points; ++gg) {
2164 auto t_dof = data.getFTensor0FieldData();
2165 int bb = 0;
2166 for (; bb != nb_dofs; ++bb) {
2167 t_data(i) += t_n_hdiv(i) * t_dof;
2168 ++t_n_hdiv;
2169 ++t_dof;
2170 }
2171 for (; bb != nb_base_functions; ++bb)
2172 ++t_n_hdiv;
2173 ++t_data;
2174 }
2175
2176 if (dataVec.use_count()) {
2177 data.getFieldData().swap(dotVector);
2178 }
2180}
2181
2182/** \brief Get vector field for H-div approximation
2183 *
2184 * \ingroup mofem_forces_and_sources_user_data_operators
2185 */
2186template <int Base_Dim, int Field_Dim = Base_Dim>
2189 Base_Dim, Field_Dim, double, ublas::row_major, DoubleAllocator> {
2191 Base_Dim, Field_Dim, double, ublas::row_major,
2192 DoubleAllocator>::OpCalculateHVecVectorField_General;
2193};
2194
2195/** \brief Get vector field for H-div approximation
2196 * \ingroup mofem_forces_and_sources_user_data_operators
2197 */
2198template <int Base_Dim, int Field_Dim = Base_Dim>
2200
2201template <int Field_Dim>
2204
2206 boost::shared_ptr<MatrixDouble> data_ptr,
2207 const EntityType zero_type = MBEDGE,
2208 const int zero_side = 0)
2211 dataPtr(data_ptr), zeroType(zero_type), zeroSide(zero_side) {
2212 if (!dataPtr)
2213 THROW_MESSAGE("Pointer is not set");
2214 }
2215
2216 /**
2217 * \brief Calculate values of vector field at integration points
2218 * @param side side entity number
2219 * @param type side entity type
2220 * @param data entity data
2221 * @return error code
2222 */
2223 MoFEMErrorCode doWork(int side, EntityType type,
2225
2226private:
2227 boost::shared_ptr<MatrixDouble> dataPtr;
2229 const int zeroSide;
2230};
2231
2232template <int Field_Dim>
2234 int side, EntityType type, EntitiesFieldData::EntData &data) {
2236
2237 const size_t nb_integration_points = this->getGaussPts().size2();
2238 if (type == zeroType && side == zeroSide) {
2239 dataPtr->resize(Field_Dim, nb_integration_points, false);
2240 dataPtr->clear();
2241 }
2242
2243 auto &local_indices = data.getIndices();
2244 const size_t nb_dofs = local_indices.size();
2245 if (nb_dofs) {
2246
2247 std::array<double, MAX_DOFS_ON_ENTITY> dot_dofs_vector;
2248 const double *array;
2249 CHKERR VecGetArrayRead(getFEMethod()->ts_u_t, &array);
2250 for (size_t i = 0; i != nb_dofs; ++i)
2251 if (local_indices[i] != -1)
2252 dot_dofs_vector[i] = array[local_indices[i]];
2253 else
2254 dot_dofs_vector[i] = 0;
2255 CHKERR VecRestoreArrayRead(getFEMethod()->ts_u_t, &array);
2256
2257 const size_t nb_base_functions = data.getN().size2() / 3;
2258 FTensor::Index<'i', Field_Dim> i;
2259 auto t_n_hdiv = data.getFTensor1N<3>();
2260 auto t_data = getFTensor1FromMat<Field_Dim>(*dataPtr);
2261 for (size_t gg = 0; gg != nb_integration_points; ++gg) {
2262 int bb = 0;
2263 for (; bb != nb_dofs; ++bb) {
2264 t_data(i) += t_n_hdiv(i) * dot_dofs_vector[bb];
2265 ++t_n_hdiv;
2266 }
2267 for (; bb != nb_base_functions; ++bb)
2268 ++t_n_hdiv;
2269 ++t_data;
2270 }
2271 }
2272
2274}
2275
2276/**
2277 * @brief Calculate divergence of vector field
2278 * @ingroup mofem_forces_and_sources_user_data_operators
2279 *
2280 * @tparam BASE_DIM
2281 * @tparam SPACE_DIM
2282 */
2283template <int BASE_DIM, int SPACE_DIM>
2286
2288 boost::shared_ptr<VectorDouble> data_ptr,
2289 const EntityType zero_type = MBEDGE,
2290 const int zero_side = 0)
2293 dataPtr(data_ptr), zeroType(zero_type), zeroSide(zero_side) {
2294 if (!dataPtr)
2295 THROW_MESSAGE("Pointer is not set");
2296 }
2297
2298 MoFEMErrorCode doWork(int side, EntityType type,
2301 const size_t nb_integration_points = getGaussPts().size2();
2302 if (type == zeroType && side == zeroSide) {
2303 dataPtr->resize(nb_integration_points, false);
2304 dataPtr->clear();
2305 }
2306 const size_t nb_dofs = data.getFieldData().size();
2307 if (!nb_dofs)
2309 const size_t nb_base_functions = data.getN().size2() / BASE_DIM;
2312 auto t_n_diff_hdiv = data.getFTensor2DiffN<BASE_DIM, SPACE_DIM>();
2313 auto t_data = getFTensor0FromVec(*dataPtr);
2314 for (size_t gg = 0; gg != nb_integration_points; ++gg) {
2315 auto t_dof = data.getFTensor0FieldData();
2316 int bb = 0;
2317 for (; bb != nb_dofs; ++bb) {
2318 t_data += t_dof * t_n_diff_hdiv(j, j);
2319 ++t_n_diff_hdiv;
2320 ++t_dof;
2321 }
2322 for (; bb != nb_base_functions; ++bb)
2323 ++t_n_diff_hdiv;
2324 ++t_data;
2325 }
2327 }
2328
2329private:
2330 boost::shared_ptr<VectorDouble> dataPtr;
2332 const int zeroSide;
2333};
2334
2335/**
2336 * @brief Calculate gradient of vector field
2337 * @ingroup mofem_forces_and_sources_user_data_operators
2338 *
2339 * @tparam BASE_DIM
2340 * @tparam SPACE_DIM
2341 */
2342template <int BASE_DIM, int SPACE_DIM>
2345
2347 boost::shared_ptr<MatrixDouble> data_ptr,
2348 const EntityType zero_type = MBEDGE,
2349 const int zero_side = 0)
2352 dataPtr(data_ptr), zeroType(zero_type), zeroSide(zero_side) {
2353 if (!dataPtr)
2354 THROW_MESSAGE("Pointer is not set");
2355 }
2356
2357 MoFEMErrorCode doWork(int side, EntityType type,
2360 const size_t nb_integration_points = getGaussPts().size2();
2361 if (type == zeroType && side == zeroSide) {
2362 dataPtr->resize(BASE_DIM * SPACE_DIM, nb_integration_points, false);
2363 dataPtr->clear();
2364 }
2365 const size_t nb_dofs = data.getFieldData().size();
2366 if (!nb_dofs)
2368 const size_t nb_base_functions = data.getN().size2() / BASE_DIM;
2371 auto t_base_diff = data.getFTensor2DiffN<BASE_DIM, SPACE_DIM>();
2372 auto t_data = getFTensor2FromMat<BASE_DIM, SPACE_DIM>(*dataPtr);
2373 for (size_t gg = 0; gg != nb_integration_points; ++gg) {
2374 auto t_dof = data.getFTensor0FieldData();
2375 int bb = 0;
2376 for (; bb != nb_dofs; ++bb) {
2377 t_data(i, j) += t_dof * t_base_diff(i, j);
2378 ++t_base_diff;
2379 ++t_dof;
2380 }
2381 for (; bb != nb_base_functions; ++bb)
2382 ++t_base_diff;
2383 ++t_data;
2384 }
2386 }
2387
2388private:
2389 boost::shared_ptr<MatrixDouble> dataPtr;
2391 const int zeroSide;
2392};
2393
2394/**
2395 * @brief Calculate gradient of vector field
2396 * @ingroup mofem_forces_and_sources_user_data_operators
2397 *
2398 * @tparam BASE_DIM
2399 * @tparam SPACE_DIM
2400 */
2401template <int BASE_DIM, int SPACE_DIM>
2404
2406 boost::shared_ptr<MatrixDouble> data_ptr,
2407 const EntityType zero_type = MBEDGE,
2408 const int zero_side = 0)
2411 dataPtr(data_ptr), zeroType(zero_type), zeroSide(zero_side) {
2412 if (!dataPtr)
2413 THROW_MESSAGE("Pointer is not set");
2414 }
2415
2416 MoFEMErrorCode doWork(int side, EntityType type,
2419 const size_t nb_integration_points = getGaussPts().size2();
2420 if (type == zeroType && side == zeroSide) {
2421 dataPtr->resize(BASE_DIM * SPACE_DIM * SPACE_DIM, nb_integration_points,
2422 false);
2423 dataPtr->clear();
2424 }
2425 const size_t nb_dofs = data.getFieldData().size();
2426 if (!nb_dofs)
2428
2429 const int nb_base_functions = data.getN().size2() / BASE_DIM;
2430
2431 #ifndef NDEBUG
2432 auto &hessian_base = data.getN(BaseDerivatives::SecondDerivative);
2433 if (hessian_base.size1() != nb_integration_points) {
2434 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
2435 "Wrong number of integration pts (%d != %d)",
2436 hessian_base.size1(), nb_integration_points);
2437 }
2438 if (hessian_base.size2() !=
2439 BASE_DIM * nb_base_functions * SPACE_DIM * SPACE_DIM) {
2440 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
2441 "Wrong number of base functions (%d != %d)",
2442 hessian_base.size2() / (BASE_DIM * SPACE_DIM * SPACE_DIM),
2443 nb_base_functions);
2444 }
2445 if (hessian_base.size2() < BASE_DIM * nb_dofs * SPACE_DIM * SPACE_DIM) {
2446 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
2447 "Wrong number of base functions (%d < %d)", hessian_base.size2(),
2448 BASE_DIM * nb_dofs * SPACE_DIM * SPACE_DIM);
2449 }
2450 #endif
2451
2455
2456 auto t_base_diff2 =
2458 auto t_data = getFTensor3FromMat<BASE_DIM, SPACE_DIM, SPACE_DIM>(*dataPtr);
2459 for (size_t gg = 0; gg != nb_integration_points; ++gg) {
2460 auto t_dof = data.getFTensor0FieldData();
2461 int bb = 0;
2462 for (; bb != nb_dofs; ++bb) {
2463 t_data(i, j, k) += t_dof * t_base_diff2(i, j, k);
2464
2465 ++t_base_diff2;
2466 ++t_dof;
2467 }
2468 for (; bb != nb_base_functions; ++bb)
2469 ++t_base_diff2;
2470 ++t_data;
2471 }
2473 }
2474
2475private:
2476 boost::shared_ptr<MatrixDouble> dataPtr;
2478 const int zeroSide;
2479};
2480
2481/**
2482 * @brief Calculate divergence of vector field dot
2483 * @ingroup mofem_forces_and_sources_user_data_operators
2484 *
2485 * @tparam Tensor_Dim dimension of space
2486 */
2487template <int Tensor_Dim1, int Tensor_Dim2>
2490
2492 boost::shared_ptr<VectorDouble> data_ptr,
2493 const EntityType zero_type = MBEDGE,
2494 const int zero_side = 0)
2497 dataPtr(data_ptr), zeroType(zero_type), zeroSide(zero_side) {
2498 if (!dataPtr)
2499 THROW_MESSAGE("Pointer is not set");
2500 }
2501
2502 MoFEMErrorCode doWork(int side, EntityType type,
2505 const size_t nb_integration_points = getGaussPts().size2();
2506 if (type == zeroType && side == zeroSide) {
2507 dataPtr->resize(nb_integration_points, false);
2508 dataPtr->clear();
2509 }
2510
2511 const auto &local_indices = data.getLocalIndices();
2512 const int nb_dofs = local_indices.size();
2513 if (nb_dofs) {
2514
2515 std::array<double, MAX_DOFS_ON_ENTITY> dot_dofs_vector;
2516 const double *array;
2517 CHKERR VecGetArrayRead(getFEMethod()->ts_u_t, &array);
2518 for (size_t i = 0; i != local_indices.size(); ++i)
2519 if (local_indices[i] != -1)
2520 dot_dofs_vector[i] = array[local_indices[i]];
2521 else
2522 dot_dofs_vector[i] = 0;
2523 CHKERR VecRestoreArrayRead(getFEMethod()->ts_u_t, &array);
2524
2525 const size_t nb_base_functions = data.getN().size2() / Tensor_Dim1;
2526 FTensor::Index<'i', Tensor_Dim1> i;
2527 auto t_n_diff_hdiv = data.getFTensor2DiffN<Tensor_Dim1, Tensor_Dim2>();
2528 auto t_data = getFTensor0FromVec(*dataPtr);
2529 for (size_t gg = 0; gg != nb_integration_points; ++gg) {
2530 int bb = 0;
2531 for (; bb != nb_dofs; ++bb) {
2532 double div = 0;
2533 for (auto ii = 0; ii != Tensor_Dim2; ++ii)
2534 div += t_n_diff_hdiv(ii, ii);
2535 t_data += dot_dofs_vector[bb] * div;
2536 ++t_n_diff_hdiv;
2537 }
2538 for (; bb != nb_base_functions; ++bb)
2539 ++t_n_diff_hdiv;
2540 ++t_data;
2541 }
2542 }
2544 }
2545
2546private:
2547 boost::shared_ptr<VectorDouble> dataPtr;
2549 const int zeroSide;
2550};
2551
2552/**
2553 * @brief Calculate curl of vector field
2554 * @ingroup mofem_forces_and_sources_user_data_operators
2555 *
2556 * @tparam Base_Dim base function dimension
2557 * @tparam Space_Dim dimension of space
2558 * @tparam Hcurl field dimension
2559 */
2560template <int Base_Dim, int Space_Dim> struct OpCalculateHcurlVectorCurl;
2561
2562/**
2563 * @brief Calculate curl of vector field
2564 * @ingroup mofem_forces_and_sources_user_data_operators
2565 *
2566 * @tparam Base_Dim base function dimension
2567 * @tparam Space_Dim dimension of space
2568 * @tparam Hcurl field dimension
2569 */
2570template <>
2573 OpCalculateHcurlVectorCurl(const std::string field_name,
2574 boost::shared_ptr<MatrixDouble> data_ptr,
2575 const EntityType zero_type = MBEDGE,
2576 const int zero_side = 0);
2577 MoFEMErrorCode doWork(int side, EntityType type,
2579
2580private:
2581 boost::shared_ptr<MatrixDouble> dataPtr;
2583 const int zeroSide;
2584};
2585
2586/**
2587 * @brief Calculate curl of vector field
2588 * @ingroup mofem_forces_and_sources_user_data_operators
2589 *
2590 * @tparam Field_Dim dimension of field
2591 * @tparam Space_Dim dimension of space
2592 */
2593template <>
2596
2597 OpCalculateHcurlVectorCurl(const std::string field_name,
2598 boost::shared_ptr<MatrixDouble> data_ptr,
2599 const EntityType zero_type = MBVERTEX,
2600 const int zero_side = 0);
2601
2602 MoFEMErrorCode doWork(int side, EntityType type,
2604
2605private:
2606 boost::shared_ptr<MatrixDouble> dataPtr;
2608 const int zeroSide;
2609};
2610
2611/**
2612 * @brief Calculate curl of vector field
2613 * @ingroup mofem_forces_and_sources_user_data_operators
2614 *
2615 * @tparam Field_Dim dimension of field
2616 * @tparam Space_Dim dimension of space
2617 */
2618template <>
2621
2622 OpCalculateHcurlVectorCurl(const std::string field_name,
2623 boost::shared_ptr<MatrixDouble> data_ptr,
2624 const EntityType zero_type = MBVERTEX,
2625 const int zero_side = 0);
2626
2627 MoFEMErrorCode doWork(int side, EntityType type,
2629
2630private:
2631 boost::shared_ptr<MatrixDouble> dataPtr;
2633 const int zeroSide;
2634};
2635
2636/**
2637 * @brief Calculate tenor field using vectorial base, i.e. Hdiv/Hcurl
2638 * \ingroup mofem_forces_and_sources_user_data_operators
2639 *
2640 * @tparam Tensor_Dim0 rank of the field
2641 * @tparam Tensor_Dim1 dimension of space
2642 */
2643template <int Tensor_Dim0, int Tensor_Dim1>
2646
2648 boost::shared_ptr<MatrixDouble> data_ptr,
2649 boost::shared_ptr<double> scale_ptr,
2651 const EntityType zero_type = MBEDGE,
2652 const int zero_side = 0)
2655 dataPtr(data_ptr), scalePtr(scale_ptr), dataVec(data_vec),
2656 zeroType(zero_type), zeroSide(zero_side) {
2657 if (!dataPtr)
2658 CHK_THROW_MESSAGE(MOFEM_DATA_INCONSISTENCY, "Pointer is not set");
2659 }
2660
2662 boost::shared_ptr<MatrixDouble> data_ptr,
2663 const EntityType zero_type = MBEDGE,
2664 const int zero_side = 0)
2665 : OpCalculateHVecTensorField(field_name, data_ptr, nullptr,
2666 SmartPetscObj<Vec>(), zero_type, zero_side) {
2667 }
2668
2669 MoFEMErrorCode doWork(int side, EntityType type,
2672 const size_t nb_integration_points = getGaussPts().size2();
2673 if (type == zeroType && side == zeroSide) {
2674 dataPtr->resize(Tensor_Dim0 * Tensor_Dim1, nb_integration_points, false);
2675 dataPtr->clear();
2676 }
2677 const size_t nb_dofs = data.getFieldData().size();
2678 if (nb_dofs) {
2679
2680 if (dataVec.use_count()) {
2681 dotVector.resize(nb_dofs, false);
2682 const double *array;
2683 CHKERR VecGetArrayRead(dataVec, &array);
2684 const auto &local_indices = data.getLocalIndices();
2685 for (int i = 0; i != local_indices.size(); ++i)
2686 if (local_indices[i] != -1)
2687 dotVector[i] = array[local_indices[i]];
2688 else
2689 dotVector[i] = 0;
2690 CHKERR VecRestoreArrayRead(dataVec, &array);
2691 data.getFieldData().swap(dotVector);
2692 }
2693
2694 double scale = (scalePtr) ? *scalePtr : 1.0;
2695 const size_t nb_base_functions = data.getN().size2() / 3;
2696 FTensor::Index<'i', Tensor_Dim0> i;
2697 FTensor::Index<'j', Tensor_Dim1> j;
2698 auto t_n_hvec = data.getFTensor1N<3>();
2699 auto t_data = getFTensor2FromMat<Tensor_Dim0, Tensor_Dim1>(*dataPtr);
2700 for (size_t gg = 0; gg != nb_integration_points; ++gg) {
2701 auto t_dof = data.getFTensor1FieldData<Tensor_Dim0>();
2702 size_t bb = 0;
2703 for (; bb != nb_dofs / Tensor_Dim0; ++bb) {
2704 t_data(i, j) += (scale * t_dof(i)) * t_n_hvec(j);
2705 ++t_n_hvec;
2706 ++t_dof;
2707 }
2708 for (; bb < nb_base_functions; ++bb)
2709 ++t_n_hvec;
2710 ++t_data;
2711 }
2712
2713 if (dataVec.use_count()) {
2714 data.getFieldData().swap(dotVector);
2715 }
2716 }
2718 }
2719
2720private:
2721 boost::shared_ptr<MatrixDouble> dataPtr;
2722 boost::shared_ptr<double> scalePtr;
2725 const int zeroSide;
2726 VectorDouble dotVector; ///< Keeps temporary values of time derivatives
2727};
2728
2729/** \brief Get tensor field for H-div approximation
2730 * \ingroup mofem_forces_and_sources_user_data_operators
2731 *
2732 * \warning This operator is not tested
2733 */
2734template <int Tensor_Dim0, int Tensor_Dim1>
2737
2739
2741 const std::string field_name, boost::shared_ptr<MatrixDouble> data_ptr,
2742 SmartPetscObj<Vec> data_vec, EntityType broken_type,
2743 boost::shared_ptr<Range> broken_range_ptr = nullptr,
2744 boost::shared_ptr<double> scale_ptr = nullptr,
2745 const EntityType zero_type = MBEDGE, const int zero_side = 0)
2747 dataPtr(data_ptr), dataVec(data_vec), brokenType(broken_type),
2748 brokenRangePtr(broken_range_ptr), zeroType(zero_type) {
2749 if (!dataPtr)
2750 CHK_THROW_MESSAGE(MOFEM_DATA_INCONSISTENCY, "Pointer is not set");
2751 }
2752
2753 /**
2754 * \brief Calculate values of vector field at integration points
2755 * @param side side entity number
2756 * @param type side entity type
2757 * @param data entity data
2758 * @return error code
2759 */
2760 MoFEMErrorCode doWork(int side, EntityType type,
2762
2763private:
2764 boost::shared_ptr<MatrixDouble> dataPtr;
2766 EntityType brokenType;
2767 boost::shared_ptr<Range> brokenRangePtr;
2768 boost::shared_ptr<double> scalePtr;
2771};
2772
2773template <int Tensor_Dim0, int Tensor_Dim1>
2776 int side, EntityType type, EntitiesFieldData::EntData &data) {
2778 const size_t nb_integration_points = OP::getGaussPts().size2();
2779 if (type == zeroType) {
2780 dataPtr->resize(Tensor_Dim0 * Tensor_Dim1, nb_integration_points, false);
2781 dataPtr->clear();
2782 }
2783 const size_t nb_dofs = data.getFieldData().size();
2784 if (!nb_dofs)
2786
2787 if (dataVec.use_count()) {
2788 dotVector.resize(nb_dofs, false);
2789 const double *array;
2790 CHKERR VecGetArrayRead(dataVec, &array);
2791 const auto &local_indices = data.getLocalIndices();
2792 for (int i = 0; i != local_indices.size(); ++i)
2793 if (local_indices[i] != -1)
2794 dotVector[i] = array[local_indices[i]];
2795 else
2796 dotVector[i] = 0;
2797 CHKERR VecRestoreArrayRead(dataVec, &array);
2798 data.getFieldData().swap(dotVector);
2799 }
2800
2801 /**
2802 * @brief Get side face dofs
2803 *
2804 * Find which base functions on borken space have adjacent given entity type
2805 * and are in the range ptr if given.
2806 *
2807 */
2808 auto get_get_side_face_dofs = [&]() {
2809 auto fe_type = OP::getFEType();
2810
2811 BaseFunction::DofsSideMap &side_dof_map =
2812 data.getFieldEntities()[0]->getDofSideMap().at(fe_type);
2813 std::vector<int> side_face_dofs;
2814 side_face_dofs.reserve(data.getIndices().size() / Tensor_Dim0);
2815
2816 for (
2817
2818 auto it = side_dof_map.get<1>().begin();
2819 it != side_dof_map.get<1>().end(); ++it
2820
2821 ) {
2822 if ((Tensor_Dim0 * it->dof) >= data.getIndices().size()) {
2823 break;
2824 }
2825 if (it->type == brokenType) {
2826 if (brokenRangePtr) {
2827 auto ent = OP::getSideEntity(it->side, brokenType);
2828 if (brokenRangePtr->find(ent) != brokenRangePtr->end()) {
2829 side_face_dofs.push_back(it->dof);
2830 }
2831 } else {
2832 side_face_dofs.push_back(it->dof);
2833 }
2834 }
2835 }
2836
2837 return side_face_dofs;
2838 };
2839
2840 auto side_face_dofs = get_get_side_face_dofs();
2841
2842 FTensor::Index<'i', Tensor_Dim0> i;
2843 FTensor::Index<'j', Tensor_Dim1> j;
2844 auto t_data = getFTensor2FromMat<Tensor_Dim0, Tensor_Dim1>(*dataPtr);
2845 for (size_t gg = 0; gg != nb_integration_points; ++gg) {
2846 for (auto b : side_face_dofs) {
2847 auto t_row_base = data.getFTensor1N<3>(gg, b);
2848 auto t_dof = getFTensor1FromPtr<Tensor_Dim0>(data.getFieldData().data() +
2849 b * Tensor_Dim0);
2850 t_data(i, j) += t_dof(i) * t_row_base(j);
2851 }
2852 ++t_data;
2853 }
2854 *dataPtr *= (scalePtr) ? *scalePtr : 1.0;
2855
2856 if (dataVec.use_count()) {
2857 data.getFieldData().swap(dotVector);
2858 }
2859
2861}
2862
2863/**
2864 * @brief Calculate tenor field using tensor base, i.e. Hdiv/Hcurl
2865 * \ingroup mofem_forces_and_sources_user_data_operators
2866 *
2867 * @tparam Tensor_Dim0 rank of the field
2868 * @tparam Tensor_Dim1 dimension of space
2869 */
2870template <int Tensor_Dim0, int Tensor_Dim1>
2873
2875 const std::string field_name, boost::shared_ptr<MatrixDouble> data_ptr,
2876 boost::shared_ptr<double> scale_ptr,
2878 const EntityType zero_type = MBEDGE, const int zero_side = 0)
2881 dataPtr(data_ptr), scalePtr(scale_ptr), dataVec(data_vec),
2882 zeroType(zero_type), zeroSide(zero_side) {
2883 if (!dataPtr)
2884 THROW_MESSAGE("Pointer is not set");
2885 }
2886
2888 boost::shared_ptr<MatrixDouble> data_ptr,
2889 const EntityType zero_type = MBEDGE,
2890 const int zero_side = 0)
2891 : OpCalculateHTensorTensorField(field_name, data_ptr, nullptr,
2892 SmartPetscObj<Vec>(), zero_type,
2893 zero_side) {}
2894
2895 MoFEMErrorCode doWork(int side, EntityType type,
2898 const size_t nb_integration_points = getGaussPts().size2();
2899 if (type == zeroType && side == zeroSide) {
2900 dataPtr->resize(Tensor_Dim0 * Tensor_Dim1, nb_integration_points, false);
2901 dataPtr->clear();
2902 }
2903 const size_t nb_dofs = data.getFieldData().size();
2904 if (!nb_dofs)
2906
2907 if (dataVec.use_count()) {
2908 dotVector.resize(nb_dofs, false);
2909 const double *array;
2910 CHKERR VecGetArrayRead(dataVec, &array);
2911 const auto &local_indices = data.getLocalIndices();
2912 for (int i = 0; i != local_indices.size(); ++i)
2913 if (local_indices[i] != -1)
2914 dotVector[i] = array[local_indices[i]];
2915 else
2916 dotVector[i] = 0;
2917 CHKERR VecRestoreArrayRead(dataVec, &array);
2918 data.getFieldData().swap(dotVector);
2919 }
2920
2921 double scale = (scalePtr) ? *scalePtr : 1.0;
2922 const size_t nb_base_functions =
2923 data.getN().size2() / (Tensor_Dim0 * Tensor_Dim1);
2924 FTensor::Index<'i', Tensor_Dim0> i;
2925 FTensor::Index<'j', Tensor_Dim1> j;
2926 auto t_n_hten = data.getFTensor2N<Tensor_Dim0, Tensor_Dim1>();
2927 auto t_data = getFTensor2FromMat<Tensor_Dim0, Tensor_Dim1>(*dataPtr);
2928 for (size_t gg = 0; gg != nb_integration_points; ++gg) {
2929 auto t_dof = data.getFTensor0FieldData();
2930 size_t bb = 0;
2931 for (; bb != nb_dofs; ++bb) {
2932 t_data(i, j) += (scale * t_dof) * t_n_hten(i, j);
2933 ++t_n_hten;
2934 ++t_dof;
2935 }
2936 for (; bb < nb_base_functions; ++bb)
2937 ++t_n_hten;
2938 ++t_data;
2939 }
2940
2941 if (dataVec.use_count()) {
2942 data.getFieldData().swap(dotVector);
2943 }
2944
2946 }
2947
2948private:
2949 boost::shared_ptr<MatrixDouble> dataPtr;
2950 boost::shared_ptr<double> scalePtr;
2953 const int zeroSide;
2954 VectorDouble dotVector; ///< Keeps temporary values of time derivatives
2955};
2956
2957/**
2958 * @brief Calculate divergence of tonsorial field using vectorial base
2959 * \ingroup mofem_forces_and_sources_user_data_operators
2960 *
2961 * @tparam Tensor_Dim0 rank of the field
2962 * @tparam Tensor_Dim1 dimension of space
2963 */
2964template <int Tensor_Dim0, int Tensor_Dim1,
2965 CoordinateTypes CoordSys = CARTESIAN>
2968
2970 boost::shared_ptr<MatrixDouble> data_ptr,
2971 SmartPetscObj<Vec> data_vec,
2972 const EntityType zero_type = MBEDGE,
2973 const int zero_side = 0)
2976 dataPtr(data_ptr), dataVec(data_vec), zeroType(zero_type),
2977 zeroSide(zero_side) {
2978 if (!dataPtr)
2979 CHK_THROW_MESSAGE(MOFEM_DATA_INCONSISTENCY, "Pointer is not set");
2980 }
2981
2983 boost::shared_ptr<MatrixDouble> data_ptr,
2984 const EntityType zero_type = MBEDGE,
2985 const int zero_side = 0)
2987 field_name, data_ptr, SmartPetscObj<Vec>(), zero_type, zero_side) {}
2988
2989 MoFEMErrorCode doWork(int side, EntityType type,
2992 const size_t nb_integration_points = getGaussPts().size2();
2993 if (type == zeroType && side == zeroSide) {
2994 dataPtr->resize(Tensor_Dim0, nb_integration_points, false);
2995 dataPtr->clear();
2996 }
2997 const size_t nb_dofs = data.getFieldData().size();
2998 if (nb_dofs) {
2999
3000 if (dataVec.use_count()) {
3001 dotVector.resize(nb_dofs, false);
3002 const double *array;
3003 CHKERR VecGetArrayRead(dataVec, &array);
3004 const auto &local_indices = data.getLocalIndices();
3005 for (int i = 0; i != local_indices.size(); ++i)
3006 if (local_indices[i] != -1)
3007 dotVector[i] = array[local_indices[i]];
3008 else
3009 dotVector[i] = 0;
3010 CHKERR VecRestoreArrayRead(dataVec, &array);
3011 data.getFieldData().swap(dotVector);
3012 }
3013
3014 const size_t nb_base_functions = data.getN().size2() / 3;
3015 FTensor::Index<'i', Tensor_Dim0> i;
3016 FTensor::Index<'j', Tensor_Dim1> j;
3017 auto t_n_diff_hvec = data.getFTensor2DiffN<3, Tensor_Dim1>();
3018 auto t_data = getFTensor1FromMat<Tensor_Dim0>(*dataPtr);
3019 auto t_base = data.getFTensor1N<3>();
3020 auto t_coords = getFTensor1CoordsAtGaussPts();
3021 for (size_t gg = 0; gg != nb_integration_points; ++gg) {
3022 auto t_dof = data.getFTensor1FieldData<Tensor_Dim0>();
3023 size_t bb = 0;
3024 for (; bb != nb_dofs / Tensor_Dim0; ++bb) {
3025 double div = t_n_diff_hvec(j, j);
3026 t_data(i) += t_dof(i) * div;
3027 if constexpr (CoordSys == CYLINDRICAL) {
3028 t_data(i) += t_base(0) * (t_dof(i) / t_coords(0));
3029 }
3030 ++t_n_diff_hvec;
3031 ++t_dof;
3032 ++t_base;
3033 }
3034 for (; bb < nb_base_functions; ++bb) {
3035 ++t_base;
3036 ++t_n_diff_hvec;
3037 }
3038 ++t_data;
3039 ++t_coords;
3040 }
3041
3042 if (dataVec.use_count()) {
3043 data.getFieldData().swap(dotVector);
3044 }
3045 }
3047 }
3048
3049private:
3050 boost::shared_ptr<MatrixDouble> dataPtr;
3053 const int zeroSide;
3054
3055 VectorDouble dotVector; ///< Keeps temporary values of time derivatives
3056};
3057
3058/**
3059 * @brief Calculate divergence of tonsorial field using vectorial base
3060 * \ingroup mofem_forces_and_sources_user_data_operators
3061 *
3062 * \warning This operator is not tested
3063 *
3064 * @tparam Tensor_Dim0 rank of the field
3065 * @tparam Tensor_Dim1 dimension of space
3066 */
3067template <int Tensor_Dim0, int Tensor_Dim1,
3068 CoordinateTypes CoordSys = CARTESIAN>
3071
3073
3075 const std::string field_name, boost::shared_ptr<MatrixDouble> data_ptr,
3076 SmartPetscObj<Vec> data_vec, EntityType broken_type,
3077 boost::shared_ptr<Range> broken_range_ptr = nullptr,
3078 boost::shared_ptr<double> scale_ptr = nullptr,
3079 const EntityType zero_type = MBEDGE)
3081 dataPtr(data_ptr), dataVec(data_vec), brokenType(broken_type),
3082 brokenRangePtr(broken_range_ptr), scalePtr(scale_ptr),
3083 zeroType(zero_type) {
3084 if (!dataPtr)
3085 CHK_THROW_MESSAGE(MOFEM_DATA_INCONSISTENCY, "Pointer is not set");
3086 }
3087
3088 MoFEMErrorCode doWork(int side, EntityType type,
3091 const size_t nb_integration_points = getGaussPts().size2();
3092 if (type == zeroType && side == 0) {
3093 dataPtr->resize(Tensor_Dim0, nb_integration_points, false);
3094 dataPtr->clear();
3095 }
3096
3097 const size_t nb_dofs = data.getFieldData().size();
3098 if (nb_dofs) {
3099
3100 if (dataVec.use_count()) {
3101 dotVector.resize(nb_dofs, false);
3102 const double *array;
3103 CHKERR VecGetArrayRead(dataVec, &array);
3104 const auto &local_indices = data.getLocalIndices();
3105 for (int i = 0; i != local_indices.size(); ++i)
3106 if (local_indices[i] != -1)
3107 dotVector[i] = array[local_indices[i]];
3108 else
3109 dotVector[i] = 0;
3110 CHKERR VecRestoreArrayRead(dataVec, &array);
3111 data.getFieldData().swap(dotVector);
3112 }
3113
3114 /**
3115 * @brief Get side face dofs
3116 *
3117 * Find which base functions on borken space have adjacent given entity
3118 * type and are in the range ptr if given.
3119 *
3120 */
3121 auto get_get_side_face_dofs = [&]() {
3122 auto fe_type = OP::getFEType();
3123
3124 BaseFunction::DofsSideMap &side_dof_map =
3125 data.getFieldEntities()[0]->getDofSideMap().at(fe_type);
3126 std::vector<int> side_face_dofs;
3127 side_face_dofs.reserve(data.getIndices().size() / Tensor_Dim0);
3128
3129 for (
3130
3131 auto it = side_dof_map.get<1>().begin();
3132 it != side_dof_map.get<1>().end(); ++it
3133
3134 ) {
3135 if ((Tensor_Dim0 * it->dof) >= data.getIndices().size()) {
3136 break;
3137 }
3138 if (it->type == brokenType) {
3139 if (brokenRangePtr) {
3140 auto ent = OP::getSideEntity(it->side, brokenType);
3141 if (brokenRangePtr->find(ent) != brokenRangePtr->end()) {
3142 side_face_dofs.push_back(it->dof);
3143 }
3144 } else {
3145 side_face_dofs.push_back(it->dof);
3146 }
3147 }
3148 }
3149
3150 return side_face_dofs;
3151 };
3152
3153 auto side_face_dofs = get_get_side_face_dofs();
3154
3155 FTensor::Index<'i', Tensor_Dim0> i;
3156 FTensor::Index<'j', Tensor_Dim1> j;
3157 auto t_data = getFTensor1FromMat<Tensor_Dim0>(*dataPtr);
3158 auto t_coords = getFTensor1CoordsAtGaussPts();
3159 for (size_t gg = 0; gg != nb_integration_points; ++gg) {
3160 for (auto b : side_face_dofs) {
3161 auto t_dof = getFTensor1FromPtr<Tensor_Dim0>(
3162 data.getFieldData().data() + b * Tensor_Dim0);
3163 auto t_base = data.getFTensor1N<3>(gg, b);
3164 auto t_diff_base = data.getFTensor2DiffN<3, Tensor_Dim1>(gg, b);
3165 double div = t_diff_base(j, j);
3166 t_data(i) += t_dof(i) * div;
3167 if constexpr (CoordSys == CYLINDRICAL) {
3168 t_data(i) += t_base(0) * (t_dof(i) / t_coords(0));
3169 }
3170 }
3171 ++t_data;
3172 ++t_coords;
3173 }
3174 }
3175
3176 if (dataVec.use_count()) {
3177 data.getFieldData().swap(dotVector);
3178 }
3179
3181 }
3182
3183private:
3184 boost::shared_ptr<MatrixDouble> dataPtr;
3186 EntityType brokenType;
3187 boost::shared_ptr<Range> brokenRangePtr;
3188 boost::shared_ptr<double> scalePtr;
3191};
3192
3193/**
3194 * @brief Calculate trace of vector (Hdiv/Hcurl) space
3195 *
3196 * @tparam Tensor_Dim
3197 * @tparam OpBase
3198 */
3199template <int Tensor_Dim, typename OpBase>
3201
3203 boost::shared_ptr<MatrixDouble> data_ptr,
3204 boost::shared_ptr<double> scale_ptr,
3205 const EntityType zero_type = MBEDGE,
3206 const int zero_side = 0)
3207 : OpBase(field_name, OpBase::OPROW), dataPtr(data_ptr),
3208 scalePtr(scale_ptr), zeroType(zero_type), zeroSide(zero_side) {
3209 if (!dataPtr)
3210 THROW_MESSAGE("Pointer is not set");
3211 }
3212
3214 boost::shared_ptr<MatrixDouble> data_ptr,
3215 const EntityType zero_type = MBEDGE,
3216 const int zero_side = 0)
3217 : OpCalculateHVecTensorTrace(field_name, data_ptr, nullptr, zero_type,
3218 zero_side) {}
3219
3220 MoFEMErrorCode doWork(int side, EntityType type,
3223 const size_t nb_integration_points = OpBase::getGaussPts().size2();
3224 if (type == zeroType && side == 0) {
3225 dataPtr->resize(Tensor_Dim, nb_integration_points, false);
3226 dataPtr->clear();
3227 }
3228 const size_t nb_dofs = data.getFieldData().size();
3229 if (nb_dofs) {
3230 double scale_val = (scalePtr) ? *scalePtr : 1.0;
3231 auto t_normal = OpBase::getFTensor1NormalsAtGaussPts();
3232 const size_t nb_base_functions = data.getN().size2() / 3;
3233 auto t_base = data.getFTensor1N<3>();
3234 auto t_data = getFTensor1FromMat<Tensor_Dim>(*dataPtr);
3235 for (size_t gg = 0; gg != nb_integration_points; ++gg) {
3236 FTensor::Tensor1<double, Tensor_Dim> t_normalized_normal;
3237 t_normalized_normal(j) = t_normal(j);
3238 t_normalized_normal.normalize();
3239 auto t_dof = data.getFTensor1FieldData<Tensor_Dim>();
3240 size_t bb = 0;
3241 for (; bb != nb_dofs / Tensor_Dim; ++bb) {
3242 t_data(i) +=
3243 (scale_val * t_dof(i)) * (t_base(j) * t_normalized_normal(j));
3244 ++t_base;
3245 ++t_dof;
3246 }
3247 for (; bb < nb_base_functions; ++bb) {
3248 ++t_base;
3249 }
3250 ++t_data;
3251 ++t_normal;
3252 }
3253 }
3255 }
3256
3257private:
3258 boost::shared_ptr<MatrixDouble> dataPtr;
3259 boost::shared_ptr<double> scalePtr;
3261 const int zeroSide;
3262 FTensor::Index<'i', Tensor_Dim> i;
3263 FTensor::Index<'j', Tensor_Dim> j;
3264};
3265
3266/**@}*/
3267
3268/** \name Other operators */
3269
3270/**@{*/
3271
3272/**@}*/
3273
3274/** \name Operators for faces */
3275
3276/**@{*/
3277
3278/** \brief Transform local reference derivatives of shape functions to global
3279derivatives
3280
3281\ingroup mofem_forces_and_sources_tri_element
3282
3283*/
3284template <int DIM, int DERIVATIVE = 1> struct OpSetInvJacSpaceForFaceImpl;
3285
3288
3290 boost::shared_ptr<MatrixDouble> inv_jac_ptr);
3291
3292protected:
3293 template <int D1, int D2, int J1, int J2>
3296
3297 static_assert(D2 == J2, "Dimension of jacobian and dimension of <out> "
3298 "directive does not match");
3299
3300 size_t nb_functions = diff_n.size2() / D1;
3301 if (nb_functions) {
3302 size_t nb_gauss_pts = diff_n.size1();
3303
3304 #ifndef NDEBUG
3305 if (nb_gauss_pts != getGaussPts().size2())
3306 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
3307 "Wrong number of Gauss Pts");
3308 if (diff_n.size2() % D1)
3309 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
3310 "Number of directives of base functions and D1 dimension does "
3311 "not match");
3312 #endif
3313
3314 diffNinvJac.resize(diff_n.size1(), D2 * nb_functions, false);
3315
3316 FTensor::Index<'i', D2> i;
3317 FTensor::Index<'k', D1> k;
3318 auto t_diff_n = getFTensor1FromPtr<D2>(&*diffNinvJac.data().begin());
3319 auto t_diff_n_ref = getFTensor1FromPtr<D1>(&*diff_n.data().begin());
3320 auto t_inv_jac = getFTensor2FromMat<J1, J2>(*invJacPtr);
3321 for (size_t gg = 0; gg != nb_gauss_pts; ++gg, ++t_inv_jac) {
3322 for (size_t dd = 0; dd != nb_functions; ++dd) {
3323 t_diff_n(i) = t_inv_jac(k, i) * t_diff_n_ref(k);
3324 ++t_diff_n;
3325 ++t_diff_n_ref;
3326 }
3327 }
3328
3329 diff_n.swap(diffNinvJac);
3330 }
3332 }
3333
3334 boost::shared_ptr<MatrixDouble> invJacPtr;
3336};
3337
3338template <>
3341
3343
3344 MoFEMErrorCode doWork(int side, EntityType type,
3346};
3347
3348template <>
3350 : public OpSetInvJacSpaceForFaceImpl<2, 1> {
3351
3352 using OpSetInvJacSpaceForFaceImpl<2, 1>::OpSetInvJacSpaceForFaceImpl;
3353
3354 MoFEMErrorCode doWork(int side, EntityType type,
3356};
3357
3358template <>
3360 : public OpSetInvJacSpaceForFaceImpl<2, 1> {
3361
3362 using OpSetInvJacSpaceForFaceImpl<2, 1>::OpSetInvJacSpaceForFaceImpl;
3363
3364 MoFEMErrorCode doWork(int side, EntityType type,
3366};
3367
3368template <int DERIVARIVE = 1>
3370 : public OpSetInvJacSpaceForFaceImpl<2, DERIVARIVE> {
3371 OpSetInvJacH1ForFace(boost::shared_ptr<MatrixDouble> inv_jac_ptr)
3372 : OpSetInvJacSpaceForFaceImpl<2, DERIVARIVE>(H1, inv_jac_ptr) {}
3373};
3374
3375template <int DERIVARIVE = 1>
3377 : public OpSetInvJacSpaceForFaceImpl<2, DERIVARIVE> {
3378 OpSetInvJacL2ForFace(boost::shared_ptr<MatrixDouble> inv_jac_ptr)
3379 : OpSetInvJacSpaceForFaceImpl<2, DERIVARIVE>(L2, inv_jac_ptr) {}
3380};
3381
3382template <int DERIVARIVE = 1>
3384 : public OpSetInvJacSpaceForFaceImpl<3, DERIVARIVE> {
3386 boost::shared_ptr<MatrixDouble> inv_jac_ptr)
3387 : OpSetInvJacSpaceForFaceImpl<3, DERIVARIVE>(H1, inv_jac_ptr) {}
3388};
3389
3390template <int DERIVARIVE = 1>
3392 : public OpSetInvJacSpaceForFaceImpl<3, DERIVARIVE> {
3394 boost::shared_ptr<MatrixDouble> inv_jac_ptr)
3395 : OpSetInvJacSpaceForFaceImpl<3, DERIVARIVE>(L2, inv_jac_ptr) {}
3396};
3397
3398/**
3399 * \brief Transform local reference derivatives of shape function to
3400 global derivatives for face
3401
3402 * \ingroup mofem_forces_and_sources_tri_element
3403 */
3404template <int DIM> struct OpSetInvJacHcurlFaceImpl;
3405
3406template <>
3409
3410 OpSetInvJacHcurlFaceImpl(boost::shared_ptr<MatrixDouble> inv_jac_ptr)
3412 invJacPtr(inv_jac_ptr) {}
3413
3414 MoFEMErrorCode doWork(int side, EntityType type,
3416
3417protected:
3418 boost::shared_ptr<MatrixDouble> invJacPtr;
3420};
3421
3422template <>
3424 using OpSetInvJacHcurlFaceImpl<2>::OpSetInvJacHcurlFaceImpl;
3425 MoFEMErrorCode doWork(int side, EntityType type,
3427};
3428
3431
3432/**
3433 * @brief Make Hdiv space from Hcurl space in 2d
3434 * @ingroup mofem_forces_and_sources_tri_element
3435 */
3445
3446/** \brief Transform Hcurl base fluxes from reference element to physical
3447 * triangle
3448 */
3450
3451/** \brief Apply contravariant (Piola) transfer to Hdiv space on face
3452 *
3453 * Covariant Piola transformation
3454 \f[
3455 \psi_i|_t = J^{-1}_{ij}\hat{\psi}_j\\
3456 \left.\frac{\partial \psi_i}{\partial \xi_j}\right|_t
3457 = J^{-1}_{ik}\frac{\partial \hat{\psi}_k}{\partial \xi_j}
3458 \f]
3459
3460 */
3461template <>
3464
3466 boost::shared_ptr<MatrixDouble> inv_jac_ptr);
3467
3468 MoFEMErrorCode doWork(int side, EntityType type,
3470
3471private:
3472 boost::shared_ptr<MatrixDouble> invJacPtr;
3473
3476};
3477
3480
3481/** \brief Apply contravariant (Piola) transfer to Hdiv space on face
3482 *
3483 * \note Hdiv space is generated by Hcurl space in 2d.
3484 *
3485 * Contravariant Piola transformation
3486 * \f[
3487 * \psi_i|_t = \frac{1}{\textrm{det}(J)}J_{ij}\hat{\psi}_j\\
3488 * \left.\frac{\partial \psi_i}{\partial \xi_j}\right|_t
3489 * =
3490 * \frac{1}{\textrm{det}(J)}J_{ik}\frac{\partial \hat{\psi}_k}{\partial \xi_j}
3491 * \f]
3492 *
3493 * \ingroup mofem_forces_and_sources
3494 *
3495 */
3497
3498template <>
3501
3503 boost::shared_ptr<MatrixDouble> jac_ptr)
3505 jacPtr(jac_ptr) {}
3506
3507 MoFEMErrorCode doWork(int side, EntityType type,
3509
3510protected:
3511 boost::shared_ptr<MatrixDouble> jacPtr;
3514};
3515
3516template <>
3520 2>::OpSetContravariantPiolaTransformOnFace2DImpl;
3521
3522 MoFEMErrorCode doWork(int side, EntityType type,
3524};
3525
3530
3531/**@}*/
3532
3533/** \name Operators for edges */
3534
3535/**@{*/
3536
3546
3547/**
3548 * @deprecated Name is deprecated and this is added for back compatibility
3549 */
3552
3553/**@}*/
3554
3555/** \name Operator for fat prisms */
3556
3557/**@{*/
3558
3559/**
3560 * @brief Operator for fat prism element updating integration weights in the
3561 * volume.
3562 *
3563 * Jacobian on the distorted element is nonconstant. This operator updates
3564 * integration weight on prism to take into account nonconstat jacobian.
3565 *
3566 * \f[
3567 * W_i = w_i \left( \frac{1}{2V} \left\| \frac{\partial \mathbf{x}}{\partial
3568 * \pmb\xi} \right\| \right)
3569 * \f]
3570 * where \f$w_i\f$ is integration weight at integration point \f$i\f$,
3571 * \f$\mathbf{x}\f$ is physical coordinate, and \f$\pmb\xi\f$ is reference
3572 * element coordinate.
3573 *
3574 */
3584
3585/** \brief Calculate inverse of jacobian for face element
3586
3587 It is assumed that face element is XY plane. Applied
3588 only for 2d problems.
3589
3590 FIXME Generalize function for arbitrary face orientation in 3d space
3591 FIXME Calculate to Jacobins for two faces
3592
3593 \ingroup mofem_forces_and_sources_prism_element
3594
3595*/
3598
3599 OpCalculateInvJacForFatPrism(boost::shared_ptr<MatrixDouble> inv_jac_ptr)
3601 invJacPtr(inv_jac_ptr), invJac(*invJacPtr) {}
3602
3606 MoFEMErrorCode doWork(int side, EntityType type,
3608
3609private:
3610 const boost::shared_ptr<MatrixDouble> invJacPtr;
3612};
3613
3614/** \brief Transform local reference derivatives of shape functions to global
3615derivatives
3616
3617FIXME Generalize to curved shapes
3618FIXME Generalize to case that top and bottom face has different shape
3619
3620\ingroup mofem_forces_and_sources_prism_element
3621
3622*/
3625
3626 OpSetInvJacH1ForFatPrism(boost::shared_ptr<MatrixDouble> inv_jac_ptr)
3628 invJacPtr(inv_jac_ptr), invJac(*invJacPtr) {}
3629
3633
3634 MoFEMErrorCode doWork(int side, EntityType type,
3636
3637private:
3638 const boost::shared_ptr<MatrixDouble> invJacPtr;
3641};
3642
3643// Flat prism
3644
3645/** \brief Calculate inverse of jacobian for face element
3646
3647 It is assumed that face element is XY plane. Applied
3648 only for 2d problems.
3649
3650 FIXME Generalize function for arbitrary face orientation in 3d space
3651 FIXME Calculate to Jacobins for two faces
3652
3653 \ingroup mofem_forces_and_sources_prism_element
3654
3655*/
3668
3669/** \brief Transform local reference derivatives of shape functions to global
3670derivatives
3671
3672FIXME Generalize to curved shapes
3673FIXME Generalize to case that top and bottom face has different shape
3674
3675\ingroup mofem_forces_and_sources_prism_element
3676
3677*/
3692
3693/**@}*/
3694
3695/** \name Operation on matrices at integration points */
3696
3697/**@{*/
3698
3699template <int DIM>
3701
3702 OpInvertMatrix(boost::shared_ptr<MatrixDouble> in_ptr,
3703 boost::shared_ptr<VectorDouble> det_ptr,
3704 boost::shared_ptr<MatrixDouble> out_ptr)
3706 outPtr(out_ptr), detPtr(det_ptr) {}
3707
3708 MoFEMErrorCode doWork(int side, EntityType type,
3710
3711private:
3712 boost::shared_ptr<MatrixDouble> inPtr;
3713 boost::shared_ptr<MatrixDouble> outPtr;
3714 boost::shared_ptr<VectorDouble> detPtr;
3715};
3716
3717template <int DIM>
3721
3722 if (!inPtr)
3723 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
3724 "Pointer for inPtr matrix not allocated");
3725 if (!detPtr)
3726 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
3727 "Pointer for detPtr matrix not allocated");
3728
3729 const auto nb_rows = inPtr->size1();
3730 const auto nb_integration_pts = inPtr->size2();
3731
3732 // Calculate determinant
3733 {
3734 detPtr->resize(nb_integration_pts, false);
3735 auto t_in = getFTensor2FromMat<DIM, DIM>(*inPtr);
3736 auto t_det = getFTensor0FromVec(*detPtr);
3737 for (size_t gg = 0; gg != nb_integration_pts; ++gg) {
3738 t_det = determinantTensor(t_in);
3739 ++t_in;
3740 ++t_det;
3741 }
3742 }
3743
3744 // Invert jacobian
3745 if (outPtr) {
3746 outPtr->resize(nb_rows, nb_integration_pts, false);
3747 auto t_in = getFTensor2FromMat<DIM, DIM>(*inPtr);
3748 auto t_out = getFTensor2FromMat<DIM, DIM>(*outPtr);
3749 auto t_det = getFTensor0FromVec(*detPtr);
3750 for (size_t gg = 0; gg != nb_integration_pts; ++gg) {
3751 CHKERR invertTensor(t_in, t_det, t_out);
3752 ++t_in;
3753 ++t_out;
3754 ++t_det;
3755 }
3756 }
3757
3759}
3760
3761/**@}*/
3762
3763/** \brief Calculates the trace of an input matrix
3764
3765\ingroup mofem_forces_and_sources
3766
3767*/
3768
3769template <int DIM>
3771
3772 OpCalculateTraceFromMat(boost::shared_ptr<MatrixDouble> in_ptr,
3773 boost::shared_ptr<VectorDouble> out_ptr)
3775 outPtr(out_ptr) {}
3776
3777 MoFEMErrorCode doWork(int side, EntityType type,
3779
3780private:
3782 boost::shared_ptr<MatrixDouble> inPtr;
3783 boost::shared_ptr<VectorDouble> outPtr;
3784};
3785
3786template <int DIM>
3791
3792 if (!inPtr)
3793 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
3794 "Pointer for inPtr matrix not allocated");
3795
3796 const auto nb_integration_pts = inPtr->size2();
3797 // Invert jacobian
3798 if (outPtr) {
3799 outPtr->resize(nb_integration_pts, false);
3800 auto t_in = getFTensor2FromMat<DIM, DIM>(*inPtr);
3801 auto t_out = getFTensor0FromVec(*outPtr);
3802
3803 for (size_t gg = 0; gg != nb_integration_pts; ++gg) {
3804 t_out = t_in(i, i);
3805 ++t_in;
3806 ++t_out;
3807 }
3808 }
3809
3811}
3812
3813/**@}*/
3814
3815/** \brief Calculates the trace of an input matrix
3816
3817\ingroup mofem_forces_and_sources
3818
3819*/
3820
3821template <int DIM>
3824
3825 OpCalculateTraceFromSymmMat(boost::shared_ptr<MatrixDouble> in_ptr,
3826 boost::shared_ptr<VectorDouble> out_ptr)
3828 outPtr(out_ptr) {}
3829
3830 MoFEMErrorCode doWork(int side, EntityType type,
3832
3833private:
3835 boost::shared_ptr<MatrixDouble> inPtr;
3836 boost::shared_ptr<VectorDouble> outPtr;
3837};
3838
3839template <int DIM>
3844
3845 if (!inPtr)
3846 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
3847 "Pointer for inPtr matrix not allocated");
3848
3849 const auto nb_integration_pts = inPtr->size2();
3850 // Invert jacobian
3851 if (outPtr) {
3852 outPtr->resize(nb_integration_pts, false);
3853 auto t_in = getFTensor2SymmetricFromMat<DIM>(*inPtr);
3854 auto t_out = getFTensor0FromVec(*outPtr);
3855
3856 for (size_t gg = 0; gg != nb_integration_pts; ++gg) {
3857 t_out = t_in(i, i);
3858 ++t_in;
3859 ++t_out;
3860 }
3861 }
3862
3864}
3865
3866} // namespace MoFEM
3867
3868#endif // __USER_DATA_OPERATORS_HPP__
3869
3870/**
3871 * \defgroup mofem_forces_and_sources_user_data_operators Users Operators
3872 *
3873 * \brief Classes and functions used to evaluate fields at integration pts,
3874 *jacobians, etc..
3875 *
3876 * \ingroup mofem_forces_and_sources
3877 **/
ForcesAndSourcesCore::UserDataOperator UserDataOperator
constexpr int SPACE_DIM
Tensor1< T, Tensor_Dim > normalize()
static const char *const CoordinateTypesNames[]
Coordinate system names.
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
FieldSpace
approximation spaces
Definition definitions.h:82
@ L2
field with C-1 continuity
Definition definitions.h:88
@ H1
continuous field
Definition definitions.h:85
@ NOSPACE
Definition definitions.h:83
@ HCURL
field with continuous tangents
Definition definitions.h:86
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ MOFEM_DATA_INCONSISTENCY
Definition definitions.h:31
@ MOFEM_NOT_IMPLEMENTED
Definition definitions.h:32
CoordinateTypes
Coodinate system.
@ CYLINDRICAL
@ POLAR
@ CARTESIAN
@ SPHERICAL
#define DEPRECATED
Definition definitions.h:17
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
#define THROW_MESSAGE(msg)
Throw MoFEM exception.
constexpr int BASE_DIM
#define MOFEM_LOG(channel, severity)
Log.
#define MOFEM_LOG_CHANNEL(channel)
Set and reset channel.
FTensor::Index< 'i', SPACE_DIM > i
FTensor::Index< 'J', DIM1 > J
Definition level_set.cpp:30
FTensor::Index< 'j', 3 > j
FTensor::Index< 'k', 3 > k
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
VecAllocator< double > DoubleAllocator
Definition Types.hpp:62
implementation of Data Operators for Forces and Sources
Definition Common.hpp:10
static MoFEMErrorCode invertTensor(FTensor::Tensor2< T1, DIM, DIM > &t, T2 &det, FTensor::Tensor2< T3, DIM, DIM > &inv_t)
static auto determinantTensor(FTensor::Tensor2< T, DIM, DIM > &t)
Calculate the determinant of a tensor of rank DIM.
static auto getFTensor0FromVec(ublas::vector< T, A > &data)
Get tensor rank 0 (scalar) form data vector.
OpCalculateScalarFieldValuesFromPetscVecImpl< PetscData::CTX_SET_X_T > OpCalculateScalarFieldValuesDot
constexpr IntegrationType I
constexpr AssemblyType A
double scale
Definition plastic.cpp:123
constexpr auto field_name
multi_index_container< DofsSideMapData, indexed_by< ordered_non_unique< tag< TypeSide_mi_tag >, composite_key< DofsSideMapData, member< DofsSideMapData, EntityType, &DofsSideMapData::type >, member< DofsSideMapData, int, &DofsSideMapData::side > > >, ordered_unique< tag< EntDofIdx_mi_tag >, member< DofsSideMapData, int, &DofsSideMapData::dof > > > > DofsSideMap
Map entity stype and side to element/entity dof index.
std::array< bool, MBMAXTYPE > doEntities
If true operator is executed for entity.
Data on single entity (This is passed as argument to DataOperator::doWork)
FTensor::Tensor2< FTensor::PackPtr< double *, Tensor_Dim0 *Tensor_Dim1 >, Tensor_Dim0, Tensor_Dim1 > getFTensor2DiffN(FieldApproximationBase base)
Get derivatives of base functions for Hdiv space.
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1DiffN(const FieldApproximationBase base)
Get derivatives of base functions.
FTensor::Tensor2< FTensor::PackPtr< double *, Tensor_Dim0 *Tensor_Dim1 >, Tensor_Dim0, Tensor_Dim1 > getFTensor2N(FieldApproximationBase base)
Get base functions for Hdiv/Hcurl spaces.
FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > getFTensor0N(const FieldApproximationBase base)
Get base function as Tensor0.
MatrixDouble & getDiffN(const FieldApproximationBase base)
get derivatives of base functions
const VectorFieldEntities & getFieldEntities() const
get field entities
MatrixDouble & getN(const FieldApproximationBase base)
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....
const VectorDouble & getFieldData() const
get dofs values
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.
const VectorInt & getLocalIndices() const
get local indices of dofs on entity
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1FieldData()
Return FTensor of rank 1, i.e. vector from field data coefficients.
FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > getFTensor0FieldData()
Resturn scalar files as a FTensor of rank 0.
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1N(FieldApproximationBase base)
Get base functions for Hdiv/Hcurl spaces.
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::Tensor2< FTensor::PackPtr< double *, Tensor_Dim0 *Tensor_Dim1 >, Tensor_Dim0, Tensor_Dim1 > getFTensor2FieldData()
Return FTensor rank 2, i.e. matrix from field data coefficients.
const VectorInt & getIndices() const
Get global indices of dofs on entity.
EntityType getFEType() const
Get dimension of finite element.
auto getFTensor1CoordsAtGaussPts()
Get coordinates at integration points assuming linear geometry.
EntityHandle getSideEntity(const int side_number, const EntityType type)
Get the side entity.
@ OPROW
operator doWork function is executed on FE rows
@ OPSPACE
operator do Work is execute on space data
const FEMethod * getFEMethod() const
Return raw pointer to Finite Element Method object.
MatrixDouble & getGaussPts()
matrix of integration (Gauss) points for Volume Element
structure to get information form mofem into EntitiesFieldData
Calculate divergence of tonsorial field using vectorial base.
OpCalculateBrokenHVecTensorDivergence(const std::string field_name, boost::shared_ptr< MatrixDouble > data_ptr, SmartPetscObj< Vec > data_vec, EntityType broken_type, boost::shared_ptr< Range > broken_range_ptr=nullptr, boost::shared_ptr< double > scale_ptr=nullptr, const EntityType zero_type=MBEDGE)
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
Get tensor field for H-div approximation.
boost::shared_ptr< MatrixDouble > dataPtr
OpCalculateBrokenHVecTensorField(const std::string field_name, boost::shared_ptr< MatrixDouble > data_ptr, SmartPetscObj< Vec > data_vec, EntityType broken_type, boost::shared_ptr< Range > broken_range_ptr=nullptr, boost::shared_ptr< double > scale_ptr=nullptr, const EntityType zero_type=MBEDGE, const int zero_side=0)
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Calculate values of vector field at integration points.
Calculate field values (template specialization) for tensor field rank 1, i.e. vector field.
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
OpCalculateDivergenceVectorFieldValues(const std::string field_name, boost::shared_ptr< VectorDouble > data_ptr, const EntityType zero_type=MBVERTEX)
Calculate tenor field using tensor base, i.e. Hdiv/Hcurl.
VectorDouble dotVector
Keeps temporary values of time derivatives.
OpCalculateHTensorTensorField(const std::string field_name, boost::shared_ptr< MatrixDouble > data_ptr, boost::shared_ptr< double > scale_ptr, SmartPetscObj< Vec > data_vec=SmartPetscObj< Vec >(), const EntityType zero_type=MBEDGE, const int zero_side=0)
boost::shared_ptr< MatrixDouble > dataPtr
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
OpCalculateHTensorTensorField(const std::string field_name, boost::shared_ptr< MatrixDouble > data_ptr, const EntityType zero_type=MBEDGE, const int zero_side=0)
Calculate divergence of tonsorial field using vectorial base.
boost::shared_ptr< MatrixDouble > dataPtr
VectorDouble dotVector
Keeps temporary values of time derivatives.
OpCalculateHVecTensorDivergence(const std::string field_name, boost::shared_ptr< MatrixDouble > data_ptr, SmartPetscObj< Vec > data_vec, const EntityType zero_type=MBEDGE, const int zero_side=0)
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
OpCalculateHVecTensorDivergence(const std::string field_name, boost::shared_ptr< MatrixDouble > data_ptr, const EntityType zero_type=MBEDGE, const int zero_side=0)
Calculate tenor field using vectorial base, i.e. Hdiv/Hcurl.
OpCalculateHVecTensorField(const std::string field_name, boost::shared_ptr< MatrixDouble > data_ptr, boost::shared_ptr< double > scale_ptr, SmartPetscObj< Vec > data_vec=SmartPetscObj< Vec >(), const EntityType zero_type=MBEDGE, const int zero_side=0)
boost::shared_ptr< MatrixDouble > dataPtr
OpCalculateHVecTensorField(const std::string field_name, boost::shared_ptr< MatrixDouble > data_ptr, const EntityType zero_type=MBEDGE, const int zero_side=0)
VectorDouble dotVector
Keeps temporary values of time derivatives.
boost::shared_ptr< double > scalePtr
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
Calculate trace of vector (Hdiv/Hcurl) space.
boost::shared_ptr< MatrixDouble > dataPtr
boost::shared_ptr< double > scalePtr
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
OpCalculateHVecTensorTrace(const std::string field_name, boost::shared_ptr< MatrixDouble > data_ptr, const EntityType zero_type=MBEDGE, const int zero_side=0)
FTensor::Index< 'j', Tensor_Dim > j
FTensor::Index< 'i', Tensor_Dim > i
OpCalculateHVecTensorTrace(const std::string field_name, boost::shared_ptr< MatrixDouble > data_ptr, boost::shared_ptr< double > scale_ptr, const EntityType zero_type=MBEDGE, const int zero_side=0)
OpCalculateHVecVectorFieldDot(const std::string field_name, boost::shared_ptr< MatrixDouble > data_ptr, const EntityType zero_type=MBEDGE, const int zero_side=0)
Get vector field for H-div approximation.
OpCalculateHVecVectorField_General(const std::string field_name, boost::shared_ptr< MatrixDouble > data_ptr, SmartPetscObj< Vec > data_vec, const EntityType zero_type=MBEDGE, const int zero_side=0)
OpCalculateHVecVectorField_General(const std::string field_name, boost::shared_ptr< MatrixDouble > data_ptr, const EntityType zero_type=MBEDGE, const int zero_side=0)
Get vector field for H-div approximation.
Get vector field for H-div approximation.
Calculate gradient of vector field.
boost::shared_ptr< MatrixDouble > dataPtr
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
OpCalculateHVecVectorGradient(const std::string field_name, boost::shared_ptr< MatrixDouble > data_ptr, const EntityType zero_type=MBEDGE, const int zero_side=0)
Calculate gradient of vector field.
boost::shared_ptr< MatrixDouble > dataPtr
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
OpCalculateHVecVectorHessian(const std::string field_name, boost::shared_ptr< MatrixDouble > data_ptr, const EntityType zero_type=MBEDGE, const int zero_side=0)
boost::shared_ptr< MatrixDouble > dataPtr
boost::shared_ptr< MatrixDouble > dataPtr
boost::shared_ptr< MatrixDouble > dataPtr
Calculate curl of vector field.
Calculate divergence of vector field dot.
boost::shared_ptr< VectorDouble > dataPtr
OpCalculateHdivVectorDivergenceDot(const std::string field_name, boost::shared_ptr< VectorDouble > data_ptr, const EntityType zero_type=MBEDGE, const int zero_side=0)
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
Calculate divergence of vector field.
boost::shared_ptr< VectorDouble > dataPtr
OpCalculateHdivVectorDivergence(const std::string field_name, boost::shared_ptr< VectorDouble > data_ptr, const EntityType zero_type=MBEDGE, const int zero_side=0)
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
Calculate inverse of jacobian for face element.
OpCalculateInvJacForFatPrism(MatrixDouble &inv_jac)
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
OpCalculateInvJacForFatPrism(boost::shared_ptr< MatrixDouble > inv_jac_ptr)
const boost::shared_ptr< MatrixDouble > invJacPtr
Calculate inverse of jacobian for face element.
OpCalculateInvJacForFlatPrism(MatrixDouble &inv_jac_f3)
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
Evaluate field gradient values for scalar field, i.e. gradient is tensor rank 1 (vector)
Get field gradients at integration pts for scalar field rank 0, i.e. vector field.
Evaluate field gradient values for scalar field, i.e. gradient is tensor rank 1 (vector),...
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
calculate gradient values of scalar field at integration points
Get rate of scalar field at integration points.
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
OpCalculateScalarFieldValuesFromPetscVecImpl(const std::string field_name, boost::shared_ptr< VectorDouble > data_ptr, const EntityType zero_at_type=MBVERTEX)
Scalar field values at integration points.
OpCalculateScalarFieldValues_General(const std::string field_name, boost::shared_ptr< ublas::vector< T, A > > data_ptr, SmartPetscObj< Vec > data_vec, const EntityType zero_type=MBVERTEX)
boost::shared_ptr< ublas::vector< T, A > > dataPtr
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
calculate values of scalar field at integration points
OpCalculateScalarFieldValues_General(const std::string field_name, boost::shared_ptr< ublas::vector< T, A > > data_ptr, const EntityType zero_type=MBVERTEX)
Get value at integration points for scalar field.
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
calculate values of scalar field at integration points
Get time direvarive values at integration pts for tensor field rank 2, i.e. matrix field.
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
OpCalculateTensor2FieldValuesDot(const std::string field_name, boost::shared_ptr< MatrixDouble > data_ptr, const EntityType zero_at_type=MBVERTEX)
EntityType zeroAtType
Zero values at Gauss point at this type.
VectorDouble dotVector
Keeps temporary values of time derivatives.
boost::shared_ptr< MatrixDouble > dataPtr
Data computed into this matrix.
OpCalculateTensor2FieldValues_General(const std::string field_name, boost::shared_ptr< ublas::matrix< double, ublas::row_major, DoubleAllocator > > data_ptr, SmartPetscObj< Vec > data_vec, const EntityType zero_type=MBVERTEX)
OpCalculateTensor2FieldValues_General(const std::string field_name, boost::shared_ptr< ublas::matrix< double, ublas::row_major, DoubleAllocator > > data_ptr, const EntityType zero_type=MBVERTEX)
Calculate field values for tenor field rank 2.
OpCalculateTensor2FieldValues_General(const std::string field_name, boost::shared_ptr< ublas::matrix< T, L, A > > data_ptr, const EntityType zero_type=MBVERTEX)
boost::shared_ptr< ublas::matrix< T, L, A > > dataPtr
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
Get values at integration pts for tensor field rank 2, i.e. matrix field.
OpCalculateTensor2SymmetricFieldGradient_General(const std::string field_name, boost::shared_ptr< MatrixDouble > data_ptr, const EntityType zero_type=MBVERTEX)
Evaluate field gradient values for symmetric 2nd order tensor field, i.e. gradient is tensor rank 3.
OpCalculateTensor2SymmetricFieldGradient_General(const std::string field_name, boost::shared_ptr< MatrixDouble > data_ptr, const EntityType zero_type=MBVERTEX)
Get field gradients at integration pts for symmetric tensorial field rank 2.
OpCalculateTensor2SymmetricFieldGradient(const std::string field_name, boost::shared_ptr< MatrixDouble > data_ptr, const EntityType zero_type=MBVERTEX)
Calculate symmetric tensor field rates ant integratio pts.
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
OpCalculateTensor2SymmetricFieldValuesDot(const std::string field_name, boost::shared_ptr< MatrixDouble > data_ptr, const EntityType zero_type=MBEDGE, const int zero_side=0)
Calculate symmetric tensor field values at integration pts.
OpCalculateTensor2SymmetricFieldValues(const std::string field_name, boost::shared_ptr< MatrixDouble > data_ptr, SmartPetscObj< Vec > data_vec, const EntityType zero_type=MBEDGE, const int zero_side=0)
OpCalculateTensor2SymmetricFieldValues(const std::string field_name, boost::shared_ptr< MatrixDouble > data_ptr, const EntityType zero_type=MBEDGE, const int zero_side=0)
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
Calculates the trace of an input matrix.
boost::shared_ptr< MatrixDouble > inPtr
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
boost::shared_ptr< VectorDouble > outPtr
OpCalculateTraceFromMat(boost::shared_ptr< MatrixDouble > in_ptr, boost::shared_ptr< VectorDouble > out_ptr)
Calculates the trace of an input matrix.
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
boost::shared_ptr< VectorDouble > outPtr
boost::shared_ptr< MatrixDouble > inPtr
OpCalculateTraceFromSymmMat(boost::shared_ptr< MatrixDouble > in_ptr, boost::shared_ptr< VectorDouble > out_ptr)
Get field gradients time derivative at integration pts for scalar field rank 0, i....
boost::shared_ptr< MatrixDouble > dataPtr
Data computed into this matrix.
OpCalculateVectorFieldGradientDot(const std::string field_name, boost::shared_ptr< MatrixDouble > data_ptr, const EntityType zero_at_type=MBVERTEX)
EntityType zeroAtType
Zero values at Gauss point at this type.
VectorDouble dotVector
Keeps temporary values of time derivatives.
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
Evaluate field gradient values for vector field, i.e. gradient is tensor rank 2.
Get field gradients at integration pts for scalar field rank 0, i.e. vector field.
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
calculate values of vector field at integration points
Approximate field values for given petsc vector.
OpCalculateVectorFieldValuesFromPetscVecImpl(const std::string field_name, boost::shared_ptr< MatrixDouble > data_ptr, const EntityType zero_at_type=MBVERTEX, bool throw_error=true)
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
OpCalculateVectorFieldValues_General(const std::string field_name, boost::shared_ptr< MatrixDouble > data_ptr, const EntityType zero_type=MBVERTEX)
OpCalculateVectorFieldValues_General(const std::string field_name, boost::shared_ptr< MatrixDouble > data_ptr, SmartPetscObj< Vec > data_vec, const EntityType zero_type=MBVERTEX)
Calculate field values for tenor field rank 1, i.e. vector field.
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
calculate values of vector field at integration points
OpCalculateVectorFieldValues_General(const std::string field_name, boost::shared_ptr< ublas::matrix< T, L, A > > data_ptr, const EntityType zero_type=MBVERTEX)
boost::shared_ptr< ublas::matrix< T, L, A > > dataPtr
Get values at integration pts for tensor field rank 1, i.e. vector field.
boost::shared_ptr< VectorDouble > detPtr
OpInvertMatrix(boost::shared_ptr< MatrixDouble > in_ptr, boost::shared_ptr< VectorDouble > det_ptr, boost::shared_ptr< MatrixDouble > out_ptr)
boost::shared_ptr< MatrixDouble > outPtr
boost::shared_ptr< MatrixDouble > inPtr
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
Make Hdiv space from Hcurl space in 2d.
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
Operator for fat prism element updating integration weights in the volume.
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
Operator for linear form, usually to calculate values on right hand side.
boost::shared_ptr< MatrixDouble > outMat
DEPRECATED OpScaleMatrix(const std::string field_name, const double scale, boost::shared_ptr< MatrixDouble > in_mat, boost::shared_ptr< MatrixDouble > out_mat)
boost::shared_ptr< double > scalePtr
OpScaleMatrix(boost::shared_ptr< double > scale_ptr, boost::shared_ptr< MatrixDouble > in_mat, boost::shared_ptr< MatrixDouble > out_mat)
boost::shared_ptr< MatrixDouble > inMat
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
OpSetContravariantPiolaTransformOnEdge2D(const FieldSpace space=HCURL)
OpSetContravariantPiolaTransformOnFace2DImpl(boost::shared_ptr< MatrixDouble > jac_ptr)
Apply contravariant (Piola) transfer to Hdiv space on face.
Apply contravariant (Piola) transfer to Hdiv space on face.
Transform Hcurl base fluxes from reference element to physical triangle.
OpSetInvJacH1ForFaceEmbeddedIn3DSpace(boost::shared_ptr< MatrixDouble > inv_jac_ptr)
OpSetInvJacH1ForFace(boost::shared_ptr< MatrixDouble > inv_jac_ptr)
Transform local reference derivatives of shape functions to global derivatives.
const boost::shared_ptr< MatrixDouble > invJacPtr
OpSetInvJacH1ForFatPrism(MatrixDouble &inv_jac)
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
OpSetInvJacH1ForFatPrism(boost::shared_ptr< MatrixDouble > inv_jac_ptr)
Transform local reference derivatives of shape functions to global derivatives.
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
OpSetInvJacH1ForFlatPrism(MatrixDouble &inv_jac_f3)
boost::shared_ptr< MatrixDouble > invJacPtr
OpSetInvJacHcurlFaceImpl(boost::shared_ptr< MatrixDouble > inv_jac_ptr)
Transform local reference derivatives of shape function to global derivatives for face.
OpSetInvJacL2ForFaceEmbeddedIn3DSpace(boost::shared_ptr< MatrixDouble > inv_jac_ptr)
OpSetInvJacL2ForFace(boost::shared_ptr< MatrixDouble > inv_jac_ptr)
Transform local reference derivatives of shape functions to global derivatives.
MoFEMErrorCode applyTransform(MatrixDouble &diff_n)
OpSetInvJacToScalarBasesBasic(FieldSpace space, boost::shared_ptr< MatrixDouble > inv_jac_ptr)
boost::shared_ptr< MatrixDouble > invJacPtr
FTensor::Index< 'i', DIM > i
boost::shared_ptr< MatrixDouble > outMat
OpSymmetrizeTensor(boost::shared_ptr< MatrixDouble > in_mat, boost::shared_ptr< MatrixDouble > out_mat)
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
Operator for linear form, usually to calculate values on right hand side.
FTensor::Index< 'j', DIM > j
DEPRECATED OpSymmetrizeTensor(const std::string field_name, boost::shared_ptr< MatrixDouble > in_mat, boost::shared_ptr< MatrixDouble > out_mat)
boost::shared_ptr< MatrixDouble > inMat
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
Operator for linear form, usually to calculate values on right hand side.
boost::shared_ptr< MatrixDouble > dMat
OpTensorTimesSymmetricTensor(boost::shared_ptr< MatrixDouble > in_mat, boost::shared_ptr< MatrixDouble > out_mat, boost::shared_ptr< MatrixDouble > d_mat)
boost::shared_ptr< MatrixDouble > inMat
DEPRECATED OpTensorTimesSymmetricTensor(const std::string field_name, boost::shared_ptr< MatrixDouble > in_mat, boost::shared_ptr< MatrixDouble > out_mat, boost::shared_ptr< MatrixDouble > d_mat)
boost::shared_ptr< MatrixDouble > outMat
static constexpr Switches CtxSetX
static constexpr Switches CtxSetX_TT
std::bitset< 8 > Switches
static constexpr Switches CtxSetX_T
intrusive_ptr for managing petsc objects