interpolate_qend Subroutine

public subroutine interpolate_qend(pf, f_lev_ptr, c_lev_ptr)

Subroutine to update the fine terminal condition from coarse increment by spatial interpolation used for adjoint solver

create local workspace restrict fine initial data to coarse get coarse level correction

interpolate correction in space

update fine inital condition

destroy local workspace

Arguments

Type IntentOptional AttributesName
type(pf_pfasst_t), intent(inout) :: pf

main pfasst structure

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

fine level

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

coarse level


Calls

proc~~interpolate_qend~~CallsGraph proc~interpolate_qend interpolate_qend proc~call_hooks call_hooks proc~interpolate_qend->proc~call_hooks proc~start_timer start_timer proc~interpolate_qend->proc~start_timer proc~end_timer end_timer proc~interpolate_qend->proc~end_timer proc~call_hooks->proc~start_timer proc~call_hooks->proc~end_timer

Contents

Source Code


Source Code

  subroutine interpolate_qend(pf, f_lev_ptr, c_lev_ptr)
  
    type(pf_pfasst_t), intent(inout) :: pf          !!  main pfasst structure
    class(pf_level_t),  intent(inout) :: f_lev_ptr  !!  fine level
    class(pf_level_t),  intent(inout) :: c_lev_ptr  !!  coarse level
    
    class(pf_encap_t), allocatable ::    c_delta    !!  coarse correction
    class(pf_encap_t), allocatable ::    f_delta    !!  fine correction

    call call_hooks(pf, f_lev_ptr%index, PF_PRE_INTERP_Q0)
    call start_timer(pf, TINTERPOLATE + f_lev_ptr%index - 1)
    !> create local workspace
    call c_lev_ptr%ulevel%factory%create_single(c_delta, c_lev_ptr%index, c_lev_ptr%shape)
    call f_lev_ptr%ulevel%factory%create_single(f_delta, f_lev_ptr%index, f_lev_ptr%shape)

    call c_delta%setval(0.0_pfdp)
    call f_delta%setval(0.0_pfdp)

    !>  restrict fine initial data to coarse
    call f_lev_ptr%ulevel%restrict(f_lev_ptr, c_lev_ptr, f_lev_ptr%qend, c_delta, pf%state%t0, 2)
    !>  get coarse level correction
    call c_delta%axpy(-1.0_pfdp, c_lev_ptr%qend, 2)    

    !>  interpolate correction in space
    call f_lev_ptr%ulevel%interpolate(f_lev_ptr, c_lev_ptr, f_delta, c_delta, pf%state%t0, 2)

    !> update fine inital condition
    call f_lev_ptr%qend%axpy(-1.0_pfdp, f_delta, 2)

    call end_timer(pf, TINTERPOLATE + f_lev_ptr%index - 1)
    call call_hooks(pf, f_lev_ptr%index, PF_POST_INTERP_Q0)

    !> destroy local workspace
    call c_lev_ptr%ulevel%factory%destroy_single(c_delta, c_lev_ptr%index, c_lev_ptr%shape)
    call f_lev_ptr%ulevel%factory%destroy_single(f_delta, f_lev_ptr%index,  f_lev_ptr%shape)

  end subroutine interpolate_qend