PRISMS-PF Manual
Loading...
Searching...
No Matches
FieldContainer< dim, degree, number > Class Template Reference

This class permits the access of a subset of indexed fields and gives an error if any non-allowed fields are requested. More...

#include <field_container.h>

Collaboration diagram for FieldContainer< dim, degree, number >:
[legend]

Classes

struct  FEEValuationDeps
 Struct to hold the relevant dealii::FEEvaluation for a given solution block index. More...
 

Public Types

using ScalarValue = dealii::VectorizedArray<number>
 Typedef for the basic value that the user manipulates.
 
using VectorValue = dealii::Tensor<1, dim, ScalarValue>
 Typedef for the basic vector that the user manipulates.
 
template<TensorRank Rank>
using Value
 
template<TensorRank Rank>
using Gradient = dealii::Tensor<int(Rank) + 1, dim, ScalarValue>
 
template<TensorRank Rank>
using Hessian = dealii::Tensor<int(Rank) + 2, dim, ScalarValue>
 
template<TensorRank Rank>
using FEEval
 

Public Member Functions

 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.
 
void reinit (unsigned int cell)
 Initialize based on cell for all dependencies.
 
void eval (const BlockVector< number > *src_solutions, bool plain)
 Read solution vector, and evaluate based on dependency flags for all dependencies.
 
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 dependencies.
 
void integrate ()
 Integrate the residuals.
 
void distribute (BlockVector< number > *dst_solutions)
 Distribute the integrated residuals.
 
void integrate_and_distribute (BlockVector< number > *dst_solutions)
 Integrate the residuals and distribute from local to global.
 
void set_q_point (unsigned int q)
 Set the current quadrature point.
 
template<TensorRank Rank, DependencyType type>
Value< Rank > get_value (Types::Index field_index) const
 Return the value of the specified field.
 
template<TensorRank Rank>
Value< Rank > get_value (Types::Index field_index, DependencyType type) const
 Return the value of the specified field.
 
template<TensorRank Rank, DependencyType type>
Gradient< Rank > get_gradient (Types::Index field_index) const
 Return the gradient of the specified field.
 
template<TensorRank Rank>
Gradient< Rank > get_gradient (Types::Index field_index, DependencyType type) const
 Return the gradient of the specified field.
 
template<TensorRank Rank, DependencyType type>
Hessian< Rank > get_hessian (Types::Index field_index) const
 Return the hessian of the specified field.
 
template<TensorRank Rank>
Hessian< Rank > get_hessian (Types::Index field_index, DependencyType type) const
 Return the hessian of the specified field.
 
template<TensorRank Rank, DependencyType type>
Gradient< Rank > get_hessian_diagonal (Types::Index field_index) const
 Return the diagonal of the hessian of the specified field.
 
template<TensorRank Rank>
Gradient< Rank > get_hessian_diagonal (Types::Index field_index, DependencyType type) const
 Return the diagonal of the hessian of the specified field.
 
template<TensorRank Rank, DependencyType type>
Value< Rank > get_laplacian (Types::Index field_index) const
 Return the laplacian of the specified field.
 
template<TensorRank Rank>
Value< Rank > get_laplacian (Types::Index field_index, DependencyType type) const
 Return the laplacian of the specified field.
 
template<TensorRank Rank, DependencyType type>
ScalarValue get_divergence (Types::Index field_index) const
 Return the divergence of the specified field.
 
template<TensorRank Rank>
ScalarValue get_divergence (Types::Index field_index, DependencyType type) const
 Return the divergence of the specified field.
 
template<TensorRank Rank, DependencyType type>
dealii::SymmetricTensor< 2, dim, ScalarValueget_symmetric_gradient (Types::Index field_index) const
 Return the symmetric gradient of the specified field.
 
template<TensorRank Rank>
dealii::SymmetricTensor< 2, dim, ScalarValueget_symmetric_gradient (Types::Index field_index, DependencyType type) const
 Return the symmetric gradient of the specified field.
 
template<TensorRank Rank, DependencyType type>
dealii::Tensor< 1,(dim==2 ? 1 :dim), ScalarValueget_curl (Types::Index field_index) const
 Return the curl of the specified field.
 
template<TensorRank Rank>
dealii::Tensor< 1,(dim==2 ? 1 :dim), ScalarValueget_curl (Types::Index field_index, DependencyType type) const
 Return the curl of the specified field.
 
