PFASST controller for block mode When starting a new block, broadcast new initial conditions to all procs For initial block, this is done when initial conditions are set Reset some flags Pack away your last solution Everyone resets their q0 Just stick qend in q0 Update the step and t0 variables for new block Call the predictor to get an initial guess on all levels and all processors
Start the loops over SDC sweeps
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
type(pf_pfasst_t), | intent(inout), | target | :: | pf | ||
class(pf_encap_t), | intent(in) | :: | q0 | |||
real(kind=pfdp), | intent(in) | :: | dt | |||
integer, | intent(in) | :: | nsteps | |||
class(pf_encap_t), | intent(inout), | optional | :: | qend | ||
integer, | intent(in), | optional | :: | flags(:) |
subroutine pf_block_run(pf, q0, dt, nsteps, qend,flags)
type(pf_pfasst_t), intent(inout), target :: pf
class(pf_encap_t), intent(in ) :: q0
real(pfdp), intent(in ) :: dt
integer, intent(in ) :: nsteps
class(pf_encap_t), intent(inout), optional :: qend
integer, intent(in ), optional :: flags(:)
class(pf_level_t), pointer :: lev_p !! pointer to the one level we are operating on
integer :: j, k
integer :: nblocks !! The number of blocks of steps to do
integer :: nproc !! The number of processors being used
integer :: level_index_c !! Coarsest leve in V-cycle
call start_timer(pf, TTOTAL)
pf%state%dt = dt
pf%state%proc = pf%rank+1
pf%state%step = pf%rank
pf%state%t0 = pf%state%step * dt
! pointer to finest level to start
lev_p => pf%levels(pf%nlevels)
! Stick the initial condition into q0 (will happen on all processors)
call lev_p%q0%copy(q0, flags=1)
nproc = pf%comm%nproc
nblocks = nsteps/nproc
! Decide what the coarsest level in the V-cycle is
level_index_c=1
if (.not. pf%Vcycle) level_index_c=pf%nlevels
do k = 1, nblocks ! Loop over blocks of time steps
! print *,'Starting step=',pf%state%step,' block k=',k
! Each block will consist of
! 1. predictor
! 2. Vcycle until max iterations, or tolerances met
! 3. Move solution to next block
! Reset some flags
!> When starting a new block, broadcast new initial conditions to all procs
!> For initial block, this is done when initial conditions are set
!> Reset some flags
pf%state%iter = -1
pf%state%itcnt = 0
pf%state%mysteps = 0
pf%state%status = PF_STATUS_PREDICTOR
pf%state%pstatus = PF_STATUS_PREDICTOR
pf%comm%statreq = -66
if (k > 1) then
if (nproc > 1) then
call lev_p%qend%pack(lev_p%send) !! Pack away your last solution
call pf_broadcast(pf, lev_p%send, lev_p%mpibuflen, pf%comm%nproc-1)
call lev_p%q0%unpack(lev_p%send) !! Everyone resets their q0
else
call lev_p%q0%copy(lev_p%qend, flags=1) !! Just stick qend in q0
end if
!> Update the step and t0 variables for new block
pf%state%step = pf%state%step + pf%comm%nproc
pf%state%t0 = pf%state%step * dt
end if
!> Call the predictor to get an initial guess on all levels and all processors
call pf_predictor(pf, pf%state%t0, dt, flags)
!> Start the loops over SDC sweeps
pf%state%iter = 0
call call_hooks(pf, -1, PF_POST_ITERATION)
call start_timer(pf, TITERATION)
do j = 1, pf%niters
call call_hooks(pf, -1, PF_PRE_ITERATION)
pf%state%iter = j
! Do a v_cycle
call pf_v_cycle(pf, k, pf%state%t0, dt,level_index_c,pf%nlevels)
! Check for convergence
call pf_check_convergence_block(pf, send_tag=1111*k+j)
if (pf%save_results) pf%results%residuals(pf%state%iter, k, lev_p%index) = lev_p%residual
! print *,pf%rank, ' post res'
call call_hooks(pf, -1, PF_POST_ITERATION)
! If we are converged, exit block
if (pf%state%status == PF_STATUS_CONVERGED) exit
end do ! Loop over the iteration in this bloc
call call_hooks(pf, -1, PF_POST_CONVERGENCE)
call end_timer(pf, TITERATION)
end do ! Loop over the blocks
call end_timer(pf, TTOTAL)
! Grab the last solution for return (if wanted)
if (present(qend)) then
call qend%copy(lev_p%qend, flags=1)
end if
end subroutine pf_block_run