pf_mod_imk Module

This module implements fully implicit Munthe-Kaas Runge Kutta methods using explicit SDC sweeping

The equation to be solved is

where is a matrix and (y)\ is a vector or matrix or if Lax_pair = true

where both and are matrices

We solve this by finding the solution to

Using PFASST


Uses

  • module~~pf_mod_imk~~UsesGraph module~pf_mod_imk pf_mod_imk module~pf_mod_utils pf_mod_utils module~pf_mod_imk->module~pf_mod_utils module~pf_mod_dtype pf_mod_dtype module~pf_mod_imk->module~pf_mod_dtype module~pf_mod_utils->module~pf_mod_dtype module~pf_mod_timer pf_mod_timer module~pf_mod_utils->module~pf_mod_timer iso_c_binding iso_c_binding module~pf_mod_dtype->iso_c_binding module~pf_mod_timer->module~pf_mod_dtype

Contents


Interfaces

interface

  • public subroutine pf_f_eval_p(this, y, t, level, f)

    Subroutine f_eval computes A(y,t)

    Arguments

    Type IntentOptional AttributesName
    class(pf_imk_t), intent(inout) :: this
    class(pf_encap_t), intent(inout) :: y
    real(kind=pfdp), intent(in) :: t
    integer(kind=c_int), intent(in) :: level
    class(pf_encap_t), intent(inout) :: f

interface

  • public subroutine pf_dexpinv_p(this, a, omega, f)

    Subroutine dexpinv computes Om'=F=dexpinv_Om(A)

    Arguments

    Type IntentOptional AttributesName
    class(pf_imk_t), intent(inout) :: this
    class(pf_encap_t), intent(inout) :: a
    class(pf_encap_t), intent(inout) :: omega
    class(pf_encap_t), intent(inout) :: f

    The resultign-level

interface

  • public subroutine pf_propagate_p(this, q0, q)

    Subroutine propagate computes y_m=expm(Om_m)y_0(expm(Om_m))-1 or (expm(Om_m))y_0 or

    Arguments

    Type IntentOptional AttributesName
    class(pf_imk_t), intent(inout) :: this
    class(pf_encap_t), intent(inout) :: q0
    class(pf_encap_t), intent(inout) :: q

interface

  • public subroutine pf_commutator_p(this, a, b, out, flags)

    Arguments

    Type IntentOptional AttributesName
    class(pf_imk_t), intent(inout) :: this
    class(pf_encap_t), intent(inout) :: a
    class(pf_encap_t), intent(inout) :: b
    class(pf_encap_t), intent(inout) :: out
    integer, intent(in), optional :: flags

Derived Types

type, public, abstract, extends(pf_sweeper_t) :: pf_imk_t

Implicit Munthe-Kaas Runge-Kutta sweeper type, extends abstract sweeper

Components

TypeVisibility AttributesNameInitial
integer, public :: npieces
logical, public :: use_LUq
class(pf_encap_t), public, allocatable:: A(:)
real(kind=pfdp), public, allocatable:: QtilE(:,:)
real(kind=pfdp), public, allocatable:: dtsdc(:)
real(kind=pfdp), public, allocatable:: tsdc(:)
real(kind=pfdp), public, allocatable:: QdiffE(:,:)

qmat-QtilE

real(kind=pfdp), public :: bernoullis(20)
real(kind=pfdp), public :: t0
real(kind=pfdp), public :: dt
integer, public :: qtype
integer, public :: nterms
logical, public :: Lax_pair
logical, public :: use_SDC
logical, public :: debug
logical, public :: mkrk
logical, public :: rk

Type-Bound Procedures

procedure(pf_destroy_p), public :: destroy
procedure, public :: sweep => imk_sweep
procedure, public :: initialize => imk_initialize
procedure, public :: evaluate => imk_evaluate
procedure, public :: integrate => imk_integrate
procedure, public :: residual => imk_residual
procedure, public :: spreadq0 => imk_spreadq0
procedure, public :: evaluate_all => imk_evaluate_all
procedure, public :: imk_destroy
procedure(pf_f_eval_p), public :: f_eval
procedure(pf_dexpinv_p), public :: dexpinv
procedure(pf_propagate_p), public :: propagate
procedure(pf_commutator_p), public :: commutator_p

Subroutines

public subroutine imk_sweep(this, pf, level_index, t0, dt, nsweeps, flags)

