![]() |
PRISMS-PF Manual v3.0-pre
|
This class permits the access of a subset of indexed fields and gives an error if any non-allowed fields are requested. More...
#include <variable_container.h>
Public Types | |
using | VectorType = dealii::LinearAlgebra::distributed::Vector< number > |
using | value_type = number |
using | size_type = dealii::VectorizedArray< number > |
Public Member Functions | |
variableContainer (const dealii::MatrixFree< dim, number > &data, const std::map< unsigned int, variableAttributes > &_subset_attributes, const std::unordered_map< std::pair< unsigned int, dependencyType >, unsigned int, pairHash > &_global_to_local_solution, const solveType &_solve_type) | |
Constructor. | |
size_type | get_scalar_value (unsigned int global_variable_index, dependencyType dependency_type=dependencyType::NORMAL) const |
Return the value of the specified scalar field. | |
dealii::Tensor< 1, dim, size_type > | get_scalar_gradient (unsigned int global_variable_index, dependencyType dependency_type=dependencyType::NORMAL) const |
Return the gradient of the specified scalar field. | |
dealii::Tensor< 2, dim, size_type > | get_scalar_hessian (unsigned int global_variable_index, dependencyType dependency_type=dependencyType::NORMAL) const |
Return the hessian of the specified scalar field. | |
dealii::Tensor< 1, dim, size_type > | get_scalar_hessian_diagonal (unsigned int global_variable_index, dependencyType dependency_type=dependencyType::NORMAL) const |
Return the diagonal of the hessian of the specified scalar field. | |
size_type | get_scalar_laplacian (unsigned int global_variable_index, dependencyType dependency_type=dependencyType::NORMAL) const |
Return the laplacian of the specified scalar field. | |
dealii::Tensor< 1, dim, size_type > | get_vector_value (unsigned int global_variable_index, dependencyType dependency_type=dependencyType::NORMAL) const |
Return the value of the specified vector field. | |
dealii::Tensor< 2, dim, size_type > | get_vector_gradient (unsigned int global_variable_index, dependencyType dependency_type=dependencyType::NORMAL) const |
Return the gradient of the specified vector field. | |
dealii::Tensor< 3, dim, size_type > | get_vector_hessian (unsigned int global_variable_index, dependencyType dependency_type=dependencyType::NORMAL) const |
Return the hessian of the specified vector field. | |
dealii::Tensor< 2, dim, size_type > | get_vector_hessian_diagonal (unsigned int global_variable_index, dependencyType dependency_type=dependencyType::NORMAL) const |
Return the diagonal of the hessian of the specified vector field. | |
dealii::Tensor< 1, dim, size_type > | get_vector_laplacian (unsigned int global_variable_index, dependencyType dependency_type=dependencyType::NORMAL) const |
Return the laplacian of the specified vector field. | |
size_type | get_vector_divergence (unsigned int global_variable_index, dependencyType dependency_type=dependencyType::NORMAL) const |
Return the divergence of the specified vector field. | |
dealii::Tensor< 2, dim, size_type > | get_vector_symmetric_gradient (unsigned int global_variable_index, dependencyType dependency_type=dependencyType::NORMAL) const |
Return the symmetric gradient of the specified vector field. | |
dealii::Tensor< 1,(dim==2 ? 1 :dim), size_type > | get_vector_curl (unsigned int global_variable_index, dependencyType dependency_type=dependencyType::NORMAL) const |
Return the curl of the specified vector field. Note that this is dealii::VectorizedArray<number> type for 2D and dealii::Tensor<1, dim, dealii::VectorizedArray<number>> type for 3D. | |
void | set_scalar_value_term (const unsigned int &global_variable_index, const size_type &val, const dependencyType &dependency_type=dependencyType::NORMAL) |
Set the residual value of the specified scalar field. | |
void | set_scalar_gradient_term (const unsigned int &global_variable_index, const dealii::Tensor< 1, dim, size_type > &grad, const dependencyType &dependency_type=dependencyType::NORMAL) |
Set the residual gradient of the specified scalar field. | |
void | set_vector_value_term (const unsigned int &global_variable_index, const dealii::Tensor< 1, dim, size_type > &val, const dependencyType &dependency_type=dependencyType::NORMAL) |
Set the residual value of the specified vector field. | |
void | set_vector_gradient_term (const unsigned int &global_variable_index, const dealii::Tensor< 2, dim, size_type > &grad, const dependencyType &dependency_type=dependencyType::NORMAL) |
Set the residual gradient of the specified vector field. | |
void | eval_local_operator (const std::function< void(variableContainer &, const dealii::Point< dim, size_type > &)> &func, std::vector< VectorType * > &dst, const std::vector< VectorType * > &src, const std::pair< unsigned int, unsigned int > &cell_range) |
Apply some operator function for a given cell range and source vector to some destination vector. | |
void | eval_local_operator (const std::function< void(variableContainer &, const dealii::Point< dim, size_type > &)> &func, VectorType &dst, const std::vector< VectorType * > &src, const std::pair< unsigned int, unsigned int > &cell_range) |
Apply some operator function for a given cell range and source vector to some destination vector. | |
void | eval_local_operator (const std::function< void(variableContainer &, const dealii::Point< dim, size_type > &)> &func, VectorType &dst, const VectorType &src, const std::vector< VectorType * > &src_subset, const std::pair< unsigned int, unsigned int > &cell_range) |
Apply some operator function for a given cell range and source vector to some destination vector. | |
void | eval_local_diagonal (const std::function< void(variableContainer &, const dealii::Point< dim, size_type > &)> &func, VectorType &dst, const std::vector< VectorType * > &src_subset, const std::pair< unsigned int, unsigned int > &cell_range) |
TODO (landinjm): Add comments. | |
This class permits the access of a subset of indexed fields and gives an error if any non-allowed fields are requested.
dim | The number of dimensions in the problem. |
degree | The polynomial degree of the shape functions. |
number | Datatype to use for dealii::VectorizedArray<number> . Either double or float. |