1c4762a1bSJed Brown static const char help[] = "Test star forest communication (PetscSF)\n\n"; 2c4762a1bSJed Brown 31bb3edfdSPatrick Sanan /* 4c4762a1bSJed Brown Description: A star is a simple tree with one root and zero or more leaves. 5c4762a1bSJed Brown A star forest is a union of disjoint stars. 6c4762a1bSJed Brown Many common communication patterns can be expressed as updates of rootdata using leafdata and vice-versa. 7c4762a1bSJed Brown This example creates a star forest, communicates values using the graph (see options for types of communication), views the graph, then destroys it. 81bb3edfdSPatrick Sanan */ 9c4762a1bSJed Brown 10c4762a1bSJed Brown /* 11c4762a1bSJed Brown Include petscsf.h so we can use PetscSF objects. Note that this automatically 12c4762a1bSJed Brown includes petscsys.h. 13c4762a1bSJed Brown */ 14c4762a1bSJed Brown #include <petscsf.h> 15c4762a1bSJed Brown #include <petscviewer.h> 16c4762a1bSJed Brown 17c4762a1bSJed Brown /* like PetscSFView() but with alternative array of local indices */ 18d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscSFViewCustomLocals_Private(PetscSF sf, const PetscInt locals[], PetscViewer viewer) 19d71ae5a4SJacob Faibussowitsch { 20c4762a1bSJed Brown const PetscSFNode *iremote; 21*6497c311SBarry Smith PetscInt i, nroots, nleaves; 22*6497c311SBarry Smith PetscMPIInt rank, nranks; 23c4762a1bSJed Brown 24c4762a1bSJed Brown PetscFunctionBeginUser; 259566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)sf), &rank)); 269566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sf, &nroots, &nleaves, NULL, &iremote)); 279566063dSJacob Faibussowitsch PetscCall(PetscSFGetRootRanks(sf, &nranks, NULL, NULL, NULL, NULL)); 289566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 299566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 30*6497c311SBarry Smith PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] Number of roots=%" PetscInt_FMT ", leaves=%" PetscInt_FMT ", remote ranks=%d\n", rank, nroots, nleaves, nranks)); 31*6497c311SBarry Smith for (i = 0; i < nleaves; i++) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] %" PetscInt_FMT " <- (%d,%" PetscInt_FMT ")\n", rank, locals[i], (PetscMPIInt)iremote[i].rank, iremote[i].index)); 329566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 339566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 349566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 36c4762a1bSJed Brown } 37c4762a1bSJed Brown 38d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 39d71ae5a4SJacob Faibussowitsch { 400dd791a8SStefano Zampini PetscInt i, bs = 2, nroots, nrootsalloc, nleaves, nleavesalloc, *mine, stride; 41c4762a1bSJed Brown PetscSFNode *remote; 42c4762a1bSJed Brown PetscMPIInt rank, size; 430dd791a8SStefano Zampini PetscSF sf, vsf; 440dd791a8SStefano Zampini PetscBool test_all, test_bcast, test_bcastop, test_reduce, test_degree, test_fetchandop, test_gather, test_scatter, test_embed, test_invert, test_sf_distribute, test_char, test_vector = PETSC_FALSE; 45c4762a1bSJed Brown MPI_Op mop = MPI_OP_NULL; /* initialize to prevent compiler warnings with cxx_quad build */ 46c4762a1bSJed Brown char opstring[256]; 47c4762a1bSJed Brown PetscBool strflg; 48c4762a1bSJed Brown 49327415f7SBarry Smith PetscFunctionBeginUser; 509566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 519566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 529566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 53c4762a1bSJed Brown 54d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, "", "PetscSF Test Options", "none"); 55c4762a1bSJed Brown test_all = PETSC_FALSE; 569566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-test_all", "Test all SF communications", "", test_all, &test_all, NULL)); 57c4762a1bSJed Brown test_bcast = test_all; 589566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-test_bcast", "Test broadcast", "", test_bcast, &test_bcast, NULL)); 59c4762a1bSJed Brown test_bcastop = test_all; 609566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-test_bcastop", "Test broadcast and reduce", "", test_bcastop, &test_bcastop, NULL)); 61c4762a1bSJed Brown test_reduce = test_all; 629566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-test_reduce", "Test reduction", "", test_reduce, &test_reduce, NULL)); 63c4762a1bSJed Brown test_char = test_all; 649566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-test_char", "Test signed char, unsigned char, and char", "", test_char, &test_char, NULL)); 65c4762a1bSJed Brown mop = MPI_SUM; 66c6a7a370SJeremy L Thompson PetscCall(PetscStrncpy(opstring, "sum", sizeof(opstring))); 679566063dSJacob Faibussowitsch PetscCall(PetscOptionsString("-test_op", "Designate which MPI_Op to use", "", opstring, opstring, sizeof(opstring), NULL)); 689566063dSJacob Faibussowitsch PetscCall(PetscStrcmp("sum", opstring, &strflg)); 69ad540459SPierre Jolivet if (strflg) mop = MPIU_SUM; 709566063dSJacob Faibussowitsch PetscCall(PetscStrcmp("prod", opstring, &strflg)); 71ad540459SPierre Jolivet if (strflg) mop = MPI_PROD; 729566063dSJacob Faibussowitsch PetscCall(PetscStrcmp("max", opstring, &strflg)); 73ad540459SPierre Jolivet if (strflg) mop = MPI_MAX; 749566063dSJacob Faibussowitsch PetscCall(PetscStrcmp("min", opstring, &strflg)); 75ad540459SPierre Jolivet if (strflg) mop = MPI_MIN; 769566063dSJacob Faibussowitsch PetscCall(PetscStrcmp("land", opstring, &strflg)); 77ad540459SPierre Jolivet if (strflg) mop = MPI_LAND; 789566063dSJacob Faibussowitsch PetscCall(PetscStrcmp("band", opstring, &strflg)); 79ad540459SPierre Jolivet if (strflg) mop = MPI_BAND; 809566063dSJacob Faibussowitsch PetscCall(PetscStrcmp("lor", opstring, &strflg)); 81ad540459SPierre Jolivet if (strflg) mop = MPI_LOR; 829566063dSJacob Faibussowitsch PetscCall(PetscStrcmp("bor", opstring, &strflg)); 83ad540459SPierre Jolivet if (strflg) mop = MPI_BOR; 849566063dSJacob Faibussowitsch PetscCall(PetscStrcmp("lxor", opstring, &strflg)); 85ad540459SPierre Jolivet if (strflg) mop = MPI_LXOR; 869566063dSJacob Faibussowitsch PetscCall(PetscStrcmp("bxor", opstring, &strflg)); 87ad540459SPierre Jolivet if (strflg) mop = MPI_BXOR; 88c4762a1bSJed Brown test_degree = test_all; 899566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-test_degree", "Test computation of vertex degree", "", test_degree, &test_degree, NULL)); 90c4762a1bSJed Brown test_fetchandop = test_all; 919566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-test_fetchandop", "Test atomic Fetch-And-Op", "", test_fetchandop, &test_fetchandop, NULL)); 92c4762a1bSJed Brown test_gather = test_all; 939566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-test_gather", "Test point gather", "", test_gather, &test_gather, NULL)); 94c4762a1bSJed Brown test_scatter = test_all; 959566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-test_scatter", "Test point scatter", "", test_scatter, &test_scatter, NULL)); 96c4762a1bSJed Brown test_embed = test_all; 979566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-test_embed", "Test point embed", "", test_embed, &test_embed, NULL)); 98c4762a1bSJed Brown test_invert = test_all; 999566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-test_invert", "Test point invert", "", test_invert, &test_invert, NULL)); 100c4762a1bSJed Brown stride = 1; 1019566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-stride", "Stride for leaf and root data", "", stride, &stride, NULL)); 102c4762a1bSJed Brown test_sf_distribute = PETSC_FALSE; 1039566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-test_sf_distribute", "Create an SF that 'distributes' to each process, like an alltoall", "", test_sf_distribute, &test_sf_distribute, NULL)); 1049566063dSJacob Faibussowitsch PetscCall(PetscOptionsString("-test_op", "Designate which MPI_Op to use", "", opstring, opstring, sizeof(opstring), NULL)); 1050dd791a8SStefano Zampini PetscCall(PetscOptionsInt("-bs", "Block size for vectorial SF", "", bs, &bs, NULL)); 1060dd791a8SStefano Zampini PetscCall(PetscOptionsBool("-test_vector", "Run tests using the vectorial SF", "", test_vector, &test_vector, NULL)); 107d0609cedSBarry Smith PetscOptionsEnd(); 108c4762a1bSJed Brown 109c4762a1bSJed Brown if (test_sf_distribute) { 110c4762a1bSJed Brown nroots = size; 111c4762a1bSJed Brown nrootsalloc = size; 112c4762a1bSJed Brown nleaves = size; 113c4762a1bSJed Brown nleavesalloc = size; 114c4762a1bSJed Brown mine = NULL; 1159566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nleaves, &remote)); 116*6497c311SBarry Smith for (PetscMPIInt ii = 0; ii < size; ii++) { 117*6497c311SBarry Smith remote[ii].rank = ii; 118*6497c311SBarry Smith remote[ii].index = rank; 119c4762a1bSJed Brown } 120c4762a1bSJed Brown } else { 121c4762a1bSJed Brown nroots = 2 + (PetscInt)(rank == 0); 122c4762a1bSJed Brown nrootsalloc = nroots * stride; 123c4762a1bSJed Brown nleaves = 2 + (PetscInt)(rank > 0); 124c4762a1bSJed Brown nleavesalloc = nleaves * stride; 125c4762a1bSJed Brown mine = NULL; 126c4762a1bSJed Brown if (stride > 1) { 127c4762a1bSJed Brown PetscInt i; 128c4762a1bSJed Brown 1299566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nleaves, &mine)); 130ad540459SPierre Jolivet for (i = 0; i < nleaves; i++) mine[i] = stride * i; 131c4762a1bSJed Brown } 1329566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nleaves, &remote)); 133c4762a1bSJed Brown /* Left periodic neighbor */ 134c4762a1bSJed Brown remote[0].rank = (rank + size - 1) % size; 135c4762a1bSJed Brown remote[0].index = 1 * stride; 136c4762a1bSJed Brown /* Right periodic neighbor */ 137c4762a1bSJed Brown remote[1].rank = (rank + 1) % size; 138c4762a1bSJed Brown remote[1].index = 0 * stride; 139c4762a1bSJed Brown if (rank > 0) { /* All processes reference rank 0, index 1 */ 140c4762a1bSJed Brown remote[2].rank = 0; 141c4762a1bSJed Brown remote[2].index = 2 * stride; 142c4762a1bSJed Brown } 143c4762a1bSJed Brown } 144c4762a1bSJed Brown 1450dd791a8SStefano Zampini /* Create a star forest for communication */ 1469566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(PETSC_COMM_WORLD, &sf)); 1479566063dSJacob Faibussowitsch PetscCall(PetscSFSetFromOptions(sf)); 1489566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(sf, nrootsalloc, nleaves, mine, PETSC_OWN_POINTER, remote, PETSC_OWN_POINTER)); 1499566063dSJacob Faibussowitsch PetscCall(PetscSFSetUp(sf)); 150c4762a1bSJed Brown 1510dd791a8SStefano Zampini /* We can also obtain a strided SF to communicate interleaved blocks of data */ 1520dd791a8SStefano Zampini PetscCall(PetscSFCreateStridedSF(sf, bs, PETSC_DECIDE, PETSC_DECIDE, &vsf)); 1530dd791a8SStefano Zampini if (test_vector) { /* perform all tests below using the vectorial SF */ 1540dd791a8SStefano Zampini PetscSF t = sf; 1550dd791a8SStefano Zampini 1560dd791a8SStefano Zampini sf = vsf; 1570dd791a8SStefano Zampini vsf = t; 1580dd791a8SStefano Zampini 1590dd791a8SStefano Zampini nroots *= bs; 1600dd791a8SStefano Zampini nleaves *= bs; 1610dd791a8SStefano Zampini nrootsalloc *= bs; 1620dd791a8SStefano Zampini nleavesalloc *= bs; 1630dd791a8SStefano Zampini } 1640dd791a8SStefano Zampini PetscCall(PetscSFDestroy(&vsf)); 1650dd791a8SStefano Zampini 166c4762a1bSJed Brown /* View graph, mostly useful for debugging purposes. */ 1679566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_ASCII_INFO_DETAIL)); 1689566063dSJacob Faibussowitsch PetscCall(PetscSFView(sf, PETSC_VIEWER_STDOUT_WORLD)); 1699566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD)); 170c4762a1bSJed Brown 171c4762a1bSJed Brown if (test_bcast) { /* broadcast rootdata into leafdata */ 172c4762a1bSJed Brown PetscInt *rootdata, *leafdata; 1732d4ee042Sprj- /* Allocate space for send and receive buffers. This example communicates PetscInt, but other types, including 174c4762a1bSJed Brown * user-defined structures, could also be used. */ 1759566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(nrootsalloc, &rootdata, nleavesalloc, &leafdata)); 176c4762a1bSJed Brown /* Set rootdata buffer to be broadcast */ 177c4762a1bSJed Brown for (i = 0; i < nrootsalloc; i++) rootdata[i] = -1; 178c4762a1bSJed Brown for (i = 0; i < nroots; i++) rootdata[i * stride] = 100 * (rank + 1) + i; 179c4762a1bSJed Brown /* Initialize local buffer, these values are never used. */ 180c4762a1bSJed Brown for (i = 0; i < nleavesalloc; i++) leafdata[i] = -1; 181c4762a1bSJed Brown /* Broadcast entries from rootdata to leafdata. Computation or other communication can be performed between the begin and end calls. */ 1828e54d7e8SToby Isaac // test persistent communication code paths by repeated bcast several times 1838e54d7e8SToby Isaac PetscCall(PetscSFRegisterPersistent(sf, MPIU_INT, rootdata, leafdata)); 1848e54d7e8SToby Isaac for (PetscInt i = 0; i < 10; i++) { 1859566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(sf, MPIU_INT, rootdata, leafdata, MPI_REPLACE)); 1869566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(sf, MPIU_INT, rootdata, leafdata, MPI_REPLACE)); 1878e54d7e8SToby Isaac } 1888e54d7e8SToby Isaac PetscCall(PetscSFDeregisterPersistent(sf, MPIU_INT, rootdata, leafdata)); 1899566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Bcast Rootdata\n")); 1909566063dSJacob Faibussowitsch PetscCall(PetscIntView(nrootsalloc, rootdata, PETSC_VIEWER_STDOUT_WORLD)); 1919566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Bcast Leafdata\n")); 1929566063dSJacob Faibussowitsch PetscCall(PetscIntView(nleavesalloc, leafdata, PETSC_VIEWER_STDOUT_WORLD)); 1939566063dSJacob Faibussowitsch PetscCall(PetscFree2(rootdata, leafdata)); 194c4762a1bSJed Brown } 195c4762a1bSJed Brown 196c4762a1bSJed Brown if (test_bcast && test_char) { /* Bcast with char */ 197c4762a1bSJed Brown PetscInt len; 198c4762a1bSJed Brown char buf[256]; 199c4762a1bSJed Brown char *rootdata, *leafdata; 2009566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(nrootsalloc, &rootdata, nleavesalloc, &leafdata)); 201c4762a1bSJed Brown /* Set rootdata buffer to be broadcast */ 202c4762a1bSJed Brown for (i = 0; i < nrootsalloc; i++) rootdata[i] = '*'; 203*6497c311SBarry Smith for (i = 0; i < nroots; i++) rootdata[i * stride] = (char)('A' + rank * 3 + i); /* rank is very small, so it is fine to compute a char */ 204c4762a1bSJed Brown /* Initialize local buffer, these values are never used. */ 205c4762a1bSJed Brown for (i = 0; i < nleavesalloc; i++) leafdata[i] = '?'; 206c4762a1bSJed Brown 2079566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(sf, MPI_CHAR, rootdata, leafdata, MPI_REPLACE)); 2089566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(sf, MPI_CHAR, rootdata, leafdata, MPI_REPLACE)); 209c4762a1bSJed Brown 2109566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Bcast Rootdata in type of char\n")); 2119371c9d4SSatish Balay len = 0; 2129371c9d4SSatish Balay PetscCall(PetscSNPrintf(buf, 256, "%4d:", rank)); 2139371c9d4SSatish Balay len += 5; 2149371c9d4SSatish Balay for (i = 0; i < nrootsalloc; i++) { 2159371c9d4SSatish Balay PetscCall(PetscSNPrintf(buf + len, 256 - len, "%5c", rootdata[i])); 2169371c9d4SSatish Balay len += 5; 2179371c9d4SSatish Balay } 2189566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%s\n", buf)); 2199566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT)); 220c4762a1bSJed Brown 2219566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Bcast Leafdata in type of char\n")); 2229371c9d4SSatish Balay len = 0; 2239371c9d4SSatish Balay PetscCall(PetscSNPrintf(buf, 256, "%4d:", rank)); 2249371c9d4SSatish Balay len += 5; 2259371c9d4SSatish Balay for (i = 0; i < nleavesalloc; i++) { 2269371c9d4SSatish Balay PetscCall(PetscSNPrintf(buf + len, 256 - len, "%5c", leafdata[i])); 2279371c9d4SSatish Balay len += 5; 2289371c9d4SSatish Balay } 2299566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%s\n", buf)); 2309566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT)); 231c4762a1bSJed Brown 2329566063dSJacob Faibussowitsch PetscCall(PetscFree2(rootdata, leafdata)); 233c4762a1bSJed Brown } 234c4762a1bSJed Brown 235c4762a1bSJed Brown if (test_bcastop) { /* Reduce rootdata into leafdata */ 236c4762a1bSJed Brown PetscInt *rootdata, *leafdata; 237c4762a1bSJed Brown /* Allocate space for send and receive buffers. This example communicates PetscInt, but other types, including 238c4762a1bSJed Brown * user-defined structures, could also be used. */ 2399566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(nrootsalloc, &rootdata, nleavesalloc, &leafdata)); 240c4762a1bSJed Brown /* Set rootdata buffer to be broadcast */ 241c4762a1bSJed Brown for (i = 0; i < nrootsalloc; i++) rootdata[i] = -1; 242c4762a1bSJed Brown for (i = 0; i < nroots; i++) rootdata[i * stride] = 100 * (rank + 1) + i; 243c4762a1bSJed Brown /* Set leaf values to reduce with */ 244c4762a1bSJed Brown for (i = 0; i < nleavesalloc; i++) leafdata[i] = -10 * (rank + 1) - i; 2459566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Pre-BcastAndOp Leafdata\n")); 2469566063dSJacob Faibussowitsch PetscCall(PetscIntView(nleavesalloc, leafdata, PETSC_VIEWER_STDOUT_WORLD)); 247c4762a1bSJed Brown /* Broadcast entries from rootdata to leafdata. Computation or other communication can be performed between the begin and end calls. */ 2489566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(sf, MPIU_INT, rootdata, leafdata, mop)); 2499566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(sf, MPIU_INT, rootdata, leafdata, mop)); 2509566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## BcastAndOp Rootdata\n")); 2519566063dSJacob Faibussowitsch PetscCall(PetscIntView(nrootsalloc, rootdata, PETSC_VIEWER_STDOUT_WORLD)); 2529566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## BcastAndOp Leafdata\n")); 2539566063dSJacob Faibussowitsch PetscCall(PetscIntView(nleavesalloc, leafdata, PETSC_VIEWER_STDOUT_WORLD)); 2549566063dSJacob Faibussowitsch PetscCall(PetscFree2(rootdata, leafdata)); 255c4762a1bSJed Brown } 256c4762a1bSJed Brown 257c4762a1bSJed Brown if (test_reduce) { /* Reduce leafdata into rootdata */ 258c4762a1bSJed Brown PetscInt *rootdata, *leafdata; 2599566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(nrootsalloc, &rootdata, nleavesalloc, &leafdata)); 260c4762a1bSJed Brown /* Initialize rootdata buffer in which the result of the reduction will appear. */ 261c4762a1bSJed Brown for (i = 0; i < nrootsalloc; i++) rootdata[i] = -1; 262c4762a1bSJed Brown for (i = 0; i < nroots; i++) rootdata[i * stride] = 100 * (rank + 1) + i; 263c4762a1bSJed Brown /* Set leaf values to reduce. */ 264c4762a1bSJed Brown for (i = 0; i < nleavesalloc; i++) leafdata[i] = -1; 265c4762a1bSJed Brown for (i = 0; i < nleaves; i++) leafdata[i * stride] = 1000 * (rank + 1) + 10 * i; 2669566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Pre-Reduce Rootdata\n")); 2679566063dSJacob Faibussowitsch PetscCall(PetscIntView(nrootsalloc, rootdata, PETSC_VIEWER_STDOUT_WORLD)); 268c4762a1bSJed Brown /* Perform reduction. Computation or other communication can be performed between the begin and end calls. 269c4762a1bSJed Brown * This example sums the values, but other MPI_Ops can be used (e.g MPI_MAX, MPI_PROD). */ 2709566063dSJacob Faibussowitsch PetscCall(PetscSFReduceBegin(sf, MPIU_INT, leafdata, rootdata, mop)); 2719566063dSJacob Faibussowitsch PetscCall(PetscSFReduceEnd(sf, MPIU_INT, leafdata, rootdata, mop)); 2729566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Reduce Leafdata\n")); 2739566063dSJacob Faibussowitsch PetscCall(PetscIntView(nleavesalloc, leafdata, PETSC_VIEWER_STDOUT_WORLD)); 2749566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Reduce Rootdata\n")); 2759566063dSJacob Faibussowitsch PetscCall(PetscIntView(nrootsalloc, rootdata, PETSC_VIEWER_STDOUT_WORLD)); 2769566063dSJacob Faibussowitsch PetscCall(PetscFree2(rootdata, leafdata)); 277c4762a1bSJed Brown } 278c4762a1bSJed Brown 279c4762a1bSJed Brown if (test_reduce && test_char) { /* Reduce with signed char */ 280c4762a1bSJed Brown PetscInt len; 281c4762a1bSJed Brown char buf[256]; 282c4762a1bSJed Brown signed char *rootdata, *leafdata; 2839566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(nrootsalloc, &rootdata, nleavesalloc, &leafdata)); 284c4762a1bSJed Brown /* Initialize rootdata buffer in which the result of the reduction will appear. */ 285c4762a1bSJed Brown for (i = 0; i < nrootsalloc; i++) rootdata[i] = -1; 286*6497c311SBarry Smith for (i = 0; i < nroots; i++) rootdata[i * stride] = (signed char)(10 * (rank + 1) + i); 287c4762a1bSJed Brown /* Set leaf values to reduce. */ 288c4762a1bSJed Brown for (i = 0; i < nleavesalloc; i++) leafdata[i] = -1; 289*6497c311SBarry Smith for (i = 0; i < nleaves; i++) leafdata[i * stride] = (signed char)(50 * (rank + 1) + 10 * i); 2909566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Pre-Reduce Rootdata in type of signed char\n")); 291c4762a1bSJed Brown 2929371c9d4SSatish Balay len = 0; 2939371c9d4SSatish Balay PetscCall(PetscSNPrintf(buf, 256, "%4d:", rank)); 2949371c9d4SSatish Balay len += 5; 2959371c9d4SSatish Balay for (i = 0; i < nrootsalloc; i++) { 2969371c9d4SSatish Balay PetscCall(PetscSNPrintf(buf + len, 256 - len, "%5d", rootdata[i])); 2979371c9d4SSatish Balay len += 5; 2989371c9d4SSatish Balay } 2999566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%s\n", buf)); 3009566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT)); 301c4762a1bSJed Brown 302c4762a1bSJed Brown /* Using MPI_CHAR should trigger an error since MPI standard does not support reduction on MPI_CHAR. 303c4762a1bSJed Brown Testing with -test_op max, one can see the sign does take effect in MPI_MAX. 304c4762a1bSJed Brown */ 3059566063dSJacob Faibussowitsch PetscCall(PetscSFReduceBegin(sf, MPI_SIGNED_CHAR, leafdata, rootdata, mop)); 3069566063dSJacob Faibussowitsch PetscCall(PetscSFReduceEnd(sf, MPI_SIGNED_CHAR, leafdata, rootdata, mop)); 307c4762a1bSJed Brown 3089566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Reduce Leafdata in type of signed char\n")); 3099371c9d4SSatish Balay len = 0; 3109371c9d4SSatish Balay PetscCall(PetscSNPrintf(buf, 256, "%4d:", rank)); 3119371c9d4SSatish Balay len += 5; 3129371c9d4SSatish Balay for (i = 0; i < nleavesalloc; i++) { 3139371c9d4SSatish Balay PetscCall(PetscSNPrintf(buf + len, 256 - len, "%5d", leafdata[i])); 3149371c9d4SSatish Balay len += 5; 3159371c9d4SSatish Balay } 3169566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%s\n", buf)); 3179566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT)); 318c4762a1bSJed Brown 3199566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Reduce Rootdata in type of signed char\n")); 3209371c9d4SSatish Balay len = 0; 3219371c9d4SSatish Balay PetscCall(PetscSNPrintf(buf, 256, "%4d:", rank)); 3229371c9d4SSatish Balay len += 5; 3239371c9d4SSatish Balay for (i = 0; i < nrootsalloc; i++) { 3249371c9d4SSatish Balay PetscCall(PetscSNPrintf(buf + len, 256 - len, "%5d", rootdata[i])); 3259371c9d4SSatish Balay len += 5; 3269371c9d4SSatish Balay } 3279566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%s\n", buf)); 3289566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT)); 329c4762a1bSJed Brown 3309566063dSJacob Faibussowitsch PetscCall(PetscFree2(rootdata, leafdata)); 331c4762a1bSJed Brown } 332c4762a1bSJed Brown 333c4762a1bSJed Brown if (test_reduce && test_char) { /* Reduce with unsigned char */ 334c4762a1bSJed Brown PetscInt len; 335c4762a1bSJed Brown char buf[256]; 336c4762a1bSJed Brown unsigned char *rootdata, *leafdata; 3379566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(nrootsalloc, &rootdata, nleavesalloc, &leafdata)); 338c4762a1bSJed Brown /* Initialize rootdata buffer in which the result of the reduction will appear. */ 339b9d03b0cSStefano Zampini for (i = 0; i < nrootsalloc; i++) rootdata[i] = 0; 340*6497c311SBarry Smith for (i = 0; i < nroots; i++) rootdata[i * stride] = (unsigned char)(10 * (rank + 1) + i); 341c4762a1bSJed Brown /* Set leaf values to reduce. */ 342b9d03b0cSStefano Zampini for (i = 0; i < nleavesalloc; i++) leafdata[i] = 0; 343*6497c311SBarry Smith for (i = 0; i < nleaves; i++) leafdata[i * stride] = (unsigned char)(50 * (rank + 1) + 10 * i); 3449566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Pre-Reduce Rootdata in type of unsigned char\n")); 345c4762a1bSJed Brown 3469371c9d4SSatish Balay len = 0; 3479371c9d4SSatish Balay PetscCall(PetscSNPrintf(buf, 256, "%4d:", rank)); 3489371c9d4SSatish Balay len += 5; 3499371c9d4SSatish Balay for (i = 0; i < nrootsalloc; i++) { 3509371c9d4SSatish Balay PetscCall(PetscSNPrintf(buf + len, 256 - len, "%5u", rootdata[i])); 3519371c9d4SSatish Balay len += 5; 3529371c9d4SSatish Balay } 3539566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%s\n", buf)); 3549566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT)); 355c4762a1bSJed Brown 356c4762a1bSJed Brown /* Using MPI_CHAR should trigger an error since MPI standard does not support reduction on MPI_CHAR. 357c4762a1bSJed Brown Testing with -test_op max, one can see the sign does take effect in MPI_MAX. 358c4762a1bSJed Brown */ 3599566063dSJacob Faibussowitsch PetscCall(PetscSFReduceBegin(sf, MPI_UNSIGNED_CHAR, leafdata, rootdata, mop)); 3609566063dSJacob Faibussowitsch PetscCall(PetscSFReduceEnd(sf, MPI_UNSIGNED_CHAR, leafdata, rootdata, mop)); 361c4762a1bSJed Brown 3629566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Reduce Leafdata in type of unsigned char\n")); 3639371c9d4SSatish Balay len = 0; 3649371c9d4SSatish Balay PetscCall(PetscSNPrintf(buf, 256, "%4d:", rank)); 3659371c9d4SSatish Balay len += 5; 3669371c9d4SSatish Balay for (i = 0; i < nleavesalloc; i++) { 3679371c9d4SSatish Balay PetscCall(PetscSNPrintf(buf + len, 256 - len, "%5u", leafdata[i])); 3689371c9d4SSatish Balay len += 5; 3699371c9d4SSatish Balay } 3709566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%s\n", buf)); 3719566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT)); 372c4762a1bSJed Brown 3739566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Reduce Rootdata in type of unsigned char\n")); 3749371c9d4SSatish Balay len = 0; 3759371c9d4SSatish Balay PetscCall(PetscSNPrintf(buf, 256, "%4d:", rank)); 3769371c9d4SSatish Balay len += 5; 3779371c9d4SSatish Balay for (i = 0; i < nrootsalloc; i++) { 3789371c9d4SSatish Balay PetscCall(PetscSNPrintf(buf + len, 256 - len, "%5u", rootdata[i])); 3799371c9d4SSatish Balay len += 5; 3809371c9d4SSatish Balay } 3819566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%s\n", buf)); 3829566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT)); 383c4762a1bSJed Brown 3849566063dSJacob Faibussowitsch PetscCall(PetscFree2(rootdata, leafdata)); 385c4762a1bSJed Brown } 386c4762a1bSJed Brown 387c4762a1bSJed Brown if (test_degree) { 388c4762a1bSJed Brown const PetscInt *degree; 3899566063dSJacob Faibussowitsch PetscCall(PetscSFComputeDegreeBegin(sf, °ree)); 3909566063dSJacob Faibussowitsch PetscCall(PetscSFComputeDegreeEnd(sf, °ree)); 3919566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Root degrees\n")); 3929566063dSJacob Faibussowitsch PetscCall(PetscIntView(nrootsalloc, degree, PETSC_VIEWER_STDOUT_WORLD)); 393c4762a1bSJed Brown } 394c4762a1bSJed Brown 395c4762a1bSJed Brown if (test_fetchandop) { 396c4762a1bSJed Brown /* Cannot use text compare here because token ordering is not deterministic */ 397c4762a1bSJed Brown PetscInt *leafdata, *leafupdate, *rootdata; 3989566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(nleavesalloc, &leafdata, nleavesalloc, &leafupdate, nrootsalloc, &rootdata)); 399c4762a1bSJed Brown for (i = 0; i < nleavesalloc; i++) leafdata[i] = -1; 400c4762a1bSJed Brown for (i = 0; i < nleaves; i++) leafdata[i * stride] = 1; 401c4762a1bSJed Brown for (i = 0; i < nrootsalloc; i++) rootdata[i] = -1; 402c4762a1bSJed Brown for (i = 0; i < nroots; i++) rootdata[i * stride] = 0; 4039566063dSJacob Faibussowitsch PetscCall(PetscSFFetchAndOpBegin(sf, MPIU_INT, rootdata, leafdata, leafupdate, mop)); 4049566063dSJacob Faibussowitsch PetscCall(PetscSFFetchAndOpEnd(sf, MPIU_INT, rootdata, leafdata, leafupdate, mop)); 4059566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Rootdata (sum of 1 from each leaf)\n")); 4069566063dSJacob Faibussowitsch PetscCall(PetscIntView(nrootsalloc, rootdata, PETSC_VIEWER_STDOUT_WORLD)); 4079566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Leafupdate (value at roots prior to my atomic update)\n")); 4089566063dSJacob Faibussowitsch PetscCall(PetscIntView(nleavesalloc, leafupdate, PETSC_VIEWER_STDOUT_WORLD)); 4099566063dSJacob Faibussowitsch PetscCall(PetscFree3(leafdata, leafupdate, rootdata)); 410c4762a1bSJed Brown } 411c4762a1bSJed Brown 412c4762a1bSJed Brown if (test_gather) { 413c4762a1bSJed Brown const PetscInt *degree; 414c4762a1bSJed Brown PetscInt inedges, *indata, *outdata; 4159566063dSJacob Faibussowitsch PetscCall(PetscSFComputeDegreeBegin(sf, °ree)); 4169566063dSJacob Faibussowitsch PetscCall(PetscSFComputeDegreeEnd(sf, °ree)); 417c4762a1bSJed Brown for (i = 0, inedges = 0; i < nrootsalloc; i++) inedges += degree[i]; 4189566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(inedges, &indata, nleavesalloc, &outdata)); 419c4762a1bSJed Brown for (i = 0; i < nleavesalloc; i++) outdata[i] = -1; 420c4762a1bSJed Brown for (i = 0; i < nleaves; i++) outdata[i * stride] = 1000 * (rank + 1) + i; 4219566063dSJacob Faibussowitsch PetscCall(PetscSFGatherBegin(sf, MPIU_INT, outdata, indata)); 4229566063dSJacob Faibussowitsch PetscCall(PetscSFGatherEnd(sf, MPIU_INT, outdata, indata)); 4239566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Gathered data at multi-roots from leaves\n")); 4249566063dSJacob Faibussowitsch PetscCall(PetscIntView(inedges, indata, PETSC_VIEWER_STDOUT_WORLD)); 4259566063dSJacob Faibussowitsch PetscCall(PetscFree2(indata, outdata)); 426c4762a1bSJed Brown } 427c4762a1bSJed Brown 428c4762a1bSJed Brown if (test_scatter) { 429c4762a1bSJed Brown const PetscInt *degree; 430c4762a1bSJed Brown PetscInt j, count, inedges, *indata, *outdata; 4319566063dSJacob Faibussowitsch PetscCall(PetscSFComputeDegreeBegin(sf, °ree)); 4329566063dSJacob Faibussowitsch PetscCall(PetscSFComputeDegreeEnd(sf, °ree)); 433c4762a1bSJed Brown for (i = 0, inedges = 0; i < nrootsalloc; i++) inedges += degree[i]; 4349566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(inedges, &indata, nleavesalloc, &outdata)); 435c4762a1bSJed Brown for (i = 0; i < nleavesalloc; i++) outdata[i] = -1; 436c4762a1bSJed Brown for (i = 0, count = 0; i < nrootsalloc; i++) { 437c4762a1bSJed Brown for (j = 0; j < degree[i]; j++) indata[count++] = 1000 * (rank + 1) + 100 * i + j; 438c4762a1bSJed Brown } 4399566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Data at multi-roots, to scatter to leaves\n")); 4409566063dSJacob Faibussowitsch PetscCall(PetscIntView(inedges, indata, PETSC_VIEWER_STDOUT_WORLD)); 441c4762a1bSJed Brown 4429566063dSJacob Faibussowitsch PetscCall(PetscSFScatterBegin(sf, MPIU_INT, indata, outdata)); 4439566063dSJacob Faibussowitsch PetscCall(PetscSFScatterEnd(sf, MPIU_INT, indata, outdata)); 4449566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Scattered data at leaves\n")); 4459566063dSJacob Faibussowitsch PetscCall(PetscIntView(nleavesalloc, outdata, PETSC_VIEWER_STDOUT_WORLD)); 4469566063dSJacob Faibussowitsch PetscCall(PetscFree2(indata, outdata)); 447c4762a1bSJed Brown } 448c4762a1bSJed Brown 449c4762a1bSJed Brown if (test_embed) { 450dd400576SPatrick Sanan const PetscInt nroots = 1 + (PetscInt)(rank == 0); 451c4762a1bSJed Brown PetscInt selected[2]; 452c4762a1bSJed Brown PetscSF esf; 453c4762a1bSJed Brown 454c4762a1bSJed Brown selected[0] = stride; 455c4762a1bSJed Brown selected[1] = 2 * stride; 4569566063dSJacob Faibussowitsch PetscCall(PetscSFCreateEmbeddedRootSF(sf, nroots, selected, &esf)); 4579566063dSJacob Faibussowitsch PetscCall(PetscSFSetUp(esf)); 4589566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Embedded PetscSF\n")); 4599566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_ASCII_INFO_DETAIL)); 4609566063dSJacob Faibussowitsch PetscCall(PetscSFView(esf, PETSC_VIEWER_STDOUT_WORLD)); 4619566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD)); 4629566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&esf)); 463c4762a1bSJed Brown } 464c4762a1bSJed Brown 465c4762a1bSJed Brown if (test_invert) { 466c4762a1bSJed Brown const PetscInt *degree; 467c4762a1bSJed Brown PetscInt *mRootsOrigNumbering; 468c4762a1bSJed Brown PetscInt inedges; 469c4762a1bSJed Brown PetscSF msf, imsf; 470c4762a1bSJed Brown 4719566063dSJacob Faibussowitsch PetscCall(PetscSFGetMultiSF(sf, &msf)); 4729566063dSJacob Faibussowitsch PetscCall(PetscSFCreateInverseSF(msf, &imsf)); 4739566063dSJacob Faibussowitsch PetscCall(PetscSFSetUp(msf)); 4749566063dSJacob Faibussowitsch PetscCall(PetscSFSetUp(imsf)); 4759566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Multi-SF\n")); 4769566063dSJacob Faibussowitsch PetscCall(PetscSFView(msf, PETSC_VIEWER_STDOUT_WORLD)); 4779566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Multi-SF roots indices in original SF roots numbering\n")); 4789566063dSJacob Faibussowitsch PetscCall(PetscSFComputeDegreeBegin(sf, °ree)); 4799566063dSJacob Faibussowitsch PetscCall(PetscSFComputeDegreeEnd(sf, °ree)); 4809566063dSJacob Faibussowitsch PetscCall(PetscSFComputeMultiRootOriginalNumbering(sf, degree, &inedges, &mRootsOrigNumbering)); 4819566063dSJacob Faibussowitsch PetscCall(PetscIntView(inedges, mRootsOrigNumbering, PETSC_VIEWER_STDOUT_WORLD)); 4829566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Inverse of Multi-SF\n")); 4839566063dSJacob Faibussowitsch PetscCall(PetscSFView(imsf, PETSC_VIEWER_STDOUT_WORLD)); 4849566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Inverse of Multi-SF, original numbering\n")); 4859566063dSJacob Faibussowitsch PetscCall(PetscSFViewCustomLocals_Private(imsf, mRootsOrigNumbering, PETSC_VIEWER_STDOUT_WORLD)); 4869566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&imsf)); 4879566063dSJacob Faibussowitsch PetscCall(PetscFree(mRootsOrigNumbering)); 488c4762a1bSJed Brown } 489c4762a1bSJed Brown 490c4762a1bSJed Brown /* Clean storage for star forest. */ 4919566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sf)); 4929566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 493b122ec5aSJacob Faibussowitsch return 0; 494c4762a1bSJed Brown } 495c4762a1bSJed Brown 496c4762a1bSJed Brown /*TEST 497c4762a1bSJed Brown 498c4762a1bSJed Brown test: 499c4762a1bSJed Brown nsize: 4 500c4762a1bSJed Brown filter: grep -v "type" | grep -v "sort" 501c4762a1bSJed Brown args: -test_bcast -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}} 502dfd57a17SPierre Jolivet requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW) 503c4762a1bSJed Brown 504c4762a1bSJed Brown test: 505c4762a1bSJed Brown suffix: 2 506c4762a1bSJed Brown nsize: 4 507c4762a1bSJed Brown filter: grep -v "type" | grep -v "sort" 508c4762a1bSJed Brown args: -test_reduce -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}} 509dfd57a17SPierre Jolivet requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW) 510c4762a1bSJed Brown 511c4762a1bSJed Brown test: 512c4762a1bSJed Brown suffix: 2_basic 513c4762a1bSJed Brown nsize: 4 514c4762a1bSJed Brown args: -test_reduce -sf_type basic 515c4762a1bSJed Brown 516c4762a1bSJed Brown test: 517c4762a1bSJed Brown suffix: 3 518c4762a1bSJed Brown nsize: 4 519c4762a1bSJed Brown filter: grep -v "type" | grep -v "sort" 520c4762a1bSJed Brown args: -test_degree -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}} 521dfd57a17SPierre Jolivet requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW) 522c4762a1bSJed Brown 523c4762a1bSJed Brown test: 524c4762a1bSJed Brown suffix: 3_basic 525c4762a1bSJed Brown nsize: 4 526c4762a1bSJed Brown args: -test_degree -sf_type basic 527c4762a1bSJed Brown 528c4762a1bSJed Brown test: 529c4762a1bSJed Brown suffix: 4 530c4762a1bSJed Brown nsize: 4 531c4762a1bSJed Brown filter: grep -v "type" | grep -v "sort" 532c4762a1bSJed Brown args: -test_gather -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}} 533dfd57a17SPierre Jolivet requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW) 534c4762a1bSJed Brown 535c4762a1bSJed Brown test: 536c4762a1bSJed Brown suffix: 4_basic 537c4762a1bSJed Brown nsize: 4 538c4762a1bSJed Brown args: -test_gather -sf_type basic 539c4762a1bSJed Brown 540c4762a1bSJed Brown test: 541c4762a1bSJed Brown suffix: 4_stride 542c4762a1bSJed Brown nsize: 4 543c4762a1bSJed Brown args: -test_gather -sf_type basic -stride 2 544c4762a1bSJed Brown 545c4762a1bSJed Brown test: 546c4762a1bSJed Brown suffix: 5 547c4762a1bSJed Brown nsize: 4 548c4762a1bSJed Brown filter: grep -v "type" | grep -v "sort" 549c4762a1bSJed Brown args: -test_scatter -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}} 550dfd57a17SPierre Jolivet requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW) 551c4762a1bSJed Brown 552c4762a1bSJed Brown test: 553c4762a1bSJed Brown suffix: 5_basic 554c4762a1bSJed Brown nsize: 4 555c4762a1bSJed Brown args: -test_scatter -sf_type basic 556c4762a1bSJed Brown 557c4762a1bSJed Brown test: 558c4762a1bSJed Brown suffix: 5_stride 559c4762a1bSJed Brown nsize: 4 560c4762a1bSJed Brown args: -test_scatter -sf_type basic -stride 2 561c4762a1bSJed Brown 562c4762a1bSJed Brown test: 563c4762a1bSJed Brown suffix: 6 564c4762a1bSJed Brown nsize: 4 565c4762a1bSJed Brown filter: grep -v "type" | grep -v "sort" 566c20d7725SJed Brown # No -sf_window_flavor dynamic due to bug https://gitlab.com/petsc/petsc/issues/555 567c20d7725SJed Brown args: -test_embed -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create allocate}} 568dfd57a17SPierre Jolivet requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW) 569c4762a1bSJed Brown 570c4762a1bSJed Brown test: 571c4762a1bSJed Brown suffix: 6_basic 572c4762a1bSJed Brown nsize: 4 573c4762a1bSJed Brown args: -test_embed -sf_type basic 574c4762a1bSJed Brown 575c4762a1bSJed Brown test: 576c4762a1bSJed Brown suffix: 7 577c4762a1bSJed Brown nsize: 4 578c4762a1bSJed Brown filter: grep -v "type" | grep -v "sort" 579c4762a1bSJed Brown args: -test_invert -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}} 580dfd57a17SPierre Jolivet requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW) 581c4762a1bSJed Brown 582c4762a1bSJed Brown test: 583c4762a1bSJed Brown suffix: 7_basic 584c4762a1bSJed Brown nsize: 4 585c4762a1bSJed Brown args: -test_invert -sf_type basic 586c4762a1bSJed Brown 587c4762a1bSJed Brown test: 588c4762a1bSJed Brown suffix: basic 589c4762a1bSJed Brown nsize: 4 590c4762a1bSJed Brown args: -test_bcast -sf_type basic 591c4762a1bSJed Brown output_file: output/ex1_1_basic.out 592c4762a1bSJed Brown 593c4762a1bSJed Brown test: 594c4762a1bSJed Brown suffix: bcastop_basic 595c4762a1bSJed Brown nsize: 4 596c4762a1bSJed Brown args: -test_bcastop -sf_type basic 597c4762a1bSJed Brown output_file: output/ex1_bcastop_basic.out 598c4762a1bSJed Brown 599c4762a1bSJed Brown test: 600c4762a1bSJed Brown suffix: 8 601c4762a1bSJed Brown nsize: 3 602c4762a1bSJed Brown filter: grep -v "type" | grep -v "sort" 603c4762a1bSJed Brown args: -test_bcast -test_sf_distribute -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}} 604dfd57a17SPierre Jolivet requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW) 605c4762a1bSJed Brown 606c4762a1bSJed Brown test: 607c4762a1bSJed Brown suffix: 8_basic 608c4762a1bSJed Brown nsize: 3 609c4762a1bSJed Brown args: -test_bcast -test_sf_distribute -sf_type basic 610c4762a1bSJed Brown 611c4762a1bSJed Brown test: 612c4762a1bSJed Brown suffix: 9_char 613c4762a1bSJed Brown nsize: 4 614c4762a1bSJed Brown args: -sf_type basic -test_bcast -test_reduce -test_op max -test_char 615c4762a1bSJed Brown 616c4762a1bSJed Brown # Here we do not test -sf_window_flavor dynamic since it is designed for repeated SFs with few different rootdata pointers 617c4762a1bSJed Brown test: 618c4762a1bSJed Brown suffix: 10 619c4762a1bSJed Brown filter: grep -v "type" | grep -v "sort" 620c4762a1bSJed Brown nsize: 4 621c4762a1bSJed Brown args: -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create allocate}} -test_all -test_bcastop 0 -test_fetchandop 0 622dfd57a17SPierre Jolivet requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW) 623c4762a1bSJed Brown 624c4762a1bSJed Brown # The nightly test suite with MPICH uses ch3:sock, which is broken when winsize == 0 in some of the processes 625c4762a1bSJed Brown test: 626c4762a1bSJed Brown suffix: 10_shared 627c4762a1bSJed Brown output_file: output/ex1_10.out 628c4762a1bSJed Brown filter: grep -v "type" | grep -v "sort" 629c4762a1bSJed Brown nsize: 4 630c4762a1bSJed Brown args: -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor shared -test_all -test_bcastop 0 -test_fetchandop 0 631e3c15826SJunchao Zhang requires: defined(PETSC_HAVE_MPI_PROCESS_SHARED_MEMORY) !defined(PETSC_HAVE_MPICH) defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW) 632c4762a1bSJed Brown 633c4762a1bSJed Brown test: 634c4762a1bSJed Brown suffix: 10_basic 635c4762a1bSJed Brown nsize: 4 636c4762a1bSJed Brown args: -sf_type basic -test_all -test_bcastop 0 -test_fetchandop 0 637c4762a1bSJed Brown 6380dd791a8SStefano Zampini test: 6390dd791a8SStefano Zampini suffix: 10_basic_vector 6400dd791a8SStefano Zampini nsize: 4 6410dd791a8SStefano Zampini args: -sf_type basic -test_all -test_bcastop 0 -test_fetchandop 0 -test_vector 6420dd791a8SStefano Zampini 643c4762a1bSJed Brown TEST*/ 644