Actual source code: petscis.h

petsc-3.4.5 2014-06-29
  1: /*
  2:    An index set is a generalization of a subset of integers.  Index sets
  3:    are used for defining scatters and gathers.
  4: */
  7: #include <petscsys.h>
  8: #include <petscviewertypes.h>
  9: #include <petscsftypes.h>

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

 14: PETSC_EXTERN PetscErrorCode ISInitializePackage(void);

 16: /*S
 17:      IS - Abstract PETSc object that allows indexing.

 19:    Level: beginner

 21:   Concepts: indexing, stride

 23: .seealso:  ISCreateGeneral(), ISCreateBlock(), ISCreateStride(), ISGetIndices(), ISDestroy()
 24: S*/
 25: typedef struct _p_IS* IS;

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

 30:    Level: beginner

 32: .seealso: ISSetType(), IS, ISCreate(), ISRegister()
 33: J*/
 34: typedef const char* ISType;
 35: #define ISGENERAL      "general"
 36: #define ISSTRIDE       "stride"
 37: #define ISBLOCK        "block"

 39: /* Dynamic creation and loading functions */
 40: PETSC_EXTERN PetscFunctionList ISList;
 41: PETSC_EXTERN PetscBool         ISRegisterAllCalled;
 42: PETSC_EXTERN PetscErrorCode ISSetType(IS, ISType);
 43: PETSC_EXTERN PetscErrorCode ISGetType(IS, ISType *);
 44: PETSC_EXTERN PetscErrorCode ISRegister(const char[],PetscErrorCode (*)(IS));
 45: PETSC_EXTERN PetscErrorCode ISRegisterAll(void);
 46: PETSC_EXTERN PetscErrorCode ISCreate(MPI_Comm,IS*);

 48: /*
 49:     Default index set data structures that PETSc provides.
 50: */
 51: PETSC_EXTERN PetscErrorCode ISCreateGeneral(MPI_Comm,PetscInt,const PetscInt[],PetscCopyMode,IS *);
 52: PETSC_EXTERN PetscErrorCode ISGeneralSetIndices(IS,PetscInt,const PetscInt[],PetscCopyMode);
 53: PETSC_EXTERN PetscErrorCode ISCreateBlock(MPI_Comm,PetscInt,PetscInt,const PetscInt[],PetscCopyMode,IS *);
 54: PETSC_EXTERN PetscErrorCode ISBlockSetIndices(IS,PetscInt,PetscInt,const PetscInt[],PetscCopyMode);
 55: PETSC_EXTERN PetscErrorCode ISCreateStride(MPI_Comm,PetscInt,PetscInt,PetscInt,IS *);
 56: PETSC_EXTERN PetscErrorCode ISStrideSetStride(IS,PetscInt,PetscInt,PetscInt);

 58: PETSC_EXTERN PetscErrorCode ISDestroy(IS*);
 59: PETSC_EXTERN PetscErrorCode ISSetPermutation(IS);
 60: PETSC_EXTERN PetscErrorCode ISPermutation(IS,PetscBool *);
 61: PETSC_EXTERN PetscErrorCode ISSetIdentity(IS);
 62: PETSC_EXTERN PetscErrorCode ISIdentity(IS,PetscBool *);
 63: PETSC_EXTERN PetscErrorCode ISContiguousLocal(IS,PetscInt,PetscInt,PetscInt*,PetscBool*);

 65: PETSC_EXTERN PetscErrorCode ISGetIndices(IS,const PetscInt *[]);
 66: PETSC_EXTERN PetscErrorCode ISRestoreIndices(IS,const PetscInt *[]);
 67: PETSC_EXTERN PetscErrorCode ISGetTotalIndices(IS,const PetscInt *[]);
 68: PETSC_EXTERN PetscErrorCode ISRestoreTotalIndices(IS,const PetscInt *[]);
 69: PETSC_EXTERN PetscErrorCode ISGetNonlocalIndices(IS,const PetscInt *[]);
 70: PETSC_EXTERN PetscErrorCode ISRestoreNonlocalIndices(IS,const PetscInt *[]);
 71: PETSC_EXTERN PetscErrorCode ISGetNonlocalIS(IS, IS *is);
 72: PETSC_EXTERN PetscErrorCode ISRestoreNonlocalIS(IS, IS *is);
 73: PETSC_EXTERN PetscErrorCode ISGetSize(IS,PetscInt *);
 74: PETSC_EXTERN PetscErrorCode ISGetLocalSize(IS,PetscInt *);
 75: PETSC_EXTERN PetscErrorCode ISInvertPermutation(IS,PetscInt,IS*);
 76: PETSC_EXTERN PetscErrorCode ISView(IS,PetscViewer);
 77: PETSC_EXTERN PetscErrorCode ISEqual(IS,IS,PetscBool  *);
 78: PETSC_EXTERN PetscErrorCode ISSort(IS);
 79: PETSC_EXTERN PetscErrorCode ISSorted(IS,PetscBool  *);
 80: PETSC_EXTERN PetscErrorCode ISDifference(IS,IS,IS*);
 81: PETSC_EXTERN PetscErrorCode ISSum(IS,IS,IS*);
 82: PETSC_EXTERN PetscErrorCode ISExpand(IS,IS,IS*);
 83: PETSC_EXTERN PetscErrorCode ISGetMinMax(IS,PetscInt*,PetscInt*);

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

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

 94: PETSC_EXTERN PetscErrorCode ISToGeneral(IS);

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

