pf_imexQ_oc_t Derived Type

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

IMEX SDC sweeper type for optimal control, extends abstract sweeper


Inherits

type~~pf_imexq_oc_t~~InheritsGraph type~pf_imexq_oc_t pf_imexQ_oc_t type~pf_sweeper_t pf_sweeper_t type~pf_imexq_oc_t->type~pf_sweeper_t type~pf_encap_t pf_encap_t type~pf_imexq_oc_t->type~pf_encap_t rhs

Contents

Source Code


Components

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

Approximate explcit quadrature rule

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

Approximate implcit 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.

Is there an explicit piece

logical, public :: implicit =.true.

Is 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, deferred :: f_eval

RHS function evaluations

  • subroutine pf_f_eval_p(this, y, t, level_index, f, piece, flags, idx, step) Prototype

    Evaluae f_piece(y), where piece is one or two

    Arguments

    Type IntentOptional AttributesName
    class(pf_imexQ_oc_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

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

procedure(pf_f_comp_p), public, deferred :: f_comp

Implicit solver

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

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

    Arguments

    Type IntentOptional AttributesName
    class(pf_imexQ_oc_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

    integer, intent(in) :: flags

procedure, public :: sweep => imexQ_oc_sweep

Set the generic functions

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

    Assign level pointer

    Read more…

    Arguments

    Type IntentOptional AttributesName
    class(pf_imexQ_oc_t), intent(inout) :: this
    type(pf_pfasst_t), intent(inout), target:: pf

    PFASST structure

    integer, intent(in) :: level_index

    which level this is

    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

procedure, public :: initialize => imexQ_oc_initialize

procedure, public :: evaluate => imexQ_oc_evaluate

  • public subroutine imexQ_oc_evaluate(this, lev, t, m, flags, step)

    Arguments

    Type IntentOptional AttributesName
    class(pf_imexQ_oc_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

procedure, public :: integrate => imexQ_oc_integrate

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

    Arguments

    Type IntentOptional AttributesName
    class(pf_imexQ_oc_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

procedure, public :: residual => imexQ_oc_residual

  • public subroutine imexQ_oc_residual(this, lev, dt, flags)

    Arguments

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

procedure, public :: evaluate_all => imexQ_oc_evaluate_all

  • public subroutine imexQ_oc_evaluate_all(this, lev, t, flags, step)

    Evaluate all function values

    Arguments

    Type IntentOptional AttributesName
    class(pf_imexQ_oc_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

procedure, public :: spreadq0 => imexQ_oc_spreadq0

  • public subroutine imexQ_oc_spreadq0(this, lev, t0, flags, step)

    Arguments

    Type IntentOptional AttributesName
    class(pf_imexQ_oc_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

procedure, public :: destroy => imexQ_oc_destroy

procedure, public :: imexQ_oc_destroy

Source Code

  type, extends(pf_sweeper_t), abstract :: pf_imexQ_oc_t
     real(pfdp), allocatable :: QtilE(:,:)   !!  Approximate explcit quadrature rule
     real(pfdp), allocatable :: QtilI(:,:)   !!  Approximate implcit quadrature rule
     real(pfdp), allocatable :: dtsdc(:)     !!  SDC step sizes
     real(pfdp), allocatable :: QdiffE(:,:)  !!  qmat-QtilE
     real(pfdp), allocatable :: QdiffI(:,:)  !!  qmat-QtilI

     logical                 :: explicit = .true. !!  Is there an explicit piece
     logical                 :: implicit = .true. !!  Is there an implicit piece

     class(pf_encap_t), allocatable :: rhs   !! holds rhs for implicit solve
     
  contains
     procedure(pf_f_eval_p), deferred :: f_eval   !!  RHS function evaluations
     procedure(pf_f_comp_p), deferred :: f_comp   !!  Implicit solver
     !>  Set the generic functions
     procedure :: sweep        => imexQ_oc_sweep
     procedure :: initialize   => imexQ_oc_initialize
     procedure :: evaluate     => imexQ_oc_evaluate
     procedure :: integrate    => imexQ_oc_integrate
     procedure :: residual     => imexQ_oc_residual
     procedure :: evaluate_all => imexQ_oc_evaluate_all
     procedure :: spreadq0     => imexQ_oc_spreadq0
     procedure :: destroy      => imexQ_oc_destroy
     procedure :: imexQ_oc_destroy
  end type pf_imexQ_oc_t