TABLE OF CONTENTS
ABINIT/xmpi_alltoall [ Functions ]
NAME
xmpi_alltoall
FUNCTION
This module contains functions that calls MPI routine, if we compile the code using the MPI CPP flags. xmpi_alltoall is the generic function.
COPYRIGHT
Copyright (C) 2001-2024 ABINIT group (AR,XG) This file is distributed under the terms of the GNU General Public License, see ~ABINIT/COPYING or http://www.gnu.org/copyleft/gpl.txt .
SOURCE
ABINIT/xmpi_alltoall_dp2d [ Functions ]
NAME
xmpi_alltoall_dp2d
FUNCTION
Sends data from all to all processes. Target: double precision two-dimensional arrays.
SOURCE
70 subroutine xmpi_alltoall_dp2d(xval, sendsize, recvbuf, recvsize, comm, ier) 71 72 !Arguments------------------------- 73 real(dp), DEV_CONTARRD intent(in) :: xval(:,:) 74 real(dp), DEV_CONTARRD intent(inout) :: recvbuf(:,:) 75 integer ,intent(in) :: sendsize, recvsize 76 integer ,intent(in) :: comm 77 integer ,intent(out) :: ier 78 79 ! ************************************************************************* 80 81 ier=0 82 #if defined HAVE_MPI 83 if (comm /= MPI_COMM_SELF .and. comm /= MPI_COMM_NULL) then 84 ! allgather xval on all proc. in comm 85 call MPI_ALLTOALL(xval, sendsize, MPI_DOUBLE_PRECISION, recvbuf, & 86 & recvsize, MPI_DOUBLE_PRECISION, comm, ier) 87 else if (comm == MPI_COMM_SELF) then 88 recvbuf=xval 89 end if 90 #else 91 recvbuf=xval 92 #endif 93 94 end subroutine xmpi_alltoall_dp2d
ABINIT/xmpi_alltoall_dp4d [ Functions ]
NAME
xmpi_alltoall_dp4d
FUNCTION
Sends data from all to all processes. Target: double precision four-dimensional arrays.
SOURCE
107 subroutine xmpi_alltoall_dp4d(xval, sendsize, recvbuf, recvsize, comm, ier) 108 109 !Arguments------------------------- 110 real(dp), DEV_CONTARRD intent(in) :: xval(:,:,:,:) 111 real(dp), DEV_CONTARRD intent(inout) :: recvbuf(:,:,:,:) 112 integer ,intent(in) :: sendsize, recvsize 113 integer ,intent(in) :: comm 114 integer ,intent(out) :: ier 115 116 ! ************************************************************************* 117 118 ier=0 119 #if defined HAVE_MPI 120 if (comm /= MPI_COMM_SELF .and. comm /= MPI_COMM_NULL) then 121 ! allgather xval on all proc. in comm 122 call MPI_ALLTOALL(xval, sendsize, MPI_DOUBLE_PRECISION, recvbuf, & 123 & recvsize, MPI_DOUBLE_PRECISION, comm, ier) 124 else if (comm == MPI_COMM_SELF) then 125 recvbuf=xval 126 end if 127 #else 128 recvbuf=xval 129 #endif 130 131 end subroutine xmpi_alltoall_dp4d
ABINIT/xmpi_alltoall_int [ Functions ]
NAME
xmpi_alltoall_int
FUNCTION
Sends data from all to all processes. Target: integer mono dimensional arrays.
SOURCE
31 subroutine xmpi_alltoall_int(xval, sendsize, recvbuf, recvsize, comm, ier) 32 33 !Arguments------------------------- 34 integer, DEV_CONTARRD intent(in) :: xval(:) 35 integer, DEV_CONTARRD intent(inout) :: recvbuf(:) 36 integer ,intent(in) :: sendsize, recvsize 37 integer ,intent(in) :: comm 38 integer ,intent(out) :: ier 39 40 !Local variables------------------- 41 42 ! ************************************************************************* 43 44 ier=0 45 #if defined HAVE_MPI 46 if (comm /= MPI_COMM_SELF .and. comm /= MPI_COMM_NULL) then 47 ! allgather xval on all proc. in comm 48 call MPI_ALLTOALL(xval, sendsize, MPI_INTEGER, recvbuf, & 49 & recvsize, MPI_INTEGER, comm, ier) 50 else if (comm == MPI_COMM_SELF) then 51 recvbuf=xval 52 end if 53 #else 54 recvbuf=xval 55 #endif 56 57 end subroutine xmpi_alltoall_int