106: /* --------------------------------------------------------------------------*/
107: PETSC_EXTERN PetscClassId IS_LTOGM_CLASSID;

109: /*S
110:    ISLocalToGlobalMapping - mappings from an arbitrary
111:       local ordering from 0 to n-1 to a global PETSc ordering
112:       used by a vector or matrix.

114:    Level: intermediate

116:    Note: mapping from Local to Global is scalable; but Global
117:   to Local may not be if the range of global values represented locally
118:   is very large.

120:    Note: the ISLocalToGlobalMapping is actually a private object; it is included
121:   here for the inline function ISLocalToGlobalMappingApply() to allow it to be inlined since
122:   it is used so often.

124: .seealso:  ISLocalToGlobalMappingCreate()
125: S*/
126: typedef struct _p_ISLocalToGlobalMapping* ISLocalToGlobalMapping;

128: /*E
129:     ISGlobalToLocalMappingType - Indicates if missing global indices are

131:    IS_GTOLM_MASK - missing global indices are replaced with -1
132:    IS_GTOLM_DROP - missing global indices are dropped

134:    Level: beginner

136: .seealso: ISGlobalToLocalMappingApply()

138: E*/
139: typedef enum {IS_GTOLM_MASK,IS_GTOLM_DROP} ISGlobalToLocalMappingType;

141: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingCreate(MPI_Comm,PetscInt,const PetscInt[],PetscCopyMode,ISLocalToGlobalMapping*);
142: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingCreateIS(IS,ISLocalToGlobalMapping *);
143: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingCreateSF(PetscSF,PetscInt,ISLocalToGlobalMapping*);
144: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingView(ISLocalToGlobalMapping,PetscViewer);
145: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingDestroy(ISLocalToGlobalMapping*);
146: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingApplyIS(ISLocalToGlobalMapping,IS,IS*);
147: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingApply(ISLocalToGlobalMapping,PetscInt,const PetscInt[],PetscInt[]);
148: PETSC_EXTERN PetscErrorCode ISGlobalToLocalMappingApply(ISLocalToGlobalMapping,ISGlobalToLocalMappingType,PetscInt,const PetscInt[],PetscInt*,PetscInt[]);
149: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingGetSize(ISLocalToGlobalMapping,PetscInt*);
150: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping,PetscInt*,PetscInt*[],PetscInt*[],PetscInt**[]);
151: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingRestoreInfo(ISLocalToGlobalMapping,PetscInt*,PetscInt*[],PetscInt*[],PetscInt**[]);
152: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingGetIndices(ISLocalToGlobalMapping,const PetscInt**);
153: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingRestoreIndices(ISLocalToGlobalMapping,const PetscInt**);
154: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingBlock(ISLocalToGlobalMapping,PetscInt,ISLocalToGlobalMapping*);
155: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingUnBlock(ISLocalToGlobalMapping,PetscInt,ISLocalToGlobalMapping*);
156: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingConcatenate(MPI_Comm,PetscInt,const ISLocalToGlobalMapping[],ISLocalToGlobalMapping*);
157: PETSC_EXTERN PetscErrorCode ISG2LMapApply(ISLocalToGlobalMapping,PetscInt,const PetscInt[],PetscInt[]);

