pf_mod_imexQ Module

IMEX Sweeper Module Module of the the derived sweeper class for doing IMEX sweeps for an equation of the form The piece is treated explicitly and implicitl Afer this sweeper is initialized (usually in main), the logical flags can be changed if desired

 explicit:  Make false if there is no explicit piece

 implicit:  Make false if there is no implicit piece

The user needs to supply the feval and fcomp routines for a given example


Uses

  • module~~pf_mod_imexq~~UsesGraph module~pf_mod_imexq pf_mod_imexQ module~pf_mod_utils pf_mod_utils module~pf_mod_imexq->module~pf_mod_utils module~pf_mod_dtype pf_mod_dtype module~pf_mod_imexq->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

Used by

  • module~~pf_mod_imexq~~UsedByGraph module~pf_mod_imexq pf_mod_imexQ module~pfasst pfasst module~pfasst->module~pf_mod_imexq

Contents


Interfaces

interface

  • public subroutine pf_f_eval_p(this, y, t, level_index, f, piece)

    This is the interface for the routine to compute the RHS function values Evaluate f_piece(y), where piece is one or two Evaluate f_piece(y), where piece is one or two

    Arguments

    Type IntentOptional AttributesName
    class(pf_imexQ_t), intent(inout) :: this
    class(pf_encap_t), intent(in) :: y

    Argument for evaluation

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

    Time at evaluation

    integer, intent(in) :: level_index

    Level index

    class(pf_encap_t), intent(inout) :: f

    RHS function value

    integer, intent(in) :: piece

    Which piece to evaluate

interface

  • public subroutine pf_f_comp_p(this, y, t, dtq, rhs, level_index, f, piece)

    Solve the equation y - dtq*f_2(y) =rhs

    Arguments

    Type IntentOptional AttributesName
    class(pf_imexQ_t), intent(inout) :: this
    class(pf_encap_t), intent(inout) :: y

    Solution of implicit solve

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

    Time of solve

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

    dt*quadrature weight

    class(pf_encap_t), intent(in) :: rhs

    RHS for solve

    integer, intent(in) :: level_index

    Level index

    class(pf_encap_t), intent(inout) :: f

    f_2 of solution y

    integer, intent(in) :: piece

    Which piece to evaluate


Derived Types

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

IMEX SDC sweeper type, extends abstract sweeper

Components

TypeVisibility AttributesNameInitial
integer, public :: npieces
logical, public :: use_LUq
real(kind=pfdp), public, allocatable:: QtilE(:,:)

Approximate explicit quadrature rule

real(kind=pfdp), public, allocatable:: QtilI(:,:)

Approximate implicit quadrature rule

real(kind=pfdp), public, allocatable:: dtsdc(:)

SDC step sizes

real(kind=pfdp), public, allocatable:: QdiffE(:,:)

qmat-QtilE

real(kind=pfdp), public, allocatable:: QdiffI(:,:)

qmat-QtilI

logical, public :: explicit =.true.

True if there is an explicit piece

logical, public :: implicit =.true.

True if there an implicit piece

class(pf_encap_t), public, allocatable:: rhs

holds rhs for implicit solve

Type-Bound Procedures

procedure(pf_f_eval_p), public :: f_eval

RHS function evaluations

procedure(pf_f_comp_p), public :: f_comp

Implicit solver

procedure, public :: sweep => imexQ_sweep

Set the generic functions

procedure, public :: initialize => imexQ_initialize
procedure, public :: evaluate => imexQ_evaluate
procedure, public :: integrate => imexQ_integrate
procedure, public :: residual => imexQ_residual
procedure, public :: spreadq0 => imexQ_spreadq0
procedure, public :: evaluate_all => imexQ_evaluate_all
procedure, public :: destroy => imexQ_destroy
procedure, public :: imexQ_destroy

Subroutines

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

Perform nsweep SDC sweeps on level level_index and set qend appropriately. Assign level pointer

Read more…

Arguments

Type IntentOptional AttributesName
class(pf_imexQ_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 imexQ_initialize(this, lev)

Subroutine to initialize matrices and space for sweeper Array of substep sizes

Read more…

Arguments

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

Current level

public subroutine imexQ_destroy(this, lev)

Subroutine to deallocate sweeper

Arguments

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

Current level

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

Subroutine to compute Picard integral of function values

Arguments

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

Current level

class(pf_encap_t), intent(in) :: qSDC(:)

Solution values

class(pf_encap_t), intent(in) :: fSDC(:,:)

RHS Function values

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

Time step

class(pf_encap_t), intent(inout) :: fintSDC(:)

Integral from t_n to t_m

integer, intent(in), optional :: flags

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

Subroutine to compute Residual

Arguments

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

Current level

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

Time step

integer, intent(in), optional :: flags

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

Arguments

Type IntentOptional AttributesName
class(pf_imexQ_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 imexQ_evaluate(this, lev, t, m, flags, step)

Subroutine to evaluate function value at node m

Arguments

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

Current level

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

Time at which to evaluate

integer, intent(in) :: m

Node at which to evaluate

integer, intent(in), optional :: flags
integer, intent(in), optional :: step

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

Subroutine to evaluate the function values at all nodes

Arguments

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

Current level

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

Array of times at each node

integer, intent(in), optional :: flags
integer, intent(in), optional :: step