dealii::Point< dim, ScalarValueget_q_point_location () const
 Return the quadrature point location.
 
ScalarValue get_element_volume () const
 Return the quadrature point location.
 
unsigned int get_n_q_points () const
 Return the number of quadrature points.
 
template<typename ValType>
void set_value_term (Types::Index field_index, const ValType &val)
 Set the residual value of the specified scalar/vector field.
 
template<typename GradType>
void set_gradient_term (Types::Index field_index, const GradType &val)
 Set the residual gradient of the specified scalar/vector field.
 
template<TensorRank Rank, DependencyType type>
DEAL_II_ALWAYS_INLINE FieldContainer< dim, degree, number >::template Value< Rank > get_value (Types::Index field_index) const
 
template<TensorRank Rank>
DEAL_II_ALWAYS_INLINE FieldContainer< dim, degree, number >::template Value< Rank > get_value (Types::Index field_index, DependencyType type) const
 
template<TensorRank Rank, DependencyType type>
DEAL_II_ALWAYS_INLINE FieldContainer< dim, degree, number >::template Gradient< Rank > get_gradient (Types::Index field_index) const
 
template<TensorRank Rank>
DEAL_II_ALWAYS_INLINE FieldContainer< dim, degree, number >::template Gradient< Rank > get_gradient (Types::Index field_index, DependencyType type) const
 
template<TensorRank Rank, DependencyType type>
DEAL_II_ALWAYS_INLINE FieldContainer< dim, degree, number >::template Hessian< Rank > get_hessian (Types::Index field_index) const
 
template<TensorRank Rank>
DEAL_II_ALWAYS_INLINE FieldContainer< dim, degree, number >::template Hessian< Rank > get_hessian (Types::Index field_index, DependencyType type) const
 
template<TensorRank Rank, DependencyType type>
DEAL_II_ALWAYS_INLINE FieldContainer< dim, degree, number >::template Gradient< Rank > get_hessian_diagonal (Types::Index field_index) const
 
template<TensorRank Rank>
DEAL_II_ALWAYS_INLINE FieldContainer< dim, degree, number >::template Gradient< Rank > get_hessian_diagonal (Types::Index field_index, DependencyType type) const
 
template<TensorRank Rank, DependencyType type>
DEAL_II_ALWAYS_INLINE FieldContainer< dim, degree, number >::template Value< Rank > get_laplacian (Types::Index field_index) const
 
template<TensorRank Rank>
DEAL_II_ALWAYS_INLINE FieldContainer< dim, degree, number >::template Value< Rank > get_laplacian (Types::Index field_index, DependencyType type) const
 
template<TensorRank Rank, DependencyType type>
DEAL_II_ALWAYS_INLINE FieldContainer< dim, degree, number >::ScalarValue get_divergence (Types::Index field_index) const
 
template<TensorRank Rank>
DEAL_II_ALWAYS_INLINE FieldContainer< dim, degree, number >::ScalarValue get_divergence (Types::Index field_index, DependencyType type) const
 
template<TensorRank Rank, DependencyType type>
DEAL_II_ALWAYS_INLINE dealii::SymmetricTensor< 2, dim, typename FieldContainer< dim, degree, number >::ScalarValueget_symmetric_gradient (Types::Index field_index) const
 
template<TensorRank Rank>
DEAL_II_ALWAYS_INLINE dealii::SymmetricTensor< 2, dim, typename FieldContainer< dim, degree, number >::ScalarValueget_symmetric_gradient (Types::Index field_index, DependencyType type) const
 
template<TensorRank Rank, DependencyType type>
DEAL_II_ALWAYS_INLINE dealii::Tensor< 1,(dim==2 ? 1 :dim), typename FieldContainer< dim, degree, number >::ScalarValueget_curl (Types::Index field_index) const
 
template<TensorRank Rank>
DEAL_II_ALWAYS_INLINE dealii::Tensor< 1,(dim==2 ? 1 :dim), typename FieldContainer< dim, degree, number >::ScalarValueget_curl (Types::Index field_index, DependencyType type) const
 
template<typename ValType>
DEAL_II_ALWAYS_INLINE void set_value_term (Types::Index field_index, const ValType &val)
 
