pf_send Subroutine

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

Subroutine to send the solution to the next processor

Arguments

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

Calls

proc~~pf_send~~CallsGraph proc~pf_send pf_send start_timer start_timer proc~pf_send->start_timer end_timer end_timer proc~pf_send->end_timer

Called by

proc~~pf_send~~CalledByGraph proc~pf_send pf_send proc~pf_predictor_oc pf_predictor_oc proc~pf_predictor_oc->proc~pf_send proc~pf_v_cycle pf_v_cycle proc~pf_v_cycle->proc~pf_send proc~pf_predictor pf_predictor proc~pf_predictor->proc~pf_send proc~pf_v_cycle_oc pf_v_cycle_oc proc~pf_v_cycle_oc->proc~pf_send 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_send(pf, level, tag, blocking, direction)
    type(pf_pfasst_t), intent(inout) :: pf
    class(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, TSEND + level%index - 1)
    if (pf%debug) print*,  'DEBUG --',pf%rank, 'begin send, tag=',tag,blocking,' pf%state%status =',pf%state%status
    if (pf%rank /= pf%comm%nproc-1  .and. dir == 1 ) then
       call pf%comm%send(pf, level, tag, blocking, ierror, dir)
    elseif (pf%rank /= 0 .and. dir == 2 ) then
       !print *, pf%rank, 'sending backward',tag,blocking,pf%rank
       call pf%comm%send(pf, level, tag, blocking, ierror, dir)
    end if
    if (ierror /= 0) then
      print *, pf%rank, 'warning: error during send', ierror
      stop "pf_parallel:pf_send"
   endif
   if (pf%debug) print*,  'DEBUG --',pf%rank, 'end send, tag=',tag,blocking                  
    call end_timer(pf, TSEND + level%index - 1)
  end subroutine pf_send