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 || dim == 1,
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
205
206 void
207 integrate();
208
209 void
210 distribute(BlockVector<number> *dst_solutions);
211
212 void
214 };
215
228 FieldContainer(const std::vector<FieldAttributes> &_field_attributes,
229 const SolutionIndexer<dim, number> &_solution_indexer,
230 unsigned int _relative_level,
231 const DependencyMap &dependency_map,
232 const SolveBlock &_solve_block,
233 const MatrixFree<dim, number> &matrix_free);
234
238 void
239 reinit(unsigned int cell);
240
245 void
246 eval(const BlockVector<number> *src_solutions, bool plain);
247
254 void
255 reinit_and_eval(unsigned int cell,
256 const BlockVector<number> *src_solutions,
257 bool plain);
258
263 void
265
269 void
270 integrate();
271
275 void
276 distribute(BlockVector<number> *dst_solutions);
277
281 void
282 distribute(unsigned int field_index, BlockVector<number> *dst_solutions);
283
289 void
291
295 void
296 set_q_point(unsigned int q);
297
301 template <TensorRank Rank, DependencyType type>
302 [[nodiscard]] Value<Rank>
303 get_value(Types::Index field_index) const;
304
308 template <TensorRank Rank>
309 [[nodiscard]] Value<Rank>
310 get_value(Types::Index field_index, DependencyType type) const;
311
315 template <TensorRank Rank, DependencyType type>
316 [[nodiscard]] Gradient<Rank>
317 get_gradient(Types::Index field_index) const;
318
322 template <TensorRank Rank>
323 [[nodiscard]] Gradient<Rank>
324 get_gradient(Types::Index field_index, DependencyType type) const;
325
329 template <TensorRank Rank, DependencyType type>
330 [[nodiscard]] Hessian<Rank>
331 get_hessian(Types::Index field_index) const;
332
336 template <TensorRank Rank>
337 [[nodiscard]] Hessian<Rank>
338 get_hessian(Types::Index field_index, DependencyType type) const;
339
343 template <TensorRank Rank, DependencyType type>
344 [[nodiscard]] Gradient<Rank>
346
350 template <TensorRank Rank>
351 [[nodiscard]] Gradient<Rank>
353
357 template <TensorRank Rank, DependencyType type>
358 [[nodiscard]] Value<Rank>
359 get_laplacian(Types::Index field_index) const;
360
364 template <TensorRank Rank>
365 [[nodiscard]] Value<Rank>
366 get_laplacian(Types::Index field_index, DependencyType type) const;
367
371 template <TensorRank Rank, DependencyType type>
372 [[nodiscard]] ScalarValue
373 get_divergence(Types::Index field_index) const;
374
378 template <TensorRank Rank>
379 [[nodiscard]] ScalarValue
381
385 template <TensorRank Rank, DependencyType type>
386 [[nodiscard]] dealii::SymmetricTensor<2, dim, ScalarValue>
388
392 template <TensorRank Rank>
393 [[nodiscard]] dealii::SymmetricTensor<2, dim, ScalarValue>
395
399 template <TensorRank Rank, DependencyType type>
400 [[nodiscard]] dealii::Tensor<1, (dim == 2 ? 1 : dim), ScalarValue>
401 get_curl(Types::Index field_index) const;
402
406 template <TensorRank Rank>
407 [[nodiscard]] dealii::Tensor<1, (dim == 2 ? 1 : dim), ScalarValue>
408 get_curl(Types::Index field_index, DependencyType type) const;
409
413 [[nodiscard]] dealii::Point<dim, ScalarValue>
414 get_q_point_location() const;
415
419 [[nodiscard]] ScalarValue
420 get_element_volume() const;
421
425 [[nodiscard]] unsigned int
426 get_n_q_points() const;
427
431 template <typename ValType>
432 void
433 set_value_term(Types::Index field_index, const ValType &val);
434
438 template <typename GradType>
439 void
440 set_gradient_term(Types::Index field_index, const GradType &val);
441
445 template <typename ValType>
446 void
447 get_dof_value_to(ValType &destination,
448 Types::Index field_index,
449 unsigned int dof_index);
450
454 template <typename ValType>
455 void
456 submit_dof_value(Types::Index field_index, const ValType &val, unsigned int dof_index);
457
458 constexpr static unsigned int dofs_per_component =
460
461private:
462 template <TensorRank Rank>
463 std::vector<FEEValuationDeps<Rank>> &
465
466 template <TensorRank Rank>
467 const std::vector<FEEValuationDeps<Rank>> &
469
474 void
475 feevaluation_size_valid(Types::Index field_index) const;
476
481 void
482 feevaluation_exists(Types::Index field_index, Types::Index dependency_index) const;
483
488 void
489 access_valid(Types::Index field_index,
490 DependencyType dependency_type,
491 dealii::EvaluationFlags::EvaluationFlags flag) const;
492
496 void
497 submission_valid(Types::Index field_index, DependencyType dependency_type) const;
498
502 const std::vector<FieldAttributes> *field_attributes_ptr;
503
508
512 std::vector<FEEValuationDeps<TensorRank::Scalar>> feeval_deps_scalar;
513
517 std::vector<FEEValuationDeps<TensorRank::Vector>> feeval_deps_vector;
518
523
534
538 unsigned int relative_level;
539
543 unsigned int q_point = 0;
544};
545
546inline static const std::map<DependencyType, std::string> dependency_type_to_string = {
547 {DependencyType::Current, "current"},
548 {DependencyType::OldOne, "old_1" },
549 {DependencyType::OldTwo, "old_2" },
550 {DependencyType::OldThree, "old_3" },
551 {DependencyType::OldFour, "old_4" },
552 {DependencyType::SRC, "src" },
553 {DependencyType::DST, "dst" }
554};
555
556template <unsigned int dim, unsigned int degree, typename number>
557template <TensorRank Rank>
559 const Dependency &dependency,
560 const std::pair<const SolutionLevel<dim, number> *, unsigned int> &mf_id_pair,
561 bool is_dst)
562 : solution_level(mf_id_pair.first)
563 , block_index(mf_id_pair.second)
564{
565 // Make an FEEvaluation if the current solution needs to be evaluated
566 if (dependency.flag)
567 {
568 fe_eval = std::make_shared<FEEDepPair>(FEEval<Rank>(solution_level->matrix_free,
569 block_index),
570 dependency.flag);
571 }
572 // Make FEEvaluations for the the old solutions
573 fe_eval_old.resize(dependency.old_flags.size(), nullptr);
574 for (unsigned int age = 0; age < dependency.old_flags.size(); ++age)
575 {
576 if (dependency.old_flags.at(age)) // kinda redundant... maybe remove?
577 {
578 fe_eval_old[age] =
579 std::make_shared<FEEDepPair>(FEEval<Rank>(solution_level->matrix_free,
581 dependency.old_flags.at(age));
582 }
583 }
584 if (dependency.src_flag || is_dst)
585 {
586 fe_eval_src_dst =
587 std::make_shared<FEEDepPair>(FEEval<Rank>(solution_level->matrix_free,
588 block_index),
589 dependency.src_flag);
590 }
591}
592
593class ExcDepNotInitialized : public dealii::ExceptionBase
594{
595public:
597 : dependency_type(_dependency_type)
598 {}
599
600 const char *
601 what() const noexcept override
602 {
603 return "Access was attempted for a field that was not declared as a "
604 "dependency for the current solve block.\n";
605 }
606
607 void
608 print_info(std::ostream &out) const override
609 {
610 out << what() << std::flush;
611 }
612
613private:
615};
616#ifdef DEBUG
617# define AssertAccessible(fe_eval_pair_ptr, dependency_type) \
618 AssertThrowDebug(block_index != -1, ExcDepNotInitialized(dependency_type)); \
619 AssertThrowDebug(fe_eval_pair_ptr, ExcDepNotInitialized(dependency_type));
620#else
621# define AssertAccessible(fe_eval_pair_ptr, dependency_type)
622#endif
623
624template <unsigned int dim, unsigned int degree, typename number>
625template <TensorRank Rank>
626template <DependencyType Type>
627inline DEAL_II_ALWAYS_INLINE const typename FieldContainer<dim, degree, number>::
628 template FEEval<Rank> &
630{
631 if constexpr (Type == DependencyType::SRC || Type == DependencyType::DST)
632 {
634 return fe_eval_src_dst->first;
635 }
636 else if constexpr (Type == DependencyType::Current)
637 {
639 return fe_eval->first;
640 }
641 else
642 {
645 return fe_eval_old[int(Type) - 1]->first;
646 }
647}
648
649template <unsigned int dim, unsigned int degree, typename number>
650template <TensorRank Rank>
651template <DependencyType Type>
652inline DEAL_II_ALWAYS_INLINE
655{
656 if constexpr (Type == DependencyType::SRC || Type == DependencyType::DST)
657 {
659 return fe_eval_src_dst->first;
660 }
661 else if constexpr (Type == DependencyType::Current)
662 {
664 return fe_eval->first;
665 }
666 else
667 {
670 return fe_eval_old[int(Type) - 1]->first;
671 }
672}
673
674template <unsigned int dim, unsigned int degree, typename number>
675template <TensorRank Rank>
676inline DEAL_II_ALWAYS_INLINE const typename FieldContainer<dim, degree, number>::
677 template FEEval<Rank> &
679 DependencyType type) const
680{
681 if (type == DependencyType::SRC || type == DependencyType::DST)
682 {
684 return fe_eval_src_dst->first;
685 }
686 if (type == DependencyType::Current)
687 {
689 return fe_eval->first;
690 }
691 {
692 AssertThrowDebug(int(type) - 1 < fe_eval_old.size(), ExcDepNotInitialized(type));
693 AssertAccessible(fe_eval_old[int(type) - 1], type);
694 return fe_eval_old[type - 1]->first;
695 }
696}
697
698template <unsigned int dim, unsigned int degree, typename number>
699template <TensorRank Rank>
700inline DEAL_II_ALWAYS_INLINE
703{
704 if (type == DependencyType::SRC || type == DependencyType::DST)
705 {
707 return fe_eval_src_dst->first;
708 }
709 if (type == DependencyType::Current)
710 {
712 return fe_eval->first;
713 }
714 {
715 AssertThrowDebug(int(type) - 1 < fe_eval_old.size(), ExcDepNotInitialized(type));
716 AssertAccessible(fe_eval_old[int(type) - 1], type);
717 return fe_eval_old[int(type) - 1]->first;
718 }
719}
720
721template <unsigned int dim, unsigned int degree, typename number>
722template <TensorRank Rank>
723inline void
725{
726 if (fe_eval)
727 {
728 fe_eval->first.reinit(cell);
729 }
730 for (auto &old_fe_eval : fe_eval_old)
731 {
732 if (old_fe_eval)
733 {
734 old_fe_eval->first.reinit(cell);
735 }
736 }
737 if (fe_eval_src_dst)
738 {
739 fe_eval_src_dst->first.reinit(cell);
740 }
741}
742
743template <unsigned int dim, unsigned int degree, typename number>
744template <TensorRank Rank>
745inline void
747 const BlockVector<number> *_src_solutions,
748 bool plain)
749{
750 // NOTE: `read_dof_values_plain` must be called here so that constraints aren't
751 // implicitly applied. This allows us to have inhomogeneous constraints.
752 if (fe_eval)
753 {
754 fe_eval->first.read_dof_values_plain(solution_level->solutions.block(block_index));
755 fe_eval->first.evaluate(fe_eval->second);
756 }
757 for (unsigned int age = 0; age < fe_eval_old.size(); ++age)
758 {
759 if (FEEDepPairPtr &old_fe_eval = fe_eval_old[age])
760 {
761 old_fe_eval->first.read_dof_values_plain(
762 solution_level->old_solutions[age].block(block_index));
763 old_fe_eval->first.evaluate(old_fe_eval->second);
764 }
765 }
766 if (fe_eval_src_dst && fe_eval_src_dst->second != EvalFlags::nothing)
767 {
768 if (plain)
769 {
770 fe_eval_src_dst->first.read_dof_values_plain(
771 _src_solutions->block(block_index));
772 }
773 else
774 {
775 fe_eval_src_dst->first.read_dof_values(_src_solutions->block(block_index));
776 }
777 fe_eval_src_dst->first.evaluate(fe_eval_src_dst->second);
778 }
779}
780
781template <unsigned int dim, unsigned int degree, typename number>
782template <TensorRank Rank>
783inline void
785 unsigned int cell,
786 const BlockVector<number> *_src_solutions,
787 bool plain)
788{
789 // NOTE: `read_dof_values_plain` must be called here so that constraints aren't
790 // implicitly applied. This allows us to have inhomogeneous constraints.
791 if (fe_eval)
792 {
793 fe_eval->first.reinit(cell);
794 fe_eval->first.read_dof_values_plain(solution_level->solutions.block(block_index));
795 fe_eval->first.evaluate(fe_eval->second);
796 }
797 for (unsigned int age = 0; age < fe_eval_old.size(); ++age)
798 {
799 if (FEEDepPairPtr &old_fe_eval = fe_eval_old[age])
800 {
801 old_fe_eval->first.reinit(cell);
802 old_fe_eval->first.read_dof_values_plain(
803 solution_level->old_solutions[age].block(block_index));
804 old_fe_eval->first.evaluate(old_fe_eval->second);
805 }
806 }
807 if (fe_eval_src_dst)
808 {
809 fe_eval_src_dst->first.reinit(cell);
810 if (fe_eval_src_dst->second != EvalFlags::nothing)
811 {
812 if (plain)
813 {
814 fe_eval_src_dst->first.read_dof_values_plain(
815 _src_solutions->block(block_index));
816 }
817 else
818 {
819 fe_eval_src_dst->first.read_dof_values(_src_solutions->block(block_index));
820 }
821 fe_eval_src_dst->first.evaluate(fe_eval_src_dst->second);
822 }
823 }
824}
825
826template <unsigned int dim, unsigned int degree, typename number>
827template <TensorRank Rank>
828inline void
830{
831 if (fe_eval_src_dst && fe_eval_src_dst->second != EvalFlags::nothing)
832 {
833 fe_eval_src_dst->first.evaluate(fe_eval_src_dst->second);
834 }
835}
836
837template <unsigned int dim, unsigned int degree, typename number>
838template <TensorRank Rank>
839inline void
847
848template <unsigned int dim, unsigned int degree, typename number>
849template <TensorRank Rank>
850inline void
852 BlockVector<number> *dst_solutions)
853{
854 if (fe_eval_src_dst)
855 {
856 fe_eval_src_dst->first.distribute_local_to_global(
857 dst_solutions->block(block_index));
858 }
859}
860
861template <unsigned int dim, unsigned int degree, typename number>
862template <TensorRank Rank>
863inline void
865 BlockVector<number> *dst_solutions)
866{
867 if (fe_eval_src_dst)
868 {
869 fe_eval_src_dst->first.integrate_scatter(integration_flags,
870 dst_solutions->block(block_index));
871 }
872}
873
874template <unsigned int dim, unsigned int degree, typename number>
875inline void
877{
878 for (auto &fe_eval : feeval_deps_scalar)
879 {
880 fe_eval.reinit(cell);
881 }
882 for (auto &fe_eval : feeval_deps_vector)
883 {
884 fe_eval.reinit(cell);
885 }
886 shared_feeval_scalar.reinit(cell);
887}
888
889template <unsigned int dim, unsigned int degree, typename number>
890inline void
892 bool plain)
893{
894 for (auto &fe_eval : feeval_deps_scalar)
895 {
896 fe_eval.eval(src_solutions, plain);
897 }
898 for (auto &fe_eval : feeval_deps_vector)
899 {
900 fe_eval.eval(src_solutions, plain);
901 }
902 // Don't eval `shared_feeval_scalar` because we only use it for information.
903}
904
905template <unsigned int dim, unsigned int degree, typename number>
906inline void
908 unsigned int cell,
909 const BlockVector<number> *src_solutions,
910 bool plain)
911{
912 for (auto &fe_eval : feeval_deps_scalar)
913 {
914 fe_eval.reinit_and_eval(cell, src_solutions, plain);
915 }
916 for (auto &fe_eval : feeval_deps_vector)
917 {
918 fe_eval.reinit_and_eval(cell, src_solutions, plain);
919 }
920 // Don't eval `shared_feeval_scalar` because we only use it for information.
921 shared_feeval_scalar.reinit(cell);
922}
923
924template <unsigned int dim, unsigned int degree, typename number>
925inline void
927{
928 for (auto &fe_eval : feeval_deps_scalar)
929 {
930 fe_eval.eval_without_read();
931 }
932 for (auto &fe_eval : feeval_deps_vector)
933 {
934 fe_eval.eval_without_read();
935 }
936}
937
938template <unsigned int dim, unsigned int degree, typename number>
939inline void
941{
942 const std::vector<FieldAttributes> &field_attributes = *field_attributes_ptr;
943 for (const Types::Index &field_index : solve_block->field_indices)
944 {
945 if (field_attributes[field_index].field_type == TensorRank::Scalar)
946 {
947 feeval_deps_scalar[field_index].integrate();
948 }
949 else /* vector */
950 {
951 feeval_deps_vector[field_index].integrate();
952 }
953 }
954 // Don't integrate `shared_feeval_scalar` because we only use it for information.
955}
956
957template <unsigned int dim, unsigned int degree, typename number>
958inline void
960{
961 const std::vector<FieldAttributes> &field_attributes = *field_attributes_ptr;
962 for (const Types::Index &field_index : solve_block->field_indices)
963 {
964 if (field_attributes[field_index].field_type == TensorRank::Scalar)
965 {
966 feeval_deps_scalar[field_index].distribute(dst_solutions);
967 }
968 else /* vector */
969 {
970 feeval_deps_vector[field_index].distribute(dst_solutions);
971 }
972 }
973}
974
975template <unsigned int dim, unsigned int degree, typename number>
976inline void
978 BlockVector<number> *dst_solutions)
979{
980 const std::vector<FieldAttributes> &field_attributes = *field_attributes_ptr;
981 if (field_attributes[field_index].field_type == TensorRank::Scalar)
982 {
983 feeval_deps_scalar[field_index].distribute(dst_solutions);
984 }
985 else /* vector */
986 {
987 feeval_deps_vector[field_index].distribute(dst_solutions);
988 }
989}
990
991template <unsigned int dim, unsigned int degree, typename number>
992inline void
994 BlockVector<number> *dst_solutions)
995{
996 const std::vector<FieldAttributes> &field_attributes = *field_attributes_ptr;
997 for (const Types::Index &field_index : solve_block->field_indices)
998 {
999 if (field_attributes[field_index].field_type == TensorRank::Scalar)
1000 {
1001 feeval_deps_scalar[field_index].integrate_and_distribute(dst_solutions);
1002 }
1003 else /* vector */
1004 {
1005 feeval_deps_vector[field_index].integrate_and_distribute(dst_solutions);
1006 }
1007 }
1008}
1009
1010template <unsigned int dim, unsigned int degree, typename number>
1011inline DEAL_II_ALWAYS_INLINE void
1013{
1014 q_point = q;
1015}
1016
1017// there are two catches we can do here.
1018// 1. Dependencies for the dependency type (current, old, src/dst) don't exist.
1019// 2. Dependency is not initialized for values/gradients.
1020// We catch these separately to give more informative error messages.
1021#define GetterTempl(dependency_type) template get<dependency_type>()
1022#define GetterNoTempl(dependency_type) get(dependency_type)
1023#define ReturnGetter(get_handle, Rank, field_index, dependency_type, getter) \
1024 try \
1025 { \
1026 AssertThrowDebug((field_index) < get_relevant_feeval_vector<Rank>().size(), \
1027 dealii::ExcMessage("Error: Field index " + \
1028 std::to_string(field_index) + \
1029 " is not associated with any field.")); \
1030 return get_relevant_feeval_vector<Rank>()[field_index] \
1031 .getter(dependency_type) \
1032 .get_handle(q_point); \
1033 } \
1034 catch (const ExcDepNotInitialized &e) \
1035 { \
1036 std::cerr << "Error when trying to access field with index " << (field_index) \
1037 << " and dependency type " \
1038 << dependency_type_to_string.at(dependency_type) << ":\n" \
1039 << e.what() \
1040 << "Ensure that each dependency is requested in the solve block.\n" \
1041 << std::flush; \
1042 throw; \
1043 } \
1044 catch (const dealii::internal::ExcAccessToUninitializedField &e) \
1045 { \
1046 std::cerr << "Error when trying to access field with index " << (field_index) \
1047 << " and dependency type " \
1048 << dependency_type_to_string.at(dependency_type) \
1049 << ":\nAccess was attempted for values or gradients that were " \
1050 "not requested.\n" \
1051 << std::flush; \
1052 throw; \
1053 }
1054
1055template <unsigned int dim, unsigned int degree, typename number>
1056template <TensorRank Rank, DependencyType type>
1057inline DEAL_II_ALWAYS_INLINE
1060{
1061 ReturnGetter(get_value, Rank, field_index, type, GetterTempl);
1062}
1063
1064template <unsigned int dim, unsigned int degree, typename number>
1065template <TensorRank Rank>
1066inline DEAL_II_ALWAYS_INLINE
1069 DependencyType type) const
1070{
1071 ReturnGetter(get_value, Rank, field_index, type, GetterNoTempl);
1072}
1073
1074template <unsigned int dim, unsigned int degree, typename number>
1075template <TensorRank Rank, DependencyType type>
1076inline DEAL_II_ALWAYS_INLINE
1079{
1080 ReturnGetter(get_gradient, Rank, field_index, type, GetterTempl);
1081}
1082
1083template <unsigned int dim, unsigned int degree, typename number>
1084template <TensorRank Rank>
1085inline DEAL_II_ALWAYS_INLINE
1088 DependencyType type) const
1089{
1090 ReturnGetter(get_gradient, Rank, field_index, type, GetterNoTempl);
1091}
1092
1093template <unsigned int dim, unsigned int degree, typename number>
1094template <TensorRank Rank, DependencyType type>
1095inline DEAL_II_ALWAYS_INLINE
1098{
1099 ReturnGetter(get_hessian, Rank, field_index, type, GetterTempl);
1100}
1101
1102template <unsigned int dim, unsigned int degree, typename number>
1103template <TensorRank Rank>
1104inline DEAL_II_ALWAYS_INLINE
1107 DependencyType type) const
1108{
1109 ReturnGetter(get_hessian, Rank, field_index, type, GetterNoTempl);
1110}
1111
1112template <unsigned int dim, unsigned int degree, typename number>
1113template <TensorRank Rank, DependencyType type>
1114inline DEAL_II_ALWAYS_INLINE
1121
1122template <unsigned int dim, unsigned int degree, typename number>
1123template <TensorRank Rank>
1124inline DEAL_II_ALWAYS_INLINE
1131
1132template <unsigned int dim, unsigned int degree, typename number>
1133template <TensorRank Rank, DependencyType type>
1134inline DEAL_II_ALWAYS_INLINE
1137{
1138 ReturnGetter(get_laplacian, Rank, field_index, type, GetterTempl);
1139}
1140
1141template <unsigned int dim, unsigned int degree, typename number>
1142template <TensorRank Rank>
1143inline DEAL_II_ALWAYS_INLINE
1146 DependencyType type) const
1147{
1148 ReturnGetter(get_laplacian, Rank, field_index, type, GetterNoTempl);
1149}
1150
1151template <unsigned int dim, unsigned int degree, typename number>
1152template <TensorRank Rank, DependencyType type>
1153inline DEAL_II_ALWAYS_INLINE typename FieldContainer<dim, degree, number>::ScalarValue
1155{
1156 static_assert(Rank == 1, "Divergences are only available for vector fields");
1157 ReturnGetter(get_divergence, Rank, field_index, type, GetterTempl);
1158}
1159
1160template <unsigned int dim, unsigned int degree, typename number>
1161template <TensorRank Rank>
1162inline DEAL_II_ALWAYS_INLINE typename FieldContainer<dim, degree, number>::ScalarValue
1164 DependencyType type) const
1165{
1166 static_assert(Rank == 1, "Divergences are only available for vector fields");
1167 ReturnGetter(get_divergence, Rank, field_index, type, GetterNoTempl);
1168}
1169
1170template <unsigned int dim, unsigned int degree, typename number>
1171template <TensorRank Rank, DependencyType type>
1172inline DEAL_II_ALWAYS_INLINE dealii::
1173 SymmetricTensor<2, dim, typename FieldContainer<dim, degree, number>::ScalarValue>
1175 Types::Index field_index) const
1176{
1177 static_assert(Rank == 1, "Symmetric gradients are only available for vector fields");
1178 ReturnGetter(get_symmetric_gradient, Rank, field_index, type, GetterTempl);
1179}
1180
1181template <unsigned int dim, unsigned int degree, typename number>
1182template <TensorRank Rank>
1183inline DEAL_II_ALWAYS_INLINE dealii::
1184 SymmetricTensor<2, dim, typename FieldContainer<dim, degree, number>::ScalarValue>
1186 DependencyType type) const
1187{
1188 static_assert(Rank == 1, "Symmetric gradients are only available for vector fields");
1189 ReturnGetter(get_symmetric_gradient, Rank, field_index, type, GetterNoTempl);
1190}
1191
1192template <unsigned int dim, unsigned int degree, typename number>
1193template <TensorRank Rank, DependencyType type>
1194inline DEAL_II_ALWAYS_INLINE
1195 dealii::Tensor<1,
1196 (dim == 2 ? 1 : dim),
1199{
1200 static_assert(Rank == 1, "Curl is only available for vector fields");
1202 dim > 1,
1203 dealii::ExcMessage(
1204 "Curl is only available for vector fields with dimension greater than 1."));
1205 ReturnGetter(get_curl, Rank, field_index, type, GetterTempl);
1206}
1207
1208template <unsigned int dim, unsigned int degree, typename number>
1209template <TensorRank Rank>
1210inline DEAL_II_ALWAYS_INLINE
1211 dealii::Tensor<1,
1212 (dim == 2 ? 1 : dim),
1215 DependencyType type) const
1216{
1217 static_assert(Rank == 1, "Curl is only available for vector fields");
1219 dim > 1,
1220 dealii::ExcMessage(
1221 "Curl is only available for vector fields with dimension greater than 1."));
1222 ReturnGetter(get_curl, Rank, field_index, type, GetterNoTempl);
1223}
1224
1225template <unsigned int dim, unsigned int degree, typename number>
1226inline DEAL_II_ALWAYS_INLINE
1227 dealii::Point<dim, typename FieldContainer<dim, degree, number>::ScalarValue>
1232
1233template <unsigned int dim, unsigned int degree, typename number>
1234inline DEAL_II_ALWAYS_INLINE typename FieldContainer<dim, degree, number>::ScalarValue
1240
1241template <unsigned int dim, unsigned int degree, typename number>
1242inline DEAL_II_ALWAYS_INLINE unsigned int
1247
1248template <unsigned int dim, unsigned int degree, typename number>
1249template <typename ValType>
1250inline DEAL_II_ALWAYS_INLINE void
1252 const ValType &val)
1253{
1254 AssertThrowDebug((field_index) <
1256 dealii::ExcMessage("Error: Field index " +
1257 std::to_string(field_index) +
1258 " is not associated with any field."));
1259 auto &relevant_feeval_vector =
1261 try
1262 {
1263 relevant_feeval_vector.template get<DependencyType::DST>().submit_value(val,
1264 q_point);
1265 relevant_feeval_vector.integration_flags |= dealii::EvaluationFlags::values;
1266 }
1267 catch (...)
1268 {
1269 std::cerr << "Error when trying to submit value for field with index "
1270 << (field_index)
1271 << ": Error: Submission for field not part of this solve block.\n"
1272 << std::flush;
1273 throw;
1274 }
1275}
1276
1277template <unsigned int dim, unsigned int degree, typename number>
1278template <typename GradType>
1279inline DEAL_II_ALWAYS_INLINE void
1281 const GradType &val)
1282{
1283 AssertThrowDebug((field_index) <
1285 dealii::ExcMessage("Error: Field index " +
1286 std::to_string(field_index) +
1287 " is not associated with any field."));
1288 auto &relevant_feeval_vector =
1290 try
1291 {
1292 relevant_feeval_vector.template get<DependencyType::DST>().submit_gradient(val,
1293 q_point);
1294 relevant_feeval_vector.integration_flags |= dealii::EvaluationFlags::gradients;
1295 }
1296 catch (...)
1297 {
1298 std::cerr << "Error when trying to submit gradient for field with index "
1299 << (field_index)
1300 << ": Error: Submission for field not part of this solve block.\n"
1301 << std::flush;
1302 throw;
1303 }
1304}
1305
1306template <unsigned int dim, unsigned int degree, typename number>
1307template <typename ValType>
1308inline DEAL_II_ALWAYS_INLINE void
1310 Types::Index field_index,
1311 unsigned int dof_index)
1312{
1313 AssertThrowDebug((field_index) <
1315 dealii::ExcMessage("Error: Field index " +
1316 std::to_string(field_index) +
1317 " is not associated with any field."));
1318 auto &relevant_feeval_vector =
1320 try
1321 {
1322 destination =
1323 relevant_feeval_vector.template get<DependencyType::DST>().get_dof_value(
1324 dof_index);
1325 }
1326 catch (...)
1327 {
1328 std::cerr << "Error when trying to get dof value for field with index "
1329 << (field_index) << ": Error: Field not part of this solve block.\n"
1330 << std::flush;
1331 throw;
1332 }
1333}
1334
1335template <unsigned int dim, unsigned int degree, typename number>
1336template <typename ValType>
1337inline DEAL_II_ALWAYS_INLINE void
1339 const ValType &val,
1340 unsigned int dof_index)
1341{
1342 AssertThrowDebug((field_index) <
1344 dealii::ExcMessage("Error: Field index " +
1345 std::to_string(field_index) +
1346 " is not associated with any field."));
1347 auto &relevant_feeval_vector =
1349 try
1350 {
1351 relevant_feeval_vector.template get<DependencyType::DST>()
1352 .submit_dof_value(val, dof_index);
1353 }
1354 catch (...)
1355 {
1356 std::cerr << "Error when trying to submit dof value for field with index "
1357 << (field_index)
1358 << ": Error: Submission for field not part of this solve block.\n"
1359 << std::flush;
1360 throw;
1361 }
1362}
1363
1364template <unsigned int dim, unsigned int degree, typename number>
1365template <TensorRank Rank>
1366inline DEAL_II_ALWAYS_INLINE std::vector<
1367 typename FieldContainer<dim, degree, number>::template FEEValuationDeps<Rank>> &
1369{
1370 if constexpr (Rank == TensorRank::Scalar || dim == 1)
1371 {
1372 return feeval_deps_scalar;
1373 }
1374 else if constexpr (Rank == TensorRank::Vector)
1375 {
1376 return feeval_deps_vector;
1377 }
1378}
1379
1380template <unsigned int dim, unsigned int degree, typename number>
1381template <TensorRank Rank>
1382inline DEAL_II_ALWAYS_INLINE const std::vector<
1383 typename FieldContainer<dim, degree, number>::template FEEValuationDeps<Rank>> &
1385{
1386 if constexpr (Rank == TensorRank::Scalar || dim == 1)
1387 {
1388 return feeval_deps_scalar;
1389 }
1390 else if constexpr (Rank == TensorRank::Vector)
1391 {
1392 return feeval_deps_vector;
1393 }
1394}
1395
1396template <unsigned int dim, unsigned int degree, typename number>
1397inline DEAL_II_ALWAYS_INLINE void
1399 [[maybe_unused]] Types::Index field_index) const
1400{
1401 // TODO
1402}
1403
1404template <unsigned int dim, unsigned int degree, typename number>
1405inline DEAL_II_ALWAYS_INLINE void
1407 [[maybe_unused]] Types::Index field_index,
1408 [[maybe_unused]] Types::Index dependency_index) const
1409{
1410 // TODO
1411}
1412
1413template <unsigned int dim, unsigned int degree, typename number>
1414inline DEAL_II_ALWAYS_INLINE void
1416 [[maybe_unused]] Types::Index field_index,
1417 [[maybe_unused]] DependencyType dependency_type,
1418 [[maybe_unused]] dealii::EvaluationFlags::EvaluationFlags flag) const
1419{
1420 // TODO
1421}
1422
1423template <unsigned int dim, unsigned int degree, typename number>
1424inline DEAL_II_ALWAYS_INLINE void
1426 [[maybe_unused]] Types::Index field_index,
1427 [[maybe_unused]] DependencyType dependency_type) const
1428{
1429 // TODO
1430}
1431
1432#undef AssertAccessible
1433#undef ReturnGetter
1434#undef GetterTempl
1435#undef GetterNoTempl
1436
1437PRISMS_PF_END_NAMESPACE
Definition field_container.h:594
ExcDepNotInitialized(DependencyType _dependency_type)
Definition field_container.h:596
const char * what() const noexcept override
Definition field_container.h:601
void print_info(std::ostream &out) const override
Definition field_container.h:608
DependencyType dependency_type
Definition field_container.h:614
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 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 gradient of the specified scalar/vector field.
void submit_dof_value(Types::Index field_index, const ValType &val, unsigned int dof_index)
Set the dof values directly at a node.
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:891
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:1415
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:533
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:538
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:1012
void submission_valid(Types::Index field_index, DependencyType dependency_type) const
Check that a value is valid for submission.
Definition field_container.h:1425
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:543
static constexpr unsigned int dofs_per_component
Definition field_container.h:458
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 values.
Definition field_container.h:940
dealii::SymmetricTensor< 2, dim, ScalarValue > get_symmetric_gradient(Types::Index field_index) const
Return the symmetric gradient of the specified field.
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:907
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:517
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 values to a solution vector.
Definition field_container.h:959
void reinit(unsigned int cell)
Initialize based on cell for all dependencies.
Definition field_container.h:876
Gradient< Rank > get_gradient(Types::Index field_index, DependencyType type) const
Return the gradient of the specified field.
std::conditional_t< Rank==TensorRank::Scalar||dim==1, ScalarValue, dealii::Tensor< int(Rank), dim, ScalarValue > > Value
Definition field_container.h:60
std::vector< FEEValuationDeps< TensorRank::Scalar > > feeval_deps_scalar
Collection of FEEvaluation dependencies for scalar fields.
Definition field_container.h:512
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:1406
void integrate_and_distribute(BlockVector< number > *dst_solutions)
Integrate the values and distribute from local to global.
Definition field_container.h:993
const std::vector< FieldAttributes > * field_attributes_ptr
Pointer to the vector of field attributes.
Definition field_container.h:502
void eval_without_read()
Evaluate src to quad pts based on dependency flags. Assume feeval has local dof values already.
Definition field_container.h:926
void get_dof_value_to(ValType &destination, Types::Index field_index, unsigned int dof_index)
Get the dof values directly from a node.
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:1243
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:1228
ScalarValue get_element_volume() const
Return the quadrature point location.
Definition field_container.h:1235
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:522
const SolutionIndexer< dim, number > * solution_indexer
Pointer to the solution indexer.
Definition field_container.h:507
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:1398
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:26
#define AssertAccessible(fe_eval_pair_ptr, dependency_type)
Definition field_container.h:621
#define GetterNoTempl(dependency_type)
Definition field_container.h:1022
#define GetterTempl(dependency_type)
Definition field_container.h:1021
static const std::map< DependencyType, std::string > dependency_type_to_string
Definition field_container.h:546
#define ReturnGetter(get_handle, Rank, field_index, dependency_type, getter)
Definition field_container.h:1023
@ 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 eval_without_read()
Definition field_container.h:829
void integrate()
Definition field_container.h:840
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:851
void reinit(unsigned int cell)
Definition field_container.h:724
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:784
void eval(const BlockVector< number > *_src_solutions, bool plain)
Definition field_container.h:746
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:864
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:94
@ Current
Definition type_enums.h:101
@ DST
Definition type_enums.h:95
@ OldFour
Definition type_enums.h:105
@ OldThree
Definition type_enums.h:104
@ SRC
Definition type_enums.h:96
@ OldOne
Definition type_enums.h:102
@ OldTwo
Definition type_enums.h:103
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