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 */ 18c4762a1bSJed Brown static PetscErrorCode PetscSFViewCustomLocals_Private(PetscSF sf,const PetscInt locals[],PetscViewer viewer) 19c4762a1bSJed Brown { 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)); 31c4762a1bSJed Brown for (i=0; i<nleaves; i++) { 329566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"[%d] %" PetscInt_FMT " <- (%" PetscInt_FMT ",%" PetscInt_FMT ")\n",rank,locals[i],iremote[i].rank,iremote[i].index)); 33c4762a1bSJed Brown } 349566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 359566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 369566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 37c4762a1bSJed Brown PetscFunctionReturn(0); 38c4762a1bSJed Brown } 39c4762a1bSJed Brown 40c4762a1bSJed Brown int main(int argc,char **argv) 41c4762a1bSJed Brown { 42c4762a1bSJed Brown PetscInt i,nroots,nrootsalloc,nleaves,nleavesalloc,*mine,stride; 43c4762a1bSJed Brown PetscSFNode *remote; 44c4762a1bSJed Brown PetscMPIInt rank,size; 45c4762a1bSJed Brown PetscSF sf; 46c4762a1bSJed Brown 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; 47c4762a1bSJed Brown MPI_Op mop=MPI_OP_NULL; /* initialize to prevent compiler warnings with cxx_quad build */ 48c4762a1bSJed Brown char opstring[256]; 49c4762a1bSJed Brown PetscBool strflg; 50c4762a1bSJed Brown 51*327415f7SBarry Smith PetscFunctionBeginUser; 529566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); 539566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 549566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 55c4762a1bSJed Brown 56d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD,"","PetscSF Test Options","none"); 57c4762a1bSJed Brown test_all = PETSC_FALSE; 589566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-test_all","Test all SF communications","",test_all,&test_all,NULL)); 59c4762a1bSJed Brown test_bcast = test_all; 609566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-test_bcast","Test broadcast","",test_bcast,&test_bcast,NULL)); 61c4762a1bSJed Brown test_bcastop = test_all; 629566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-test_bcastop","Test broadcast and reduce","",test_bcastop,&test_bcastop,NULL)); 63c4762a1bSJed Brown test_reduce = test_all; 649566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-test_reduce","Test reduction","",test_reduce,&test_reduce,NULL)); 65c4762a1bSJed Brown test_char = test_all; 669566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-test_char","Test signed char, unsigned char, and char","",test_char,&test_char,NULL)); 67c4762a1bSJed Brown mop = MPI_SUM; 689566063dSJacob Faibussowitsch PetscCall(PetscStrcpy(opstring,"sum")); 699566063dSJacob Faibussowitsch PetscCall(PetscOptionsString("-test_op","Designate which MPI_Op to use","",opstring,opstring,sizeof(opstring),NULL)); 709566063dSJacob Faibussowitsch PetscCall(PetscStrcmp("sum",opstring,&strflg)); 71c4762a1bSJed Brown if (strflg) { 72c4762a1bSJed Brown mop = MPIU_SUM; 73c4762a1bSJed Brown } 749566063dSJacob Faibussowitsch PetscCall(PetscStrcmp("prod",opstring,&strflg)); 75c4762a1bSJed Brown if (strflg) { 76c4762a1bSJed Brown mop = MPI_PROD; 77c4762a1bSJed Brown } 789566063dSJacob Faibussowitsch PetscCall(PetscStrcmp("max",opstring,&strflg)); 79c4762a1bSJed Brown if (strflg) { 80c4762a1bSJed Brown mop = MPI_MAX; 81c4762a1bSJed Brown } 829566063dSJacob Faibussowitsch PetscCall(PetscStrcmp("min",opstring,&strflg)); 83c4762a1bSJed Brown if (strflg) { 84c4762a1bSJed Brown mop = MPI_MIN; 85c4762a1bSJed Brown } 869566063dSJacob Faibussowitsch PetscCall(PetscStrcmp("land",opstring,&strflg)); 87c4762a1bSJed Brown if (strflg) { 88c4762a1bSJed Brown mop = MPI_LAND; 89c4762a1bSJed Brown } 909566063dSJacob Faibussowitsch PetscCall(PetscStrcmp("band",opstring,&strflg)); 91c4762a1bSJed Brown if (strflg) { 92c4762a1bSJed Brown mop = MPI_BAND; 93c4762a1bSJed Brown } 949566063dSJacob Faibussowitsch PetscCall(PetscStrcmp("lor",opstring,&strflg)); 95c4762a1bSJed Brown if (strflg) { 96c4762a1bSJed Brown mop = MPI_LOR; 97c4762a1bSJed Brown } 989566063dSJacob Faibussowitsch PetscCall(PetscStrcmp("bor",opstring,&strflg)); 99c4762a1bSJed Brown if (strflg) { 100c4762a1bSJed Brown mop = MPI_BOR; 101c4762a1bSJed Brown } 1029566063dSJacob Faibussowitsch PetscCall(PetscStrcmp("lxor",opstring,&strflg)); 103c4762a1bSJed Brown if (strflg) { 104c4762a1bSJed Brown mop = MPI_LXOR; 105c4762a1bSJed Brown } 1069566063dSJacob Faibussowitsch PetscCall(PetscStrcmp("bxor",opstring,&strflg)); 107c4762a1bSJed Brown if (strflg) { 108c4762a1bSJed Brown mop = MPI_BXOR; 109c4762a1bSJed Brown } 110c4762a1bSJed Brown test_degree = test_all; 1119566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-test_degree","Test computation of vertex degree","",test_degree,&test_degree,NULL)); 112c4762a1bSJed Brown test_fetchandop = test_all; 1139566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-test_fetchandop","Test atomic Fetch-And-Op","",test_fetchandop,&test_fetchandop,NULL)); 114c4762a1bSJed Brown test_gather = test_all; 1159566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-test_gather","Test point gather","",test_gather,&test_gather,NULL)); 116c4762a1bSJed Brown test_scatter = test_all; 1179566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-test_scatter","Test point scatter","",test_scatter,&test_scatter,NULL)); 118c4762a1bSJed Brown test_embed = test_all; 1199566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-test_embed","Test point embed","",test_embed,&test_embed,NULL)); 120c4762a1bSJed Brown test_invert = test_all; 1219566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-test_invert","Test point invert","",test_invert,&test_invert,NULL)); 122c4762a1bSJed Brown stride = 1; 1239566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-stride","Stride for leaf and root data","",stride,&stride,NULL)); 124c4762a1bSJed Brown test_sf_distribute = PETSC_FALSE; 1259566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-test_sf_distribute","Create an SF that 'distributes' to each process, like an alltoall","",test_sf_distribute,&test_sf_distribute,NULL)); 1269566063dSJacob Faibussowitsch PetscCall(PetscOptionsString("-test_op","Designate which MPI_Op to use","",opstring,opstring,sizeof(opstring),NULL)); 127d0609cedSBarry Smith PetscOptionsEnd(); 128c4762a1bSJed Brown 129c4762a1bSJed Brown if (test_sf_distribute) { 130c4762a1bSJed Brown nroots = size; 131c4762a1bSJed Brown nrootsalloc = size; 132c4762a1bSJed Brown nleaves = size; 133c4762a1bSJed Brown nleavesalloc = size; 134c4762a1bSJed Brown mine = NULL; 1359566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nleaves,&remote)); 136c4762a1bSJed Brown for (i=0; i<size; i++) { 137c4762a1bSJed Brown remote[i].rank = i; 138c4762a1bSJed Brown remote[i].index = rank; 139c4762a1bSJed Brown } 140c4762a1bSJed Brown } else { 141c4762a1bSJed Brown nroots = 2 + (PetscInt)(rank == 0); 142c4762a1bSJed Brown nrootsalloc = nroots * stride; 143c4762a1bSJed Brown nleaves = 2 + (PetscInt)(rank > 0); 144c4762a1bSJed Brown nleavesalloc = nleaves * stride; 145c4762a1bSJed Brown mine = NULL; 146c4762a1bSJed Brown if (stride > 1) { 147c4762a1bSJed Brown PetscInt i; 148c4762a1bSJed Brown 1499566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nleaves,&mine)); 150c4762a1bSJed Brown for (i = 0; i < nleaves; i++) { 151c4762a1bSJed Brown mine[i] = stride * i; 152c4762a1bSJed Brown } 153c4762a1bSJed Brown } 1549566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nleaves,&remote)); 155c4762a1bSJed Brown /* Left periodic neighbor */ 156c4762a1bSJed Brown remote[0].rank = (rank+size-1)%size; 157c4762a1bSJed Brown remote[0].index = 1 * stride; 158c4762a1bSJed Brown /* Right periodic neighbor */ 159c4762a1bSJed Brown remote[1].rank = (rank+1)%size; 160c4762a1bSJed Brown remote[1].index = 0 * stride; 161c4762a1bSJed Brown if (rank > 0) { /* All processes reference rank 0, index 1 */ 162c4762a1bSJed Brown remote[2].rank = 0; 163c4762a1bSJed Brown remote[2].index = 2 * stride; 164c4762a1bSJed Brown } 165c4762a1bSJed Brown } 166c4762a1bSJed Brown 167c4762a1bSJed Brown /* Create a star forest for communication. In this example, the leaf space is dense, so we pass NULL. */ 1689566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(PETSC_COMM_WORLD,&sf)); 1699566063dSJacob Faibussowitsch PetscCall(PetscSFSetFromOptions(sf)); 1709566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(sf,nrootsalloc,nleaves,mine,PETSC_OWN_POINTER,remote,PETSC_OWN_POINTER)); 1719566063dSJacob Faibussowitsch PetscCall(PetscSFSetUp(sf)); 172c4762a1bSJed Brown 173c4762a1bSJed Brown /* View graph, mostly useful for debugging purposes. */ 1749566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL)); 1759566063dSJacob Faibussowitsch PetscCall(PetscSFView(sf,PETSC_VIEWER_STDOUT_WORLD)); 1769566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD)); 177c4762a1bSJed Brown 178c4762a1bSJed Brown if (test_bcast) { /* broadcast rootdata into leafdata */ 179c4762a1bSJed Brown PetscInt *rootdata,*leafdata; 1802d4ee042Sprj- /* Allocate space for send and receive buffers. This example communicates PetscInt, but other types, including 181c4762a1bSJed Brown * user-defined structures, could also be used. */ 1829566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(nrootsalloc,&rootdata,nleavesalloc,&leafdata)); 183c4762a1bSJed Brown /* Set rootdata buffer to be broadcast */ 184c4762a1bSJed Brown for (i=0; i<nrootsalloc; i++) rootdata[i] = -1; 185c4762a1bSJed Brown for (i=0; i<nroots; i++) rootdata[i*stride] = 100*(rank+1) + i; 186c4762a1bSJed Brown /* Initialize local buffer, these values are never used. */ 187c4762a1bSJed Brown for (i=0; i<nleavesalloc; i++) leafdata[i] = -1; 188c4762a1bSJed Brown /* Broadcast entries from rootdata to leafdata. Computation or other communication can be performed between the begin and end calls. */ 1899566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(sf,MPIU_INT,rootdata,leafdata,MPI_REPLACE)); 1909566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(sf,MPIU_INT,rootdata,leafdata,MPI_REPLACE)); 1919566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Bcast Rootdata\n")); 1929566063dSJacob Faibussowitsch PetscCall(PetscIntView(nrootsalloc,rootdata,PETSC_VIEWER_STDOUT_WORLD)); 1939566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Bcast Leafdata\n")); 1949566063dSJacob Faibussowitsch PetscCall(PetscIntView(nleavesalloc,leafdata,PETSC_VIEWER_STDOUT_WORLD)); 1959566063dSJacob Faibussowitsch PetscCall(PetscFree2(rootdata,leafdata)); 196c4762a1bSJed Brown } 197c4762a1bSJed Brown 198c4762a1bSJed Brown if (test_bcast && test_char) { /* Bcast with char */ 199c4762a1bSJed Brown PetscInt len; 200c4762a1bSJed Brown char buf[256]; 201c4762a1bSJed Brown char *rootdata,*leafdata; 2029566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(nrootsalloc,&rootdata,nleavesalloc,&leafdata)); 203c4762a1bSJed Brown /* Set rootdata buffer to be broadcast */ 204c4762a1bSJed Brown for (i=0; i<nrootsalloc; i++) rootdata[i] = '*'; 205c4762a1bSJed 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 */ 206c4762a1bSJed Brown /* Initialize local buffer, these values are never used. */ 207c4762a1bSJed Brown for (i=0; i<nleavesalloc; i++) leafdata[i] = '?'; 208c4762a1bSJed Brown 2099566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(sf,MPI_CHAR,rootdata,leafdata,MPI_REPLACE)); 2109566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(sf,MPI_CHAR,rootdata,leafdata,MPI_REPLACE)); 211c4762a1bSJed Brown 2129566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Bcast Rootdata in type of char\n")); 2139566063dSJacob Faibussowitsch len = 0; PetscCall(PetscSNPrintf(buf,256,"%4d:",rank)); len += 5; 2149566063dSJacob Faibussowitsch for (i=0; i<nrootsalloc; i++) {PetscCall(PetscSNPrintf(buf+len,256-len,"%5c",rootdata[i])); len += 5;} 2159566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"%s\n",buf)); 2169566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT)); 217c4762a1bSJed Brown 2189566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Bcast Leafdata in type of char\n")); 2199566063dSJacob Faibussowitsch len = 0; PetscCall(PetscSNPrintf(buf,256,"%4d:",rank)); len += 5; 2209566063dSJacob Faibussowitsch for (i=0; i<nleavesalloc; i++) {PetscCall(PetscSNPrintf(buf+len,256-len,"%5c",leafdata[i])); len += 5;} 2219566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"%s\n",buf)); 2229566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT)); 223c4762a1bSJed Brown 2249566063dSJacob Faibussowitsch PetscCall(PetscFree2(rootdata,leafdata)); 225c4762a1bSJed Brown } 226c4762a1bSJed Brown 227c4762a1bSJed Brown if (test_bcastop) { /* Reduce rootdata into leafdata */ 228c4762a1bSJed Brown PetscInt *rootdata,*leafdata; 229c4762a1bSJed Brown /* Allocate space for send and receive buffers. This example communicates PetscInt, but other types, including 230c4762a1bSJed Brown * user-defined structures, could also be used. */ 2319566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(nrootsalloc,&rootdata,nleavesalloc,&leafdata)); 232c4762a1bSJed Brown /* Set rootdata buffer to be broadcast */ 233c4762a1bSJed Brown for (i=0; i<nrootsalloc; i++) rootdata[i] = -1; 234c4762a1bSJed Brown for (i=0; i<nroots; i++) rootdata[i*stride] = 100*(rank+1) + i; 235c4762a1bSJed Brown /* Set leaf values to reduce with */ 236c4762a1bSJed Brown for (i=0; i<nleavesalloc; i++) leafdata[i] = -10*(rank+1) - i; 2379566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Pre-BcastAndOp Leafdata\n")); 2389566063dSJacob Faibussowitsch PetscCall(PetscIntView(nleavesalloc,leafdata,PETSC_VIEWER_STDOUT_WORLD)); 239c4762a1bSJed Brown /* Broadcast entries from rootdata to leafdata. Computation or other communication can be performed between the begin and end calls. */ 2409566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(sf,MPIU_INT,rootdata,leafdata,mop)); 2419566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(sf,MPIU_INT,rootdata,leafdata,mop)); 2429566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## BcastAndOp Rootdata\n")); 2439566063dSJacob Faibussowitsch PetscCall(PetscIntView(nrootsalloc,rootdata,PETSC_VIEWER_STDOUT_WORLD)); 2449566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## BcastAndOp Leafdata\n")); 2459566063dSJacob Faibussowitsch PetscCall(PetscIntView(nleavesalloc,leafdata,PETSC_VIEWER_STDOUT_WORLD)); 2469566063dSJacob Faibussowitsch PetscCall(PetscFree2(rootdata,leafdata)); 247c4762a1bSJed Brown } 248c4762a1bSJed Brown 249c4762a1bSJed Brown if (test_reduce) { /* Reduce leafdata into rootdata */ 250c4762a1bSJed Brown PetscInt *rootdata,*leafdata; 2519566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(nrootsalloc,&rootdata,nleavesalloc,&leafdata)); 252c4762a1bSJed Brown /* Initialize rootdata buffer in which the result of the reduction will appear. */ 253c4762a1bSJed Brown for (i=0; i<nrootsalloc; i++) rootdata[i] = -1; 254c4762a1bSJed Brown for (i=0; i<nroots; i++) rootdata[i*stride] = 100*(rank+1) + i; 255c4762a1bSJed Brown /* Set leaf values to reduce. */ 256c4762a1bSJed Brown for (i=0; i<nleavesalloc; i++) leafdata[i] = -1; 257c4762a1bSJed Brown for (i=0; i<nleaves; i++) leafdata[i*stride] = 1000*(rank+1) + 10*i; 2589566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Pre-Reduce Rootdata\n")); 2599566063dSJacob Faibussowitsch PetscCall(PetscIntView(nrootsalloc,rootdata,PETSC_VIEWER_STDOUT_WORLD)); 260c4762a1bSJed Brown /* Perform reduction. Computation or other communication can be performed between the begin and end calls. 261c4762a1bSJed Brown * This example sums the values, but other MPI_Ops can be used (e.g MPI_MAX, MPI_PROD). */ 2629566063dSJacob Faibussowitsch PetscCall(PetscSFReduceBegin(sf,MPIU_INT,leafdata,rootdata,mop)); 2639566063dSJacob Faibussowitsch PetscCall(PetscSFReduceEnd(sf,MPIU_INT,leafdata,rootdata,mop)); 2649566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Reduce Leafdata\n")); 2659566063dSJacob Faibussowitsch PetscCall(PetscIntView(nleavesalloc,leafdata,PETSC_VIEWER_STDOUT_WORLD)); 2669566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Reduce Rootdata\n")); 2679566063dSJacob Faibussowitsch PetscCall(PetscIntView(nrootsalloc,rootdata,PETSC_VIEWER_STDOUT_WORLD)); 2689566063dSJacob Faibussowitsch PetscCall(PetscFree2(rootdata,leafdata)); 269c4762a1bSJed Brown } 270c4762a1bSJed Brown 271c4762a1bSJed Brown if (test_reduce && test_char) { /* Reduce with signed char */ 272c4762a1bSJed Brown PetscInt len; 273c4762a1bSJed Brown char buf[256]; 274c4762a1bSJed Brown signed char *rootdata,*leafdata; 2759566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(nrootsalloc,&rootdata,nleavesalloc,&leafdata)); 276c4762a1bSJed Brown /* Initialize rootdata buffer in which the result of the reduction will appear. */ 277c4762a1bSJed Brown for (i=0; i<nrootsalloc; i++) rootdata[i] = -1; 278c4762a1bSJed Brown for (i=0; i<nroots; i++) rootdata[i*stride] = 10*(rank+1) + i; 279c4762a1bSJed Brown /* Set leaf values to reduce. */ 280c4762a1bSJed Brown for (i=0; i<nleavesalloc; i++) leafdata[i] = -1; 281c4762a1bSJed Brown for (i=0; i<nleaves; i++) leafdata[i*stride] = 50*(rank+1) + 10*i; 2829566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Pre-Reduce Rootdata in type of signed char\n")); 283c4762a1bSJed Brown 2849566063dSJacob Faibussowitsch len = 0; PetscCall(PetscSNPrintf(buf,256,"%4d:",rank)); len += 5; 2859566063dSJacob Faibussowitsch for (i=0; i<nrootsalloc; i++) {PetscCall(PetscSNPrintf(buf+len,256-len,"%5d",rootdata[i])); len += 5;} 2869566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"%s\n",buf)); 2879566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT)); 288c4762a1bSJed Brown 289c4762a1bSJed Brown /* Using MPI_CHAR should trigger an error since MPI standard does not support reduction on MPI_CHAR. 290c4762a1bSJed Brown Testing with -test_op max, one can see the sign does take effect in MPI_MAX. 291c4762a1bSJed Brown */ 2929566063dSJacob Faibussowitsch PetscCall(PetscSFReduceBegin(sf,MPI_SIGNED_CHAR,leafdata,rootdata,mop)); 2939566063dSJacob Faibussowitsch PetscCall(PetscSFReduceEnd(sf,MPI_SIGNED_CHAR,leafdata,rootdata,mop)); 294c4762a1bSJed Brown 2959566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Reduce Leafdata in type of signed char\n")); 2969566063dSJacob Faibussowitsch len = 0; PetscCall(PetscSNPrintf(buf,256,"%4d:",rank)); len += 5; 2979566063dSJacob Faibussowitsch for (i=0; i<nleavesalloc; i++) {PetscCall(PetscSNPrintf(buf+len,256-len,"%5d",leafdata[i])); len += 5;} 2989566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"%s\n",buf)); 2999566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT)); 300c4762a1bSJed Brown 3019566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Reduce Rootdata in type of signed char\n")); 3029566063dSJacob Faibussowitsch len = 0; PetscCall(PetscSNPrintf(buf,256,"%4d:",rank)); len += 5; 3039566063dSJacob Faibussowitsch for (i=0; i<nrootsalloc; i++) {PetscCall(PetscSNPrintf(buf+len,256-len,"%5d",rootdata[i])); len += 5;} 3049566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"%s\n",buf)); 3059566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT)); 306c4762a1bSJed Brown 3079566063dSJacob Faibussowitsch PetscCall(PetscFree2(rootdata,leafdata)); 308c4762a1bSJed Brown } 309c4762a1bSJed Brown 310c4762a1bSJed Brown if (test_reduce && test_char) { /* Reduce with unsigned char */ 311c4762a1bSJed Brown PetscInt len; 312c4762a1bSJed Brown char buf[256]; 313c4762a1bSJed Brown unsigned char *rootdata,*leafdata; 3149566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(nrootsalloc,&rootdata,nleavesalloc,&leafdata)); 315c4762a1bSJed Brown /* Initialize rootdata buffer in which the result of the reduction will appear. */ 316b9d03b0cSStefano Zampini for (i=0; i<nrootsalloc; i++) rootdata[i] = 0; 317c4762a1bSJed Brown for (i=0; i<nroots; i++) rootdata[i*stride] = 10*(rank+1) + i; 318c4762a1bSJed Brown /* Set leaf values to reduce. */ 319b9d03b0cSStefano Zampini for (i=0; i<nleavesalloc; i++) leafdata[i] = 0; 320c4762a1bSJed Brown for (i=0; i<nleaves; i++) leafdata[i*stride] = 50*(rank+1) + 10*i; 3219566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Pre-Reduce Rootdata in type of unsigned char\n")); 322c4762a1bSJed Brown 3239566063dSJacob Faibussowitsch len = 0; PetscCall(PetscSNPrintf(buf,256,"%4d:",rank)); len += 5; 3249566063dSJacob Faibussowitsch for (i=0; i<nrootsalloc; i++) {PetscCall(PetscSNPrintf(buf+len,256-len,"%5u",rootdata[i])); len += 5;} 3259566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"%s\n",buf)); 3269566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT)); 327c4762a1bSJed Brown 328c4762a1bSJed Brown /* Using MPI_CHAR should trigger an error since MPI standard does not support reduction on MPI_CHAR. 329c4762a1bSJed Brown Testing with -test_op max, one can see the sign does take effect in MPI_MAX. 330c4762a1bSJed Brown */ 3319566063dSJacob Faibussowitsch PetscCall(PetscSFReduceBegin(sf,MPI_UNSIGNED_CHAR,leafdata,rootdata,mop)); 3329566063dSJacob Faibussowitsch PetscCall(PetscSFReduceEnd(sf,MPI_UNSIGNED_CHAR,leafdata,rootdata,mop)); 333c4762a1bSJed Brown 3349566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Reduce Leafdata in type of unsigned char\n")); 3359566063dSJacob Faibussowitsch len = 0; PetscCall(PetscSNPrintf(buf,256,"%4d:",rank)); len += 5; 3369566063dSJacob Faibussowitsch for (i=0; i<nleavesalloc; i++) {PetscCall(PetscSNPrintf(buf+len,256-len,"%5u",leafdata[i])); len += 5;} 3379566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"%s\n",buf)); 3389566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT)); 339c4762a1bSJed Brown 3409566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Reduce Rootdata in type of unsigned char\n")); 3419566063dSJacob Faibussowitsch len = 0; PetscCall(PetscSNPrintf(buf,256,"%4d:",rank)); len += 5; 3429566063dSJacob Faibussowitsch for (i=0; i<nrootsalloc; i++) {PetscCall(PetscSNPrintf(buf+len,256-len,"%5u",rootdata[i])); len += 5;} 3439566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"%s\n",buf)); 3449566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT)); 345c4762a1bSJed Brown 3469566063dSJacob Faibussowitsch PetscCall(PetscFree2(rootdata,leafdata)); 347c4762a1bSJed Brown } 348c4762a1bSJed Brown 349c4762a1bSJed Brown if (test_degree) { 350c4762a1bSJed Brown const PetscInt *degree; 3519566063dSJacob Faibussowitsch PetscCall(PetscSFComputeDegreeBegin(sf,°ree)); 3529566063dSJacob Faibussowitsch PetscCall(PetscSFComputeDegreeEnd(sf,°ree)); 3539566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Root degrees\n")); 3549566063dSJacob Faibussowitsch PetscCall(PetscIntView(nrootsalloc,degree,PETSC_VIEWER_STDOUT_WORLD)); 355c4762a1bSJed Brown } 356c4762a1bSJed Brown 357c4762a1bSJed Brown if (test_fetchandop) { 358c4762a1bSJed Brown /* Cannot use text compare here because token ordering is not deterministic */ 359c4762a1bSJed Brown PetscInt *leafdata,*leafupdate,*rootdata; 3609566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(nleavesalloc,&leafdata,nleavesalloc,&leafupdate,nrootsalloc,&rootdata)); 361c4762a1bSJed Brown for (i=0; i<nleavesalloc; i++) leafdata[i] = -1; 362c4762a1bSJed Brown for (i=0; i<nleaves; i++) leafdata[i*stride] = 1; 363c4762a1bSJed Brown for (i=0; i<nrootsalloc; i++) rootdata[i] = -1; 364c4762a1bSJed Brown for (i=0; i<nroots; i++) rootdata[i*stride] = 0; 3659566063dSJacob Faibussowitsch PetscCall(PetscSFFetchAndOpBegin(sf,MPIU_INT,rootdata,leafdata,leafupdate,mop)); 3669566063dSJacob Faibussowitsch PetscCall(PetscSFFetchAndOpEnd(sf,MPIU_INT,rootdata,leafdata,leafupdate,mop)); 3679566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Rootdata (sum of 1 from each leaf)\n")); 3689566063dSJacob Faibussowitsch PetscCall(PetscIntView(nrootsalloc,rootdata,PETSC_VIEWER_STDOUT_WORLD)); 3699566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Leafupdate (value at roots prior to my atomic update)\n")); 3709566063dSJacob Faibussowitsch PetscCall(PetscIntView(nleavesalloc,leafupdate,PETSC_VIEWER_STDOUT_WORLD)); 3719566063dSJacob Faibussowitsch PetscCall(PetscFree3(leafdata,leafupdate,rootdata)); 372c4762a1bSJed Brown } 373c4762a1bSJed Brown 374c4762a1bSJed Brown if (test_gather) { 375c4762a1bSJed Brown const PetscInt *degree; 376c4762a1bSJed Brown PetscInt inedges,*indata,*outdata; 3779566063dSJacob Faibussowitsch PetscCall(PetscSFComputeDegreeBegin(sf,°ree)); 3789566063dSJacob Faibussowitsch PetscCall(PetscSFComputeDegreeEnd(sf,°ree)); 379c4762a1bSJed Brown for (i=0,inedges=0; i<nrootsalloc; i++) inedges += degree[i]; 3809566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(inedges,&indata,nleavesalloc,&outdata)); 381c4762a1bSJed Brown for (i=0; i<nleavesalloc; i++) outdata[i] = -1; 382c4762a1bSJed Brown for (i=0; i<nleaves; i++) outdata[i*stride] = 1000*(rank+1) + i; 3839566063dSJacob Faibussowitsch PetscCall(PetscSFGatherBegin(sf,MPIU_INT,outdata,indata)); 3849566063dSJacob Faibussowitsch PetscCall(PetscSFGatherEnd(sf,MPIU_INT,outdata,indata)); 3859566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Gathered data at multi-roots from leaves\n")); 3869566063dSJacob Faibussowitsch PetscCall(PetscIntView(inedges,indata,PETSC_VIEWER_STDOUT_WORLD)); 3879566063dSJacob Faibussowitsch PetscCall(PetscFree2(indata,outdata)); 388c4762a1bSJed Brown } 389c4762a1bSJed Brown 390c4762a1bSJed Brown if (test_scatter) { 391c4762a1bSJed Brown const PetscInt *degree; 392c4762a1bSJed Brown PetscInt j,count,inedges,*indata,*outdata; 3939566063dSJacob Faibussowitsch PetscCall(PetscSFComputeDegreeBegin(sf,°ree)); 3949566063dSJacob Faibussowitsch PetscCall(PetscSFComputeDegreeEnd(sf,°ree)); 395c4762a1bSJed Brown for (i=0,inedges=0; i<nrootsalloc; i++) inedges += degree[i]; 3969566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(inedges,&indata,nleavesalloc,&outdata)); 397c4762a1bSJed Brown for (i=0; i<nleavesalloc; i++) outdata[i] = -1; 398c4762a1bSJed Brown for (i=0,count=0; i<nrootsalloc; i++) { 399c4762a1bSJed Brown for (j=0; j<degree[i]; j++) indata[count++] = 1000*(rank+1) + 100*i + j; 400c4762a1bSJed Brown } 4019566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Data at multi-roots, to scatter to leaves\n")); 4029566063dSJacob Faibussowitsch PetscCall(PetscIntView(inedges,indata,PETSC_VIEWER_STDOUT_WORLD)); 403c4762a1bSJed Brown 4049566063dSJacob Faibussowitsch PetscCall(PetscSFScatterBegin(sf,MPIU_INT,indata,outdata)); 4059566063dSJacob Faibussowitsch PetscCall(PetscSFScatterEnd(sf,MPIU_INT,indata,outdata)); 4069566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Scattered data at leaves\n")); 4079566063dSJacob Faibussowitsch PetscCall(PetscIntView(nleavesalloc,outdata,PETSC_VIEWER_STDOUT_WORLD)); 4089566063dSJacob Faibussowitsch PetscCall(PetscFree2(indata,outdata)); 409c4762a1bSJed Brown } 410c4762a1bSJed Brown 411c4762a1bSJed Brown if (test_embed) { 412dd400576SPatrick Sanan const PetscInt nroots = 1 + (PetscInt) (rank == 0); 413c4762a1bSJed Brown PetscInt selected[2]; 414c4762a1bSJed Brown PetscSF esf; 415c4762a1bSJed Brown 416c4762a1bSJed Brown selected[0] = stride; 417c4762a1bSJed Brown selected[1] = 2*stride; 4189566063dSJacob Faibussowitsch PetscCall(PetscSFCreateEmbeddedRootSF(sf,nroots,selected,&esf)); 4199566063dSJacob Faibussowitsch PetscCall(PetscSFSetUp(esf)); 4209566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Embedded PetscSF\n")); 4219566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL)); 4229566063dSJacob Faibussowitsch PetscCall(PetscSFView(esf,PETSC_VIEWER_STDOUT_WORLD)); 4239566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD)); 4249566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&esf)); 425c4762a1bSJed Brown } 426c4762a1bSJed Brown 427c4762a1bSJed Brown if (test_invert) { 428c4762a1bSJed Brown const PetscInt *degree; 429c4762a1bSJed Brown PetscInt *mRootsOrigNumbering; 430c4762a1bSJed Brown PetscInt inedges; 431c4762a1bSJed Brown PetscSF msf,imsf; 432c4762a1bSJed Brown 4339566063dSJacob Faibussowitsch PetscCall(PetscSFGetMultiSF(sf,&msf)); 4349566063dSJacob Faibussowitsch PetscCall(PetscSFCreateInverseSF(msf,&imsf)); 4359566063dSJacob Faibussowitsch PetscCall(PetscSFSetUp(msf)); 4369566063dSJacob Faibussowitsch PetscCall(PetscSFSetUp(imsf)); 4379566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Multi-SF\n")); 4389566063dSJacob Faibussowitsch PetscCall(PetscSFView(msf,PETSC_VIEWER_STDOUT_WORLD)); 4399566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Multi-SF roots indices in original SF roots numbering\n")); 4409566063dSJacob Faibussowitsch PetscCall(PetscSFComputeDegreeBegin(sf,°ree)); 4419566063dSJacob Faibussowitsch PetscCall(PetscSFComputeDegreeEnd(sf,°ree)); 4429566063dSJacob Faibussowitsch PetscCall(PetscSFComputeMultiRootOriginalNumbering(sf,degree,&inedges,&mRootsOrigNumbering)); 4439566063dSJacob Faibussowitsch PetscCall(PetscIntView(inedges,mRootsOrigNumbering,PETSC_VIEWER_STDOUT_WORLD)); 4449566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Inverse of Multi-SF\n")); 4459566063dSJacob Faibussowitsch PetscCall(PetscSFView(imsf,PETSC_VIEWER_STDOUT_WORLD)); 4469566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Inverse of Multi-SF, original numbering\n")); 4479566063dSJacob Faibussowitsch PetscCall(PetscSFViewCustomLocals_Private(imsf,mRootsOrigNumbering,PETSC_VIEWER_STDOUT_WORLD)); 4489566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&imsf)); 4499566063dSJacob Faibussowitsch PetscCall(PetscFree(mRootsOrigNumbering)); 450c4762a1bSJed Brown } 451c4762a1bSJed Brown 452c4762a1bSJed Brown /* Clean storage for star forest. */ 4539566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sf)); 4549566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 455b122ec5aSJacob Faibussowitsch return 0; 456c4762a1bSJed Brown } 457c4762a1bSJed Brown 458c4762a1bSJed Brown /*TEST 459c4762a1bSJed Brown 460c4762a1bSJed Brown test: 461c4762a1bSJed Brown nsize: 4 462c4762a1bSJed Brown filter: grep -v "type" | grep -v "sort" 463c4762a1bSJed Brown args: -test_bcast -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}} 464dfd57a17SPierre Jolivet requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW) 465c4762a1bSJed Brown 466c4762a1bSJed Brown test: 467c4762a1bSJed Brown suffix: 2 468c4762a1bSJed Brown nsize: 4 469c4762a1bSJed Brown filter: grep -v "type" | grep -v "sort" 470c4762a1bSJed Brown args: -test_reduce -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}} 471dfd57a17SPierre Jolivet requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW) 472c4762a1bSJed Brown 473c4762a1bSJed Brown test: 474c4762a1bSJed Brown suffix: 2_basic 475c4762a1bSJed Brown nsize: 4 476c4762a1bSJed Brown args: -test_reduce -sf_type basic 477c4762a1bSJed Brown 478c4762a1bSJed Brown test: 479c4762a1bSJed Brown suffix: 3 480c4762a1bSJed Brown nsize: 4 481c4762a1bSJed Brown filter: grep -v "type" | grep -v "sort" 482c4762a1bSJed Brown args: -test_degree -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}} 483dfd57a17SPierre Jolivet requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW) 484c4762a1bSJed Brown 485c4762a1bSJed Brown test: 486c4762a1bSJed Brown suffix: 3_basic 487c4762a1bSJed Brown nsize: 4 488c4762a1bSJed Brown args: -test_degree -sf_type basic 489c4762a1bSJed Brown 490c4762a1bSJed Brown test: 491c4762a1bSJed Brown suffix: 4 492c4762a1bSJed Brown nsize: 4 493c4762a1bSJed Brown filter: grep -v "type" | grep -v "sort" 494c4762a1bSJed Brown args: -test_gather -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}} 495dfd57a17SPierre Jolivet requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW) 496c4762a1bSJed Brown 497c4762a1bSJed Brown test: 498c4762a1bSJed Brown suffix: 4_basic 499c4762a1bSJed Brown nsize: 4 500c4762a1bSJed Brown args: -test_gather -sf_type basic 501c4762a1bSJed Brown 502c4762a1bSJed Brown test: 503c4762a1bSJed Brown suffix: 4_stride 504c4762a1bSJed Brown nsize: 4 505c4762a1bSJed Brown args: -test_gather -sf_type basic -stride 2 506c4762a1bSJed Brown 507c4762a1bSJed Brown test: 508c4762a1bSJed Brown suffix: 5 509c4762a1bSJed Brown nsize: 4 510c4762a1bSJed Brown filter: grep -v "type" | grep -v "sort" 511c4762a1bSJed Brown args: -test_scatter -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}} 512dfd57a17SPierre Jolivet requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW) 513c4762a1bSJed Brown 514c4762a1bSJed Brown test: 515c4762a1bSJed Brown suffix: 5_basic 516c4762a1bSJed Brown nsize: 4 517c4762a1bSJed Brown args: -test_scatter -sf_type basic 518c4762a1bSJed Brown 519c4762a1bSJed Brown test: 520c4762a1bSJed Brown suffix: 5_stride 521c4762a1bSJed Brown nsize: 4 522c4762a1bSJed Brown args: -test_scatter -sf_type basic -stride 2 523c4762a1bSJed Brown 524c4762a1bSJed Brown test: 525c4762a1bSJed Brown suffix: 6 526c4762a1bSJed Brown nsize: 4 527c4762a1bSJed Brown filter: grep -v "type" | grep -v "sort" 528c20d7725SJed Brown # No -sf_window_flavor dynamic due to bug https://gitlab.com/petsc/petsc/issues/555 529c20d7725SJed Brown args: -test_embed -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create allocate}} 530dfd57a17SPierre Jolivet requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW) 531c4762a1bSJed Brown 532c4762a1bSJed Brown test: 533c4762a1bSJed Brown suffix: 6_basic 534c4762a1bSJed Brown nsize: 4 535c4762a1bSJed Brown args: -test_embed -sf_type basic 536c4762a1bSJed Brown 537c4762a1bSJed Brown test: 538c4762a1bSJed Brown suffix: 7 539c4762a1bSJed Brown nsize: 4 540c4762a1bSJed Brown filter: grep -v "type" | grep -v "sort" 541c4762a1bSJed Brown args: -test_invert -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}} 542dfd57a17SPierre Jolivet requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW) 543c4762a1bSJed Brown 544c4762a1bSJed Brown test: 545c4762a1bSJed Brown suffix: 7_basic 546c4762a1bSJed Brown nsize: 4 547c4762a1bSJed Brown args: -test_invert -sf_type basic 548c4762a1bSJed Brown 549c4762a1bSJed Brown test: 550c4762a1bSJed Brown suffix: basic 551c4762a1bSJed Brown nsize: 4 552c4762a1bSJed Brown args: -test_bcast -sf_type basic 553c4762a1bSJed Brown output_file: output/ex1_1_basic.out 554c4762a1bSJed Brown 555c4762a1bSJed Brown test: 556c4762a1bSJed Brown suffix: bcastop_basic 557c4762a1bSJed Brown nsize: 4 558c4762a1bSJed Brown args: -test_bcastop -sf_type basic 559c4762a1bSJed Brown output_file: output/ex1_bcastop_basic.out 560c4762a1bSJed Brown 561c4762a1bSJed Brown test: 562c4762a1bSJed Brown suffix: 8 563c4762a1bSJed Brown nsize: 3 564c4762a1bSJed Brown filter: grep -v "type" | grep -v "sort" 565c4762a1bSJed Brown args: -test_bcast -test_sf_distribute -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}} 566dfd57a17SPierre Jolivet requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW) 567c4762a1bSJed Brown 568c4762a1bSJed Brown test: 569c4762a1bSJed Brown suffix: 8_basic 570c4762a1bSJed Brown nsize: 3 571c4762a1bSJed Brown args: -test_bcast -test_sf_distribute -sf_type basic 572c4762a1bSJed Brown 573c4762a1bSJed Brown test: 574c4762a1bSJed Brown suffix: 9_char 575c4762a1bSJed Brown nsize: 4 576c4762a1bSJed Brown args: -sf_type basic -test_bcast -test_reduce -test_op max -test_char 577c4762a1bSJed Brown 578c4762a1bSJed Brown # Here we do not test -sf_window_flavor dynamic since it is designed for repeated SFs with few different rootdata pointers 579c4762a1bSJed Brown test: 580c4762a1bSJed Brown suffix: 10 581c4762a1bSJed Brown filter: grep -v "type" | grep -v "sort" 582c4762a1bSJed Brown nsize: 4 583c4762a1bSJed Brown args: -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create allocate}} -test_all -test_bcastop 0 -test_fetchandop 0 584dfd57a17SPierre Jolivet requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW) 585c4762a1bSJed Brown 586c4762a1bSJed Brown # The nightly test suite with MPICH uses ch3:sock, which is broken when winsize == 0 in some of the processes 587c4762a1bSJed Brown test: 588c4762a1bSJed Brown suffix: 10_shared 589c4762a1bSJed Brown output_file: output/ex1_10.out 590c4762a1bSJed Brown filter: grep -v "type" | grep -v "sort" 591c4762a1bSJed Brown nsize: 4 592c4762a1bSJed Brown args: -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor shared -test_all -test_bcastop 0 -test_fetchandop 0 593dfd57a17SPierre 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) 594c4762a1bSJed Brown 595c4762a1bSJed Brown test: 596c4762a1bSJed Brown suffix: 10_basic 597c4762a1bSJed Brown nsize: 4 598c4762a1bSJed Brown args: -sf_type basic -test_all -test_bcastop 0 -test_fetchandop 0 599c4762a1bSJed Brown 600c4762a1bSJed Brown TEST*/ 601