On 01/24/2011 11:28 PM, Harald Anlauf wrote: > Hi, > > MPI_Allreduce works for me with MPI_INTEGER8 for all OpenMPI > versions up to 1.4.3. However, with OpenMPI 1.5.1 I get a > failure at runtime: > > [proton:23642] *** An error occurred in MPI_Allreduce: the reduction > operation MPI_SUM is not defined on the MPI_INTEGER8 datatype > [proton:23642] *** on communicator MPI_COMM_WORLD > [proton:23642] *** MPI_ERR_OP: invalid reduce operation > [proton:23642] *** MPI_ERRORS_ARE_FATAL (your MPI job will now abort)
Since I got no reply yet, I have attached an enhanced test case. With openmpi-1.5.1 and np=1, also tested with gfortran: Real kind, digits: 8 53 Integer kind, bits: 8 64 Default Integer : 4 32 Sum[real(8)]: 1.0000000000000000 2.0000000000000000 3.0000000000000000 Sum[integer(4)]: 1 2 3 [proton:16920] *** An error occurred in MPI_Allreduce: the reduction operation MPI_SUM is not defined on the MPI_INTEGER8 datatype [proton:16920] *** on communicator MPI_COMM_WORLD [proton:16920] *** MPI_ERR_OP: invalid reduce operation [proton:16920] *** MPI_ERRORS_ARE_FATAL (your MPI job will now abort) With openmpi-1.4.3: Real kind, digits: 8 53 Integer kind, bits: 8 64 Default Integer : 4 32 Sum[real(8)]: 1.0000000000000000 2.0000000000000000 3.0000000000000000 Sum[integer(4)]: 1 2 3 Sum[integer(8)]: 1 2 3 That's clearly a regression. Harald
program test use mpi implicit none integer, parameter :: i8 = selected_int_kind (15) integer, parameter :: r8 = selected_real_kind (15,90) integer, parameter :: N = 3 integer :: i4i(N), i4s(N) integer(i8) :: i8i(N), i8s(N) real(r8) :: r8i(N), r8s(N) integer :: ierr, nproc, myrank, i i4i = (/ (i, i=1,N) /) i8i = (/ (i, i=1,N) /) r8i = (/ (i, i=1,N) /) call MPI_INIT (ierr) call MPI_COMM_SIZE (MPI_COMM_WORLD, nproc, ierr) call MPI_COMM_RANK (MPI_COMM_WORLD, myrank, ierr) if (myrank == 0) then print *, "Real kind, digits:", r8, digits (1._r8) print *, "Integer kind, bits:", i8, bit_size (1_i8) print *, "Default Integer :", kind (1), bit_size (1) end if call MPI_ALLREDUCE (r8i, r8s, N, MPI_REAL8, MPI_SUM, MPI_COMM_WORLD, ierr) if (myrank == 0) print *, "Sum[real(8)]:", r8s call MPI_ALLREDUCE (i4i, i4s, N, MPI_INTEGER4, MPI_SUM, MPI_COMM_WORLD, ierr) if (myrank == 0) print *, "Sum[integer(4)]:", i4s call MPI_ALLREDUCE (i8i, i8s, N, MPI_INTEGER8, MPI_SUM, MPI_COMM_WORLD, ierr) if (myrank == 0) print *, "Sum[integer(8)]:", i8s call MPI_FINALIZE (ierr) end program test