4.11 Time Integration

This directory consists of the following files:



TimeIntegrationScheme.h TimeIntegrationScheme.cpp


TimeIntegrationWrapper.h TimeIntegrationWrapper.cpp


The original time-stepping classes (found in TimeIntegrationScheme) where implemented by Vos et al. [62] and a detailed explanation of the mathematics and computer science concepts are provided there. After the original implementation, Cantwell et al. extended Vos’ original work by adding the time-stepping routines into a factory pattern (found in TimeIntegrationWrapper).

4.11.1 General Linear Methods

General linear methods (GLMs) can be considered as the generalization of a broad range of different numerical methods for ordinary differential equations. They were introduced by Butcher and they provide a unified formulation for traditional methods such as the Runge-Kutta methods and the linear multi-step methods such as Adams-Bashforth. From an implementation point of view, this means that all these numerical methods can be abstracted in a similar way. As this allows a high level of generality, it is chosen in Nektar++ to cast all time integration schemes in the framework of general linear methods.

For background information about general linear methods, please consult [14].

4.11.2 Introduction

The standard initial value problem can written in the form

dy-
dt = f(t,y), y(t0) = y0

where y is a vector containing the variable (or an array of array containing the variables).

In the formulation of general linear methods, it is more convenient to consider the ODE in autonomous form, i.e.

d^y-
dt = ^ f (^y ), ^ y (t0) = ^ y 0.

4.11.3 Formulation

Suppose the governing differential equation is given in autonomous form, the nth step of the GLM comprising

is formulated as:

Y i = Δt∑ j=0s-1a ijFj + ∑ j=0r-1u ij^ y j[n-1], i = 0,1,…,s - 1
^y i[n] = Δt∑ j=0s-1b ijFj + ∑ j=0r-1v ij^ y j[n-1], i = 0,1,…,r - 1

where Y i are referred to as the stage values and Fj as the stage derivatives. Both quantities are related by the differential equation:

Fi = ^ f (Y i).

The matrices A = [aij], U = [uij], B = [bij], V = [vij] are characteristic of a specific method. Each scheme can then be uniquely defined by a partitioned (s + r) × (s + r) matrix

[       ]
  A  U
  B  V

4.11.4 Matrix notation

Adopting the notation:

^y [n-1] = ⌊   [n-1] ⌋
| y^0[n-1] |
|| y^1    ||
||    ..   ||
⌈   [n.-1] ⌉
  y^r-1, ^ y [n] = ⌊   [n] ⌋
|  ^y0[n] |
||  ^y1  ||
||   ..  ||
⌈  [.n] ⌉
  ^yr-1, Y = ⌊       ⌋
|  Y 0  |
||  Y 1  ||
|⌈   ...   |⌉

  Y s-1, F = ⌊      ⌋
|  F 0 |
||  F 1 ||
|⌈   ...  |⌉

  F s- 1

the general linear method can be written more compactly in the following form:

[     ]
  Y
   [n]
  ^y = [                  ]
  A ⊗ IN   U ⊗ IN
  B ⊗  I   V ⊗ I
       N        N [       ]
   ΔtF
    [n-1]
  y^

where IN is the identity matrix of dimension N × N.

4.11.5 General Linear Methods in Nektar++

Although the GLM is essentially presented for ODE’s in its autonomous form, in Nektar++ it will be used to solve ODE’s formulated in non-autonomous form. Given the ODE,

dy-
dt = f(t,y), y(t0) = y0

a single step of GLM can then be evaluated in the following way:

For a detailed description of the formulation and a deeper insight of the numerical method see [62].

4.11.6 Types of time Integration Schemes

Nektar++ contains various classes and methods which implement the concept of GLMs. This toolbox is capable of numerically solving the generalised ODE using a broad range of different time-stepping methods. We distinguish the following types of general linear methods:

The aim in Nektar++ is to fully support the first three types of GLMs. Fully implicit methods are currently not implemented.

4.11.7 Usage

