[petsc-users] Issue with window SF type using derived datatypes for reduction

Lawrence Mitchell lawrence.mitchell at imperial.ac.uk
Mon Feb 9 12:22:01 CST 2015


Hi all,

I'm trying to use an SFReduce to update some data using a contiguous derived type:

sketch:

MPI_Type_contiguous(bs, MPI_DOUBLE, &dtype);
MPI_Type_commit(&dtype);

PetscSFReduceBegin(sf, dtype, rootdata, leafdata, MPI_SUM);
PetscSFReduceEnd(sf, dtype, rootdata, leafdata, MPI_SUM);

All is well if the block size (bs) is 1.  However, for larger block sizes I get unexpected behaviour (the reduction appears not to SUM, but rather overwrite the local data).  This is using -sf_type window, with OpenMPI 1.8.3.  -sf_type basic works fine.  The attached code illustrates the problem.

Run with:

$ mpiexec -n 2 ./petsc_sftest -data_bs 2 -sf_type basic
PetscSF Object: 2 MPI processes
  type: basic
    sort=rank-order
  [0] Number of roots=3, leaves=8, remote ranks=2
  [0] 0 <- (0,0)
  [0] 1 <- (0,1)
  [0] 2 <- (0,2)
  [0] 3 <- (1,1)
  [0] 7 <- (1,5)
  [0] 5 <- (1,3)
  [0] 4 <- (1,2)
  [0] 6 <- (1,4)
  [1] Number of roots=6, leaves=8, remote ranks=2
  [1] 1 <- (1,1)
  [1] 5 <- (1,5)
  [1] 0 <- (1,0)
  [1] 3 <- (1,3)
  [1] 2 <- (1,2)
  [1] 4 <- (1,4)
  [1] 6 <- (0,1)
  [1] 7 <- (0,2)
Vec Object: 2 MPI processes
  type: mpi
Process [0]
10
10
10
10
10
10
0
0
0
0
0
0
0
0
0
0
Process [1]
12
12
12
12
12
12
12
12
12
12
12
12
0
0
0
0
Vec Object: 2 MPI processes
  type: mpi
Process [0]
10
10
10
10
10
10
Process [1]
12
12
12
12
12
12
12
12
12
12
12
12


With the window type:

$ mpiexec -n 2 ./petsc_sftest -data_bs 2 -sf_type window
PetscSF Object: 2 MPI processes
  type: window
    synchronization=FENCE sort=rank-order
  [0] Number of roots=3, leaves=8, remote ranks=2
  [0] 0 <- (0,0)
  [0] 1 <- (0,1)
  [0] 2 <- (0,2)
  [0] 3 <- (1,1)
  [0] 7 <- (1,5)
  [0] 5 <- (1,3)
  [0] 4 <- (1,2)
  [0] 6 <- (1,4)
  [1] Number of roots=6, leaves=8, remote ranks=2
  [1] 1 <- (1,1)
  [1] 5 <- (1,5)
  [1] 0 <- (1,0)
  [1] 3 <- (1,3)
  [1] 2 <- (1,2)
  [1] 4 <- (1,4)
  [1] 6 <- (0,1)
  [1] 7 <- (0,2)
Vec Object: 2 MPI processes
  type: mpi
Process [0]
10
10
10
10
10
10
0
0
0
0
0
0
0
0
0
0
Process [1]
12
12
12
12
12
12
12
12
12
12
12
12
0
0
0
0
Vec Object: 2 MPI processes
  type: mpi
Process [0]
10
10
0
0
0
0
Process [1]
12
12
12
12
12
12
12
12
12
12
12
12

Note how Vec B on process 0 has entries 10, 10, 0, 0, 0, 0 (rather than all 10).

Having just tried a build with --download-mpich, I notice this problem does not occur.  So should I shout at the OpenMPI team?

Cheers,

Lawrence

static const char help[] = "Test overlapped communication on a single star forest (PetscSF)\n\n";

#include <petscvec.h>
#include <petscsf.h>
#include <petscviewer.h>

