xref: /petsc/src/vec/is/sf/tests/ex24.c (revision 2ff79c18c26c94ed8cb599682f680f231dca6444) !
1 const char help[] = "Test overlapping PetscSF communication with empty roots and leaves";
2 
3 #include <petscsf.h>
4 
5 static PetscErrorCode testOverlappingCommunication(PetscSF sf)
6 {
7   PetscInt  nroots, maxleaf;
8   PetscInt *leafa, *leafb, *roota, *rootb;
9 
10   PetscFunctionBegin;
11   PetscCall(PetscSFSetUp(sf));
12   PetscCall(PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL));
13   PetscCall(PetscSFGetLeafRange(sf, NULL, &maxleaf));
14   PetscCall(PetscMalloc4(nroots, &roota, nroots, &rootb, maxleaf + 1, &leafa, maxleaf + 1, &leafb));
15 
16   // test reduce
17   for (PetscInt i = 0; i < nroots; i++) roota[i] = 0;
18   for (PetscInt i = 0; i < nroots; i++) rootb[i] = 0;
19   for (PetscInt i = 0; i < maxleaf + 1; i++) leafa[i] = (i + 1);
20   for (PetscInt i = 0; i < maxleaf + 1; i++) leafb[i] = -(i + 1);
21 
22   PetscCall(PetscSFReduceBegin(sf, MPIU_INT, leafa, roota, MPI_REPLACE));
23   PetscCall(PetscSFReduceBegin(sf, MPIU_INT, leafb, rootb, MPI_REPLACE));
24   PetscCall(PetscSFReduceEnd(sf, MPIU_INT, leafa, roota, MPI_REPLACE));
25   PetscCall(PetscSFReduceEnd(sf, MPIU_INT, leafb, rootb, MPI_REPLACE));
26   for (PetscInt i = 0; i < nroots; i++) PetscCheck(roota[i] >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PetscSFReduce in (A,B,A,B) order crosses separate reductions");
27   for (PetscInt i = 0; i < nroots; i++) PetscCheck(rootb[i] <= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PetscSFReduce in (A,B,A,B) order crosses separate reductions");
28 
29   for (PetscInt i = 0; i < nroots; i++) roota[i] = 0;
30   for (PetscInt i = 0; i < nroots; i++) rootb[i] = 0;
31 
32   PetscCall(PetscSFReduceBegin(sf, MPIU_INT, leafa, roota, MPI_REPLACE));
33   PetscCall(PetscSFReduceBegin(sf, MPIU_INT, leafb, rootb, MPI_REPLACE));
34   PetscCall(PetscSFReduceEnd(sf, MPIU_INT, leafb, rootb, MPI_REPLACE));
35   PetscCall(PetscSFReduceEnd(sf, MPIU_INT, leafa, roota, MPI_REPLACE));
36 
37   for (PetscInt i = 0; i < nroots; i++) PetscCheck(roota[i] >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PetscSFReduce in (A,B,B,A) order crosses separate reductions");
38   for (PetscInt i = 0; i < nroots; i++) PetscCheck(rootb[i] <= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PetscSFReduce in (A,B,B,A) order crosses separate reductions");
39 
40   // test bcast
41   for (PetscInt i = 0; i < nroots; i++) roota[i] = (i + 1);
42   for (PetscInt i = 0; i < nroots; i++) rootb[i] = -(i + 1);
43   for (PetscInt i = 0; i < maxleaf + 1; i++) leafa[i] = 0;
44   for (PetscInt i = 0; i < maxleaf + 1; i++) leafb[i] = 0;
45 
46   PetscCall(PetscSFBcastBegin(sf, MPIU_INT, roota, leafa, MPI_REPLACE));
47   PetscCall(PetscSFBcastBegin(sf, MPIU_INT, rootb, leafb, MPI_REPLACE));
48   PetscCall(PetscSFBcastEnd(sf, MPIU_INT, roota, leafa, MPI_REPLACE));
49   PetscCall(PetscSFBcastEnd(sf, MPIU_INT, rootb, leafb, MPI_REPLACE));
50 
51   for (PetscInt i = 0; i < maxleaf + 1; i++) PetscCheck(leafa[i] >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PetscSFBcast in (A,B,A,B) order crosses separate broadcasts");
52   for (PetscInt i = 0; i < maxleaf + 1; i++) PetscCheck(leafb[i] <= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PetscSFBcast in (A,B,A,B) order crosses separate broadcasts");
53 
54   for (PetscInt i = 0; i < maxleaf + 1; i++) leafa[i] = 0;
55   for (PetscInt i = 0; i < maxleaf + 1; i++) leafb[i] = 0;
56 
57   PetscCall(PetscSFBcastBegin(sf, MPIU_INT, roota, leafa, MPI_REPLACE));
58   PetscCall(PetscSFBcastBegin(sf, MPIU_INT, rootb, leafb, MPI_REPLACE));
59   PetscCall(PetscSFBcastEnd(sf, MPIU_INT, rootb, leafb, MPI_REPLACE));
60   PetscCall(PetscSFBcastEnd(sf, MPIU_INT, roota, leafa, MPI_REPLACE));
61 
62   for (PetscInt i = 0; i < maxleaf + 1; i++) PetscCheck(leafa[i] >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PetscSFBcast in (A,B,B,A) order crosses separate broadcasts");
63   for (PetscInt i = 0; i < maxleaf + 1; i++) PetscCheck(leafb[i] <= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PetscSFBcast in (A,B,B,A) order crosses separate broadcasts");
64 
65   PetscCall(PetscFree4(roota, rootb, leafa, leafb));
66   PetscFunctionReturn(PETSC_SUCCESS);
67 }
68 
69 static PetscErrorCode createSparseSF(MPI_Comm comm, PetscSF *sf)
70 {
71   PetscMPIInt rank;
72   PetscInt    nroots, nleaves;
73   PetscSFNode remote;
74 
75   PetscFunctionBegin;
76   PetscCallMPI(MPI_Comm_rank(comm, &rank));
77   PetscCall(PetscSFCreate(comm, sf));
78   nroots       = (rank & 1) ? 1 : 0;
79   nleaves      = (rank & 2) ? 1 : 0;
80   remote.rank  = -1;
81   remote.index = 0;
82   if (nleaves == 1) remote.rank = (rank & 1) ? (rank ^ 2) : (rank ^ 1);
83   PetscCall(PetscSFSetGraph(*sf, nroots, nleaves, NULL, PETSC_COPY_VALUES, &remote, PETSC_COPY_VALUES));
84   PetscFunctionReturn(PETSC_SUCCESS);
85 }
86 
87 int main(int argc, char **argv)
88 {
89   PetscSF sf;
90 
91   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
92   PetscCall(createSparseSF(PETSC_COMM_WORLD, &sf));
93   PetscCall(PetscSFSetFromOptions(sf));
94   PetscCall(testOverlappingCommunication(sf));
95   PetscCall(PetscSFDestroy(&sf));
96   PetscCall(PetscFinalize());
97   return 0;
98 }
99 
100 /*TEST
101 
102   test:
103     nsize: 4
104     suffix: 0
105     output_file: output/empty.out
106 
107   test:
108     TODO: frequent timeout with the CI job linux-hip-cmplx
109     nsize: 4
110     suffix: 0_window
111     output_file: output/empty.out
112     args: -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}}
113     requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW)
114 
115 TEST*/
116