PRISMS-PF Manual
Loading...
Searching...
No Matches
field_container.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: © 2025 PRISMS Center at the University of Michigan
2// SPDX-License-Identifier: GNU Lesser General Public Version 2.1
3
4#pragma once
5
6#include <deal.II/base/config.h>
7#include <deal.II/base/exceptions.h>
8#include <deal.II/base/point.h>
9#include <deal.II/base/symmetric_tensor.h>
10#include <deal.II/base/tensor.h>
11#include <deal.II/lac/vector.h>
12#include <deal.II/matrix_free/evaluation_flags.h>
13#include <deal.II/matrix_free/fe_evaluation.h>
14#include <deal.II/matrix_free/matrix_free.h>
15
22#include <prismspf/core/types.h>
23
24#include <prismspf/config.h>
25
26#include <utility>
27
29
45template <unsigned int dim, unsigned int degree, typename number>
47{
48public:
52 using ScalarValue = dealii::VectorizedArray<number>;
53
57 using VectorValue = dealii::Tensor<1, dim, ScalarValue>;
58
59 template <TensorRank Rank>
60 using Value = std::conditional_t<Rank == TensorRank::Scalar,
62 dealii::Tensor<int(Rank), dim, ScalarValue>>;
63
64 template <TensorRank Rank>
65 using Gradient = dealii::Tensor<int(Rank) + 1, dim, ScalarValue>;
66
67 template <TensorRank Rank>
68 using Hessian = dealii::Tensor<int(Rank) + 2, dim, ScalarValue>;
69
70 template <TensorRank Rank>
71 using FEEval =
72 dealii::FEEvaluation<dim,
73 degree,
74 degree + 1,
75 dealii::Tensor<int(Rank), dim>::n_independent_components,
76 number,
78
82 template <typename ValType>
83 static constexpr TensorRank RankFromVal = []() constexpr
84 {
85 if constexpr (std::is_same_v<ValType, ScalarValue> ||
86 std::is_same_v<ValType, typename ScalarValue::value_type>)
87 {
88 return TensorRank::Scalar;
89 }
90 else
91 {
92 return TensorRank(ValType::rank);
93 }
94 }();
95
99 template <typename GradType>
100 static constexpr TensorRank RankFromGrad = []() constexpr
101 {
102 if constexpr (std::is_same_v<GradType, ScalarValue>)
103 {
104 return TensorRank::Scalar;
105 }
106 else
107 {
108 return TensorRank(GradType::rank - 1);
109 }
110 }();
111
116 template <TensorRank Rank>
118 {
122 using FEEDepPair = std::pair<FEEval<Rank>, EvalFlags>;
123
128 using FEEDepPairPtr = std::shared_ptr<FEEDepPair>;
129
134
139
143 std::vector<FEEDepPairPtr> fe_eval_old;
144
149 EvalFlags integration_flags = EvalFlags::nothing;
150
158
166 unsigned int block_index = -1;
167
171 FEEValuationDeps() = default;
172
174 const Dependency &dependency,
175 const std::pair<const SolutionLevel<dim, number> *, unsigned int> &mf_id_pair,
176 bool is_dst);
177
178 template <DependencyType type>
179 const FEEval<Rank> &
180 get() const;
181
182 template <DependencyType type>
185
186 const FEEval<Rank> &
187 get(DependencyType type) const;
188
190 get(DependencyType type);
191
192 void
193 reinit(unsigned int cell);
194
195 void
196 eval(const BlockVector<number> *_src_solutions, bool plain);
197
198 void
199 reinit_and_eval(unsigned int cell,
200 const BlockVector<number> *_src_solutions,
201 bool plain);
202
203 void
204 integrate();
205
206 void
207 distribute(BlockVector<number> *dst_solutions);
208
209 void
211 };
212
225 FieldContainer(const std::vector<FieldAttributes> &_field_attributes,
226 const SolutionIndexer<dim, number> &_solution_indexer,
227 unsigned int _relative_level,
228 const DependencyMap &dependency_map,
229 const SolveBlock &_solve_block,
230 const MatrixFree<dim, number> &matrix_free);
231
235 void
236 reinit(unsigned int cell);
237
242 void
243 eval(const BlockVector<number> *src_solutions, bool plain);
244
251 void
252 reinit_and_eval(unsigned int cell,
253 const BlockVector<number> *src_solutions,
254 bool plain);
255
259 void
260 integrate();
261
265 void
266 distribute(BlockVector<number> *dst_solutions);
267
273 void
275
279 void
280 set_q_point(unsigned int q);
281
285 template <TensorRank Rank, DependencyType type>
286 [[nodiscard]] Value<Rank>
287 get_value(Types::Index field_index) const;
288
292 template <TensorRank Rank>
293 [[nodiscard]] Value<Rank>
294 get_value(Types::Index field_index, DependencyType type) const;
295
299 template <TensorRank Rank, DependencyType type>
300 [[nodiscard]] Gradient<Rank>
301 get_gradient(Types::Index field_index) const;
302
306 template <TensorRank Rank>
307 [[nodiscard]] Gradient<Rank>
308 get_gradient(Types::Index field_index, DependencyType type) const;
309
313 template <TensorRank Rank, DependencyType type>
314 [[nodiscard]] Hessian<Rank>
315 get_hessian(Types::Index field_index) const;
316
320 template <TensorRank Rank>
321 [[nodiscard]] Hessian<Rank>
322 get_hessian(Types::Index field_index, DependencyType type) const;
323
327 template <TensorRank Rank, DependencyType type>
328 [[nodiscard]] Gradient<Rank>
330
334 template <TensorRank Rank>
335 [[nodiscard]] Gradient<Rank>
337
341 template <TensorRank Rank, DependencyType type>
342 [[nodiscard]] Value<Rank>
343 get_laplacian(Types::Index field_index) const;
344
348 template <TensorRank Rank>
349 [[nodiscard]] Value<Rank>
350 get_laplacian(Types::Index field_index, DependencyType type) const;
351
355 template <TensorRank Rank, DependencyType type>
356 [[nodiscard]] ScalarValue
357 get_divergence(Types::Index field_index) const;
358
362 template <TensorRank Rank>
363 [[nodiscard]] ScalarValue
365
369 template <TensorRank Rank, DependencyType type>
370 [[nodiscard]] dealii::SymmetricTensor<2, dim, ScalarValue>
372
376 template <TensorRank Rank>
377 [[nodiscard]] dealii::SymmetricTensor<2, dim, ScalarValue>
379
383 template <TensorRank Rank, DependencyType type>
384 [[nodiscard]] dealii::Tensor<1, (dim == 2 ? 1 : dim), ScalarValue>
385 get_curl(Types::Index field_index) const;
386
390 template <TensorRank Rank>
391 [[nodiscard]] dealii::Tensor<1, (dim == 2 ? 1 : dim), ScalarValue>
392 get_curl(Types::Index field_index, DependencyType type) const;
393
397 [[nodiscard]] dealii::Point<dim, ScalarValue>
398 get_q_point_location() const;
399
403 [[nodiscard]] ScalarValue
404 get_element_volume() const;
405
409 [[nodiscard]] unsigned int
410 get_n_q_points() const;
411
415 template <typename ValType>
416 void
417 set_value_term(Types::Index field_index, const ValType &val);
418
422 template <typename GradType>
423 void
424 set_gradient_term(Types::Index field_index, const GradType &val);
425
426private:
427 template <TensorRank Rank>
428 std::vector<FEEValuationDeps<Rank>> &
430
431 template <TensorRank Rank>
432 const std::vector<FEEValuationDeps<Rank>> &
434
439 void
440 feevaluation_size_valid(Types::Index field_index) const;
441
446 void
447 feevaluation_exists(Types::Index field_index, Types::Index dependency_index) const;
448
453 void
454 access_valid(Types::Index field_index,
455 DependencyType dependency_type,
456 dealii::EvaluationFlags::EvaluationFlags flag) const;
457
461 void
462 submission_valid(Types::Index field_index, DependencyType dependency_type) const;
463
467 const std::vector<FieldAttributes> *field_attributes_ptr;
468
473
477 std::vector<FEEValuationDeps<TensorRank::Scalar>> feeval_deps_scalar;
478
482 std::vector<FEEValuationDeps<TensorRank::Vector>> feeval_deps_vector;
483
488
499
503 unsigned int relative_level;
504
508 unsigned int q_point = 0;
509};
510
511inline static const std::map<DependencyType, std::string> dependency_type_to_string = {
512 {DependencyType::Current, "current"},
513 {DependencyType::OldOne, "old_1" },
514 {DependencyType::OldTwo, "old_2" },
515 {DependencyType::OldThree, "old_3" },
516 {DependencyType::OldFour, "old_4" },
517 {DependencyType::SRC, "src" },
518 {DependencyType::DST, "dst" }
519};
520
521template <unsigned int dim, unsigned int degree, typename number>
522template <TensorRank Rank>
524 const Dependency &dependency,
525 const std::pair<const SolutionLevel<dim, number> *, unsigned int> &mf_id_pair,
526 bool is_dst)
527 : solution_level(mf_id_pair.first)
528 , block_index(mf_id_pair.second)
529{
530 // Make an FEEvaluation if the current solution needs to be evaluated
531 if (dependency.flag)
532 {
533 fe_eval = std::make_shared<FEEDepPair>(FEEval<Rank>(solution_level->matrix_free,
534 block_index),
535 dependency.flag);
536 }
537 // Make FEEvaluations for the the old solutions
538 fe_eval_old.resize(dependency.old_flags.size(), nullptr);
539 for (unsigned int age = 0; age < dependency.old_flags.size(); ++age)
540 {
541 if (dependency.old_flags.at(age)) // kinda redundant... maybe remove?
542 {
543 fe_eval_old[age] =
544 std::make_shared<FEEDepPair>(FEEval<Rank>(solution_level->matrix_free,
546 dependency.old_flags.at(age));
547 }
548 }
549 if (dependency.src_flag || is_dst)
550 {
551 fe_eval_src_dst =
552 std::make_shared<FEEDepPair>(FEEval<Rank>(solution_level->matrix_free,
553 block_index),
554 dependency.src_flag);
555 }
556}
557
558class ExcDepNotInitialized : public dealii::ExceptionBase
559{
560public:
562 : dependency_type(_dependency_type)
563 {}
564
565 const char *
566 what() const noexcept override
567 {
568 return "Access was attempted for a field that was not declared as a "
569 "dependency for the current solve block.\n";
570 }
571
572 void
573 print_info(std::ostream &out) const override
574 {
575 out << what() << std::flush;
576 }
577
578private:
580};
581#ifdef DEBUG
582# define AssertAccessible(fe_eval_pair_ptr, dependency_type) \
583 AssertThrowDebug(block_index != -1, ExcDepNotInitialized(dependency_type)); \
584 AssertThrowDebug(fe_eval_pair_ptr, ExcDepNotInitialized(dependency_type));
585#else
586# define AssertAccessible(fe_eval_pair_ptr, dependency_type)
587#endif
588
589template <unsigned int dim, unsigned int degree, typename number>
590template <TensorRank Rank>
591template <DependencyType Type>
592inline DEAL_II_ALWAYS_INLINE const typename FieldContainer<dim, degree, number>::
593 template FEEval<Rank> &
595{
596 if constexpr (Type == DependencyType::SRC || Type == DependencyType::DST)
597 {
599 return fe_eval_src_dst->first;
600 }
601 else if constexpr (Type == DependencyType::Current)
602 {
604 return fe_eval->first;
605 }
606 else
607 {
610 return fe_eval_old[int(Type) - 1]->first;
611 }
612}
613
614template <unsigned int dim, unsigned int degree, typename number>
615template <TensorRank Rank>
616template <DependencyType Type>
617inline DEAL_II_ALWAYS_INLINE
620{
621 if constexpr (Type == DependencyType::SRC || Type == DependencyType::DST)
622 {
624 return fe_eval_src_dst->first;
625 }
626 else if constexpr (Type == DependencyType::Current)
627 {
629 return fe_eval->first;
630 }
631 else
632 {
635 return fe_eval_old[int(Type) - 1]->first;
636 }
637}
638
639template <unsigned int dim, unsigned int degree, typename number>
640template <TensorRank Rank>
641inline DEAL_II_ALWAYS_INLINE const typename FieldContainer<dim, degree, number>::
642 template FEEval<Rank> &
644 DependencyType type) const
645{
646 if (type == DependencyType::SRC || type == DependencyType::DST)
647 {
649 return fe_eval_src_dst->first;
650 }
651 if (type == DependencyType::Current)
652 {
654 return fe_eval->first;
655 }
656 {
657 AssertThrowDebug(int(type) - 1 < fe_eval_old.size(), ExcDepNotInitialized(type));
658 AssertAccessible(fe_eval_old[int(type) - 1], type);
659 return fe_eval_old[type - 1]->first;
660 }
661}
662
663template <unsigned int dim, unsigned int degree, typename number>
664template <TensorRank Rank>
665inline DEAL_II_ALWAYS_INLINE
668{
669 if (type == DependencyType::SRC || type == DependencyType::DST)
670 {
672 return fe_eval_src_dst->first;
673 }
674 if (type == DependencyType::Current)
675 {
677 return fe_eval->first;
678 }
679 {
680 AssertThrowDebug(int(type) - 1 < fe_eval_old.size(), ExcDepNotInitialized(type));
681 AssertAccessible(fe_eval_old[int(type) - 1], type);
682 return fe_eval_old[int(type) - 1]->first;
683 }
684}
685
686template <unsigned int dim, unsigned int degree, typename number>
687template <TensorRank Rank>
688inline void
690{
691 if (fe_eval)
692 {
693 fe_eval->first.reinit(cell);
694 }
695 for (auto &old_fe_eval : fe_eval_old)
696 {
697 if (old_fe_eval)
698 {
699 old_fe_eval->first.reinit(cell);
700 }
701 }
702 if (fe_eval_src_dst)
703 {
704 fe_eval_src_dst->first.reinit(cell);
705 }
706}
707
708template <unsigned int dim, unsigned int degree, typename number>
709template <TensorRank Rank>
710inline void
712 const BlockVector<number> *_src_solutions,
713 bool plain)
714{
715 // NOTE: `read_dof_values_plain` must be called here so that constraints aren't
716 // implicitly applied. This allows us to have inhomogeneous constraints.
717 if (fe_eval)
718 {
719 fe_eval->first.read_dof_values_plain(solution_level->solutions.block(block_index));
720 fe_eval->first.evaluate(fe_eval->second);
721 }
722 for (unsigned int age = 0; age < fe_eval_old.size(); ++age)
723 {
724 if (FEEDepPairPtr &old_fe_eval = fe_eval_old[age])
725 {
726 old_fe_eval->first.read_dof_values_plain(
727 solution_level->old_solutions[age].block(block_index));
728 old_fe_eval->first.evaluate(old_fe_eval->second);
729 }
730 }
731 if (fe_eval_src_dst && fe_eval_src_dst->second != EvalFlags::nothing)
732 {
733 if (plain)
734 {
735 fe_eval_src_dst->first.read_dof_values_plain(
736 _src_solutions->block(block_index));
737 }
738 else
739 {
740 fe_eval_src_dst->first.read_dof_values(_src_solutions->block(block_index));
741 }
742 fe_eval_src_dst->first.evaluate(fe_eval_src_dst->second);
743 }
744}
745
746template <unsigned int dim, unsigned int degree, typename number>
747template <TensorRank Rank>
748inline void
750 unsigned int cell,
751 const BlockVector<number> *_src_solutions,
752 bool plain)
753{
754 // NOTE: `read_dof_values_plain` must be called here so that constraints aren't
755 // implicitly applied. This allows us to have inhomogeneous constraints.
756 if (fe_eval)
757 {
758 fe_eval->first.reinit(cell);
759 fe_eval->first.read_dof_values_plain(solution_level->solutions.block(block_index));
760 fe_eval->first.evaluate(fe_eval->second);
761 }
762 for (unsigned int age = 0; age < fe_eval_old.size(); ++age)
763 {
764 if (FEEDepPairPtr &old_fe_eval = fe_eval_old[age])
765 {
766 old_fe_eval->first.reinit(cell);
767 old_fe_eval->first.read_dof_values_plain(
768 solution_level->old_solutions[age].block(block_index));
769 old_fe_eval->first.evaluate(old_fe_eval->second);
770 }
771 }
772 if (fe_eval_src_dst)
773 {
774 fe_eval_src_dst->first.reinit(cell);
775 if (fe_eval_src_dst->second != EvalFlags::nothing)
776 {
777 if (plain)
778 {
779 fe_eval_src_dst->first.read_dof_values_plain(
780 _src_solutions->block(block_index));
781 }
782 else
783 {
784 fe_eval_src_dst->first.read_dof_values(_src_solutions->block(block_index));
785 }
786 fe_eval_src_dst->first.evaluate(fe_eval_src_dst->second);
787 }
788 }
789}
790
791template <unsigned int dim, unsigned int degree, typename number>
792template <TensorRank Rank>
793inline void
801
802template <unsigned int dim, unsigned int degree, typename number>
803template <TensorRank Rank>
804inline void
806 BlockVector<number> *dst_solutions)
807{
808 if (fe_eval_src_dst)
809 {
810 fe_eval_src_dst->first.distribute_local_to_global(
811 dst_solutions->block(block_index));
812 }
813}
814
815template <unsigned int dim, unsigned int degree, typename number>
816template <TensorRank Rank>
817inline void
819 BlockVector<number> *dst_solutions)
820{
821 if (fe_eval_src_dst)
822 {
823 fe_eval_src_dst->first.integrate_scatter(integration_flags,
824 dst_solutions->block(block_index));
825 }
826}
827
828template <unsigned int dim, unsigned int degree, typename number>
829inline void
831{
832 for (auto &fe_eval : feeval_deps_scalar)
833 {
834 fe_eval.reinit(cell);
835 }
836 for (auto &fe_eval : feeval_deps_vector)
837 {
838 fe_eval.reinit(cell);
839 }
840 shared_feeval_scalar.reinit(cell);
841}
842
843template <unsigned int dim, unsigned int degree, typename number>
844inline void
846 bool plain)
847{
848 for (auto &fe_eval : feeval_deps_scalar)
849 {
850 fe_eval.eval(src_solutions, plain);
851 }
852 for (auto &fe_eval : feeval_deps_vector)
853 {
854 fe_eval.eval(src_solutions, plain);
855 }
856 // Don't eval `shared_feeval_scalar` because we only use it for information.
857}
858
859template <unsigned int dim, unsigned int degree, typename number>
860inline void
862 unsigned int cell,
863 const BlockVector<number> *src_solutions,
864 bool plain)
865{
866 for (auto &fe_eval : feeval_deps_scalar)
867 {
868 fe_eval.reinit_and_eval(cell, src_solutions, plain);
869 }
870 for (auto &fe_eval : feeval_deps_vector)
871 {
872 fe_eval.reinit_and_eval(cell, src_solutions, plain);
873 }
874 // Don't eval `shared_feeval_scalar` because we only use it for information.
875 shared_feeval_scalar.reinit(cell);
876}
877
878template <unsigned int dim, unsigned int degree, typename number>
879inline void
881{
882 const std::vector<FieldAttributes> &field_attributes = *field_attributes_ptr;
883 for (const Types::Index &field_index : solve_block->field_indices)
884 {
885 if (field_attributes[field_index].field_type == TensorRank::Scalar)
886 {
887 feeval_deps_scalar[field_index].integrate();
888 }
889 else /* vector */
890 {
891 feeval_deps_vector[field_index].integrate();
892 }
893 }
894 // Don't integrate `shared_feeval_scalar` because we only use it for information.
895}
896
897template <unsigned int dim, unsigned int degree, typename number>
898inline void
900{
901 const std::vector<FieldAttributes> &field_attributes = *field_attributes_ptr;
902 for (const Types::Index &field_index : solve_block->field_indices)
903 {
904 if (field_attributes[field_index].field_type == TensorRank::Scalar)
905 {
906 feeval_deps_scalar[field_index].distribute(dst_solutions);
907 }
908 else /* vector */
909 {
910 feeval_deps_vector[field_index].distribute(dst_solutions);
911 }
912 }
913 // Don't distribute `shared_feeval_scalar` because we only use it for information.
914}
915
916template <unsigned int dim, unsigned int degree, typename number>
917inline void
919 BlockVector<number> *dst_solutions)
920{
921 const std::vector<FieldAttributes> &field_attributes = *field_attributes_ptr;
922 for (const Types::Index &field_index : solve_block->field_indices)
923 {
924 if (field_attributes[field_index].field_type == TensorRank::Scalar)
925 {
926 feeval_deps_scalar[field_index].integrate_and_distribute(dst_solutions);
927 }
928 else /* vector */
929 {
930 feeval_deps_vector[field_index].integrate_and_distribute(dst_solutions);
931 }
932 }
933 // Don't integrate and distribute `shared_feeval_scalar` because we only use it for
934 // information.
935}
936
937template <unsigned int dim, unsigned int degree, typename number>
938inline DEAL_II_ALWAYS_INLINE void
943
944// there are two catches we can do here.
945// 1. Dependencies for the dependency type (current, old, src/dst) don't exist.
946// 2. Dependency is not initialized for values/gradients.
947// We catch these separately to give more informative error messages.
948#define GetterTempl(dependency_type) template get<dependency_type>()
949#define GetterNoTempl(dependency_type) get(dependency_type)
950#define ReturnGetter(get_handle, Rank, field_index, dependency_type, getter) \
951 try \
952 { \
953 AssertThrowDebug((field_index) < get_relevant_feeval_vector<Rank>().size(), \
954 dealii::ExcMessage("Error: Field index " + \
955 std::to_string(field_index) + \
956 " is not associated with any field.")); \
957 return get_relevant_feeval_vector<Rank>()[field_index] \
958 .getter(dependency_type) \
959 .get_handle(q_point); \
960 } \
961 catch (const ExcDepNotInitialized &e) \
962 { \
963 std::cerr << "Error when trying to access field with index " << (field_index) \
964 << " and dependency type " \
965 << dependency_type_to_string.at(dependency_type) << ":\n" \
966 << e.what() \
967 << "Ensure that each dependency is requested in the solve block.\n" \
968 << std::flush; \
969 throw; \
970 } \
971 catch (const dealii::internal::ExcAccessToUninitializedField &e) \
972 { \
973 std::cerr << "Error when trying to access field with index " << (field_index) \
974 << " and dependency type " \
975 << dependency_type_to_string.at(dependency_type) \
976 << ":\nAccess was attempted for values or gradients that were " \
977 "not requested.\n" \
978 << std::flush; \
979 throw; \
980 }
981
982template <unsigned int dim, unsigned int degree, typename number>
983template <TensorRank Rank, DependencyType type>
984inline DEAL_II_ALWAYS_INLINE
987{
988 ReturnGetter(get_value, Rank, field_index, type, GetterTempl);
989}
990
991template <unsigned int dim, unsigned int degree, typename number>
992template <TensorRank Rank>
993inline DEAL_II_ALWAYS_INLINE
996 DependencyType type) const
997{
998 ReturnGetter(get_value, Rank, field_index, type, GetterNoTempl);
999}
1000
1001template <unsigned int dim, unsigned int degree, typename number>
1002template <TensorRank Rank, DependencyType type>
1003inline DEAL_II_ALWAYS_INLINE
1006{
1007 ReturnGetter(get_gradient, Rank, field_index, type, GetterTempl);
1008}
1009
1010template <unsigned int dim, unsigned int degree, typename number>
1011template <TensorRank Rank>
1012inline DEAL_II_ALWAYS_INLINE
1015 DependencyType type) const
1016{
1017 ReturnGetter(get_gradient, Rank, field_index, type, GetterNoTempl);
1018}
1019
1020template <unsigned int dim, unsigned int degree, typename number>
1021template <TensorRank Rank, DependencyType type>
1022inline DEAL_II_ALWAYS_INLINE
1025{
1026 ReturnGetter(get_hessian, Rank, field_index, type, GetterTempl);
1027}
1028
1029template <unsigned int dim, unsigned int degree, typename number>
1030template <TensorRank Rank>
1031inline DEAL_II_ALWAYS_INLINE
1034 DependencyType type) const
1035{
1036 ReturnGetter(get_hessian, Rank, field_index, type, GetterNoTempl);
1037}
1038
1039template <unsigned int dim, unsigned int degree, typename number>
1040template <TensorRank Rank, DependencyType type>
1041inline DEAL_II_ALWAYS_INLINE
1048
1049template <unsigned int dim, unsigned int degree, typename number>
1050template <TensorRank Rank>
1051inline DEAL_II_ALWAYS_INLINE
1058
1059template <unsigned int dim, unsigned int degree, typename number>
1060template <TensorRank Rank, DependencyType type>
1061inline DEAL_II_ALWAYS_INLINE
1064{
1065 ReturnGetter(get_laplacian, Rank, field_index, type, GetterTempl);
1066}
1067
1068template <unsigned int dim, unsigned int degree, typename number>
1069template <TensorRank Rank>
1070inline DEAL_II_ALWAYS_INLINE
1073 DependencyType type) const
1074{
1075 ReturnGetter(get_laplacian, Rank, field_index, type, GetterNoTempl);
1076}
1077
1078template <unsigned int dim, unsigned int degree, typename number>
1079template <TensorRank Rank, DependencyType type>
1080inline DEAL_II_ALWAYS_INLINE typename FieldContainer<dim, degree, number>::ScalarValue
1082{
1083 static_assert(Rank == 1, "Divergences are only available for vector fields");
1084 ReturnGetter(get_divergence, Rank, field_index, type, GetterTempl);
1085}
1086
1087template <unsigned int dim, unsigned int degree, typename number>
1088template <TensorRank Rank>
1089inline DEAL_II_ALWAYS_INLINE typename FieldContainer<dim, degree, number>::ScalarValue
1091 DependencyType type) const
1092{
1093 static_assert(Rank == 1, "Divergences are only available for vector fields");
1094 ReturnGetter(get_divergence, Rank, field_index, type, GetterNoTempl);
1095}
1096
1097template <unsigned int dim, unsigned int degree, typename number>
1098template <TensorRank Rank, DependencyType type>
1099inline DEAL_II_ALWAYS_INLINE dealii::
1100 SymmetricTensor<2, dim, typename FieldContainer<dim, degree, number>::ScalarValue>
1102 Types::Index field_index) const
1103{
1104 static_assert(Rank == 1, "Symmetric gradients are only available for vector fields");
1105 ReturnGetter(get_symmetric_gradient, Rank, field_index, type, GetterTempl);
1106}
1107
1108template <unsigned int dim, unsigned int degree, typename number>
1109template <TensorRank Rank>
1110inline DEAL_II_ALWAYS_INLINE dealii::
1111 SymmetricTensor<2, dim, typename FieldContainer<dim, degree, number>::ScalarValue>
1113 DependencyType type) const
1114{
1115 static_assert(Rank == 1, "Symmetric gradients are only available for vector fields");
1116 ReturnGetter(get_symmetric_gradient, Rank, field_index, type, GetterNoTempl);
1117}
1118
1119template <unsigned int dim, unsigned int degree, typename number>
1120template <TensorRank Rank, DependencyType type>
1121inline DEAL_II_ALWAYS_INLINE
1122 dealii::Tensor<1,
1123 (dim == 2 ? 1 : dim),
1126{
1127 static_assert(Rank == 1, "Curl is only available for vector fields");
1129 dim > 1,
1130 dealii::ExcMessage(
1131 "Curl is only available for vector fields with dimension greater than 1."));
1132 ReturnGetter(get_curl, Rank, field_index, type, GetterTempl);
1133}
1134
1135template <unsigned int dim, unsigned int degree, typename number>
1136template <TensorRank Rank>
1137inline DEAL_II_ALWAYS_INLINE
1138 dealii::Tensor<1,
1139 (dim == 2 ? 1 : dim),
1142 DependencyType type) const
1143{
1144 static_assert(Rank == 1, "Curl is only available for vector fields");
1146 dim > 1,
1147 dealii::ExcMessage(
1148 "Curl is only available for vector fields with dimension greater than 1."));
1149 ReturnGetter(get_curl, Rank, field_index, type, GetterNoTempl);
1150}
1151
1152template <unsigned int dim, unsigned int degree, typename number>
1153inline DEAL_II_ALWAYS_INLINE
1154 dealii::Point<dim, typename FieldContainer<dim, degree, number>::ScalarValue>
1159
1160template <unsigned int dim, unsigned int degree, typename number>
1161inline DEAL_II_ALWAYS_INLINE typename FieldContainer<dim, degree, number>::ScalarValue
1167
1168template <unsigned int dim, unsigned int degree, typename number>
1169inline DEAL_II_ALWAYS_INLINE unsigned int
1174
1175template <unsigned int dim, unsigned int degree, typename number>
1176template <typename ValType>
1177inline DEAL_II_ALWAYS_INLINE void
1179 const ValType &val)
1180{
1181 AssertThrowDebug((field_index) <
1183 dealii::ExcMessage("Error: Field index " +
1184 std::to_string(field_index) +
1185 " is not associated with any field."));
1186 auto &relevant_feeval_vector =
1188 try
1189 {
1190 relevant_feeval_vector.template get<DependencyType::DST>().submit_value(val,
1191 q_point);
1192 relevant_feeval_vector.integration_flags |= dealii::EvaluationFlags::values;
1193 }
1194 catch (...)
1195 {
1196 std::cerr << "Error when trying to submit value for field with index "
1197 << (field_index)
1198 << ": Error: Submission for field not part of this solve block.\n"
1199 << std::flush;
1200 throw;
1201 }
1202}
1203
1204template <unsigned int dim, unsigned int degree, typename number>
1205template <typename GradType>
1206inline DEAL_II_ALWAYS_INLINE void
1208 const GradType &val)
1209{
1210 AssertThrowDebug((field_index) <
1212 dealii::ExcMessage("Error: Field index " +
1213 std::to_string(field_index) +
1214 " is not associated with any field."));
1215 auto &relevant_feeval_vector =
1217 try
1218 {
1219 relevant_feeval_vector.template get<DependencyType::DST>().submit_gradient(val,
1220 q_point);
1221 relevant_feeval_vector.integration_flags |= dealii::EvaluationFlags::gradients;
1222 }
1223 catch (...)
1224 {
1225 std::cerr << "Error when trying to submit gradient for field with index "
1226 << (field_index)
1227 << ": Error: Submission for field not part of this solve block.\n"
1228 << std::flush;
1229 throw;
1230 }
1231}
1232
1233template <unsigned int dim, unsigned int degree, typename number>
1234template <TensorRank Rank>
1235inline DEAL_II_ALWAYS_INLINE std::vector<
1236 typename FieldContainer<dim, degree, number>::template FEEValuationDeps<Rank>> &
1238{
1239 if constexpr (Rank == TensorRank::Scalar)
1240 {
1241 return feeval_deps_scalar;
1242 }
1243 else if constexpr (Rank == TensorRank::Vector)
1244 {
1245 return feeval_deps_vector;
1246 }
1247}
1248
1249template <unsigned int dim, unsigned int degree, typename number>
1250template <TensorRank Rank>
1251inline DEAL_II_ALWAYS_INLINE const std::vector<
1252 typename FieldContainer<dim, degree, number>::template FEEValuationDeps<Rank>> &
1254{
1255 if constexpr (Rank == TensorRank::Scalar)
1256 {
1257 return feeval_deps_scalar;
1258 }
1259 else if constexpr (Rank == TensorRank::Vector)
1260 {
1261 return feeval_deps_vector;
1262 }
1263}
1264
1265template <unsigned int dim, unsigned int degree, typename number>
1266inline DEAL_II_ALWAYS_INLINE void
1268 [[maybe_unused]] Types::Index field_index) const
1269{
1270 // TODO
1271}
1272
1273template <unsigned int dim, unsigned int degree, typename number>
1274inline DEAL_II_ALWAYS_INLINE void
1276 [[maybe_unused]] Types::Index field_index,
1277 [[maybe_unused]] Types::Index dependency_index) const
1278{
1279 // TODO
1280}
1281
1282template <unsigned int dim, unsigned int degree, typename number>
1283inline DEAL_II_ALWAYS_INLINE void
1285 [[maybe_unused]] Types::Index field_index,
1286 [[maybe_unused]] DependencyType dependency_type,
1287 [[maybe_unused]] dealii::EvaluationFlags::EvaluationFlags flag) const
1288{
1289 // TODO
1290}
1291
1292template <unsigned int dim, unsigned int degree, typename number>
1293inline DEAL_II_ALWAYS_INLINE void
1295 [[maybe_unused]] Types::Index field_index,
1296 [[maybe_unused]] DependencyType dependency_type) const
1297{
1298 // TODO
1299}
1300
1301#undef AssertAccessible
1302#undef ReturnGetter
1303#undef GetterTempl
1304#undef GetterNoTempl
1305
1306PRISMS_PF_END_NAMESPACE
Definition field_container.h:559
ExcDepNotInitialized(DependencyType _dependency_type)
Definition field_container.h:561
const char * what() const noexcept override
Definition field_container.h:566
void print_info(std::ostream &out) const override
Definition field_container.h:573
DependencyType dependency_type
Definition field_container.h:579
This class permits the access of a subset of indexed fields and gives an error if any non-allowed fie...
Definition field_container.h:47
std::vector< FEEValuationDeps< Rank > > & get_relevant_feeval_vector()
void set_value_term(Types::Index field_index, const ValType &val)
Set the residual value of the specified scalar/vector field.
Value< Rank > get_value(Types::Index field_index, DependencyType type) const
Return the value of the specified field.
void set_gradient_term(Types::Index field_index, const GradType &val)
Set the residual gradient of the specified scalar/vector field.
Gradient< Rank > get_hessian_diagonal(Types::Index field_index, DependencyType type) const
Return the diagonal of the hessian of the specified field.
dealii::Tensor< int(Rank)+2, dim, ScalarValue > Hessian
Definition field_container.h:68
void eval(const BlockVector< number > *src_solutions, bool plain)
Read solution vector, and evaluate based on dependency flags for all dependencies.
Definition field_container.h:845
void access_valid(Types::Index field_index, DependencyType dependency_type, dealii::EvaluationFlags::EvaluationFlags flag) const
Check that a variable value/gradient/hessians was marked as needed and thus properly initialized.
Definition field_container.h:1284
Hessian< Rank > get_hessian(Types::Index field_index) const
Return the hessian of the specified field.
Value< Rank > get_laplacian(Types::Index field_index, DependencyType type) const
Return the laplacian of the specified field.
FEEval< TensorRank::Scalar > shared_feeval_scalar
FEEvaluation object for generic cell operations.
Definition field_container.h:498
dealii::Tensor< 1, dim, ScalarValue > VectorValue
Typedef for the basic vector that the user manipulates.
Definition field_container.h:57
unsigned int relative_level
Multigrid level.
Definition field_container.h:503
FieldContainer(const std::vector< FieldAttributes > &_field_attributes, const SolutionIndexer< dim, number > &_solution_indexer, unsigned int _relative_level, const DependencyMap &dependency_map, const SolveBlock &_solve_block, const MatrixFree< dim, number > &matrix_free)
Constructor.
Definition field_container.cc:14
dealii::Tensor< 1,(dim==2 ? 1 :dim), ScalarValue > get_curl(Types::Index field_index) const
Return the curl of the specified field.
ScalarValue get_divergence(Types::Index field_index) const
Return the divergence of the specified field.
Hessian< Rank > get_hessian(Types::Index field_index, DependencyType type) const
Return the hessian of the specified field.
void set_q_point(unsigned int q)
Set the current quadrature point.
Definition field_container.h:939
void submission_valid(Types::Index field_index, DependencyType dependency_type) const
Check that a value is valid for submission.
Definition field_container.h:1294
dealii::FEEvaluation< dim, degree, degree+1, dealii::Tensor< int(Rank), dim >::n_independent_components, number, ScalarValue > FEEval
Definition field_container.h:71
unsigned int q_point
The quadrature point.
Definition field_container.h:508
dealii::Tensor< int(Rank)+1, dim, ScalarValue > Gradient
Definition field_container.h:65
dealii::VectorizedArray< number > ScalarValue
Typedef for the basic value that the user manipulates.
Definition field_container.h:52
void integrate()
Integrate the residuals.
Definition field_container.h:880
dealii::SymmetricTensor< 2, dim, ScalarValue > get_symmetric_gradient(Types::Index field_index) const
Return the symmetric gradient of the specified field.
std::conditional_t< Rank==TensorRank::Scalar, ScalarValue, dealii::Tensor< int(Rank), dim, ScalarValue > > Value
Definition field_container.h:60
void reinit_and_eval(unsigned int cell, const BlockVector< number > *src_solutions, bool plain)
Initialize based on cell, read solution vector, and evaluate based on dependency flags for all depend...
Definition field_container.h:861
Value< Rank > get_laplacian(Types::Index field_index) const
Return the laplacian of the specified field.
ScalarValue get_divergence(Types::Index field_index, DependencyType type) const
Return the divergence of the specified field.
std::vector< FEEValuationDeps< TensorRank::Vector > > feeval_deps_vector
Collection of FEEvaluation dependencies for vector fields.
Definition field_container.h:482
const std::vector< FEEValuationDeps< Rank > > & get_relevant_feeval_vector() const
static constexpr TensorRank RankFromVal
Return the tensor rank from the specified template value.
Definition field_container.h:83
void distribute(BlockVector< number > *dst_solutions)
Distribute the integrated residuals.
Definition field_container.h:899
void reinit(unsigned int cell)
Initialize based on cell for all dependencies.
Definition field_container.h:830
Gradient< Rank > get_gradient(Types::Index field_index, DependencyType type) const
Return the gradient of the specified field.
std::vector< FEEValuationDeps< TensorRank::Scalar > > feeval_deps_scalar
Collection of FEEvaluation dependencies for scalar fields.
Definition field_container.h:477
void feevaluation_exists(Types::Index field_index, Types::Index dependency_index) const
Check whether the entry for the FEEvaluation is within the bounds of the vector and not a nullptr.
Definition field_container.h:1275
void integrate_and_distribute(BlockVector< number > *dst_solutions)
Integrate the residuals and distribute from local to global.
Definition field_container.h:918
const std::vector< FieldAttributes > * field_attributes_ptr
Pointer to the vector of field attributes.
Definition field_container.h:467
dealii::SymmetricTensor< 2, dim, ScalarValue > get_symmetric_gradient(Types::Index field_index, DependencyType type) const
Return the symmetric gradient of the specified field.
Gradient< Rank > get_gradient(Types::Index field_index) const
Return the gradient of the specified field.
unsigned int get_n_q_points() const
Return the number of quadrature points.
Definition field_container.h:1170
Gradient< Rank > get_hessian_diagonal(Types::Index field_index) const
Return the diagonal of the hessian of the specified field.
dealii::Point< dim, ScalarValue > get_q_point_location() const
Return the quadrature point location.
Definition field_container.h:1155
ScalarValue get_element_volume() const
Return the quadrature point location.
Definition field_container.h:1162
static constexpr TensorRank RankFromGrad
Return the tensor rank from the specified template gradient.
Definition field_container.h:100
const SolveBlock * solve_block
Solve block information.
Definition field_container.h:487
const SolutionIndexer< dim, number > * solution_indexer
Pointer to the solution indexer.
Definition field_container.h:472
dealii::Tensor< 1,(dim==2 ? 1 :dim), ScalarValue > get_curl(Types::Index field_index, DependencyType type) const
Return the curl of the specified field.
void feevaluation_size_valid(Types::Index field_index) const
Check whether the entry for the FEEvaluation is within the bounds of the vector.
Definition field_container.h:1267
Value< Rank > get_value(Types::Index field_index) const
Return the value of the specified field.
Class that provides access to solution vectors spread across different groups.
Definition solution_indexer.h:20
Structure to hold the attributes of a solve-block.
Definition solve_block.h:56
static const dealii::QGaussLobatto< dim > quadrature
Quadrature rule.
Definition system_wide.h:41
std::map< Types::Index, Dependency > DependencyMap
Definition dependencies.h:129
#define AssertThrowDebug(cond, exc)
Definition exceptions.h:24
#define AssertAccessible(fe_eval_pair_ptr, dependency_type)
Definition field_container.h:586
#define GetterNoTempl(dependency_type)
Definition field_container.h:949
#define GetterTempl(dependency_type)
Definition field_container.h:948
static const std::map< DependencyType, std::string > dependency_type_to_string
Definition field_container.h:511
#define ReturnGetter(get_handle, Rank, field_index, dependency_type, getter)
Definition field_container.h:950
@ Value
Use value of the variable as a criterion for refinement.
Definition grid_refiner_criterion.h:31
@ Gradient
Use gradient of the variable as a criterion for refinement.
Definition grid_refiner_criterion.h:36
dealii::LinearAlgebra::distributed::BlockVector< number > BlockVector
Typedef for solution block vector.
Definition group_solution_handler.h:29
Definition conditional_ostreams.cc:20
unsigned int Index
Type for field indices.
Definition types.h:19
Type
Factory function to create appropriate reader based on input file type not a member of ReadFieldBase ...
Definition read_field_factory.h:30
Dependency struct containing evaluation flags for each field.
Definition dependencies.h:30
EvalFlags src_flag
Evaluation flags for the DependencyType::SRC solution.
Definition dependencies.h:42
std::vector< EvalFlags > old_flags
Collection of evaluation flags for the old solutions.
Definition dependencies.h:47
EvalFlags flag
Evaluation flags for the current solution.
Definition dependencies.h:36
FEEDepPairPtr fe_eval_src_dst
FEEvaluation for the src and dst solutions.
Definition field_container.h:138
void integrate()
Definition field_container.h:794
std::shared_ptr< FEEDepPair > FEEDepPairPtr
Pointer the the FEEDepPair. We have to use ptrs because we don't want to construct fields we don't ne...
Definition field_container.h:128
const FEEval< Rank > & get() const
EvalFlags integration_flags
Destination integration flags. We use this to integrate value and gradient terms, if they exist.
Definition field_container.h:149
void distribute(BlockVector< number > *dst_solutions)
Definition field_container.h:805
void reinit(unsigned int cell)
Definition field_container.h:689
FEEDepPairPtr fe_eval
FEEvaluation for the current solution.
Definition field_container.h:133
void reinit_and_eval(unsigned int cell, const BlockVector< number > *_src_solutions, bool plain)
Definition field_container.h:749
void eval(const BlockVector< number > *_src_solutions, bool plain)
Definition field_container.h:711
FEEValuationDeps()=default
Default constructor.
std::vector< FEEDepPairPtr > fe_eval_old
Collection of FEEvaluation for old solutions (< n-1 state).
Definition field_container.h:143
void integrate_and_distribute(BlockVector< number > *dst_solutions)
Definition field_container.h:818
const SolutionLevel< dim, number > * solution_level
The solution block.
Definition field_container.h:157
std::pair< FEEval< Rank >, EvalFlags > FEEDepPair
dealii::FEEvaluation object and the evaluation flags.
Definition field_container.h:122
unsigned int block_index
The solution block index.
Definition field_container.h:166
The solution vectors and their corresponding matrix free object with respect to some multigrid level.
Definition group_solution_handler.h:61
DependencyType
Internal classification for types of variable dependencies.
Definition type_enums.h:73
@ Current
Definition type_enums.h:80
@ DST
Definition type_enums.h:74
@ OldFour
Definition type_enums.h:84
@ OldThree
Definition type_enums.h:83
@ SRC
Definition type_enums.h:75
@ OldOne
Definition type_enums.h:81
@ OldTwo
Definition type_enums.h:82
TensorRank
Tensor rank of the field.
Definition type_enums.h:52
@ Vector
Definition type_enums.h:55
@ Scalar
Definition type_enums.h:54
dealii::EvaluationFlags::EvaluationFlags EvalFlags
Definition types.h:86