TABLE OF CONTENTS


ABINIT/xmpi_exch_dp1d [ Functions ]

[ Top ] [ Functions ]

NAME

  xmpi_exch_dp1d

FUNCTION

  Sends and receives data.
  Target: double precision one-dimensional arrays.

INPUTS

  mtag= message tag
  nt= vector length
  vsend= send buffer
  sender= node sending the data
  recever= node receiving the data
  comm= MPI communicator

OUTPUT

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

SIDE EFFECTS

  vrecv= receive buffer

SOURCE

145 subroutine xmpi_exch_dp1d(vsend,nt,sender,vrecv,recever,comm,mtag,ier)
146 
147 !Arguments----------------
148  integer,intent(in) :: mtag,nt
149  real(dp), DEV_CONTARRD intent(in) :: vsend(:)
150  real(dp), DEV_CONTARRD intent(inout) :: vrecv(:)
151  integer,intent(in) :: sender,recever,comm
152  integer,intent(out) :: ier
153 
154 !Local variables--------------
155 #if defined HAVE_MPI
156  integer :: statux(MPI_STATUS_SIZE)
157  integer :: me,my_dt,my_op,tag
158  integer(kind=int64) :: ntot
159 #endif
160 
161 ! *************************************************************************
162 
163  ier=0
164 #if defined HAVE_MPI
165  if (sender==recever.or.comm==MPI_COMM_NULL.or.nt==0) return
166  call MPI_COMM_RANK(comm,me,ier)
167  tag = MOD(mtag, xmpi_tag_ub)
168 
169  !The dimension can be greater than a 32bit integer
170  !We use a INT64 to store it. If it is too large, we switch to an
171  !alternate routine because MPI<4 doesnt handle 64 bit counts.
172  ntot=int(nt,kind=int64)
173 
174  if (recever==me) then
175    if (ntot<=xmpi_maxint32_64) then
176      call MPI_RECV(vrecv,nt,MPI_DOUBLE_PRECISION,sender,tag,comm,statux,ier)
177    else
178      call xmpi_largetype_create(ntot,MPI_DOUBLE_PRECISION,my_dt,my_op,MPI_OP_NULL)
179      call MPI_RECV(vrecv,1,my_dt,sender,tag,comm,statux,ier)
180      call xmpi_largetype_free(my_dt,my_op)
181    end if
182  end if
183  if (sender==me) then
184    if (ntot<=xmpi_maxint32_64) then
185      call MPI_SEND(vsend,nt,MPI_DOUBLE_PRECISION,recever,tag,comm,ier)
186    else
187      call xmpi_largetype_create(ntot,MPI_DOUBLE_PRECISION,my_dt,my_op,MPI_OP_NULL)
188      call MPI_SEND(vsend,1,my_dt,recever,tag,comm,ier)
189      call xmpi_largetype_free(my_dt,my_op)
190    end if
191  end if
192 #endif
193 
194 end subroutine xmpi_exch_dp1d

ABINIT/xmpi_exch_dp2d [ Functions ]

[ Top ] [ Functions ]

NAME

  xmpi_exch_dp2d

FUNCTION

  Sends and receives data.
  Target: double precision two-dimensional arrays.

INPUTS

  mtag= message tag
  nt= vector length
  vsend= send buffer
  sender= node sending the data
  recever= node receiving the data
  comm= MPI communicator

OUTPUT

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

SIDE EFFECTS

  vrecv= receive buffer

SOURCE

