pf_mod_magnus_picard Module

This module implements fully implicit Magnus method using explicit Picard sweeping


Uses

  • module~~pf_mod_magnus_picard~~UsesGraph module~pf_mod_magnus_picard pf_mod_magnus_picard module~pf_mod_utils pf_mod_utils module~pf_mod_magnus_picard->module~pf_mod_utils module~pf_mod_dtype pf_mod_dtype module~pf_mod_magnus_picard->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)

    Arguments

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

interface

  • public subroutine pf_compute_single_commutators_p(this, f)

    Arguments

    Type IntentOptional AttributesName
    class(pf_magpicard_t), intent(inout) :: this
    class(pf_encap_t), intent(inout) :: f(:,:)

interface

  • public subroutine pf_compute_omega_p(this, omega, integrals, f, nodes, qmat, dt, this_node, coefs)

    Arguments

    Type IntentOptional AttributesName
    class(pf_magpicard_t), intent(inout) :: this
    class(pf_encap_t), intent(inout) :: omega
    class(pf_encap_t), intent(inout) :: integrals(:)
    class(pf_encap_t), intent(inout) :: f(:,:)
    real(kind=pfdp), intent(in) :: nodes(:)
    real(kind=pfdp), intent(in) :: qmat(:,:)
    real(kind=pfdp), intent(in) :: dt
    integer, intent(in) :: this_node
    real(kind=pfdp), intent(in) :: coefs(:,:)

interface

  • public subroutine pf_propagate_solution_p(this, sol_t0, sol_tn, omega, level)

    Arguments

    Type IntentOptional AttributesName
    class(pf_magpicard_t), intent(inout) :: this
    class(pf_encap_t), intent(inout) :: sol_t0
    class(pf_encap_t), intent(inout) :: sol_tn
    class(pf_encap_t), intent(inout) :: omega
    integer, intent(in) :: level

interface

  • public subroutine pf_destroy_magpicard_p(this, Lev)

    Arguments

    Type IntentOptional AttributesName
    class(pf_magpicard_t), intent(inout) :: this
    class(pf_level_t), intent(inout) :: Lev

Derived Types

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

Magnus Picard sweeper type, extends abstract sweeper

Components

TypeVisibility AttributesNameInitial
integer, public :: npieces
logical, public :: use_LUq
real(kind=pfdp), public, allocatable:: dtsdc(:)
integer, public :: magnus_order
integer, public :: qtype
real(kind=pfdp), public :: dt
real(kind=pfdp), public :: commutator_coefs(9,3,4)
complex(kind=pfdp), public, allocatable:: commutators(:,:,:)
class(pf_encap_t), public, allocatable:: omega(:)
class(pf_encap_t), public, allocatable:: time_ev_op(:)

Type-Bound Procedures

procedure, public :: sweep => magpicard_sweep
procedure, public :: initialize => magpicard_initialize
procedure, public :: evaluate => magpicard_evaluate
procedure, public :: integrate => magpicard_integrate
procedure, public :: residual => magpicard_residual
procedure, public :: spreadq0 => magpicard_spreadq0
procedure, public :: evaluate_all => magpicard_evaluate_all
procedure(pf_f_eval_p), public :: f_eval
procedure(pf_compute_single_commutators_p), public :: compute_single_commutators
procedure(pf_compute_omega_p), public :: compute_omega
procedure(pf_propagate_solution_p), public :: propagate_solution
procedure(pf_destroy_magpicard_p), public :: destroy
procedure, public :: magpicard_destroy

Subroutines

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

this loop not OMP'd because the deferred procs are OMP'd

Arguments

Type IntentOptional AttributesName
class(pf_magpicard_t), intent(inout) :: this
type(pf_pfasst_t), intent(inout), target:: pf
integer, intent(in) :: level_index
real(kind=pfdp), intent(in) :: t0
real(kind=pfdp), intent(in) :: dt
integer, intent(in) :: nsweeps
integer, intent(in), optional :: flags

public subroutine magpicard_initialize(this, lev)

Arguments

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

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

Compute SDC integral fintSDC = \int_{t_n}^{t_m} fSDC dt

Arguments

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

Arguments

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

Arguments

Type IntentOptional AttributesName
class(pf_magpicard_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 magpicard_residual(this, lev, dt, flags)

Arguments

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

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

Arguments

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

Arguments

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

public subroutine get_commutator_coefs(qtype, nnodes, dt, coefs)

Arguments

Type IntentOptional AttributesName
integer, intent(in) :: qtype
integer, intent(in) :: nnodes
real(kind=pfdp), intent(in) :: dt
real(kind=pfdp), intent(inout) :: coefs(:,:,:)