MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
mfem::ODESolver Class Referenceabstract

Abstract class for solving systems of ODEs: dx/dt = f(x,t) More...

#include <ode.hpp>

Inheritance diagram for mfem::ODESolver:
[legend]
Collaboration diagram for mfem::ODESolver:
[legend]

Public Member Functions

 ODESolver ()
 
virtual void Init (TimeDependentOperator &f_)
 Associate a TimeDependentOperator with the ODE solver.
 
virtual void Step (Vector &x, real_t &t, real_t &dt)=0
 Perform a time step from time t [in] to time t [out] based on the requested step size dt [in].
 
virtual void Run (Vector &x, real_t &t, real_t &dt, real_t tf)
 Perform time integration from time t [in] to time tf [in].
 
virtual int GetMaxStateSize ()
 Function for getting and setting the state vectors.
 
virtual int GetStateSize ()
 
virtual const VectorGetStateVector (int i)
 
virtual void GetStateVector (int i, Vector &state)
 
virtual void SetStateVector (int i, Vector &state)
 
virtual ~ODESolver ()
 

Protected Attributes

TimeDependentOperatorf
 Pointer to the associated TimeDependentOperator.
 
MemoryType mem_type
 

Detailed Description

Abstract class for solving systems of ODEs: dx/dt = f(x,t)

Definition at line 23 of file ode.hpp.

Constructor & Destructor Documentation

◆ ODESolver()

mfem::ODESolver::ODESolver ( )
inline

Definition at line 31 of file ode.hpp.

◆ ~ODESolver()

virtual mfem::ODESolver::~ODESolver ( )
inlinevirtual

Definition at line 112 of file ode.hpp.

Member Function Documentation

◆ GetMaxStateSize()

virtual int mfem::ODESolver::GetMaxStateSize ( )
inlinevirtual

Function for getting and setting the state vectors.

Reimplemented in mfem::AdamsBashforthSolver, mfem::AdamsMoultonSolver, and mfem::GeneralizedAlphaSolver.

Definition at line 96 of file ode.hpp.

◆ GetStateSize()

virtual int mfem::ODESolver::GetStateSize ( )
inlinevirtual

Reimplemented in mfem::AdamsBashforthSolver, mfem::AdamsMoultonSolver, and mfem::GeneralizedAlphaSolver.

Definition at line 97 of file ode.hpp.

◆ GetStateVector() [1/2]

virtual const Vector & mfem::ODESolver::GetStateVector ( int i)
inlinevirtual

Reimplemented in mfem::AdamsBashforthSolver, mfem::AdamsMoultonSolver, and mfem::GeneralizedAlphaSolver.

Definition at line 98 of file ode.hpp.

◆ GetStateVector() [2/2]

virtual void mfem::ODESolver::GetStateVector ( int i,
Vector & state )
inlinevirtual

◆ Init()

void mfem::ODESolver::Init ( TimeDependentOperator & f_)
virtual

◆ Run()

virtual void mfem::ODESolver::Run ( Vector & x,
real_t & t,
real_t & dt,
real_t tf )
inlinevirtual

Perform time integration from time t [in] to time tf [in].

Parameters
[in,out]xApproximate solution.
[in,out]tTime associated with the approximate solution x.
[in,out]dtTime step size.
[in]tfRequested final time.

The default implementation makes consecutive calls to Step() until reaching tf. The following rules describe the common behavior of the method:

  • The input x [in] is the approximate solution for the input time t [in].
  • The input dt [in] is the initial time step size.
  • The output dt [out] is the last time step taken by the method which may be smaller or larger than the input dt [in] value, e.g. because of time step control.
  • The output value of t [out] is not smaller than tf [in].

Reimplemented in mfem::PetscODESolver.

Definition at line 90 of file ode.hpp.

◆ SetStateVector()

virtual void mfem::ODESolver::SetStateVector ( int i,
Vector & state )
inlinevirtual

◆ Step()

virtual void mfem::ODESolver::Step ( Vector & x,
real_t & t,
real_t & dt )
pure virtual

Perform a time step from time t [in] to time t [out] based on the requested step size dt [in].

Parameters
[in,out]xApproximate solution.
[in,out]tTime associated with the approximate solution x.
[in,out]dtTime step size.

The following rules describe the common behavior of the method:

  • The input x [in] is the approximate solution for the input time t [in].
  • The input dt [in] is the desired time step size, defining the desired target time: t [target] = t [in] + dt [in].
  • The output x [out] is the approximate solution for the output time t [out].
  • The output dt [out] is the last time step taken by the method which may be smaller or larger than the input dt [in] value, e.g. because of time step control.
  • The method may perform more than one time step internally; in this case dt [out] is the last internal time step size.
  • The output value of t [out] may be smaller or larger than t [target], however, it is not smaller than t [in] + dt [out], if at least one internal time step was performed.
  • The value x [out] may be obtained by interpolation using internally stored data.
  • In some cases, the contents of x [in] may not be used, e.g. when x [out] from a previous Step() call was obtained by interpolation.
  • In consecutive calls to this method, the output t [out] of one Step() call has to be the same as the input t [in] to the next Step() call.
  • If the previous rule has to be broken, e.g. to restart a time stepping sequence, then the ODE solver must be re-initialized by calling Init() between the two Step() calls.

Implemented in mfem::AdamsBashforthSolver, mfem::AdamsMoultonSolver, mfem::BackwardEulerSolver, mfem::ESDIRK32Solver, mfem::ESDIRK33Solver, mfem::ExplicitRKSolver, mfem::ForwardEulerSolver, mfem::GeneralizedAlphaSolver, mfem::ImplicitMidpointSolver, mfem::PetscODESolver, mfem::RK2Solver, mfem::RK3SSPSolver, mfem::RK4Solver, mfem::SDIRK23Solver, mfem::SDIRK33Solver, mfem::SDIRK34Solver, and mfem::TrapezoidalRuleSolver.

Member Data Documentation

◆ f

TimeDependentOperator* mfem::ODESolver::f
protected

Pointer to the associated TimeDependentOperator.

Definition at line 27 of file ode.hpp.

◆ mem_type

MemoryType mfem::ODESolver::mem_type
protected

Definition at line 28 of file ode.hpp.


The documentation for this class was generated from the following files: