Skip to content

Commit

Permalink
Merge pull request #36 from terminationshock/implement_split_type
Browse files Browse the repository at this point in the history
Implemented mpifx_comm_split_type
  • Loading branch information
aradi authored Feb 21, 2022
2 parents b32d4cc + aa15bed commit d8682dc
Show file tree
Hide file tree
Showing 5 changed files with 67 additions and 0 deletions.
51 changes: 51 additions & 0 deletions lib/mpifx_comm.fpp
Original file line number Diff line number Diff line change
Expand Up @@ -21,6 +21,9 @@ module mpifx_comm_module
!> Creates a new communicator by splitting the old one.
procedure :: split => mpifx_comm_split

!> Creates a new communicator by splitting the old one given a split type.
procedure :: split_type => mpifx_comm_split_type

!> Frees the communicator. The communicator should not be used after this.
procedure :: free => mpifx_comm_free

Expand Down Expand Up @@ -111,6 +114,54 @@ contains
end subroutine mpifx_comm_split


!> Creates a new communicator by splitting the old one applying a given split type.
!!
!! \param self Communicator instance.
!! \param splittype Determines which ranks to be grouped together. In MPI 3.0,
!! this can only be MPI_COMM_TYPE_SHARED grouping all MPI ranks together
!! that can share memory (usually on a node).
!! \param rankkey Is used to determine the rank of the process in its new
!! communicator. Processes calling the routine with a higher value will
!! have a higher rank in the new communicator.
!! \param newcomm New communicator for the given process.
!! \param error Optional error code on return.
!!
!! Example:
!!
!! program test_split_type
!! use libmpifx_module
!! implicit none
!!
!! type(mpifx_comm) :: allproc, splitproc
!!
!! call mpifx_init()
!! call allproc%init()
!! call allproc%split_type(MPI_COMM_TYPE_SHARED, allproc%rank, splitproc)
!! write(*, "(2(A,1X,I0,1X))") "ID:", allproc%rank, "SPLIT ID", splitproc%rank
!! call mpifx_finalize()
!!
!! end program test_split_type
!!
!! \see MPI documentation (\c MPI_COMM_SPLIT_TYPE)
!!
subroutine mpifx_comm_split_type(self, splittype, rankkey, newcomm, error)
class(mpifx_comm), intent(inout) :: self
integer, intent(in) :: splittype, rankkey
class(mpifx_comm), intent(out) :: newcomm
integer, intent(out), optional :: error

integer :: error0, newcommid

call mpi_comm_split_type(self%id, splittype, rankkey, MPI_INFO_NULL, newcommid, error0)
call handle_errorflag(error0, "mpi_comm_split_type() in mpifx_comm_split_type()", error)
if (error0 /= 0) then
return
end if
call newcomm%init(newcommid, error)

end subroutine mpifx_comm_split_type


!> Frees the MPI communicator.
!>
!> After this call, the passed communicator should not be used any more.
Expand Down
1 change: 1 addition & 0 deletions lib/mpifx_constants.fpp
Original file line number Diff line number Diff line change
Expand Up @@ -8,6 +8,7 @@ module mpifx_constants_module
public :: MPI_LAND, MPI_BAND, MPI_LOR, MPI_BOR, MPI_LXOR ,MPI_BXOR
public :: MPI_MAXLOC, MPI_MINLOC
public :: MPI_THREAD_SINGLE, MPI_THREAD_FUNNELED, MPI_THREAD_SERIALIZED, MPI_THREAD_MULTIPLE
public :: MPI_COMM_TYPE_SHARED
public :: MPIFX_UNHANDLED_ERROR, MPIFX_ASSERT_FAILED


Expand Down
1 change: 1 addition & 0 deletions test/CMakeLists.txt
Original file line number Diff line number Diff line change
Expand Up @@ -4,6 +4,7 @@ set(targets
test_allreduce
test_bcast
test_comm_split
test_comm_split_type
test_gather
test_gatherv
test_reduce
Expand Down
1 change: 1 addition & 0 deletions test/meson.build
Original file line number Diff line number Diff line change
Expand Up @@ -6,6 +6,7 @@ tests = [
'allreduce',
'bcast',
'comm_split',
'comm_split_type',
'gather',
'gatherv',
'reduce',
Expand Down
13 changes: 13 additions & 0 deletions test/test_comm_split_type.f90
Original file line number Diff line number Diff line change
@@ -0,0 +1,13 @@
program test_split_type
use libmpifx_module
implicit none

type(mpifx_comm) :: allproc, splitproc

call mpifx_init()
call allproc%init()
call allproc%split_type(MPI_COMM_TYPE_SHARED, allproc%rank, splitproc)
write(*, "(2(A,1X,I0,1X))") "ID:", allproc%rank, "SPLIT ID", splitproc%rank
call mpifx_finalize()

end program test_split_type

0 comments on commit d8682dc

Please sign in to comment.