TABLE OF CONTENTS
- ABINIT/xmpi_gather
- ABINIT/xmpi_gather_dp
- ABINIT/xmpi_gather_dp2d
- ABINIT/xmpi_gather_dp3d
- ABINIT/xmpi_gather_dp4d
- ABINIT/xmpi_gather_int
- ABINIT/xmpi_gather_int2d
ABINIT/xmpi_gather [ Functions ]
NAME
xmpi_gather
FUNCTION
This module contains functions that calls MPI routine, if we compile the code using the MPI CPP flags. xmpi_gather is the generic function.
COPYRIGHT
Copyright (C) 2001-2024 ABINIT group (MT) 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_gather_dp [ Functions ]
NAME
xmpi_gather_dp
FUNCTION
Gathers data from all tasks and delivers it to all. Target: one-dimensional real arrays.
INPUTS
xval= buffer array sendcont= number of sent elements recvcount= number of received elements root= rank of receiving process comm= MPI communicator
OUTPUT
ier= exit status, a non-zero value meaning there is an error
SIDE EFFECTS
recvbuf= received buffer
SOURCE
138 subroutine xmpi_gather_dp(xval,sendcount,recvbuf,recvcount,root,comm,ier) 139 140 !Arguments------------------------- 141 integer,intent(in) :: sendcount,recvcount 142 real(dp), DEV_CONTARRD intent(in) :: xval(:) 143 real(dp), DEV_CONTARRD intent(inout) :: recvbuf(:) 144 integer,intent(in) :: root,comm 145 integer,intent(out) :: ier 146 147 ! ************************************************************************* 148 149 ier=0 150 #if defined HAVE_MPI 151 if (comm /= MPI_COMM_SELF .and. comm /= MPI_COMM_NULL) then 152 call MPI_gather(xval,sendcount,MPI_DOUBLE_PRECISION,recvbuf,recvcount,MPI_DOUBLE_PRECISION,& 153 & root,comm,ier) 154 else if (comm == MPI_COMM_SELF) then 155 recvbuf=xval 156 end if 157 #else 158 recvbuf=xval 159 #endif 160 end subroutine xmpi_gather_dp
ABINIT/xmpi_gather_dp2d [ Functions ]
NAME
xmpi_gather_dp2d
FUNCTION
Gathers data from all tasks and delivers it to all. Target: two-dimensional real arrays.
INPUTS
xval= buffer array sendcont= number of sent elements recvcount= number of received elements root= rank of receiving process comm= MPI communicator
OUTPUT
ier= exit status, a non-zero value meaning there is an error
SIDE EFFECTS
recvbuf= received buffer
SOURCE
186 subroutine xmpi_gather_dp2d(xval,sendcount,recvbuf,recvcount,root,comm,ier) 187 188 !Arguments------------------------- 189 integer,intent(in) :: sendcount,recvcount 190 real(dp), DEV_CONTARRD intent(in) :: xval(:,:) 191 real(dp), DEV_CONTARRD intent(inout) :: recvbuf(:,:) 192 integer,intent(in) :: root,comm 193 integer,intent(out) :: ier 194 195 ! ************************************************************************* 196 197 ier=0 198 #if defined HAVE_MPI 199 if (comm /= MPI_COMM_SELF .and. comm /= MPI_COMM_NULL) then 200 call MPI_gather(xval,sendcount,MPI_DOUBLE_PRECISION,recvbuf,recvcount,MPI_DOUBLE_PRECISION,& 201 & root,comm,ier) 202 else if (comm == MPI_COMM_SELF) then 203 recvbuf=xval 204 end if 205 #else 206 recvbuf=xval 207 #endif 208 end subroutine xmpi_gather_dp2d
ABINIT/xmpi_gather_dp3d [ Functions ]
NAME
xmpi_gather_dp3d
FUNCTION
Gathers data from all tasks and delivers it to all. Target: three-dimensional real arrays.
INPUTS
xval= buffer array sendcont= number of sent elements recvcount= number of received elements root= rank of receiving process comm= MPI communicator
OUTPUT
ier= exit status, a non-zero value meaning there is an error
SIDE EFFECTS
recvbuf= received buffer
SOURCE
234 subroutine xmpi_gather_dp3d(xval,sendcount,recvbuf,recvcount,root,comm,ier) 235 236 !Arguments------------------------- 237 integer,intent(in) :: sendcount,recvcount 238 real(dp), DEV_CONTARRD intent(in) :: xval(:,:,:) 239 real(dp), DEV_CONTARRD intent(inout) :: recvbuf(:,:,:) 240 integer,intent(in) :: root,comm 241 integer,intent(out) :: ier 242 243 ! ************************************************************************* 244 245 ier=0 246 #if defined HAVE_MPI 247 if (comm /= MPI_COMM_SELF .and. comm /= MPI_COMM_NULL) then 248 call MPI_gather(xval,sendcount,MPI_DOUBLE_PRECISION,recvbuf,recvcount,MPI_DOUBLE_PRECISION,& 249 & root,comm,ier) 250 else if (comm == MPI_COMM_SELF) then 251 recvbuf=xval 252 end if 253 #else 254 recvbuf=xval 255 #endif 256 257 end subroutine xmpi_gather_dp3d
ABINIT/xmpi_gather_dp4d [ Functions ]
NAME
xmpi_gather_dp4d
FUNCTION
Gathers data from all tasks and delivers it to all. Target: four-dimensional real arrays.
INPUTS
xval= buffer array sendcont= number of sent elements recvcount= number of received elements root= rank of receiving process comm= MPI communicator
OUTPUT
ier= exit status, a non-zero value meaning there is an error
SIDE EFFECTS
recvbuf= received buffer
SOURCE
283 subroutine xmpi_gather_dp4d(xval,sendcount,recvbuf,recvcount,root,comm,ier) 284 285 !Arguments------------------------- 286 integer,intent(in) :: sendcount,recvcount 287 real(dp), DEV_CONTARRD intent(in) :: xval(:,:,:,:) 288 real(dp), DEV_CONTARRD intent(inout) :: recvbuf(:,:,:,:) 289 integer,intent(in) :: root,comm 290 integer,intent(out) :: ier 291 292 ! ************************************************************************* 293 294 ier=0 295 #if defined HAVE_MPI 296 if (comm /= MPI_COMM_SELF .and. comm /= MPI_COMM_NULL) then 297 call MPI_gather(xval,sendcount,MPI_DOUBLE_PRECISION,recvbuf,recvcount,MPI_DOUBLE_PRECISION,& 298 & root,comm,ier) 299 else if (comm == MPI_COMM_SELF) then 300 recvbuf=xval 301 end if 302 #else 303 recvbuf=xval 304 #endif 305 306 end subroutine xmpi_gather_dp4d
ABINIT/xmpi_gather_int [ Functions ]
NAME
xmpi_gather_int
FUNCTION
Gathers data from all tasks and delivers it to all. Target: one-dimensional integer arrays.
INPUTS
xval= buffer array sendcont= number of sent elements recvcount= number of received elements root= rank of receiving process comm= MPI communicator
OUTPUT
ier= exit status, a non-zero value meaning there is an error
SIDE EFFECTS
recvbuf= received buffer
SOURCE
44 subroutine xmpi_gather_int(xval,sendcount,recvbuf,recvcount,root,comm,ier) 45 46 !Arguments------------------------- 47 integer,intent(in) :: sendcount,recvcount 48 integer, DEV_CONTARRD intent(in) :: xval(:) 49 integer, DEV_CONTARRD intent(inout) :: recvbuf(:) 50 integer,intent(in) :: root,comm 51 integer,intent(out) :: ier 52 53 ! ************************************************************************* 54 55 ier=0 56 #if defined HAVE_MPI 57 if (comm /= MPI_COMM_SELF .and. comm /= MPI_COMM_NULL) then 58 call MPI_gather(xval,sendcount,MPI_INTEGER,recvbuf,recvcount,MPI_INTEGER,root,comm,ier) 59 else if (comm == MPI_COMM_SELF) then 60 recvbuf=xval 61 end if 62 #else 63 recvbuf=xval 64 #endif 65 end subroutine xmpi_gather_int
ABINIT/xmpi_gather_int2d [ Functions ]
NAME
xmpi_gather_int2d
FUNCTION
Gathers data from all tasks and delivers it to all. Target: two-dimensional integer arrays.
INPUTS
xval= buffer array sendcont= number of sent elements recvcount= number of received elements root= rank of receiving process comm= MPI communicator
OUTPUT
ier= exit status, a non-zero value meaning there is an error
SIDE EFFECTS
recvbuf= received buffer
SOURCE
91 subroutine xmpi_gather_int2d(xval,sendcount,recvbuf,recvcount,root,comm,ier) 92 93 !Arguments------------------------- 94 integer,intent(in) :: sendcount,recvcount 95 integer, DEV_CONTARRD intent(in) :: xval(:,:) 96 integer, DEV_CONTARRD intent(inout) :: recvbuf(:,:) 97 integer,intent(in) :: root,comm 98 integer,intent(out) :: ier 99 100 ! ************************************************************************* 101 102 ier=0 103 #if defined HAVE_MPI 104 if (comm /= MPI_COMM_SELF .and. comm /= MPI_COMM_NULL) then 105 call MPI_gather(xval,sendcount,MPI_INTEGER,recvbuf,recvcount,MPI_INTEGER,root,comm,ier) 106 else if (comm == MPI_COMM_SELF) then 107 recvbuf=xval 108 end if 109 #else 110 recvbuf=xval 111 #endif 112 end subroutine xmpi_gather_int2d