template<typename GradType>
DEAL_II_ALWAYS_INLINE void set_gradient_term (Types::Index field_index, const GradType &val)
 
template<TensorRank Rank>
DEAL_II_ALWAYS_INLINE std::vector< typename FieldContainer< dim, degree, number >::template FEEValuationDeps< Rank > > & get_relevant_feeval_vector ()
 
template<TensorRank Rank>
DEAL_II_ALWAYS_INLINE const std::vector< typename FieldContainer< dim, degree, number >::template FEEValuationDeps< Rank > > & get_relevant_feeval_vector () const
 

Static Public Attributes

template<typename ValType>
static constexpr TensorRank RankFromVal
 Return the tensor rank from the specified template value.
 
template<typename GradType>
static constexpr TensorRank RankFromGrad
 Return the tensor rank from the specified template gradient.
 

Private Member Functions

template<TensorRank Rank>
std::vector< FEEValuationDeps< Rank > > & get_relevant_feeval_vector ()
 
template<TensorRank Rank>
const std::vector< FEEValuationDeps< Rank > > & get_relevant_feeval_vector () const
 
void feevaluation_size_valid (Types::Index field_index) const
 Check whether the entry for the FEEvaluation is within the bounds of the vector.
 
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.
 
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.
 
void submission_valid (Types::Index field_index, DependencyType dependency_type) const
 Check that a value is valid for submission.
 

Private Attributes

const std::vector< FieldAttributes > * field_attributes_ptr
 Pointer to the vector of field attributes.
 
const SolutionIndexer< dim, number > * solution_indexer
 Pointer to the solution indexer.
 
std::vector< FEEValuationDeps< TensorRank::Scalar > > feeval_deps_scalar
 Collection of FEEvaluation dependencies for scalar fields.
 
std::vector< FEEValuationDeps< TensorRank::Vector > > feeval_deps_vector
 Collection of FEEvaluation dependencies for vector fields.
 
const SolveBlocksolve_block
 Solve block information.
 
FEEval< TensorRank::Scalarshared_feeval_scalar
 FEEvaluation object for generic cell operations.
 
unsigned int relative_level
 Multigrid level.
 
unsigned int q_point = 0
 The quadrature point.
 

Detailed Description

template<unsigned int dim, unsigned int degree, typename number>
class FieldContainer< dim, degree, number >

This class permits the access of a subset of indexed fields and gives an error if any non-allowed fields are requested.

Template Parameters
dimThe number of dimensions in the problem.
degreeThe polynomial degree of the shape functions.
numberDatatype to use for dealii::VectorizedArray<number>. Either double or float.

Importantly, this class is mostly a wrapper for the dealii::FEEvaluation class, which allows for the evaluation of functions at quadrature points and cell integrations (basically the backbone of the FEM solver). As such, it is one of the most important classes when it comes to performance because these functions are called millions of times.

Member Typedef Documentation

◆ FEEval

template<unsigned int dim, unsigned int degree, typename number>
template<TensorRank Rank>
using FieldContainer< dim, degree, number >::FEEval
Initial value:
dealii::FEEvaluation<dim,
degree,
degree + 1,
dealii::Tensor<int(Rank), dim>::n_independent_components,
number,
dealii::VectorizedArray< number > ScalarValue
Typedef for the basic value that the user manipulates.
Definition field_container.h:52

◆ Gradient

template<unsigned int dim, unsigned int degree, typename number>
template<TensorRank Rank>
using FieldContainer< dim, degree, number >::Gradient = dealii::Tensor<int(Rank) + 1, dim, ScalarValue>

◆ Hessian

template<unsigned int dim, unsigned int degree, typename number>
template<TensorRank Rank>
using FieldContainer< dim, degree, number >::Hessian = dealii::Tensor<int(Rank) + 2, dim, ScalarValue>

◆ ScalarValue

template<unsigned int dim, unsigned int degree, typename number>
using FieldContainer< dim, degree, number >::ScalarValue = dealii::VectorizedArray<number>

Typedef for the basic value that the user manipulates.

◆ Value

template<unsigned int dim, unsigned int degree, typename number>
template<TensorRank Rank>
using FieldContainer< dim, degree, number >::Value
Initial value:
std::conditional_t<Rank == TensorRank::Scalar,
dealii::Tensor<int(Rank), dim, ScalarValue>>
@ Scalar
Definition type_enums.h:54

