restrict_time_space_fas Subroutine

public subroutine restrict_time_space_fas(pf, t0, dt, level_index, flags)

Restrict (in time and space) fine level to coarse and set coarse level FAS correction.

The coarse function values are re-evaluated after restriction. Note that even if the number of variables and nodes is the same, we should still compute the FAS correction since the function evaluations may be different. create workspaces restrict q's and recompute f's Recompute the functions

Compute FAS correction Clean up

Arguments

Type IntentOptional AttributesName
type(pf_pfasst_t), intent(inout), target:: pf
real(kind=pfdp), intent(in) :: t0

time at beginning of step

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

time step

integer, intent(in) :: level_index

defines which level to restrict

integer, intent(in), optional :: flags

Calls

proc~~restrict_time_space_fas~~CallsGraph proc~restrict_time_space_fas restrict_time_space_fas proc~call_hooks call_hooks proc~restrict_time_space_fas->proc~call_hooks proc~restrict_sdc restrict_sdc proc~restrict_time_space_fas->proc~restrict_sdc proc~start_timer start_timer proc~restrict_time_space_fas->proc~start_timer proc~end_timer end_timer proc~restrict_time_space_fas->proc~end_timer proc~call_hooks->proc~start_timer proc~call_hooks->proc~end_timer proc~pf_apply_mat pf_apply_mat proc~restrict_sdc->proc~pf_apply_mat

Called by

proc~~restrict_time_space_fas~~CalledByGraph proc~restrict_time_space_fas restrict_time_space_fas proc~pf_predictor_oc pf_predictor_oc proc~pf_predictor_oc->proc~restrict_time_space_fas proc~pf_v_cycle pf_v_cycle proc~pf_v_cycle->proc~restrict_time_space_fas proc~pf_predictor pf_predictor proc~pf_predictor->proc~restrict_time_space_fas proc~pf_v_cycle_oc pf_v_cycle_oc proc~pf_v_cycle_oc->proc~restrict_time_space_fas proc~pf_block_run pf_block_run proc~pf_block_run->proc~pf_v_cycle proc~pf_block_run->proc~pf_predictor proc~pf_pfasst_block_oc pf_pfasst_block_oc proc~pf_pfasst_block_oc->proc~pf_predictor_oc proc~pf_pfasst_block_oc->proc~pf_v_cycle_oc proc~pf_pfasst_run pf_pfasst_run proc~pf_pfasst_run->proc~pf_block_run

Contents


Source Code

  subroutine restrict_time_space_fas(pf, t0, dt, level_index, flags)
    !! Restrict (in time and space) fine level to coarse and set coarse level FAS correction.
    !!
    !! The coarse function values are re-evaluated after restriction.
    !! Note that even if the number of variables and nodes is the same,
    !! we should still compute the FAS correction since the function
    !! evaluations may be different.
    type(pf_pfasst_t), intent(inout),target :: pf
    real(pfdp),        intent(in)    :: t0            !!  time at beginning of step
    real(pfdp),        intent(in)    :: dt            !!  time step
    integer,           intent(in)    :: level_index   !! defines which level to restrict
    integer, optional, intent(in)    :: flags    

    !>  Local variables
    class(pf_level_t), pointer :: c_lev_ptr    
    class(pf_level_t), pointer :: f_lev_ptr

    integer    :: m, step

    real(pfdp), allocatable :: c_times(:)  !!  Simulation time at coarse nodes  
    real(pfdp), allocatable :: f_times(:)  !!  Simulation time at fine nodes
    class(pf_encap_t), allocatable :: &
         c_tmp_array(:), &    ! coarse integral of coarse function values
         f_int_array(:), &    ! fine integral of fine function values
         f_int_arrayr(:)      ! coarse integral of restricted fine function values
    
    f_lev_ptr => pf%levels(level_index);
    c_lev_ptr => pf%levels(level_index-1)

    
    step = pf%state%step+1
    
    
    call call_hooks(pf, level_index, PF_PRE_RESTRICT_ALL)
    call start_timer(pf, TRESTRICT + level_index - 1)
    
    !> create workspaces
    call c_lev_ptr%ulevel%factory%create_array(c_tmp_array, c_lev_ptr%nnodes, &
      c_lev_ptr%index,   c_lev_ptr%shape)
    call c_lev_ptr%ulevel%factory%create_array(f_int_arrayr, c_lev_ptr%nnodes, &
      c_lev_ptr%index,   c_lev_ptr%shape)
    call c_lev_ptr%ulevel%factory%create_array(f_int_array, f_lev_ptr%nnodes, &
      f_lev_ptr%index,   f_lev_ptr%shape)
    allocate(c_times(c_lev_ptr%nnodes))
    allocate(f_times(f_lev_ptr%nnodes))

    !> restrict q's and recompute f's
    c_times = t0 + dt*c_lev_ptr%nodes
    f_times = t0 + dt*f_lev_ptr%nodes

    call restrict_sdc(f_lev_ptr, c_lev_ptr, f_lev_ptr%Q, c_lev_ptr%Q, .false., f_times, flags)

    !>  Recompute the functions
    call c_lev_ptr%ulevel%sweeper%evaluate_all(c_lev_ptr, c_times, flags=flags, step=step)


    !>  Compute  FAS correction
    do m = 1, c_lev_ptr%nnodes-1
       call c_lev_ptr%tauQ(m)%setval(0.0_pfdp, flags)
    end do
    if (pf%state%iter >= pf%taui0)  then
       ! compute '0 to node' integral on the coarse level
      call c_lev_ptr%ulevel%sweeper%integrate(c_lev_ptr, c_lev_ptr%Q, &
        c_lev_ptr%F, dt, c_tmp_array, flags)
       ! compute '0 to node' integral on the fine level
      call f_lev_ptr%ulevel%sweeper%integrate(f_lev_ptr, f_lev_ptr%Q, &
        f_lev_ptr%F, dt, f_lev_ptr%I, flags)
       !  put tau in on fine level
       if (allocated(f_lev_ptr%tauQ)) then
          do m = 1, f_lev_ptr%nnodes-1
             call f_lev_ptr%I(m)%axpy(1.0_pfdp, f_lev_ptr%tauQ(m), flags)
          end do
       end if
  
       ! restrict '0 to node' integral on the fine level  in time and space
       call restrict_sdc(f_lev_ptr, c_lev_ptr, f_lev_ptr%I, f_int_arrayr, .true.,f_times, flags)

      ! compute '0 to node' tau correction
       do m = 1, c_lev_ptr%nnodes-1
          call c_lev_ptr%tauQ(m)%axpy(1.0_pfdp, f_int_arrayr(m), flags)
          call c_lev_ptr%tauQ(m)%axpy(-1.0_pfdp, c_tmp_array(m), flags)
       end do
    end if

    call end_timer(pf, TRESTRICT + level_index - 1)
    call call_hooks(pf, level_index, PF_POST_RESTRICT_ALL)

    !>  Clean up
    call c_lev_ptr%ulevel%factory%destroy_array(c_tmp_array, c_lev_ptr%nnodes, &
      c_lev_ptr%index,   c_lev_ptr%shape)
    call c_lev_ptr%ulevel%factory%destroy_array(f_int_arrayr, c_lev_ptr%nnodes, &
      c_lev_ptr%index,  c_lev_ptr%shape)
    call f_lev_ptr%ulevel%factory%destroy_array(f_int_array, f_lev_ptr%nnodes, &
      f_lev_ptr%index,   f_lev_ptr%shape)

    deallocate(c_times)
    deallocate(f_times)
  end subroutine restrict_time_space_fas