Perform nsweep sweeps on level and set qend appropriately.

Arguments

Type IntentOptional AttributesName
class(pf_imk_t), intent(inout) :: this

Inputs

type(pf_pfasst_t), intent(inout), target:: pf

PFASST structure

integer, intent(in) :: level_index

which level to sweep on

real(kind=pfdp), intent(in) :: t0

Time at beginning of time step

real(kind=pfdp), intent(in) :: dt

time step size

integer, intent(in) :: nsweeps

number of sweeps to do

integer, intent(in), optional :: flags

public subroutine rk_step(this, pf, t0, dt)

Arguments

Type IntentOptional AttributesName
class(pf_imk_t), intent(inout) :: this

Inputs

type(pf_pfasst_t), intent(inout), target:: pf

PFASST structure

real(kind=pfdp), intent(in) :: t0

Time at beginning of time step

real(kind=pfdp), intent(in) :: dt

time step size

public subroutine mkrk_step(this, pf, t0, dt)

Arguments

Type IntentOptional AttributesName
class(pf_imk_t), intent(inout) :: this

Inputs

type(pf_pfasst_t), intent(inout), target:: pf

PFASST structure

real(kind=pfdp), intent(in) :: t0

Time at beginning of time step

real(kind=pfdp), intent(in) :: dt

time step size

public subroutine imk_actually_sweep(this, pf, level_index, t0, dt, nsweeps)

Assign level pointer

Read more…

Arguments

Type IntentOptional AttributesName
class(pf_imk_t), intent(inout) :: this

Inputs

type(pf_pfasst_t), intent(inout), target:: pf

PFASST structure

integer, intent(in) :: level_index

which level to sweep on

real(kind=pfdp), intent(in) :: t0

Time at beginning of time step

real(kind=pfdp), intent(in) :: dt

time step size

integer, intent(in) :: nsweeps

number of sweeps to do

public subroutine imk_initialize(this, lev)

Assign explicit approximate quadrature rule Make space for temporary variables

Arguments

Type IntentOptional AttributesName
class(pf_imk_t), intent(inout) :: this
class(pf_level_t), intent(inout) :: lev

public subroutine imk_integrate(this, lev, qSDC, fSDC, dt, fintSDC, flags)

Arguments

Type IntentOptional AttributesName
class(pf_imk_t), intent(inout) :: this
class(pf_level_t), intent(in) :: lev
class(pf_encap_t), intent(in) :: qSDC(:)
class(pf_encap_t), intent(in) :: fSDC(:,:)
real(kind=pfdp), intent(in) :: dt
class(pf_encap_t), intent(inout) :: fintSDC(:)
integer, intent(in), optional :: flags

public subroutine imk_evaluate(this, lev, t, m, flags, step)

Arguments

Type IntentOptional AttributesName
class(pf_imk_t), intent(inout) :: this
class(pf_level_t), intent(inout) :: lev
real(kind=pfdp), intent(in) :: t
integer, intent(in) :: m
integer, intent(in), optional :: flags
integer, intent(in), optional :: step

public subroutine imk_evaluate_all(this, lev, t, flags, step)

Arguments

Type IntentOptional AttributesName
class(pf_imk_t), intent(inout) :: this
class(pf_level_t), intent(inout) :: lev
real(kind=pfdp), intent(in) :: t(:)
integer, intent(in), optional :: flags
integer, intent(in), optional :: step

public subroutine imk_residual(this, lev, dt, flags)

Arguments

Type IntentOptional AttributesName
class(pf_imk_t), intent(inout) :: this
class(pf_level_t), intent(inout) :: lev
real(kind=pfdp), intent(in) :: dt
integer, intent(in), optional :: flags

public subroutine imk_spreadq0(this, lev, t0, flags, step)

Arguments

Type IntentOptional AttributesName
class(pf_imk_t), intent(inout) :: this
class(pf_level_t), intent(inout) :: lev
real(kind=pfdp), intent(in) :: t0
integer, intent(in), optional :: flags
integer, intent(in), optional :: step

public subroutine imk_save(lev)

Save function values so that difference can be computed

Arguments

Type IntentOptional AttributesName
class(pf_level_t), intent(inout) :: lev

Level to save on

public subroutine imk_destroy(this, lev)

Arguments

Type IntentOptional AttributesName
class(pf_imk_t), intent(inout) :: this
class(pf_level_t), intent(inout) :: lev