open-mpi / ompi

Open MPI main development repository
https://www.open-mpi.org
Other
2.12k stars 857 forks source link

opal_unpack_general_function might upack too much data #2535

Closed ggouaillardet closed 7 years ago

ggouaillardet commented 7 years ago

@bosilca can you please take care of that ?

the following program works fine if both task run on homogeneous nodes, but fail otherwise.

i strongly suspect the root cause is opal_unpack_general_function that does not exit when the unpacked data reaches *max_data

#include <assert.h>
#include <mpi.h>

/* simple bug reproducer
 * task 0 MPI_Recv 2 MPI_INT but task 1 only MPI_Send 1 MPI_INT
 * ok on homogeneous node(s), ko on heterogeneous nodes */
int main(int argc, char *argv[]) {
    int b[2];
    int rank;
    MPI_Status status;
    int count;
    int err;
    MPI_Init(&argc, &argv);
    MPI_Comm_rank(MPI_COMM_WORLD, &rank);
    if (0 == rank) {
        b[0] = -1;
        b[1] = -1;
        err = MPI_Recv(b, 2, MPI_INT, 1, 0, MPI_COMM_WORLD, &status);
        assert(MPI_SUCCESS == err);
        err = MPI_Get_elements(&status, MPI_INT, &count);
        assert(MPI_SUCCESS == err);
        assert(1 == count);
        assert(0 == b[0]);
        assert(-1 == b[1]); // currently fails on heterogeneous nodes
    } else if (1 == rank) {
        b[0] = 0;
        b[1] = 2;
        MPI_Send(b, 1, MPI_INT, 0, 0, MPI_COMM_WORLD);
    }

    MPI_Finalize();
    return 0;
}
bosilca commented 7 years ago

runs on open-mpi/ompi@f95f11d without issues.

ggouaillardet commented 7 years ago

i did not make it very clear ... the issue only occurs on heterogeneous clusters (e.g. one x86_64 node and one sparcv9 node, using for example btl/tcp) i guess you might not have any hardware to test this, so i suggested you review the code of opal_unpack_general_function and try to figure out what goes wrong in this case (less bytes received than expected). meanwhile, i ll try to write an opal only reproducer

ggouaillardet commented 7 years ago

@bosilca here is the standalone test you need to configure with --enable-heterogeneous and then

cd test/datatype
make unpack_hetero
./unpack_hetero
commit a236ecc39f75bcb6d6623e9fb3a9dabed3eedae3
Author: Gilles Gouaillardet <gilles@rist.or.jp>
Date:   Thu Dec 8 10:23:42 2016 +0900

    test/datatype: add test for short unpack on heteregeneous cluster

    Signed-off-by: Gilles Gouaillardet <gilles@rist.or.jp>

diff --git a/test/datatype/Makefile.am b/test/datatype/Makefile.am
index cf54cb3..ba65499 100644
--- a/test/datatype/Makefile.am
+++ b/test/datatype/Makefile.am
@@ -18,7 +18,7 @@ if PROJECT_OMPI
     MPI_TESTS = checksum position position_noncontig ddt_test ddt_raw unpack_ooo ddt_pack external32
     MPI_CHECKS = to_self ddt_pack_hetero
 endif
-TESTS = opal_datatype_test $(MPI_TESTS)
+TESTS = opal_datatype_test unpack_hetero $(MPI_TESTS)

 check_PROGRAMS = $(TESTS) $(MPI_CHECKS)

@@ -83,5 +83,10 @@ external32_LDADD = \
         $(top_builddir)/ompi/lib@OMPI_LIBMPI_NAME@.la \
         $(top_builddir)/opal/lib@OPAL_LIB_PREFIX@open-pal.la

+unpack_hetero_SOURCES = unpack_hetero.c
+unpack_hetero_LDFLAGS = $(OMPI_PKG_CONFIG_LDFLAGS)
+unpack_hetero_LDADD = \
+        $(top_builddir)/opal/lib@OPAL_LIB_PREFIX@open-pal.la
+
 distclean:
        rm -rf *.dSYM .deps .libs *.log *.o *.trs $(check_PROGRAMS) Makefile
