[P]arallel [Hi]gh-order [Li]brary for [P]DEs  Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
ode_solver_factory.h
1 #ifndef __ODE_SOLVER_FACTORY__
2 #define __ODE_SOLVER_FACTORY__
3 
4 #include "dg/dg_base.hpp"
5 #include "ode_solver_base.h"
6 #include "parameters/all_parameters.h"
7 #include "reduced_order/pod_basis_base.h"
8 #include "runge_kutta_methods/rk_tableau_base.h"
9 #include "runge_kutta_methods/low_storage_rk_tableau_base.h"
10 #include "relaxation_runge_kutta/empty_RRK_base.h"
11 
12 namespace PHiLiP {
13 namespace ODE {
15 
17 #if PHILIP_DIM==1
18 template <int dim, typename real, typename MeshType = dealii::Triangulation<dim>>
19 #else
20 template <int dim, typename real, typename MeshType = dealii::parallel::distributed::Triangulation<dim>>
21 #endif
23 {
24 public:
26  static std::shared_ptr<ODESolverBase<dim,real,MeshType>> create_ODESolver(std::shared_ptr< DGBase<dim, real, MeshType> > dg_input);
27 
29  static std::shared_ptr<ODESolverBase<dim,real,MeshType>> create_ODESolver(std::shared_ptr< DGBase<dim, real, MeshType> > dg_input, std::shared_ptr<ProperOrthogonalDecomposition::PODBase<dim>> pod);
30 
32  static std::shared_ptr<ODESolverBase<dim,real,MeshType>> create_ODESolver(std::shared_ptr< DGBase<dim, real, MeshType> > dg_input, std::shared_ptr<ProperOrthogonalDecomposition::PODBase<dim>> pod, Epetra_Vector weights);
33 
35  static std::shared_ptr<ODESolverBase<dim,real,MeshType>> create_ODESolver_manual(Parameters::ODESolverParam::ODESolverEnum ode_solver_type, std::shared_ptr< DGBase<dim, real, MeshType> > dg_input);
36 
38  static std::shared_ptr<ODESolverBase<dim,real,MeshType>> create_ODESolver_manual(Parameters::ODESolverParam::ODESolverEnum ode_solver_type, std::shared_ptr< DGBase<dim, real, MeshType> > dg_input, std::shared_ptr<ProperOrthogonalDecomposition::PODBase<dim>> pod);
39 
41  static std::shared_ptr<ODESolverBase<dim,real,MeshType>> create_ODESolver_manual(Parameters::ODESolverParam::ODESolverEnum ode_solver_type, std::shared_ptr< DGBase<dim, real, MeshType> > dg_input, std::shared_ptr<ProperOrthogonalDecomposition::PODBase<dim>> pod, Epetra_Vector weights);
42 
44  static void display_error_ode_solver_factory(Parameters::ODESolverParam::ODESolverEnum ode_solver_type, bool reduced_order);
45 
47  static std::shared_ptr<ODESolverBase<dim,real,MeshType>> create_RungeKuttaODESolver(std::shared_ptr< DGBase<dim, real, MeshType> > dg_input);
48 
50  static std::shared_ptr<ODESolverBase<dim,real,MeshType>> create_RungeKuttaODESolver(std::shared_ptr< DGBase<dim, real, MeshType> > dg_input, std::shared_ptr<ProperOrthogonalDecomposition::PODBase<dim>> pod);
51 
53  static std::shared_ptr<RKTableauBase<dim,real,MeshType>> create_RKTableau(std::shared_ptr< DGBase<dim,real,MeshType> > dg_input);
54 
56  static std::shared_ptr<LowStorageRKTableauBase<dim,real,MeshType>> create_LowStorageRKTableau(std::shared_ptr< DGBase<dim,real,MeshType> > dg_input);
57 
59  static std::shared_ptr<EmptyRRKBase<dim,real,MeshType>> create_RRKObject(std::shared_ptr< DGBase<dim,real,MeshType> > dg_input,
60  std::shared_ptr<RKTableauBase<dim,real,MeshType>> rk_tableau);
61 
62 };
63 
64 } // ODE namespace
65 } // PHiLiP namespace
66 
67 #endif
static std::shared_ptr< ODESolverBase< dim, real, MeshType > > create_ODESolver_manual(Parameters::ODESolverParam::ODESolverEnum ode_solver_type, std::shared_ptr< DGBase< dim, real, MeshType > > dg_input)
Creates either implicit or explicit ODE solver based on manual input (no POD basis given) ...
static std::shared_ptr< LowStorageRKTableauBase< dim, real, MeshType > > create_LowStorageRKTableau(std::shared_ptr< DGBase< dim, real, MeshType > > dg_input)
Creates an RKTableau object based on the specified RK method.
Files for the baseline physics.
Definition: ADTypes.hpp:10
Create specified ODE solver as ODESolverBase object.
static std::shared_ptr< ODESolverBase< dim, real, MeshType > > create_RungeKuttaODESolver(std::shared_ptr< DGBase< dim, real, MeshType > > dg_input)
Creates an ODESolver object based on the specified RK method, including derived classes.
static void display_error_ode_solver_factory(Parameters::ODESolverParam::ODESolverEnum ode_solver_type, bool reduced_order)
Output error message for Implicit and Explicit solver.
static std::shared_ptr< RKTableauBase< dim, real, MeshType > > create_RKTableau(std::shared_ptr< DGBase< dim, real, MeshType > > dg_input)
Creates an RKTableau object based on the specified RK method.
Base class for storing the RK method.
DGBase is independent of the number of state variables.
Definition: dg_base.hpp:82
static std::shared_ptr< ODESolverBase< dim, real, MeshType > > create_ODESolver(std::shared_ptr< DGBase< dim, real, MeshType > > dg_input)
Creates either implicit or explicit ODE solver based on parameter value(no POD basis given) ...
static std::shared_ptr< EmptyRRKBase< dim, real, MeshType > > create_RRKObject(std::shared_ptr< DGBase< dim, real, MeshType > > dg_input, std::shared_ptr< RKTableauBase< dim, real, MeshType >> rk_tableau)
Creates an RRK object with specified RRK type; if no RRK is being used, creates an RRK object with em...