pf_pfasst_block_oc Subroutine

public subroutine pf_pfasst_block_oc(pf, dt, nsteps, predict, flags, step)

Routine to do the pfasst iterations for optimal control problems on one block of processors until completion. Each processor will do either a fixed number of iterations, or iterate until a tolerance is met On calling, it is assumed that the levels are already loaded with the initial guesses

Arguments

Type IntentOptional AttributesName
type(pf_pfasst_t), intent(inout), target:: pf
real(kind=pfdp), intent(in) :: dt
integer, intent(in) :: nsteps
logical, intent(in) :: predict
integer, intent(in), optional :: flags
integer, intent(in), optional :: step

Calls

proc~~pf_pfasst_block_oc~~CallsGraph proc~pf_pfasst_block_oc pf_pfasst_block_oc proc~end_timer end_timer proc~pf_pfasst_block_oc->proc~end_timer proc~pf_check_convergence_oc pf_check_convergence_oc proc~pf_pfasst_block_oc->proc~pf_check_convergence_oc proc~pf_v_cycle_oc pf_v_cycle_oc proc~pf_pfasst_block_oc->proc~pf_v_cycle_oc proc~call_hooks call_hooks proc~pf_pfasst_block_oc->proc~call_hooks proc~pf_predictor_oc pf_predictor_oc proc~pf_pfasst_block_oc->proc~pf_predictor_oc proc~start_timer start_timer proc~pf_pfasst_block_oc->proc~start_timer proc~pf_check_convergence_oc->proc~call_hooks proc~pf_check_residual_oc pf_check_residual_oc proc~pf_check_convergence_oc->proc~pf_check_residual_oc proc~pf_send_status pf_send_status proc~pf_check_convergence_oc->proc~pf_send_status 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~interpolate_time_space interpolate_time_space proc~pf_v_cycle_oc->proc~interpolate_time_space 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~call_hooks->proc~end_timer proc~call_hooks->proc~start_timer proc~pf_predictor_oc->proc~end_timer proc~pf_predictor_oc->proc~call_hooks proc~pf_predictor_oc->proc~start_timer proc~pf_predictor_oc->proc~restrict_time_space_fas proc~pf_predictor_oc->proc~pf_residual proc~pf_predictor_oc->proc~interpolate_time_space proc~pf_predictor_oc->proc~pf_send proc~pf_predictor_oc->proc~pf_recv 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 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 start_timer start_timer proc~pf_send->start_timer end_timer end_timer proc~pf_send->end_timer proc~pf_recv->start_timer proc~pf_recv->end_timer proc~restrict_sdc->proc~pf_apply_mat

Contents

Source Code


Source Code

  subroutine pf_pfasst_block_oc(pf, dt, nsteps, predict, flags, step)
    type(pf_pfasst_t), intent(inout), target :: pf
    real(pfdp),        intent(in)    :: dt
    integer,           intent(in)    :: nsteps 
    logical,           intent(in)    :: predict
    integer, optional, intent(in)    :: flags    !0 (default): sweep on y and p, 1: just y, 2: just p
    integer, optional, intent(in)    :: step
    ! not yet clear how to handle send and receive for forward and backward combined 

    type(pf_level_t), pointer :: fine_lev_p, coarse_lev_p
    integer                   :: k, j, l, which, pred_flags(1), dir, ierror !dir to choose forward or backward send 
    real(pfdp)                :: residual

    logical :: converged, qbroadcast
    logical :: did_post_step_hook

    call start_timer(pf, TTOTAL)

    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
    pred_flags(1) = which
    
    if( present(step) ) then
      pf%state%step    = step
    else
      pf%state%step    = pf%rank
    end if
      
    pf%state%dt      = dt
    pf%state%proc    = pf%rank+1
    pf%state%t0      = pf%state%step * dt
    pf%state%iter    = -1
    pf%state%cycle   = -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
    pf%state%nsteps  = nsteps
    

    residual = -1
    did_post_step_hook = .false.
    
!    call pf%results%initialize(nsteps, pf%niters, pf%comm%nproc, pf%nlevels)
    