The goal of abstracting the concept of general linear methods is to provide users with a single interface for time-stepping, independent of the chosen method. The classes tree allows the user to numerically integrate ODE’s using high-order complex schemes, as if it was done using the Forward-Euler method. Switching between time-stepping schemes is as easy as changing a parameter in an input file.

In the already implemented solvers the time-integration schemes have been set up according to the nature of the equations. For example the incompressible Navier-Stokes equations solver allows the use of three different Implicit-Explicit time-schemes if solving the equations with a splitting-scheme. This is because this kind of scheme has an explicit and an implicit operator that combined solve the ODE’s system.

Once aware of the problem’s nature and implementation, the user can easily switch between some (depending on the problem) of the following time-integration schemes:

AdamsBashforthOrder1

Adams-Bashforth Forward multi-step scheme of order 1

AdamsBashforthOrder2

Adams-Bashforth Forward multi-step scheme of order 2

AdamsBashforthOrder3

Adams-Bashforth Forward multi-step scheme of order 3

AdamsMoultonOrder1

Adams-Moulton Forward multi-step scheme of order 1

AdamsMoultonOrder2

Adams-Moulton Forward multi-step scheme of order 2

BDFImplicitOrder1

BDF multi-step scheme of order 1 (implicit)

BDFImplicitOrder2

BDF multi-step scheme of order 2 (implicit)

ClassicalRungeKutta4

Runge-Kutta multi-stage scheme 4th order explicit (old name)

RungeKutta4

Runge-Kutta multi-stage scheme 4th order explicit (new name for ClassicalRungeKutta4)

RungeKutta3_SSP

Nonlinear SSP RungeKutta3 explicit

RungeKutta2_ImprovedEuler

Improved RungeKutta2 explicit (old name meaning Heun’s method)

RungeKutta2_SSP

Nonlinear SSP RungeKutta2 explicit (surrogate for RungeKutta2_ImprovedEuler)

ForwardEuler

Forward-Euler scheme

BackwardEuler

Backward Euler scheme

IMEXOrder1

IMEX 1st order scheme using Euler Backwards Euler Forwards

IMEXOrder2

IMEX 2nd order scheme using Backward Different Formula & Extrapolation

IMEXOrder3

IMEX 3rd order scheme using Backward Different Formula & Extrapolation

Midpoint

Midpoint method (old name)

RungeKutta2

Classical RungeKutta2 method (new name for Midpoint)

DIRKOrder2

Diagonally Implicit Runge-Kutta scheme of order 2

DIRKOrder3

Diagonally Implicit Runge-Kutta scheme of order 3

CNAB

Crank-Nicolson-Adams-Bashforth Order 2 (CNAB)

IMEXGear

IMEX Gear Order 2

MCNAB

Modified Crank-Nicolson-Adams-Bashforth Order 2 (MCNAB)

IMEXdirk_1_1_1

Forward-Backward Euler IMEX DIRK(1,1,1)

IMEXdirk_1_2_1

Forward-Backward Euler IMEX DIRK(1,2,1)

IMEXdirk_1_2_2

Implicit-Explicit Midpoint IMEX DIRK(1,2,2)

IMEXdirk_2_2_2

L-stable, two stage, second order IMEX DIRK(2,2,2)

IMEXdirk_2_3_2

L-stable, three stage, third order IMEX DIRK(3,4,3)

IMEXdirk_2_3_3

L-stable, two stage, third order IMEX DIRK(2,3,3)

IMEXdirk_3_4_3

L-stable, three stage, third order IMEX DIRK(3,4,3)

IMEXdirk_4_4_3

L-stable, four stage, third order IMEX DIRK(4,4,3)

Nektar++ input file for your problem will ask you just the string corresponding the time-stepping scheme you want to use (between quotation marks in the previous list), and few parameters to define your integration in time (time-step and number of steps or final time). For example:

