pf_v_cycle_oc Subroutine

public subroutine pf_v_cycle_oc(pf, iteration, t0, dt, level_index_c, level_index_f, flags)

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

Arguments

Type IntentOptional AttributesName
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

Calls

proc~~pf_v_cycle_oc~~CallsGraph proc~pf_v_cycle_oc pf_v_cycle_oc proc~restrict_time_space_fas restrict_time_space_fas proc~pf_v_cycle_oc->proc~restrict_time_space_fas proc~pf_post pf_post proc~pf_v_cycle_oc->proc~pf_post proc~pf_residual pf_residual proc~pf_v_cycle_oc->proc~pf_residual proc~pf_send pf_send proc~pf_v_cycle_oc->proc~pf_send proc~pf_recv pf_recv proc~pf_v_cycle_oc->proc~pf_recv proc~interpolate_time_space interpolate_time_space proc~pf_v_cycle_oc->proc~interpolate_time_space proc~restrict_sdc restrict_sdc proc~restrict_time_space_fas->proc~restrict_sdc proc~end_timer end_timer proc~restrict_time_space_fas->proc~end_timer proc~call_hooks call_hooks proc~restrict_time_space_fas->proc~call_hooks proc~start_timer start_timer proc~restrict_time_space_fas->proc~start_timer proc~pf_residual->proc~end_timer proc~pf_residual->proc~start_timer end_timer end_timer proc~pf_send->end_timer start_timer start_timer proc~pf_send->start_timer proc~pf_recv->end_timer proc~pf_recv->start_timer proc~interpolate_time_space->proc~end_timer proc~interpolate_time_space->proc~call_hooks proc~pf_apply_mat pf_apply_mat proc~interpolate_time_space->proc~pf_apply_mat proc~interpolate_time_space->proc~start_timer proc~restrict_sdc->proc~pf_apply_mat proc~call_hooks->proc~end_timer proc~call_hooks->proc~start_timer

Called by

proc~~pf_v_cycle_oc~~CalledByGraph proc~pf_v_cycle_oc pf_v_cycle_oc proc~pf_pfasst_block_oc pf_pfasst_block_oc proc~pf_pfasst_block_oc->proc~pf_v_cycle_oc

Contents

Source Code


Source Code

  subroutine pf_v_cycle_oc(pf, iteration, t0, dt, level_index_c,level_index_f, flags)
  ! Execute a V-cycle between levels nfine and ncoarse

    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, which, dir

    which = 1
    if(present(flags)) which = flags
    ! send forward by default, even if sweeping on both components; send backwards if sweeping on p only
    dir = 1
    if(which == 2) dir = 2 !

!     !  For a single level, just get new initial conditions and return
!     if (pf%nlevels == 1) then
!        f_lev_p => pf%levels(1)
!        call pf_recv(pf, f_lev_p, f_lev_p%index*10000+iteration, .true., dir)
!        return
!     end if

    !>  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, dir)
    end do
    
!     !
!     ! down (fine to coarse)
!     !
!     do level_index = pf%nlevels-1, 2, -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, which)
!       call pf_send(pf, f_lev_p, level_index*10000+iteration, .false., dir)
!       call restrict_time_space_fas(pf, t0, dt, level_index, which)
!       call save(c_lev_p, which)
!     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, which)
       call pf_send(pf, f_lev_p, level_index*10000+iteration, .false., dir)
       call restrict_time_space_fas(pf, t0, dt, level_index, which)
       call save(c_lev_p, which)
    end do

!     !
!     ! bottom  (coarsest level)
!     !
!     level_index=1
!     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., dir)
!           call f_lev_p%ulevel%sweeper%sweep(pf, level_index, t0, dt, 1, which)
!           call pf_send(pf, f_lev_p, f_lev_p%index*10000+iteration+j, .false., dir)
!        end do
!     else
! !       if (which == 0) then
! !         call pf_recv(pf, f_lev_p, f_lev_p%index*10000+iteration, .true., dir)
! !         call f_lev_p%ulevel%sweeper%sweep(pf, level_index, t0, dt, f_lev_p%nsweeps, 1)
! !         call pf_send(pf, f_lev_p, level_index*10000+iteration, .false., dir)
! !         call f_lev_p%ulevel%sweeper%sweep(pf, level_index, t0, dt, f_lev_p%nsweeps, 2) ! this interferes with skipping y sweeps: have to check 
! !                                                                                        ! state residual in case of which==1 in sweeper as well
! !       else
!         call pf_recv(pf, f_lev_p, f_lev_p%index*10000+iteration, .true., dir)
!         call f_lev_p%ulevel%sweeper%sweep(pf, level_index, t0, dt, f_lev_p%nsweeps, which)
!         call pf_send(pf, f_lev_p, level_index*10000+iteration, .false., dir)
! !       endif
!     endif
    ! 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., dir)
          call f_lev_p%ulevel%sweeper%sweep(pf, level_index, t0, dt, 1, which)
          call pf_send(pf, f_lev_p, f_lev_p%index*10000+iteration+j, .false., dir)
       end do
    else
       call pf_recv(pf, f_lev_p, f_lev_p%index*10000+iteration, .true., dir)
       call f_lev_p%ulevel%sweeper%sweep(pf, level_index, t0, dt, f_lev_p%nsweeps, which)
       call pf_send(pf, f_lev_p, level_index*10000+iteration, .false., dir)
    endif
    
!     !
!     ! up  (coarse to fine)
!     !
!     do level_index = 2, pf%nlevels
!       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, which)
!       call pf_recv(pf, f_lev_p, level_index*10000+iteration, .false., dir)
! 
!        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
!           if ((which .eq. 0) .or. (which .eq. 1)) call interpolate_q0(pf, f_lev_p, c_lev_p, 1)
!        end if
!        if (pf%rank /= pf%comm%nproc-1) then
!           if (which .eq. 2) call interpolate_qend(pf, f_lev_p, c_lev_p)
!        end if
! 
!        if (level_index < pf%nlevels) then
!           call call_hooks(pf, level_index, PF_PRE_SWEEP)
!           ! compute residual
!           ! do while residual > tol and j < nswps
!           ! assuming residual computed at end of sweep 
!           call f_lev_p%ulevel%sweeper%sweep(pf, level_index, t0, dt, f_lev_p%nsweeps, which)
!        end if
!     end do
    
    ! 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, which)
       call pf_recv(pf, f_lev_p, level_index*10000+iteration, .false., dir)
       
       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
          if ((which .eq. 0) .or. (which .eq. 1)) call interpolate_q0(pf, f_lev_p, c_lev_p, 1)
       end if
       if (pf%rank /= pf%comm%nproc-1) then
          if (which .eq. 2)                       call interpolate_qend(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, which)
       else  !  compute residual for diagnostics since we didn't sweep
          call pf_residual(pf, f_lev_p, dt, which)
       end if
    end do


  end subroutine pf_v_cycle_oc