TransientProblemInterface
-
class TransientProblemInterface
Interface for transient simulations.
Subclassed by godzilla::ExplicitProblemInterface, godzilla::ImplicitFENonlinearProblem
Public Types
-
enum class TimeScheme
Time stepping schemes.
Values:
-
enumerator BEULER
Backward Euler.
-
enumerator CN
Crank-Nicolson.
-
enumerator EULER
Forward Euler.
-
enumerator RK_2
Runge-Kutta 2 (midpoint)
-
enumerator HEUN
Heun’s.
-
enumerator SSP_RK_2
Strong stability preserving RK-2.
-
enumerator SSP_RK_3
Strong stability preserving RK-3.
-
enumerator BEULER
Public Functions
-
void set_time_stepping_adaptor(TimeSteppingAdaptor *adaptor)
Set time stepping adaptivity using an adaptor class.
- Parameters:
adaptor – Time stepping adaptor object
-
TimeSteppingAdaptor *get_time_stepping_adaptor() const
Get time stepping adaptor.
- Returns:
Time stepping adaptor
-
TS get_ts() const
Get TS object.
- Returns:
PETSc TS object
-
Real get_time_step() const
Get the current timestep size.
- Returns:
the current timestep size
-
void set_time_step(Real dt) const
Allows one to reset the timestep at any time, useful for simple pseudo-timestepping codes.
- Parameters:
dt – The size of the timestep
-
void set_max_time(Real max_time)
Sets the maximum (or final) time for time-stepping.
- Parameters:
max_time – Final time to step to
-
Real get_max_time() const
Gets the maximum (or final) time for time-stepping.
- Returns:
Final time to step to
-
void set_converged_reason(TSConvergedReason reason)
Sets the reason for handling the convergence.
- Parameters:
reason – Converged reason
-
TSConvergedReason get_converged_reason() const
Gets the reason the time iteration was stopped.
NOTE: Can only be called after the call to solve() is complete.
- Returns:
Converged reason
-
void set_scheme(TimeScheme scheme)
Set time-stepping scheme.
-
void set_scheme(const std::string &scheme_name)
Set time-stepping scheme.
-
std::string get_scheme() const
Get the name of time stepping scheme.
-
void set_time(Real t)
Set simulation time.
- Parameters:
t – New simulation time
-
virtual void pre_step()
Called before the time step solve.
-
virtual void post_step()
Called after the time step is done solving.
-
virtual void pre_stage(Real stage_time)
Runs the user-defined pre-stage function.
- Parameters:
stage_time – The absolute time of the current stage
-
virtual void post_stage(Real stage_time, Int stage_index, const std::vector<Vector> &Y)
Runs the user-defined post-stage function.
- Parameters:
stage_time – The absolute time of the current stage
stage_index – Stage number
Y – Array of vectors (of size = total number of stages) with the stage solutions
-
void compute_rhs(Real time, const Vector &x, Vector &F)
Compute right-hand side.
- Parameters:
time – Current time
x – Solution at time
time
F – Right-hand side vector
-
enum class TimeScheme