#undef __FUNCT__
#define __FUNCT__ "main"
int main(int argc, char **argv)
{
    PetscInt ierr;
    PetscSF sf;
    Vec A;
    Vec B;
    double *bufA;
    double *bufB;
    MPI_Comm c;
    PetscMPIInt rank, size;
    PetscInt nroots, nleaves;
    PetscInt i;
    PetscInt *ilocal;
    PetscSFNode *iremote;
    MPI_Datatype dtype;
    PetscInt bs = 1;
    PetscInitialize(&argc,&argv,NULL,help);

    ierr = PetscOptionsGetInt(NULL, "-data_bs", &bs, NULL);CHKERRQ(ierr);
    c = PETSC_COMM_WORLD;

    ierr = MPI_Comm_rank(c,&rank);CHKERRQ(ierr);
    ierr = MPI_Comm_size(c,&size);CHKERRQ(ierr);

    if (size != 2) {
        SETERRQ(c, PETSC_ERR_USER, "Only coded for two MPI processes\n");
    }

    ierr = PetscSFCreate(c,&sf);CHKERRQ(ierr);
    ierr = PetscSFSetFromOptions(sf);CHKERRQ(ierr);

    nleaves = 8;
    nroots = 3 * (rank + 1);

    ierr = PetscMalloc1(nleaves,&ilocal);CHKERRQ(ierr);

    ierr = PetscMalloc1(nleaves,&iremote);CHKERRQ(ierr);
    if ( rank == 0 ) {
        ilocal[0] = 0; ilocal[1] = 1; ilocal[2] = 2; ilocal[3] = 3;
        ilocal[4] = 7; ilocal[5] = 5; ilocal[6] = 4; ilocal[7] = 6;
        iremote[0].rank = 0;
        iremote[0].index = 0;
        iremote[1].rank = 0;
        iremote[1].index = 1;
        iremote[2].rank = 0;
        iremote[2].index = 2;
        iremote[3].rank = 1;
        iremote[3].index = 1;
        iremote[4].rank = 1;
        iremote[4].index = 5;
        iremote[5].rank = 1;
        iremote[5].index = 3;
        iremote[6].rank = 1;
        iremote[6].index = 2;
        iremote[7].rank = 1;
        iremote[7].index = 4;
    } else {
        ilocal[0] = 1; ilocal[1] = 5; ilocal[2] = 0; ilocal[3] = 3;
        ilocal[4] = 2; ilocal[5] = 4; ilocal[6] = 6; ilocal[7] = 7;
        iremote[0].rank = 1;
        iremote[0].index = 1;
        iremote[1].rank = 1;
        iremote[1].index = 5;
        iremote[2].rank = 1;
        iremote[2].index = 0;
        iremote[3].rank = 1;
        iremote[3].index = 3;
        iremote[4].rank = 1;
        iremote[4].index = 2;
        iremote[5].rank = 1;
        iremote[5].index = 4;
        iremote[6].rank = 0;
        iremote[6].index = 1;
        iremote[7].rank = 0;
        iremote[7].index = 2;
    }
    ierr = PetscSFSetGraph(sf,nroots,nleaves,ilocal,PETSC_OWN_POINTER,
                           iremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
    ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
    ierr = PetscSFView(sf,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
    ierr = VecCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
    ierr = VecSetSizes(A,nleaves*bs,PETSC_DETERMINE);CHKERRQ(ierr);
    ierr = VecSetFromOptions(A);CHKERRQ(ierr);
    ierr = VecSetUp(A);CHKERRQ(ierr);
    ierr = VecCreate(PETSC_COMM_WORLD,&B);CHKERRQ(ierr);
    ierr = VecSetSizes(B,nroots*bs,PETSC_DETERMINE);CHKERRQ(ierr);
    ierr = VecSetFromOptions(B);CHKERRQ(ierr);
    ierr = VecSetUp(B);CHKERRQ(ierr);
    ierr = VecGetArray(A,&bufA);CHKERRQ(ierr);
    for (i=0; i < nroots*bs; i++) {
        bufA[i] = 10.0 + 2*rank;
    }

    ierr = VecRestoreArray(A,&bufA);CHKERRQ(ierr);
    ierr = VecGetArray(A,&bufA);CHKERRQ(ierr);
    ierr = VecGetArray(B,&bufB);CHKERRQ(ierr);
    ierr = MPI_Type_contiguous(bs, MPI_DOUBLE, &dtype); CHKERRQ(ierr);
    ierr = MPI_Type_commit(&dtype); CHKERRQ(ierr);
    ierr = PetscSFReduceBegin(sf,dtype,(const void*)bufA,(void *)bufB, MPI_SUM);CHKERRQ(ierr);
    ierr = PetscSFReduceEnd(sf,dtype,(const void*)bufA,(void *)bufB, MPI_SUM);CHKERRQ(ierr);

    ierr = VecRestoreArray(A,&bufA);CHKERRQ(ierr);
    ierr = VecRestoreArray(B,&bufB);CHKERRQ(ierr);

    ierr = VecView(A,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
    ierr = VecView(B,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
    ierr = VecDestroy(&A);CHKERRQ(ierr);
    ierr = VecDestroy(&B);CHKERRQ(ierr);
    ierr = PetscSFDestroy(&sf);CHKERRQ(ierr);

    ierr = MPI_Type_free(&dtype); CHKERRQ(ierr);
    PetscFinalize();

    return 0;
}

-------------- next part --------------
A non-text attachment was scrubbed...
Name: signature.asc
Type: application/pgp-signature
Size: 455 bytes
Desc: Message signed with OpenPGP using GPGMail
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20150209/848608fb/attachment.pgp>


More information about the petsc-users mailing list