1<SOLVERINFO> 
2   <I PROPERTY="EQTYPE" VALUE="UnsteadyNavierStokes" /> 
3   <I PROPERTY="SolverType" VALUE="VelocityCorrectionScheme" /> 
4   <I PROPERTY="AdvectionForm" VALUE="Convective" /> 
5   <I PROPERTY="Projection" VALUE="Galerkin" /> 
6   <I PROPERTY="TimeIntegrationMethod" VALUE="IMEXOrder1" /> 
7 </SOLVERINFO> 
8 
9 <PARAMETERS> 
10   <P> TimeStep = 0.001     </P> 
11   <P> NumSteps = 1000      </P> 
12   <P> Kinvis = 1           </P> 
13</PARAMETERS>

4.11.8 Implementation of a time-dependent problem

In order to implement a new solver which takes advantage of the time-integration class in Nektar++, two main ingredients are required:

Your pseudo-main file, where you actually loop over the time steps, will look like

1NekDouble timestep    = 0.1; 
2NekDouble time        = 0.0; 
3int NumSteps          = 1000; 
4 
5YourClass solver(Input); 
6 
7LibUtilities::TimeIntegrationMethod  TIME_SCHEME; 
8LibUtilities::TimeIntegrationSchemeOperators ODE; 
9 
10ODE.DefineOdeRhs(&YourClass::YourExplicitOperatorFunction,solver); 
11ODE.DefineProjection(&YourClass::YourProjectionFunction,solver); 
12ODE.DefineImplicitSolve(&YourClass::YourImplicitOperatorFunction,solver); 
13 
14Array<OneD, LibUtilities::TimeIntegrationSchemeSharedPtr> IntScheme; 
15 
16LibUtilities::TimeIntegrationSolutionSharedPtr ode_solution; 
17 
18Array<OneD, Array<OneD, NekDouble> > U; 
19 
20TIME_SCHEME = LibUtilities::eForwardEuler; 
21int numMultiSteps=1; 
22IntScheme = Array<OneD,LibUtilities::TimeIntegrationSchemeSharedPtr>(numMultiSteps); 
23LibUtilities::TimeIntegrationSchemeKey IntKey(TIME_SCHEME); 
24IntScheme[0] = LibUtilities::TimeIntegrationSchemeManager()[IntKey]; 
25ode_solution = IntScheme[0]->InitializeScheme(timestep,U,time,ODE); 
26 
27for(int n = 0; n < NumSteps; ++n) { 
28    U = IntScheme[0]->TimeIntegrate(timestep,ode_solution,ODE); 
29}

We can distinguish three different sections in the code above

Definitions
1NekDouble timestep    = 0.1; 
2NekDouble time        = 0.0; 
3int NumSteps          = 1000; 
4 
5YourClass equation(Input); 
6 
7LibUtilities::TimeIntegrationMethod  TIME_SCHEME; 
8LibUtilities::TimeIntegrationSchemeOperators ODE; 
9 
10ODE.DefineOdeRhs(&YourClass::YourExplicitOperatorFunction,equation); 
11ODE.DefineProjection(&YourClass::YourProjectionFunction, equation); 
12ODE.DefineImplicitSolve(&YourClass::YourImplicitOperatorFunction, equation);

In this section you define the basic parameters (like time-step, initial time, etc.) and the time-integration objects. The operators are not all required, it depends on the nature of your problem and on the type of time integration schemes you want to use. In this case, the problem has been set up to work just with Forward-Euler, then for sure you will not need the implicit operator. An object named equation has been initialized, is an object of type YourClass, where your spatial discretization and the functions which actually represent your operators are implemented. An example of this class will be shown later in this page.

Initialisations
1Array<OneD,LibUtilities::TimeIntegrationSchemeSharedPtr> IntScheme; 
2 
3LibUtilities::TimeIntegrationSolutionSharedPtr ode_solution; 
4 
5Array<OneD, Array<OneD, NekDouble> > U; 
6 
7TIME_SCHEME = LibUtilities::eForwardEuler; 
8int numMultiSteps=1; 
9IntScheme = Array<OneD,LibUtilities::TimeIntegrationSchemeSharedPtr>(numMultiSteps); 
10LibUtilities::TimeIntegrationSchemeKey IntKey(TIME_SCHEME); 
11IntScheme[0] = LibUtilities::TimeIntegrationSchemeManager()[IntKey]; 
12ode_solution = IntScheme[0]->InitializeScheme(timestep,U,time,ODE);

