Execute a V-cycle between levels nfine and ncoarse
Post the nonblocking receives on the all the levels that will be recieving later (for single level this will be skipped) move from fine to coarse doing sweeps
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
type(pf_pfasst_t), | intent(inout), | target | :: | pf | ||
integer, | intent(in) | :: | iteration | |||
real(kind=pfdp), | intent(in) | :: | t0 | |||
real(kind=pfdp), | intent(in) | :: | dt | |||
integer, | intent(in) | :: | level_index_c | Coarsest level of V-cycle |
||
integer, | intent(in) | :: | level_index_f | Finest level of V-cycle |
||
integer, | intent(in), | optional | :: | flags |
subroutine pf_v_cycle(pf, iteration, t0, dt,level_index_c,level_index_f, flags)
type(pf_pfasst_t), intent(inout), target :: pf
real(pfdp), intent(in) :: t0, dt
integer, intent(in) :: iteration
integer, intent(in) :: level_index_c !! Coarsest level of V-cycle
integer, intent(in) :: level_index_f !! Finest level of V-cycle
integer, optional, intent(in) :: flags
type(pf_level_t), pointer :: f_lev_p, c_lev_p
integer :: level_index, j
!> Post the nonblocking receives on the all the levels that will be recieving later
!> (for single level this will be skipped)
do level_index = level_index_c+1, level_index_f
f_lev_p => pf%levels(level_index)
call pf_post(pf, f_lev_p, f_lev_p%index*10000+iteration)
end do
!> move from fine to coarse doing sweeps
do level_index = level_index_f, level_index_c+1, -1
f_lev_p => pf%levels(level_index);
c_lev_p => pf%levels(level_index-1)
call f_lev_p%ulevel%sweeper%sweep(pf, level_index, t0, dt, f_lev_p%nsweeps)
call pf_send(pf, f_lev_p, level_index*10000+iteration, .false.)
call restrict_time_space_fas(pf, t0, dt, level_index)
call save(c_lev_p)
end do
! Do the coarsest level
level_index=level_index_c
f_lev_p => pf%levels(level_index)
if (pf%pipeline_pred) then
do j = 1, f_lev_p%nsweeps
call pf_recv(pf, f_lev_p, f_lev_p%index*10000+iteration+j, .true.)
call f_lev_p%ulevel%sweeper%sweep(pf, level_index, t0, dt, 1)
call pf_send(pf, f_lev_p, f_lev_p%index*10000+iteration+j, .false.)
end do
else
call pf_recv(pf, f_lev_p, f_lev_p%index*10000+iteration, .true.)
call f_lev_p%ulevel%sweeper%sweep(pf, level_index, t0, dt, f_lev_p%nsweeps)
call pf_send(pf, f_lev_p, level_index*10000+iteration, .false.)
endif
! Now move coarse to fine interpolating and sweeping
do level_index = level_index_c+1,level_index_f
f_lev_p => pf%levels(level_index);
c_lev_p => pf%levels(level_index-1)
call interpolate_time_space(pf, t0, dt, level_index, c_lev_p%Finterp)
call pf_recv(pf, f_lev_p, level_index*10000+iteration, .false.)
if (pf%rank /= 0) then
! interpolate increment to q0 -- the fine initial condition
! needs the same increment that Q(1) got, but applied to the
! new fine initial condition
call interpolate_q0(pf, f_lev_p, c_lev_p)
end if
! don't sweep on the finest level since that is only done at beginning
if (level_index < level_index_f) then
call f_lev_p%ulevel%sweeper%sweep(pf, level_index, t0, dt, f_lev_p%nsweeps)
else ! compute residual for diagnostics since we didn't sweep
call pf_residual(pf, f_lev_p, dt)
end if
end do
end subroutine pf_v_cycle