159: /* --------------------------------------------------------------------------*/
160: /*E
161:     ISColoringType - determines if the coloring is for the entire parallel grid/graph/matrix
162:                      or for just the local ghosted portion

164:     Level: beginner

166: $   IS_COLORING_GLOBAL - does not include the colors for ghost points, this is used when the function
167: $                        is called synchronously in parallel. This requires generating a "parallel coloring".
168: $   IS_COLORING_GHOSTED - includes colors for ghost points, this is used when the function can be called
169: $                         seperately on individual processes with the ghost points already filled in. Does not
170: $                         require a "parallel coloring", rather each process colors its local + ghost part.
171: $                         Using this can result in much less parallel communication. In the paradigm of
172: $                         DMGetLocalVector() and DMGetGlobalVector() this could be called IS_COLORING_LOCAL

174: .seealso: DMCreateColoring()
175: E*/
176: typedef enum {IS_COLORING_GLOBAL,IS_COLORING_GHOSTED} ISColoringType;
177: PETSC_EXTERN const char *const ISColoringTypes[];
178: typedef unsigned PETSC_IS_COLOR_VALUE_TYPE ISColoringValue;
179: PETSC_EXTERN PetscErrorCode ISAllGatherColors(MPI_Comm,PetscInt,ISColoringValue*,PetscInt*,ISColoringValue*[]);

181: /*S
182:      ISColoring - sets of IS's that define a coloring
183:               of the underlying indices

185:    Level: intermediate

187:     Notes:
188:         One should not access the *is records below directly because they may not yet
189:     have been created. One should use ISColoringGetIS() to make sure they are
190:     created when needed.

192:     Developer Note: this is not a PetscObject

194: .seealso:  ISColoringCreate(), ISColoringGetIS(), ISColoringView(), ISColoringGetIS()
195: S*/
196: struct _n_ISColoring {
197:   PetscInt        refct;
198:   PetscInt        n;                /* number of colors */
199:   IS              *is;              /* for each color indicates columns */
200:   MPI_Comm        comm;
201:   ISColoringValue *colors;          /* for each column indicates color */
202:   PetscInt        N;                /* number of columns */
203:   ISColoringType  ctype;
204: };
205: typedef struct _n_ISColoring* ISColoring;

207: PETSC_EXTERN PetscErrorCode ISColoringCreate(MPI_Comm,PetscInt,PetscInt,const ISColoringValue[],ISColoring*);
208: PETSC_EXTERN PetscErrorCode ISColoringDestroy(ISColoring*);
209: PETSC_EXTERN PetscErrorCode ISColoringView(ISColoring,PetscViewer);
210: PETSC_EXTERN PetscErrorCode ISColoringGetIS(ISColoring,PetscInt*,IS*[]);
211: PETSC_EXTERN PetscErrorCode ISColoringRestoreIS(ISColoring,IS*[]);
212: PETSC_EXTERN PetscErrorCode ISColoringReference(ISColoring);
213: PETSC_EXTERN PetscErrorCode ISColoringSetType(ISColoring,ISColoringType);


216: /* --------------------------------------------------------------------------*/

218: PETSC_EXTERN PetscErrorCode ISPartitioningToNumbering(IS,IS*);
219: PETSC_EXTERN PetscErrorCode ISPartitioningCount(IS,PetscInt,PetscInt[]);

221: PETSC_EXTERN PetscErrorCode ISCompressIndicesGeneral(PetscInt,PetscInt,PetscInt,PetscInt,const IS[],IS[]);
222: PETSC_EXTERN PetscErrorCode ISCompressIndicesSorted(PetscInt,PetscInt,PetscInt,const IS[],IS[]);
223: PETSC_EXTERN PetscErrorCode ISExpandIndicesGeneral(PetscInt,PetscInt,PetscInt,PetscInt,const IS[],IS[]);

