Subroutine to update the fine initial condition from coarse increment by spatial interpolation
create local workspace restrict fine initial data to coarse get coarse level correction interpolate correction in space update fine inital condition destroy local workspace
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
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 |
||
integer, | intent(in), | optional | :: | flags |
subroutine interpolate_q0(pf, f_lev_ptr, c_lev_ptr, flags)
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
integer, optional, intent(in) :: flags !! optional: specify component on which to operate
! here flags more or less is logical, if it is present we operate on component 1
! of the ndarray-type
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%q0, c_delta, pf%state%t0, 1)
!> get coarse level correction
call c_delta%axpy(-1.0_pfdp, c_lev_ptr%q0, 1)
!> interpolate correction in space
call f_lev_ptr%ulevel%interpolate(f_lev_ptr, c_lev_ptr, f_delta, c_delta, pf%state%t0, 1)
!> update fine inital condition
call f_lev_ptr%q0%axpy(-1.0_pfdp, f_delta, 1)
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_q0