static const char help[] = "Test star forest communication (PetscSF)\n\n"; /*T Description: A star is a simple tree with one root and zero or more leaves. A star forest is a union of disjoint stars. Many common communication patterns can be expressed as updates of rootdata using leafdata and vice-versa. This example creates a star forest, communicates values using the graph (see options for types of communication), views the graph, then destroys it. T*/ /* Include petscsf.h so we can use PetscSF objects. Note that this automatically includes petscsys.h. */ #include #include /* like PetscSFView() but with alternative array of local indices */ static PetscErrorCode PetscSFViewCustomLocals_Private(PetscSF sf,const PetscInt locals[],PetscViewer viewer) { const PetscSFNode *iremote; PetscInt i,nroots,nleaves,nranks; PetscMPIInt rank; PetscErrorCode ierr; PetscFunctionBeginUser; ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr); ierr = PetscSFGetGraph(sf,&nroots,&nleaves,NULL,&iremote);CHKERRQ(ierr); ierr = PetscSFGetRootRanks(sf,&nranks,NULL,NULL,NULL,NULL);CHKERRQ(ierr); ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Number of roots=%D, leaves=%D, remote ranks=%D\n",rank,nroots,nleaves,nranks);CHKERRQ(ierr); for (i=0; i 0); nleavesalloc = nleaves * stride; mine = NULL; if (stride > 1) { PetscInt i; ierr = PetscMalloc1(nleaves,&mine);CHKERRQ(ierr); for (i = 0; i < nleaves; i++) { mine[i] = stride * i; } } ierr = PetscMalloc1(nleaves,&remote);CHKERRQ(ierr); /* Left periodic neighbor */ remote[0].rank = (rank+size-1)%size; remote[0].index = 1 * stride; /* Right periodic neighbor */ remote[1].rank = (rank+1)%size; remote[1].index = 0 * stride; if (rank > 0) { /* All processes reference rank 0, index 1 */ remote[2].rank = 0; remote[2].index = 2 * stride; } } /* Create a star forest for communication. In this example, the leaf space is dense, so we pass NULL. */ ierr = PetscSFCreate(PETSC_COMM_WORLD,&sf);CHKERRQ(ierr); ierr = PetscSFSetFromOptions(sf);CHKERRQ(ierr); ierr = PetscSFSetGraph(sf,nrootsalloc,nleaves,mine,PETSC_OWN_POINTER,remote,PETSC_OWN_POINTER);CHKERRQ(ierr); ierr = PetscSFSetUp(sf);CHKERRQ(ierr); /* View graph, mostly useful for debugging purposes. */ ierr = PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);CHKERRQ(ierr); ierr = PetscSFView(sf,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); ierr = PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); if (test_bcast) { /* broadcast rootdata into leafdata */ PetscInt *rootdata,*leafdata; /* Allocate space for send and receive buffers. This example communicates PetscInt, but other types, including * user-defined structures, could also be used. */ ierr = PetscMalloc2(nrootsalloc,&rootdata,nleavesalloc,&leafdata);CHKERRQ(ierr); /* Set rootdata buffer to be broadcast */ for (i=0; i