Actual source code: petscis.h
petsc-3.5.4 2015-05-23
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 <petscsftypes.h>
10: #define IS_FILE_CLASSID 1211218
11: PETSC_EXTERN PetscClassId IS_CLASSID;
13: PETSC_EXTERN PetscErrorCode ISInitializePackage(void);
15: /*S
16: IS - Abstract PETSc object that allows indexing.
18: Level: beginner
20: Concepts: indexing, stride
22: .seealso: ISCreateGeneral(), ISCreateBlock(), ISCreateStride(), ISGetIndices(), ISDestroy()
23: S*/
24: typedef struct _p_IS* IS;
26: /*J
27: ISType - String with the name of a PETSc index set type
29: Level: beginner
31: .seealso: ISSetType(), IS, ISCreate(), ISRegister()
32: J*/
33: typedef const char* ISType;
34: #define ISGENERAL "general"
35: #define ISSTRIDE "stride"
36: #define ISBLOCK "block"
38: /* Dynamic creation and loading functions */
39: PETSC_EXTERN PetscFunctionList ISList;
40: PETSC_EXTERN PetscBool ISRegisterAllCalled;
41: PETSC_EXTERN PetscErrorCode ISSetType(IS, ISType);
42: PETSC_EXTERN PetscErrorCode ISGetType(IS, ISType *);
43: PETSC_EXTERN PetscErrorCode ISRegister(const char[],PetscErrorCode (*)(IS));
44: PETSC_EXTERN PetscErrorCode ISRegisterAll(void);
45: PETSC_EXTERN PetscErrorCode ISCreate(MPI_Comm,IS*);
47: /*
48: Default index set data structures that PETSc provides.
49: */
50: PETSC_EXTERN PetscErrorCode ISCreateGeneral(MPI_Comm,PetscInt,const PetscInt[],PetscCopyMode,IS *);
51: PETSC_EXTERN PetscErrorCode ISGeneralSetIndices(IS,PetscInt,const PetscInt[],PetscCopyMode);
52: PETSC_EXTERN PetscErrorCode ISCreateBlock(MPI_Comm,PetscInt,PetscInt,const PetscInt[],PetscCopyMode,IS *);
53: PETSC_EXTERN PetscErrorCode ISBlockSetIndices(IS,PetscInt,PetscInt,const PetscInt[],PetscCopyMode);
54: PETSC_EXTERN PetscErrorCode ISCreateStride(MPI_Comm,PetscInt,PetscInt,PetscInt,IS *);
55: PETSC_EXTERN PetscErrorCode ISStrideSetStride(IS,PetscInt,PetscInt,PetscInt);
57: PETSC_EXTERN PetscErrorCode ISDestroy(IS*);
58: PETSC_EXTERN PetscErrorCode ISSetPermutation(IS);
59: PETSC_EXTERN PetscErrorCode ISPermutation(IS,PetscBool *);
60: PETSC_EXTERN PetscErrorCode ISSetIdentity(IS);
61: PETSC_EXTERN PetscErrorCode ISIdentity(IS,PetscBool *);
62: PETSC_EXTERN PetscErrorCode ISContiguousLocal(IS,PetscInt,PetscInt,PetscInt*,PetscBool*);
64: PETSC_EXTERN PetscErrorCode ISGetIndices(IS,const PetscInt *[]);
65: PETSC_EXTERN PetscErrorCode ISRestoreIndices(IS,const PetscInt *[]);
66: PETSC_EXTERN PetscErrorCode ISGetTotalIndices(IS,const PetscInt *[]);
67: PETSC_EXTERN PetscErrorCode ISRestoreTotalIndices(IS,const PetscInt *[]);
68: PETSC_EXTERN PetscErrorCode ISGetNonlocalIndices(IS,const PetscInt *[]);
69: PETSC_EXTERN PetscErrorCode ISRestoreNonlocalIndices(IS,const PetscInt *[]);
70: PETSC_EXTERN PetscErrorCode ISGetNonlocalIS(IS, IS *is);
71: PETSC_EXTERN PetscErrorCode ISRestoreNonlocalIS(IS, IS *is);
72: PETSC_EXTERN PetscErrorCode ISGetSize(IS,PetscInt *);
73: PETSC_EXTERN PetscErrorCode ISGetLocalSize(IS,PetscInt *);
74: PETSC_EXTERN PetscErrorCode ISInvertPermutation(IS,PetscInt,IS*);
75: PETSC_EXTERN PetscErrorCode ISView(IS,PetscViewer);
76: PETSC_STATIC_INLINE PetscErrorCode ISViewFromOptions(IS A,const char prefix[],const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,prefix,name);}
77: PETSC_EXTERN PetscErrorCode ISLoad(IS,PetscViewer);
78: PETSC_EXTERN PetscErrorCode ISEqual(IS,IS,PetscBool *);
79: PETSC_EXTERN PetscErrorCode ISSort(IS);
80: PETSC_EXTERN PetscErrorCode ISSortRemoveDups(IS);
81: PETSC_EXTERN PetscErrorCode ISSorted(IS,PetscBool *);
82: PETSC_EXTERN PetscErrorCode ISDifference(IS,IS,IS*);
83: PETSC_EXTERN PetscErrorCode ISSum(IS,IS,IS*);
84: PETSC_EXTERN PetscErrorCode ISExpand(IS,IS,IS*);
85: PETSC_EXTERN PetscErrorCode ISGetMinMax(IS,PetscInt*,PetscInt*);
87: PETSC_EXTERN PetscErrorCode ISBlockGetIndices(IS,const PetscInt *[]);
88: PETSC_EXTERN PetscErrorCode ISBlockRestoreIndices(IS,const PetscInt *[]);
89: PETSC_EXTERN PetscErrorCode ISBlockGetLocalSize(IS,PetscInt *);
90: PETSC_EXTERN PetscErrorCode ISBlockGetSize(IS,PetscInt *);
91: PETSC_EXTERN PetscErrorCode ISGetBlockSize(IS,PetscInt*);
92: PETSC_EXTERN PetscErrorCode ISSetBlockSize(IS,PetscInt);
94: PETSC_EXTERN PetscErrorCode ISStrideGetInfo(IS,PetscInt *,PetscInt*);
96: PETSC_EXTERN PetscErrorCode ISToGeneral(IS);
98: PETSC_EXTERN PetscErrorCode ISDuplicate(IS,IS*);
99: PETSC_EXTERN PetscErrorCode ISCopy(IS,IS);
100: PETSC_EXTERN PetscErrorCode ISAllGather(IS,IS*);
101: PETSC_EXTERN PetscErrorCode ISComplement(IS,PetscInt,PetscInt,IS*);
102: PETSC_EXTERN PetscErrorCode ISConcatenate(MPI_Comm,PetscInt,const IS[],IS*);
103: PETSC_EXTERN PetscErrorCode ISListToPair(MPI_Comm,PetscInt, IS[],IS*,IS*);
104: PETSC_EXTERN PetscErrorCode ISPairToList(IS,IS,PetscInt*, IS *[]);
105: PETSC_EXTERN PetscErrorCode ISEmbed(IS,IS,PetscBool,IS*);
106: PETSC_EXTERN PetscErrorCode ISOnComm(IS,MPI_Comm,PetscCopyMode,IS*);
108: /* --------------------------------------------------------------------------*/
109: PETSC_EXTERN PetscClassId IS_LTOGM_CLASSID;
111: /*S
112: ISLocalToGlobalMapping - mappings from an arbitrary
113: local ordering from 0 to n-1 to a global PETSc ordering
114: used by a vector or matrix.
116: Level: intermediate
118: Note: mapping from Local to Global is scalable; but Global
119: to Local may not be if the range of global values represented locally
120: is very large.
122: Note: the ISLocalToGlobalMapping is actually a private object; it is included
123: here for the inline function ISLocalToGlobalMappingApply() to allow it to be inlined since
124: it is used so often.
126: .seealso: ISLocalToGlobalMappingCreate()
127: S*/
128: typedef struct _p_ISLocalToGlobalMapping* ISLocalToGlobalMapping;
130: /*E
131: ISGlobalToLocalMappingType - Indicates if missing global indices are
133: IS_GTOLM_MASK - missing global indices are replaced with -1
134: IS_GTOLM_DROP - missing global indices are dropped
136: Level: beginner
138: .seealso: ISGlobalToLocalMappingApplyBlock()
140: E*/
141: typedef enum {IS_GTOLM_MASK,IS_GTOLM_DROP} ISGlobalToLocalMappingType;
143: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingCreate(MPI_Comm,PetscInt,PetscInt,const PetscInt[],PetscCopyMode,ISLocalToGlobalMapping*);
144: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingCreateIS(IS,ISLocalToGlobalMapping *);
145: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingCreateSF(PetscSF,PetscInt,ISLocalToGlobalMapping*);
146: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingView(ISLocalToGlobalMapping,PetscViewer);
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,ISGlobalToLocalMappingType,PetscInt,const PetscInt[],PetscInt*,PetscInt[]);
152: PETSC_EXTERN PetscErrorCode ISGlobalToLocalMappingApplyBlock(ISLocalToGlobalMapping,ISGlobalToLocalMappingType,PetscInt,const PetscInt[],PetscInt*,PetscInt[]);
153: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingGetSize(ISLocalToGlobalMapping,PetscInt*);
154: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping,PetscInt*,PetscInt*[],PetscInt*[],PetscInt**[]);
155: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingRestoreInfo(ISLocalToGlobalMapping,PetscInt*,PetscInt*[],PetscInt*[],PetscInt**[]);
156: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingGetBlockInfo(ISLocalToGlobalMapping,PetscInt*,PetscInt*[],PetscInt*[],PetscInt**[]);
157: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingRestoreBlockInfo(ISLocalToGlobalMapping,PetscInt*,PetscInt*[],PetscInt*[],PetscInt**[]);
158: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingGetIndices(ISLocalToGlobalMapping,const PetscInt**);
159: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingRestoreIndices(ISLocalToGlobalMapping,const PetscInt**);
160: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingGetBlockIndices(ISLocalToGlobalMapping,const PetscInt**);
161: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingRestoreBlockIndices(ISLocalToGlobalMapping,const PetscInt**);
162: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingConcatenate(MPI_Comm,PetscInt,const ISLocalToGlobalMapping[],ISLocalToGlobalMapping*);
163: PETSC_EXTERN PetscErrorCode ISG2LMapApply(ISLocalToGlobalMapping,PetscInt,const PetscInt[],PetscInt[]);
164: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingGetBlockSize(ISLocalToGlobalMapping,PetscInt*);
167: /* --------------------------------------------------------------------------*/
168: /*E
169: ISColoringType - determines if the coloring is for the entire parallel grid/graph/matrix
170: or for just the local ghosted portion
172: Level: beginner
174: $ IS_COLORING_GLOBAL - does not include the colors for ghost points, this is used when the function
175: $ is called synchronously in parallel. This requires generating a "parallel coloring".
176: $ IS_COLORING_GHOSTED - includes colors for ghost points, this is used when the function can be called
177: $ seperately on individual processes with the ghost points already filled in. Does not
178: $ require a "parallel coloring", rather each process colors its local + ghost part.
179: $ Using this can result in much less parallel communication. In the paradigm of
180: $ DMGetLocalVector() and DMGetGlobalVector() this could be called IS_COLORING_LOCAL
182: .seealso: DMCreateColoring()
183: E*/
184: typedef enum {IS_COLORING_GLOBAL,IS_COLORING_GHOSTED} ISColoringType;
185: PETSC_EXTERN const char *const ISColoringTypes[];
186: typedef unsigned PETSC_IS_COLOR_VALUE_TYPE ISColoringValue;
187: PETSC_EXTERN PetscErrorCode ISAllGatherColors(MPI_Comm,PetscInt,ISColoringValue*,PetscInt*,ISColoringValue*[]);
189: /*S
190: ISColoring - sets of IS's that define a coloring
191: of the underlying indices
193: Level: intermediate
195: Notes:
196: One should not access the *is records below directly because they may not yet
197: have been created. One should use ISColoringGetIS() to make sure they are
198: created when needed.
200: Developer Note: this is not a PetscObject
202: .seealso: ISColoringCreate(), ISColoringGetIS(), ISColoringView(), ISColoringGetIS()
203: S*/
204: struct _n_ISColoring {
205: PetscInt refct;
206: PetscInt n; /* number of colors */
207: IS *is; /* for each color indicates columns */
208: MPI_Comm comm;
209: ISColoringValue *colors; /* for each column indicates color */
210: PetscInt N; /* number of columns */
211: ISColoringType ctype;
212: };
213: typedef struct _n_ISColoring* ISColoring;
215: PETSC_EXTERN PetscErrorCode ISColoringCreate(MPI_Comm,PetscInt,PetscInt,const ISColoringValue[],ISColoring*);
216: PETSC_EXTERN PetscErrorCode ISColoringDestroy(ISColoring*);
217: PETSC_EXTERN PetscErrorCode ISColoringView(ISColoring,PetscViewer);
218: PETSC_EXTERN PetscErrorCode ISColoringViewFromOptions(ISColoring,const char[],const char[]);
219: PETSC_EXTERN PetscErrorCode ISColoringGetIS(ISColoring,PetscInt*,IS*[]);
220: PETSC_EXTERN PetscErrorCode ISColoringRestoreIS(ISColoring,IS*[]);
221: PETSC_EXTERN PetscErrorCode ISColoringReference(ISColoring);
222: PETSC_EXTERN PetscErrorCode ISColoringSetType(ISColoring,ISColoringType);
225: /* --------------------------------------------------------------------------*/
227: PETSC_EXTERN PetscErrorCode ISPartitioningToNumbering(IS,IS*);
228: PETSC_EXTERN PetscErrorCode ISPartitioningCount(IS,PetscInt,PetscInt[]);
230: PETSC_EXTERN PetscErrorCode ISCompressIndicesGeneral(PetscInt,PetscInt,PetscInt,PetscInt,const IS[],IS[]);
231: PETSC_EXTERN PetscErrorCode ISCompressIndicesSorted(PetscInt,PetscInt,PetscInt,const IS[],IS[]);
232: PETSC_EXTERN PetscErrorCode ISExpandIndicesGeneral(PetscInt,PetscInt,PetscInt,PetscInt,const IS[],IS[]);
234: /*S
235: PetscLayout - defines layout of vectors and matrices across processes (which rows are owned by which processes)
237: Level: developer
240: .seealso: PetscLayoutCreate(), PetscLayoutDestroy()
241: S*/
242: typedef struct _n_PetscLayout* PetscLayout;
243: struct _n_PetscLayout{
244: MPI_Comm comm;
245: PetscInt n,N; /* local, global vector size */
246: PetscInt rstart,rend; /* local start, local end + 1 */
247: PetscInt *range; /* the offset of each processor */
248: PetscInt bs; /* number of elements in each block (generally for multi-component
249: * problems). Defaults to -1 and can be arbitrarily lazy so always use
250: * PetscAbs(map->bs) when accessing directly and expecting result to be
251: * positive. Do NOT multiply above numbers by bs */
252: PetscInt refcnt; /* MPI Vecs obtained with VecDuplicate() and from MatGetVecs() reuse map of input object */
253: ISLocalToGlobalMapping mapping; /* mapping used in Vec/MatSetValuesLocal() */
254: PetscInt *trstarts; /* local start for each thread */
255: };
257: PETSC_EXTERN PetscErrorCode PetscLayoutCreate(MPI_Comm,PetscLayout*);
258: PETSC_EXTERN PetscErrorCode PetscLayoutSetUp(PetscLayout);
259: PETSC_EXTERN PetscErrorCode PetscLayoutDestroy(PetscLayout*);
260: PETSC_EXTERN PetscErrorCode PetscLayoutDuplicate(PetscLayout,PetscLayout*);
261: PETSC_EXTERN PetscErrorCode PetscLayoutReference(PetscLayout,PetscLayout*);
262: PETSC_EXTERN PetscErrorCode PetscLayoutSetLocalSize(PetscLayout,PetscInt);
263: PETSC_EXTERN PetscErrorCode PetscLayoutGetLocalSize(PetscLayout,PetscInt *);
264: PETSC_EXTERN PetscErrorCode PetscLayoutSetSize(PetscLayout,PetscInt);
265: PETSC_EXTERN PetscErrorCode PetscLayoutGetSize(PetscLayout,PetscInt *);
266: PETSC_EXTERN PetscErrorCode PetscLayoutSetBlockSize(PetscLayout,PetscInt);
267: PETSC_EXTERN PetscErrorCode PetscLayoutGetBlockSize(PetscLayout,PetscInt*);
268: PETSC_EXTERN PetscErrorCode PetscLayoutGetRange(PetscLayout,PetscInt *,PetscInt *);
269: PETSC_EXTERN PetscErrorCode PetscLayoutGetRanges(PetscLayout,const PetscInt *[]);
270: PETSC_EXTERN PetscErrorCode PetscLayoutSetISLocalToGlobalMapping(PetscLayout,ISLocalToGlobalMapping);
271: PETSC_EXTERN PetscErrorCode PetscSFSetGraphLayout(PetscSF,PetscLayout,PetscInt,const PetscInt*,PetscCopyMode,const PetscInt*);
275: /*@C
276: PetscLayoutFindOwner - Find the owning rank for a global index
278: Not Collective
280: Input Parameters:
281: + map - the layout
282: - idx - global index to find the owner of
284: Output Parameter:
285: . owner - the owning rank
287: Level: developer
289: Fortran Notes:
290: Not available from Fortran
292: @*/
293: PETSC_STATIC_INLINE PetscErrorCode PetscLayoutFindOwner(PetscLayout map,PetscInt idx,PetscInt *owner)
294: {
296: PetscMPIInt lo = 0,hi,t;
299: *owner = -1; /* GCC erroneously issues warning about possibly uninitialized use when error condition */
300: if (!((map->n >= 0) && (map->N >= 0) && (map->range))) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"PetscLayoutSetUp() must be called first");
301: if (idx < 0 || idx > map->N) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Index %D is out of range",idx);
302: MPI_Comm_size(map->comm,&hi);
303: while (hi - lo > 1) {
304: t = lo + (hi - lo) / 2;
305: if (idx < map->range[t]) hi = t;
306: else lo = t;
307: }
308: *owner = lo;
309: return(0);
310: }
314: /*@C
315: PetscLayoutFindOwnerIndex - Find the owning rank and the local index for a global index
317: Not Collective
319: Input Parameters:
320: + map - the layout
321: - idx - global index to find the owner of
323: Output Parameter:
324: + owner - the owning rank
325: - lidx - local index used by the owner for idx
327: Level: developer
329: Fortran Notes:
330: Not available from Fortran
332: @*/
333: PETSC_STATIC_INLINE PetscErrorCode PetscLayoutFindOwnerIndex(PetscLayout map,PetscInt idx,PetscInt *owner, PetscInt *lidx)
334: {
336: PetscMPIInt lo = 0,hi,t;
339: if (!((map->n >= 0) && (map->N >= 0) && (map->range))) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"PetscLayoutSetUp() must be called first");
340: if (idx < 0 || idx > map->N) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Index %D is out of range",idx);
341: MPI_Comm_size(map->comm,&hi);
342: while (hi - lo > 1) {
343: t = lo + (hi - lo) / 2;
344: if (idx < map->range[t]) hi = t;
345: else lo = t;
346: }
347: *owner = lo;
348: *lidx = idx-map->range[lo];
349: return(0);
350: }
352: PETSC_EXTERN PetscClassId PETSC_SECTION_CLASSID;
354: /*S
355: PetscSection - Mapping from integers in a designated range to contiguous sets of integers.
357: In contrast to IS, which maps from integers to single integers, the range of a PetscSection is in the space of
358: contiguous sets of integers. These ranges are frequently interpreted as domains of other array-like objects,
359: especially other PetscSections, Vecs, and ISs. The domain is set with PetscSectionSetChart() and does not need to
360: start at 0. For each point in the domain of a PetscSection, the output set is represented through an offset and a
361: count, which are set using PetscSectionSetOffset() and PetscSectionSetDof() respectively. Lookup is typically using
362: accessors or routines like VecGetValuesSection().
364: Level: developer
366: .seealso: PetscSectionCreate(), PetscSectionDestroy()
367: S*/
368: typedef struct _p_PetscSection *PetscSection;
369: PETSC_EXTERN PetscErrorCode PetscSectionCreate(MPI_Comm,PetscSection*);
370: PETSC_EXTERN PetscErrorCode PetscSectionClone(PetscSection, PetscSection*);
371: PETSC_EXTERN PetscErrorCode PetscSectionGetNumFields(PetscSection, PetscInt *);
372: PETSC_EXTERN PetscErrorCode PetscSectionSetNumFields(PetscSection, PetscInt);
373: PETSC_EXTERN PetscErrorCode PetscSectionGetFieldName(PetscSection, PetscInt, const char *[]);
374: PETSC_EXTERN PetscErrorCode PetscSectionSetFieldName(PetscSection, PetscInt, const char []);
375: PETSC_EXTERN PetscErrorCode PetscSectionGetFieldComponents(PetscSection, PetscInt, PetscInt *);
376: PETSC_EXTERN PetscErrorCode PetscSectionSetFieldComponents(PetscSection, PetscInt, PetscInt);
377: PETSC_EXTERN PetscErrorCode PetscSectionGetChart(PetscSection, PetscInt *, PetscInt *);
378: PETSC_EXTERN PetscErrorCode PetscSectionSetChart(PetscSection, PetscInt, PetscInt);
379: PETSC_EXTERN PetscErrorCode PetscSectionGetPermutation(PetscSection, IS *);
380: PETSC_EXTERN PetscErrorCode PetscSectionSetPermutation(PetscSection, IS);
381: PETSC_EXTERN PetscErrorCode PetscSectionGetDof(PetscSection, PetscInt, PetscInt*);
382: PETSC_EXTERN PetscErrorCode PetscSectionSetDof(PetscSection, PetscInt, PetscInt);
383: PETSC_EXTERN PetscErrorCode PetscSectionAddDof(PetscSection, PetscInt, PetscInt);
384: PETSC_EXTERN PetscErrorCode PetscSectionGetFieldDof(PetscSection, PetscInt, PetscInt, PetscInt*);
385: PETSC_EXTERN PetscErrorCode PetscSectionSetFieldDof(PetscSection, PetscInt, PetscInt, PetscInt);
386: PETSC_EXTERN PetscErrorCode PetscSectionAddFieldDof(PetscSection, PetscInt, PetscInt, PetscInt);
387: PETSC_EXTERN PetscErrorCode PetscSectionHasConstraints(PetscSection, PetscBool *);
388: PETSC_EXTERN PetscErrorCode PetscSectionGetConstraintDof(PetscSection, PetscInt, PetscInt*);
389: PETSC_EXTERN PetscErrorCode PetscSectionSetConstraintDof(PetscSection, PetscInt, PetscInt);
390: PETSC_EXTERN PetscErrorCode PetscSectionAddConstraintDof(PetscSection, PetscInt, PetscInt);
391: PETSC_EXTERN PetscErrorCode PetscSectionGetFieldConstraintDof(PetscSection, PetscInt, PetscInt, PetscInt*);
392: PETSC_EXTERN PetscErrorCode PetscSectionSetFieldConstraintDof(PetscSection, PetscInt, PetscInt, PetscInt);
393: PETSC_EXTERN PetscErrorCode PetscSectionAddFieldConstraintDof(PetscSection, PetscInt, PetscInt, PetscInt);
394: PETSC_EXTERN PetscErrorCode PetscSectionGetConstraintIndices(PetscSection, PetscInt, const PetscInt**);
395: PETSC_EXTERN PetscErrorCode PetscSectionSetConstraintIndices(PetscSection, PetscInt, const PetscInt*);
396: PETSC_EXTERN PetscErrorCode PetscSectionGetFieldConstraintIndices(PetscSection, PetscInt, PetscInt, const PetscInt**);
397: PETSC_EXTERN PetscErrorCode PetscSectionSetFieldConstraintIndices(PetscSection, PetscInt, PetscInt, const PetscInt*);
398: PETSC_EXTERN PetscErrorCode PetscSectionSetUpBC(PetscSection);
399: PETSC_EXTERN PetscErrorCode PetscSectionSetUp(PetscSection);
400: PETSC_EXTERN PetscErrorCode PetscSectionGetMaxDof(PetscSection, PetscInt*);
401: PETSC_EXTERN PetscErrorCode PetscSectionGetStorageSize(PetscSection, PetscInt*);
402: PETSC_EXTERN PetscErrorCode PetscSectionGetConstrainedStorageSize(PetscSection, PetscInt*);
403: PETSC_EXTERN PetscErrorCode PetscSectionGetOffset(PetscSection, PetscInt, PetscInt*);
404: PETSC_EXTERN PetscErrorCode PetscSectionSetOffset(PetscSection, PetscInt, PetscInt);
405: PETSC_EXTERN PetscErrorCode PetscSectionGetFieldOffset(PetscSection, PetscInt, PetscInt, PetscInt*);
406: PETSC_EXTERN PetscErrorCode PetscSectionSetFieldOffset(PetscSection, PetscInt, PetscInt, PetscInt);
407: PETSC_EXTERN PetscErrorCode PetscSectionGetOffsetRange(PetscSection, PetscInt *, PetscInt *);
408: PETSC_EXTERN PetscErrorCode PetscSectionView(PetscSection, PetscViewer);
409: PETSC_STATIC_INLINE PetscErrorCode PetscSectionViewFromOptions(PetscSection A,const char prefix[],const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,prefix,name);}
410: PETSC_EXTERN PetscErrorCode PetscSectionReset(PetscSection);
411: PETSC_EXTERN PetscErrorCode PetscSectionDestroy(PetscSection*);
412: PETSC_EXTERN PetscErrorCode PetscSectionCreateGlobalSection(PetscSection, PetscSF, PetscBool, PetscSection *);
413: PETSC_EXTERN PetscErrorCode PetscSectionCreateGlobalSectionCensored(PetscSection, PetscSF, PetscBool, PetscInt, const PetscInt [], PetscSection *);
414: PETSC_EXTERN PetscErrorCode PetscSectionCreateSubsection(PetscSection, PetscInt, PetscInt [], PetscSection *);
415: PETSC_EXTERN PetscErrorCode PetscSectionCreateSubmeshSection(PetscSection, IS, PetscSection *);
416: PETSC_EXTERN PetscErrorCode PetscSectionGetPointLayout(MPI_Comm, PetscSection, PetscLayout *);
417: PETSC_EXTERN PetscErrorCode PetscSectionGetValueLayout(MPI_Comm, PetscSection, PetscLayout *);
418: PETSC_EXTERN PetscErrorCode PetscSectionPermute(PetscSection, IS, PetscSection *);
419: PETSC_EXTERN PetscErrorCode PetscSectionGetField(PetscSection, PetscInt, PetscSection *);
421: PETSC_EXTERN PetscErrorCode PetscSectionSetClosureIndex(PetscSection, PetscObject, PetscSection, IS);
422: PETSC_EXTERN PetscErrorCode PetscSectionGetClosureIndex(PetscSection, PetscObject, PetscSection *, IS *);
424: /* PetscSF support */
425: PETSC_EXTERN PetscErrorCode PetscSFConvertPartition(PetscSF, PetscSection, IS, ISLocalToGlobalMapping *, PetscSF *);
426: PETSC_EXTERN PetscErrorCode PetscSFCreateRemoteOffsets(PetscSF, PetscSection, PetscSection, PetscInt **);
427: PETSC_EXTERN PetscErrorCode PetscSFDistributeSection(PetscSF, PetscSection, PetscInt **, PetscSection);
428: PETSC_EXTERN PetscErrorCode PetscSFCreateSectionSF(PetscSF, PetscSection, PetscInt [], PetscSection, PetscSF *);
430: /* Reset __FUNCT__ in case the user does not define it themselves */
434: #endif