15 #include "../config/config.hpp"
71 virtual void Step(
Vector &x,
double &t,
double &dt) = 0;
89 virtual void Run(
Vector &x,
double &t,
double &dt,
double tf)
91 while (t < tf) {
Step(x, t, dt); }
107 virtual void Step(
Vector &x,
double &t,
double &dt);
127 virtual void Step(
Vector &x,
double &t,
double &dt);
140 virtual void Step(
Vector &x,
double &t,
double &dt);
153 virtual void Step(
Vector &x,
double &t,
double &dt);
170 const double *a, *b, *c;
179 virtual void Step(
Vector &x,
double &t,
double &dt);
190 static const double a[28], b[8], c[7];
202 static const double a[66], b[12], c[11];
224 virtual void Step(
Vector &x,
double &t,
double &dt);
228 if (RKsolver) {
delete RKsolver; }
237 static const double a[1];
247 static const double a[2];
257 static const double a[3];
267 static const double a[4];
277 static const double a[5];
299 virtual void Step(
Vector &x,
double &t,
double &dt);
303 if (RKsolver) {
delete RKsolver; }
311 static const double a[1];
322 static const double a[2];
332 static const double a[3];
342 static const double a[4];
352 static const double a[5];
368 virtual void Step(
Vector &x,
double &t,
double &dt);
381 virtual void Step(
Vector &x,
double &t,
double &dt);
402 virtual void Step(
Vector &x,
double &t,
double &dt);
416 virtual void Step(
Vector &x,
double &t,
double &dt);
430 virtual void Step(
Vector &x,
double &t,
double &dt);
452 virtual void Step(
Vector &x,
double &t,
double &dt);
482 while (t < tf) {
Step(q, p, t, dt); }
606 while (t < tf) {
Step(x, dxdt, t, dt); }
626 NewmarkSolver(
double beta_ = 0.25,
double gamma_ = 0.5) { beta = beta_; gamma = gamma_; };
669 rho_inf = (rho_inf > 1.0) ? 1.0 : rho_inf;
670 rho_inf = (rho_inf < 0.0) ? 0.0 : rho_inf;
672 alpha_m = (2.0 - rho_inf)/(1.0 + rho_inf);
730 rho_inf = (rho_inf > 1.0) ? 1.0 : rho_inf;
731 rho_inf = (rho_inf < 0.0) ? 0.0 : rho_inf;
virtual ~SecondOrderODESolver()
virtual void Init(TimeDependentOperator &_f)
Associate a TimeDependentOperator with the ODE solver.
void PrintProperties(std::ostream &out=mfem::out)
virtual void Init(Operator &P, TimeDependentOperator &F)
AverageAccelerationSolver()
void SetRhoInf(double rho_inf)
void Step(Vector &q, Vector &p, double &t, double &dt)
virtual void Step(Vector &x, Vector &dxdt, double &t, double &dt)
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in]...
virtual void Init(TimeDependentOperator &_f)
Associate a TimeDependentOperator with the ODE solver.
virtual void Step(Vector &x, double &t, double &dt)
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in]...
Base abstract class for first order time dependent operators.
virtual void Step(Vector &x, double &t, double &dt)=0
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in]...
virtual void Init(TimeDependentOperator &f)
Associate a TimeDependentOperator with the ODE solver.
Abstract class for solving systems of ODEs: dx/dt = f(x,t)
virtual void Step(Vector &x, Vector &dxdt, double &t, double &dt)
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in]...
virtual void Init(SecondOrderTimeDependentOperator &_f)
Associate a TimeDependentOperator with the ODE solver.
virtual void Step(Vector &x, double &t, double &dt)
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in]...
virtual void Step(Vector &x, double &t, double &dt)
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in]...
Abstract class for solving systems of ODEs: d2x/dt2 = f(x,dx/dt,t)
Backward Euler ODE solver. L-stable.
CentralDifferenceSolver()
virtual void Step(Vector &x, double &t, double &dt)
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in]...
LinearAccelerationSolver()
virtual void PrintProperties(std::ostream &out=mfem::out)
AdamsMoultonSolver(int _s, const double *_a)
virtual void Step(Vector &x, double &t, double &dt)
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in]...
virtual void Step(Vector &x, Vector &dxdt, double &t, double &dt)=0
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in]...
virtual void Init(TimeDependentOperator &_f)
Associate a TimeDependentOperator with the ODE solver.
GeneralizedAlpha2Solver(double rho_inf=1.0)
virtual void Init(TimeDependentOperator &_f)
Associate a TimeDependentOperator with the ODE solver.
virtual void Init(TimeDependentOperator &_f)
Associate a TimeDependentOperator with the ODE solver.
virtual ~ExplicitRKSolver()
The classical midpoint method.
virtual void Step(Vector &x, double &t, double &dt)
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in]...
virtual void Init(SecondOrderTimeDependentOperator &f)
Associate a TimeDependentOperator with the ODE solver.
virtual void Init(SecondOrderTimeDependentOperator &_f)
Associate a TimeDependentOperator with the ODE solver.
virtual void Init(TimeDependentOperator &_f)
Associate a TimeDependentOperator with the ODE solver.
virtual void Step(Vector &x, double &t, double &dt)
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in]...
ExplicitRKSolver(int _s, const double *_a, const double *_b, const double *_c)
WBZAlphaSolver(double rho_inf=1.0)
virtual void Step(Vector &x, double &t, double &dt)
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in]...
virtual void Init(TimeDependentOperator &_f)
Associate a TimeDependentOperator with the ODE solver.
virtual void Init(TimeDependentOperator &_f)
Associate a TimeDependentOperator with the ODE solver.
virtual void Init(TimeDependentOperator &_f)
Associate a TimeDependentOperator with the ODE solver.
The classical explicit forth-order Runge-Kutta method, RK4.
static MemoryType GetHostMemoryType()
Get the current Host MemoryType. This is the MemoryType used by most MFEM classes when allocating mem...
virtual void Init(TimeDependentOperator &_f)
Associate a TimeDependentOperator with the ODE solver.
TimeDependentOperator * F_
MemoryType
Memory types supported by MFEM.
virtual void Init(TimeDependentOperator &_f)
Associate a TimeDependentOperator with the ODE solver.
Third-order, strong stability preserving (SSP) Runge-Kutta method.
Base abstract class for second order time dependent operators.
NewmarkSolver(double beta_=0.25, double gamma_=0.5)
virtual void PrintProperties(std::ostream &out=mfem::out)
RK2Solver(const double _a=2./3.)
Implicit midpoint method. A-stable, not L-stable.
Host memory; using new[] and delete[].
void Step(Vector &q, Vector &p, double &t, double &dt)
void Step(Vector &q, Vector &p, double &t, double &dt)
virtual void Run(Vector &x, Vector &dxdt, double &t, double &dt, double tf)
Perform time integration from time t [in] to time tf [in].
virtual void Step(Vector &x, double &t, double &dt)
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in]...
virtual void Step(Vector &x, double &t, double &dt)
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in]...
AdamsBashforthSolver(int _s, const double *_a)
virtual void Step(Vector &q, Vector &p, double &t, double &dt)=0
GeneralizedAlphaSolver(double rho=1.0)
virtual void Run(Vector &q, Vector &p, double &t, double &dt, double tf)
virtual void Init(TimeDependentOperator &_f)
Associate a TimeDependentOperator with the ODE solver.
virtual void Run(Vector &x, double &t, double &dt, double tf)
Perform time integration from time t [in] to time tf [in].
virtual void Step(Vector &x, double &t, double &dt)
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in]...
TimeDependentOperator * f
Pointer to the associated TimeDependentOperator.
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
The classical forward Euler method.
virtual void Step(Vector &x, double &t, double &dt)
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in]...
SDIRK23Solver(int gamma_opt=1)
HHTAlphaSolver(double alpha=1.0)
virtual void Step(Vector &x, double &t, double &dt)
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in]...
SecondOrderTimeDependentOperator * f
Pointer to the associated TimeDependentOperator.
virtual void Init(TimeDependentOperator &_f)
Associate a TimeDependentOperator with the ODE solver.