!     do k = 1, 666666666 
! 
!        qbroadcast = .false.
! 
!        if (pf%state%status == PF_STATUS_CONVERGED .and. .not. did_post_step_hook) then
!          call call_hooks(pf, -1, PF_POST_STEP)
!          did_post_step_hook = .true.
!          pf%state%itcnt = pf%state%itcnt + pf%state%iter
!          pf%state%mysteps = pf%state%mysteps + 1
!          exit
!        end if
! 
!        ! jump to next block if we've reached the max iteration count
!        if (pf%state%iter >= pf%niters) then
! !           print *, pf%rank, 'pf%state%iter >= pf%niters'
!           if (.not. did_post_step_hook) then
!             call call_hooks(pf, -1, PF_POST_STEP)
!             pf%state%itcnt = pf%state%itcnt + pf%state%iter
!             pf%state%mysteps = pf%state%mysteps + 1
!           end if
!           did_post_step_hook = .false.
! 
!           pf%state%step = pf%state%step + pf%comm%nproc
!           pf%state%t0   = pf%state%step * dt
!           
!           if (pf%state%step >= pf%state%nsteps) exit  ! for optimal control this exit should always happen
!           
!           pf%state%status = PF_STATUS_PREDICTOR
!           !pf%state%block  = pf%state%block + 1
!           residual = -1
!           qbroadcast = .true.
!        end if
! 
!        if (k > 1 .and. qbroadcast) then
!           if (pf%comm%nproc > 1) then
!              stop "broadcast not supported" 
!              !fine_lev_p => pf%levels(pf%nlevels)
!              !call pf%comm%wait(pf, pf%nlevels)
!              !call fine_lev_p%encap%pack(fine_lev_p%send, fine_lev_p%qend)
!              !call pf_broadcast(pf, fine_lev_p%send, fine_lev_p%nvars, pf%comm%nproc-1)
!              !call fine_lev_p%encap%unpack(fine_lev_p%q0,fine_lev_p%send)
!           else
!              stop "we should not be here I guess"
!              ! for sequential optimal control, we need to save the Q(m) values for state solution
!              ! and load them when solving the adjoint
!              ! additionally, state solution is needed for objective, adjoint for gradient
!              
!              !print *, 'copying initial/terminal value'
!              fine_lev_p => pf%levels(pf%nlevels)
!              if ((which .eq. 0) .or. (which .eq. 1)) call fine_lev_p%q0%copy(fine_lev_p%qend, 1)
!              if (which .eq. 2) call fine_lev_p%qend%copy(fine_lev_p%q0, 2)
!           end if
!        end if

!       if (pf%state%status == PF_STATUS_PREDICTOR) then
!         !print *, 'pf%state%status == PF_STATUS_PREDICTOR', pf%state%t0, dt, which
        if (predict) then
          !print *, 'calling predictor'
           call pf_predictor_oc(pf, pf%state%t0, dt, pred_flags)
        else
           pf%state%iter = 0
           pf%state%status  = PF_STATUS_ITERATING
           pf%state%pstatus = PF_STATUS_ITERATING
        end if
!       end if    
        call call_hooks(pf, -1, PF_POST_ITERATION)

!       pf%state%iter  = pf%state%iter + 1
! 
! !       exit! just do predictor
!       
!       call start_timer(pf, TITERATION)
!       call call_hooks(pf, -1, PF_PRE_ITERATION)
!       
!       if (pf%state%status /= PF_STATUS_CONVERGED) then
!           fine_lev_p => pf%levels(pf%nlevels)
!           call fine_lev_p%ulevel%sweeper%sweep(pf, pf%nlevels, pf%state%t0, dt, fine_lev_p%nsweeps, which)
!        end if
!       
!       ! check convergence  (should always be not converged)
!       call pf_check_convergence_oc(pf, k,  residual, converged, dir)
!       
!       if (pf%state%step >= pf%state%nsteps) exit
!       
!       if (.not. converged) then
!         !   non-blocking receive at all but the coarsest level
!         do l = 2, pf%nlevels
!           fine_lev_p => pf%levels(l)
!           call pf_post(pf, fine_lev_p, fine_lev_p%index*10000+k, dir)
!         end do
! 
!         if (pf%state%status /= PF_STATUS_CONVERGED) then
!           fine_lev_p => pf%levels(pf%nlevels)
!           call pf_send(pf, fine_lev_p, fine_lev_p%index*10000+k, .false., dir)
!           if (pf%nlevels > 1) then
!             coarse_lev_p => pf%levels(pf%nlevels-1)
!             call restrict_time_space_fas(pf, pf%state%t0, dt, pf%nlevels, which)
!             call save(coarse_lev_p, which)
!           end if             
!         end if
!         
!         call pf_v_cycle_oc(pf, k, pf%state%t0, dt, which)
!         call call_hooks(pf, -1, PF_POST_ITERATION)
!         call end_timer(pf, TITERATION)
!       end if
!     end do  !  Niter loop
! 
!     pf%state%iter = -1
!     call end_timer(pf, TTOTAL)
           
    k = 1 ! always one block
    do j = 1, pf%niters
      call start_timer(pf, TITERATION)
      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, 1 ,pf%nlevels)
      call pf_v_cycle_oc(pf, j, pf%state%t0, dt, 1, pf%nlevels,  which)

      !  Check for convergence
      call pf_check_convergence_oc(pf, send_tag=1111*k+j, flags=dir)
      if (pf%save_results) pf%results%residuals(pf%state%iter, k, pf%nlevels) = pf%levels(pf%nlevels)%residual

      call call_hooks(pf, -1, PF_POST_ITERATION)
      call end_timer(pf, TITERATION)  
      
      !  If we are converged, exit block
      if (pf%state%status == PF_STATUS_CONVERGED)  exit
    end do  !  Loop over the iteration in this block
    call call_hooks(pf, -1, PF_POST_CONVERGENCE)
    call call_hooks(pf, -1, PF_POST_STEP)
    pf%state%itcnt = pf%state%itcnt + pf%state%iter
    
    if (pf%save_results) call pf%results%dump(pf%results)
    
    call end_timer(pf, TTOTAL)
  end subroutine pf_pfasst_block_oc