Actual source code: petscis.h

petsc-3.12.5 2020-03-29
Report Typos and Errors
  1: /*
  2:    An index set is a generalization of a subset of integers.  Index sets
  3:    are used for defining scatters and gathers.
  4: */
  5: #if !defined(PETSCIS_H)
  6: #define PETSCIS_H
  7:  #include <petscsys.h>
  8:  #include <petscsftypes.h>
  9:  #include <petscsectiontypes.h>
 10:  #include <petscistypes.h>

 12: #define IS_FILE_CLASSID 1211218
 13: PETSC_EXTERN PetscClassId IS_CLASSID;

 15: PETSC_EXTERN PetscErrorCode ISInitializePackage(void);

 17: /*J
 18:     ISType - String with the name of a PETSc index set type

 20:    Level: beginner

 22: .seealso: ISSetType(), IS, ISCreate(), ISRegister()
 23: J*/
 24: typedef const char* ISType;
 25: #define ISGENERAL      "general"
 26: #define ISSTRIDE       "stride"
 27: #define ISBLOCK        "block"

 29: /* Dynamic creation and loading functions */
 30: PETSC_EXTERN PetscFunctionList ISList;
 31: PETSC_EXTERN PetscErrorCode ISSetType(IS, ISType);
 32: PETSC_EXTERN PetscErrorCode ISGetType(IS, ISType *);
 33: PETSC_EXTERN PetscErrorCode ISRegister(const char[],PetscErrorCode (*)(IS));
 34: PETSC_EXTERN PetscErrorCode ISCreate(MPI_Comm,IS*);

 36: /*
 37:     Default index set data structures that PETSc provides.
 38: */
 39: PETSC_EXTERN PetscErrorCode ISCreateGeneral(MPI_Comm,PetscInt,const PetscInt[],PetscCopyMode,IS *);
 40: PETSC_EXTERN PetscErrorCode ISGeneralSetIndices(IS,PetscInt,const PetscInt[],PetscCopyMode);
 41: PETSC_EXTERN PetscErrorCode ISCreateBlock(MPI_Comm,PetscInt,PetscInt,const PetscInt[],PetscCopyMode,IS *);
 42: PETSC_EXTERN PetscErrorCode ISBlockSetIndices(IS,PetscInt,PetscInt,const PetscInt[],PetscCopyMode);
 43: PETSC_EXTERN PetscErrorCode ISCreateStride(MPI_Comm,PetscInt,PetscInt,PetscInt,IS *);
 44: PETSC_EXTERN PetscErrorCode ISStrideSetStride(IS,PetscInt,PetscInt,PetscInt);

 46: PETSC_EXTERN PetscErrorCode ISDestroy(IS*);
 47: PETSC_EXTERN PetscErrorCode ISSetPermutation(IS);
 48: PETSC_EXTERN PetscErrorCode ISPermutation(IS,PetscBool *);
 49: PETSC_EXTERN PetscErrorCode ISSetIdentity(IS);
 50: PETSC_EXTERN PetscErrorCode ISIdentity(IS,PetscBool *);
 51: PETSC_EXTERN PetscErrorCode ISContiguousLocal(IS,PetscInt,PetscInt,PetscInt*,PetscBool*);

 53: PETSC_EXTERN PetscErrorCode ISGetIndices(IS,const PetscInt *[]);
 54: PETSC_EXTERN PetscErrorCode ISRestoreIndices(IS,const PetscInt *[]);
 55: PETSC_EXTERN PetscErrorCode ISGetTotalIndices(IS,const PetscInt *[]);
 56: PETSC_EXTERN PetscErrorCode ISRestoreTotalIndices(IS,const PetscInt *[]);
 57: PETSC_EXTERN PetscErrorCode ISGetNonlocalIndices(IS,const PetscInt *[]);
 58: PETSC_EXTERN PetscErrorCode ISRestoreNonlocalIndices(IS,const PetscInt *[]);
 59: PETSC_EXTERN PetscErrorCode ISGetNonlocalIS(IS, IS *is);
 60: PETSC_EXTERN PetscErrorCode ISRestoreNonlocalIS(IS, IS *is);
 61: PETSC_EXTERN PetscErrorCode ISGetSize(IS,PetscInt *);
 62: PETSC_EXTERN PetscErrorCode ISGetLocalSize(IS,PetscInt *);
 63: PETSC_EXTERN PetscErrorCode ISInvertPermutation(IS,PetscInt,IS*);
 64: PETSC_EXTERN PetscErrorCode ISView(IS,PetscViewer);
 65: PETSC_STATIC_INLINE PetscErrorCode ISViewFromOptions(IS A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);}
 66: PETSC_EXTERN PetscErrorCode ISLoad(IS,PetscViewer);
 67: PETSC_EXTERN PetscErrorCode ISEqual(IS,IS,PetscBool  *);
 68: PETSC_EXTERN PetscErrorCode ISEqualUnsorted(IS,IS,PetscBool *);
 69: PETSC_EXTERN PetscErrorCode ISSort(IS);
 70: PETSC_EXTERN PetscErrorCode ISSortRemoveDups(IS);
 71: PETSC_EXTERN PetscErrorCode ISSorted(IS,PetscBool  *);
 72: PETSC_EXTERN PetscErrorCode ISDifference(IS,IS,IS*);
 73: PETSC_EXTERN PetscErrorCode ISSum(IS,IS,IS*);
 74: PETSC_EXTERN PetscErrorCode ISExpand(IS,IS,IS*);
 75: PETSC_EXTERN PetscErrorCode ISIntersect(IS,IS,IS*);
 76: PETSC_EXTERN PetscErrorCode ISGetMinMax(IS,PetscInt*,PetscInt*);

 78: PETSC_EXTERN PetscErrorCode ISLocate(IS,PetscInt,PetscInt*);
 79: PETSC_EXTERN PetscErrorCode ISGetPointRange(IS,PetscInt*,PetscInt*,const PetscInt**);
 80: PETSC_EXTERN PetscErrorCode ISRestorePointRange(IS,PetscInt*,PetscInt*,const PetscInt**);
 81: PETSC_EXTERN PetscErrorCode ISGetPointSubrange(IS,PetscInt,PetscInt,const PetscInt*);

 83: PETSC_EXTERN PetscErrorCode ISBlockGetIndices(IS,const PetscInt *[]);
 84: PETSC_EXTERN PetscErrorCode ISBlockRestoreIndices(IS,const PetscInt *[]);
 85: PETSC_EXTERN PetscErrorCode ISBlockGetLocalSize(IS,PetscInt *);
 86: PETSC_EXTERN PetscErrorCode ISBlockGetSize(IS,PetscInt *);
 87: PETSC_EXTERN PetscErrorCode ISGetBlockSize(IS,PetscInt*);
 88: PETSC_EXTERN PetscErrorCode ISSetBlockSize(IS,PetscInt);

 90: PETSC_EXTERN PetscErrorCode ISStrideGetInfo(IS,PetscInt *,PetscInt*);

 92: PETSC_EXTERN PetscErrorCode ISToGeneral(IS);

 94: PETSC_EXTERN PetscErrorCode ISDuplicate(IS,IS*);
 95: PETSC_EXTERN PetscErrorCode ISCopy(IS,IS);
 96: PETSC_EXTERN PetscErrorCode ISAllGather(IS,IS*);
 97: PETSC_EXTERN PetscErrorCode ISComplement(IS,PetscInt,PetscInt,IS*);
 98: PETSC_EXTERN PetscErrorCode ISConcatenate(MPI_Comm,PetscInt,const IS[],IS*);
 99: PETSC_EXTERN PetscErrorCode ISListToPair(MPI_Comm,PetscInt, IS[],IS*,IS*);
