pf_predictor Subroutine

public subroutine pf_predictor(pf, t0, dt, flags)

PFASST Predictor. Subroutine to initialize the solution on each processor The goal is to have a solution at each level and each node set to a consistent value When this is called, the value of q0 at the fine level on each processor has been set somehow (see q0_style below)

This can be broken down into four substeps 1. Get the initial condition on the finest level at each node 2. Coarsen the initial condition to each coarser level with tau corrections 3. Do the "Burn in" step on the coarse level to make the coarse values consistent (this is skipped if the fine initial conditions are already consistent) 4. Do some coarse grid sweeps to improve initial solutions on coarsest nodes 5. Interpolating coarse correction back to finer levels sweeping along the way.

There are several parameters or flags that determine how things are done: integer q0_style: can take 3 values 0: Only the q0 at t=0 is valid (default) 1: The q0 at each processor is valid 2: q0 and all nodes at each processor is valid logical PFASST_pred: If true, the burn-in step uses the "PFASST predictor" trick integer nsweeps_burn: Determines how many sweeps are done on the coarse level during burn in integer nsweeps_pred: Determines how many sweeps are done at the coarse level after burn in logical Pipeline_burn: True if coarse sweeps during burn in are pipelined (meaningless if nsweeps_burn>1 on coarse level) logical Pipeline_pred: True if coarse sweeps after burn in are pipelined (meaningless if nsweeps_pred>1 on coarse level) Pipeline variables do nothing if there is only one processor logical RK_pred: If true, the coarse level is initialized with Runge-Kutta instead of the PFASST burn in. We will still do coarse sweeps after and correct finer levels

The user defined flags(:) parameter is used to determine whether we are in a (standard) forward-in-time run (flags(1) == 1) or backward-in-time (for the adjoint) with a given terminal condition qend instead of initial condition q0 (flags(1) == 2). In the latter case, e.g., sweeper%spreadq0 has to do the correct thing (i.e., spread qend instead of q0).

No time communication is performed during the predictor since all procesors can do the work themselves

The iteration count is reset to 0, and the status is reset to ITERATING.

Step 1. Getting the initial condition on the finest level at each processor If we are doing multiple levels, then we need to coarsen to fine level

Step 2: Proceed fine to coarse levels coarsening the fine solution and computing tau correction

Step 3. Do the "Burn in" step on the coarse level to make the coarse values consistent (this is skipped if the fine initial conditions are already consistent) The first processor does nothing, the second does one set of sweeps, the third two, etc Hence, this is skipped completely if nprocs=1 If RK_pred is true, just do some RK_steps

Step 4: Now we have everyone burned in, so do some coarse sweeps

Step 5: Return to fine level sweeping on any level in between coarsest and finest

Arguments

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

PFASST main data structure

real(kind=pfdp), intent(in) :: t0

Initial time of this processor

real(kind=pfdp), intent(in) :: dt

time step

integer, intent(in), optional :: flags(:)

User defined flags


Calls

proc~~pf_predictor~~CallsGraph proc~pf_predictor pf_predictor proc~end_timer end_timer proc~pf_predictor->proc~end_timer proc~restrict_time_space_fas restrict_time_space_fas proc~pf_predictor->proc~restrict_time_space_fas proc~pf_residual pf_residual proc~pf_predictor->proc~pf_residual proc~pf_send pf_send proc~pf_predictor->proc~pf_send proc~interpolate_q0 interpolate_q0 proc~pf_predictor->proc~interpolate_q0 proc~call_hooks call_hooks proc~pf_predictor->proc~call_hooks proc~pf_recv pf_recv proc~pf_predictor->proc~pf_recv proc~start_timer start_timer proc~pf_predictor->proc~start_timer proc~interpolate_time_space interpolate_time_space proc~pf_predictor->proc~interpolate_time_space proc~restrict_time_space_fas->proc~end_timer proc~restrict_time_space_fas->proc~call_hooks proc~restrict_time_space_fas->proc~start_timer proc~restrict_sdc restrict_sdc proc~restrict_time_space_fas->proc~restrict_sdc proc~pf_residual->proc~end_timer proc~pf_residual->proc~start_timer start_timer start_timer proc~pf_send->start_timer end_timer end_timer proc~pf_send->end_timer proc~interpolate_q0->proc~end_timer proc~interpolate_q0->proc~call_hooks proc~interpolate_q0->proc~start_timer proc~call_hooks->proc~end_timer proc~call_hooks->proc~start_timer proc~pf_recv->start_timer proc~pf_recv->end_timer proc~interpolate_time_space->proc~end_timer proc~interpolate_time_space->proc~call_hooks proc~interpolate_time_space->proc~start_timer proc~pf_apply_mat pf_apply_mat proc~interpolate_time_space->proc~pf_apply_mat proc~restrict_sdc->proc~pf_apply_mat

