xref: /petsc/src/vec/is/sf/tests/ex5.c (revision 8112c1cbf372cb53bf7c5aca994a84a6a303db4d)
114357523SPierre Jolivet static char help[] = "Test PetscSFCompose() and PetscSFCreateStridedSF() when the ilocal arrays are not identity nor dense\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petsc.h>
4c4762a1bSJed Brown #include <petscsf.h>
5c4762a1bSJed Brown 
TestVector(PetscSF sf,const char * sfname)60dd791a8SStefano Zampini static PetscErrorCode TestVector(PetscSF sf, const char *sfname)
70dd791a8SStefano Zampini {
80dd791a8SStefano Zampini   PetscInt mr, ml;
90dd791a8SStefano Zampini   MPI_Comm comm;
100dd791a8SStefano Zampini 
110dd791a8SStefano Zampini   PetscFunctionBeginUser;
120dd791a8SStefano Zampini   comm = PetscObjectComm((PetscObject)sf);
130dd791a8SStefano Zampini   PetscCall(PetscSFGetGraph(sf, &mr, NULL, NULL, NULL));
140dd791a8SStefano Zampini   PetscCall(PetscSFGetLeafRange(sf, NULL, &ml));
150dd791a8SStefano Zampini   for (PetscInt bs = 1; bs < 3; bs++) {
160dd791a8SStefano Zampini     for (PetscInt r = 0; r < 2; r++) {
170dd791a8SStefano Zampini       for (PetscInt l = 0; l < 2; l++) {
180dd791a8SStefano Zampini         PetscSF   vsf;
190dd791a8SStefano Zampini         PetscInt *rdata, *ldata, *rdatav, *ldatav;
200dd791a8SStefano Zampini         PetscInt  ldr = PETSC_DECIDE;
210dd791a8SStefano Zampini         PetscInt  ldl = PETSC_DECIDE;
220dd791a8SStefano Zampini         PetscBool flg;
230dd791a8SStefano Zampini 
240dd791a8SStefano Zampini         if (r == 1) ldr = mr;
250dd791a8SStefano Zampini         if (r == 2) ldr = mr + 7;
260dd791a8SStefano Zampini         if (l == 1) ldl = ml + 1;
270dd791a8SStefano Zampini         if (l == 2) ldl = ml + 5;
280dd791a8SStefano Zampini 
290dd791a8SStefano Zampini         PetscCall(PetscSFCreateStridedSF(sf, bs, ldr, ldl, &vsf));
300dd791a8SStefano Zampini         if (ldr == PETSC_DECIDE) ldr = mr;
310dd791a8SStefano Zampini         if (ldl == PETSC_DECIDE) ldl = ml + 1;
320dd791a8SStefano Zampini 
330dd791a8SStefano Zampini         PetscCall(PetscCalloc4(bs * ldr, &rdata, bs * ldl, &ldata, bs * ldr, &rdatav, bs * ldl, &ldatav));
340dd791a8SStefano Zampini         for (PetscInt i = 0; i < bs * ldr; i++) rdata[i] = i + 1;
350dd791a8SStefano Zampini 
360dd791a8SStefano Zampini         for (PetscInt i = 0; i < bs; i++) {
370dd791a8SStefano Zampini           PetscCall(PetscSFBcastBegin(sf, MPIU_INT, PetscSafePointerPlusOffset(rdata, i * ldr), PetscSafePointerPlusOffset(ldata, i * ldl), MPI_REPLACE));
380dd791a8SStefano Zampini           PetscCall(PetscSFBcastEnd(sf, MPIU_INT, PetscSafePointerPlusOffset(rdata, i * ldr), PetscSafePointerPlusOffset(ldata, i * ldl), MPI_REPLACE));
390dd791a8SStefano Zampini         }
400dd791a8SStefano Zampini         PetscCall(PetscSFBcastBegin(vsf, MPIU_INT, rdata, ldatav, MPI_REPLACE));
410dd791a8SStefano Zampini         PetscCall(PetscSFBcastEnd(vsf, MPIU_INT, rdata, ldatav, MPI_REPLACE));
420dd791a8SStefano Zampini         PetscCall(PetscArraycmp(ldata, ldatav, bs * ldl, &flg));
430dd791a8SStefano Zampini 
44*5440e5dcSBarry Smith         PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &flg, 1, MPI_C_BOOL, MPI_LAND, comm));
450dd791a8SStefano Zampini         if (!flg) {
460dd791a8SStefano Zampini           PetscCall(PetscPrintf(comm, "Error with Bcast on %s: block size %" PetscInt_FMT ", ldr %" PetscInt_FMT ", ldl %" PetscInt_FMT "\n", sfname, bs, ldr, ldl));
470dd791a8SStefano Zampini           PetscCall(PetscPrintf(comm, "Single SF\n"));
480dd791a8SStefano Zampini           PetscCall(PetscIntView(bs * ldl, ldata, PETSC_VIEWER_STDOUT_(comm)));
490dd791a8SStefano Zampini           PetscCall(PetscPrintf(comm, "Vector SF\n"));
500dd791a8SStefano Zampini           PetscCall(PetscIntView(bs * ldl, ldatav, PETSC_VIEWER_STDOUT_(comm)));
510dd791a8SStefano Zampini         }
520dd791a8SStefano Zampini         PetscCall(PetscArrayzero(rdata, bs * ldr));
530dd791a8SStefano Zampini         PetscCall(PetscArrayzero(rdatav, bs * ldr));
540dd791a8SStefano Zampini 
550dd791a8SStefano Zampini         for (PetscInt i = 0; i < bs; i++) {
560dd791a8SStefano Zampini           PetscCall(PetscSFReduceBegin(sf, MPIU_INT, PetscSafePointerPlusOffset(ldata, i * ldl), PetscSafePointerPlusOffset(rdata, i * ldr), MPI_SUM));
570dd791a8SStefano Zampini           PetscCall(PetscSFReduceEnd(sf, MPIU_INT, PetscSafePointerPlusOffset(ldata, i * ldl), PetscSafePointerPlusOffset(rdata, i * ldr), MPI_SUM));
580dd791a8SStefano Zampini         }
590dd791a8SStefano Zampini         PetscCall(PetscSFReduceBegin(vsf, MPIU_INT, ldata, rdatav, MPI_SUM));
600dd791a8SStefano Zampini         PetscCall(PetscSFReduceEnd(vsf, MPIU_INT, ldata, rdatav, MPI_SUM));
610dd791a8SStefano Zampini         PetscCall(PetscArraycmp(rdata, rdatav, bs * ldr, &flg));
62*5440e5dcSBarry Smith         PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &flg, 1, MPI_C_BOOL, MPI_LAND, comm));
630dd791a8SStefano Zampini         if (!flg) {
640dd791a8SStefano Zampini           PetscCall(PetscPrintf(comm, "Error with Reduce on %s: block size %" PetscInt_FMT ", ldr %" PetscInt_FMT ", ldl %" PetscInt_FMT "\n", sfname, bs, ldr, ldl));
650dd791a8SStefano Zampini           PetscCall(PetscPrintf(comm, "Single SF\n"));
660dd791a8SStefano Zampini           PetscCall(PetscIntView(bs * ldr, rdata, PETSC_VIEWER_STDOUT_(comm)));
670dd791a8SStefano Zampini           PetscCall(PetscPrintf(comm, "Vector SF\n"));
680dd791a8SStefano Zampini           PetscCall(PetscIntView(bs * ldr, rdatav, PETSC_VIEWER_STDOUT_(comm)));
690dd791a8SStefano Zampini         }
700dd791a8SStefano Zampini         PetscCall(PetscFree4(rdata, ldata, rdatav, ldatav));
710dd791a8SStefano Zampini         PetscCall(PetscSFDestroy(&vsf));
720dd791a8SStefano Zampini       }
730dd791a8SStefano Zampini     }
740dd791a8SStefano Zampini   }
750dd791a8SStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
760dd791a8SStefano Zampini }
770dd791a8SStefano Zampini 
main(int argc,char ** argv)78d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
79d71ae5a4SJacob Faibussowitsch {
80c4762a1bSJed Brown   PetscSF      sfA, sfB, sfBA, sfAAm, sfBBm, sfAm, sfBm;
81c4762a1bSJed Brown   PetscInt     nrootsA, nleavesA, nrootsB, nleavesB;
82c4762a1bSJed Brown   PetscInt    *ilocalA, *ilocalB;
83c4762a1bSJed Brown   PetscSFNode *iremoteA, *iremoteB;
84c4762a1bSJed Brown   PetscMPIInt  rank, size;
85c4762a1bSJed Brown   PetscInt     i, m, n, k, nl = 2, mA, mB, nldataA, nldataB;
86c4762a1bSJed Brown   PetscInt    *rdA, *rdB, *ldA, *ldB;
870dd791a8SStefano Zampini   PetscBool    inverse = PETSC_FALSE, test_vector = PETSC_TRUE;
88c4762a1bSJed Brown 
89327415f7SBarry Smith   PetscFunctionBeginUser;
909566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
919566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-nl", &nl, NULL));
929566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-explicit_inverse", &inverse, NULL));
939566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
949566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
95c4762a1bSJed Brown 
969566063dSJacob Faibussowitsch   PetscCall(PetscSFCreate(PETSC_COMM_WORLD, &sfA));
979566063dSJacob Faibussowitsch   PetscCall(PetscSFCreate(PETSC_COMM_WORLD, &sfB));
989566063dSJacob Faibussowitsch   PetscCall(PetscSFSetFromOptions(sfA));
999566063dSJacob Faibussowitsch   PetscCall(PetscSFSetFromOptions(sfB));
100c4762a1bSJed Brown 
1010dd791a8SStefano Zampini   // disable vector tests with linux-misc-32bit and sftype window
1020dd791a8SStefano Zampini #if (PETSC_SIZEOF_SIZE_T == 4)
1030dd791a8SStefano Zampini   {
1040dd791a8SStefano Zampini     PetscBool iswin;
1050dd791a8SStefano Zampini 
1060dd791a8SStefano Zampini     PetscCall(PetscObjectTypeCompare((PetscObject)sfA, PETSCSFWINDOW, &iswin));
1070dd791a8SStefano Zampini     if (iswin) test_vector = PETSC_FALSE;
1080dd791a8SStefano Zampini   }
1090dd791a8SStefano Zampini #endif
1100dd791a8SStefano Zampini   PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_vector", &test_vector, NULL));
1110dd791a8SStefano Zampini 
112c4762a1bSJed Brown   n = 4 * nl * size;
113c4762a1bSJed Brown   m = 2 * nl;
114c4762a1bSJed Brown   k = nl;
115c4762a1bSJed Brown 
116dd400576SPatrick Sanan   nldataA = rank == 0 ? n : 0;
117c4762a1bSJed Brown   nldataB = 3 * nl;
118c4762a1bSJed Brown 
119c4762a1bSJed Brown   nrootsA  = m;
120dd400576SPatrick Sanan   nleavesA = rank == 0 ? size * m : 0;
121dd400576SPatrick Sanan   nrootsB  = rank == 0 ? n : 0;
122c4762a1bSJed Brown   nleavesB = k;
123c4762a1bSJed Brown 
1249566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nleavesA, &ilocalA));
1259566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nleavesA, &iremoteA));
1269566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nleavesB, &ilocalB));
1279566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nleavesB, &iremoteB));
128c4762a1bSJed Brown 
129c4762a1bSJed Brown   /* sf A bcast is equivalent to a sparse gather on process 0
130c4762a1bSJed Brown      process 0 receives data in the middle [nl,3*nl] of the leaf data array for A */
131c4762a1bSJed Brown   for (i = 0; i < nleavesA; i++) {
132300f1712SStefano Zampini     iremoteA[i].rank  = i / m;
133c4762a1bSJed Brown     iremoteA[i].index = i % m;
134c4762a1bSJed Brown     ilocalA[i]        = nl + i / m * 4 * nl + i % m;
135c4762a1bSJed Brown   }
136c4762a1bSJed Brown 
137c4762a1bSJed Brown   /* sf B bcast is equivalent to a sparse scatter from process 0
138c4762a1bSJed Brown      process 0 sends data from [nl,2*nl] of the leaf data array for A
139c4762a1bSJed Brown      each process receives, in reverse order, in the middle [nl,2*nl] of the leaf data array for B */
140c4762a1bSJed Brown   for (i = 0; i < nleavesB; i++) {
141c4762a1bSJed Brown     iremoteB[i].rank  = 0;
142c4762a1bSJed Brown     iremoteB[i].index = rank * 4 * nl + nl + i % m;
143c4762a1bSJed Brown     ilocalB[i]        = 2 * nl - i - 1;
144c4762a1bSJed Brown   }
1459566063dSJacob Faibussowitsch   PetscCall(PetscSFSetGraph(sfA, nrootsA, nleavesA, ilocalA, PETSC_OWN_POINTER, iremoteA, PETSC_OWN_POINTER));
1469566063dSJacob Faibussowitsch   PetscCall(PetscSFSetGraph(sfB, nrootsB, nleavesB, ilocalB, PETSC_OWN_POINTER, iremoteB, PETSC_OWN_POINTER));
1479566063dSJacob Faibussowitsch   PetscCall(PetscSFSetUp(sfA));
1489566063dSJacob Faibussowitsch   PetscCall(PetscSFSetUp(sfB));
1499566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)sfA, "sfA"));
1509566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)sfB, "sfB"));
1519566063dSJacob Faibussowitsch   PetscCall(PetscSFViewFromOptions(sfA, NULL, "-view"));
1529566063dSJacob Faibussowitsch   PetscCall(PetscSFViewFromOptions(sfB, NULL, "-view"));
153c4762a1bSJed Brown 
1540dd791a8SStefano Zampini   if (test_vector) PetscCall(TestVector(sfA, "sfA"));
1550dd791a8SStefano Zampini   if (test_vector) PetscCall(TestVector(sfB, "sfB"));
1560dd791a8SStefano Zampini 
1579566063dSJacob Faibussowitsch   PetscCall(PetscSFGetLeafRange(sfA, NULL, &mA));
1589566063dSJacob Faibussowitsch   PetscCall(PetscSFGetLeafRange(sfB, NULL, &mB));
1599566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(nrootsA, &rdA, nldataA, &ldA));
1609566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(nrootsB, &rdB, nldataB, &ldB));
161c4762a1bSJed Brown   for (i = 0; i < nrootsA; i++) rdA[i] = m * rank + i;
162c4762a1bSJed Brown   for (i = 0; i < nldataA; i++) ldA[i] = -1;
163c4762a1bSJed Brown   for (i = 0; i < nldataB; i++) ldB[i] = -1;
164c4762a1bSJed Brown 
1659566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "BcastB(BcastA)\n"));
1669566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "A: root data\n"));
1679566063dSJacob Faibussowitsch   PetscCall(PetscIntView(nrootsA, rdA, PETSC_VIEWER_STDOUT_WORLD));
1689566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastBegin(sfA, MPIU_INT, rdA, ldA, MPI_REPLACE));
1699566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastEnd(sfA, MPIU_INT, rdA, ldA, MPI_REPLACE));
1709566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "A: leaf data (all)\n"));
1719566063dSJacob Faibussowitsch   PetscCall(PetscIntView(nldataA, ldA, PETSC_VIEWER_STDOUT_WORLD));
1729566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastBegin(sfB, MPIU_INT, ldA, ldB, MPI_REPLACE));
1739566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastEnd(sfB, MPIU_INT, ldA, ldB, MPI_REPLACE));
1749566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "B: leaf data (all)\n"));
1759566063dSJacob Faibussowitsch   PetscCall(PetscIntView(nldataB, ldB, PETSC_VIEWER_STDOUT_WORLD));
176c4762a1bSJed Brown 
1779566063dSJacob Faibussowitsch   PetscCall(PetscSFCompose(sfA, sfB, &sfBA));
1789566063dSJacob Faibussowitsch   PetscCall(PetscSFSetFromOptions(sfBA));
1799566063dSJacob Faibussowitsch   PetscCall(PetscSFSetUp(sfBA));
1809566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)sfBA, "sfBA"));
1819566063dSJacob Faibussowitsch   PetscCall(PetscSFViewFromOptions(sfBA, NULL, "-view"));
1820dd791a8SStefano Zampini   if (test_vector) PetscCall(TestVector(sfBA, "sfBA"));
183c4762a1bSJed Brown 
184c4762a1bSJed Brown   for (i = 0; i < nldataB; i++) ldB[i] = -1;
1859566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "BcastBA\n"));
1869566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "BA: root data\n"));
1879566063dSJacob Faibussowitsch   PetscCall(PetscIntView(nrootsA, rdA, PETSC_VIEWER_STDOUT_WORLD));
1889566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastBegin(sfBA, MPIU_INT, rdA, ldB, MPI_REPLACE));
1899566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastEnd(sfBA, MPIU_INT, rdA, ldB, MPI_REPLACE));
1909566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "BA: leaf data (all)\n"));
1919566063dSJacob Faibussowitsch   PetscCall(PetscIntView(nldataB, ldB, PETSC_VIEWER_STDOUT_WORLD));
192c4762a1bSJed Brown 
1939566063dSJacob Faibussowitsch   PetscCall(PetscSFCreateInverseSF(sfA, &sfAm));
1949566063dSJacob Faibussowitsch   PetscCall(PetscSFSetFromOptions(sfAm));
1959566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)sfAm, "sfAm"));
1969566063dSJacob Faibussowitsch   PetscCall(PetscSFViewFromOptions(sfAm, NULL, "-view"));
1970dd791a8SStefano Zampini   if (test_vector) PetscCall(TestVector(sfAm, "sfAm"));
198c4762a1bSJed Brown 
199c4762a1bSJed Brown   if (!inverse) {
2009566063dSJacob Faibussowitsch     PetscCall(PetscSFComposeInverse(sfA, sfA, &sfAAm));
201c4762a1bSJed Brown   } else {
2029566063dSJacob Faibussowitsch     PetscCall(PetscSFCompose(sfA, sfAm, &sfAAm));
203c4762a1bSJed Brown   }
2049566063dSJacob Faibussowitsch   PetscCall(PetscSFSetFromOptions(sfAAm));
2059566063dSJacob Faibussowitsch   PetscCall(PetscSFSetUp(sfAAm));
2069566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)sfAAm, "sfAAm"));
2079566063dSJacob Faibussowitsch   PetscCall(PetscSFViewFromOptions(sfAAm, NULL, "-view"));
2080dd791a8SStefano Zampini   if (test_vector) PetscCall(TestVector(sfAAm, "sfAAm"));
209c4762a1bSJed Brown 
2109566063dSJacob Faibussowitsch   PetscCall(PetscSFCreateInverseSF(sfB, &sfBm));
2119566063dSJacob Faibussowitsch   PetscCall(PetscSFSetFromOptions(sfBm));
2129566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)sfBm, "sfBm"));
2139566063dSJacob Faibussowitsch   PetscCall(PetscSFViewFromOptions(sfBm, NULL, "-view"));
2140dd791a8SStefano Zampini   if (test_vector) PetscCall(TestVector(sfBm, "sfBm"));
215c4762a1bSJed Brown 
216c4762a1bSJed Brown   if (!inverse) {
2179566063dSJacob Faibussowitsch     PetscCall(PetscSFComposeInverse(sfB, sfB, &sfBBm));
218c4762a1bSJed Brown   } else {
2199566063dSJacob Faibussowitsch     PetscCall(PetscSFCompose(sfB, sfBm, &sfBBm));
220c4762a1bSJed Brown   }
2219566063dSJacob Faibussowitsch   PetscCall(PetscSFSetFromOptions(sfBBm));
2229566063dSJacob Faibussowitsch   PetscCall(PetscSFSetUp(sfBBm));
2239566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)sfBBm, "sfBBm"));
2249566063dSJacob Faibussowitsch   PetscCall(PetscSFViewFromOptions(sfBBm, NULL, "-view"));
2250dd791a8SStefano Zampini   if (test_vector) PetscCall(TestVector(sfBBm, "sfBBm"));
226c4762a1bSJed Brown 
2279566063dSJacob Faibussowitsch   PetscCall(PetscFree2(rdA, ldA));
2289566063dSJacob Faibussowitsch   PetscCall(PetscFree2(rdB, ldB));
229c4762a1bSJed Brown 
2309566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&sfA));
2319566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&sfB));
2329566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&sfBA));
2339566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&sfAm));
2349566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&sfBm));
2359566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&sfAAm));
2369566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&sfBBm));
237c4762a1bSJed Brown 
2389566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
239b122ec5aSJacob Faibussowitsch   return 0;
240c4762a1bSJed Brown }
241c4762a1bSJed Brown 
242c4762a1bSJed Brown /*TEST
243c4762a1bSJed Brown 
244c4762a1bSJed Brown    test:
245c4762a1bSJed Brown      suffix: 1
246c4762a1bSJed Brown      args: -view -explicit_inverse {{0 1}}
247c4762a1bSJed Brown 
248c4762a1bSJed Brown    test:
249c4762a1bSJed Brown      nsize: 7
250c4762a1bSJed Brown      filter: grep -v "type" | grep -v "sort"
251c4762a1bSJed Brown      suffix: 2
252c4762a1bSJed Brown      args: -view -nl 5 -explicit_inverse {{0 1}}
253c4762a1bSJed Brown 
254c4762a1bSJed Brown    # we cannot test for -sf_window_flavor dynamic because SFCompose with sparse leaves may change the root data pointer only locally, and this is not supported by the dynamic case
255c4762a1bSJed Brown    test:
256c204e120SJunchao Zhang      TODO: frequent timeout with the CI job linux-hip-cmplx
257c4762a1bSJed Brown      nsize: 7
258c4762a1bSJed Brown      suffix: 2_window
259c4762a1bSJed Brown      filter: grep -v "type" | grep -v "sort"
260c4762a1bSJed Brown      output_file: output/ex5_2.out
261c4762a1bSJed Brown      args: -view -nl 5 -explicit_inverse {{0 1}} -sf_type window -sf_window_sync {{fence lock active}} -sf_window_flavor {{create allocate}}
262dfd57a17SPierre Jolivet      requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW)
263c4762a1bSJed Brown 
264c4762a1bSJed Brown    # The nightly test suite with MPICH uses ch3:sock, which is broken when winsize == 0 in some of the processes
265c4762a1bSJed Brown    test:
266c204e120SJunchao Zhang      TODO: frequent timeout with the CI job linux-hip-cmplx
267c4762a1bSJed Brown      nsize: 7
268c4762a1bSJed Brown      suffix: 2_window_shared
269c4762a1bSJed Brown      filter: grep -v "type" | grep -v "sort"
270c4762a1bSJed Brown      output_file: output/ex5_2.out
271c4762a1bSJed Brown      args: -view -nl 5 -explicit_inverse {{0 1}} -sf_type window -sf_window_sync {{fence lock active}} -sf_window_flavor shared
272e3c15826SJunchao Zhang      requires: defined(PETSC_HAVE_MPI_PROCESS_SHARED_MEMORY) !defined(PETSC_HAVE_MPICH) defined(PETSC_HAVE_MPI_ONE_SIDED)
273c4762a1bSJed Brown 
274c4762a1bSJed Brown TEST*/
275