100: PETSC_EXTERN PetscErrorCode ISPairToList(IS,IS,PetscInt*, IS *[]);
101: PETSC_EXTERN PetscErrorCode ISEmbed(IS,IS,PetscBool,IS*);
102: PETSC_EXTERN PetscErrorCode ISSortPermutation(IS,PetscBool,IS*);
103: PETSC_EXTERN PetscErrorCode ISOnComm(IS,MPI_Comm,PetscCopyMode,IS*);
104: PETSC_EXTERN PetscErrorCode ISRenumber(IS,IS,PetscInt*,IS*);
105: PETSC_EXTERN PetscErrorCode ISCreateSubIS(IS,IS,IS*);

107: PETSC_EXTERN PetscErrorCode ISGeneralFilter(IS,PetscInt,PetscInt);

109: /* --------------------------------------------------------------------------*/
110: PETSC_EXTERN PetscClassId IS_LTOGM_CLASSID;

112: /*E
113:     ISGlobalToLocalMappingMode - Indicates mapping behavior if global indices are missing

115:    IS_GTOLM_MASK - missing global indices are masked by mapping them to a local index of -1
116:    IS_GTOLM_DROP - missing global indices are dropped

118:    Level: beginner

120: .seealso: ISGlobalToLocalMappingApplyBlock(), ISGlobalToLocalMappingApply()

122: E*/
123: typedef enum {IS_GTOLM_MASK,IS_GTOLM_DROP} ISGlobalToLocalMappingMode;