◆ VectorValue

template<unsigned int dim, unsigned int degree, typename number>
using FieldContainer< dim, degree, number >::VectorValue = dealii::Tensor<1, dim, ScalarValue>

Typedef for the basic vector that the user manipulates.

Constructor & Destructor Documentation

◆ FieldContainer()

template<unsigned int dim, unsigned int degree, typename number>
PRISMS_PF_BEGIN_NAMESPACE FieldContainer< dim, degree, number >::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.

field_attributes is the collection of field attributes. We only use this for the rank of field (Scalar/Vector). solution_indexer is the solution indexer object, which allows us to get access to solution vectors from in and outside the current solve block. relative_level is the multigrid level. dependency_map is the map of field indices and their dependency flags. solve_block is the object that contains attributes about the current solve block. matrix_free is the dealii::MatrixFree object.

Member Function Documentation

◆ access_valid()

template<unsigned int dim, unsigned int degree, typename number>
DEAL_II_ALWAYS_INLINE void FieldContainer< dim, degree, number >::access_valid ( Types::Index field_index,
DependencyType dependency_type,
dealii::EvaluationFlags::EvaluationFlags flag ) const
inlineprivate

Check that a variable value/gradient/hessians was marked as needed and thus properly initialized.

◆ distribute()

template<unsigned int dim, unsigned int degree, typename number>
void FieldContainer< dim, degree, number >::distribute ( BlockVector< number > * dst_solutions)
inline

Distribute the integrated residuals.

◆ eval()

template<unsigned int dim, unsigned int degree, typename number>
void FieldContainer< dim, degree, number >::eval ( const BlockVector< number > * src_solutions,
bool plain )
inline

Read solution vector, and evaluate based on dependency flags for all dependencies.

◆ feevaluation_exists()

template<unsigned int dim, unsigned int degree, typename number>
DEAL_II_ALWAYS_INLINE void FieldContainer< dim, degree, number >::feevaluation_exists ( Types::Index field_index,
Types::Index dependency_index ) const
inlineprivate

Check whether the entry for the FEEvaluation is within the bounds of the vector and not a nullptr.

◆ feevaluation_size_valid()

template<unsigned int dim, unsigned int degree, typename number>
DEAL_II_ALWAYS_INLINE void FieldContainer< dim, degree, number >::feevaluation_size_valid ( Types::Index field_index) const
inlineprivate

Check whether the entry for the FEEvaluation is within the bounds of the vector.

◆ get_curl() [1/4]

template<unsigned int dim, unsigned int degree, typename number>
template<TensorRank Rank, DependencyType type>
dealii::Tensor< 1,(dim==2 ? 1 :dim), ScalarValue > FieldContainer< dim, degree, number >::get_curl ( Types::Index field_index) const
nodiscard

Return the curl of the specified field.

◆ get_curl() [2/4]

template<unsigned int dim, unsigned int degree, typename number>
template<TensorRank Rank, DependencyType type>
DEAL_II_ALWAYS_INLINE dealii::Tensor< 1,(dim==2 ? 1 :dim), typename FieldContainer< dim, degree, number >::ScalarValue > FieldContainer< dim, degree, number >::get_curl ( Types::Index field_index) const
inline

◆ get_curl() [3/4]

template<unsigned int dim, unsigned int degree, typename number>
template<TensorRank Rank>
dealii::Tensor< 1,(dim==2 ? 1 :dim), ScalarValue > FieldContainer< dim, degree, number >::get_curl ( Types::Index field_index,
DependencyType type ) const
nodiscard

Return the curl of the specified field.

◆ get_curl() [4/4]

template<unsigned int dim, unsigned int degree, typename number>
template<TensorRank Rank>
DEAL_II_ALWAYS_INLINE dealii::Tensor< 1,(dim==2 ? 1 :dim), typename FieldContainer< dim, degree, number >::ScalarValue > FieldContainer< dim, degree, number >::get_curl ( Types::Index field_index,
DependencyType type ) const
inline

◆ get_divergence() [1/4]

template<unsigned int dim, unsigned int degree, typename number>
template<TensorRank Rank, DependencyType type>
ScalarValue FieldContainer< dim, degree, number >::get_divergence ( Types::Index field_index) const
nodiscard

Return the divergence of the specified field.

