[P]arallel [Hi]gh-order [Li]brary for [P]DEs
Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
|
Create specified ODE solver as ODESolverBase object. More...
#include <ode_solver_factory.h>
Static Public Member Functions | |
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< ODESolverBase< dim, real, MeshType > > | create_ODESolver (std::shared_ptr< DGBase< dim, real, MeshType > > dg_input, std::shared_ptr< ProperOrthogonalDecomposition::PODBase< dim >> pod) |
Creates either POD-Galerkin or POD-Petrov-Galerkin ODE solver based on parameter value (POD basis given) | |
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) |
Creates Hyper-reduced POD-Petrov-Galerkin ODE solver based on parameter value (POD basis and ECSW weights given) | |
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< 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) |
Creates either POD-Galerkin or POD-Petrov-Galerkin ODE solver based on manual input (POD basis given) | |
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) |
Creates Hyper-reduced POD-Petrov-Galerkin ODE solver based on manual input (POD basis and ECSW weights given) | |
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< 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 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) |
Creates an ODESolver object based on the specified RK method, with a POD Reduced Order Basis. | |
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. | |
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. | |
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 empty functions. | |
Create specified ODE solver as ODESolverBase object.
Factory design pattern whose job is to create the correct ODE solver
Definition at line 22 of file ode_solver_factory.h.