TABLE OF CONTENTS


ABINIT/xmpi_gather [ Functions ]

[ Top ] [ 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-2018 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 ]

[ Top ] [ 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
  spaceComm= MPI communicator

OUTPUT

  ier= exit status, a non-zero value meaning there is an error

SIDE EFFECTS

  recvbuf= received buffer

SOURCE

154 subroutine xmpi_gather_dp(xval,sendcount,recvbuf,recvcount,root,spaceComm,ier)
155 
156 
157 
158 !This section has been created automatically by the script Abilint (TD).
159 !Do not modify the following lines by hand.
160 #undef ABI_FUNC
161 #define ABI_FUNC 'xmpi_gather_dp'
162 !End of the abilint section
163 
164  implicit none
165 
166 
167 !Arguments-------------------------
168  integer,intent(in) :: sendcount,recvcount
169  real(dp), DEV_CONTARRD intent(in) :: xval(:)
170  real(dp), DEV_CONTARRD intent(inout)   :: recvbuf(:)
171  integer,intent(in) :: root,spaceComm
172  integer,intent(out) :: ier
173 
174 ! *************************************************************************
175 
176  ier=0
177 #if defined HAVE_MPI
178  if (spaceComm /= MPI_COMM_SELF .and. spaceComm /= MPI_COMM_NULL) then
179    call MPI_gather(xval,sendcount,MPI_DOUBLE_PRECISION,recvbuf,recvcount,MPI_DOUBLE_PRECISION,&
180 &   root,spaceComm,ier)
181  else if (spaceComm == MPI_COMM_SELF) then
182    recvbuf=xval
183  end if
184 #else
185  recvbuf=xval
186 #endif
187 end subroutine xmpi_gather_dp

ABINIT/xmpi_gather_dp2d [ Functions ]

[ Top ] [ 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
  spaceComm= MPI communicator

OUTPUT

  ier= exit status, a non-zero value meaning there is an error

SIDE EFFECTS

  recvbuf= received buffer

SOURCE

212 subroutine xmpi_gather_dp2d(xval,sendcount,recvbuf,recvcount,root,spaceComm,ier)
213 
214 !This section has been created automatically by the script Abilint (TD).
215 !Do not modify the following lines by hand.
216 #undef ABI_FUNC
217 #define ABI_FUNC 'xmpi_gather_dp2d'
218 !End of the abilint section
219 
220  implicit none
221 
222 !Arguments-------------------------
223  integer,intent(in) :: sendcount,recvcount
224  real(dp), DEV_CONTARRD intent(in) :: xval(:,:)
225  real(dp), DEV_CONTARRD intent(inout) :: recvbuf(:,:)
226  integer,intent(in) :: root,spaceComm
227  integer,intent(out)   :: ier
228 
229 ! *************************************************************************
230 
231  ier=0
232 #if defined HAVE_MPI
233  if (spaceComm /= MPI_COMM_SELF .and. spaceComm /= MPI_COMM_NULL) then
234    call MPI_gather(xval,sendcount,MPI_DOUBLE_PRECISION,recvbuf,recvcount,MPI_DOUBLE_PRECISION,&
235 &   root,spaceComm,ier)
236  else if (spaceComm == MPI_COMM_SELF) then
237    recvbuf=xval
238  end if
239 #else
240  recvbuf=xval
241 #endif
242 end subroutine xmpi_gather_dp2d

ABINIT/xmpi_gather_dp3d [ Functions ]

[ Top ] [ 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
  spaceComm= MPI communicator

OUTPUT

  ier= exit status, a non-zero value meaning there is an error

SIDE EFFECTS

  recvbuf= received buffer

SOURCE

267 subroutine xmpi_gather_dp3d(xval,sendcount,recvbuf,recvcount,root,spaceComm,ier)
268 
269 
270 
271 !This section has been created automatically by the script Abilint (TD).
272 !Do not modify the following lines by hand.
273 #undef ABI_FUNC
274 #define ABI_FUNC 'xmpi_gather_dp3d'
275 !End of the abilint section
276 
277  implicit none
278 
279 !Arguments-------------------------
280  integer,intent(in) :: sendcount,recvcount
281  real(dp), DEV_CONTARRD intent(in) :: xval(:,:,:)
282  real(dp), DEV_CONTARRD intent(inout) :: recvbuf(:,:,:)
283  integer,intent(in) :: root,spaceComm
284  integer,intent(out) :: ier
285 
286 ! *************************************************************************
287 
288  ier=0
289 #if defined HAVE_MPI
290  if (spaceComm /= MPI_COMM_SELF .and. spaceComm /= MPI_COMM_NULL) then
291    call MPI_gather(xval,sendcount,MPI_DOUBLE_PRECISION,recvbuf,recvcount,MPI_DOUBLE_PRECISION,&
292 &   root,spaceComm,ier)
293  else if (spaceComm == MPI_COMM_SELF) then
294    recvbuf=xval
295  end if
296 #else
297  recvbuf=xval
298 #endif
299 
300 end subroutine xmpi_gather_dp3d

ABINIT/xmpi_gather_dp4d [ Functions ]

[ Top ] [ 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
  spaceComm= MPI communicator

OUTPUT

  ier= exit status, a non-zero value meaning there is an error

SIDE EFFECTS

  recvbuf= received buffer

SOURCE

325 subroutine xmpi_gather_dp4d(xval,sendcount,recvbuf,recvcount,root,spaceComm,ier)
326 
327 
328 !This section has been created automatically by the script Abilint (TD).
329 !Do not modify the following lines by hand.
330 #undef ABI_FUNC
331 #define ABI_FUNC 'xmpi_gather_dp4d'
332 !End of the abilint section
333 
334  implicit none
335 
336 
337 !Arguments-------------------------
338  integer,intent(in) :: sendcount,recvcount
339  real(dp), DEV_CONTARRD intent(in) :: xval(:,:,:,:)
340  real(dp), DEV_CONTARRD intent(inout) :: recvbuf(:,:,:,:)
341  integer,intent(in) :: root,spaceComm
342  integer,intent(out) :: ier
343 
344 ! *************************************************************************
345 
346  ier=0
347 #if defined HAVE_MPI
348  if (spaceComm /= MPI_COMM_SELF .and. spaceComm /= MPI_COMM_NULL) then
349    call MPI_gather(xval,sendcount,MPI_DOUBLE_PRECISION,recvbuf,recvcount,MPI_DOUBLE_PRECISION,&
350 &   root,spaceComm,ier)
351  else if (spaceComm == MPI_COMM_SELF) then
352    recvbuf=xval
353  end if
354 #else
355  recvbuf=xval
356 #endif
357 
358 end subroutine xmpi_gather_dp4d

ABINIT/xmpi_gather_int [ Functions ]

[ Top ] [ 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
  spaceComm= 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,spaceComm,ier)
45 
46 !This section has been created automatically by the script Abilint (TD).
47 !Do not modify the following lines by hand.
48 #undef ABI_FUNC
49 #define ABI_FUNC 'xmpi_gather_int'
50 !End of the abilint section
51 
52  implicit none
53 
54 !Arguments-------------------------
55  integer,intent(in) :: sendcount,recvcount
56  integer, DEV_CONTARRD intent(in) :: xval(:)
57  integer, DEV_CONTARRD intent(inout) :: recvbuf(:)
58  integer,intent(in) :: root,spaceComm
59  integer,intent(out) :: ier
60 
61 ! *************************************************************************
62 
63  ier=0
64 #if defined HAVE_MPI
65  if (spaceComm /= MPI_COMM_SELF .and. spaceComm /= MPI_COMM_NULL) then
66    call MPI_gather(xval,sendcount,MPI_INTEGER,recvbuf,recvcount,MPI_INTEGER,root,spaceComm,ier)
67  else if (spaceComm == MPI_COMM_SELF) then
68    recvbuf=xval
69  end if
70 #else
71  recvbuf=xval
72 #endif
73 end subroutine xmpi_gather_int

ABINIT/xmpi_gather_int2d [ Functions ]

[ Top ] [ 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
  spaceComm= MPI communicator

OUTPUT

  ier= exit status, a non-zero value meaning there is an error

SIDE EFFECTS

  recvbuf= received buffer

SOURCE

 98 subroutine xmpi_gather_int2d(xval,sendcount,recvbuf,recvcount,root,spaceComm,ier)
 99 
100 
101 
102 !This section has been created automatically by the script Abilint (TD).
103 !Do not modify the following lines by hand.
104 #undef ABI_FUNC
105 #define ABI_FUNC 'xmpi_gather_int2d'
106 !End of the abilint section
107 
108  implicit none
109 
110 !Arguments-------------------------
111  integer,intent(in) :: sendcount,recvcount
112  integer, DEV_CONTARRD intent(in) :: xval(:,:)
113  integer, DEV_CONTARRD intent(inout) :: recvbuf(:,:)
114  integer,intent(in) :: root,spaceComm
115  integer,intent(out) :: ier
116 
117 ! *************************************************************************
118 
119  ier=0
120 #if defined HAVE_MPI
121  if (spaceComm /= MPI_COMM_SELF .and. spaceComm /= MPI_COMM_NULL) then
122    call MPI_gather(xval,sendcount,MPI_INTEGER,recvbuf,recvcount,MPI_INTEGER,root,spaceComm,ier)
123  else if (spaceComm == MPI_COMM_SELF) then
124    recvbuf=xval
125  end if
126 #else
127  recvbuf=xval
128 #endif
129 end subroutine xmpi_gather_int2d