◆ get_divergence() [2/4]

template<unsigned int dim, unsigned int degree, typename number>
template<TensorRank Rank, DependencyType type>
DEAL_II_ALWAYS_INLINE FieldContainer< dim, degree, number >::ScalarValue FieldContainer< dim, degree, number >::get_divergence ( Types::Index field_index) const
inline

◆ get_divergence() [3/4]

template<unsigned int dim, unsigned int degree, typename number>
template<TensorRank Rank>
ScalarValue FieldContainer< dim, degree, number >::get_divergence ( Types::Index field_index,
DependencyType type ) const
nodiscard

Return the divergence of the specified field.

◆ get_divergence() [4/4]

template<unsigned int dim, unsigned int degree, typename number>
template<TensorRank Rank>
DEAL_II_ALWAYS_INLINE FieldContainer< dim, degree, number >::ScalarValue FieldContainer< dim, degree, number >::get_divergence ( Types::Index field_index,
DependencyType type ) const
inline

◆ get_element_volume()

template<unsigned int dim, unsigned int degree, typename number>
DEAL_II_ALWAYS_INLINE FieldContainer< dim, degree, number >::ScalarValue FieldContainer< dim, degree, number >::get_element_volume ( ) const
inlinenodiscard

Return the quadrature point location.

◆ get_gradient() [1/4]

template<unsigned int dim, unsigned int degree, typename number>
template<TensorRank Rank, DependencyType type>
Gradient< Rank > FieldContainer< dim, degree, number >::get_gradient ( Types::Index field_index) const
nodiscard

Return the gradient of the specified field.

◆ get_gradient() [2/4]

template<unsigned int dim, unsigned int degree, typename number>
template<TensorRank Rank, DependencyType type>
DEAL_II_ALWAYS_INLINE FieldContainer< dim, degree, number >::template Gradient< Rank > FieldContainer< dim, degree, number >::get_gradient ( Types::Index field_index) const
inline

◆ get_gradient() [3/4]

template<unsigned int dim, unsigned int degree, typename number>
template<TensorRank Rank>
Gradient< Rank > FieldContainer< dim, degree, number >::get_gradient ( Types::Index field_index,
DependencyType type ) const
nodiscard

Return the gradient of the specified field.

◆ get_gradient() [4/4]

template<unsigned int dim, unsigned int degree, typename number>
template<TensorRank Rank>
DEAL_II_ALWAYS_INLINE FieldContainer< dim, degree, number >::template Gradient< Rank > FieldContainer< dim, degree, number >::get_gradient ( Types::Index field_index,
DependencyType type ) const
inline

◆ get_hessian() [1/4]

template<unsigned int dim, unsigned int degree, typename number>
template<TensorRank Rank, DependencyType type>
Hessian< Rank > FieldContainer< dim, degree, number >::get_hessian ( Types::Index field_index) const
nodiscard

Return the hessian of the specified field.

◆ get_hessian() [2/4]

template<unsigned int dim, unsigned int degree, typename number>
template<TensorRank Rank, DependencyType type>
DEAL_II_ALWAYS_INLINE FieldContainer< dim, degree, number >::template Hessian< Rank > FieldContainer< dim, degree, number >::get_hessian ( Types::Index field_index) const
inline

◆ get_hessian() [3/4]

template<unsigned int dim, unsigned int degree, typename number>
template<TensorRank Rank>
Hessian< Rank > FieldContainer< dim, degree, number >::get_hessian ( Types::Index field_index,
DependencyType type ) const
nodiscard

Return the hessian of the specified field.

◆ get_hessian() [4/4]

template<unsigned int dim, unsigned int degree, typename number>
template<TensorRank Rank>
DEAL_II_ALWAYS_INLINE FieldContainer< dim, degree, number >::template Hessian< Rank > FieldContainer< dim, degree, number >::get_hessian ( Types::Index field_index,
DependencyType type ) const
inline

◆ get_hessian_diagonal() [1/4]

template<unsigned int dim, unsigned int degree, typename number>
template<TensorRank Rank, DependencyType type>
Gradient< Rank > FieldContainer< dim, degree, number >::get_hessian_diagonal ( Types::Index field_index) const
nodiscard

Return the diagonal of the hessian of the specified field.

◆ get_hessian_diagonal() [2/4]