The second part consists in the scheme initialization. In this example we set up just Forward-Euler, but we can set up more then one time-integration scheme and quickly switch between them from the input file. Forward-Euler does not require any other scheme for the start-up procedure. High order multi-step schemes may need lower-order schemes for the start up.

Integration
1for(int n = 0; n < NumSteps; ++n) { 
2    U = IntScheme[0]->TimeIntegrate(timestep,ode_solution,ODE); 
3}

The last step is the typical time-loop, where you iterate in time to get your new solution at each time-level. The solution at time tn+1 is stored into vector U (you need to properly initialize this vector). U is an Array of Arrays, where the first dimension corresponds to the number of variables (eg. u,v,w) and the second dimension corresponds to the variables size (e.g. the number of modes or the number of physical points).

The variable ODE is an object which contains the methods. A class representing a PDE equation (or a system of equations) must have a series of functions representing the implicit/explicit part of the method, which represents the reduction of the PDE’s to a system of ODE’s. The spatial discretization and the definition of this method should be implemented in YourClass. &YourClass::YourExplicitOperatorFunction is a functor, i.e. a pointer to a function where the method is implemented. equation is a pointer to the object, i.e. the class, where the function/method is implemented. Here a pseudo-example of the .h file of your hypothetical class representing the set of equations. The implementation of the functions is meant to be in the related .cpp file.

1class YourCalss 
2
3public: 
4    YourClass(INPUT); 
5 
6    ~YourClass(void); 
7 
8    void YourExplicitOperatorFunction( 
9        const Array<OneD, Array<OneD, NekDouble> > & inarray, 
10              Array<OneD, Array<OneD, NekDouble> > & outarray, 
11        const NekDouble time); 
12 
13    void YourProjectionFunction( 
14        const Array<OneD, Array<OneD, NekDouble> > & inarray, 
15              Array<OneD, Array<OneD, NekDouble> > & outarray, const 
16              NekDouble time); 
17 
18    void YourImplicitOperatorFunction( 
19        const Array<OneD, Array<OneD, NekDouble> > & inarray, 
20              Array<OneD, Array<OneD, NekDouble> > & outarray, const 
21              NekDouble time, 
22        const NekDouble lambda); 
23 
24    void InternalMethod1(...); 
25 
26    NekDouble internalvariale1; 
27 
28protected: 
29    ... 
30 
31private: 
32    ... 
33};

4.11.9 Strongly imposed essential boundary conditions

Dirichlet boundary conditions can be strongly imposed by lifting the known Dirichlet solution. This is equivalent to decompose the approximate solution y into an known part, yD, which satisfies the Dirichlet boundary conditions, and an unknown part, yH, which is zero on the Dirichlet boundaries, i.e.

y = yD + yH

In a Finite Element discretisation, this corresponds to splitting the solution vector of coefficients y into the known Dirichlet degrees of freedom yD and the unknown homogeneous degrees of freedom yH. Ordering the known coefficients first, this corresponds to:

y = [     ]
  yD
  yH

The generalised formulation of the general linear method (i.e. the introduction of a left hand side operator) allows for an easier treatment of these types of boundary conditions. To better appreciate this, consider the equation for the stage values for an explicit general linear method where both the left and right hand side operator are linear operators, i.e. they can be represented by a matrix.

MY i = Δt∑ j=0i-1a ijLY j + ∑ j=0r-1u ijyj[n-1], i = 0,1,…,s - 1

In case of a lifted known solution, this can be written as:

[ M DD   M  DH ]
    HD      HH
  M      M[ Y D  ]
    iH
  Y i = Δt∑ j=0i-1a ij[ LDD   LDH   ]
    HD    HH
  L     L [ Y D ]
    jH
   Yj + ∑ j=0r-1u ij[  yD[n-1]]
    jH[n-1]
  y j,
i = 0,1,…,s - 1

In order to calculate the stage values correctly, the explicit operator should now be implemented to do the following:

