TABLE OF CONTENTS


ABINIT/xmpi_isum [ Functions ]

[ Top ] [ Functions ]

NAME

  xmpi_isum

FUNCTION

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

COPYRIGHT

  Copyright (C) 2001-2018 ABINIT group (MG)
  This file is distributed under the terms of the
  GNU General Public License, see ~ABINIT/COPYING
  or http://www.gnu.org/copyleft/gpl.txt .

NOTES

    MPI_IALLREDUCE(SENDBUF, RECVBUF, COUNT, DATATYPE, OP, COMM, REQUEST, IERROR)
                   <type>    SENDBUF(*), RECVBUF(*)
                   INTEGER    COUNT, DATATYPE, OP, COMM, REQUEST, IERROR

SOURCE


ABINIT/xmpi_isum_int0d [ Functions ]

[ Top ] [ Functions ]

NAME

  xmpi_isum_int0d

FUNCTION

  Combines values from all processes and distribute
  the result back to all processes.
  Target: scalar integers.
  Non-blocking version.

INPUTS

OUTPUT

SOURCE

42 subroutine xmpi_isum_int0d(xval, xsum, comm, request, ierr)
43 
44 !This section has been created automatically by the script Abilint (TD).
45 !Do not modify the following lines by hand.
46 #undef ABI_FUNC
47 #define ABI_FUNC 'xmpi_alltoall_dp4d'
48 !End of the abilint section
49 
50  implicit none
51 
52 !Arguments-------------------------
53  integer, intent(in) :: xval
54  integer ABI_ASYNC, intent(out) :: xsum
55  integer,intent(in) :: comm
56  integer,intent(out) :: ierr, request
57 
58 !Local variables-------------------
59  integer :: itmp
60 
61 ! *************************************************************************
62 
63 #ifdef HAVE_MPI_IALLREDUCE
64  if (comm /= MPI_COMM_SELF .and. comm /= MPI_COMM_NULL) then
65    call MPI_IALLREDUCE(xval, xsum, 1, MPI_INTEGER, MPI_SUM, comm, request,ierr)
66    !write(std_out,*)"After MPI_IALLREDUCE"
67    return
68  end if
69  xsum = xval; request = xmpi_request_null
70  return
71 #endif
72 
73  ! Call the blocking version and return null request.
74  !write(std_out,*)"will block here and return fake request"
75  itmp = xval
76  call xmpi_sum(itmp, comm, ierr)
77  xsum = itmp; request = xmpi_request_null
78 
79 end subroutine xmpi_isum_int0d