pf_recv Subroutine

public subroutine pf_recv(pf, level, tag, blocking, direction)

Subroutine to recieve the solution from the previous processor

Arguments

Type IntentOptional AttributesName
type(pf_pfasst_t), intent(inout) :: pf
type(pf_level_t), intent(inout) :: level
integer, intent(in) :: tag
logical, intent(in) :: blocking
integer, intent(in), optional :: direction

Calls

proc~~pf_recv~~CallsGraph proc~pf_recv pf_recv start_timer start_timer proc~pf_recv->start_timer end_timer end_timer proc~pf_recv->end_timer

Called by

proc~~pf_recv~~CalledByGraph proc~pf_recv pf_recv proc~pf_predictor_oc pf_predictor_oc proc~pf_predictor_oc->proc~pf_recv proc~pf_v_cycle pf_v_cycle proc~pf_v_cycle->proc~pf_recv proc~pf_predictor pf_predictor proc~pf_predictor->proc~pf_recv proc~pf_v_cycle_oc pf_v_cycle_oc proc~pf_v_cycle_oc->proc~pf_recv 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 pf_recv(pf, level, tag, blocking, direction)
    type(pf_pfasst_t), intent(inout) :: pf
    type(pf_level_t),  intent(inout) :: level
    integer,           intent(in)    :: tag
    logical,           intent(in)    :: blocking
    integer, optional, intent(in)    :: direction
    integer                          :: dir, ierror
    
    dir = 1 ! default: send forward
    if(present(direction)) dir = direction

    ierror = 0
    call start_timer(pf, TRECEIVE + level%index - 1)
    if (pf%debug) print*,  'DEBUG --',pf%rank, 'begin recv, tag=',tag,blocking, "pf%state%pstatus=",pf%state%pstatus
    if (pf%rank /= 0 .and. pf%state%pstatus == PF_STATUS_ITERATING &
                                  .and. dir == 1) then
       call pf%comm%recv(pf, level,tag, blocking, ierror, dir)
       if (ierror .eq. 0) then
          if (present(direction)) then
             call level%q0%unpack(level%recv, 1)
          else 
             call level%q0%unpack(level%recv)
          end if
        end if
    elseif (pf%rank /= pf%comm%nproc-1 .and. pf%state%pstatus == PF_STATUS_ITERATING &
                                     .and. dir == 2) then
       call pf%comm%recv(pf, level, tag, blocking, ierror, dir)
       if (ierror .eq. 0) then 
         if (present(direction)) then
           call level%qend%unpack(level%recv, 2)
         else
           call level%qend%unpack(level%recv)
        end if
  
       end if
    end if
    if (pf%debug) print*,  'DEBUG --',pf%rank, 'end recv, tag=',tag,blocking    
    if(ierror .ne. 0) then
      print *, pf%rank, 'warning: mpi error during receive', ierror
      stop "pf_parallel:pf_recv"     
    end if
    
    call end_timer(pf, TRECEIVE + level%index - 1)
  end subroutine pf_recv