223 subroutine xmpi_exch_dp2d(vsend,nt,sender,vrecv,recever,comm,mtag,ier)
224 
225 !Arguments----------------
226  integer,intent(in) :: mtag,nt
227  real(dp), DEV_CONTARRD intent(in) :: vsend(:,:)
228  real(dp), DEV_CONTARRD intent(inout) :: vrecv(:,:)
229  integer,intent(in) :: sender,recever,comm
230  integer,intent(out) :: ier
231 
232 !Local variables--------------
233 #if defined HAVE_MPI
234  integer :: statux(MPI_STATUS_SIZE)
235  integer :: me,my_dt,my_op,tag
236  integer(kind=int64) :: ntot
237 #endif
238 
239 ! *************************************************************************
240 
241  ier=0
242 #if defined HAVE_MPI
243  if (sender==recever.or.comm==MPI_COMM_NULL.or.nt==0) return
244  call MPI_COMM_RANK(comm,me,ier)
245  tag = MOD(mtag, xmpi_tag_ub)
246 
247  !The dimension can be greater than a 32bit integer
248  !We use a INT64 to store it. If it is too large, we switch to an
249  !alternate routine because MPI<4 doesnt handle 64 bit counts.
250  ntot=int(nt,kind=int64)
251 
252  if (recever==me) then
253    if (ntot<=xmpi_maxint32_64) then
254      call MPI_RECV(vrecv,nt,MPI_DOUBLE_PRECISION,sender,tag,comm,statux,ier)
255    else
256      call xmpi_largetype_create(ntot,MPI_DOUBLE_PRECISION,my_dt,my_op,MPI_OP_NULL)
257      call MPI_RECV(vrecv,1,my_dt,sender,tag,comm,statux,ier)
258      call xmpi_largetype_free(my_dt,my_op)
259    end if
260  end if
261  if (sender==me) then
262    if (ntot<=xmpi_maxint32_64) then
263      call MPI_SEND(vsend,nt,MPI_DOUBLE_PRECISION,recever,tag,comm,ier)
264    else
265      call xmpi_largetype_create(ntot,MPI_DOUBLE_PRECISION,my_dt,my_op,MPI_OP_NULL)
266      call MPI_SEND(vsend,1,my_dt,recever,tag,comm,ier)
267      call xmpi_largetype_free(my_dt,my_op)
268    end if
269  end if
270 #endif
271 
272 end subroutine xmpi_exch_dp2d

ABINIT/xmpi_exch_dp3d [ Functions ]

[ Top ] [ Functions ]

NAME

  xmpi_exch_dp3d

FUNCTION

  Sends and receives data.
  Target: double precision three-dimensional arrays.

INPUTS

  mtag= message tag
  nt= vector length
  vsend= send buffer
  sender= node sending the data
  recever= node receiving the data
  comm= MPI communicator

OUTPUT

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

SIDE EFFECTS

  vrecv= receive buffer

SOURCE

301 subroutine xmpi_exch_dp3d(vsend,nt,sender,vrecv,recever,comm,mtag,ier)
302 
303 !Arguments----------------
304  integer,intent(in) :: mtag,nt
305  real(dp), DEV_CONTARRD intent(in) :: vsend(:,:,:)
306  real(dp), DEV_CONTARRD intent(inout) :: vrecv(:,:,:)
307  integer,intent(in) :: sender,recever,comm
308  integer,intent(out) :: ier
309 
310 !Local variables--------------
311 #if defined HAVE_MPI
312  integer :: statux(MPI_STATUS_SIZE)
313  integer :: me,my_dt,my_op,tag
314  integer(kind=int64) :: ntot
315 #endif
316 
317 ! *************************************************************************
318 
319  ier=0
320 #if defined HAVE_MPI
321  if (sender==recever.or.comm==MPI_COMM_NULL.or.nt==0) return
322  call MPI_COMM_RANK(comm,me,ier)
323  tag = MOD(mtag, xmpi_tag_ub)
324 
325  !The dimension can be greater than a 32bit integer
326  !We use a INT64 to store it. If it is too large, we switch to an
327  !alternate routine because MPI<4 doesnt handle 64 bit counts.
328  ntot=int(nt,kind=int64)
329 
330  if (recever==me) then
331    if (ntot<=xmpi_maxint32_64) then
332      call MPI_RECV(vrecv,nt,MPI_DOUBLE_PRECISION,sender,tag,comm,statux,ier)
333    else
334      call xmpi_largetype_create(ntot,MPI_DOUBLE_PRECISION,my_dt,my_op,MPI_OP_NULL)
335      call MPI_RECV(vrecv,1,my_dt,sender,tag,comm,statux,ier)
336      call xmpi_largetype_free(my_dt,my_op)
337    end if
338  end if
339  if (sender==me) then
340    if (ntot<=xmpi_maxint32_64) then
341      call MPI_SEND(vsend,nt,MPI_DOUBLE_PRECISION,recever,tag,comm,ier)
342    else
343      call xmpi_largetype_create(ntot,MPI_DOUBLE_PRECISION,my_dt,my_op,MPI_OP_NULL)
344      call MPI_SEND(vsend,1,my_dt,recever,tag,comm,ier)
345      call xmpi_largetype_free(my_dt,my_op)
346    end if
347  end if
348 #endif
349 
350 end subroutine xmpi_exch_dp3d

