xref: /petsc/src/vec/is/sf/tests/ex25.c (revision e91c04dfc8a52dee1965211bb1cc8e5bf775178f)
1 static const char help[] = "Test PetscSF with derived data types created with MPI large count\n\n";
2 
3 #include <petscsys.h>
4 #include <petscsf.h>
5 
6 int main(int argc, char **argv)
7 {
8   PetscSF        sf;
9   PetscInt       i, nroots, nleaves;
10   const PetscInt m = 4, n = 64;
11   PetscSFNode   *iremote = NULL;
12   PetscMPIInt    rank, size;
13   int           *rootdata = NULL, *leafdata = NULL;
14   MPI_Datatype   newtype;
15 
16   PetscFunctionBeginUser;
17   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
18   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
19   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
20   PetscCheck(size == 2, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "The test can only run with two MPI ranks");
21 
22   PetscCall(PetscSFCreate(PETSC_COMM_WORLD, &sf));
23   PetscCall(PetscSFSetFromOptions(sf));
24 
25   if (rank == 0) {
26     nroots  = n;
27     nleaves = 0;
28   } else {
29     nroots  = 0;
30     nleaves = n;
31     PetscCall(PetscMalloc1(nleaves, &iremote));
32     for (i = 0; i < nleaves; i++) {
33       iremote[i].rank  = 0;
34       iremote[i].index = i;
35     }
36   }
37   PetscCall(PetscSFSetGraph(sf, nroots, nleaves, NULL, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER));
38 
39   PetscCall(PetscCalloc2(nroots * m, &rootdata, nleaves * m, &leafdata)); // allocate fat nodes to apply a derived data type of m MPI_INTs
40 
41   if (rank == 0) rootdata[nroots * m - 1] = 123; // set the last integer in rootdata and then check on leaves
42 
43 #if defined(PETSC_HAVE_MPI_LARGE_COUNT)
44   PetscCallMPI(MPI_Type_contiguous_c(m, MPI_INT, &newtype));
45 #else
46   PetscCallMPI(MPI_Type_contiguous(m, MPI_INT, &newtype));
47 #endif
48 
49   PetscCallMPI(MPI_Type_commit(&newtype));
50 
51   PetscCall(PetscSFBcastBegin(sf, newtype, rootdata, leafdata, MPI_REPLACE)); //  bcast rootdata to leafdata
52   PetscCall(PetscSFBcastEnd(sf, newtype, rootdata, leafdata, MPI_REPLACE));
53 
54   if (rank == 1) PetscCheck(leafdata[nleaves * m - 1] == 123, PETSC_COMM_SELF, PETSC_ERR_PLIB, "SF: wrong results");
55 
56   PetscCallMPI(MPI_Type_free(&newtype));
57   PetscCall(PetscFree2(rootdata, leafdata));
58   PetscCall(PetscSFDestroy(&sf));
59   PetscCall(PetscFinalize());
60   return 0;
61 }
62 
63 /**TEST
64    test:
65      nsize: 2
66      output_file: output/empty.out
67      args: -sf_type {{basic neighbor}}
68 
69 TEST**/
70