restrict_sdc Subroutine

public subroutine restrict_sdc(f_lev_ptr, c_lev_ptr, f_encap_array, c_encap_array, IS_INTEGRAL, f_time, flags)

Restrict (in time and space) f_sol_array to c_sol_array Depending on the flag INTEGRAL, we may be restricting solutions, or integrals of F

do the restriction

Arguments

Type IntentOptional AttributesName
class(pf_level_t), intent(inout) :: f_lev_ptr

pointer to fine level

class(pf_level_t), intent(inout) :: c_lev_ptr

pointer to coarse level

class(pf_encap_t), intent(inout) :: f_encap_array(:)

array of fine level data to be restricted

class(pf_encap_t), intent(inout) :: c_encap_array(:)

array of coarse level data to be computed

logical, intent(in) :: IS_INTEGRAL

flag determines if it is integral data being restricted

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

time at the fine nodes

integer, intent(in), optional :: flags

Calls

proc~~restrict_sdc~~CallsGraph proc~restrict_sdc restrict_sdc proc~pf_apply_mat pf_apply_mat proc~restrict_sdc->proc~pf_apply_mat

Called by

proc~~restrict_sdc~~CalledByGraph proc~restrict_sdc restrict_sdc proc~restrict_time_space_fas restrict_time_space_fas proc~restrict_time_space_fas->proc~restrict_sdc 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


Source Code

  subroutine restrict_sdc(f_lev_ptr, c_lev_ptr, f_encap_array, c_encap_array, IS_INTEGRAL,f_time, flags)

    !! Restrict (in time and space) f_sol_array  to c_sol_array
    !! Depending on the flag INTEGRAL, we may be restricting solutions, or integrals of F
    
    class(pf_level_t),  intent(inout) :: f_lev_ptr   !!   pointer to fine level
    class(pf_level_t),  intent(inout) :: c_lev_ptr   !!   pointer to coarse level
    class(pf_encap_t),  intent(inout) :: f_encap_array(:)   !! array of fine level data to be restricted
    class(pf_encap_t),  intent(inout) :: c_encap_array(:)   !! array of coarse level data to be computed
    logical,            intent(in)    :: IS_INTEGRAL       !! flag determines if it is integral data being restricted
    real(pfdp),         intent(in) :: f_time(:)             !! time at the fine nodes
    integer, optional, intent(in)    :: flags    

    class(pf_encap_t), allocatable :: f_encap_array_c(:)  !!  fine solution restricted in space only
    integer :: m
    integer :: f_nnodes


    f_nnodes = f_lev_ptr%nnodes

    !!  do the restriction
    if (IS_INTEGRAL) then   ! Restriction of integrals
       call c_lev_ptr%ulevel%factory%create_array(f_encap_array_c, f_nnodes-1, c_lev_ptr%index, c_lev_ptr%shape)
       !  spatial restriction
       do m = 1, f_nnodes-1
          call f_lev_ptr%ulevel%restrict(f_lev_ptr, c_lev_ptr, f_encap_array(m), f_encap_array_c(m), f_time(m), flags)
       end do

       ! temporal restriction
       ! when restricting '0 to node' integral terms, skip the first entry since it is zero
       if (present(flags)) then
          if ((flags .eq. 0) .or. (flags .eq. 1)) &
            call pf_apply_mat(c_encap_array, 1.0_pfdp, f_lev_ptr%rmat(2:,2:), f_encap_array_c, .true., 1)
          if ((flags .eq. 0) .or. (flags .eq. 2)) &
            call pf_apply_mat_backward(c_encap_array, 1.0_pfdp, f_lev_ptr%rmat(2:,2:), f_encap_array_c, .true., 2)
       else
          call pf_apply_mat(c_encap_array, 1.0_pfdp, f_lev_ptr%rmat(2:,2:), f_encap_array_c, .true.)
       end if
       call c_lev_ptr%ulevel%factory%destroy_array(f_encap_array_c, f_nnodes-1, c_lev_ptr%index, c_lev_ptr%shape)
    else
       call c_lev_ptr%ulevel%factory%create_array(f_encap_array_c, f_nnodes, c_lev_ptr%index, c_lev_ptr%shape)
       !  spatial restriction
       do m = 1, f_nnodes
          call f_lev_ptr%ulevel%restrict(f_lev_ptr, c_lev_ptr, f_encap_array(m), f_encap_array_c(m), f_time(m), flags)
       end do! temporal restriction
       if (present(flags)) then
          if ((flags .eq. 0) .or. (flags .eq. 1)) &
            call pf_apply_mat(c_encap_array, 1.0_pfdp, f_lev_ptr%rmat, f_encap_array_c, .true., 1)
          if ((flags .eq. 0) .or. (flags .eq. 2)) &
            call pf_apply_mat_backward(c_encap_array, 1.0_pfdp, f_lev_ptr%rmat, f_encap_array_c, .true., 2)
        else
           call pf_apply_mat(c_encap_array, 1.0_pfdp, f_lev_ptr%rmat, f_encap_array_c, .true.)
        end if
       call c_lev_ptr%ulevel%factory%destroy_array(f_encap_array_c, f_nnodes, c_lev_ptr%index, c_lev_ptr%shape)
    end if

  end subroutine restrict_sdc