TABLE OF CONTENTS
ABINIT/xmpi_min_dp [ Functions ]
NAME
xmpi_min_dp
FUNCTION
Combines values from all processes and distribute the result back to all processes. Target: scalar double precisions, in-place version
INPUTS
comm= MPI communicator
OUTPUT
ier= exit status, a non-zero value meaning there is an error
SIDE EFFECTS
xval= buffer array
SOURCE
179 subroutine xmpi_min_dp(xval, comm, ier) 180 181 !Arguments------------------------- 182 real(dp),intent(inout) :: xval 183 integer ,intent(in) :: comm 184 integer ,intent(out) :: ier 185 186 !Local variables------------------- 187 #if defined HAVE_MPI 188 real(dp) :: arr_xval(1), arr_xmin(1) 189 #endif 190 191 ! ************************************************************************* 192 193 ier=0 194 #if defined HAVE_MPI 195 if (comm /= MPI_COMM_SELF .and. comm /= MPI_COMM_NULL) then 196 arr_xval(1) = xval 197 call MPI_ALLREDUCE(arr_xval,arr_xmin,1,MPI_DOUBLE_PRECISION,MPI_MIN,comm,ier) 198 xval = arr_xmin(1) 199 end if 200 #endif 201 202 end subroutine xmpi_min_dp
ABINIT/xmpi_min_dpv [ Functions ]
NAME
xmpi_min_dpv
FUNCTION
Combines values from all processes and distribute the result back to all processes. Target: scalar double precisions.
INPUTS
comm= MPI communicator
OUTPUT
ier= exit status, a non-zero value meaning there is an error
SIDE EFFECTS
xval= buffer array xmin= number of elements in send buffer
PARENTS
CHILDREN
mpi_allreduce
SOURCE
81 subroutine xmpi_min_dpv(xval,xmin,comm,ier) 82 83 !Arguments------------------------- 84 real(dp),intent(in) :: xval 85 real(dp),intent(out) :: xmin 86 integer ,intent(in) :: comm 87 integer ,intent(out) :: ier 88 89 !Local variables------------------- 90 #if defined HAVE_MPI 91 real(dp) :: arr_xval(1),arr_xmin(1) 92 #endif 93 94 ! ************************************************************************* 95 96 ier=0 97 #if defined HAVE_MPI 98 if (comm /= MPI_COMM_SELF .and. comm /= MPI_COMM_NULL) then 99 arr_xval(1) = xval 100 call MPI_ALLREDUCE(arr_xval,arr_xmin,1,MPI_DOUBLE_PRECISION,MPI_MIN,comm,ier) 101 xmin = arr_xmin(1) 102 else 103 #endif 104 xmin=xval 105 #if defined HAVE_MPI 106 end if 107 #endif 108 109 end subroutine xmpi_min_dpv
ABINIT/xmpi_min_int1d [ Functions ]
NAME
xmpi_min_int1d
FUNCTION
Combines values from all processes and distribute the result back to all processes. Target: 1d int array, in-place version
INPUTS
comm= MPI communicator
OUTPUT
ier= exit status, a non-zero value meaning there is an error
SIDE EFFECTS
xval= buffer array
SOURCE
133 subroutine xmpi_min_int1d(xval, comm, ier) 134 135 !Arguments------------------------- 136 integer,intent(inout) :: xval(:) 137 integer,intent(in) :: comm 138 integer,intent(out) :: ier 139 140 !Local variables------------------- 141 #if defined HAVE_MPI 142 integer :: arr_xmin(size(xval)) 143 #endif 144 145 ! ************************************************************************* 146 147 ier=0 148 #if defined HAVE_MPI 149 if (comm /= MPI_COMM_SELF .and. comm /= MPI_COMM_NULL) then 150 call MPI_ALLREDUCE(xval, arr_xmin, size(xval), MPI_INTEGER, MPI_MIN, comm, ier) 151 xval = arr_xmin 152 end if 153 #endif 154 155 end subroutine xmpi_min_int1d
ABINIT/xmpi_min_intv [ Functions ]
NAME
xmpi_min_intv
FUNCTION
This module contains functions that calls MPI routine, if we compile the code using the MPI CPP flags. xmpi_min is the generic function.
COPYRIGHT
Copyright (C) 2001-2022 ABINIT group (AR,XG,MB) This file is distributed under the terms of the GNU General Public License, see ~ABINIT/COPYING or http://www.gnu.org/copyleft/gpl.txt .
PARENTS
CHILDREN
mpi_allreduce
SOURCE
24 subroutine xmpi_min_intv(xval,xmin,comm,ier) 25 26 !Arguments------------------------- 27 integer ,intent(in) :: xval 28 integer ,intent(out) :: xmin 29 integer ,intent(in) :: comm 30 integer ,intent(out) :: ier 31 32 !Local variables------------------- 33 #if defined HAVE_MPI 34 integer :: arr_xval(1),arr_xmin(1) 35 #endif 36 37 ! ************************************************************************* 38 39 ier=0 40 #if defined HAVE_MPI 41 if (comm /= MPI_COMM_SELF .and. comm /= MPI_COMM_NULL) then 42 arr_xval(1) = xval 43 call MPI_ALLREDUCE(arr_xval,arr_xmin,1,MPI_INTEGER,MPI_MIN,comm,ier) 44 xmin = arr_xmin(1) 45 else 46 #endif 47 xmin=xval 48 #if defined HAVE_MPI 49 end if 50 #endif 51 52 end subroutine xmpi_min_intv