[P]arallel [Hi]gh-order [Li]brary for [P]DEs
Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
|
Base class ODE solver. More...
#include <ode_solver_base.h>
Public Member Functions | |
ODESolverBase (std::shared_ptr< DGBase< dim, real, MeshType > > dg_input, std::shared_ptr< ProperOrthogonalDecomposition::PODBase< dim >> pod) | |
Default constructor that will set the constants. More... | |
ODESolverBase (std::shared_ptr< DGBase< dim, real, MeshType > > dg_input) | |
Situational Constructor that will call the default with no POD. | |
virtual | ~ODESolverBase ()=default |
Destructor. | |
void | write_ode_solver_steady_state_convergence_data_to_table (const unsigned int current_iteration, const double current_residual, const std::shared_ptr< dealii::TableHandler > data_table) const |
Writes the ode solver steady state convergence data to a file. | |
virtual int | steady_state () |
Evaluate steady state solution. | |
void | initialize_steady_polynomial_ramping (const unsigned int global_final_poly_degree) |
Ramps up the solution from p0 all the way up to the given global_final_poly_degree. More... | |
void | valid_initial_conditions () const |
Checks whether the DG vector has valid values. More... | |
int | advance_solution_time (double time_advance) |
Function to advance solution to time+dt. | |
virtual void | step_in_time (real dt, const bool pseudotime)=0 |
Virtual function to evaluate solution update. | |
virtual double | get_automatic_error_adaptive_step_size (real dt, const bool pseudotime) |
Virtual function to evaluate automatic error adaptive time step. | |
virtual double | get_automatic_initial_step_size (real dt, const bool pseudotime) |
Virtual function to evaluate initial automatic error adaptive time step. | |
virtual void | allocate_ode_system ()=0 |
Virtual function to allocate the ODE system. | |
double | get_original_time_step () const |
Getter for original_time_step. | |
double | get_modified_time_step () const |
Getter for modified_time_step. | |
Public Attributes | |
dealii::TableHandler | solutions_table |
Table used to output solution vector at each time step. | |
double | residual_norm |
Current residual norm. Only makes sense for steady state. | |
double | residual_norm_decrease |
Current residual norm normalized by initial residual. Only makes sense for steady state. | |
std::shared_ptr< DGBase< dim, real, MeshType > > | dg |
Smart pointer to DGBase. | |
std::shared_ptr< ProperOrthogonalDecomposition::PODBase< dim > > | pod |
Smart pointer to PODBasis. | |
std::unique_ptr< BoundPreservingLimiter< dim, real > > | limiter |
Pointer to BoundPreservingLimiter. | |
double | current_time |
Useful for accurate time-stepping. More... | |
unsigned int | current_iteration |
Current iteration. More... | |
double | current_desired_time_for_output_solution_every_dt_time_intervals |
Current desired time for output solution every dt time intervals. More... | |
double | FR_entropy_contribution_RRK_solver = 0 |
Entropy FR correction at the current timestep. More... | |
double | relaxation_parameter_RRK_solver =1 |
Relaxation parameter. More... | |
Protected Attributes | |
double | CFL_factor |
CFL factor for (un)successful linesearches. More... | |
double | update_norm |
Norm of the solution update. | |
double | initial_residual_norm |
Initial residual norm. | |
dealii::LinearAlgebra::distributed::Vector< double > | solution_update |
Solution update given by the ODE solver. | |
const Parameters::AllParameters *const | all_parameters |
Input parameters. | |
const Parameters::ODESolverParam | ode_param |
Input ODE solver parameters. | |
double | original_time_step |
Original time step before calling step_in_time. | |
double | modified_time_step |
Modified time step after calling step_in_time. | |
const MPI_Comm | mpi_communicator |
MPI communicator. | |
const int | mpi_rank |
MPI rank. | |
dealii::ConditionalOStream | pcout |
Parallel std::cout that only outputs on mpi_rank==0. | |
Base class ODE solver.
Definition at line 26 of file ode_solver_base.h.
|
explicit |
Default constructor that will set the constants.
Default Constructor.
Definition at line 8 of file ode_solver_base.cpp.
void PHiLiP::ODE::ODESolverBase< dim, real, MeshType >::initialize_steady_polynomial_ramping | ( | const unsigned int | global_final_poly_degree | ) |
Ramps up the solution from p0 all the way up to the given global_final_poly_degree.
This first interpolates the current solution to the P0 space as an initial solution. The P0 is then solved, interpolated to P1, and the process is repeated until the desired polynomial is solved.
This is mainly usely for grid studies.
Definition at line 57 of file ode_solver_base.cpp.
void PHiLiP::ODE::ODESolverBase< dim, real, MeshType >::valid_initial_conditions | ( | ) | const |
Checks whether the DG vector has valid values.
By default, the DG solution vector is initialized with the lowest possible value.
Definition at line 86 of file ode_solver_base.cpp.
|
protected |
CFL factor for (un)successful linesearches.
When the linesearch succeeds on its first try, double the CFL on top of the CFL ramping. If the linesearch fails and needs to look at the other direction or accept a higher residual, halve the CFL on top of the residual (de)ramping
Definition at line 90 of file ode_solver_base.h.
double PHiLiP::ODE::ODESolverBase< dim, real, MeshType >::current_desired_time_for_output_solution_every_dt_time_intervals |
Current desired time for output solution every dt time intervals.
This variable will change when advance_solution_time() is called if ODE parameter "output_solution_every_dt_time_intervals" > 0.
Definition at line 127 of file ode_solver_base.h.
unsigned int PHiLiP::ODE::ODESolverBase< dim, real, MeshType >::current_iteration |
Current iteration.
This variable will change when step_in_time() is called.
Definition at line 122 of file ode_solver_base.h.
double PHiLiP::ODE::ODESolverBase< dim, real, MeshType >::current_time |
Useful for accurate time-stepping.
This variable will change when step_in_time() is called.
Definition at line 118 of file ode_solver_base.h.
double PHiLiP::ODE::ODESolverBase< dim, real, MeshType >::FR_entropy_contribution_RRK_solver = 0 |
Entropy FR correction at the current timestep.
Used in entropy-RRK ODE solver. This is stored in ode_solver_base such that both flow solver case and ode solver can access it.
Definition at line 143 of file ode_solver_base.h.
double PHiLiP::ODE::ODESolverBase< dim, real, MeshType >::relaxation_parameter_RRK_solver =1 |
Relaxation parameter.
Used in RRK ODE solver. This is stored in ode_solver_base such that both flow solver case and ode solver can access it.
Definition at line 148 of file ode_solver_base.h.