template<unsigned int dim, unsigned int degree, typename number>
template<TensorRank Rank, DependencyType type>
DEAL_II_ALWAYS_INLINE FieldContainer< dim, degree, number >::template Gradient< Rank > FieldContainer< dim, degree, number >::get_hessian_diagonal ( Types::Index field_index) const
inline

◆ get_hessian_diagonal() [3/4]

template<unsigned int dim, unsigned int degree, typename number>
template<TensorRank Rank>
Gradient< Rank > FieldContainer< dim, degree, number >::get_hessian_diagonal ( Types::Index field_index,
DependencyType type ) const
nodiscard

Return the diagonal of the hessian of the specified field.

◆ get_hessian_diagonal() [4/4]

template<unsigned int dim, unsigned int degree, typename number>
template<TensorRank Rank>
DEAL_II_ALWAYS_INLINE FieldContainer< dim, degree, number >::template Gradient< Rank > FieldContainer< dim, degree, number >::get_hessian_diagonal ( Types::Index field_index,
DependencyType type ) const
inline

◆ get_laplacian() [1/4]

template<unsigned int dim, unsigned int degree, typename number>
template<TensorRank Rank, DependencyType type>
Value< Rank > FieldContainer< dim, degree, number >::get_laplacian ( Types::Index field_index) const
nodiscard

Return the laplacian of the specified field.

◆ get_laplacian() [2/4]

template<unsigned int dim, unsigned int degree, typename number>
template<TensorRank Rank, DependencyType type>
DEAL_II_ALWAYS_INLINE FieldContainer< dim, degree, number >::template Value< Rank > FieldContainer< dim, degree, number >::get_laplacian ( Types::Index field_index) const
inline

◆ get_laplacian() [3/4]

template<unsigned int dim, unsigned int degree, typename number>
template<TensorRank Rank>
Value< Rank > FieldContainer< dim, degree, number >::get_laplacian ( Types::Index field_index,
DependencyType type ) const
nodiscard

Return the laplacian of the specified field.

◆ get_laplacian() [4/4]

template<unsigned int dim, unsigned int degree, typename number>
template<TensorRank Rank>
DEAL_II_ALWAYS_INLINE FieldContainer< dim, degree, number >::template Value< Rank > FieldContainer< dim, degree, number >::get_laplacian ( Types::Index field_index,
DependencyType type ) const
inline

◆ get_n_q_points()

template<unsigned int dim, unsigned int degree, typename number>
DEAL_II_ALWAYS_INLINE unsigned int FieldContainer< dim, degree, number >::get_n_q_points ( ) const
inlinenodiscard

Return the number of quadrature points.

◆ get_q_point_location()

template<unsigned int dim, unsigned int degree, typename number>
DEAL_II_ALWAYS_INLINE dealii::Point< dim, typename FieldContainer< dim, degree, number >::ScalarValue > FieldContainer< dim, degree, number >::get_q_point_location ( ) const
inlinenodiscard

Return the quadrature point location.

◆ get_relevant_feeval_vector() [1/4]

template<unsigned int dim, unsigned int degree, typename number>
template<TensorRank Rank>
std::vector< FEEValuationDeps< Rank > > & FieldContainer< dim, degree, number >::get_relevant_feeval_vector ( )
private

◆ get_relevant_feeval_vector() [2/4]

template<unsigned int dim, unsigned int degree, typename number>
template<TensorRank Rank>
DEAL_II_ALWAYS_INLINE std::vector< typename FieldContainer< dim, degree, number >::template FEEValuationDeps< Rank > > & FieldContainer< dim, degree, number >::get_relevant_feeval_vector ( )
inline

◆ get_relevant_feeval_vector() [3/4]

template<unsigned int dim, unsigned int degree, typename number>
template<TensorRank Rank>
const std::vector< FEEValuationDeps< Rank > > & FieldContainer< dim, degree, number >::get_relevant_feeval_vector ( ) const
private

◆ get_relevant_feeval_vector() [4/4]

template<unsigned int dim, unsigned int degree, typename number>
template<TensorRank Rank>
DEAL_II_ALWAYS_INLINE const std::vector< typename FieldContainer< dim, degree, number >::template FEEValuationDeps< Rank > > & FieldContainer< dim, degree, number >::get_relevant_feeval_vector ( ) const
inline