225: /*S
226:      PetscLayout - defines layout of vectors and matrices across processes (which rows are owned by which processes)

228:    Level: developer


231: .seealso:  PetscLayoutCreate(), PetscLayoutDestroy()
232: S*/
233: typedef struct _n_PetscLayout* PetscLayout;
234: struct _n_PetscLayout{
235:   MPI_Comm               comm;
236:   PetscInt               n,N;         /* local, global vector size */
237:   PetscInt               rstart,rend; /* local start, local end + 1 */
238:   PetscInt               *range;      /* the offset of each processor */
239:   PetscInt               bs;          /* number of elements in each block (generally for multi-component problems) Do NOT multiply above numbers by bs */
240:   PetscInt               refcnt;      /* MPI Vecs obtained with VecDuplicate() and from MatGetVecs() reuse map of input object */
241:   ISLocalToGlobalMapping mapping;     /* mapping used in Vec/MatSetValuesLocal() */
242:   ISLocalToGlobalMapping bmapping;    /* mapping used in Vec/MatSetValuesBlockedLocal() */
243:   PetscInt               *trstarts;   /* local start for each thread */
244: };

246: PETSC_EXTERN PetscErrorCode PetscLayoutCreate(MPI_Comm,PetscLayout*);
247: PETSC_EXTERN PetscErrorCode PetscLayoutSetUp(PetscLayout);
248: PETSC_EXTERN PetscErrorCode PetscLayoutDestroy(PetscLayout*);
249: PETSC_EXTERN PetscErrorCode PetscLayoutDuplicate(PetscLayout,PetscLayout*);
250: PETSC_EXTERN PetscErrorCode PetscLayoutReference(PetscLayout,PetscLayout*);
251: PETSC_EXTERN PetscErrorCode PetscLayoutSetLocalSize(PetscLayout,PetscInt);
252: PETSC_EXTERN PetscErrorCode PetscLayoutGetLocalSize(PetscLayout,PetscInt *);
253: PETSC_EXTERN PetscErrorCode PetscLayoutSetSize(PetscLayout,PetscInt);
254: PETSC_EXTERN PetscErrorCode PetscLayoutGetSize(PetscLayout,PetscInt *);
255: PETSC_EXTERN PetscErrorCode PetscLayoutSetBlockSize(PetscLayout,PetscInt);
256: PETSC_EXTERN PetscErrorCode PetscLayoutGetBlockSize(PetscLayout,PetscInt*);
257: PETSC_EXTERN PetscErrorCode PetscLayoutGetRange(PetscLayout,PetscInt *,PetscInt *);
258: PETSC_EXTERN PetscErrorCode PetscLayoutGetRanges(PetscLayout,const PetscInt *[]);
259: PETSC_EXTERN PetscErrorCode PetscLayoutSetISLocalToGlobalMapping(PetscLayout,ISLocalToGlobalMapping);
260: PETSC_EXTERN PetscErrorCode PetscLayoutSetISLocalToGlobalMappingBlock(PetscLayout,ISLocalToGlobalMapping);
261: PETSC_EXTERN PetscErrorCode PetscSFSetGraphLayout(PetscSF,PetscLayout,PetscInt,const PetscInt*,PetscCopyMode,const PetscInt*);

265: /*@C
266:      PetscLayoutFindOwner - Find the owning rank for a global index

268:     Not Collective

270:    Input Parameters:
271: +    map - the layout
272: -    idx - global index to find the owner of

274:    Output Parameter:
275: .    owner - the owning rank

277:    Level: developer

279:     Fortran Notes:
280:       Not available from Fortran

282: @*/
283: PETSC_STATIC_INLINE PetscErrorCode PetscLayoutFindOwner(PetscLayout map,PetscInt idx,PetscInt *owner)
284: {
286:   PetscMPIInt    lo = 0,hi,t;

289:   *owner = -1;                  /* GCC erroneously issues warning about possibly uninitialized use when error condition */
290:   if (!((map->n >= 0) && (map->N >= 0) && (map->range))) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"PetscLayoutSetUp() must be called first");
291:   if (idx < 0 || idx > map->N) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Index %D is out of range",idx);
292:   MPI_Comm_size(map->comm,&hi);
293:   while (hi - lo > 1) {
294:     t = lo + (hi - lo) / 2;
295:     if (idx < map->range[t]) hi = t;
296:     else                     lo = t;
297:   }
298:   *owner = lo;
299:   return(0);
300: }

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

307:     Not Collective

309:    Input Parameters:
310: +    map   - the layout
311: -    idx   - global index to find the owner of

313:    Output Parameter:
314: +    owner - the owning rank
315: -    lidx  - local index used by the owner for idx

317:    Level: developer

