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>
24#include <prismspf/config.h>
45template <
unsigned int dim,
unsigned int degree,
typename number>
59 template <TensorRank Rank>
64 template <TensorRank Rank>
67 template <TensorRank Rank>
70 template <TensorRank Rank>
72 dealii::FEEvaluation<dim,
75 dealii::Tensor<int(Rank), dim>::n_independent_components,
82 template <
typename ValType>
85 if constexpr (std::is_same_v<ValType, ScalarValue> ||
86 std::is_same_v<ValType, typename ScalarValue::value_type>)
99 template <
typename GradType>
102 if constexpr (std::is_same_v<GradType, ScalarValue>)
116 template <TensorRank Rank>
178 template <DependencyType type>
182 template <DependencyType type>
193 reinit(
unsigned int cell);
225 FieldContainer(
const std::vector<FieldAttributes> &_field_attributes,
227 unsigned int _relative_level,
230 const MatrixFree<dim, number> &matrix_free);
236 reinit(
unsigned int cell);
285 template <TensorRank Rank, DependencyType type>
292 template <TensorRank Rank>
299 template <TensorRank Rank, DependencyType type>
306 template <TensorRank Rank>
313 template <TensorRank Rank, DependencyType type>
320 template <TensorRank Rank>
327 template <TensorRank Rank, DependencyType type>
334 template <TensorRank Rank>
341 template <TensorRank Rank, DependencyType type>
348 template <TensorRank Rank>
355 template <TensorRank Rank, DependencyType type>
362 template <TensorRank Rank>
369 template <TensorRank Rank, DependencyType type>
370 [[nodiscard]] dealii::SymmetricTensor<2, dim, ScalarValue>
376 template <TensorRank Rank>
377 [[nodiscard]] dealii::SymmetricTensor<2, dim, ScalarValue>
383 template <TensorRank Rank, DependencyType type>
384 [[nodiscard]] dealii::Tensor<1, (dim == 2 ? 1 : dim),
ScalarValue>
390 template <TensorRank Rank>
391 [[nodiscard]] dealii::Tensor<1, (dim == 2 ? 1 : dim),
ScalarValue>
397 [[nodiscard]] dealii::Point<dim, ScalarValue>
409 [[nodiscard]]
unsigned int
415 template <
typename ValType>
422 template <
typename GradType>
427 template <TensorRank Rank>
428 std::vector<FEEValuationDeps<Rank>> &
431 template <TensorRank Rank>
432 const std::vector<FEEValuationDeps<Rank>> &
456 dealii::EvaluationFlags::EvaluationFlags flag)
const;
521template <
unsigned int dim,
unsigned int degree,
typename number>
522template <TensorRank Rank>
533 fe_eval = std::make_shared<FEEDepPair>(FEEval<Rank>(solution_level->matrix_free,
539 for (
unsigned int age = 0; age < dependency.
old_flags.size(); ++age)
552 std::make_shared<FEEDepPair>(FEEval<Rank>(solution_level->matrix_free,
554 dependency.src_flag);
568 return "Access was attempted for a field that was not declared as a "
569 "dependency for the current solve block.\n";
575 out <<
what() << std::flush;
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));
586# define AssertAccessible(fe_eval_pair_ptr, dependency_type)
589template <
unsigned int dim,
unsigned int degree,
typename number>
590template <TensorRank Rank>
591template <DependencyType Type>
614template <
unsigned int dim,
unsigned int degree,
typename number>
615template <TensorRank Rank>
616template <DependencyType Type>
617inline DEAL_II_ALWAYS_INLINE
639template <
unsigned int dim,
unsigned int degree,
typename number>
640template <TensorRank Rank>
663template <
unsigned int dim,
unsigned int degree,
typename number>
664template <TensorRank Rank>
665inline DEAL_II_ALWAYS_INLINE
686template <
unsigned int dim,
unsigned int degree,
typename number>
687template <TensorRank Rank>
699 old_fe_eval->first.reinit(cell);
708template <
unsigned int dim,
unsigned int degree,
typename number>
709template <TensorRank Rank>
722 for (
unsigned int age = 0; age <
fe_eval_old.size(); ++age)
726 old_fe_eval->first.read_dof_values_plain(
728 old_fe_eval->first.evaluate(old_fe_eval->second);
746template <
unsigned int dim,
unsigned int degree,
typename number>
747template <TensorRank Rank>
762 for (
unsigned int age = 0; age <
fe_eval_old.size(); ++age)
766 old_fe_eval->first.reinit(cell);
767 old_fe_eval->first.read_dof_values_plain(
769 old_fe_eval->first.evaluate(old_fe_eval->second);
791template <
unsigned int dim,
unsigned int degree,
typename number>
792template <TensorRank Rank>
802template <
unsigned int dim,
unsigned int degree,
typename number>
803template <TensorRank Rank>
815template <
unsigned int dim,
unsigned int degree,
typename number>
816template <TensorRank Rank>
828template <
unsigned int dim,
unsigned int degree,
typename number>
834 fe_eval.reinit(cell);
838 fe_eval.reinit(cell);
843template <
unsigned int dim,
unsigned int degree,
typename number>
850 fe_eval.eval(src_solutions, plain);
854 fe_eval.eval(src_solutions, plain);
859template <
unsigned int dim,
unsigned int degree,
typename number>
868 fe_eval.reinit_and_eval(cell, src_solutions, plain);
872 fe_eval.reinit_and_eval(cell, src_solutions, plain);
878template <
unsigned int dim,
unsigned int degree,
typename number>
897template <
unsigned int dim,
unsigned int degree,
typename number>
916template <
unsigned int dim,
unsigned int degree,
typename number>
937template <
unsigned int dim,
unsigned int degree,
typename number>
938inline DEAL_II_ALWAYS_INLINE
void
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) \
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); \
961 catch (const ExcDepNotInitialized &e) \
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" \
967 << "Ensure that each dependency is requested in the solve block.\n" \
971 catch (const dealii::internal::ExcAccessToUninitializedField &e) \
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 " \
982template <
unsigned int dim,
unsigned int degree,
typename number>
983template <TensorRank Rank, DependencyType type>
984inline DEAL_II_ALWAYS_INLINE
991template <
unsigned int dim,
unsigned int degree,
typename number>
992template <TensorRank Rank>
993inline DEAL_II_ALWAYS_INLINE
1001template <
unsigned int dim,
unsigned int degree,
typename number>
1002template <TensorRank Rank, DependencyType type>
1003inline DEAL_II_ALWAYS_INLINE
1010template <
unsigned int dim,
unsigned int degree,
typename number>
1011template <TensorRank Rank>
1012inline DEAL_II_ALWAYS_INLINE
1020template <
unsigned int dim,
unsigned int degree,
typename number>
1021template <TensorRank Rank, DependencyType type>
1022inline DEAL_II_ALWAYS_INLINE
1029template <
unsigned int dim,
unsigned int degree,
typename number>
1030template <TensorRank Rank>
1031inline DEAL_II_ALWAYS_INLINE
1039template <
unsigned int dim,
unsigned int degree,
typename number>
1040template <TensorRank Rank, DependencyType type>
1041inline DEAL_II_ALWAYS_INLINE
1049template <
unsigned int dim,
unsigned int degree,
typename number>
1050template <TensorRank Rank>
1051inline DEAL_II_ALWAYS_INLINE
1059template <
unsigned int dim,
unsigned int degree,
typename number>
1060template <TensorRank Rank, DependencyType type>
1061inline DEAL_II_ALWAYS_INLINE
1068template <
unsigned int dim,
unsigned int degree,
typename number>
1069template <TensorRank Rank>
1070inline DEAL_II_ALWAYS_INLINE
1078template <
unsigned int dim,
unsigned int degree,
typename number>
1079template <TensorRank Rank, DependencyType type>
1083 static_assert(Rank == 1,
"Divergences are only available for vector fields");
1087template <
unsigned int dim,
unsigned int degree,
typename number>
1088template <TensorRank Rank>
1093 static_assert(Rank == 1,
"Divergences are only available for vector fields");
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>
1104 static_assert(Rank == 1,
"Symmetric gradients are only available for vector fields");
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>
1115 static_assert(Rank == 1,
"Symmetric gradients are only available for vector fields");
1119template <
unsigned int dim,
unsigned int degree,
typename number>
1120template <TensorRank Rank, DependencyType type>
1121inline DEAL_II_ALWAYS_INLINE
1123 (dim == 2 ? 1 : dim),
1127 static_assert(Rank == 1,
"Curl is only available for vector fields");
1131 "Curl is only available for vector fields with dimension greater than 1."));
1135template <
unsigned int dim,
unsigned int degree,
typename number>
1136template <TensorRank Rank>
1137inline DEAL_II_ALWAYS_INLINE
1139 (dim == 2 ? 1 : dim),
1144 static_assert(Rank == 1,
"Curl is only available for vector fields");
1148 "Curl is only available for vector fields with dimension greater than 1."));
1152template <
unsigned int dim,
unsigned int degree,
typename number>
1153inline DEAL_II_ALWAYS_INLINE
1154 dealii::Point<dim, typename FieldContainer<dim, degree, number>::ScalarValue>
1160template <
unsigned int dim,
unsigned int degree,
typename number>
1168template <
unsigned int dim,
unsigned int degree,
typename number>
1169inline DEAL_II_ALWAYS_INLINE
unsigned int
1175template <
unsigned int dim,
unsigned int degree,
typename number>
1176template <
typename ValType>
1177inline DEAL_II_ALWAYS_INLINE
void
1183 dealii::ExcMessage(
"Error: Field index " +
1184 std::to_string(field_index) +
1185 " is not associated with any field."));
1186 auto &relevant_feeval_vector =
1190 relevant_feeval_vector.template get<DependencyType::DST>().submit_value(val,
1192 relevant_feeval_vector.integration_flags |= dealii::EvaluationFlags::values;
1196 std::cerr <<
"Error when trying to submit value for field with index "
1198 <<
": Error: Submission for field not part of this solve block.\n"
1204template <
unsigned int dim,
unsigned int degree,
typename number>
1205template <
typename GradType>
1206inline DEAL_II_ALWAYS_INLINE
void
1208 const GradType &val)
1212 dealii::ExcMessage(
"Error: Field index " +
1213 std::to_string(field_index) +
1214 " is not associated with any field."));
1215 auto &relevant_feeval_vector =
1219 relevant_feeval_vector.template get<DependencyType::DST>().submit_gradient(val,
1221 relevant_feeval_vector.integration_flags |= dealii::EvaluationFlags::gradients;
1225 std::cerr <<
"Error when trying to submit gradient for field with index "
1227 <<
": Error: Submission for field not part of this solve block.\n"
1233template <
unsigned int dim,
unsigned int degree,
typename number>
1234template <TensorRank Rank>
1235inline DEAL_II_ALWAYS_INLINE std::vector<
1249template <
unsigned int dim,
unsigned int degree,
typename number>
1250template <TensorRank Rank>
1251inline DEAL_II_ALWAYS_INLINE
const std::vector<
1265template <
unsigned int dim,
unsigned int degree,
typename number>
1266inline DEAL_II_ALWAYS_INLINE
void
1273template <
unsigned int dim,
unsigned int degree,
typename number>
1274inline DEAL_II_ALWAYS_INLINE
void
1282template <
unsigned int dim,
unsigned int degree,
typename number>
1283inline DEAL_II_ALWAYS_INLINE
void
1287 [[maybe_unused]] dealii::EvaluationFlags::EvaluationFlags flag)
const
1292template <
unsigned int dim,
unsigned int degree,
typename number>
1293inline DEAL_II_ALWAYS_INLINE
void
1301#undef AssertAccessible
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