1c4762a1bSJed Brown static const char help[] = "Test star forest communication (PetscSF)\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown /*T 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. 8c4762a1bSJed Brown T*/ 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 PetscErrorCode ierr; 24c4762a1bSJed Brown 25c4762a1bSJed Brown PetscFunctionBeginUser; 26ffc4695bSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRMPI(ierr); 27c4762a1bSJed Brown ierr = PetscSFGetGraph(sf,&nroots,&nleaves,NULL,&iremote);CHKERRQ(ierr); 28c4762a1bSJed Brown ierr = PetscSFGetRootRanks(sf,&nranks,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 29c4762a1bSJed Brown ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 30c4762a1bSJed Brown ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 31c4762a1bSJed Brown ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Number of roots=%D, leaves=%D, remote ranks=%D\n",rank,nroots,nleaves,nranks);CHKERRQ(ierr); 32c4762a1bSJed Brown for (i=0; i<nleaves; i++) { 33c4762a1bSJed Brown ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] %D <- (%D,%D)\n",rank,locals[i],iremote[i].rank,iremote[i].index);CHKERRQ(ierr); 34c4762a1bSJed Brown } 35c4762a1bSJed Brown ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 36c4762a1bSJed Brown ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 37c4762a1bSJed Brown ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 38c4762a1bSJed Brown PetscFunctionReturn(0); 39c4762a1bSJed Brown } 40c4762a1bSJed Brown 41c4762a1bSJed Brown int main(int argc,char **argv) 42c4762a1bSJed Brown { 43c4762a1bSJed Brown PetscErrorCode ierr; 44c4762a1bSJed Brown PetscInt i,nroots,nrootsalloc,nleaves,nleavesalloc,*mine,stride; 45c4762a1bSJed Brown PetscSFNode *remote; 46c4762a1bSJed Brown PetscMPIInt rank,size; 47c4762a1bSJed Brown PetscSF sf; 48c4762a1bSJed 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; 49c4762a1bSJed Brown MPI_Op mop=MPI_OP_NULL; /* initialize to prevent compiler warnings with cxx_quad build */ 50c4762a1bSJed Brown char opstring[256]; 51c4762a1bSJed Brown PetscBool strflg; 52c4762a1bSJed Brown 53c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 54ffc4695bSBarry Smith ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr); 55ffc4695bSBarry Smith ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); 56c4762a1bSJed Brown 57c4762a1bSJed Brown ierr = PetscOptionsBegin(PETSC_COMM_WORLD,"","PetscSF Test Options","none");CHKERRQ(ierr); 58c4762a1bSJed Brown test_all = PETSC_FALSE; 59c4762a1bSJed Brown ierr = PetscOptionsBool("-test_all","Test all SF communications","",test_all,&test_all,NULL);CHKERRQ(ierr); 60c4762a1bSJed Brown test_bcast = test_all; 61c4762a1bSJed Brown ierr = PetscOptionsBool("-test_bcast","Test broadcast","",test_bcast,&test_bcast,NULL);CHKERRQ(ierr); 62c4762a1bSJed Brown test_bcastop = test_all; 63c4762a1bSJed Brown ierr = PetscOptionsBool("-test_bcastop","Test broadcast and reduce","",test_bcastop,&test_bcastop,NULL);CHKERRQ(ierr); 64c4762a1bSJed Brown test_reduce = test_all; 65c4762a1bSJed Brown ierr = PetscOptionsBool("-test_reduce","Test reduction","",test_reduce,&test_reduce,NULL);CHKERRQ(ierr); 66c4762a1bSJed Brown test_char = test_all; 67c4762a1bSJed Brown ierr = PetscOptionsBool("-test_char","Test signed char, unsigned char, and char","",test_char,&test_char,NULL);CHKERRQ(ierr); 68c4762a1bSJed Brown mop = MPI_SUM; 69c4762a1bSJed Brown ierr = PetscStrcpy(opstring,"sum");CHKERRQ(ierr); 70589a23caSBarry Smith ierr = PetscOptionsString("-test_op","Designate which MPI_Op to use","",opstring,opstring,sizeof(opstring),NULL);CHKERRQ(ierr); 71c4762a1bSJed Brown ierr = PetscStrcmp("sum",opstring,&strflg);CHKERRQ(ierr); 72c4762a1bSJed Brown if (strflg) { 73c4762a1bSJed Brown mop = MPIU_SUM; 74c4762a1bSJed Brown } 75c4762a1bSJed Brown ierr = PetscStrcmp("prod",opstring,&strflg);CHKERRQ(ierr); 76c4762a1bSJed Brown if (strflg) { 77c4762a1bSJed Brown mop = MPI_PROD; 78c4762a1bSJed Brown } 79c4762a1bSJed Brown ierr = PetscStrcmp("max",opstring,&strflg);CHKERRQ(ierr); 80c4762a1bSJed Brown if (strflg) { 81c4762a1bSJed Brown mop = MPI_MAX; 82c4762a1bSJed Brown } 83c4762a1bSJed Brown ierr = PetscStrcmp("min",opstring,&strflg);CHKERRQ(ierr); 84c4762a1bSJed Brown if (strflg) { 85c4762a1bSJed Brown mop = MPI_MIN; 86c4762a1bSJed Brown } 87c4762a1bSJed Brown ierr = PetscStrcmp("land",opstring,&strflg);CHKERRQ(ierr); 88c4762a1bSJed Brown if (strflg) { 89c4762a1bSJed Brown mop = MPI_LAND; 90c4762a1bSJed Brown } 91c4762a1bSJed Brown ierr = PetscStrcmp("band",opstring,&strflg);CHKERRQ(ierr); 92c4762a1bSJed Brown if (strflg) { 93c4762a1bSJed Brown mop = MPI_BAND; 94c4762a1bSJed Brown } 95c4762a1bSJed Brown ierr = PetscStrcmp("lor",opstring,&strflg);CHKERRQ(ierr); 96c4762a1bSJed Brown if (strflg) { 97c4762a1bSJed Brown mop = MPI_LOR; 98c4762a1bSJed Brown } 99c4762a1bSJed Brown ierr = PetscStrcmp("bor",opstring,&strflg);CHKERRQ(ierr); 100c4762a1bSJed Brown if (strflg) { 101c4762a1bSJed Brown mop = MPI_BOR; 102c4762a1bSJed Brown } 103c4762a1bSJed Brown ierr = PetscStrcmp("lxor",opstring,&strflg);CHKERRQ(ierr); 104c4762a1bSJed Brown if (strflg) { 105c4762a1bSJed Brown mop = MPI_LXOR; 106c4762a1bSJed Brown } 107c4762a1bSJed Brown ierr = PetscStrcmp("bxor",opstring,&strflg);CHKERRQ(ierr); 108c4762a1bSJed Brown if (strflg) { 109c4762a1bSJed Brown mop = MPI_BXOR; 110c4762a1bSJed Brown } 111c4762a1bSJed Brown test_degree = test_all; 112c4762a1bSJed Brown ierr = PetscOptionsBool("-test_degree","Test computation of vertex degree","",test_degree,&test_degree,NULL);CHKERRQ(ierr); 113c4762a1bSJed Brown test_fetchandop = test_all; 114c4762a1bSJed Brown ierr = PetscOptionsBool("-test_fetchandop","Test atomic Fetch-And-Op","",test_fetchandop,&test_fetchandop,NULL);CHKERRQ(ierr); 115c4762a1bSJed Brown test_gather = test_all; 116c4762a1bSJed Brown ierr = PetscOptionsBool("-test_gather","Test point gather","",test_gather,&test_gather,NULL);CHKERRQ(ierr); 117c4762a1bSJed Brown test_scatter = test_all; 118c4762a1bSJed Brown ierr = PetscOptionsBool("-test_scatter","Test point scatter","",test_scatter,&test_scatter,NULL);CHKERRQ(ierr); 119c4762a1bSJed Brown test_embed = test_all; 120c4762a1bSJed Brown ierr = PetscOptionsBool("-test_embed","Test point embed","",test_embed,&test_embed,NULL);CHKERRQ(ierr); 121c4762a1bSJed Brown test_invert = test_all; 122c4762a1bSJed Brown ierr = PetscOptionsBool("-test_invert","Test point invert","",test_invert,&test_invert,NULL);CHKERRQ(ierr); 123c4762a1bSJed Brown stride = 1; 124c4762a1bSJed Brown ierr = PetscOptionsInt("-stride","Stride for leaf and root data","",stride,&stride,NULL);CHKERRQ(ierr); 125c4762a1bSJed Brown test_sf_distribute = PETSC_FALSE; 126c4762a1bSJed Brown ierr = PetscOptionsBool("-test_sf_distribute","Create an SF that 'distributes' to each process, like an alltoall","",test_sf_distribute,&test_sf_distribute,NULL);CHKERRQ(ierr); 127589a23caSBarry Smith ierr = PetscOptionsString("-test_op","Designate which MPI_Op to use","",opstring,opstring,sizeof(opstring),NULL);CHKERRQ(ierr); 128c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 129c4762a1bSJed Brown 130c4762a1bSJed Brown if (test_sf_distribute) { 131c4762a1bSJed Brown nroots = size; 132c4762a1bSJed Brown nrootsalloc = size; 133c4762a1bSJed Brown nleaves = size; 134c4762a1bSJed Brown nleavesalloc = size; 135c4762a1bSJed Brown mine = NULL; 136c4762a1bSJed Brown ierr = PetscMalloc1(nleaves,&remote);CHKERRQ(ierr); 137c4762a1bSJed Brown for (i=0; i<size; i++) { 138c4762a1bSJed Brown remote[i].rank = i; 139c4762a1bSJed Brown remote[i].index = rank; 140c4762a1bSJed Brown } 141c4762a1bSJed Brown } else { 142c4762a1bSJed Brown nroots = 2 + (PetscInt)(rank == 0); 143c4762a1bSJed Brown nrootsalloc = nroots * stride; 144c4762a1bSJed Brown nleaves = 2 + (PetscInt)(rank > 0); 145c4762a1bSJed Brown nleavesalloc = nleaves * stride; 146c4762a1bSJed Brown mine = NULL; 147c4762a1bSJed Brown if (stride > 1) { 148c4762a1bSJed Brown PetscInt i; 149c4762a1bSJed Brown 150c4762a1bSJed Brown ierr = PetscMalloc1(nleaves,&mine);CHKERRQ(ierr); 151c4762a1bSJed Brown for (i = 0; i < nleaves; i++) { 152c4762a1bSJed Brown mine[i] = stride * i; 153c4762a1bSJed Brown } 154c4762a1bSJed Brown } 155c4762a1bSJed Brown ierr = PetscMalloc1(nleaves,&remote);CHKERRQ(ierr); 156c4762a1bSJed Brown /* Left periodic neighbor */ 157c4762a1bSJed Brown remote[0].rank = (rank+size-1)%size; 158c4762a1bSJed Brown remote[0].index = 1 * stride; 159c4762a1bSJed Brown /* Right periodic neighbor */ 160c4762a1bSJed Brown remote[1].rank = (rank+1)%size; 161c4762a1bSJed Brown remote[1].index = 0 * stride; 162c4762a1bSJed Brown if (rank > 0) { /* All processes reference rank 0, index 1 */ 163c4762a1bSJed Brown remote[2].rank = 0; 164c4762a1bSJed Brown remote[2].index = 2 * stride; 165c4762a1bSJed Brown } 166c4762a1bSJed Brown } 167c4762a1bSJed Brown 168c4762a1bSJed Brown /* Create a star forest for communication. In this example, the leaf space is dense, so we pass NULL. */ 169c4762a1bSJed Brown ierr = PetscSFCreate(PETSC_COMM_WORLD,&sf);CHKERRQ(ierr); 170c4762a1bSJed Brown ierr = PetscSFSetFromOptions(sf);CHKERRQ(ierr); 171c4762a1bSJed Brown ierr = PetscSFSetGraph(sf,nrootsalloc,nleaves,mine,PETSC_OWN_POINTER,remote,PETSC_OWN_POINTER);CHKERRQ(ierr); 172c4762a1bSJed Brown ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 173c4762a1bSJed Brown 174c4762a1bSJed Brown /* View graph, mostly useful for debugging purposes. */ 175c4762a1bSJed Brown ierr = PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);CHKERRQ(ierr); 176c4762a1bSJed Brown ierr = PetscSFView(sf,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 177c4762a1bSJed Brown ierr = PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 178c4762a1bSJed Brown 179c4762a1bSJed Brown if (test_bcast) { /* broadcast rootdata into leafdata */ 180c4762a1bSJed Brown PetscInt *rootdata,*leafdata; 1812d4ee042Sprj- /* Allocate space for send and receive buffers. This example communicates PetscInt, but other types, including 182c4762a1bSJed Brown * user-defined structures, could also be used. */ 183c4762a1bSJed Brown ierr = PetscMalloc2(nrootsalloc,&rootdata,nleavesalloc,&leafdata);CHKERRQ(ierr); 184c4762a1bSJed Brown /* Set rootdata buffer to be broadcast */ 185c4762a1bSJed Brown for (i=0; i<nrootsalloc; i++) rootdata[i] = -1; 186c4762a1bSJed Brown for (i=0; i<nroots; i++) rootdata[i*stride] = 100*(rank+1) + i; 187c4762a1bSJed Brown /* Initialize local buffer, these values are never used. */ 188c4762a1bSJed Brown for (i=0; i<nleavesalloc; i++) leafdata[i] = -1; 189c4762a1bSJed Brown /* Broadcast entries from rootdata to leafdata. Computation or other communication can be performed between the begin and end calls. */ 190ad227feaSJunchao Zhang ierr = PetscSFBcastBegin(sf,MPIU_INT,rootdata,leafdata,MPI_REPLACE);CHKERRQ(ierr); 191ad227feaSJunchao Zhang ierr = PetscSFBcastEnd(sf,MPIU_INT,rootdata,leafdata,MPI_REPLACE);CHKERRQ(ierr); 192c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Bcast Rootdata\n");CHKERRQ(ierr); 193c4762a1bSJed Brown ierr = PetscIntView(nrootsalloc,rootdata,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 194c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Bcast Leafdata\n");CHKERRQ(ierr); 195c4762a1bSJed Brown ierr = PetscIntView(nleavesalloc,leafdata,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 196c4762a1bSJed Brown ierr = PetscFree2(rootdata,leafdata);CHKERRQ(ierr); 197c4762a1bSJed Brown } 198c4762a1bSJed Brown 199c4762a1bSJed Brown if (test_bcast && test_char) { /* Bcast with char */ 200c4762a1bSJed Brown PetscInt len; 201c4762a1bSJed Brown char buf[256]; 202c4762a1bSJed Brown char *rootdata,*leafdata; 203c4762a1bSJed Brown ierr = PetscMalloc2(nrootsalloc,&rootdata,nleavesalloc,&leafdata);CHKERRQ(ierr); 204c4762a1bSJed Brown /* Set rootdata buffer to be broadcast */ 205c4762a1bSJed Brown for (i=0; i<nrootsalloc; i++) rootdata[i] = '*'; 206c4762a1bSJed 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 */ 207c4762a1bSJed Brown /* Initialize local buffer, these values are never used. */ 208c4762a1bSJed Brown for (i=0; i<nleavesalloc; i++) leafdata[i] = '?'; 209c4762a1bSJed Brown 210ad227feaSJunchao Zhang ierr = PetscSFBcastBegin(sf,MPI_CHAR,rootdata,leafdata,MPI_REPLACE);CHKERRQ(ierr); 211ad227feaSJunchao Zhang ierr = PetscSFBcastEnd(sf,MPI_CHAR,rootdata,leafdata,MPI_REPLACE);CHKERRQ(ierr); 212c4762a1bSJed Brown 213c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Bcast Rootdata in type of char\n");CHKERRQ(ierr); 214c4762a1bSJed Brown len = 0; ierr = PetscSNPrintf(buf,256,"%4d:",rank);CHKERRQ(ierr); len += 5; 215c4762a1bSJed Brown for (i=0; i<nrootsalloc; i++) {ierr = PetscSNPrintf(buf+len,256-len,"%5c",rootdata[i]);CHKERRQ(ierr); len += 5;} 216c4762a1bSJed Brown ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD,"%s\n",buf);CHKERRQ(ierr); 217c4762a1bSJed Brown ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT);CHKERRQ(ierr); 218c4762a1bSJed Brown 219c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Bcast Leafdata in type of char\n");CHKERRQ(ierr); 220c4762a1bSJed Brown len = 0; ierr = PetscSNPrintf(buf,256,"%4d:",rank);CHKERRQ(ierr); len += 5; 221c4762a1bSJed Brown for (i=0; i<nleavesalloc; i++) {ierr = PetscSNPrintf(buf+len,256-len,"%5c",leafdata[i]);CHKERRQ(ierr); len += 5;} 222c4762a1bSJed Brown ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD,"%s\n",buf);CHKERRQ(ierr); 223c4762a1bSJed Brown ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT);CHKERRQ(ierr); 224c4762a1bSJed Brown 225c4762a1bSJed Brown ierr = PetscFree2(rootdata,leafdata);CHKERRQ(ierr); 226c4762a1bSJed Brown } 227c4762a1bSJed Brown 228c4762a1bSJed Brown if (test_bcastop) { /* Reduce rootdata into leafdata */ 229c4762a1bSJed Brown PetscInt *rootdata,*leafdata; 230c4762a1bSJed Brown /* Allocate space for send and receive buffers. This example communicates PetscInt, but other types, including 231c4762a1bSJed Brown * user-defined structures, could also be used. */ 232c4762a1bSJed Brown ierr = PetscMalloc2(nrootsalloc,&rootdata,nleavesalloc,&leafdata);CHKERRQ(ierr); 233c4762a1bSJed Brown /* Set rootdata buffer to be broadcast */ 234c4762a1bSJed Brown for (i=0; i<nrootsalloc; i++) rootdata[i] = -1; 235c4762a1bSJed Brown for (i=0; i<nroots; i++) rootdata[i*stride] = 100*(rank+1) + i; 236c4762a1bSJed Brown /* Set leaf values to reduce with */ 237c4762a1bSJed Brown for (i=0; i<nleavesalloc; i++) leafdata[i] = -10*(rank+1) - i; 238c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Pre-BcastAndOp Leafdata\n");CHKERRQ(ierr); 239c4762a1bSJed Brown ierr = PetscIntView(nleavesalloc,leafdata,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 240c4762a1bSJed Brown /* Broadcast entries from rootdata to leafdata. Computation or other communication can be performed between the begin and end calls. */ 241ad227feaSJunchao Zhang ierr = PetscSFBcastBegin(sf,MPIU_INT,rootdata,leafdata,mop);CHKERRQ(ierr); 242ad227feaSJunchao Zhang ierr = PetscSFBcastEnd(sf,MPIU_INT,rootdata,leafdata,mop);CHKERRQ(ierr); 243c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## BcastAndOp Rootdata\n");CHKERRQ(ierr); 244c4762a1bSJed Brown ierr = PetscIntView(nrootsalloc,rootdata,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 245c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## BcastAndOp Leafdata\n");CHKERRQ(ierr); 246c4762a1bSJed Brown ierr = PetscIntView(nleavesalloc,leafdata,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 247c4762a1bSJed Brown ierr = PetscFree2(rootdata,leafdata);CHKERRQ(ierr); 248c4762a1bSJed Brown } 249c4762a1bSJed Brown 250c4762a1bSJed Brown if (test_reduce) { /* Reduce leafdata into rootdata */ 251c4762a1bSJed Brown PetscInt *rootdata,*leafdata; 252c4762a1bSJed Brown ierr = PetscMalloc2(nrootsalloc,&rootdata,nleavesalloc,&leafdata);CHKERRQ(ierr); 253c4762a1bSJed Brown /* Initialize rootdata buffer in which the result of the reduction will appear. */ 254c4762a1bSJed Brown for (i=0; i<nrootsalloc; i++) rootdata[i] = -1; 255c4762a1bSJed Brown for (i=0; i<nroots; i++) rootdata[i*stride] = 100*(rank+1) + i; 256c4762a1bSJed Brown /* Set leaf values to reduce. */ 257c4762a1bSJed Brown for (i=0; i<nleavesalloc; i++) leafdata[i] = -1; 258c4762a1bSJed Brown for (i=0; i<nleaves; i++) leafdata[i*stride] = 1000*(rank+1) + 10*i; 259c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Pre-Reduce Rootdata\n");CHKERRQ(ierr); 260c4762a1bSJed Brown ierr = PetscIntView(nrootsalloc,rootdata,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 261c4762a1bSJed Brown /* Perform reduction. Computation or other communication can be performed between the begin and end calls. 262c4762a1bSJed Brown * This example sums the values, but other MPI_Ops can be used (e.g MPI_MAX, MPI_PROD). */ 263c4762a1bSJed Brown ierr = PetscSFReduceBegin(sf,MPIU_INT,leafdata,rootdata,mop);CHKERRQ(ierr); 264c4762a1bSJed Brown ierr = PetscSFReduceEnd(sf,MPIU_INT,leafdata,rootdata,mop);CHKERRQ(ierr); 265c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Reduce Leafdata\n");CHKERRQ(ierr); 266c4762a1bSJed Brown ierr = PetscIntView(nleavesalloc,leafdata,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 267c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Reduce Rootdata\n");CHKERRQ(ierr); 268c4762a1bSJed Brown ierr = PetscIntView(nrootsalloc,rootdata,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 269c4762a1bSJed Brown ierr = PetscFree2(rootdata,leafdata);CHKERRQ(ierr); 270c4762a1bSJed Brown } 271c4762a1bSJed Brown 272c4762a1bSJed Brown if (test_reduce && test_char) { /* Reduce with signed char */ 273c4762a1bSJed Brown PetscInt len; 274c4762a1bSJed Brown char buf[256]; 275c4762a1bSJed Brown signed char *rootdata,*leafdata; 276c4762a1bSJed Brown ierr = PetscMalloc2(nrootsalloc,&rootdata,nleavesalloc,&leafdata);CHKERRQ(ierr); 277c4762a1bSJed Brown /* Initialize rootdata buffer in which the result of the reduction will appear. */ 278c4762a1bSJed Brown for (i=0; i<nrootsalloc; i++) rootdata[i] = -1; 279c4762a1bSJed Brown for (i=0; i<nroots; i++) rootdata[i*stride] = 10*(rank+1) + i; 280c4762a1bSJed Brown /* Set leaf values to reduce. */ 281c4762a1bSJed Brown for (i=0; i<nleavesalloc; i++) leafdata[i] = -1; 282c4762a1bSJed Brown for (i=0; i<nleaves; i++) leafdata[i*stride] = 50*(rank+1) + 10*i; 283c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Pre-Reduce Rootdata in type of signed char\n");CHKERRQ(ierr); 284c4762a1bSJed Brown 285c4762a1bSJed Brown len = 0; ierr = PetscSNPrintf(buf,256,"%4d:",rank);CHKERRQ(ierr); len += 5; 286c4762a1bSJed Brown for (i=0; i<nrootsalloc; i++) {ierr = PetscSNPrintf(buf+len,256-len,"%5d",rootdata[i]);CHKERRQ(ierr); len += 5;} 287c4762a1bSJed Brown ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD,"%s\n",buf);CHKERRQ(ierr); 288c4762a1bSJed Brown ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT);CHKERRQ(ierr); 289c4762a1bSJed Brown 290c4762a1bSJed Brown /* Using MPI_CHAR should trigger an error since MPI standard does not support reduction on MPI_CHAR. 291c4762a1bSJed Brown Testing with -test_op max, one can see the sign does take effect in MPI_MAX. 292c4762a1bSJed Brown */ 293c4762a1bSJed Brown ierr = PetscSFReduceBegin(sf,MPI_SIGNED_CHAR,leafdata,rootdata,mop);CHKERRQ(ierr); 294c4762a1bSJed Brown ierr = PetscSFReduceEnd(sf,MPI_SIGNED_CHAR,leafdata,rootdata,mop);CHKERRQ(ierr); 295c4762a1bSJed Brown 296c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Reduce Leafdata in type of signed char\n");CHKERRQ(ierr); 297c4762a1bSJed Brown len = 0; ierr = PetscSNPrintf(buf,256,"%4d:",rank);CHKERRQ(ierr); len += 5; 298c4762a1bSJed Brown for (i=0; i<nleavesalloc; i++) {ierr = PetscSNPrintf(buf+len,256-len,"%5d",leafdata[i]);CHKERRQ(ierr); len += 5;} 299c4762a1bSJed Brown ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD,"%s\n",buf);CHKERRQ(ierr); 300c4762a1bSJed Brown ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT);CHKERRQ(ierr); 301c4762a1bSJed Brown 302c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Reduce Rootdata in type of signed char\n");CHKERRQ(ierr); 303c4762a1bSJed Brown len = 0; ierr = PetscSNPrintf(buf,256,"%4d:",rank);CHKERRQ(ierr); len += 5; 304c4762a1bSJed Brown for (i=0; i<nrootsalloc; i++) {ierr = PetscSNPrintf(buf+len,256-len,"%5d",rootdata[i]);CHKERRQ(ierr); len += 5;} 305c4762a1bSJed Brown ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD,"%s\n",buf);CHKERRQ(ierr); 306c4762a1bSJed Brown ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT);CHKERRQ(ierr); 307c4762a1bSJed Brown 308c4762a1bSJed Brown ierr = PetscFree2(rootdata,leafdata);CHKERRQ(ierr); 309c4762a1bSJed Brown } 310c4762a1bSJed Brown 311c4762a1bSJed Brown if (test_reduce && test_char) { /* Reduce with unsigned char */ 312c4762a1bSJed Brown PetscInt len; 313c4762a1bSJed Brown char buf[256]; 314c4762a1bSJed Brown unsigned char *rootdata,*leafdata; 315c4762a1bSJed Brown ierr = PetscMalloc2(nrootsalloc,&rootdata,nleavesalloc,&leafdata);CHKERRQ(ierr); 316c4762a1bSJed Brown /* Initialize rootdata buffer in which the result of the reduction will appear. */ 317c4762a1bSJed Brown for (i=0; i<nrootsalloc; i++) rootdata[i] = -1; 318c4762a1bSJed Brown for (i=0; i<nroots; i++) rootdata[i*stride] = 10*(rank+1) + i; 319c4762a1bSJed Brown /* Set leaf values to reduce. */ 320c4762a1bSJed Brown for (i=0; i<nleavesalloc; i++) leafdata[i] = -1; 321c4762a1bSJed Brown for (i=0; i<nleaves; i++) leafdata[i*stride] = 50*(rank+1) + 10*i; 322c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Pre-Reduce Rootdata in type of unsigned char\n");CHKERRQ(ierr); 323c4762a1bSJed Brown 324c4762a1bSJed Brown len = 0; ierr = PetscSNPrintf(buf,256,"%4d:",rank);CHKERRQ(ierr); len += 5; 325c4762a1bSJed Brown for (i=0; i<nrootsalloc; i++) {ierr = PetscSNPrintf(buf+len,256-len,"%5u",rootdata[i]);CHKERRQ(ierr); len += 5;} 326c4762a1bSJed Brown ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD,"%s\n",buf);CHKERRQ(ierr); 327c4762a1bSJed Brown ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT);CHKERRQ(ierr); 328c4762a1bSJed Brown 329c4762a1bSJed Brown /* Using MPI_CHAR should trigger an error since MPI standard does not support reduction on MPI_CHAR. 330c4762a1bSJed Brown Testing with -test_op max, one can see the sign does take effect in MPI_MAX. 331c4762a1bSJed Brown */ 332c4762a1bSJed Brown ierr = PetscSFReduceBegin(sf,MPI_UNSIGNED_CHAR,leafdata,rootdata,mop);CHKERRQ(ierr); 333c4762a1bSJed Brown ierr = PetscSFReduceEnd(sf,MPI_UNSIGNED_CHAR,leafdata,rootdata,mop);CHKERRQ(ierr); 334c4762a1bSJed Brown 335c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Reduce Leafdata in type of unsigned char\n");CHKERRQ(ierr); 336c4762a1bSJed Brown len = 0; ierr = PetscSNPrintf(buf,256,"%4d:",rank);CHKERRQ(ierr); len += 5; 337c4762a1bSJed Brown for (i=0; i<nleavesalloc; i++) {ierr = PetscSNPrintf(buf+len,256-len,"%5u",leafdata[i]);CHKERRQ(ierr); len += 5;} 338c4762a1bSJed Brown ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD,"%s\n",buf);CHKERRQ(ierr); 339c4762a1bSJed Brown ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT);CHKERRQ(ierr); 340c4762a1bSJed Brown 341c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Reduce Rootdata in type of unsigned char\n");CHKERRQ(ierr); 342c4762a1bSJed Brown len = 0; ierr = PetscSNPrintf(buf,256,"%4d:",rank);CHKERRQ(ierr); len += 5; 343c4762a1bSJed Brown for (i=0; i<nrootsalloc; i++) {ierr = PetscSNPrintf(buf+len,256-len,"%5u",rootdata[i]);CHKERRQ(ierr); len += 5;} 344c4762a1bSJed Brown ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD,"%s\n",buf);CHKERRQ(ierr); 345c4762a1bSJed Brown ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT);CHKERRQ(ierr); 346c4762a1bSJed Brown 347c4762a1bSJed Brown ierr = PetscFree2(rootdata,leafdata);CHKERRQ(ierr); 348c4762a1bSJed Brown } 349c4762a1bSJed Brown 350c4762a1bSJed Brown if (test_degree) { 351c4762a1bSJed Brown const PetscInt *degree; 352c4762a1bSJed Brown ierr = PetscSFComputeDegreeBegin(sf,°ree);CHKERRQ(ierr); 353c4762a1bSJed Brown ierr = PetscSFComputeDegreeEnd(sf,°ree);CHKERRQ(ierr); 354c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Root degrees\n");CHKERRQ(ierr); 355c4762a1bSJed Brown ierr = PetscIntView(nrootsalloc,degree,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 356c4762a1bSJed Brown } 357c4762a1bSJed Brown 358c4762a1bSJed Brown if (test_fetchandop) { 359c4762a1bSJed Brown /* Cannot use text compare here because token ordering is not deterministic */ 360c4762a1bSJed Brown PetscInt *leafdata,*leafupdate,*rootdata; 361c4762a1bSJed Brown ierr = PetscMalloc3(nleavesalloc,&leafdata,nleavesalloc,&leafupdate,nrootsalloc,&rootdata);CHKERRQ(ierr); 362c4762a1bSJed Brown for (i=0; i<nleavesalloc; i++) leafdata[i] = -1; 363c4762a1bSJed Brown for (i=0; i<nleaves; i++) leafdata[i*stride] = 1; 364c4762a1bSJed Brown for (i=0; i<nrootsalloc; i++) rootdata[i] = -1; 365c4762a1bSJed Brown for (i=0; i<nroots; i++) rootdata[i*stride] = 0; 366c4762a1bSJed Brown ierr = PetscSFFetchAndOpBegin(sf,MPIU_INT,rootdata,leafdata,leafupdate,mop);CHKERRQ(ierr); 367c4762a1bSJed Brown ierr = PetscSFFetchAndOpEnd(sf,MPIU_INT,rootdata,leafdata,leafupdate,mop);CHKERRQ(ierr); 368c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Rootdata (sum of 1 from each leaf)\n");CHKERRQ(ierr); 369c4762a1bSJed Brown ierr = PetscIntView(nrootsalloc,rootdata,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 370c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Leafupdate (value at roots prior to my atomic update)\n");CHKERRQ(ierr); 371c4762a1bSJed Brown ierr = PetscIntView(nleavesalloc,leafupdate,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 372c4762a1bSJed Brown ierr = PetscFree3(leafdata,leafupdate,rootdata);CHKERRQ(ierr); 373c4762a1bSJed Brown } 374c4762a1bSJed Brown 375c4762a1bSJed Brown if (test_gather) { 376c4762a1bSJed Brown const PetscInt *degree; 377c4762a1bSJed Brown PetscInt inedges,*indata,*outdata; 378c4762a1bSJed Brown ierr = PetscSFComputeDegreeBegin(sf,°ree);CHKERRQ(ierr); 379c4762a1bSJed Brown ierr = PetscSFComputeDegreeEnd(sf,°ree);CHKERRQ(ierr); 380c4762a1bSJed Brown for (i=0,inedges=0; i<nrootsalloc; i++) inedges += degree[i]; 381c4762a1bSJed Brown ierr = PetscMalloc2(inedges,&indata,nleavesalloc,&outdata);CHKERRQ(ierr); 382c4762a1bSJed Brown for (i=0; i<nleavesalloc; i++) outdata[i] = -1; 383c4762a1bSJed Brown for (i=0; i<nleaves; i++) outdata[i*stride] = 1000*(rank+1) + i; 384c4762a1bSJed Brown ierr = PetscSFGatherBegin(sf,MPIU_INT,outdata,indata);CHKERRQ(ierr); 385c4762a1bSJed Brown ierr = PetscSFGatherEnd(sf,MPIU_INT,outdata,indata);CHKERRQ(ierr); 386c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Gathered data at multi-roots from leaves\n");CHKERRQ(ierr); 387c4762a1bSJed Brown ierr = PetscIntView(inedges,indata,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 388c4762a1bSJed Brown ierr = PetscFree2(indata,outdata);CHKERRQ(ierr); 389c4762a1bSJed Brown } 390c4762a1bSJed Brown 391c4762a1bSJed Brown if (test_scatter) { 392c4762a1bSJed Brown const PetscInt *degree; 393c4762a1bSJed Brown PetscInt j,count,inedges,*indata,*outdata; 394c4762a1bSJed Brown ierr = PetscSFComputeDegreeBegin(sf,°ree);CHKERRQ(ierr); 395c4762a1bSJed Brown ierr = PetscSFComputeDegreeEnd(sf,°ree);CHKERRQ(ierr); 396c4762a1bSJed Brown for (i=0,inedges=0; i<nrootsalloc; i++) inedges += degree[i]; 397c4762a1bSJed Brown ierr = PetscMalloc2(inedges,&indata,nleavesalloc,&outdata);CHKERRQ(ierr); 398c4762a1bSJed Brown for (i=0; i<nleavesalloc; i++) outdata[i] = -1; 399c4762a1bSJed Brown for (i=0,count=0; i<nrootsalloc; i++) { 400c4762a1bSJed Brown for (j=0; j<degree[i]; j++) indata[count++] = 1000*(rank+1) + 100*i + j; 401c4762a1bSJed Brown } 402c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Data at multi-roots, to scatter to leaves\n");CHKERRQ(ierr); 403c4762a1bSJed Brown ierr = PetscIntView(inedges,indata,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 404c4762a1bSJed Brown 405c4762a1bSJed Brown ierr = PetscSFScatterBegin(sf,MPIU_INT,indata,outdata);CHKERRQ(ierr); 406c4762a1bSJed Brown ierr = PetscSFScatterEnd(sf,MPIU_INT,indata,outdata);CHKERRQ(ierr); 407c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Scattered data at leaves\n");CHKERRQ(ierr); 408c4762a1bSJed Brown ierr = PetscIntView(nleavesalloc,outdata,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 409c4762a1bSJed Brown ierr = PetscFree2(indata,outdata);CHKERRQ(ierr); 410c4762a1bSJed Brown } 411c4762a1bSJed Brown 412c4762a1bSJed Brown if (test_embed) { 413c4762a1bSJed Brown const PetscInt nroots = 1 + (PetscInt) !rank; 414c4762a1bSJed Brown PetscInt selected[2]; 415c4762a1bSJed Brown PetscSF esf; 416c4762a1bSJed Brown 417c4762a1bSJed Brown selected[0] = stride; 418c4762a1bSJed Brown selected[1] = 2*stride; 41972502a1fSJunchao Zhang ierr = PetscSFCreateEmbeddedRootSF(sf,nroots,selected,&esf);CHKERRQ(ierr); 420c4762a1bSJed Brown ierr = PetscSFSetUp(esf);CHKERRQ(ierr); 421c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Embedded PetscSF\n");CHKERRQ(ierr); 422c4762a1bSJed Brown ierr = PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);CHKERRQ(ierr); 423c4762a1bSJed Brown ierr = PetscSFView(esf,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 424c4762a1bSJed Brown ierr = PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 425c4762a1bSJed Brown ierr = PetscSFDestroy(&esf);CHKERRQ(ierr); 426c4762a1bSJed Brown } 427c4762a1bSJed Brown 428c4762a1bSJed Brown if (test_invert) { 429c4762a1bSJed Brown const PetscInt *degree; 430c4762a1bSJed Brown PetscInt *mRootsOrigNumbering; 431c4762a1bSJed Brown PetscInt inedges; 432c4762a1bSJed Brown PetscSF msf,imsf; 433c4762a1bSJed Brown 434c4762a1bSJed Brown ierr = PetscSFGetMultiSF(sf,&msf);CHKERRQ(ierr); 435c4762a1bSJed Brown ierr = PetscSFCreateInverseSF(msf,&imsf);CHKERRQ(ierr); 436c4762a1bSJed Brown ierr = PetscSFSetUp(msf);CHKERRQ(ierr); 437c4762a1bSJed Brown ierr = PetscSFSetUp(imsf);CHKERRQ(ierr); 438c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Multi-SF\n");CHKERRQ(ierr); 439c4762a1bSJed Brown ierr = PetscSFView(msf,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 440c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Multi-SF roots indices in original SF roots numbering\n");CHKERRQ(ierr); 441c4762a1bSJed Brown ierr = PetscSFComputeDegreeBegin(sf,°ree);CHKERRQ(ierr); 442c4762a1bSJed Brown ierr = PetscSFComputeDegreeEnd(sf,°ree);CHKERRQ(ierr); 443c4762a1bSJed Brown ierr = PetscSFComputeMultiRootOriginalNumbering(sf,degree,&inedges,&mRootsOrigNumbering);CHKERRQ(ierr); 444c4762a1bSJed Brown ierr = PetscIntView(inedges,mRootsOrigNumbering,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 445c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Inverse of Multi-SF\n");CHKERRQ(ierr); 446c4762a1bSJed Brown ierr = PetscSFView(imsf,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 447c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Inverse of Multi-SF, original numbering\n");CHKERRQ(ierr); 448c4762a1bSJed Brown ierr = PetscSFViewCustomLocals_Private(imsf,mRootsOrigNumbering,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 449c4762a1bSJed Brown ierr = PetscSFDestroy(&imsf);CHKERRQ(ierr); 450c4762a1bSJed Brown ierr = PetscFree(mRootsOrigNumbering);CHKERRQ(ierr); 451c4762a1bSJed Brown } 452c4762a1bSJed Brown 453c4762a1bSJed Brown /* Clean storage for star forest. */ 454c4762a1bSJed Brown ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); 455c4762a1bSJed Brown ierr = PetscFinalize(); 456c4762a1bSJed Brown return ierr; 457c4762a1bSJed Brown } 458c4762a1bSJed Brown 459c4762a1bSJed Brown /*TEST 460c4762a1bSJed Brown 461c4762a1bSJed Brown test: 462c4762a1bSJed Brown nsize: 4 463c4762a1bSJed Brown filter: grep -v "type" | grep -v "sort" 464c4762a1bSJed Brown args: -test_bcast -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}} 465*dfd57a17SPierre Jolivet requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW) 466c4762a1bSJed Brown 467c4762a1bSJed Brown test: 468c4762a1bSJed Brown suffix: 2 469c4762a1bSJed Brown nsize: 4 470c4762a1bSJed Brown filter: grep -v "type" | grep -v "sort" 471c4762a1bSJed Brown args: -test_reduce -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}} 472*dfd57a17SPierre Jolivet requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW) 473c4762a1bSJed Brown 474c4762a1bSJed Brown test: 475c4762a1bSJed Brown suffix: 2_basic 476c4762a1bSJed Brown nsize: 4 477c4762a1bSJed Brown args: -test_reduce -sf_type basic 478c4762a1bSJed Brown 479c4762a1bSJed Brown test: 480c4762a1bSJed Brown suffix: 3 481c4762a1bSJed Brown nsize: 4 482c4762a1bSJed Brown filter: grep -v "type" | grep -v "sort" 483c4762a1bSJed Brown args: -test_degree -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}} 484*dfd57a17SPierre Jolivet requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW) 485c4762a1bSJed Brown 486c4762a1bSJed Brown test: 487c4762a1bSJed Brown suffix: 3_basic 488c4762a1bSJed Brown nsize: 4 489c4762a1bSJed Brown args: -test_degree -sf_type basic 490c4762a1bSJed Brown 491c4762a1bSJed Brown test: 492c4762a1bSJed Brown suffix: 4 493c4762a1bSJed Brown nsize: 4 494c4762a1bSJed Brown filter: grep -v "type" | grep -v "sort" 495c4762a1bSJed Brown args: -test_gather -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}} 496*dfd57a17SPierre Jolivet requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW) 497c4762a1bSJed Brown 498c4762a1bSJed Brown test: 499c4762a1bSJed Brown suffix: 4_basic 500c4762a1bSJed Brown nsize: 4 501c4762a1bSJed Brown args: -test_gather -sf_type basic 502c4762a1bSJed Brown 503c4762a1bSJed Brown test: 504c4762a1bSJed Brown suffix: 4_stride 505c4762a1bSJed Brown nsize: 4 506c4762a1bSJed Brown args: -test_gather -sf_type basic -stride 2 507c4762a1bSJed Brown 508c4762a1bSJed Brown test: 509c4762a1bSJed Brown suffix: 5 510c4762a1bSJed Brown nsize: 4 511c4762a1bSJed Brown filter: grep -v "type" | grep -v "sort" 512c4762a1bSJed Brown args: -test_scatter -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}} 513*dfd57a17SPierre Jolivet requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW) 514c4762a1bSJed Brown 515c4762a1bSJed Brown test: 516c4762a1bSJed Brown suffix: 5_basic 517c4762a1bSJed Brown nsize: 4 518c4762a1bSJed Brown args: -test_scatter -sf_type basic 519c4762a1bSJed Brown 520c4762a1bSJed Brown test: 521c4762a1bSJed Brown suffix: 5_stride 522c4762a1bSJed Brown nsize: 4 523c4762a1bSJed Brown args: -test_scatter -sf_type basic -stride 2 524c4762a1bSJed Brown 525c4762a1bSJed Brown test: 526c4762a1bSJed Brown suffix: 6 527c4762a1bSJed Brown nsize: 4 528c4762a1bSJed Brown filter: grep -v "type" | grep -v "sort" 529c20d7725SJed Brown # No -sf_window_flavor dynamic due to bug https://gitlab.com/petsc/petsc/issues/555 530c20d7725SJed Brown args: -test_embed -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create allocate}} 531*dfd57a17SPierre Jolivet requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW) 532c4762a1bSJed Brown 533c4762a1bSJed Brown test: 534c4762a1bSJed Brown suffix: 6_basic 535c4762a1bSJed Brown nsize: 4 536c4762a1bSJed Brown args: -test_embed -sf_type basic 537c4762a1bSJed Brown 538c4762a1bSJed Brown test: 539c4762a1bSJed Brown suffix: 7 540c4762a1bSJed Brown nsize: 4 541c4762a1bSJed Brown filter: grep -v "type" | grep -v "sort" 542c4762a1bSJed Brown args: -test_invert -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}} 543*dfd57a17SPierre Jolivet requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW) 544c4762a1bSJed Brown 545c4762a1bSJed Brown test: 546c4762a1bSJed Brown suffix: 7_basic 547c4762a1bSJed Brown nsize: 4 548c4762a1bSJed Brown args: -test_invert -sf_type basic 549c4762a1bSJed Brown 550c4762a1bSJed Brown test: 551c4762a1bSJed Brown suffix: basic 552c4762a1bSJed Brown nsize: 4 553c4762a1bSJed Brown args: -test_bcast -sf_type basic 554c4762a1bSJed Brown output_file: output/ex1_1_basic.out 555c4762a1bSJed Brown 556c4762a1bSJed Brown test: 557c4762a1bSJed Brown suffix: bcastop_basic 558c4762a1bSJed Brown nsize: 4 559c4762a1bSJed Brown args: -test_bcastop -sf_type basic 560c4762a1bSJed Brown output_file: output/ex1_bcastop_basic.out 561c4762a1bSJed Brown 562c4762a1bSJed Brown test: 563c4762a1bSJed Brown suffix: 8 564c4762a1bSJed Brown nsize: 3 565c4762a1bSJed Brown filter: grep -v "type" | grep -v "sort" 566c4762a1bSJed Brown args: -test_bcast -test_sf_distribute -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}} 567*dfd57a17SPierre Jolivet requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW) 568c4762a1bSJed Brown 569c4762a1bSJed Brown test: 570c4762a1bSJed Brown suffix: 8_basic 571c4762a1bSJed Brown nsize: 3 572c4762a1bSJed Brown args: -test_bcast -test_sf_distribute -sf_type basic 573c4762a1bSJed Brown 574c4762a1bSJed Brown test: 575c4762a1bSJed Brown suffix: 9_char 576c4762a1bSJed Brown nsize: 4 577c4762a1bSJed Brown args: -sf_type basic -test_bcast -test_reduce -test_op max -test_char 578c4762a1bSJed Brown 579c4762a1bSJed Brown # Here we do not test -sf_window_flavor dynamic since it is designed for repeated SFs with few different rootdata pointers 580c4762a1bSJed Brown test: 581c4762a1bSJed Brown suffix: 10 582c4762a1bSJed Brown filter: grep -v "type" | grep -v "sort" 583c4762a1bSJed Brown nsize: 4 584c4762a1bSJed Brown args: -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create allocate}} -test_all -test_bcastop 0 -test_fetchandop 0 585*dfd57a17SPierre Jolivet requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW) 586c4762a1bSJed Brown 587c4762a1bSJed Brown # The nightly test suite with MPICH uses ch3:sock, which is broken when winsize == 0 in some of the processes 588c4762a1bSJed Brown test: 589c4762a1bSJed Brown suffix: 10_shared 590c4762a1bSJed Brown output_file: output/ex1_10.out 591c4762a1bSJed Brown filter: grep -v "type" | grep -v "sort" 592c4762a1bSJed Brown nsize: 4 593c4762a1bSJed Brown args: -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor shared -test_all -test_bcastop 0 -test_fetchandop 0 594*dfd57a17SPierre 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) 595c4762a1bSJed Brown 596c4762a1bSJed Brown test: 597c4762a1bSJed Brown suffix: 10_basic 598c4762a1bSJed Brown nsize: 4 599c4762a1bSJed Brown args: -sf_type basic -test_all -test_bcastop 0 -test_fetchandop 0 600c4762a1bSJed Brown 601c4762a1bSJed Brown TEST*/ 602