[  D ]
  bH
  b = [   DD    DH  ]
   LHD   LHH
   L     L[  D  ]
  yH
  y

Note that only the homogeneous part bH will be used to calculate the stage values. This means essentially that only the bottom part of the operation above, i.e. LHDyD + LHHyH is required. However, sometimes it might be more convenient to use/implement routines for the explicit operator that also calculate bD.

An implicit method should solve the system:

([               ]    [             ])
   M  DD   M DH         LDD   LDH
   M  HD  M  HH   -  λ  LHD   LHH[     ]
  yD
  yH = [              ]
  HDD    HDH
  HHD    HHH[     ]
  yD
  yH = [     ]
  bD
  bH

for the unknown vector y. This can be done in three steps:

4.11.10 How to add a new time-stepping method

To add a new time integration scheme, follow the steps below:

4.11.11 Examples of already implemented time stepping schemes

Here we show some examples time-stepping schemes implemented in Nektar++, to give an idea of what is required to add one of them.

Forward Euler

[       ]
  A  |U
 -B--|V--
     | = [      ]
  0 |1
 -1-|1-
    |, y[n] = [     ]
  y[n0] = [           ]
  m (tn,yn)

Backward Euler

[    |  ]
  A  |U
 -B--|V--
     | = [   |  ]
  1 |1
 -1-|1-
    |, y[n] = [  [n] ]
  y0 = [     n  n  ]
  m (t ,y )

2nd order Adams-Bashforth

[   |   ]
  A |U
 -B-|V---
    | = ⌊    |      ⌋
| -03-|1--0-1-|
⌈  2 |1  -2-⌉
   1 |0  0, y[n] = [      ]
  y [n]
    0[n]
  y 1 = [                ]
    m  (tn,yn )
  Δtl(tn-1,yn-1)

1st order IMEX Euler-Backward/ Euler-Forward

[            |  ]
 -AIM--AEM---|U--
  BIM  BEM   |V = ⌊ [   ]  [   ] |[      ] ⌋
    1      0   |  1  1
||-[---]--[---]-|[------]-||
⌈   1      0   |  1  1   ⌉
    0      1   |  0  0 with y[n] = [      ]
  y [n0]
  y [n]
    1 = [            ]
   m (tn,yn)
  Δtl (tn,yn )

2nd order IMEX Backward Different Formula & Extrapolation

[            |   ]
 -AIM---AEM--|U--
  BIM   BEM  |V = ⌊  [   ]  [   ] | [                ] ⌋
     2      0   |   4  - 1  4  - 2
||-⌊--32-⌋--⌊---⌋-|⌊--34----31--34----32-⌋-||
|| |  3 |  | 0 | ||  3  - 3  3  - 3 | ||
|| ||  0 ||  || 0 || |||  1   0   0   0  || ||
⌈ ⌈  0 ⌉  ⌈ 1 ⌉ |⌈  0   0   0   0  ⌉ ⌉
     0      0   |   0   0   1   0

with

y[n] = ⌊     ⌋
  y[0n]
|| y[n]||
||  1[n]||
⌈ y2  ⌉
  y[3n] = ⌊ m (tn,yn)      ⌋
| m (tn-1,yn-1)  |
||      n  n      ||
⌈ Δtl((tn,-y1 )n- 1) ⌉
  Δtl t   ,y

Note: The first two rows are normalised so the coefficient on yn[n+1] is one. In the standard formulation it is 3∕2.

3rdorder IMEX Backward Different Formula & Extrapolation

[            |   ]
  AIM   AEM  |U
 -BIM---BEM--|V--
             | =    [    ]  [   ] | [                           ]
⌊    6-      0   |   18- - -9  -2  18  - 18  6-   ⌋
|-⌊--161-⌋--⌊---⌋-|⌊--1118----119--112--1118----1118--161-⌋-|
||    11      0   |   11- - 11  11  11  - 11  11   ||
|| ||  0  ||  || 0 || |||  1    0    0   0   0    0  || ||
|| ||  0  ||  || 0 || |||  0    1    0   0   0    0  || ||
|| ||  0  ||  || 1 || |||  0    0    0   0   0    0  || ||
|⌈ |⌈  0  |⌉  |⌈ 0 |⌉ ||⌈  0    0    0   1   0    0  |⌉ |⌉
     0       0   |   0    0    0   0   1    0