ABINIT/xmpi_exch_dp4d_tag [ Functions ]

[ Top ] [ Functions ]

NAME

  xmpi_exch_dp4d_tag

FUNCTION

  Sends and receives data.
  Target: double precision four-dimensional arrays.

INPUTS

  mtag= message tag
  nt= vector length
  vsend= send buffer
  sender= node sending the data
  recever= node receiving the data
  comm= MPI communicator

OUTPUT

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

SIDE EFFECTS

  vrecv= receive buffer

SOURCE

379 subroutine xmpi_exch_dp4d(vsend,nt,sender,vrecv,recever,comm,mtag,ier)
380 
381 !Arguments----------------
382  integer,intent(in) :: mtag,nt
383  real(dp), DEV_CONTARRD intent(in) :: vsend(:,:,:,:)
384  real(dp), DEV_CONTARRD intent(inout) :: vrecv(:,:,:,:)
385  integer,intent(in) :: sender,recever,comm
386  integer,intent(out) :: ier
387 
388 !Local variables--------------
389 #if defined HAVE_MPI
390  integer :: statux(MPI_STATUS_SIZE)
391  integer :: me,my_dt,my_op,tag
392  integer(kind=int64) :: ntot
393 #endif
394 
395 ! *************************************************************************
396 
397  ier=0
398 #if defined HAVE_MPI
399  if (sender==recever.or.comm==MPI_COMM_NULL.or.nt==0) return
400  call MPI_COMM_RANK(comm,me,ier)
401  tag = MOD(mtag, xmpi_tag_ub)
402 
403  !The dimension can be greater than a 32bit integer
404  !We use a INT64 to store it. If it is too large, we switch to an
405  !alternate routine because MPI<4 doesnt handle 64 bit counts.
406  ntot=int(nt,kind=int64)
407 
408  if (recever==me) then
409    if (ntot<=xmpi_maxint32_64) then
410      call MPI_RECV(vrecv,nt,MPI_DOUBLE_PRECISION,sender,tag,comm,statux,ier)
411    else
412      call xmpi_largetype_create(ntot,MPI_DOUBLE_PRECISION,my_dt,my_op,MPI_OP_NULL)
413      call MPI_RECV(vrecv,1,my_dt,sender,tag,comm,statux,ier)
414      call xmpi_largetype_free(my_dt,my_op)
415    end if
416  end if
417  if (sender==me) then
418    if (ntot<=xmpi_maxint32_64) then
419      call MPI_SEND(vsend,nt,MPI_DOUBLE_PRECISION,recever,tag,comm,ier)
420    else
421      call xmpi_largetype_create(ntot,MPI_DOUBLE_PRECISION,my_dt,my_op,MPI_OP_NULL)
422      call MPI_SEND(vsend,1,my_dt,recever,tag,comm,ier)
423      call xmpi_largetype_free(my_dt,my_op)
424    end if
425  end if
426 #endif
427 
428 end subroutine xmpi_exch_dp4d

