[P]arallel [Hi]gh-order [Li]brary for [P]DEs  Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
local_solution.hpp
1 #ifndef PHILIP_LOCAL_SOLUTION_HPP
2 #define PHILIP_LOCAL_SOLUTION_HPP
3 
4 #include <deal.II/fe/fe_system.h>
5 #include <deal.II/base/types.h>
6 
7 
8 namespace PHiLiP {
9 
11 
19 template <typename real, int dim, int n_components>
21  public:
23  std::vector<real> coefficients;
25  const dealii::FESystem<dim, dim> &finite_element;
27  const unsigned int n_dofs;
28 
30  explicit LocalSolution(const dealii::FESystem<dim, dim> &finite_element);
31 
33  std::vector<std::array<real, n_components>> evaluate_values(const std::vector<dealii::Point<dim>> &unit_points) const;
34 
37  std::vector<std::array<dealii::Tensor<1, dim, real>, n_components>> evaluate_reference_gradients(
38  const std::vector<dealii::Point<dim>> &unit_points) const;
39 };
40 
41 } // namespace PHiLiP
42 #endif // PHILIP_LOCAL_SOLUTION_HPP
LocalSolution(const dealii::FESystem< dim, dim > &finite_element)
Constructor.
const unsigned int n_dofs
Number of degrees of freedom.
Class to store local solution coefficients and provide evaluation functions.
Files for the baseline physics.
Definition: ADTypes.hpp:10
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
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.