125: /*J
126:     ISLocalToGlobalMappingType - String with the name of a mapping method

128:    Level: beginner

130: .seealso: ISLocalToGlobalMappingSetType(), ISLocalToGlobalSetFromOptions()
131: J*/
132: typedef const char* ISLocalToGlobalMappingType;
133: #define ISLOCALTOGLOBALMAPPINGBASIC "basic"
134: #define ISLOCALTOGLOBALMAPPINGHASH  "hash"

136: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingSetType(ISLocalToGlobalMapping,ISLocalToGlobalMappingType);
137: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingRegister(const char[],PetscErrorCode (*)(ISLocalToGlobalMapping));
138: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingRegisterAll(void);
139: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingCreate(MPI_Comm,PetscInt,PetscInt,const PetscInt[],PetscCopyMode,ISLocalToGlobalMapping*);
140: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingCreateIS(IS,ISLocalToGlobalMapping *);
141: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingCreateSF(PetscSF,PetscInt,ISLocalToGlobalMapping*);
142: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingSetFromOptions(ISLocalToGlobalMapping);
143: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingSetUp(ISLocalToGlobalMapping);
144: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingView(ISLocalToGlobalMapping,PetscViewer);
145: PETSC_STATIC_INLINE PetscErrorCode ISLocalToGlobalMappingViewFromOptions(ISLocalToGlobalMapping A,PetscObject B,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,B,name);}

