pf_mpi_send_status Subroutine

public subroutine pf_mpi_send_status(pf, tag, istatus, ierror, direction)

Uses

  • proc~~pf_mpi_send_status~~UsesGraph proc~pf_mpi_send_status pf_mpi_send_status module~pf_mod_mpi pf_mod_mpi proc~pf_mpi_send_status->module~pf_mod_mpi

Subroutine to send convergence status information

Arguments

Type IntentOptional AttributesName
type(pf_pfasst_t), intent(inout) :: pf

main pfasst structure

integer, intent(in) :: tag

message tag

integer, intent(in) :: istatus

status flag to send

integer, intent(inout) :: ierror

error flag

integer, intent(in), optional :: direction

Calls

proc~~pf_mpi_send_status~~CallsGraph proc~pf_mpi_send_status pf_mpi_send_status mpi_issend mpi_issend proc~pf_mpi_send_status->mpi_issend mpi_wait mpi_wait proc~pf_mpi_send_status->mpi_wait

Contents

Source Code


Source Code

  subroutine pf_mpi_send_status(pf, tag,istatus,ierror, direction)
    use pf_mod_mpi, only: MPI_INTEGER, MPI_STATUS_SIZE, MPI_REQUEST_NULL

    type(pf_pfasst_t), intent(inout) :: pf        !!  main pfasst structure
    integer,           intent(in)    :: tag       !!  message tag
    integer,           intent(in) :: istatus      !!  status flag to send
    integer,           intent(inout) :: ierror    !!  error flag
    integer, optional, intent(in)    :: direction
    integer                          :: dest
    integer    ::  stat(MPI_STATUS_SIZE), dir

    integer :: message
    message = istatus
    
    dir = 1 ! default 1: send forward; set to 2 for send backwards
    if(present(direction)) dir = direction
    
    if( dir == 2) then
       dest = modulo(pf%rank-1, pf%comm%nproc) 
    else
       dest = modulo(pf%rank+1, pf%comm%nproc) 
    end if

    if (pf%comm%statreq /= -66) then
       if (pf%debug) print*, 'DEBUG --',pf%rank, 'waiting in send_status with statreq',pf%comm%statreq
       call mpi_wait(pf%comm%statreq, stat, ierror)
       if (pf%debug) print*, 'DEBUG --',pf%rank, 'done waiting in send_status'
    end if

    call mpi_issend(message, 1, MPI_INTEGER, &
                    dest, tag, pf%comm%comm, pf%comm%statreq, ierror)

  end subroutine pf_mpi_send_status