pf_mod_restrict Module

Module to restrict solutions between pfasst levels and create the FAS tau correction


Uses

  • module~~pf_mod_restrict~~UsesGraph module~pf_mod_restrict pf_mod_restrict module~pf_mod_dtype pf_mod_dtype module~pf_mod_restrict->module~pf_mod_dtype module~pf_mod_timer pf_mod_timer module~pf_mod_restrict->module~pf_mod_timer module~pf_mod_hooks pf_mod_hooks module~pf_mod_restrict->module~pf_mod_hooks iso_c_binding iso_c_binding module~pf_mod_dtype->iso_c_binding module~pf_mod_timer->module~pf_mod_dtype module~pf_mod_hooks->module~pf_mod_dtype

Used by

  • module~~pf_mod_restrict~~UsedByGraph module~pf_mod_restrict pf_mod_restrict module~pf_mod_interpolate pf_mod_interpolate module~pf_mod_interpolate->module~pf_mod_restrict module~pf_mod_parallel_oc pf_mod_parallel_oc module~pf_mod_parallel_oc->module~pf_mod_restrict module~pf_mod_parallel_oc->module~pf_mod_interpolate module~pf_mod_parallel pf_mod_parallel module~pf_mod_parallel->module~pf_mod_restrict module~pf_mod_parallel->module~pf_mod_interpolate module~pfasst pfasst module~pfasst->module~pf_mod_parallel

Contents


Subroutines

public subroutine restrict_time_space_fas(pf, t0, dt, level_index, flags)

Restrict (in time and space) fine level to coarse and set coarse level FAS correction.

Read more…

Arguments

Type IntentOptional AttributesName
type(pf_pfasst_t), intent(inout), target:: pf
real(kind=pfdp), intent(in) :: t0

time at beginning of step

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

time step

integer, intent(in) :: level_index

defines which level to restrict

integer, intent(in), optional :: flags

public subroutine restrict_sdc(f_lev_ptr, c_lev_ptr, f_encap_array, c_encap_array, IS_INTEGRAL, f_time, flags)

Restrict (in time and space) f_sol_array to c_sol_array Depending on the flag INTEGRAL, we may be restricting solutions, or integrals of F

Read more…

Arguments

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

pointer to fine level

class(pf_level_t), intent(inout) :: c_lev_ptr

pointer to coarse level

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

array of fine level data to be restricted

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

array of coarse level data to be computed

logical, intent(in) :: IS_INTEGRAL

flag determines if it is integral data being restricted

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

time at the fine nodes

integer, intent(in), optional :: flags

public subroutine pf_apply_mat(dst, a, mat, src, zero, flags)

Apply a matrix (tmat or rmat) to src and add to dst. Mathematically this is dst= dst + amatsrc Where dst and src are vectors, mat is a matrix, and a is a scalar If the optional variable "zero" is provided and is true, then we compute dst= amatsrc

Arguments

Type IntentOptional AttributesName
class(pf_encap_t), intent(inout) :: dst(:)

destination vector

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

scalar

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

matrix

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

src vector

logical, intent(in), optional :: zero

If false, don't zero out the the dst variable before computing

integer, intent(in), optional :: flags

Used for choosing which variable to operate on

Read more…

public subroutine pf_apply_mat_backward(dst, a, mat, src, zero, flags)

Apply a matrix (tmat or rmat) to src and add to dst.

Arguments

Type IntentOptional AttributesName
class(pf_encap_t), intent(inout) :: dst(:)

destination vector

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

scalar

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

matrix

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

src vector

logical, intent(in), optional :: zero

If false, don't zero out the the dst variable before computing

integer, intent(in), optional :: flags

Used for choosing which variable to operate on

Read more…