Called by

proc~~pf_predictor~~CalledByGraph proc~pf_predictor pf_predictor proc~pf_block_run pf_block_run proc~pf_block_run->proc~pf_predictor proc~pf_pfasst_run pf_pfasst_run proc~pf_pfasst_run->proc~pf_block_run

Contents

Source Code


Source Code

  subroutine pf_predictor(pf, t0, dt, flags)
    type(pf_pfasst_t), intent(inout), target :: pf     !! PFASST main data structure
    real(pfdp),        intent(in   )         :: t0     !! Initial time of this processor
    real(pfdp),        intent(in   )         :: dt     !! time step
    integer,           intent(in   ), optional :: flags(:)  !!  User defined flags

    class(pf_level_t), pointer :: c_lev_p
    class(pf_level_t), pointer :: f_lev_p     !!
    integer                   :: k               !!  Loop indices
    integer                   :: level_index     !!  Local variable for looping over levels
    real(pfdp)                :: t0k             !!  Initial time at time step k

    call call_hooks(pf, 1, PF_PRE_PREDICTOR)
    call start_timer(pf, TPREDICTOR)

    if (pf%debug) print*, 'DEBUG --', pf%rank, 'beginning predictor'
    !!
    !! Step 1. Getting the  initial condition on the finest level at each processor
    !!         If we are doing multiple levels, then we need to coarsen to fine level
    f_lev_p => pf%levels(pf%nlevels)
    if (pf%q0_style < 2) then  !  Spread q0 to all the nodes
       call f_lev_p%ulevel%sweeper%spreadq0(f_lev_p, t0)
    endif
    !!
    !!  Step 2:   Proceed fine to coarse levels coarsening the fine solution and computing tau correction
    if (pf%debug) print*,  'DEBUG --', pf%rank, 'do coarsen  in predictor'    
    if (pf%nlevels > 1) then  
       do level_index = pf%nlevels, 2, -1
          f_lev_p => pf%levels(level_index);
          c_lev_p => pf%levels(level_index-1)
          call pf_residual(pf, f_lev_p, dt)  
          call f_lev_p%ulevel%restrict(f_lev_p, c_lev_p, f_lev_p%q0, c_lev_p%q0, t0)          
          call restrict_time_space_fas(pf, t0, dt, level_index)  !  Restrict
          call save(c_lev_p)
       end do  !  level_index = pf%nlevels, 2, -1
    else
      level_index = 1
      c_lev_p => pf%levels(1)
    end if

    !!
    !! Step 3. Do the "Burn in" step on the coarse level to make the coarse values consistent
    !!         (this is skipped if the fine initial conditions are already consistent)
    !! The first processor does nothing, the second does one set of sweeps, the third two, etc
    !! Hence, this is skipped completely if nprocs=1
    if (pf%debug) print*,  'DEBUG --', pf%rank, 'do burnin  in predictor'    
    if (pf%q0_style .eq. 0) then  !  The coarse level needs burn in
       !! If RK_pred is true, just do some RK_steps
       if (pf%RK_pred) then  !  Use Runge-Kutta to get the coarse initial data
          !  Get new initial conditions
          call pf_recv(pf, c_lev_p, 100000+pf%rank, .true.)

          !  Do a RK_step
          call c_lev_p%ulevel%sweeper%sweep(pf, level_index, t0, dt, 1)
          !  Send forward
          call pf_send(pf, c_lev_p,  100000+pf%rank+1, .false.)
       else  !  Normal PFASST burn in
          level_index=1
          c_lev_p => pf%levels(level_index)
          do k = 1, pf%rank + 1
             pf%state%iter = -k
             t0k = t0-(pf%rank)*dt + (k-1)*dt   ! Remember t0=pf%rank*dt is the beginning of this time slice so t0-(pf%rank)*dt is 0
                                                ! and we iterate up to the correct time step.
                                                ! for optimal control problem t, t0k has no influence on f_eval, so there this does something else

             ! Get new initial value (skip on first iteration)
             if (k > 1) then
                call c_lev_p%q0%copy(c_lev_p%qend,flags=1)
                ! If we are doing PFASST_pred, we use the old values at nodes, otherwise spread q0
                if (.not. pf%PFASST_pred) then
                   call c_lev_p%ulevel%sweeper%spreadq0(c_lev_p, t0k)
                end if
             end if
             !  Do some sweeps
             call c_lev_p%ulevel%sweeper%sweep(pf, level_index, t0k, dt,pf%nsweeps_burn)
          end do
       endif  !  RK_pred
    end if  ! (q0_style .eq. 0)

    !!
    !! Step 4: Now we have everyone burned in, so do some coarse sweeps
    if (pf%nlevels > 1) then
       if (pf%debug) print*,  'DEBUG --', pf%rank, 'do sweeps  in predictor', 'Pipeline_pred',pf%Pipeline_pred    
       pf%state%pstatus = PF_STATUS_ITERATING
       pf%state%status = PF_STATUS_ITERATING
       if (pf%Pipeline_pred) then
          do k = 1, c_lev_p%nsweeps_pred
             pf%state%iter =-(pf%rank + 1) -k
             
             !  Get new initial conditions
             call pf_recv(pf, c_lev_p, c_lev_p%index*110000+pf%rank+k, .true.)
          
             !  Do a sweep
             call c_lev_p%ulevel%sweeper%sweep(pf, level_index, t0, dt, 1)
             !  Send forward
             call pf_send(pf, c_lev_p,  c_lev_p%index*110000+pf%rank+1+k, .false.)
          end do ! k = 1, c_lev_p%nsweeps_pred-1
       else  !  Don't pipeline
          !  Get new initial conditions
          call pf_recv(pf, c_lev_p, c_lev_p%index*110000+pf%rank, .true.)
          
          !  Do a sweeps
          call c_lev_p%ulevel%sweeper%sweep(pf, level_index, t0, dt, c_lev_p%nsweeps_pred)
          !  Send forward
          call pf_send(pf, c_lev_p,  c_lev_p%index*110000+pf%rank+1, .false.)
       endif  ! (Pipeline_pred .eq. .true) then
    end if
    
    if (pf%debug) print*,  'DEBUG --', pf%rank, 'returning to fine level in predictor'
    !!
    !!  Step 5:  Return to fine level sweeping on any level in between coarsest and finest
    do level_index = 2, pf%nlevels  !  Will do nothing with one level
       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 interpolate_q0(pf, f_lev_p, c_lev_p)

       !  Do a sweep on level unless we are at the finest level
       if (level_index < pf%nlevels) then
          call f_lev_p%ulevel%sweeper%sweep(pf, level_index, t0, dt, f_lev_p%nsweeps_pred)
       end if
    end do

    call end_timer(pf, TPREDICTOR)
    call call_hooks(pf, -1, PF_POST_PREDICTOR)

    pf%state%iter   = 0
    pf%state%status = PF_STATUS_ITERATING
    pf%state%pstatus = PF_STATUS_ITERATING
    if (pf%debug) print*,  'DEBUG --', pf%rank, 'ending predictor'    
  end subroutine pf_predictor