[P]arallel [Hi]gh-order [Li]brary for [P]DEs
Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
|
Class to store local solution coefficients and provide evaluation functions. More...
#include <local_solution.hpp>
Public Member Functions | |
LocalSolution (const dealii::FESystem< dim, dim > &finite_element) | |
Constructor. | |
std::vector< std::array< real, n_components > > | evaluate_values (const std::vector< dealii::Point< dim >> &unit_points) const |
Obtain values at unit points. | |
std::vector< std::array< dealii::Tensor< 1, dim, real >, n_components > > | evaluate_reference_gradients (const std::vector< dealii::Point< dim >> &unit_points) const |
Public Attributes | |
std::vector< real > | coefficients |
Solution coefficients in the finite element basis. | |
const dealii::FESystem< dim, dim > & | finite_element |
Reference to the finite element system used to represent the solution. | |
const unsigned int | n_dofs |
Number of degrees of freedom. | |
Class to store local solution coefficients and provide evaluation functions.
This class is used to store the solution coefficients in the finite element basis and provide functions to evaluate the solution value and gradients at arbitrary points. It can be used for both state and metric solutions since they are both represented by a finite element discretization. The template parameters are:
Definition at line 20 of file local_solution.hpp.
std::vector< std::array< dealii::Tensor< 1, dim, real >, n_components > > PHiLiP::LocalSolution< real, dim, n_components >::evaluate_reference_gradients | ( | const std::vector< dealii::Point< dim >> & | unit_points | ) | const |
Obtain reference gradients at unit points. Note that the gradients are not physical since they do not include metric terms.
Definition at line 37 of file local_solution.cpp.