◆ get_symmetric_gradient() [1/4]

template<unsigned int dim, unsigned int degree, typename number>
template<TensorRank Rank, DependencyType type>
dealii::SymmetricTensor< 2, dim, ScalarValue > FieldContainer< dim, degree, number >::get_symmetric_gradient ( Types::Index field_index) const
nodiscard

Return the symmetric gradient of the specified field.

◆ get_symmetric_gradient() [2/4]

template<unsigned int dim, unsigned int degree, typename number>
template<TensorRank Rank, DependencyType type>
DEAL_II_ALWAYS_INLINE dealii::SymmetricTensor< 2, dim, typename FieldContainer< dim, degree, number >::ScalarValue > FieldContainer< dim, degree, number >::get_symmetric_gradient ( Types::Index field_index) const
inline

◆ get_symmetric_gradient() [3/4]

template<unsigned int dim, unsigned int degree, typename number>
template<TensorRank Rank>
dealii::SymmetricTensor< 2, dim, ScalarValue > FieldContainer< dim, degree, number >::get_symmetric_gradient ( Types::Index field_index,
DependencyType type ) const
nodiscard

Return the symmetric gradient of the specified field.

◆ get_symmetric_gradient() [4/4]

template<unsigned int dim, unsigned int degree, typename number>
template<TensorRank Rank>
DEAL_II_ALWAYS_INLINE dealii::SymmetricTensor< 2, dim, typename FieldContainer< dim, degree, number >::ScalarValue > FieldContainer< dim, degree, number >::get_symmetric_gradient ( Types::Index field_index,
DependencyType type ) const
inline

◆ get_value() [1/4]

template<unsigned int dim, unsigned int degree, typename number>
template<TensorRank Rank, DependencyType type>
Value< Rank > FieldContainer< dim, degree, number >::get_value ( Types::Index field_index) const
nodiscard

Return the value of the specified field.

◆ get_value() [2/4]

template<unsigned int dim, unsigned int degree, typename number>
template<TensorRank Rank, DependencyType type>
DEAL_II_ALWAYS_INLINE FieldContainer< dim, degree, number >::template Value< Rank > FieldContainer< dim, degree, number >::get_value ( Types::Index field_index) const
inline

◆ get_value() [3/4]

template<unsigned int dim, unsigned int degree, typename number>
template<TensorRank Rank>
Value< Rank > FieldContainer< dim, degree, number >::get_value ( Types::Index field_index,
DependencyType type ) const
nodiscard

Return the value of the specified field.

◆ get_value() [4/4]

template<unsigned int dim, unsigned int degree, typename number>
template<TensorRank Rank>
DEAL_II_ALWAYS_INLINE FieldContainer< dim, degree, number >::template Value< Rank > FieldContainer< dim, degree, number >::get_value ( Types::Index field_index,
DependencyType type ) const
inline

◆ integrate()

template<unsigned int dim, unsigned int degree, typename number>
void FieldContainer< dim, degree, number >::integrate ( )
inline

Integrate the residuals.

◆ integrate_and_distribute()

template<unsigned int dim, unsigned int degree, typename number>
void FieldContainer< dim, degree, number >::integrate_and_distribute ( BlockVector< number > * dst_solutions)
inline

Integrate the residuals and distribute from local to global.

This is more efficient that calling integrate and distribute individually.

◆ reinit()

template<unsigned int dim, unsigned int degree, typename number>
void FieldContainer< dim, degree, number >::reinit ( unsigned int cell)
inline

Initialize based on cell for all dependencies.

◆ reinit_and_eval()

template<unsigned int dim, unsigned int degree, typename number>
void FieldContainer< dim, degree, number >::reinit_and_eval ( unsigned int cell,
const BlockVector< number > * src_solutions,
bool plain )
inline

Initialize based on cell, read solution vector, and evaluate based on dependency flags for all dependencies.

This is more efficient than calling reinit and eval individually.

◆ set_gradient_term() [1/2]

template<unsigned int dim, unsigned int degree, typename number>
template<typename GradType>
void FieldContainer< dim, degree, number >::set_gradient_term ( Types::Index field_index,
const GradType & val )

Set the residual gradient of the specified scalar/vector field.

◆ set_gradient_term() [2/2]

