Tempus
Version of the Day
Time Integration
|
Runge-Kutta methods. More...
#include <Tempus_RKButcherTableau.hpp>
Public Member Functions | |
RKButcherTableau (std::string stepperType, const Teuchos::SerialDenseMatrix< int, Scalar > &A, const Teuchos::SerialDenseVector< int, Scalar > &b, const Teuchos::SerialDenseVector< int, Scalar > &c, const int order, const int orderMin, const int orderMax, const Teuchos::SerialDenseVector< int, Scalar > &bstar=Teuchos::SerialDenseVector< int, Scalar >(), bool checkC=true) | |
virtual std::size_t | numStages () const |
Return the number of stages. More... | |
virtual const Teuchos::SerialDenseMatrix < int, Scalar > & | A () const |
Return the matrix coefficients. More... | |
virtual const Teuchos::SerialDenseVector < int, Scalar > & | b () const |
Return the vector of quadrature weights. More... | |
virtual const Teuchos::SerialDenseVector < int, Scalar > & | bstar () const |
Return the vector of quadrature weights for embedded methods. More... | |
virtual const Teuchos::SerialDenseVector < int, Scalar > & | c () const |
Return the vector of stage positions. More... | |
virtual int | order () const |
Return the order. More... | |
virtual int | orderMin () const |
Return the minimum order. More... | |
virtual int | orderMax () const |
Return the maximum order. More... | |
virtual bool | isImplicit () const |
Return true if the RK method is implicit. More... | |
virtual bool | isDIRK () const |
Return true if the RK method is Diagonally Implicit. More... | |
virtual bool | isEmbedded () const |
Return true if the RK method has embedded capabilities. More... | |
virtual std::string | description () const |
virtual void | describe (Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const |
Protected Member Functions | |
void | set_isImplicit () |
void | set_isDIRK () |
DIRK is defined as if a_ij = 0 for j>i and a_ii != 0 for at least one i. More... | |
Protected Attributes | |
std::string | description_ |
Teuchos::SerialDenseMatrix < int, Scalar > | A_ |
Teuchos::SerialDenseVector < int, Scalar > | b_ |
Teuchos::SerialDenseVector < int, Scalar > | c_ |
int | order_ |
int | orderMin_ |
int | orderMax_ |
bool | isImplicit_ |
bool | isDIRK_ |
bool | isEmbedded_ |
Teuchos::SerialDenseVector < int, Scalar > | bstar_ |
Private Member Functions | |
RKButcherTableau () | |
Runge-Kutta methods.
This base class specifies the Butcher tableau which defines the Runge-Kutta (RK) method. Both explicit and implicit RK methods can be specified here, and of arbitrary number of stages and orders. Embedded methods are also supported.
Since this is a generic RK class, no low-storage methods are incorporated here, however any RK method with a Butcher tableau can be created with the base class.
There are over 40 derived RK methods that have been implemented, ranging from first order and eight order, and from single stage to 5 stages.
This class was taken and modified from Rythmos' RKButcherTableau class.
Definition at line 50 of file Tempus_RKButcherTableau.hpp.
|
inline |
Definition at line 56 of file Tempus_RKButcherTableau.hpp.
|
private |
|
inlinevirtual |
Return the matrix coefficients.
Definition at line 126 of file Tempus_RKButcherTableau.hpp.
|
inlinevirtual |
Return the vector of quadrature weights.
Definition at line 129 of file Tempus_RKButcherTableau.hpp.
|
inlinevirtual |
Return the vector of quadrature weights for embedded methods.
Definition at line 132 of file Tempus_RKButcherTableau.hpp.
|
inlinevirtual |
Return the vector of stage positions.
Definition at line 135 of file Tempus_RKButcherTableau.hpp.
|
inlinevirtual |
Definition at line 154 of file Tempus_RKButcherTableau.hpp.
|
inlinevirtual |
Definition at line 152 of file Tempus_RKButcherTableau.hpp.
|
inlinevirtual |
Return true if the RK method is Diagonally Implicit.
Definition at line 146 of file Tempus_RKButcherTableau.hpp.
|
inlinevirtual |
Return true if the RK method has embedded capabilities.
Definition at line 148 of file Tempus_RKButcherTableau.hpp.
|
inlinevirtual |
Return true if the RK method is implicit.
Definition at line 144 of file Tempus_RKButcherTableau.hpp.
|
inlinevirtual |
Return the number of stages.
Definition at line 124 of file Tempus_RKButcherTableau.hpp.
|
inlinevirtual |
Return the order.
Definition at line 138 of file Tempus_RKButcherTableau.hpp.
|
inlinevirtual |
Return the maximum order.
Definition at line 142 of file Tempus_RKButcherTableau.hpp.
|
inlinevirtual |
Return the minimum order.
Definition at line 140 of file Tempus_RKButcherTableau.hpp.
|
inlineprotected |
DIRK is defined as if a_ij = 0 for j>i and a_ii != 0 for at least one i.
Definition at line 188 of file Tempus_RKButcherTableau.hpp.
|
inlineprotected |
Definition at line 180 of file Tempus_RKButcherTableau.hpp.
|
protected |
Definition at line 201 of file Tempus_RKButcherTableau.hpp.
|
protected |
Definition at line 202 of file Tempus_RKButcherTableau.hpp.
|
protected |
Definition at line 210 of file Tempus_RKButcherTableau.hpp.
|
protected |
Definition at line 203 of file Tempus_RKButcherTableau.hpp.
|
protected |
Definition at line 199 of file Tempus_RKButcherTableau.hpp.
|
protected |
Definition at line 208 of file Tempus_RKButcherTableau.hpp.
|
protected |
Definition at line 209 of file Tempus_RKButcherTableau.hpp.
|
protected |
Definition at line 207 of file Tempus_RKButcherTableau.hpp.
|
protected |
Definition at line 204 of file Tempus_RKButcherTableau.hpp.
|
protected |
Definition at line 206 of file Tempus_RKButcherTableau.hpp.
|
protected |
Definition at line 205 of file Tempus_RKButcherTableau.hpp.