147: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingDestroy(ISLocalToGlobalMapping*);
148: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingApply(ISLocalToGlobalMapping,PetscInt,const PetscInt[],PetscInt[]);
149: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingApplyBlock(ISLocalToGlobalMapping,PetscInt,const PetscInt[],PetscInt[]);
150: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingApplyIS(ISLocalToGlobalMapping,IS,IS*);
151: PETSC_EXTERN PetscErrorCode ISGlobalToLocalMappingApply(ISLocalToGlobalMapping,ISGlobalToLocalMappingMode,PetscInt,const PetscInt[],PetscInt*,PetscInt[]);
152: PETSC_EXTERN PetscErrorCode ISGlobalToLocalMappingApplyBlock(ISLocalToGlobalMapping,ISGlobalToLocalMappingMode,PetscInt,const PetscInt[],PetscInt*,PetscInt[]);
153: PETSC_EXTERN PetscErrorCode ISGlobalToLocalMappingApplyIS(ISLocalToGlobalMapping,ISGlobalToLocalMappingMode,IS,IS*);
154: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingGetSize(ISLocalToGlobalMapping,PetscInt*);
155: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingGetNodeInfo(ISLocalToGlobalMapping,PetscInt*,PetscInt*[],PetscInt**[]);
156: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingRestoreNodeInfo(ISLocalToGlobalMapping,PetscInt*,PetscInt*[],PetscInt**[]);
157: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping,PetscInt*,PetscInt*[],PetscInt*[],PetscInt**[]);
158: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingRestoreInfo(ISLocalToGlobalMapping,PetscInt*,PetscInt*[],PetscInt*[],PetscInt**[]);
159: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingGetBlockInfo(ISLocalToGlobalMapping,PetscInt*,PetscInt*[],PetscInt*[],PetscInt**[]);
160: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingRestoreBlockInfo(ISLocalToGlobalMapping,PetscInt*,PetscInt*[],PetscInt*[],PetscInt**[]);
161: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingGetIndices(ISLocalToGlobalMapping,const PetscInt**);
162: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingRestoreIndices(ISLocalToGlobalMapping,const PetscInt**);
163: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingGetBlockIndices(ISLocalToGlobalMapping,const PetscInt**);
164: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingRestoreBlockIndices(ISLocalToGlobalMapping,const PetscInt**);
165: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingConcatenate(MPI_Comm,PetscInt,const ISLocalToGlobalMapping[],ISLocalToGlobalMapping*);
166: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingGetBlockSize(ISLocalToGlobalMapping,PetscInt*);
167: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingSetBlockSize(ISLocalToGlobalMapping,PetscInt);
168: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingDuplicate(ISLocalToGlobalMapping,ISLocalToGlobalMapping*);


171: /* --------------------------------------------------------------------------*/
172: /*E
173:     ISColoringType - determines if the coloring is for the entire parallel grid/graph/matrix
174:                      or for just the local ghosted portion

176:     Level: beginner

178: $   IS_COLORING_GLOBAL - does not include the colors for ghost points, this is used when the function
179: $                        is called synchronously in parallel. This requires generating a "parallel coloring".
180: $   IS_COLORING_LOCAL - includes colors for ghost points, this is used when the function can be called
181: $                         separately on individual processes with the ghost points already filled in. Does not
182: $                         require a "parallel coloring", rather each process colors its local + ghost part.
183: $                         Using this can result in much less parallel communication. Currently only works
184: $                         with DMDA and if you call MatFDColoringSetFunction() with the local function.

186: .seealso: DMCreateColoring()
187: E*/
188: typedef enum {IS_COLORING_GLOBAL,IS_COLORING_LOCAL} ISColoringType;
189: PETSC_EXTERN const char *const ISColoringTypes[];
190: typedef unsigned PETSC_IS_COLOR_VALUE_TYPE ISColoringValue;
191: PETSC_EXTERN PetscErrorCode ISAllGatherColors(MPI_Comm,PetscInt,ISColoringValue*,PetscInt*,ISColoringValue*[]);

193: PETSC_EXTERN PetscErrorCode ISColoringCreate(MPI_Comm,PetscInt,PetscInt,const ISColoringValue[],PetscCopyMode,ISColoring*);
194: PETSC_EXTERN PetscErrorCode ISColoringDestroy(ISColoring*);
195: PETSC_EXTERN PetscErrorCode ISColoringView(ISColoring,PetscViewer);
196: PETSC_EXTERN PetscErrorCode ISColoringViewFromOptions(ISColoring,PetscObject,const char[]);
197: PETSC_EXTERN PetscErrorCode ISColoringGetIS(ISColoring,PetscCopyMode,PetscInt*,IS*[]);
198: PETSC_EXTERN PetscErrorCode ISColoringRestoreIS(ISColoring,PetscCopyMode,IS*[]);
199: PETSC_EXTERN PetscErrorCode ISColoringReference(ISColoring);
200: PETSC_EXTERN PetscErrorCode ISColoringSetType(ISColoring,ISColoringType);
201: PETSC_EXTERN PetscErrorCode ISColoringGetType(ISColoring,ISColoringType*);
202: PETSC_EXTERN PetscErrorCode ISColoringGetColors(ISColoring,PetscInt*,PetscInt*,const ISColoringValue**);