ABINIT/xmpi_exch_dp5d_tag [ Functions ]

[ Top ] [ Functions ]

NAME

  xmpi_exch_dp5d_tag

FUNCTION

  Sends and receives data.
  Target: double precision five-dimensional arrays.

INPUTS

  mtag= message tag
  nt= vector length
  vsend= send buffer
  sender= node sending the data
  recever= node receiving the data
  comm= MPI communicator

OUTPUT

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

SIDE EFFECTS

  vrecv= receive buffer

SOURCE

457 subroutine xmpi_exch_dp5d(vsend,nt,sender,vrecv,recever,comm,mtag,ier)
458 
459 !Arguments----------------
460  integer,intent(in) :: mtag,nt
461  real(dp), DEV_CONTARRD intent(in) :: vsend(:,:,:,:,:)
462  real(dp), DEV_CONTARRD intent(inout) :: vrecv(:,:,:,:,:)
463  integer,intent(in) :: sender,recever,comm
464  integer,intent(out) :: ier
465 
466 !Local variables--------------
467 #if defined HAVE_MPI
468  integer :: statux(MPI_STATUS_SIZE)
469  integer :: me,my_dt,my_op,tag
470  integer(kind=int64) :: ntot
471 #endif
472 
473 ! *************************************************************************
474 
475  ier=0
476 #if defined HAVE_MPI
477  if (sender==recever.or.comm==MPI_COMM_NULL.or.nt==0) return
478  call MPI_COMM_RANK(comm,me,ier)
479  tag = MOD(mtag,xmpi_tag_ub)
480 
481  !The dimension can be greater than a 32bit integer
482  !We use a INT64 to store it. If it is too large, we switch to an
483  !alternate routine because MPI<4 doesnt handle 64 bit counts.
484  ntot=int(nt,kind=int64)
485 
486  if (recever==me) then
487    if (ntot<=xmpi_maxint32_64) then
488      call MPI_RECV(vrecv,nt,MPI_DOUBLE_PRECISION,sender,tag,comm,statux,ier)
489    else
490      call xmpi_largetype_create(ntot,MPI_DOUBLE_PRECISION,my_dt,my_op,MPI_OP_NULL)
491      call MPI_RECV(vrecv,1,my_dt,sender,tag,comm,statux,ier)
492      call xmpi_largetype_free(my_dt,my_op)
493    end if
494  end if
495  if (sender==me) then
496    if (ntot<=xmpi_maxint32_64) then
497      call MPI_SEND(vsend,nt,MPI_DOUBLE_PRECISION,recever,tag,comm,ier)
498    else
499      call xmpi_largetype_create(ntot,MPI_DOUBLE_PRECISION,my_dt,my_op,MPI_OP_NULL)
500      call MPI_SEND(vsend,1,my_dt,recever,tag,comm,ier)
501      call xmpi_largetype_free(my_dt,my_op)
502    end if
503  end if
504 #endif
505 
506 end subroutine xmpi_exch_dp5d

ABINIT/xmpi_exch_dpc1d [ Functions ]

[ Top ] [ Functions ]

NAME

  xmpi_exch_dpc1d

FUNCTION

  Sends and receives data.
  Target: one-dimensional double precision complex arrays.

SOURCE

