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); }
100 Vector *s = NULL;
return *s;
124 void Step(
Vector &x,
double &t,
double &dt)
override;
144 void Step(
Vector &x,
double &t,
double &dt)
override;
157 void Step(
Vector &x,
double &t,
double &dt)
override;
170 void Step(
Vector &x,
double &t,
double &dt)
override;
187 const double *a, *b, *c;
196 void Step(
Vector &x,
double &t,
double &dt)
override;
207 static const double a[28], b[8], c[7];
219 static const double a[66], b[12], c[11];
241 void Step(
Vector &x,
double &t,
double &dt)
override;
251 if (RKsolver) {
delete RKsolver; }
260 static const double a[1];
270 static const double a[2];
280 static const double a[3];
290 static const double a[4];
300 static const double a[5];
322 void Step(
Vector &x,
double &t,
double &dt)
override;
332 if (RKsolver) {
delete RKsolver; }
340 static const double a[1];
351 static const double a[2];
361 static const double a[3];
371 static const double a[4];
381 static const double a[5];
397 void Step(
Vector &x,
double &t,
double &dt)
override;
410 void Step(
Vector &x,
double &t,
double &dt)
override;
431 void Step(
Vector &x,
double &t,
double &dt)
override;
445 void Step(
Vector &x,
double &t,
double &dt)
override;
459 void Step(
Vector &x,
double &t,
double &dt)
override;
481 void Step(
Vector &x,
double &t,
double &dt)
override;
517 while (t < tf) {
Step(q, p, t, dt); }
641 while (t < tf) {
Step(x, dxdt, t, dt); }
650 Vector *s = NULL;
return *s;
676 NewmarkSolver(
double beta_ = 0.25,
double gamma_ = 0.5) { beta = beta_; gamma = gamma_; };
719 rho_inf = (rho_inf > 1.0) ? 1.0 : rho_inf;
720 rho_inf = (rho_inf < 0.0) ? 0.0 : rho_inf;
722 alpha_m = (2.0 - rho_inf)/(1.0 + rho_inf);
786 rho_inf = (rho_inf > 1.0) ? 1.0 : rho_inf;
787 rho_inf = (rho_inf < 0.0) ? 0.0 : rho_inf;
void Step(Vector &x, double &t, double &dt) override
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in]...
void Init(TimeDependentOperator &_f) override
Associate a TimeDependentOperator with the ODE solver.
virtual ~SecondOrderODESolver()
void Step(Vector &x, double &t, double &dt) override
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in]...
void PrintProperties(std::ostream &out=mfem::out)
virtual void Init(Operator &P, TimeDependentOperator &F)
void Step(Vector &x, double &t, double &dt) override
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in]...
void Step(Vector &x, double &t, double &dt) override
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in]...
int GetMaxStateSize() override
Function for getting and setting the state vectors.
AverageAccelerationSolver()
int GetStateSize() override
void SetRhoInf(double rho_inf)
virtual void SetStateVector(int i, Vector &state)
void Init(TimeDependentOperator &_f) override
Associate a TimeDependentOperator with the ODE solver.
Base abstract class for first order time dependent operators.
virtual int GetMaxStateSize()
Function for getting and setting the state vectors.
void Init(TimeDependentOperator &_f) override
Associate a TimeDependentOperator with the ODE solver.
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 int GetStateSize()
void Step(Vector &x, Vector &dxdt, double &t, double &dt) override
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.
void Step(Vector &q, Vector &p, double &t, double &dt) override
void Init(TimeDependentOperator &_f) override
Associate a TimeDependentOperator with the ODE solver.
int GetStateSize() override
Abstract class for solving systems of ODEs: dx/dt = f(x,t)
void Init(SecondOrderTimeDependentOperator &_f) override
Associate a TimeDependentOperator with the ODE solver.
const Vector & GetStateVector(int i) override
Abstract class for solving systems of ODEs: d2x/dt2 = f(x,dx/dt,t)
void Init(TimeDependentOperator &_f) override
Associate a TimeDependentOperator with the ODE solver.
Backward Euler ODE solver. L-stable.
CentralDifferenceSolver()
Second Order Symplectic Integration Algorithm.
LinearAccelerationSolver()
void Step(Vector &q, Vector &p, double &t, double &dt) override
virtual const Vector & GetStateVector(int i)
void Step(Vector &x, Vector &dxdt, double &t, double &dt) override
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in]...
void Init(TimeDependentOperator &_f) override
Associate a TimeDependentOperator with the ODE solver.
int GetMaxStateSize() override
Function for getting and setting the state vectors.
void PrintProperties(std::ostream &out=mfem::out)
void Init(TimeDependentOperator &_f) override
Associate a TimeDependentOperator with the ODE solver.
AdamsMoultonSolver(int _s, const double *_a)
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]...
void SetStateVector(int i, Vector &state) override
GeneralizedAlpha2Solver(double rho_inf=1.0)
void Init(TimeDependentOperator &_f) override
Associate a TimeDependentOperator with the ODE solver.
virtual ~ExplicitRKSolver()
The classical midpoint method.
void mfem_error(const char *msg)
Function called when an error is encountered. Used by the macros MFEM_ABORT, MFEM_ASSERT, MFEM_VERIFY.
virtual const Vector & GetStateVector(int i)
virtual void Init(SecondOrderTimeDependentOperator &f)
Associate a TimeDependentOperator with the ODE solver.
void SetStateVector(int i, Vector &state) override
void Step(Vector &x, double &t, double &dt) override
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in]...
void Init(TimeDependentOperator &_f) override
Associate a TimeDependentOperator with the ODE solver.
void Step(Vector &x, double &t, double &dt) override
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in]...
First Order Symplectic Integration Algorithm.
ExplicitRKSolver(int _s, const double *_a, const double *_b, const double *_c)
WBZAlphaSolver(double rho_inf=1.0)
void Step(Vector &x, double &t, double &dt) override
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in]...
const Vector & GetStateVector(int i) override
void Step(Vector &x, double &t, double &dt) override
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in]...
virtual int GetMaxStateSize()
Function for getting and setting the state vectors.
The classical explicit forth-order Runge-Kutta method, RK4.
double p(const Vector &x, double t)
static MemoryType GetHostMemoryType()
Get the current Host MemoryType. This is the MemoryType used by most MFEM classes when allocating mem...
TimeDependentOperator * F_
MemoryType
Memory types supported by MFEM.
Third-order, strong stability preserving (SSP) Runge-Kutta method.
Base abstract class for second order time dependent operators.
void Init(TimeDependentOperator &_f) override
Associate a TimeDependentOperator with the ODE solver.
NewmarkSolver(double beta_=0.25, double gamma_=0.5)
void PrintProperties(std::ostream &out=mfem::out)
void Init(TimeDependentOperator &_f) override
Associate a TimeDependentOperator with the ODE solver.
void Init(SecondOrderTimeDependentOperator &_f) override
Associate a TimeDependentOperator with the ODE solver.
void SetStateVector(int i, Vector &state) override
RK2Solver(const double _a=2./3.)
Implicit midpoint method. A-stable, not L-stable.
Host memory; using new[] and delete[].
const Vector & GetStateVector(int i) override
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 int GetStateSize()
void Step(Vector &x, double &t, double &dt) override
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in]...
const Vector & GetStateVector(int i) override
AdamsBashforthSolver(int _s, const double *_a)
virtual void Step(Vector &q, Vector &p, double &t, double &dt)=0
void Step(Vector &x, double &t, double &dt) override
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in]...
int GetStateSize() override
void Step(Vector &x, double &t, double &dt) override
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in]...
GeneralizedAlphaSolver(double rho=1.0)
virtual void Run(Vector &q, Vector &p, double &t, double &dt, double tf)
virtual void GetStateVector(int i, Vector &state)
void Init(TimeDependentOperator &_f) override
Associate a TimeDependentOperator with the ODE solver.
int GetMaxStateSize() override
Function for getting and setting the state vectors.
virtual void GetStateVector(int i, Vector &state)
void Step(Vector &x, double &t, double &dt) override
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in]...
virtual void SetStateVector(int i, Vector &state)
virtual void Run(Vector &x, double &t, double &dt, double tf)
Perform time integration from time t [in] to time tf [in].
TimeDependentOperator * f
Pointer to the associated TimeDependentOperator.
void SetStateVector(int i, Vector &state) override
void Step(Vector &q, Vector &p, double &t, double &dt) override
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.
SDIRK23Solver(int gamma_opt=1)
void Step(Vector &x, double &t, double &dt) override
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in]...
HHTAlphaSolver(double alpha=1.0)
void Init(TimeDependentOperator &_f) override
Associate a TimeDependentOperator with the ODE solver.
int GetStateSize() override
int GetMaxStateSize() override
Function for getting and setting the state vectors.
Variable order Symplectic Integration Algorithm (orders 1-4)
SecondOrderTimeDependentOperator * f
Pointer to the associated TimeDependentOperator.