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; 21c4762a1bSJed Brown PetscInt i, nroots, nleaves, nranks; 22c4762a1bSJed Brown PetscMPIInt rank; 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)); 309566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] Number of roots=%" PetscInt_FMT ", leaves=%" PetscInt_FMT ", remote ranks=%" PetscInt_FMT "\n", rank, nroots, nleaves, nranks)); 3148a46eb9SPierre Jolivet for (i = 0; i < nleaves; i++) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] %" PetscInt_FMT " <- (%" PetscInt_FMT ",%" PetscInt_FMT ")\n", rank, locals[i], 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 { 40*0dd791a8SStefano Zampini PetscInt i, bs = 2, nroots, nrootsalloc, nleaves, nleavesalloc, *mine, stride; 41c4762a1bSJed Brown PetscSFNode *remote; 42c4762a1bSJed Brown PetscMPIInt rank, size; 43*0dd791a8SStefano Zampini PetscSF sf, vsf; 44*0dd791a8SStefano 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)); 105*0dd791a8SStefano Zampini PetscCall(PetscOptionsInt("-bs", "Block size for vectorial SF", "", bs, &bs, NULL)); 106*0dd791a8SStefano 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)); 116c4762a1bSJed Brown for (i = 0; i < size; i++) { 117c4762a1bSJed Brown remote[i].rank = i; 118c4762a1bSJed Brown remote[i].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 145*0dd791a8SStefano 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 151*0dd791a8SStefano Zampini /* We can also obtain a strided SF to communicate interleaved blocks of data */ 152*0dd791a8SStefano Zampini PetscCall(PetscSFCreateStridedSF(sf, bs, PETSC_DECIDE, PETSC_DECIDE, &vsf)); 153*0dd791a8SStefano Zampini if (test_vector) { /* perform all tests below using the vectorial SF */ 154*0dd791a8SStefano Zampini PetscSF t = sf; 155*0dd791a8SStefano Zampini 156*0dd791a8SStefano Zampini sf = vsf; 157*0dd791a8SStefano Zampini vsf = t; 158*0dd791a8SStefano Zampini 159*0dd791a8SStefano Zampini nroots *= bs; 160*0dd791a8SStefano Zampini nleaves *= bs; 161*0dd791a8SStefano Zampini nrootsalloc *= bs; 162*0dd791a8SStefano Zampini nleavesalloc *= bs; 163*0dd791a8SStefano Zampini } 164*0dd791a8SStefano Zampini PetscCall(PetscSFDestroy(&vsf)); 165*0dd791a8SStefano 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. */ 1829566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(sf, MPIU_INT, rootdata, leafdata, MPI_REPLACE)); 1839566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(sf, MPIU_INT, rootdata, leafdata, MPI_REPLACE)); 1849566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Bcast Rootdata\n")); 1859566063dSJacob Faibussowitsch PetscCall(PetscIntView(nrootsalloc, rootdata, PETSC_VIEWER_STDOUT_WORLD)); 1869566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Bcast Leafdata\n")); 1879566063dSJacob Faibussowitsch PetscCall(PetscIntView(nleavesalloc, leafdata, PETSC_VIEWER_STDOUT_WORLD)); 1889566063dSJacob Faibussowitsch PetscCall(PetscFree2(rootdata, leafdata)); 189c4762a1bSJed Brown } 190c4762a1bSJed Brown 191c4762a1bSJed Brown if (test_bcast && test_char) { /* Bcast with char */ 192c4762a1bSJed Brown PetscInt len; 193c4762a1bSJed Brown char buf[256]; 194c4762a1bSJed Brown char *rootdata, *leafdata; 1959566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(nrootsalloc, &rootdata, nleavesalloc, &leafdata)); 196c4762a1bSJed Brown /* Set rootdata buffer to be broadcast */ 197c4762a1bSJed Brown for (i = 0; i < nrootsalloc; i++) rootdata[i] = '*'; 198c4762a1bSJed Brown for (i = 0; i < nroots; i++) rootdata[i * stride] = 'A' + rank * 3 + i; /* rank is very small, so it is fine to compute a char */ 199c4762a1bSJed Brown /* Initialize local buffer, these values are never used. */ 200c4762a1bSJed Brown for (i = 0; i < nleavesalloc; i++) leafdata[i] = '?'; 201c4762a1bSJed Brown 2029566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(sf, MPI_CHAR, rootdata, leafdata, MPI_REPLACE)); 2039566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(sf, MPI_CHAR, rootdata, leafdata, MPI_REPLACE)); 204c4762a1bSJed Brown 2059566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Bcast Rootdata in type of char\n")); 2069371c9d4SSatish Balay len = 0; 2079371c9d4SSatish Balay PetscCall(PetscSNPrintf(buf, 256, "%4d:", rank)); 2089371c9d4SSatish Balay len += 5; 2099371c9d4SSatish Balay for (i = 0; i < nrootsalloc; i++) { 2109371c9d4SSatish Balay PetscCall(PetscSNPrintf(buf + len, 256 - len, "%5c", rootdata[i])); 2119371c9d4SSatish Balay len += 5; 2129371c9d4SSatish Balay } 2139566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%s\n", buf)); 2149566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT)); 215c4762a1bSJed Brown 2169566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Bcast Leafdata in type of char\n")); 2179371c9d4SSatish Balay len = 0; 2189371c9d4SSatish Balay PetscCall(PetscSNPrintf(buf, 256, "%4d:", rank)); 2199371c9d4SSatish Balay len += 5; 2209371c9d4SSatish Balay for (i = 0; i < nleavesalloc; i++) { 2219371c9d4SSatish Balay PetscCall(PetscSNPrintf(buf + len, 256 - len, "%5c", leafdata[i])); 2229371c9d4SSatish Balay len += 5; 2239371c9d4SSatish Balay } 2249566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%s\n", buf)); 2259566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT)); 226c4762a1bSJed Brown 2279566063dSJacob Faibussowitsch PetscCall(PetscFree2(rootdata, leafdata)); 228c4762a1bSJed Brown } 229c4762a1bSJed Brown 230c4762a1bSJed Brown if (test_bcastop) { /* Reduce rootdata into leafdata */ 231c4762a1bSJed Brown PetscInt *rootdata, *leafdata; 232c4762a1bSJed Brown /* Allocate space for send and receive buffers. This example communicates PetscInt, but other types, including 233c4762a1bSJed Brown * user-defined structures, could also be used. */ 2349566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(nrootsalloc, &rootdata, nleavesalloc, &leafdata)); 235c4762a1bSJed Brown /* Set rootdata buffer to be broadcast */ 236c4762a1bSJed Brown for (i = 0; i < nrootsalloc; i++) rootdata[i] = -1; 237c4762a1bSJed Brown for (i = 0; i < nroots; i++) rootdata[i * stride] = 100 * (rank + 1) + i; 238c4762a1bSJed Brown /* Set leaf values to reduce with */ 239c4762a1bSJed Brown for (i = 0; i < nleavesalloc; i++) leafdata[i] = -10 * (rank + 1) - i; 2409566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Pre-BcastAndOp Leafdata\n")); 2419566063dSJacob Faibussowitsch PetscCall(PetscIntView(nleavesalloc, leafdata, PETSC_VIEWER_STDOUT_WORLD)); 242c4762a1bSJed Brown /* Broadcast entries from rootdata to leafdata. Computation or other communication can be performed between the begin and end calls. */ 2439566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(sf, MPIU_INT, rootdata, leafdata, mop)); 2449566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(sf, MPIU_INT, rootdata, leafdata, mop)); 2459566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## BcastAndOp Rootdata\n")); 2469566063dSJacob Faibussowitsch PetscCall(PetscIntView(nrootsalloc, rootdata, PETSC_VIEWER_STDOUT_WORLD)); 2479566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## BcastAndOp Leafdata\n")); 2489566063dSJacob Faibussowitsch PetscCall(PetscIntView(nleavesalloc, leafdata, PETSC_VIEWER_STDOUT_WORLD)); 2499566063dSJacob Faibussowitsch PetscCall(PetscFree2(rootdata, leafdata)); 250c4762a1bSJed Brown } 251c4762a1bSJed Brown 252c4762a1bSJed Brown if (test_reduce) { /* Reduce leafdata into rootdata */ 253c4762a1bSJed Brown PetscInt *rootdata, *leafdata; 2549566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(nrootsalloc, &rootdata, nleavesalloc, &leafdata)); 255c4762a1bSJed Brown /* Initialize rootdata buffer in which the result of the reduction will appear. */ 256c4762a1bSJed Brown for (i = 0; i < nrootsalloc; i++) rootdata[i] = -1; 257c4762a1bSJed Brown for (i = 0; i < nroots; i++) rootdata[i * stride] = 100 * (rank + 1) + i; 258c4762a1bSJed Brown /* Set leaf values to reduce. */ 259c4762a1bSJed Brown for (i = 0; i < nleavesalloc; i++) leafdata[i] = -1; 260c4762a1bSJed Brown for (i = 0; i < nleaves; i++) leafdata[i * stride] = 1000 * (rank + 1) + 10 * i; 2619566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Pre-Reduce Rootdata\n")); 2629566063dSJacob Faibussowitsch PetscCall(PetscIntView(nrootsalloc, rootdata, PETSC_VIEWER_STDOUT_WORLD)); 263c4762a1bSJed Brown /* Perform reduction. Computation or other communication can be performed between the begin and end calls. 264c4762a1bSJed Brown * This example sums the values, but other MPI_Ops can be used (e.g MPI_MAX, MPI_PROD). */ 2659566063dSJacob Faibussowitsch PetscCall(PetscSFReduceBegin(sf, MPIU_INT, leafdata, rootdata, mop)); 2669566063dSJacob Faibussowitsch PetscCall(PetscSFReduceEnd(sf, MPIU_INT, leafdata, rootdata, mop)); 2679566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Reduce Leafdata\n")); 2689566063dSJacob Faibussowitsch PetscCall(PetscIntView(nleavesalloc, leafdata, PETSC_VIEWER_STDOUT_WORLD)); 2699566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Reduce Rootdata\n")); 2709566063dSJacob Faibussowitsch PetscCall(PetscIntView(nrootsalloc, rootdata, PETSC_VIEWER_STDOUT_WORLD)); 2719566063dSJacob Faibussowitsch PetscCall(PetscFree2(rootdata, leafdata)); 272c4762a1bSJed Brown } 273c4762a1bSJed Brown 274c4762a1bSJed Brown if (test_reduce && test_char) { /* Reduce with signed char */ 275c4762a1bSJed Brown PetscInt len; 276c4762a1bSJed Brown char buf[256]; 277c4762a1bSJed Brown signed char *rootdata, *leafdata; 2789566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(nrootsalloc, &rootdata, nleavesalloc, &leafdata)); 279c4762a1bSJed Brown /* Initialize rootdata buffer in which the result of the reduction will appear. */ 280c4762a1bSJed Brown for (i = 0; i < nrootsalloc; i++) rootdata[i] = -1; 281c4762a1bSJed Brown for (i = 0; i < nroots; i++) rootdata[i * stride] = 10 * (rank + 1) + i; 282c4762a1bSJed Brown /* Set leaf values to reduce. */ 283c4762a1bSJed Brown for (i = 0; i < nleavesalloc; i++) leafdata[i] = -1; 284c4762a1bSJed Brown for (i = 0; i < nleaves; i++) leafdata[i * stride] = 50 * (rank + 1) + 10 * i; 2859566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Pre-Reduce Rootdata in type of signed char\n")); 286c4762a1bSJed Brown 2879371c9d4SSatish Balay len = 0; 2889371c9d4SSatish Balay PetscCall(PetscSNPrintf(buf, 256, "%4d:", rank)); 2899371c9d4SSatish Balay len += 5; 2909371c9d4SSatish Balay for (i = 0; i < nrootsalloc; i++) { 2919371c9d4SSatish Balay PetscCall(PetscSNPrintf(buf + len, 256 - len, "%5d", rootdata[i])); 2929371c9d4SSatish Balay len += 5; 2939371c9d4SSatish Balay } 2949566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%s\n", buf)); 2959566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT)); 296c4762a1bSJed Brown 297c4762a1bSJed Brown /* Using MPI_CHAR should trigger an error since MPI standard does not support reduction on MPI_CHAR. 298c4762a1bSJed Brown Testing with -test_op max, one can see the sign does take effect in MPI_MAX. 299c4762a1bSJed Brown */ 3009566063dSJacob Faibussowitsch PetscCall(PetscSFReduceBegin(sf, MPI_SIGNED_CHAR, leafdata, rootdata, mop)); 3019566063dSJacob Faibussowitsch PetscCall(PetscSFReduceEnd(sf, MPI_SIGNED_CHAR, leafdata, rootdata, mop)); 302c4762a1bSJed Brown 3039566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Reduce Leafdata in type of signed char\n")); 3049371c9d4SSatish Balay len = 0; 3059371c9d4SSatish Balay PetscCall(PetscSNPrintf(buf, 256, "%4d:", rank)); 3069371c9d4SSatish Balay len += 5; 3079371c9d4SSatish Balay for (i = 0; i < nleavesalloc; i++) { 3089371c9d4SSatish Balay PetscCall(PetscSNPrintf(buf + len, 256 - len, "%5d", leafdata[i])); 3099371c9d4SSatish Balay len += 5; 3109371c9d4SSatish Balay } 3119566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%s\n", buf)); 3129566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT)); 313c4762a1bSJed Brown 3149566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Reduce Rootdata in type of signed char\n")); 3159371c9d4SSatish Balay len = 0; 3169371c9d4SSatish Balay PetscCall(PetscSNPrintf(buf, 256, "%4d:", rank)); 3179371c9d4SSatish Balay len += 5; 3189371c9d4SSatish Balay for (i = 0; i < nrootsalloc; i++) { 3199371c9d4SSatish Balay PetscCall(PetscSNPrintf(buf + len, 256 - len, "%5d", rootdata[i])); 3209371c9d4SSatish Balay len += 5; 3219371c9d4SSatish Balay } 3229566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%s\n", buf)); 3239566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT)); 324c4762a1bSJed Brown 3259566063dSJacob Faibussowitsch PetscCall(PetscFree2(rootdata, leafdata)); 326c4762a1bSJed Brown } 327c4762a1bSJed Brown 328c4762a1bSJed Brown if (test_reduce && test_char) { /* Reduce with unsigned char */ 329c4762a1bSJed Brown PetscInt len; 330c4762a1bSJed Brown char buf[256]; 331c4762a1bSJed Brown unsigned char *rootdata, *leafdata; 3329566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(nrootsalloc, &rootdata, nleavesalloc, &leafdata)); 333c4762a1bSJed Brown /* Initialize rootdata buffer in which the result of the reduction will appear. */ 334b9d03b0cSStefano Zampini for (i = 0; i < nrootsalloc; i++) rootdata[i] = 0; 335c4762a1bSJed Brown for (i = 0; i < nroots; i++) rootdata[i * stride] = 10 * (rank + 1) + i; 336c4762a1bSJed Brown /* Set leaf values to reduce. */ 337b9d03b0cSStefano Zampini for (i = 0; i < nleavesalloc; i++) leafdata[i] = 0; 338c4762a1bSJed Brown for (i = 0; i < nleaves; i++) leafdata[i * stride] = 50 * (rank + 1) + 10 * i; 3399566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Pre-Reduce Rootdata in type of unsigned char\n")); 340c4762a1bSJed Brown 3419371c9d4SSatish Balay len = 0; 3429371c9d4SSatish Balay PetscCall(PetscSNPrintf(buf, 256, "%4d:", rank)); 3439371c9d4SSatish Balay len += 5; 3449371c9d4SSatish Balay for (i = 0; i < nrootsalloc; i++) { 3459371c9d4SSatish Balay PetscCall(PetscSNPrintf(buf + len, 256 - len, "%5u", rootdata[i])); 3469371c9d4SSatish Balay len += 5; 3479371c9d4SSatish Balay } 3489566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%s\n", buf)); 3499566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT)); 350c4762a1bSJed Brown 351c4762a1bSJed Brown /* Using MPI_CHAR should trigger an error since MPI standard does not support reduction on MPI_CHAR. 352c4762a1bSJed Brown Testing with -test_op max, one can see the sign does take effect in MPI_MAX. 353c4762a1bSJed Brown */ 3549566063dSJacob Faibussowitsch PetscCall(PetscSFReduceBegin(sf, MPI_UNSIGNED_CHAR, leafdata, rootdata, mop)); 3559566063dSJacob Faibussowitsch PetscCall(PetscSFReduceEnd(sf, MPI_UNSIGNED_CHAR, leafdata, rootdata, mop)); 356c4762a1bSJed Brown 3579566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Reduce Leafdata in type of unsigned char\n")); 3589371c9d4SSatish Balay len = 0; 3599371c9d4SSatish Balay PetscCall(PetscSNPrintf(buf, 256, "%4d:", rank)); 3609371c9d4SSatish Balay len += 5; 3619371c9d4SSatish Balay for (i = 0; i < nleavesalloc; i++) { 3629371c9d4SSatish Balay PetscCall(PetscSNPrintf(buf + len, 256 - len, "%5u", leafdata[i])); 3639371c9d4SSatish Balay len += 5; 3649371c9d4SSatish Balay } 3659566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%s\n", buf)); 3669566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT)); 367c4762a1bSJed Brown 3689566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Reduce Rootdata in type of unsigned char\n")); 3699371c9d4SSatish Balay len = 0; 3709371c9d4SSatish Balay PetscCall(PetscSNPrintf(buf, 256, "%4d:", rank)); 3719371c9d4SSatish Balay len += 5; 3729371c9d4SSatish Balay for (i = 0; i < nrootsalloc; i++) { 3739371c9d4SSatish Balay PetscCall(PetscSNPrintf(buf + len, 256 - len, "%5u", rootdata[i])); 3749371c9d4SSatish Balay len += 5; 3759371c9d4SSatish Balay } 3769566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%s\n", buf)); 3779566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT)); 378c4762a1bSJed Brown 3799566063dSJacob Faibussowitsch PetscCall(PetscFree2(rootdata, leafdata)); 380c4762a1bSJed Brown } 381c4762a1bSJed Brown 382c4762a1bSJed Brown if (test_degree) { 383c4762a1bSJed Brown const PetscInt *degree; 3849566063dSJacob Faibussowitsch PetscCall(PetscSFComputeDegreeBegin(sf, °ree)); 3859566063dSJacob Faibussowitsch PetscCall(PetscSFComputeDegreeEnd(sf, °ree)); 3869566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Root degrees\n")); 3879566063dSJacob Faibussowitsch PetscCall(PetscIntView(nrootsalloc, degree, PETSC_VIEWER_STDOUT_WORLD)); 388c4762a1bSJed Brown } 389c4762a1bSJed Brown 390c4762a1bSJed Brown if (test_fetchandop) { 391c4762a1bSJed Brown /* Cannot use text compare here because token ordering is not deterministic */ 392c4762a1bSJed Brown PetscInt *leafdata, *leafupdate, *rootdata; 3939566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(nleavesalloc, &leafdata, nleavesalloc, &leafupdate, nrootsalloc, &rootdata)); 394c4762a1bSJed Brown for (i = 0; i < nleavesalloc; i++) leafdata[i] = -1; 395c4762a1bSJed Brown for (i = 0; i < nleaves; i++) leafdata[i * stride] = 1; 396c4762a1bSJed Brown for (i = 0; i < nrootsalloc; i++) rootdata[i] = -1; 397c4762a1bSJed Brown for (i = 0; i < nroots; i++) rootdata[i * stride] = 0; 3989566063dSJacob Faibussowitsch PetscCall(PetscSFFetchAndOpBegin(sf, MPIU_INT, rootdata, leafdata, leafupdate, mop)); 3999566063dSJacob Faibussowitsch PetscCall(PetscSFFetchAndOpEnd(sf, MPIU_INT, rootdata, leafdata, leafupdate, mop)); 4009566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Rootdata (sum of 1 from each leaf)\n")); 4019566063dSJacob Faibussowitsch PetscCall(PetscIntView(nrootsalloc, rootdata, PETSC_VIEWER_STDOUT_WORLD)); 4029566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Leafupdate (value at roots prior to my atomic update)\n")); 4039566063dSJacob Faibussowitsch PetscCall(PetscIntView(nleavesalloc, leafupdate, PETSC_VIEWER_STDOUT_WORLD)); 4049566063dSJacob Faibussowitsch PetscCall(PetscFree3(leafdata, leafupdate, rootdata)); 405c4762a1bSJed Brown } 406c4762a1bSJed Brown 407c4762a1bSJed Brown if (test_gather) { 408c4762a1bSJed Brown const PetscInt *degree; 409c4762a1bSJed Brown PetscInt inedges, *indata, *outdata; 4109566063dSJacob Faibussowitsch PetscCall(PetscSFComputeDegreeBegin(sf, °ree)); 4119566063dSJacob Faibussowitsch PetscCall(PetscSFComputeDegreeEnd(sf, °ree)); 412c4762a1bSJed Brown for (i = 0, inedges = 0; i < nrootsalloc; i++) inedges += degree[i]; 4139566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(inedges, &indata, nleavesalloc, &outdata)); 414c4762a1bSJed Brown for (i = 0; i < nleavesalloc; i++) outdata[i] = -1; 415c4762a1bSJed Brown for (i = 0; i < nleaves; i++) outdata[i * stride] = 1000 * (rank + 1) + i; 4169566063dSJacob Faibussowitsch PetscCall(PetscSFGatherBegin(sf, MPIU_INT, outdata, indata)); 4179566063dSJacob Faibussowitsch PetscCall(PetscSFGatherEnd(sf, MPIU_INT, outdata, indata)); 4189566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Gathered data at multi-roots from leaves\n")); 4199566063dSJacob Faibussowitsch PetscCall(PetscIntView(inedges, indata, PETSC_VIEWER_STDOUT_WORLD)); 4209566063dSJacob Faibussowitsch PetscCall(PetscFree2(indata, outdata)); 421c4762a1bSJed Brown } 422c4762a1bSJed Brown 423c4762a1bSJed Brown if (test_scatter) { 424c4762a1bSJed Brown const PetscInt *degree; 425c4762a1bSJed Brown PetscInt j, count, inedges, *indata, *outdata; 4269566063dSJacob Faibussowitsch PetscCall(PetscSFComputeDegreeBegin(sf, °ree)); 4279566063dSJacob Faibussowitsch PetscCall(PetscSFComputeDegreeEnd(sf, °ree)); 428c4762a1bSJed Brown for (i = 0, inedges = 0; i < nrootsalloc; i++) inedges += degree[i]; 4299566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(inedges, &indata, nleavesalloc, &outdata)); 430c4762a1bSJed Brown for (i = 0; i < nleavesalloc; i++) outdata[i] = -1; 431c4762a1bSJed Brown for (i = 0, count = 0; i < nrootsalloc; i++) { 432c4762a1bSJed Brown for (j = 0; j < degree[i]; j++) indata[count++] = 1000 * (rank + 1) + 100 * i + j; 433c4762a1bSJed Brown } 4349566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Data at multi-roots, to scatter to leaves\n")); 4359566063dSJacob Faibussowitsch PetscCall(PetscIntView(inedges, indata, PETSC_VIEWER_STDOUT_WORLD)); 436c4762a1bSJed Brown 4379566063dSJacob Faibussowitsch PetscCall(PetscSFScatterBegin(sf, MPIU_INT, indata, outdata)); 4389566063dSJacob Faibussowitsch PetscCall(PetscSFScatterEnd(sf, MPIU_INT, indata, outdata)); 4399566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Scattered data at leaves\n")); 4409566063dSJacob Faibussowitsch PetscCall(PetscIntView(nleavesalloc, outdata, PETSC_VIEWER_STDOUT_WORLD)); 4419566063dSJacob Faibussowitsch PetscCall(PetscFree2(indata, outdata)); 442c4762a1bSJed Brown } 443c4762a1bSJed Brown 444c4762a1bSJed Brown if (test_embed) { 445dd400576SPatrick Sanan const PetscInt nroots = 1 + (PetscInt)(rank == 0); 446c4762a1bSJed Brown PetscInt selected[2]; 447c4762a1bSJed Brown PetscSF esf; 448c4762a1bSJed Brown 449c4762a1bSJed Brown selected[0] = stride; 450c4762a1bSJed Brown selected[1] = 2 * stride; 4519566063dSJacob Faibussowitsch PetscCall(PetscSFCreateEmbeddedRootSF(sf, nroots, selected, &esf)); 4529566063dSJacob Faibussowitsch PetscCall(PetscSFSetUp(esf)); 4539566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Embedded PetscSF\n")); 4549566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_ASCII_INFO_DETAIL)); 4559566063dSJacob Faibussowitsch PetscCall(PetscSFView(esf, PETSC_VIEWER_STDOUT_WORLD)); 4569566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD)); 4579566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&esf)); 458c4762a1bSJed Brown } 459c4762a1bSJed Brown 460c4762a1bSJed Brown if (test_invert) { 461c4762a1bSJed Brown const PetscInt *degree; 462c4762a1bSJed Brown PetscInt *mRootsOrigNumbering; 463c4762a1bSJed Brown PetscInt inedges; 464c4762a1bSJed Brown PetscSF msf, imsf; 465c4762a1bSJed Brown 4669566063dSJacob Faibussowitsch PetscCall(PetscSFGetMultiSF(sf, &msf)); 4679566063dSJacob Faibussowitsch PetscCall(PetscSFCreateInverseSF(msf, &imsf)); 4689566063dSJacob Faibussowitsch PetscCall(PetscSFSetUp(msf)); 4699566063dSJacob Faibussowitsch PetscCall(PetscSFSetUp(imsf)); 4709566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Multi-SF\n")); 4719566063dSJacob Faibussowitsch PetscCall(PetscSFView(msf, PETSC_VIEWER_STDOUT_WORLD)); 4729566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Multi-SF roots indices in original SF roots numbering\n")); 4739566063dSJacob Faibussowitsch PetscCall(PetscSFComputeDegreeBegin(sf, °ree)); 4749566063dSJacob Faibussowitsch PetscCall(PetscSFComputeDegreeEnd(sf, °ree)); 4759566063dSJacob Faibussowitsch PetscCall(PetscSFComputeMultiRootOriginalNumbering(sf, degree, &inedges, &mRootsOrigNumbering)); 4769566063dSJacob Faibussowitsch PetscCall(PetscIntView(inedges, mRootsOrigNumbering, PETSC_VIEWER_STDOUT_WORLD)); 4779566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Inverse of Multi-SF\n")); 4789566063dSJacob Faibussowitsch PetscCall(PetscSFView(imsf, PETSC_VIEWER_STDOUT_WORLD)); 4799566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Inverse of Multi-SF, original numbering\n")); 4809566063dSJacob Faibussowitsch PetscCall(PetscSFViewCustomLocals_Private(imsf, mRootsOrigNumbering, PETSC_VIEWER_STDOUT_WORLD)); 4819566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&imsf)); 4829566063dSJacob Faibussowitsch PetscCall(PetscFree(mRootsOrigNumbering)); 483c4762a1bSJed Brown } 484c4762a1bSJed Brown 485c4762a1bSJed Brown /* Clean storage for star forest. */ 4869566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sf)); 4879566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 488b122ec5aSJacob Faibussowitsch return 0; 489c4762a1bSJed Brown } 490c4762a1bSJed Brown 491c4762a1bSJed Brown /*TEST 492c4762a1bSJed Brown 493c4762a1bSJed Brown test: 494c4762a1bSJed Brown nsize: 4 495c4762a1bSJed Brown filter: grep -v "type" | grep -v "sort" 496c4762a1bSJed Brown args: -test_bcast -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}} 497dfd57a17SPierre Jolivet requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW) 498c4762a1bSJed Brown 499c4762a1bSJed Brown test: 500c4762a1bSJed Brown suffix: 2 501c4762a1bSJed Brown nsize: 4 502c4762a1bSJed Brown filter: grep -v "type" | grep -v "sort" 503c4762a1bSJed Brown args: -test_reduce -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}} 504dfd57a17SPierre Jolivet requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW) 505c4762a1bSJed Brown 506c4762a1bSJed Brown test: 507c4762a1bSJed Brown suffix: 2_basic 508c4762a1bSJed Brown nsize: 4 509c4762a1bSJed Brown args: -test_reduce -sf_type basic 510c4762a1bSJed Brown 511c4762a1bSJed Brown test: 512c4762a1bSJed Brown suffix: 3 513c4762a1bSJed Brown nsize: 4 514c4762a1bSJed Brown filter: grep -v "type" | grep -v "sort" 515c4762a1bSJed Brown args: -test_degree -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}} 516dfd57a17SPierre Jolivet requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW) 517c4762a1bSJed Brown 518c4762a1bSJed Brown test: 519c4762a1bSJed Brown suffix: 3_basic 520c4762a1bSJed Brown nsize: 4 521c4762a1bSJed Brown args: -test_degree -sf_type basic 522c4762a1bSJed Brown 523c4762a1bSJed Brown test: 524c4762a1bSJed Brown suffix: 4 525c4762a1bSJed Brown nsize: 4 526c4762a1bSJed Brown filter: grep -v "type" | grep -v "sort" 527c4762a1bSJed Brown args: -test_gather -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}} 528dfd57a17SPierre Jolivet requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW) 529c4762a1bSJed Brown 530c4762a1bSJed Brown test: 531c4762a1bSJed Brown suffix: 4_basic 532c4762a1bSJed Brown nsize: 4 533c4762a1bSJed Brown args: -test_gather -sf_type basic 534c4762a1bSJed Brown 535c4762a1bSJed Brown test: 536c4762a1bSJed Brown suffix: 4_stride 537c4762a1bSJed Brown nsize: 4 538c4762a1bSJed Brown args: -test_gather -sf_type basic -stride 2 539c4762a1bSJed Brown 540c4762a1bSJed Brown test: 541c4762a1bSJed Brown suffix: 5 542c4762a1bSJed Brown nsize: 4 543c4762a1bSJed Brown filter: grep -v "type" | grep -v "sort" 544c4762a1bSJed Brown args: -test_scatter -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}} 545dfd57a17SPierre Jolivet requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW) 546c4762a1bSJed Brown 547c4762a1bSJed Brown test: 548c4762a1bSJed Brown suffix: 5_basic 549c4762a1bSJed Brown nsize: 4 550c4762a1bSJed Brown args: -test_scatter -sf_type basic 551c4762a1bSJed Brown 552c4762a1bSJed Brown test: 553c4762a1bSJed Brown suffix: 5_stride 554c4762a1bSJed Brown nsize: 4 555c4762a1bSJed Brown args: -test_scatter -sf_type basic -stride 2 556c4762a1bSJed Brown 557c4762a1bSJed Brown test: 558c4762a1bSJed Brown suffix: 6 559c4762a1bSJed Brown nsize: 4 560c4762a1bSJed Brown filter: grep -v "type" | grep -v "sort" 561c20d7725SJed Brown # No -sf_window_flavor dynamic due to bug https://gitlab.com/petsc/petsc/issues/555 562c20d7725SJed Brown args: -test_embed -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create allocate}} 563dfd57a17SPierre Jolivet requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW) 564c4762a1bSJed Brown 565c4762a1bSJed Brown test: 566c4762a1bSJed Brown suffix: 6_basic 567c4762a1bSJed Brown nsize: 4 568c4762a1bSJed Brown args: -test_embed -sf_type basic 569c4762a1bSJed Brown 570c4762a1bSJed Brown test: 571c4762a1bSJed Brown suffix: 7 572c4762a1bSJed Brown nsize: 4 573c4762a1bSJed Brown filter: grep -v "type" | grep -v "sort" 574c4762a1bSJed Brown args: -test_invert -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}} 575dfd57a17SPierre Jolivet requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW) 576c4762a1bSJed Brown 577c4762a1bSJed Brown test: 578c4762a1bSJed Brown suffix: 7_basic 579c4762a1bSJed Brown nsize: 4 580c4762a1bSJed Brown args: -test_invert -sf_type basic 581c4762a1bSJed Brown 582c4762a1bSJed Brown test: 583c4762a1bSJed Brown suffix: basic 584c4762a1bSJed Brown nsize: 4 585c4762a1bSJed Brown args: -test_bcast -sf_type basic 586c4762a1bSJed Brown output_file: output/ex1_1_basic.out 587c4762a1bSJed Brown 588c4762a1bSJed Brown test: 589c4762a1bSJed Brown suffix: bcastop_basic 590c4762a1bSJed Brown nsize: 4 591c4762a1bSJed Brown args: -test_bcastop -sf_type basic 592c4762a1bSJed Brown output_file: output/ex1_bcastop_basic.out 593c4762a1bSJed Brown 594c4762a1bSJed Brown test: 595c4762a1bSJed Brown suffix: 8 596c4762a1bSJed Brown nsize: 3 597c4762a1bSJed Brown filter: grep -v "type" | grep -v "sort" 598c4762a1bSJed Brown args: -test_bcast -test_sf_distribute -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}} 599dfd57a17SPierre Jolivet requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW) 600c4762a1bSJed Brown 601c4762a1bSJed Brown test: 602c4762a1bSJed Brown suffix: 8_basic 603c4762a1bSJed Brown nsize: 3 604c4762a1bSJed Brown args: -test_bcast -test_sf_distribute -sf_type basic 605c4762a1bSJed Brown 606c4762a1bSJed Brown test: 607c4762a1bSJed Brown suffix: 9_char 608c4762a1bSJed Brown nsize: 4 609c4762a1bSJed Brown args: -sf_type basic -test_bcast -test_reduce -test_op max -test_char 610c4762a1bSJed Brown 611c4762a1bSJed Brown # Here we do not test -sf_window_flavor dynamic since it is designed for repeated SFs with few different rootdata pointers 612c4762a1bSJed Brown test: 613c4762a1bSJed Brown suffix: 10 614c4762a1bSJed Brown filter: grep -v "type" | grep -v "sort" 615c4762a1bSJed Brown nsize: 4 616c4762a1bSJed Brown args: -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create allocate}} -test_all -test_bcastop 0 -test_fetchandop 0 617dfd57a17SPierre Jolivet requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW) 618c4762a1bSJed Brown 619c4762a1bSJed Brown # The nightly test suite with MPICH uses ch3:sock, which is broken when winsize == 0 in some of the processes 620c4762a1bSJed Brown test: 621c4762a1bSJed Brown suffix: 10_shared 622c4762a1bSJed Brown output_file: output/ex1_10.out 623c4762a1bSJed Brown filter: grep -v "type" | grep -v "sort" 624c4762a1bSJed Brown nsize: 4 625c4762a1bSJed Brown args: -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor shared -test_all -test_bcastop 0 -test_fetchandop 0 626dfd57a17SPierre Jolivet requires: defined(PETSC_HAVE_MPI_PROCESS_SHARED_MEMORY) !defined(PETSC_HAVE_MPICH_NUMVERSION) defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW) 627c4762a1bSJed Brown 628c4762a1bSJed Brown test: 629c4762a1bSJed Brown suffix: 10_basic 630c4762a1bSJed Brown nsize: 4 631c4762a1bSJed Brown args: -sf_type basic -test_all -test_bcastop 0 -test_fetchandop 0 632c4762a1bSJed Brown 633*0dd791a8SStefano Zampini test: 634*0dd791a8SStefano Zampini suffix: 10_basic_vector 635*0dd791a8SStefano Zampini nsize: 4 636*0dd791a8SStefano Zampini args: -sf_type basic -test_all -test_bcastop 0 -test_fetchandop 0 -test_vector 637*0dd791a8SStefano Zampini 638c4762a1bSJed Brown TEST*/ 639