204: /* --------------------------------------------------------------------------*/
205: PETSC_EXTERN PetscErrorCode ISBuildTwoSided(IS,IS,IS*);
206: PETSC_EXTERN PetscErrorCode ISPartitioningToNumbering(IS,IS*);
207: PETSC_EXTERN PetscErrorCode ISPartitioningCount(IS,PetscInt,PetscInt[]);

209: PETSC_EXTERN PetscErrorCode ISCompressIndicesGeneral(PetscInt,PetscInt,PetscInt,PetscInt,const IS[],IS[]);
210: PETSC_EXTERN PetscErrorCode ISCompressIndicesSorted(PetscInt,PetscInt,PetscInt,const IS[],IS[]);
211: PETSC_EXTERN PetscErrorCode ISExpandIndicesGeneral(PetscInt,PetscInt,PetscInt,PetscInt,const IS[],IS[]);


214: struct _n_PetscLayout{
215:   MPI_Comm               comm;
216:   PetscInt               n,N;         /* local, global vector size */
217:   PetscInt               rstart,rend; /* local start, local end + 1 */
218:   PetscInt               *range;      /* the offset of each processor */
219:   PetscBool              range_alloc; /* should range be freed in Destroy? */
220:   PetscInt               bs;          /* number of elements in each block (generally for multi-component
221:                                        * problems). Defaults to -1 and can be arbitrarily lazy so always use
222:                                        * PetscAbs(map->bs) when accessing directly and expecting result to be
223:                                        * positive. Do NOT multiply above numbers by bs */
224:   PetscInt               refcnt;      /* MPI Vecs obtained with VecDuplicate() and from MatCreateVecs() reuse map of input object */
225:   ISLocalToGlobalMapping mapping;     /* mapping used in Vec/MatSetValuesLocal() */
226:   PetscBool              setupcalled; /* Forbid setup more than once */
227:   PetscInt               oldn,oldN;   /* Checking if setup is allowed */
228:   PetscInt               oldbs;       /* And again */
229: };

231: /*@C
232:      PetscLayoutFindOwner - Find the owning rank for a global index

234:     Not Collective

236:    Input Parameters:
237: +    map - the layout
238: -    idx - global index to find the owner of

240:    Output Parameter:
241: .    owner - the owning rank

243:    Level: developer

245:     Fortran Notes:
246:       Not available from Fortran

248: @*/
249: PETSC_STATIC_INLINE PetscErrorCode PetscLayoutFindOwner(PetscLayout map,PetscInt idx,PetscInt *owner)
250: {
252:   PetscMPIInt    lo = 0,hi,t;

255:   *owner = -1;                  /* GCC erroneously issues warning about possibly uninitialized use when error condition */
256:   if (!((map->n >= 0) && (map->N >= 0) && (map->range))) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"PetscLayoutSetUp() must be called first");
257:   if (idx < 0 || idx > map->N) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Index %D is out of range",idx);
258:   MPI_Comm_size(map->comm,&hi);
259:   while (hi - lo > 1) {
260:     t = lo + (hi - lo) / 2;
261:     if (idx < map->range[t]) hi = t;
262:     else                     lo = t;
263:   }
264:   *owner = lo;
265:   return(0);
266: }

