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);
228 FieldContainer(
const std::vector<FieldAttributes> &_field_attributes,
230 unsigned int _relative_level,
233 const MatrixFree<dim, number> &matrix_free);
239 reinit(
unsigned int cell);
301 template <TensorRank Rank, DependencyType type>
308 template <TensorRank Rank>
315 template <TensorRank Rank, DependencyType type>
322 template <TensorRank Rank>
329 template <TensorRank Rank, DependencyType type>
336 template <TensorRank Rank>
343 template <TensorRank Rank, DependencyType type>
350 template <TensorRank Rank>
357 template <TensorRank Rank, DependencyType type>
364 template <TensorRank Rank>
371 template <TensorRank Rank, DependencyType type>
378 template <TensorRank Rank>
385 template <TensorRank Rank, DependencyType type>
386 [[nodiscard]] dealii::SymmetricTensor<2, dim, ScalarValue>
392 template <TensorRank Rank>
393 [[nodiscard]] dealii::SymmetricTensor<2, dim, ScalarValue>
399 template <TensorRank Rank, DependencyType type>
400 [[nodiscard]] dealii::Tensor<1, (dim == 2 ? 1 : dim),
ScalarValue>
406 template <TensorRank Rank>
407 [[nodiscard]] dealii::Tensor<1, (dim == 2 ? 1 : dim),
ScalarValue>
413 [[nodiscard]] dealii::Point<dim, ScalarValue>
425 [[nodiscard]]
unsigned int
431 template <
typename ValType>
438 template <
typename GradType>
445 template <
typename ValType>
449 unsigned int dof_index);
454 template <
typename ValType>
462 template <TensorRank Rank>
463 std::vector<FEEValuationDeps<Rank>> &
466 template <TensorRank Rank>
467 const std::vector<FEEValuationDeps<Rank>> &
491 dealii::EvaluationFlags::EvaluationFlags flag)
const;
556template <
unsigned int dim,
unsigned int degree,
typename number>
557template <TensorRank Rank>
568 fe_eval = std::make_shared<FEEDepPair>(FEEval<Rank>(solution_level->matrix_free,
574 for (
unsigned int age = 0; age < dependency.
old_flags.size(); ++age)
587 std::make_shared<FEEDepPair>(FEEval<Rank>(solution_level->matrix_free,
589 dependency.src_flag);
603 return "Access was attempted for a field that was not declared as a "
604 "dependency for the current solve block.\n";
610 out <<
what() << std::flush;
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));
621# define AssertAccessible(fe_eval_pair_ptr, dependency_type)
624template <
unsigned int dim,
unsigned int degree,
typename number>
625template <TensorRank Rank>
626template <DependencyType Type>
649template <
unsigned int dim,
unsigned int degree,
typename number>
650template <TensorRank Rank>
651template <DependencyType Type>
652inline DEAL_II_ALWAYS_INLINE
674template <
unsigned int dim,
unsigned int degree,
typename number>
675template <TensorRank Rank>
698template <
unsigned int dim,
unsigned int degree,
typename number>
699template <TensorRank Rank>
700inline DEAL_II_ALWAYS_INLINE
721template <
unsigned int dim,
unsigned int degree,
typename number>
722template <TensorRank Rank>
734 old_fe_eval->first.reinit(cell);
743template <
unsigned int dim,
unsigned int degree,
typename number>
744template <TensorRank Rank>
757 for (
unsigned int age = 0; age <
fe_eval_old.size(); ++age)
761 old_fe_eval->first.read_dof_values_plain(
763 old_fe_eval->first.evaluate(old_fe_eval->second);
781template <
unsigned int dim,
unsigned int degree,
typename number>
782template <TensorRank Rank>
797 for (
unsigned int age = 0; age <
fe_eval_old.size(); ++age)
801 old_fe_eval->first.reinit(cell);
802 old_fe_eval->first.read_dof_values_plain(
804 old_fe_eval->first.evaluate(old_fe_eval->second);
826template <
unsigned int dim,
unsigned int degree,
typename number>
827template <TensorRank Rank>
837template <
unsigned int dim,
unsigned int degree,
typename number>
838template <TensorRank Rank>
848template <
unsigned int dim,
unsigned int degree,
typename number>
849template <TensorRank Rank>
861template <
unsigned int dim,
unsigned int degree,
typename number>
862template <TensorRank Rank>
874template <
unsigned int dim,
unsigned int degree,
typename number>
880 fe_eval.reinit(cell);
884 fe_eval.reinit(cell);
889template <
unsigned int dim,
unsigned int degree,
typename number>
896 fe_eval.eval(src_solutions, plain);
900 fe_eval.eval(src_solutions, plain);
905template <
unsigned int dim,
unsigned int degree,
typename number>
914 fe_eval.reinit_and_eval(cell, src_solutions, plain);
918 fe_eval.reinit_and_eval(cell, src_solutions, plain);
924template <
unsigned int dim,
unsigned int degree,
typename number>
930 fe_eval.eval_without_read();
934 fe_eval.eval_without_read();
938template <
unsigned int dim,
unsigned int degree,
typename number>
957template <
unsigned int dim,
unsigned int degree,
typename number>
975template <
unsigned int dim,
unsigned int degree,
typename number>
991template <
unsigned int dim,
unsigned int degree,
typename number>
1010template <
unsigned int dim,
unsigned int degree,
typename number>
1011inline DEAL_II_ALWAYS_INLINE
void
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) \
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); \
1034 catch (const ExcDepNotInitialized &e) \
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" \
1040 << "Ensure that each dependency is requested in the solve block.\n" \
1044 catch (const dealii::internal::ExcAccessToUninitializedField &e) \
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" \
1055template <
unsigned int dim,
unsigned int degree,
typename number>
1056template <TensorRank Rank, DependencyType type>
1057inline DEAL_II_ALWAYS_INLINE
1064template <
unsigned int dim,
unsigned int degree,
typename number>
1065template <TensorRank Rank>
1066inline DEAL_II_ALWAYS_INLINE
1074template <
unsigned int dim,
unsigned int degree,
typename number>
1075template <TensorRank Rank, DependencyType type>
1076inline DEAL_II_ALWAYS_INLINE
1083template <
unsigned int dim,
unsigned int degree,
typename number>
1084template <TensorRank Rank>
1085inline DEAL_II_ALWAYS_INLINE
1093template <
unsigned int dim,
unsigned int degree,
typename number>
1094template <TensorRank Rank, DependencyType type>
1095inline DEAL_II_ALWAYS_INLINE
1102template <
unsigned int dim,
unsigned int degree,
typename number>
1103template <TensorRank Rank>
1104inline DEAL_II_ALWAYS_INLINE
1112template <
unsigned int dim,
unsigned int degree,
typename number>
1113template <TensorRank Rank, DependencyType type>
1114inline DEAL_II_ALWAYS_INLINE
1122template <
unsigned int dim,
unsigned int degree,
typename number>
1123template <TensorRank Rank>
1124inline DEAL_II_ALWAYS_INLINE
1132template <
unsigned int dim,
unsigned int degree,
typename number>
1133template <TensorRank Rank, DependencyType type>
1134inline DEAL_II_ALWAYS_INLINE
1141template <
unsigned int dim,
unsigned int degree,
typename number>
1142template <TensorRank Rank>
1143inline DEAL_II_ALWAYS_INLINE
1151template <
unsigned int dim,
unsigned int degree,
typename number>
1152template <TensorRank Rank, DependencyType type>
1156 static_assert(Rank == 1,
"Divergences are only available for vector fields");
1160template <
unsigned int dim,
unsigned int degree,
typename number>
1161template <TensorRank Rank>
1166 static_assert(Rank == 1,
"Divergences are only available for vector fields");
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>
1177 static_assert(Rank == 1,
"Symmetric gradients are only available for vector fields");
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>
1188 static_assert(Rank == 1,
"Symmetric gradients are only available for vector fields");
1192template <
unsigned int dim,
unsigned int degree,
typename number>
1193template <TensorRank Rank, DependencyType type>
1194inline DEAL_II_ALWAYS_INLINE
1196 (dim == 2 ? 1 : dim),
1200 static_assert(Rank == 1,
"Curl is only available for vector fields");
1204 "Curl is only available for vector fields with dimension greater than 1."));
1208template <
unsigned int dim,
unsigned int degree,
typename number>
1209template <TensorRank Rank>
1210inline DEAL_II_ALWAYS_INLINE
1212 (dim == 2 ? 1 : dim),
1217 static_assert(Rank == 1,
"Curl is only available for vector fields");
1221 "Curl is only available for vector fields with dimension greater than 1."));
1225template <
unsigned int dim,
unsigned int degree,
typename number>
1226inline DEAL_II_ALWAYS_INLINE
1227 dealii::Point<dim, typename FieldContainer<dim, degree, number>::ScalarValue>
1233template <
unsigned int dim,
unsigned int degree,
typename number>
1241template <
unsigned int dim,
unsigned int degree,
typename number>
1242inline DEAL_II_ALWAYS_INLINE
unsigned int
1248template <
unsigned int dim,
unsigned int degree,
typename number>
1249template <
typename ValType>
1250inline DEAL_II_ALWAYS_INLINE
void
1256 dealii::ExcMessage(
"Error: Field index " +
1257 std::to_string(field_index) +
1258 " is not associated with any field."));
1259 auto &relevant_feeval_vector =
1263 relevant_feeval_vector.template get<DependencyType::DST>().submit_value(val,
1265 relevant_feeval_vector.integration_flags |= dealii::EvaluationFlags::values;
1269 std::cerr <<
"Error when trying to submit value for field with index "
1271 <<
": Error: Submission for field not part of this solve block.\n"
1277template <
unsigned int dim,
unsigned int degree,
typename number>
1278template <
typename GradType>
1279inline DEAL_II_ALWAYS_INLINE
void
1281 const GradType &val)
1285 dealii::ExcMessage(
"Error: Field index " +
1286 std::to_string(field_index) +
1287 " is not associated with any field."));
1288 auto &relevant_feeval_vector =
1292 relevant_feeval_vector.template get<DependencyType::DST>().submit_gradient(val,
1294 relevant_feeval_vector.integration_flags |= dealii::EvaluationFlags::gradients;
1298 std::cerr <<
"Error when trying to submit gradient for field with index "
1300 <<
": Error: Submission for field not part of this solve block.\n"
1306template <
unsigned int dim,
unsigned int degree,
typename number>
1307template <
typename ValType>
1308inline DEAL_II_ALWAYS_INLINE
void
1311 unsigned int dof_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 =
1323 relevant_feeval_vector.template get<DependencyType::DST>().get_dof_value(
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"
1335template <
unsigned int dim,
unsigned int degree,
typename number>
1336template <
typename ValType>
1337inline DEAL_II_ALWAYS_INLINE
void
1340 unsigned int dof_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 =
1351 relevant_feeval_vector.template get<DependencyType::DST>()
1352 .submit_dof_value(val, dof_index);
1356 std::cerr <<
"Error when trying to submit dof value for field with index "
1358 <<
": Error: Submission for field not part of this solve block.\n"
1364template <
unsigned int dim,
unsigned int degree,
typename number>
1365template <TensorRank Rank>
1366inline DEAL_II_ALWAYS_INLINE std::vector<
1380template <
unsigned int dim,
unsigned int degree,
typename number>
1381template <TensorRank Rank>
1382inline DEAL_II_ALWAYS_INLINE
const std::vector<
1396template <
unsigned int dim,
unsigned int degree,
typename number>
1397inline DEAL_II_ALWAYS_INLINE
void
1404template <
unsigned int dim,
unsigned int degree,
typename number>
1405inline DEAL_II_ALWAYS_INLINE
void
1413template <
unsigned int dim,
unsigned int degree,
typename number>
1414inline DEAL_II_ALWAYS_INLINE
void
1418 [[maybe_unused]] dealii::EvaluationFlags::EvaluationFlags flag)
const
1423template <
unsigned int dim,
unsigned int degree,
typename number>
1424inline DEAL_II_ALWAYS_INLINE
void
1432#undef AssertAccessible
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