567 subroutine xmpi_exch_dpc1d(vsend,n1,sender,vrecv,recever,comm,mtag,ier)
568 
569 !Arguments----------------
570  integer,intent(in) :: mtag,n1
571  complex(dpc), DEV_CONTARRD intent(in) :: vsend(:)
572  complex(dpc), DEV_CONTARRD intent(inout) :: vrecv(:)
573  integer,intent(in) :: sender,recever,comm
574  integer,intent(out) :: ier
575 
576 !Local variables--------------
577 #if defined HAVE_MPI
578  integer :: statux(MPI_STATUS_SIZE)
579  integer :: tag,me
580 #endif
581 
582 ! *************************************************************************
583 
584  ier=0
585 #if defined HAVE_MPI
586  if (sender==recever.or.comm==MPI_COMM_NULL.or.n1==0) return
587  call MPI_COMM_RANK(comm,me,ier)
588  tag = MOD(mtag, xmpi_tag_ub)
589 
590  if (recever==me) then
591    call MPI_RECV(vrecv,n1,MPI_DOUBLE_COMPLEX,sender, tag,comm,statux,ier)
592  end if
593  if (sender==me) then
594    call MPI_SEND(vsend,n1,MPI_DOUBLE_COMPLEX,recever,tag,comm,ier)
595  end if
596 #endif
597 
598 end subroutine xmpi_exch_dpc1d

ABINIT/xmpi_exch_dpc2d [ Functions ]

[ Top ] [ Functions ]

NAME

  xmpi_exch_dpc2d

FUNCTION

  Sends and receives data.
  Target: two-dimensional double precision complex arrays.

SOURCE

613 subroutine xmpi_exch_dpc2d(vsend,nt,sender,vrecv,recever,comm,mtag,ier)
614 
615 !Arguments----------------
616  integer,intent(in) :: mtag,nt
617  complex(dpc), DEV_CONTARRD intent(in) :: vsend(:,:)
618  complex(dpc), DEV_CONTARRD intent(inout) :: vrecv(:,:)
619  integer,intent(in) :: sender,recever,comm
620  integer,intent(out) :: ier
621 
622 !Local variables--------------
623 #if defined HAVE_MPI
624  integer :: statux(MPI_STATUS_SIZE)
625  integer :: tag,me
626 #endif
627 ! *************************************************************************
628 
629  ier=0
630 #if defined HAVE_MPI
631  if (sender==recever.or.comm==MPI_COMM_NULL.or.(nt==0)) return
632  call MPI_COMM_RANK(comm,me,ier)
633  tag = MOD(mtag, xmpi_tag_ub)
634 
635  if (recever==me) then
636    call MPI_RECV(vrecv,nt,MPI_DOUBLE_COMPLEX,sender, tag,comm,statux,ier)
637  end if
638  if (sender==me) then
639    call MPI_SEND(vsend,nt,MPI_DOUBLE_COMPLEX,recever,tag,comm,ier)
640  end if
641 #endif
642 
643 end subroutine xmpi_exch_dpc2d

ABINIT/xmpi_exch_int1d [ Functions ]

[ Top ] [ Functions ]

NAME

  xmpi_exch_int1d

FUNCTION

  This module contains functions that calls MPI routine,
  if we compile the code using the MPI CPP flags.
  xmpi_exch is the generic function.

