xref: /petsc/src/vec/is/sf/tests/ex24.c (revision 834855d6effb0d027771461c8e947ee1ce5a1e17)
133e3ecb2SToby Isaac const char help[] = "Test overlapping PetscSF communication with empty roots and leaves";
233e3ecb2SToby Isaac 
333e3ecb2SToby Isaac #include <petscsf.h>
433e3ecb2SToby Isaac 
testOverlappingCommunication(PetscSF sf)533e3ecb2SToby Isaac static PetscErrorCode testOverlappingCommunication(PetscSF sf)
633e3ecb2SToby Isaac {
733e3ecb2SToby Isaac   PetscInt  nroots, maxleaf;
833e3ecb2SToby Isaac   PetscInt *leafa, *leafb, *roota, *rootb;
933e3ecb2SToby Isaac 
1033e3ecb2SToby Isaac   PetscFunctionBegin;
1133e3ecb2SToby Isaac   PetscCall(PetscSFSetUp(sf));
1233e3ecb2SToby Isaac   PetscCall(PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL));
1333e3ecb2SToby Isaac   PetscCall(PetscSFGetLeafRange(sf, NULL, &maxleaf));
1433e3ecb2SToby Isaac   PetscCall(PetscMalloc4(nroots, &roota, nroots, &rootb, maxleaf + 1, &leafa, maxleaf + 1, &leafb));
1533e3ecb2SToby Isaac 
1633e3ecb2SToby Isaac   // test reduce
1733e3ecb2SToby Isaac   for (PetscInt i = 0; i < nroots; i++) roota[i] = 0;
1833e3ecb2SToby Isaac   for (PetscInt i = 0; i < nroots; i++) rootb[i] = 0;
1933e3ecb2SToby Isaac   for (PetscInt i = 0; i < maxleaf + 1; i++) leafa[i] = (i + 1);
2033e3ecb2SToby Isaac   for (PetscInt i = 0; i < maxleaf + 1; i++) leafb[i] = -(i + 1);
2133e3ecb2SToby Isaac 
2233e3ecb2SToby Isaac   PetscCall(PetscSFReduceBegin(sf, MPIU_INT, leafa, roota, MPI_REPLACE));
2333e3ecb2SToby Isaac   PetscCall(PetscSFReduceBegin(sf, MPIU_INT, leafb, rootb, MPI_REPLACE));
2433e3ecb2SToby Isaac   PetscCall(PetscSFReduceEnd(sf, MPIU_INT, leafa, roota, MPI_REPLACE));
2533e3ecb2SToby Isaac   PetscCall(PetscSFReduceEnd(sf, MPIU_INT, leafb, rootb, MPI_REPLACE));
2633e3ecb2SToby Isaac   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");
2733e3ecb2SToby Isaac   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");
2833e3ecb2SToby Isaac 
2933e3ecb2SToby Isaac   for (PetscInt i = 0; i < nroots; i++) roota[i] = 0;
3033e3ecb2SToby Isaac   for (PetscInt i = 0; i < nroots; i++) rootb[i] = 0;
3133e3ecb2SToby Isaac 
3233e3ecb2SToby Isaac   PetscCall(PetscSFReduceBegin(sf, MPIU_INT, leafa, roota, MPI_REPLACE));
3333e3ecb2SToby Isaac   PetscCall(PetscSFReduceBegin(sf, MPIU_INT, leafb, rootb, MPI_REPLACE));
3433e3ecb2SToby Isaac   PetscCall(PetscSFReduceEnd(sf, MPIU_INT, leafb, rootb, MPI_REPLACE));
3533e3ecb2SToby Isaac   PetscCall(PetscSFReduceEnd(sf, MPIU_INT, leafa, roota, MPI_REPLACE));
3633e3ecb2SToby Isaac 
3733e3ecb2SToby Isaac   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");
3833e3ecb2SToby Isaac   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");
3933e3ecb2SToby Isaac 
4033e3ecb2SToby Isaac   // test bcast
4133e3ecb2SToby Isaac   for (PetscInt i = 0; i < nroots; i++) roota[i] = (i + 1);
4233e3ecb2SToby Isaac   for (PetscInt i = 0; i < nroots; i++) rootb[i] = -(i + 1);
4333e3ecb2SToby Isaac   for (PetscInt i = 0; i < maxleaf + 1; i++) leafa[i] = 0;
4433e3ecb2SToby Isaac   for (PetscInt i = 0; i < maxleaf + 1; i++) leafb[i] = 0;
4533e3ecb2SToby Isaac 
4633e3ecb2SToby Isaac   PetscCall(PetscSFBcastBegin(sf, MPIU_INT, roota, leafa, MPI_REPLACE));
4733e3ecb2SToby Isaac   PetscCall(PetscSFBcastBegin(sf, MPIU_INT, rootb, leafb, MPI_REPLACE));
4833e3ecb2SToby Isaac   PetscCall(PetscSFBcastEnd(sf, MPIU_INT, roota, leafa, MPI_REPLACE));
4933e3ecb2SToby Isaac   PetscCall(PetscSFBcastEnd(sf, MPIU_INT, rootb, leafb, MPI_REPLACE));
5033e3ecb2SToby Isaac 
5133e3ecb2SToby Isaac   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");
5233e3ecb2SToby Isaac   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");
5333e3ecb2SToby Isaac 
5433e3ecb2SToby Isaac   for (PetscInt i = 0; i < maxleaf + 1; i++) leafa[i] = 0;
5533e3ecb2SToby Isaac   for (PetscInt i = 0; i < maxleaf + 1; i++) leafb[i] = 0;
5633e3ecb2SToby Isaac 
5733e3ecb2SToby Isaac   PetscCall(PetscSFBcastBegin(sf, MPIU_INT, roota, leafa, MPI_REPLACE));
5833e3ecb2SToby Isaac   PetscCall(PetscSFBcastBegin(sf, MPIU_INT, rootb, leafb, MPI_REPLACE));
5933e3ecb2SToby Isaac   PetscCall(PetscSFBcastEnd(sf, MPIU_INT, rootb, leafb, MPI_REPLACE));
6033e3ecb2SToby Isaac   PetscCall(PetscSFBcastEnd(sf, MPIU_INT, roota, leafa, MPI_REPLACE));
6133e3ecb2SToby Isaac 
6233e3ecb2SToby Isaac   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");
6333e3ecb2SToby Isaac   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");
6433e3ecb2SToby Isaac 
6533e3ecb2SToby Isaac   PetscCall(PetscFree4(roota, rootb, leafa, leafb));
6633e3ecb2SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
6733e3ecb2SToby Isaac }
6833e3ecb2SToby Isaac 
createSparseSF(MPI_Comm comm,PetscSF * sf)6933e3ecb2SToby Isaac static PetscErrorCode createSparseSF(MPI_Comm comm, PetscSF *sf)
7033e3ecb2SToby Isaac {
7133e3ecb2SToby Isaac   PetscMPIInt rank;
7233e3ecb2SToby Isaac   PetscInt    nroots, nleaves;
7333e3ecb2SToby Isaac   PetscSFNode remote;
7433e3ecb2SToby Isaac 
7533e3ecb2SToby Isaac   PetscFunctionBegin;
7633e3ecb2SToby Isaac   PetscCallMPI(MPI_Comm_rank(comm, &rank));
7733e3ecb2SToby Isaac   PetscCall(PetscSFCreate(comm, sf));
7833e3ecb2SToby Isaac   nroots       = (rank & 1) ? 1 : 0;
7933e3ecb2SToby Isaac   nleaves      = (rank & 2) ? 1 : 0;
8033e3ecb2SToby Isaac   remote.rank  = -1;
8133e3ecb2SToby Isaac   remote.index = 0;
82*ac530a7eSPierre Jolivet   if (nleaves == 1) remote.rank = (rank & 1) ? (rank ^ 2) : (rank ^ 1);
8333e3ecb2SToby Isaac   PetscCall(PetscSFSetGraph(*sf, nroots, nleaves, NULL, PETSC_COPY_VALUES, &remote, PETSC_COPY_VALUES));
8433e3ecb2SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
8533e3ecb2SToby Isaac }
8633e3ecb2SToby Isaac 
main(int argc,char ** argv)8733e3ecb2SToby Isaac int main(int argc, char **argv)
8833e3ecb2SToby Isaac {
8933e3ecb2SToby Isaac   PetscSF sf;
9033e3ecb2SToby Isaac 
9133e3ecb2SToby Isaac   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
9233e3ecb2SToby Isaac   PetscCall(createSparseSF(PETSC_COMM_WORLD, &sf));
9333e3ecb2SToby Isaac   PetscCall(PetscSFSetFromOptions(sf));
9433e3ecb2SToby Isaac   PetscCall(testOverlappingCommunication(sf));
9533e3ecb2SToby Isaac   PetscCall(PetscSFDestroy(&sf));
9633e3ecb2SToby Isaac   PetscCall(PetscFinalize());
9733e3ecb2SToby Isaac   return 0;
9833e3ecb2SToby Isaac }
9933e3ecb2SToby Isaac 
10033e3ecb2SToby Isaac /*TEST
10133e3ecb2SToby Isaac 
10233e3ecb2SToby Isaac   test:
10333e3ecb2SToby Isaac     nsize: 4
10433e3ecb2SToby Isaac     suffix: 0
1053886731fSPierre Jolivet     output_file: output/empty.out
10633e3ecb2SToby Isaac 
10733e3ecb2SToby Isaac   test:
108c204e120SJunchao Zhang     TODO: frequent timeout with the CI job linux-hip-cmplx
10933e3ecb2SToby Isaac     nsize: 4
11033e3ecb2SToby Isaac     suffix: 0_window
1113886731fSPierre Jolivet     output_file: output/empty.out
11233e3ecb2SToby Isaac     args: -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}}
11333e3ecb2SToby Isaac     requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW)
11433e3ecb2SToby Isaac 
11533e3ecb2SToby Isaac TEST*/
116