1936c5a86SJed Brown /* 2936c5a86SJed Brown A star forest (SF) describes a communication pattern 3936c5a86SJed Brown */ 4936c5a86SJed Brown #if !defined(__PETSCSF_H) 5936c5a86SJed Brown #define __PETSCSF_H 62c8e378dSBarry Smith #include <petscsys.h> 70c312b8eSJed Brown #include <petscsftypes.h> 8936c5a86SJed Brown 9014dd563SJed Brown PETSC_EXTERN PetscClassId PETSCSF_CLASSID; 10936c5a86SJed Brown 115af141bcSJed Brown /*J 125af141bcSJed Brown PetscSFType - String with the name of a PetscSF method or the creation function 135af141bcSJed Brown with an optional dynamic library name, for example 145af141bcSJed Brown http://www.mcs.anl.gov/petsc/lib.so:mysfcreate() 155af141bcSJed Brown 165af141bcSJed Brown Level: beginner 175af141bcSJed Brown 18ed658588SBarry Smith Notes: The two approaches provided are 19ed658588SBarry Smith $ PETSCSFBASIC which uses MPI 1 message passing to perform the communication and 20ed658588SBarry Smith $ PETSCSFWINDOW which uses MPI 2 one-sided operations to perform the communication, this may be more efficient, 21ed658588SBarry Smith $ but may not be available for all MPI distributions. In particular OpenMPI has bugs in its one-sided 22ed658588SBarry Smith $ operations that prevent its use. 23ed658588SBarry Smith 245af141bcSJed Brown .seealso: PetscSFSetType(), PetscSF 255af141bcSJed Brown J*/ 265af141bcSJed Brown typedef const char *PetscSFType; 27ac762476SJed Brown #define PETSCSFBASIC "basic" 28ed658588SBarry Smith #define PETSCSFWINDOW "window" 295af141bcSJed Brown 30936c5a86SJed Brown /*S 31936c5a86SJed Brown PetscSFNode - specifier of owner and index 32936c5a86SJed Brown 33936c5a86SJed Brown Level: beginner 34936c5a86SJed Brown 35936c5a86SJed Brown Concepts: indexing, stride, distribution 36936c5a86SJed Brown 37936c5a86SJed Brown .seealso: PetscSFSetGraph() 38936c5a86SJed Brown S*/ 39936c5a86SJed Brown typedef struct { 40936c5a86SJed Brown PetscInt rank; /* Rank of owner */ 41936c5a86SJed Brown PetscInt index; /* Index of node on rank */ 42936c5a86SJed Brown } PetscSFNode; 43936c5a86SJed Brown 44936c5a86SJed Brown /*E 455af141bcSJed Brown PetscSFWindowSyncType - Type of synchronization for PETSCSFWINDOW 46936c5a86SJed Brown 475af141bcSJed Brown $ PETSCSF_WINDOW_SYNC_FENCE - simplest model, synchronizing across communicator 485af141bcSJed Brown $ PETSCSF_WINDOW_SYNC_LOCK - passive model, less synchronous, requires less setup than PETSCSF_WINDOW_SYNC_ACTIVE, but may require more handshakes 495af141bcSJed Brown $ PETSCSF_WINDOW_SYNC_ACTIVE - active model, provides most information to MPI implementation, needs to construct 2-way process groups (more setup than PETSCSF_WINDOW_SYNC_LOCK) 50936c5a86SJed Brown 51e84a5f06SJed Brown Level: advanced 52936c5a86SJed Brown 53e84a5f06SJed Brown .seealso: PetscSFWindowSetSyncType(), PetscSFWindowGetSyncType() 54936c5a86SJed Brown E*/ 555af141bcSJed Brown typedef enum {PETSCSF_WINDOW_SYNC_FENCE,PETSCSF_WINDOW_SYNC_LOCK,PETSCSF_WINDOW_SYNC_ACTIVE} PetscSFWindowSyncType; 565af141bcSJed Brown PETSC_EXTERN const char *const PetscSFWindowSyncTypes[]; 57936c5a86SJed Brown 58e84a5f06SJed Brown /*E 59e84a5f06SJed Brown PetscSFDuplicateOption - Aspects to preserve when duplicating a PetscSF 60e84a5f06SJed Brown 61e84a5f06SJed Brown $ PETSCSF_DUPLICATE_CONFONLY - configuration only, user must call PetscSFSetGraph() 62e84a5f06SJed Brown $ PETSCSF_DUPLICATE_RANKS - communication ranks preserved, but different graph (allows simpler setup after calling PetscSFSetGraph()) 63e84a5f06SJed Brown $ PETSCSF_DUPLICATE_GRAPH - entire graph duplicated 64e84a5f06SJed Brown 65e84a5f06SJed Brown Level: beginner 66e84a5f06SJed Brown 67e84a5f06SJed Brown .seealso: PetscSFDuplicate() 68e84a5f06SJed Brown E*/ 69e84a5f06SJed Brown typedef enum {PETSCSF_DUPLICATE_CONFONLY,PETSCSF_DUPLICATE_RANKS,PETSCSF_DUPLICATE_GRAPH} PetscSFDuplicateOption; 70e84a5f06SJed Brown PETSC_EXTERN const char *const PetscSFDuplicateOptions[]; 71090c6444SJed Brown 72adc40e5bSBarry Smith PETSC_EXTERN PetscFunctionList PetscSFList; 73607a6623SBarry Smith PETSC_EXTERN PetscErrorCode PetscSFRegisterAll(void); 74bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode PetscSFRegister(const char[],PetscErrorCode (*)(PetscSF)); 755af141bcSJed Brown 76607a6623SBarry Smith PETSC_EXTERN PetscErrorCode PetscSFInitializePackage(void); 77014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscSFFinalizePackage(void); 78014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscSFCreate(MPI_Comm comm,PetscSF*); 79014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscSFDestroy(PetscSF*); 805af141bcSJed Brown PETSC_EXTERN PetscErrorCode PetscSFSetType(PetscSF,PetscSFType); 81014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscSFView(PetscSF,PetscViewer); 825af141bcSJed Brown PETSC_EXTERN PetscErrorCode PetscSFSetUp(PetscSF); 83014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscSFSetFromOptions(PetscSF); 84e84a5f06SJed Brown PETSC_EXTERN PetscErrorCode PetscSFDuplicate(PetscSF,PetscSFDuplicateOption,PetscSF*); 855af141bcSJed Brown PETSC_EXTERN PetscErrorCode PetscSFWindowSetSyncType(PetscSF,PetscSFWindowSyncType); 865af141bcSJed Brown PETSC_EXTERN PetscErrorCode PetscSFWindowGetSyncType(PetscSF,PetscSFWindowSyncType*); 87014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscSFSetRankOrder(PetscSF,PetscBool); 8863f4a732SJed Brown PETSC_EXTERN PetscErrorCode PetscSFSetGraph(PetscSF,PetscInt,PetscInt,const PetscInt*,PetscCopyMode,const PetscSFNode*,PetscCopyMode); 89014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscSFGetGraph(PetscSF,PetscInt *nroots,PetscInt *nleaves,const PetscInt **ilocal,const PetscSFNode **iremote); 90f723732fSJed Brown PETSC_EXTERN PetscErrorCode PetscSFGetLeafRange(PetscSF,PetscInt*,PetscInt*); 91014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscSFCreateEmbeddedSF(PetscSF,PetscInt nroots,const PetscInt *selected,PetscSF *newsf); 92014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscSFReset(PetscSF); 935b6cb184SJed Brown PETSC_EXTERN PetscErrorCode PetscSFGetRanks(PetscSF,PetscInt*,const PetscMPIInt**,const PetscInt**,const PetscInt**,const PetscInt**); 94014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscSFGetGroups(PetscSF,MPI_Group*,MPI_Group*); 95014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscSFGetMultiSF(PetscSF,PetscSF*); 96014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscSFCreateInverseSF(PetscSF,PetscSF*); 97936c5a86SJed Brown 98936c5a86SJed Brown /* broadcasts rootdata to leafdata */ 99894dd566SJed Brown PETSC_EXTERN PetscErrorCode PetscSFBcastBegin(PetscSF,MPI_Datatype,const void *rootdata,void *leafdata) 100894dd566SJed Brown PetscAttrMPIPointerWithType(3,2) PetscAttrMPIPointerWithType(4,2); 101894dd566SJed Brown PETSC_EXTERN PetscErrorCode PetscSFBcastEnd(PetscSF,MPI_Datatype,const void *rootdata,void *leafdata) 102894dd566SJed Brown PetscAttrMPIPointerWithType(3,2) PetscAttrMPIPointerWithType(4,2); 103936c5a86SJed Brown /* Reduce leafdata into rootdata using provided operation */ 104894dd566SJed Brown PETSC_EXTERN PetscErrorCode PetscSFReduceBegin(PetscSF,MPI_Datatype,const void *leafdata,void *rootdata,MPI_Op) 10519436ca2SJed Brown PetscAttrMPIPointerWithType(3,2) PetscAttrMPIPointerWithType(4,2); 106894dd566SJed Brown PETSC_EXTERN PetscErrorCode PetscSFReduceEnd(PetscSF,MPI_Datatype,const void *leafdata,void *rootdata,MPI_Op) 10719436ca2SJed Brown PetscAttrMPIPointerWithType(3,2) PetscAttrMPIPointerWithType(4,2); 108936c5a86SJed Brown /* Atomically modifies (using provided operation) rootdata using leafdata from each leaf, value at root at time of modification is returned in leafupdate. */ 109894dd566SJed Brown PETSC_EXTERN PetscErrorCode PetscSFFetchAndOpBegin(PetscSF,MPI_Datatype,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op) 110894dd566SJed Brown PetscAttrMPIPointerWithType(3,2) PetscAttrMPIPointerWithType(4,2) PetscAttrMPIPointerWithType(5,2); 111894dd566SJed Brown PETSC_EXTERN PetscErrorCode PetscSFFetchAndOpEnd(PetscSF,MPI_Datatype,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op) 112894dd566SJed Brown PetscAttrMPIPointerWithType(3,2) PetscAttrMPIPointerWithType(4,2) PetscAttrMPIPointerWithType(5,2); 113936c5a86SJed Brown /* Compute the degree of every root vertex (number of leaves in its star) */ 114014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscSFComputeDegreeBegin(PetscSF,const PetscInt **degree); 115014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscSFComputeDegreeEnd(PetscSF,const PetscInt **degree); 116936c5a86SJed Brown /* Concatenate data from all leaves into roots */ 117894dd566SJed Brown PETSC_EXTERN PetscErrorCode PetscSFGatherBegin(PetscSF,MPI_Datatype,const void *leafdata,void *multirootdata) 118894dd566SJed Brown PetscAttrMPIPointerWithType(3,2) PetscAttrMPIPointerWithType(4,2); 119894dd566SJed Brown PETSC_EXTERN PetscErrorCode PetscSFGatherEnd(PetscSF,MPI_Datatype,const void *leafdata,void *multirootdata) 120894dd566SJed Brown PetscAttrMPIPointerWithType(3,2) PetscAttrMPIPointerWithType(4,2); 121936c5a86SJed Brown /* Distribute distinct values to each leaf from roots */ 122894dd566SJed Brown PETSC_EXTERN PetscErrorCode PetscSFScatterBegin(PetscSF,MPI_Datatype,const void *multirootdata,void *leafdata) 123894dd566SJed Brown PetscAttrMPIPointerWithType(3,2) PetscAttrMPIPointerWithType(4,2); 124894dd566SJed Brown PETSC_EXTERN PetscErrorCode PetscSFScatterEnd(PetscSF,MPI_Datatype,const void *multirootdata,void *leafdata) 125894dd566SJed Brown PetscAttrMPIPointerWithType(3,2) PetscAttrMPIPointerWithType(4,2); 126936c5a86SJed Brown 127*8bfbc91cSJed Brown #if defined(MPI_REPLACE) 128*8bfbc91cSJed Brown # define MPIU_REPLACE MPI_REPLACE 129*8bfbc91cSJed Brown #else 130*8bfbc91cSJed Brown /* When using an old MPI such that MPI_REPLACE is not defined, we do not pass MPI_REPLACE to MPI at all. Instead, we 131*8bfbc91cSJed Brown * use it as a flag for our own reducer in the PETSCSFBASIC implementation. This could be any unique value unlikely to 132*8bfbc91cSJed Brown * collide with another MPI_Op so we'll just use the value that has been used by every version of MPICH since 133*8bfbc91cSJed Brown * MPICH2-1.0.6. */ 134*8bfbc91cSJed Brown # define MPIU_REPLACE (MPI_Op)(0x5800000d) 135*8bfbc91cSJed Brown #endif 136*8bfbc91cSJed Brown 137936c5a86SJed Brown #endif 138