319:     Fortran Notes:
320:       Not available from Fortran

322: @*/
323: PETSC_STATIC_INLINE PetscErrorCode PetscLayoutFindOwnerIndex(PetscLayout map,PetscInt idx,PetscInt *owner, PetscInt *lidx)
324: {
326:   PetscMPIInt    lo = 0,hi,t;

329:   if (!((map->n >= 0) && (map->N >= 0) && (map->range))) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"PetscLayoutSetUp() must be called first");
330:   if (idx < 0 || idx > map->N) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Index %D is out of range",idx);
331:   MPI_Comm_size(map->comm,&hi);
332:   while (hi - lo > 1) {
333:     t = lo + (hi - lo) / 2;
334:     if (idx < map->range[t]) hi = t;
335:     else                     lo = t;
336:   }
337:   *owner = lo;
338:   *lidx  = idx-map->range[lo];
339:   return(0);
340: }

342: PETSC_EXTERN PetscClassId PETSC_SECTION_CLASSID;

344: /*S
345:   PetscSection - Mapping from integers in a designated range to contiguous sets of integers.

347:   In contrast to IS, which maps from integers to single integers, the range of a PetscSection is in the space of
348:   contiguous sets of integers. These ranges are frequently interpreted as domains of other array-like objects,
349:   especially other PetscSections, Vecs, and ISs. The domain is set with PetscSectionSetChart() and does not need to
350:   start at 0. For each point in the domain of a PetscSection, the output set is represented through an offset and a
351:   count, which are set using PetscSectionSetOffset() and PetscSectionSetDof() respectively. Lookup is typically using
352:   accessors or routines like VecGetValuesSection().

354:   Level: developer

356: .seealso:  PetscSectionCreate(), PetscSectionDestroy()
357: S*/
358: typedef struct _p_PetscSection *PetscSection;
359: PETSC_EXTERN PetscErrorCode PetscSectionCreate(MPI_Comm,PetscSection*);
360: PETSC_EXTERN PetscErrorCode PetscSectionClone(PetscSection, PetscSection*);
361: PETSC_EXTERN PetscErrorCode PetscSectionGetNumFields(PetscSection, PetscInt *);
362: PETSC_EXTERN PetscErrorCode PetscSectionSetNumFields(PetscSection, PetscInt);
363: PETSC_EXTERN PetscErrorCode PetscSectionGetFieldName(PetscSection, PetscInt, const char *[]);
364: PETSC_EXTERN PetscErrorCode PetscSectionSetFieldName(PetscSection, PetscInt, const char []);
365: PETSC_EXTERN PetscErrorCode PetscSectionGetFieldComponents(PetscSection, PetscInt, PetscInt *);
366: PETSC_EXTERN PetscErrorCode PetscSectionSetFieldComponents(PetscSection, PetscInt, PetscInt);
367: PETSC_EXTERN PetscErrorCode PetscSectionGetChart(PetscSection, PetscInt *, PetscInt *);
368: PETSC_EXTERN PetscErrorCode PetscSectionSetChart(PetscSection, PetscInt, PetscInt);
369: PETSC_EXTERN PetscErrorCode PetscSectionGetDof(PetscSection, PetscInt, PetscInt*);
370: PETSC_EXTERN PetscErrorCode PetscSectionSetDof(PetscSection, PetscInt, PetscInt);
371: PETSC_EXTERN PetscErrorCode PetscSectionAddDof(PetscSection, PetscInt, PetscInt);
372: PETSC_EXTERN PetscErrorCode PetscSectionGetFieldDof(PetscSection, PetscInt, PetscInt, PetscInt*);
373: PETSC_EXTERN PetscErrorCode PetscSectionSetFieldDof(PetscSection, PetscInt, PetscInt, PetscInt);
374: PETSC_EXTERN PetscErrorCode PetscSectionAddFieldDof(PetscSection, PetscInt, PetscInt, PetscInt);
375: PETSC_EXTERN PetscErrorCode PetscSectionGetConstraintDof(PetscSection, PetscInt, PetscInt*);
376: PETSC_EXTERN PetscErrorCode PetscSectionSetConstraintDof(PetscSection, PetscInt, PetscInt);
377: PETSC_EXTERN PetscErrorCode PetscSectionAddConstraintDof(PetscSection, PetscInt, PetscInt);
378: PETSC_EXTERN PetscErrorCode PetscSectionGetFieldConstraintDof(PetscSection, PetscInt, PetscInt, PetscInt*);
379: PETSC_EXTERN PetscErrorCode PetscSectionSetFieldConstraintDof(PetscSection, PetscInt, PetscInt, PetscInt);
380: PETSC_EXTERN PetscErrorCode PetscSectionAddFieldConstraintDof(PetscSection, PetscInt, PetscInt, PetscInt);
381: PETSC_EXTERN PetscErrorCode PetscSectionGetConstraintIndices(PetscSection, PetscInt, const PetscInt**);
382: PETSC_EXTERN PetscErrorCode PetscSectionSetConstraintIndices(PetscSection, PetscInt, const PetscInt*);
383: PETSC_EXTERN PetscErrorCode PetscSectionGetFieldConstraintIndices(PetscSection, PetscInt, PetscInt, const PetscInt**);
384: PETSC_EXTERN PetscErrorCode PetscSectionSetFieldConstraintIndices(PetscSection, PetscInt, PetscInt, const PetscInt*);
385: PETSC_EXTERN PetscErrorCode PetscSectionSetUpBC(PetscSection);
386: PETSC_EXTERN PetscErrorCode PetscSectionSetUp(PetscSection);
387: PETSC_EXTERN PetscErrorCode PetscSectionGetMaxDof(PetscSection, PetscInt*);
388: PETSC_EXTERN PetscErrorCode PetscSectionGetStorageSize(PetscSection, PetscInt*);
389: PETSC_EXTERN PetscErrorCode PetscSectionGetConstrainedStorageSize(PetscSection, PetscInt*);
390: PETSC_EXTERN PetscErrorCode PetscSectionGetOffset(PetscSection, PetscInt, PetscInt*);
391: PETSC_EXTERN PetscErrorCode PetscSectionSetOffset(PetscSection, PetscInt, PetscInt);
392: PETSC_EXTERN PetscErrorCode PetscSectionGetFieldOffset(PetscSection, PetscInt, PetscInt, PetscInt*);
393: PETSC_EXTERN PetscErrorCode PetscSectionSetFieldOffset(PetscSection, PetscInt, PetscInt, PetscInt);
394: PETSC_EXTERN PetscErrorCode PetscSectionGetOffsetRange(PetscSection, PetscInt *, PetscInt *);
395: PETSC_EXTERN PetscErrorCode PetscSectionView(PetscSection, PetscViewer);
396: PETSC_EXTERN PetscErrorCode PetscSectionReset(PetscSection);
397: PETSC_EXTERN PetscErrorCode PetscSectionDestroy(PetscSection*);
398: PETSC_EXTERN PetscErrorCode PetscSectionCreateGlobalSection(PetscSection, PetscSF, PetscBool, PetscSection *);
399: PETSC_EXTERN PetscErrorCode PetscSectionCreateGlobalSectionCensored(PetscSection, PetscSF, PetscBool, PetscInt, const PetscInt [], PetscSection *);
400: PETSC_EXTERN PetscErrorCode PetscSectionCreateSubsection(PetscSection, PetscInt, PetscInt [], PetscSection *);
401: PETSC_EXTERN PetscErrorCode PetscSectionCreateSubmeshSection(PetscSection, IS, PetscSection *);
402: PETSC_EXTERN PetscErrorCode PetscSectionGetPointLayout(MPI_Comm, PetscSection, PetscLayout *);
403: PETSC_EXTERN PetscErrorCode PetscSectionGetValueLayout(MPI_Comm, PetscSection, PetscLayout *);

405: /* PetscSF support */
406: PETSC_EXTERN PetscErrorCode PetscSFConvertPartition(PetscSF, PetscSection, IS, ISLocalToGlobalMapping *, PetscSF *);
407: PETSC_EXTERN PetscErrorCode PetscSFCreateRemoteOffsets(PetscSF, PetscSection, PetscSection, PetscInt **);
408: PETSC_EXTERN PetscErrorCode PetscSFDistributeSection(PetscSF, PetscSection, PetscInt **, PetscSection);
409: PETSC_EXTERN PetscErrorCode PetscSFCreateSectionSF(PetscSF, PetscSection, PetscInt [], PetscSection, PetscSF *);

411: /* Reset __FUNCT__ in case the user does not define it themselves */

415: #endif