[P]arallel [Hi]gh-order [Li]brary for [P]DEs  Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
exact_solution.h
1 #ifndef __EXACT_SOLUTION_H__
2 #define __EXACT_SOLUTION_H__
3 
4 // for the exact_solution function
5 #include <deal.II/lac/vector.h>
6 #include <deal.II/base/function.h>
7 #include "parameters/all_parameters.h"
8 
9 namespace PHiLiP {
10 
12 template <int dim, int nstate, typename real>
13 class ExactSolutionFunction : public dealii::Function<dim,real>
14 {
15 protected:
16  using dealii::Function<dim,real>::value;
17 
18 public:
21 
23  virtual real value (const dealii::Point<dim,real> &point, const unsigned int istate = 0) const = 0;
24 };
25 
27 template <int dim, int nstate, typename real>
29  : public ExactSolutionFunction<dim,nstate,real>
30 {
31 protected:
32  using dealii::Function<dim,real>::value;
33 
34 public:
36  explicit ExactSolutionFunction_Zero (double time_compare);
37 
39  const double t;
40 
42  real value (const dealii::Point<dim,real> &point, const unsigned int istate = 0) const override;
43 };
44 
46 template <int dim, int nstate, typename real>
48  : public ExactSolutionFunction<dim,nstate,real>
49 {
50 protected:
51  using dealii::Function<dim,real>::value;
52 
53 public:
55 
56  explicit ExactSolutionFunction_1DSine (double time_compare);
57 
59  const double t;
60 
62  real value (const dealii::Point<dim,real> &point, const unsigned int istate = 0) const override;
63 };
64 
66 template <int dim, int nstate, typename real>
68  : public ExactSolutionFunction<dim,nstate,real>
69 {
70 protected:
71  using dealii::Function<dim,real>::value;
72 
73 public:
75 
76  explicit ExactSolutionFunction_IsentropicVortex (double time_compare);
77 
79  const double t;
80 
82  real value (const dealii::Point<dim,real> &point, const unsigned int istate = 0) const override;
83 };
84 
86 template <int dim, int nstate, typename real>
88 {
89 protected:
92 
93 public:
95  static std::shared_ptr<ExactSolutionFunction<dim,nstate,real>>
96  create_ExactSolutionFunction(const Parameters::FlowSolverParam& flow_solver_parameters, const double time_compare);
97 };
98 
99 } // PHiLiP namespace
100 #endif
Exact Solution Function: 1D Sine Function; used for temporal convergence.
FlowCaseType
Selects the flow case to be simulated.
const double t
Time at which to compute the exact solution.
Files for the baseline physics.
Definition: ADTypes.hpp:10
Exact Solution Function: Zero Function; used as a placeholder when there is no exact solution...
ExactSolutionFunction()
< dealii::Function we are templating on
Parameters related to the flow solver.
Exact solution function factory.
virtual real value(const dealii::Point< dim, real > &point, const unsigned int istate=0) const =0
Value of the exact solution at a point.
const double t
Time at which to compute the exact solution.
Exact solution function used for a particular flow setup/case.
const double t
Time at which to compute the exact solution.
Exact Solution Function: Isentropic vortex.