xref: /petsc/include/petscsf.h (revision 014dd563d73e9fc78d056590fa6cf997782bf92d)
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
6936c5a86SJed Brown #include "petscsys.h"
7936c5a86SJed Brown 
8*014dd563SJed Brown PETSC_EXTERN PetscClassId PETSCSF_CLASSID;
9936c5a86SJed Brown 
10936c5a86SJed Brown /*S
11378b0e49SJed Brown    PetscSF - PETSc object for communication using star forests
12936c5a86SJed Brown 
13936c5a86SJed Brown    Level: intermediate
14936c5a86SJed Brown 
15936c5a86SJed Brown   Concepts: star forest
16936c5a86SJed Brown 
17936c5a86SJed Brown .seealso: PetscSFCreate()
18936c5a86SJed Brown S*/
19936c5a86SJed Brown typedef struct _p_PetscSF* PetscSF;
20936c5a86SJed Brown 
21936c5a86SJed Brown /*S
22936c5a86SJed Brown    PetscSFNode - specifier of owner and index
23936c5a86SJed Brown 
24936c5a86SJed Brown    Level: beginner
25936c5a86SJed Brown 
26936c5a86SJed Brown   Concepts: indexing, stride, distribution
27936c5a86SJed Brown 
28936c5a86SJed Brown .seealso: PetscSFSetGraph()
29936c5a86SJed Brown S*/
30936c5a86SJed Brown typedef struct {
31936c5a86SJed Brown   PetscInt rank;                /* Rank of owner */
32936c5a86SJed Brown   PetscInt index;               /* Index of node on rank */
33936c5a86SJed Brown } PetscSFNode;
34936c5a86SJed Brown 
35936c5a86SJed Brown /*E
36936c5a86SJed Brown     PetscSFSynchronizationType - Type of synchronization for PetscSF
37936c5a86SJed Brown 
38936c5a86SJed Brown $  PETSCSF_SYNCHRONIZATION_FENCE - simplest model, synchronizing across communicator
39936c5a86SJed Brown $  PETSCSF_SYNCHRONIZATION_LOCK - passive model, less synchronous, requires less setup than PETSCSF_SYNCHRONIZATION_ACTIVE, but may require more handshakes
40936c5a86SJed Brown $  PETSCSF_SYNCHRONIZATION_ACTIVE - active model, provides most information to MPI implementation, needs to construct 2-way process groups (more setup than PETSCSF_SYNCHRONIZATION_LOCK)
41936c5a86SJed Brown 
42936c5a86SJed Brown    Level: beginner
43936c5a86SJed Brown 
44936c5a86SJed Brown .seealso: PetscSFSetSynchronizationType()
45936c5a86SJed Brown E*/
46936c5a86SJed Brown typedef enum {PETSCSF_SYNCHRONIZATION_FENCE,PETSCSF_SYNCHRONIZATION_LOCK,PETSCSF_SYNCHRONIZATION_ACTIVE} PetscSFSynchronizationType;
47*014dd563SJed Brown PETSC_EXTERN const char *const PetscSFSynchronizationTypes[];
48936c5a86SJed Brown 
49090c6444SJed Brown #if !defined(PETSC_HAVE_MPI_WIN_CREATE) /* The intent here is to be able to compile even without a complete MPI. */
50090c6444SJed Brown typedef struct MPI_Win_MISSING *MPI_Win;
51090c6444SJed Brown #endif
52090c6444SJed Brown 
53*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscSFInitializePackage(const char*);
54*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscSFFinalizePackage(void);
55*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscSFCreate(MPI_Comm comm,PetscSF*);
56*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscSFDestroy(PetscSF*);
57*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscSFView(PetscSF,PetscViewer);
58*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscSFSetFromOptions(PetscSF);
59*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscSFSetSynchronizationType(PetscSF,PetscSFSynchronizationType);
60*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscSFSetRankOrder(PetscSF,PetscBool);
61*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscSFSetGraph(PetscSF,PetscInt nroots,PetscInt nleaves,const PetscInt *ilocal,PetscCopyMode modelocal,const PetscSFNode *remote,PetscCopyMode moderemote);
62*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscSFGetGraph(PetscSF,PetscInt *nroots,PetscInt *nleaves,const PetscInt **ilocal,const PetscSFNode **iremote);
63*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscSFCreateEmbeddedSF(PetscSF,PetscInt nroots,const PetscInt *selected,PetscSF *newsf);
64*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscSFCreateArray(PetscSF,MPI_Datatype,void*,void*);
65*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscSFDestroyArray(PetscSF,MPI_Datatype,void*,void*);
66*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscSFReset(PetscSF);
67*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscSFGetRanks(PetscSF,PetscInt*,const PetscMPIInt**,const PetscInt**,const PetscMPIInt**,const PetscMPIInt**);
68*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscSFGetDataTypes(PetscSF,MPI_Datatype,const MPI_Datatype**,const MPI_Datatype**);
69*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscSFGetWindow(PetscSF,MPI_Datatype,void*,PetscBool,PetscMPIInt,PetscMPIInt,PetscMPIInt,MPI_Win*);
70*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscSFFindWindow(PetscSF,MPI_Datatype,const void*,MPI_Win*);
71*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscSFRestoreWindow(PetscSF,MPI_Datatype,const void*,PetscBool,PetscMPIInt,MPI_Win*);
72*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscSFGetGroups(PetscSF,MPI_Group*,MPI_Group*);
73*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscSFGetMultiSF(PetscSF,PetscSF*);
74*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscSFCreateInverseSF(PetscSF,PetscSF*);
75936c5a86SJed Brown 
76936c5a86SJed Brown /* broadcasts rootdata to leafdata */
77*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscSFBcastBegin(PetscSF,MPI_Datatype,const void *rootdata,void *leafdata);
78*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscSFBcastEnd(PetscSF,MPI_Datatype,const void *rootdata,void *leafdata);
79936c5a86SJed Brown /* Reduce leafdata into rootdata using provided operation */
80*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscSFReduceBegin(PetscSF,MPI_Datatype,const void *leafdata,void *rootdata,MPI_Op);
81*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscSFReduceEnd(PetscSF,MPI_Datatype,const void *leafdata,void *rootdata,MPI_Op);
82936c5a86SJed Brown /* Atomically modifies (using provided operation) rootdata using leafdata from each leaf, value at root at time of modification is returned in leafupdate. */
83*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscSFFetchAndOpBegin(PetscSF,MPI_Datatype,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op);
84*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscSFFetchAndOpEnd(PetscSF,MPI_Datatype,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op);
85936c5a86SJed Brown /* Compute the degree of every root vertex (number of leaves in its star) */
86*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscSFComputeDegreeBegin(PetscSF,const PetscInt **degree);
87*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscSFComputeDegreeEnd(PetscSF,const PetscInt **degree);
88936c5a86SJed Brown /* Concatenate data from all leaves into roots */
89*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscSFGatherBegin(PetscSF,MPI_Datatype,const void *leafdata,void *multirootdata);
90*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscSFGatherEnd(PetscSF,MPI_Datatype,const void *leafdata,void *multirootdata);
91936c5a86SJed Brown /* Distribute distinct values to each leaf from roots */
92*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscSFScatterBegin(PetscSF,MPI_Datatype,const void *multirootdata,void *leafdata);
93*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscSFScatterEnd(PetscSF,MPI_Datatype,const void *multirootdata,void *leafdata);
94936c5a86SJed Brown 
95936c5a86SJed Brown 
96936c5a86SJed Brown #endif
97