template<unsigned int dim, unsigned int degree, typename number>
template<typename GradType>
DEAL_II_ALWAYS_INLINE void FieldContainer< dim, degree, number >::set_gradient_term ( Types::Index field_index,
const GradType & val )
inline

◆ set_q_point()

template<unsigned int dim, unsigned int degree, typename number>
DEAL_II_ALWAYS_INLINE void FieldContainer< dim, degree, number >::set_q_point ( unsigned int q)
inline

Set the current quadrature point.

◆ set_value_term() [1/2]

template<unsigned int dim, unsigned int degree, typename number>
template<typename ValType>
void FieldContainer< dim, degree, number >::set_value_term ( Types::Index field_index,
const ValType & val )

Set the residual value of the specified scalar/vector field.

◆ set_value_term() [2/2]

template<unsigned int dim, unsigned int degree, typename number>
template<typename ValType>
DEAL_II_ALWAYS_INLINE void FieldContainer< dim, degree, number >::set_value_term ( Types::Index field_index,
const ValType & val )
inline

◆ submission_valid()

template<unsigned int dim, unsigned int degree, typename number>
DEAL_II_ALWAYS_INLINE void FieldContainer< dim, degree, number >::submission_valid ( Types::Index field_index,
DependencyType dependency_type ) const
inlineprivate

Check that a value is valid for submission.

Member Data Documentation

◆ feeval_deps_scalar

template<unsigned int dim, unsigned int degree, typename number>
std::vector<FEEValuationDeps<TensorRank::Scalar> > FieldContainer< dim, degree, number >::feeval_deps_scalar
private

Collection of FEEvaluation dependencies for scalar fields.

◆ feeval_deps_vector

template<unsigned int dim, unsigned int degree, typename number>
std::vector<FEEValuationDeps<TensorRank::Vector> > FieldContainer< dim, degree, number >::feeval_deps_vector
private

Collection of FEEvaluation dependencies for vector fields.

◆ field_attributes_ptr

template<unsigned int dim, unsigned int degree, typename number>
const std::vector<FieldAttributes>* FieldContainer< dim, degree, number >::field_attributes_ptr
private

Pointer to the vector of field attributes.

◆ q_point

template<unsigned int dim, unsigned int degree, typename number>
unsigned int FieldContainer< dim, degree, number >::q_point = 0
private

The quadrature point.

◆ RankFromGrad

template<unsigned int dim, unsigned int degree, typename number>
template<typename GradType>
TensorRank FieldContainer< dim, degree, number >::RankFromGrad
staticconstexpr
Initial value:
= []() constexpr
{
if constexpr (std::is_same_v<GradType, ScalarValue>)
{
}
else
{
return TensorRank(GradType::rank - 1);
}
}()
TensorRank
Tensor rank of the field.
Definition type_enums.h:52

Return the tensor rank from the specified template gradient.

◆ RankFromVal

template<unsigned int dim, unsigned int degree, typename number>
template<typename ValType>
TensorRank FieldContainer< dim, degree, number >::RankFromVal
staticconstexpr
Initial value:
= []() constexpr
{
if constexpr (std::is_same_v<ValType, ScalarValue> ||
std::is_same_v<ValType, typename ScalarValue::value_type>)
{
}
else
{
return TensorRank(ValType::rank);
}
}()

Return the tensor rank from the specified template value.

◆ relative_level

template<unsigned int dim, unsigned int degree, typename number>
unsigned int FieldContainer< dim, degree, number >::relative_level
private

Multigrid level.

◆ shared_feeval_scalar

template<unsigned int dim, unsigned int degree, typename number>
FEEval<TensorRank::Scalar> FieldContainer< dim, degree, number >::shared_feeval_scalar
private

FEEvaluation object for generic cell operations.

This is what we use when we call get_q_point_location, get_element_volume, and get_n_q_points().

Note
This is constructed with the 0th index of the MatrixFree object. If that is not a scalar you'll run into issues.

◆ solution_indexer

template<unsigned int dim, unsigned int degree, typename number>
const SolutionIndexer<dim, number>* FieldContainer< dim, degree, number >::solution_indexer
private

Pointer to the solution indexer.

◆ solve_block

template<unsigned int dim, unsigned int degree, typename number>
const SolveBlock* FieldContainer< dim, degree, number >::solve_block
private

Solve block information.


The documentation for this class was generated from the following files: