Hi,

I'm currently working on a bug occuring at the client site with openmpi when calling MPI_Sendreceive() on datatypes built by the application. I think I've found where the bug comes from (it is located in opal_generic_simple_pack_function() - file opal/datatype/opal_datatype_pack.c). But this code is so complicated that I'm more than unsure of my fix. What I can say is that it fixes things for me, but I need some advices from the datatypes specialists.

---------------

You will find in attachment the reproducer provided by the client, as well as the resulting output.
datatypes.c : reproducer
to run the binary: salloc --exclusive -p B510 -N 1 -n 1 mpirun ./datatypes
trc_ko: traces got without the patch applied
trc_ok: traces got with the patch applied.

---------------

The proposed patch is the following: (Note that the very first change in this patch was enough in my case, but I thought all the "source_base" settings should follow this model.)

-------------------------
opal_generic_simple_pack_function: add the datatype lb when progressing in the input buffer

diff -r cb23c2f07e1f opal/datatype/opal_datatype_pack.c
--- a/opal/datatype/opal_datatype_pack.c Sun Nov 24 17:06:51 2013 +0000 +++ b/opal/datatype/opal_datatype_pack.c Mon Nov 25 10:48:00 2013 +0100
@@ -301,7 +301,7 @@ opal_generic_simple_pack_function( opal_
                 PACK_PREDEFINED_DATATYPE( pConvertor, pElem, count_desc,
source_base, destination, iov_len_local );
                 if( 0 == count_desc ) {  /* completed */
-                    source_base = pConvertor->pBaseBuf + pStack->disp;
+ source_base = pConvertor->pBaseBuf + pStack->disp + pData->lb;
                     pos_desc++;  /* advance to the next data */
UPDATE_INTERNAL_COUNTERS( description, pos_desc, pElem, count_desc );
                     continue;
@@ -333,7 +333,7 @@ opal_generic_simple_pack_function( opal_
pStack->disp += description[pStack->index].loop.extent;
                     }
                 }
-                source_base = pConvertor->pBaseBuf + pStack->disp;
+ source_base = pConvertor->pBaseBuf + pStack->disp + pData->lb; UPDATE_INTERNAL_COUNTERS( description, pos_desc, pElem, count_desc ); DO_DEBUG( opal_output( 0, "pack new_loop count %d stack_pos %d pos_desc %d disp %ld space %lu\n", (int)pStack->count, pConvertor->stack_pos, pos_desc, (long)pStack->disp, (unsigned long)iov_len_local ); );
@@ -354,7 +354,7 @@ opal_generic_simple_pack_function( opal_
                             pStack->disp + local_disp);
                 pos_desc++;
             update_loop_description:  /* update the current state */
-                source_base = pConvertor->pBaseBuf + pStack->disp;
+ source_base = pConvertor->pBaseBuf + pStack->disp + pData->lb; UPDATE_INTERNAL_COUNTERS( description, pos_desc, pElem, count_desc ); DDT_DUMP_STACK( pConvertor->pStack, pConvertor->stack_pos, pElem, "advance loop" );
                 continue;
@@ -374,7 +374,7 @@ opal_generic_simple_pack_function( opal_
     }
     /* I complete an element, next step I should go to the next one */
PUSH_STACK( pStack, pConvertor->stack_pos, pos_desc, OPAL_DATATYPE_INT8, count_desc,
-                source_base - pStack->disp - pConvertor->pBaseBuf );
+ source_base - pStack->disp - pConvertor->pBaseBuf - pData->lb ); DO_DEBUG( opal_output( 0, "pack save stack stack_pos %d pos_desc %d count_desc %d disp %ld\n", pConvertor->stack_pos, pStack->index, (int)pStack->count, (long)pStack->disp ); );
     return 0;

-------------------------------

Regards,
Nadia

--
Nadia Derbey
Bull, Architect of an Open World
http://www.bull.com

#include <mpi.h>
#include <stdio.h>

/**
 * expected output for defective OpenMPI versions:
results_2[6]     = 8
ref_results_2[6] = 12
results_2[7]     = 9
ref_results_2[7] = 13

*/

static int
do_test(MPI_Datatype * recvs, MPI_Datatype * sends, int *raw_inputs);




int main(void) {

    int rank, size;
    int ierror = 0;
    MPI_Datatype sends[2], recvs[2];

    MPI_Init(NULL, NULL);
    MPI_Comm_rank(MPI_COMM_WORLD, &rank);
    MPI_Comm_size(MPI_COMM_WORLD, &size);


    {
        int count = 2, blocklen = 2, stride = 4;

        MPI_Type_vector(count, blocklen, stride, MPI_INT, &recvs[0]);
        MPI_Type_commit(&recvs[0]);

        MPI_Type_vector(count, blocklen, stride, MPI_INT, &sends[1]);
        MPI_Type_commit(&sends[1]);
    }

    {
        int count = 1;
        int blocklength = 4;
        int array_of_displacements[] = {4};

        MPI_Type_create_indexed_block(count, blocklength, array_of_displacements, MPI_INT, &sends[0]);
        MPI_Type_commit(&sends[0]);

        MPI_Type_create_indexed_block(count, blocklength, array_of_displacements, MPI_INT, &recvs[1]);
        MPI_Type_commit(&recvs[1]);
    }

    {
        int raw_input[16] = {0,1,2,3,4,5,6,7,8,9,10,11,12,13,14,15};
        ierror = do_test(recvs, sends, raw_input);
    }

    MPI_Type_free(&sends[1]);
    MPI_Type_free(&recvs[1]);
    MPI_Type_free(&sends[0]);
    MPI_Type_free(&recvs[0]);

    MPI_Finalize();

  return ierror;
}

static int
do_test(MPI_Datatype * recvs, MPI_Datatype * sends, int *raw_inputs)
{

    int results_buffer[16] = {-1,-1,-1,-1,-1,-1,-1,-1,-1,-1,-1,-1,-1,-1,-1,-1};
    int i;
    MPI_Datatype send_dt, recv_dt;
    MPI_Status status;
    int ref_results[16] = {4,5,-1,-1,6,7,-1,-1,-1,-1,-1,-1,8,9,12,13};
    int retcode = 0;

    {
        int count = 2;
        int array_of_blocklengths[] = {1, 1};
        MPI_Aint array_of_displacements[] = {0, 32};

        MPI_Type_create_struct(count, array_of_blocklengths, array_of_displacements, recvs, &recv_dt);
        MPI_Type_commit(&recv_dt);

        MPI_Type_create_struct(count, array_of_blocklengths, array_of_displacements, sends, &send_dt);
        MPI_Type_commit(&send_dt);
    }

    MPI_Sendrecv(raw_inputs, 1, send_dt, 0, 3, results_buffer, 1, recv_dt, 0, 3, MPI_COMM_WORLD, MPI_STATUS_IGNORE);

    for (i = 0; i < 16; i++) {
        retcode |= (results_buffer[i] != ref_results[i]);
        if (results_buffer[i] != ref_results[i])
            fprintf(stdout, "   results_buffer[%d]     = %2d\t!=\tref_results[%d] = %d\n", i, results_buffer[i], i, ref_results[i]);
        else
            fprintf(stdout, "   results_buffer[%d]     = %2d\tEQ\tref_results[%d] = %d\n", i, results_buffer[i], i, ref_results[i]);
    }

    MPI_Type_free(&send_dt);
    MPI_Type_free(&recv_dt);

    return retcode;
}
   results_buffer[0]     =  4   EQ      ref_results[0] = 4
   results_buffer[1]     =  5   EQ      ref_results[1] = 5
   results_buffer[2]     = -1   EQ      ref_results[2] = -1
   results_buffer[3]     = -1   EQ      ref_results[3] = -1
   results_buffer[4]     =  6   EQ      ref_results[4] = 6
   results_buffer[5]     =  7   EQ      ref_results[5] = 7
   results_buffer[6]     = -1   EQ      ref_results[6] = -1
   results_buffer[7]     = -1   EQ      ref_results[7] = -1
   results_buffer[8]     = -1   EQ      ref_results[8] = -1
   results_buffer[9]     = -1   EQ      ref_results[9] = -1
   results_buffer[10]     = -1  EQ      ref_results[10] = -1
   results_buffer[11]     = -1  EQ      ref_results[11] = -1
   results_buffer[12]     =  8  EQ      ref_results[12] = 8
   results_buffer[13]     =  9  EQ      ref_results[13] = 9
   results_buffer[14]     =  8  !=      ref_results[14] = 12
   results_buffer[15]     =  9  !=      ref_results[15] = 13
   results_buffer[0]     =  4   EQ      ref_results[0] = 4
   results_buffer[1]     =  5   EQ      ref_results[1] = 5
   results_buffer[2]     = -1   EQ      ref_results[2] = -1
   results_buffer[3]     = -1   EQ      ref_results[3] = -1
   results_buffer[4]     =  6   EQ      ref_results[4] = 6
   results_buffer[5]     =  7   EQ      ref_results[5] = 7
   results_buffer[6]     = -1   EQ      ref_results[6] = -1
   results_buffer[7]     = -1   EQ      ref_results[7] = -1
   results_buffer[8]     = -1   EQ      ref_results[8] = -1
   results_buffer[9]     = -1   EQ      ref_results[9] = -1
   results_buffer[10]     = -1  EQ      ref_results[10] = -1
   results_buffer[11]     = -1  EQ      ref_results[11] = -1
   results_buffer[12]     =  8  EQ      ref_results[12] = 8
   results_buffer[13]     =  9  EQ      ref_results[13] = 9
   results_buffer[14]     = 12  EQ      ref_results[14] = 12
   results_buffer[15]     = 13  EQ      ref_results[15] = 13

Reply via email to