with

y[n] = ⌊     ⌋
  y[0n]
||  [n]||
|| y1[n]||
|| y2  ||
|| y[3n]||
|| y[n]||
⌈  4[n]⌉
  y5 = ⌊     n  n       ⌋
| m ((tn,-y1 )n-1)  |
|| m (t   ,y   )  ||
|| m  tn-2,yn-2   ||
|| Δtl((tn,yn)   ) ||
⌈ Δtl tn-1,yn- 1 ⌉
  Δtl(tn-2,yn- 2)

Note: The first two rows are normalised so the coefficient on yn[n+1] is one. In the standard formulation it is 11∕6.

2nd order Adams-Moulton

[  A |U  ]
  ---|----
  B  |V = ⌊   |      ⌋
 -12-|1---12-
|⌈ 12 |1   12 |⌉
  1 |0   0, y[n] = [ y[n]]
   0[n]
  y1 = [  m (tn,yn) ]
       n   n
  Δtl (t ,y  )

Midpoint Method

[    |   ]
 -A--|U--
  B  |V = ⌊      |  ⌋
| 01  0 |1 |
⌈-2--0-|1-⌉
  0  1 |1, y[n] = [   [n] ]
  y 0 = [    n   n  ]
 m  (t ,y )

RK4: the standard fourth-order Runge-Kutta scheme

[    |  ]
 -A--|U--
  B  |V =               |
⌊ 0  0  0  0  |1 ⌋
| 1  0  0  0  |1 |
|| 20  1  0  0  |1 ||
||    2        |  ||
⌈-01--01--11--01-|1-⌉
  6  3  3   6 |1, y[n] = [  [n] ]
  y0 = [     n  n  ]
  m (t ,y )

2nd order Diagonally Implicit Runge-Kutta (DIRK)

[   |   ]
 -A--U---
  B |V
    | = ⌊            |  ⌋
|    λ     0 |1 |
⌈-(1---λ)--λ-|1-⌉
  (1 - λ)  λ |1 with λ =     √--
2----2-
   2, y[n] = [  [n]]
  y0 = [     n   n ]
  m  (t ,y  )

3rd order Diagonally Implicit Runge-Kutta (DIRK)

[    |   ]
 -A--|U--
  B  |V = ⌊                                         |  ⌋
|      1  λ                  0          0 |1 |
|| 1 (  2 (21 - λ)   )  1(  2  λ       )  0 |1 ||
⌈-4-(- 6λ-+-16λ---1)--4(6λ----20λ-+-5)-λ--|1-⌉
  14  - 6λ2 + 16λ - 1  14 6λ2 - 20λ + 5  λ  |1

with

λ = 0.4358665215, y[n] = [  [n] ]
 y 0 = [ m (tn,yn) ]

3rd order L-stable, three stage IMEX DIRK(3,4,3)

[ AIM-AEM-|U--]
 BIM  BEM  V=
⌊ ⌊                                   ⌋                                         |      ⌋
    00        0λ              00        00    [     0          0         0      0 ] |[ 1 ]
|| ⌈ 0   (  12(1- λ)  )    (   λ     )  0⌉     0.3212λ788860  0.39660543747     00      00   |  11   ||
⌈-[-0--14(-6λ2+16λ--1)--14-6(λ2--20λ+-5)--λ]---[--0.1058(58296--0.55292)91479(-0.5529291479)--0-]-|[-1-]-⌉
    0  14 -6λ2+16λ- 1  14 6λ2- 20λ+ 5   λ      0 14 - 6λ2+ 16λ - 1  14 6λ2- 20λ+ 5 λ   |  1

with

λ = 0.4358665215, y[n] = [  [n] ]
 y 0 = [           ]
  m (tn,yn)