diff --git a/test/datatype/unpack_hetero.c b/test/datatype/unpack_hetero.c
new file mode 100644
index 0000000..48c9c1c
--- /dev/null
+++ b/test/datatype/unpack_hetero.c
@@ -0,0 +1,99 @@
+/* -*- Mode: C; c-basic-offset:4 ; -*- */
+/*
+ * Copyright (c) 2014-2016 Research Organization for Information Science
+ *                         and Technology (RIST). All rights reserved.
+ * $COPYRIGHT$
+ *
+ * Additional copyrights may follow
+ *
+ * $HEADER$
+ */
+
+#include "opal_config.h"
+#include "opal/runtime/opal.h"
+#include "opal/datatype/opal_datatype.h"
+#include "opal/datatype/opal_datatype_internal.h"
+#include "opal/datatype/opal_convertor.h"
+#include "opal/datatype/opal_datatype_prototypes.h"
+#include "opal/util/arch.h"
+#include <time.h>
+#include <stdlib.h>
+#ifdef HAVE_SYS_TIME_H
+#include <sys/time.h>
+#endif
+#include <stdio.h>
+#include <string.h>
+
+/* Compile with:
+gcc -DHAVE_CONFIG_H -I. -I../../include -I../.. -I../../include -I../../../ompi-trunk/opal -I../../../ompi-trunk/orte -g opal_datatype_test.c -o opal_datatype_test
+*/
+
+uint32_t remote_arch = 0xffffffff;
+
+/**
+ * Main function. Call several tests and print-out the results. It try to stress the convertor
+ * using difficult data-type constructions as well as strange segment sizes for the conversion.
+ * Usually, it is able to detect most of the data-type and convertor problems. Any modifications
+ * on the data-type engine should first pass all the tests from this file, before going into other
+ * tests.
+ */
+int main( int argc, char* argv[] )
+{
+    opal_datatype_init();
+
+    /**
+     * By default simulate homogeneous architectures.
+     */
+    remote_arch = opal_local_arch ^ OPAL_ARCH_ISBIGENDIAN;
+
+    opal_convertor_t * pConv;
+    int sbuf[2], rbuf[2];
+    size_t max_data;
+    struct iovec a;
+    uint32_t iov_count;
+
+    sbuf[0] = 0x01000000; sbuf[1] = 0x02000000;
+
+    printf( "\n\n#\n * TEST UNPACKING 1 int out of 1\n#\n\n" );
+
+    pConv = opal_convertor_create( remote_arch, 0 );
+    rbuf[0] = -1; rbuf[1] = -1;
+    if( OPAL_SUCCESS != opal_convertor_prepare_for_recv( pConv, &opal_datatype_int4, 1, rbuf ) ) {
+        printf( "Cannot attach the datatype to a convertor\n" );
+        return OPAL_ERROR;
+    }
+    
+    a.iov_base = sbuf;
+    a.iov_len = 4;
+    iov_count = 1;
+    max_data = 4;
+    opal_unpack_general( pConv, &a, &iov_count, &max_data );
+
+    assert(1 == rbuf[0]);
+    assert(-1 == rbuf[1]);
+    OBJ_RELEASE(pConv);
+
+    printf( "\n\n#\n * TEST UNPACKING 1 int out of 2\n#\n\n" );
+    pConv = opal_convertor_create( remote_arch, 0 );
+    rbuf[0] = -1; rbuf[1] = -1;
+    if( OPAL_SUCCESS != opal_convertor_prepare_for_recv( pConv, &opal_datatype_int4, 2, rbuf ) ) {
+        printf( "Cannot attach the datatype to a convertor\n" );
+        return OPAL_ERROR;
+    }
+    
+
+    a.iov_base = sbuf;
+    a.iov_len = 4;
+    iov_count = 1;
+    max_data = 4;
+    opal_unpack_general( pConv, &a, &iov_count, &max_data );
+
+    assert(1 == rbuf[0]);
+    assert(-1 == rbuf[1]);
+    OBJ_RELEASE(pConv);
+
+    /* clean-ups all data allocations */
+    opal_datatype_finalize();
+    opal_finalize();
+    return OPAL_SUCCESS;
+}
ggouaillardet commented 7 years ago

@bosilca, surprisingly it seems you do not even need to configure with --enable-heterogeneous in order to reproduce the issue

bosilca commented 7 years ago

A fix has been added to #3441.