268: /*@C
269:      PetscLayoutFindOwnerIndex - Find the owning rank and the local index for a global index

271:     Not Collective

273:    Input Parameters:
274: +    map   - the layout
275: -    idx   - global index to find the owner of

277:    Output Parameter:
278: +    owner - the owning rank
279: -    lidx  - local index used by the owner for idx

281:    Level: developer

283:     Fortran Notes:
284:       Not available from Fortran

286: @*/
287: PETSC_STATIC_INLINE PetscErrorCode PetscLayoutFindOwnerIndex(PetscLayout map,PetscInt idx,PetscInt *owner, PetscInt *lidx)
288: {
290:   PetscMPIInt    lo = 0,hi,t;

293:   if (!((map->n >= 0) && (map->N >= 0) && (map->range))) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"PetscLayoutSetUp() must be called first");
294:   if (idx < 0 || idx > map->N) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Index %D is out of range",idx);
295:   MPI_Comm_size(map->comm,&hi);
296:   while (hi - lo > 1) {
297:     t = lo + (hi - lo) / 2;
298:     if (idx < map->range[t]) hi = t;
299:     else                     lo = t;
300:   }
301:   if (owner) *owner = lo;
302:   if (lidx) *lidx  = idx-map->range[lo];
303:   return(0);
304: }

306: PETSC_EXTERN PetscErrorCode PetscLayoutCreate(MPI_Comm,PetscLayout*);
307: PETSC_EXTERN PetscErrorCode PetscLayoutCreateFromSizes(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscLayout*);
308: PETSC_EXTERN PetscErrorCode PetscLayoutCreateFromRanges(MPI_Comm,const PetscInt[],PetscCopyMode,PetscInt,PetscLayout*);
309: PETSC_EXTERN PetscErrorCode PetscLayoutSetUp(PetscLayout);
310: PETSC_EXTERN PetscErrorCode PetscLayoutDestroy(PetscLayout*);
311: PETSC_EXTERN PetscErrorCode PetscLayoutDuplicate(PetscLayout,PetscLayout*);
312: PETSC_EXTERN PetscErrorCode PetscLayoutReference(PetscLayout,PetscLayout*);
313: PETSC_EXTERN PetscErrorCode PetscLayoutSetLocalSize(PetscLayout,PetscInt);
314: PETSC_EXTERN PetscErrorCode PetscLayoutGetLocalSize(PetscLayout,PetscInt *);
315: PETSC_EXTERN PetscErrorCode PetscLayoutSetSize(PetscLayout,PetscInt);
316: PETSC_EXTERN PetscErrorCode PetscLayoutGetSize(PetscLayout,PetscInt *);
317: PETSC_EXTERN PetscErrorCode PetscLayoutSetBlockSize(PetscLayout,PetscInt);
318: PETSC_EXTERN PetscErrorCode PetscLayoutGetBlockSize(PetscLayout,PetscInt*);
319: PETSC_EXTERN PetscErrorCode PetscLayoutGetRange(PetscLayout,PetscInt *,PetscInt *);
320: PETSC_EXTERN PetscErrorCode PetscLayoutGetRanges(PetscLayout,const PetscInt *[]);
321: PETSC_EXTERN PetscErrorCode PetscLayoutCompare(PetscLayout,PetscLayout,PetscBool*);
322: PETSC_EXTERN PetscErrorCode PetscLayoutSetISLocalToGlobalMapping(PetscLayout,ISLocalToGlobalMapping);
323: PETSC_EXTERN PetscErrorCode PetscLayoutMapLocal(PetscLayout,PetscInt,const PetscInt[],PetscInt*,PetscInt**,PetscInt**);
324: PETSC_EXTERN PetscErrorCode PetscSFSetGraphLayout(PetscSF,PetscLayout,PetscInt,const PetscInt*,PetscCopyMode,const PetscInt*);

326: /* PetscSF support */
327: PETSC_EXTERN PetscErrorCode PetscSFConvertPartition(PetscSF, PetscSection, IS, ISLocalToGlobalMapping *, PetscSF *);
328: PETSC_EXTERN PetscErrorCode PetscSFCreateRemoteOffsets(PetscSF, PetscSection, PetscSection, PetscInt **);
329: PETSC_EXTERN PetscErrorCode PetscSFDistributeSection(PetscSF, PetscSection, PetscInt **, PetscSection);
330: PETSC_EXTERN PetscErrorCode PetscSFCreateSectionSF(PetscSF, PetscSection, PetscInt [], PetscSection, PetscSF *);

332: #endif