pf_mpi_send Subroutine

public subroutine pf_mpi_send(pf, level, tag, blocking, ierror, direction)

Uses

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

Subroutine to send solutions

Arguments

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

main pfasst structure

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

level to send from

integer, intent(in) :: tag

message tag

logical, intent(in) :: blocking

true if send is blocking

integer, intent(inout) :: ierror

error flag

integer, intent(in), optional :: direction

Calls

proc~~pf_mpi_send~~CallsGraph proc~pf_mpi_send pf_mpi_send mpi_send mpi_send proc~pf_mpi_send->mpi_send mpi_isend mpi_isend proc~pf_mpi_send->mpi_isend mpi_wait mpi_wait proc~pf_mpi_send->mpi_wait

Contents

Source Code


Source Code

  subroutine pf_mpi_send(pf, level, tag, blocking,ierror, direction)
    use pf_mod_mpi, only:  MPI_STATUS_SIZE

    type(pf_pfasst_t), intent(inout) :: pf       !!  main pfasst structure
    class(pf_level_t), intent(inout) :: level    !!  level to send from
    integer,           intent(in   ) :: tag      !!  message tag
    logical,           intent(in   ) :: blocking !!  true if send is blocking
    integer,           intent(inout) :: ierror   !!  error flag
    integer, optional, intent(in)    :: direction
    integer                          :: dest, dir
    integer ::  stat(MPI_STATUS_SIZE)

    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 (blocking) then
      if(dir == 2) then
          call level%q0%pack(level%send, 2)
       else
          if(present(direction)) then
            call level%qend%pack(level%send, 1)
          else 
            call level%qend%pack(level%send)
          end if
       end if       
       call mpi_send(level%send, level%mpibuflen, myMPI_Datatype, &
                     dest, tag, pf%comm%comm, stat, ierror)
    else
       call mpi_wait(pf%comm%sendreq(level%index), stat, ierror)
!        call level%qend%pack(level%send)
       if(dir == 2) then
          call level%q0%pack(level%send, 2)
       else
          if(present(direction)) then
             call level%qend%pack(level%send, 1)
          else
            call level%qend%pack(level%send)
          end if
       end if
       call mpi_isend(level%send, level%mpibuflen, myMPI_Datatype, &
                      dest, tag, pf%comm%comm, pf%comm%sendreq(level%index), ierror)
    end if
  end subroutine pf_mpi_send