COPYRIGHT

  Copyright (C) 2001-2024 ABINIT group (MB,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 .

 NOTE
  The tag conforms to the MPI request that the tag is lower than 32768,
    by using a modulo call.

SOURCE

23 !--------------------------------------------------------------------
24 
25 subroutine xmpi_exch_int1d(vsend,n1,sender,vrecv,recever,comm,mtag,ier)
26 
27 !Arguments----------------
28  integer,intent(in) :: mtag,n1
29  integer, DEV_CONTARRD intent(in) :: vsend(:)
30  integer, DEV_CONTARRD intent(inout) :: vrecv(:)
31  integer,intent(in) :: sender,recever,comm
32  integer,intent(out) :: ier
33 
34 !Local variables--------------
35 #if defined HAVE_MPI
36  integer :: statux(MPI_STATUS_SIZE)
37  integer :: tag,me
38 #endif
39 
40 ! *************************************************************************
41 
42  ier=0
43 #if defined HAVE_MPI
44  if (sender==recever.or.comm==MPI_COMM_NULL.or.n1==0) return
45  call MPI_COMM_RANK(comm,me,ier)
46  tag = MOD(mtag,xmpi_tag_ub)
47 
48  if (recever==me) then
49    call MPI_RECV(vrecv,n1,MPI_INTEGER,sender,tag,comm,statux,ier)
50  end if
51  if (sender==me) then
52    call MPI_SEND(vsend,n1,MPI_INTEGER,recever,tag,comm,ier)
53  end if
54 #endif
55 
56 end subroutine xmpi_exch_int1d

ABINIT/xmpi_exch_int2d [ Functions ]

[ Top ] [ Functions ]

NAME

  xmpi_exch_int2d

FUNCTION

  Sends and receives data.
  Target: two-dimensional integer arrays.

INPUTS

  mtag= message tag
  nt= vector length
  vsend= send buffer
  sender= node sending the data
  recever= node receiving the data
  comm= MPI communicator

OUTPUT

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

SIDE EFFECTS

  vrecv= receive buffer

SOURCE

 85 subroutine xmpi_exch_int2d(vsend,nt,sender,vrecv,recever,comm,mtag,ier)
 86 
 87 !Arguments----------------
 88  integer,intent(in) :: mtag,nt
 89  integer, DEV_CONTARRD intent(in) :: vsend(:,:)
 90  integer, DEV_CONTARRD intent(inout) :: vrecv(:,:)
 91  integer,intent(in) :: sender,recever,comm
 92  integer,intent(out) :: ier
 93 
 94 !Local variables--------------
 95 #if defined HAVE_MPI
 96  integer :: statux(MPI_STATUS_SIZE)
 97  integer :: tag,me
 98 #endif
 99 
100 ! *************************************************************************
101 
102  ier=0
103 #if defined HAVE_MPI
104  if (sender==recever.or.comm==MPI_COMM_NULL.or.nt==0) return
105  call MPI_COMM_RANK(comm,me,ier)
106  tag = MOD(mtag, xmpi_tag_ub)
107 
108  if (recever==me) then
109    call MPI_RECV(vrecv,nt,MPI_INTEGER,sender,tag,comm,statux,ier)
110  end if
111  if (sender==me) then
112    call MPI_SEND(vsend,nt,MPI_INTEGER,recever,tag,comm,ier)
113  end if
114 #endif
115 
116 end subroutine xmpi_exch_int2d

ABINIT/xmpi_exch_spc1d [ Functions ]

[ Top ] [ Functions ]

NAME

  xmpi_exch_spc1d

FUNCTION

  Sends and receives data.
  Target: one-dimensional single precision complex arrays.

SOURCE

521 subroutine xmpi_exch_spc1d(vsend,n1,sender,vrecv,recever,comm,mtag,ier)
522 
523 !Arguments----------------
524  integer,intent(in) :: mtag,n1
525  complex(spc), DEV_CONTARRD intent(in) :: vsend(:)
526  complex(spc), DEV_CONTARRD intent(inout) :: vrecv(:)
527  integer,intent(in) :: sender,recever,comm
528  integer,intent(out) :: ier
529 
530 !Local variables--------------
531 #if defined HAVE_MPI
532  integer :: statux(MPI_STATUS_SIZE)
533  integer :: tag,me
534 #endif
535 
536 ! *************************************************************************
537 
538  ier=0
539 #if defined HAVE_MPI
540  if (sender==recever.or.comm==MPI_COMM_NULL.or.(n1==0)) return
541  call MPI_COMM_RANK(comm,me,ier)
542  tag = MOD(mtag, xmpi_tag_ub)
543 
544  if (recever==me) then
545    call MPI_RECV(vrecv,n1,MPI_COMPLEX,sender, tag,comm,statux,ier)
546  end if
547  if (sender==me) then
548    call MPI_SEND(vsend,n1,MPI_COMPLEX,recever,tag,comm,ier)
549  end if
550 #endif
551 
552 end subroutine xmpi_exch_spc1d