Actual source code: petscis.h

petsc-3.8.4 2018-03-24
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: */
  7:  #include <petscsys.h>
  8:  #include <petscsftypes.h>
  9:  #include <petscistypes.h>

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

 14: PETSC_EXTERN PetscErrorCode ISInitializePackage(void);

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

 19:    Level: beginner

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

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

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

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

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

 76: PETSC_EXTERN PetscErrorCode ISLocate(IS,PetscInt,PetscInt*);

 78: PETSC_EXTERN PetscErrorCode ISBlockGetIndices(IS,const PetscInt *[]);
 79: PETSC_EXTERN PetscErrorCode ISBlockRestoreIndices(IS,const PetscInt *[]);
 80: PETSC_EXTERN PetscErrorCode ISBlockGetLocalSize(IS,PetscInt *);
 81: PETSC_EXTERN PetscErrorCode ISBlockGetSize(IS,PetscInt *);
 82: PETSC_EXTERN PetscErrorCode ISGetBlockSize(IS,PetscInt*);
 83: PETSC_EXTERN PetscErrorCode ISSetBlockSize(IS,PetscInt);

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

 87: PETSC_EXTERN PetscErrorCode ISToGeneral(IS);

 89: PETSC_EXTERN PetscErrorCode ISDuplicate(IS,IS*);
 90: PETSC_EXTERN PetscErrorCode ISCopy(IS,IS);
 91: PETSC_EXTERN PetscErrorCode ISAllGather(IS,IS*);
 92: PETSC_EXTERN PetscErrorCode ISComplement(IS,PetscInt,PetscInt,IS*);
 93: PETSC_EXTERN PetscErrorCode ISConcatenate(MPI_Comm,PetscInt,const IS[],IS*);
 94: PETSC_EXTERN PetscErrorCode ISListToPair(MPI_Comm,PetscInt, IS[],IS*,IS*);
 95: PETSC_EXTERN PetscErrorCode ISPairToList(IS,IS,PetscInt*, IS *[]);
 96: PETSC_EXTERN PetscErrorCode ISEmbed(IS,IS,PetscBool,IS*);
 97: PETSC_EXTERN PetscErrorCode ISSortPermutation(IS,PetscBool,IS*);
 98: PETSC_EXTERN PetscErrorCode ISOnComm(IS,MPI_Comm,PetscCopyMode,IS*);
 99: PETSC_EXTERN PetscErrorCode ISRenumber(IS,IS,PetscInt*,IS*);

101: /* --------------------------------------------------------------------------*/
102: PETSC_EXTERN PetscClassId IS_LTOGM_CLASSID;

104: /*E
105:     ISGlobalToLocalMappingMode - Indicates if missing global indices are

107:    IS_GTOLM_MASK - missing global indices are replaced with -1
108:    IS_GTOLM_DROP - missing global indices are dropped

110:    Level: beginner

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

114: E*/
115: typedef enum {IS_GTOLM_MASK,IS_GTOLM_DROP} ISGlobalToLocalMappingMode;

117: /*J
118:     ISLocalToGlobalMappingType - String with the name of a mapping method

120:    Level: beginner

122: .seealso: ISLocalToGlobalMappingSetType(), ISLocalToGlobalSetFromOptions()
123: J*/
124: typedef const char* ISLocalToGlobalMappingType;
125: #define ISLOCALTOGLOBALMAPPINGBASIC "basic"
126: #define ISLOCALTOGLOBALMAPPINGHASH  "hash"

128: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingSetType(ISLocalToGlobalMapping,ISLocalToGlobalMappingType);
129: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingRegisterAll(void);
130: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingCreate(MPI_Comm,PetscInt,PetscInt,const PetscInt[],PetscCopyMode,ISLocalToGlobalMapping*);
131: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingCreateIS(IS,ISLocalToGlobalMapping *);
132: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingCreateSF(PetscSF,PetscInt,ISLocalToGlobalMapping*);
133: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingSetFromOptions(ISLocalToGlobalMapping);
134: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingSetUp(ISLocalToGlobalMapping);
135: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingView(ISLocalToGlobalMapping,PetscViewer);
136: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingDestroy(ISLocalToGlobalMapping*);
137: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingApply(ISLocalToGlobalMapping,PetscInt,const PetscInt[],PetscInt[]);
138: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingApplyBlock(ISLocalToGlobalMapping,PetscInt,const PetscInt[],PetscInt[]);
139: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingApplyIS(ISLocalToGlobalMapping,IS,IS*);
140: PETSC_EXTERN PetscErrorCode ISGlobalToLocalMappingApply(ISLocalToGlobalMapping,ISGlobalToLocalMappingMode,PetscInt,const PetscInt[],PetscInt*,PetscInt[]);
141: PETSC_EXTERN PetscErrorCode ISGlobalToLocalMappingApplyBlock(ISLocalToGlobalMapping,ISGlobalToLocalMappingMode,PetscInt,const PetscInt[],PetscInt*,PetscInt[]);
142: PETSC_EXTERN PetscErrorCode ISGlobalToLocalMappingApplyIS(ISLocalToGlobalMapping,ISGlobalToLocalMappingMode,IS,IS*);
143: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingGetSize(ISLocalToGlobalMapping,PetscInt*);
144: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping,PetscInt*,PetscInt*[],PetscInt*[],PetscInt**[]);
145: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingRestoreInfo(ISLocalToGlobalMapping,PetscInt*,PetscInt*[],PetscInt*[],PetscInt**[]);
146: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingGetBlockInfo(ISLocalToGlobalMapping,PetscInt*,PetscInt*[],PetscInt*[],PetscInt**[]);
147: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingRestoreBlockInfo(ISLocalToGlobalMapping,PetscInt*,PetscInt*[],PetscInt*[],PetscInt**[]);
148: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingGetIndices(ISLocalToGlobalMapping,const PetscInt**);
149: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingRestoreIndices(ISLocalToGlobalMapping,const PetscInt**);
150: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingGetBlockIndices(ISLocalToGlobalMapping,const PetscInt**);
151: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingRestoreBlockIndices(ISLocalToGlobalMapping,const PetscInt**);
152: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingConcatenate(MPI_Comm,PetscInt,const ISLocalToGlobalMapping[],ISLocalToGlobalMapping*);
153: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingGetBlockSize(ISLocalToGlobalMapping,PetscInt*);
154: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingSetBlockSize(ISLocalToGlobalMapping,PetscInt);
155: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingDuplicate(ISLocalToGlobalMapping,ISLocalToGlobalMapping*);


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

163:     Level: beginner

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

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

180: PETSC_EXTERN PetscErrorCode ISColoringCreate(MPI_Comm,PetscInt,PetscInt,const ISColoringValue[],PetscCopyMode,ISColoring*);
181: PETSC_EXTERN PetscErrorCode ISColoringDestroy(ISColoring*);
182: PETSC_EXTERN PetscErrorCode ISColoringView(ISColoring,PetscViewer);
183: PETSC_EXTERN PetscErrorCode ISColoringViewFromOptions(ISColoring,PetscObject,const char[]);
184: PETSC_EXTERN PetscErrorCode ISColoringGetIS(ISColoring,PetscInt*,IS*[]);
185: PETSC_EXTERN PetscErrorCode ISColoringRestoreIS(ISColoring,IS*[]);
186: PETSC_EXTERN PetscErrorCode ISColoringReference(ISColoring);
187: PETSC_EXTERN PetscErrorCode ISColoringSetType(ISColoring,ISColoringType);


190: /* --------------------------------------------------------------------------*/
191: PETSC_EXTERN PetscErrorCode ISBuildTwoSided(IS,IS,IS*);
192: PETSC_EXTERN PetscErrorCode ISPartitioningToNumbering(IS,IS*);
193: PETSC_EXTERN PetscErrorCode ISPartitioningCount(IS,PetscInt,PetscInt[]);

195: PETSC_EXTERN PetscErrorCode ISCompressIndicesGeneral(PetscInt,PetscInt,PetscInt,PetscInt,const IS[],IS[]);
196: PETSC_EXTERN PetscErrorCode ISCompressIndicesSorted(PetscInt,PetscInt,PetscInt,const IS[],IS[]);
197: PETSC_EXTERN PetscErrorCode ISExpandIndicesGeneral(PetscInt,PetscInt,PetscInt,PetscInt,const IS[],IS[]);


200: struct _n_PetscLayout{
201:   MPI_Comm               comm;
202:   PetscInt               n,N;         /* local, global vector size */
203:   PetscInt               rstart,rend; /* local start, local end + 1 */
204:   PetscInt               *range;      /* the offset of each processor */
205:   PetscInt               bs;          /* number of elements in each block (generally for multi-component
206:                                        * problems). Defaults to -1 and can be arbitrarily lazy so always use
207:                                        * PetscAbs(map->bs) when accessing directly and expecting result to be
208:                                        * positive. Do NOT multiply above numbers by bs */
209:   PetscInt               refcnt;      /* MPI Vecs obtained with VecDuplicate() and from MatCreateVecs() reuse map of input object */
210:   ISLocalToGlobalMapping mapping;     /* mapping used in Vec/MatSetValuesLocal() */
211: };

213: /*@C
214:      PetscLayoutFindOwner - Find the owning rank for a global index

216:     Not Collective

218:    Input Parameters:
219: +    map - the layout
220: -    idx - global index to find the owner of

222:    Output Parameter:
223: .    owner - the owning rank

225:    Level: developer

227:     Fortran Notes:
228:       Not available from Fortran

230: @*/
231: PETSC_STATIC_INLINE PetscErrorCode PetscLayoutFindOwner(PetscLayout map,PetscInt idx,PetscInt *owner)
232: {
234:   PetscMPIInt    lo = 0,hi,t;

237:   *owner = -1;                  /* GCC erroneously issues warning about possibly uninitialized use when error condition */
238:   if (!((map->n >= 0) && (map->N >= 0) && (map->range))) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"PetscLayoutSetUp() must be called first");
239:   if (idx < 0 || idx > map->N) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Index %D is out of range",idx);
240:   MPI_Comm_size(map->comm,&hi);
241:   while (hi - lo > 1) {
242:     t = lo + (hi - lo) / 2;
243:     if (idx < map->range[t]) hi = t;
244:     else                     lo = t;
245:   }
246:   *owner = lo;
247:   return(0);
248: }

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

253:     Not Collective

255:    Input Parameters:
256: +    map   - the layout
257: -    idx   - global index to find the owner of

259:    Output Parameter:
260: +    owner - the owning rank
261: -    lidx  - local index used by the owner for idx

263:    Level: developer

265:     Fortran Notes:
266:       Not available from Fortran

268: @*/
269: PETSC_STATIC_INLINE PetscErrorCode PetscLayoutFindOwnerIndex(PetscLayout map,PetscInt idx,PetscInt *owner, PetscInt *lidx)
270: {
272:   PetscMPIInt    lo = 0,hi,t;

275:   if (!((map->n >= 0) && (map->N >= 0) && (map->range))) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"PetscLayoutSetUp() must be called first");
276:   if (idx < 0 || idx > map->N) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Index %D is out of range",idx);
277:   MPI_Comm_size(map->comm,&hi);
278:   while (hi - lo > 1) {
279:     t = lo + (hi - lo) / 2;
280:     if (idx < map->range[t]) hi = t;
281:     else                     lo = t;
282:   }
283:   if (owner) *owner = lo;
284:   if (lidx) *lidx  = idx-map->range[lo];
285:   return(0);
286: }

288: PETSC_EXTERN PetscErrorCode PetscLayoutCreate(MPI_Comm,PetscLayout*);
289: PETSC_EXTERN PetscErrorCode PetscLayoutSetUp(PetscLayout);
290: PETSC_EXTERN PetscErrorCode PetscLayoutDestroy(PetscLayout*);
291: PETSC_EXTERN PetscErrorCode PetscLayoutDuplicate(PetscLayout,PetscLayout*);
292: PETSC_EXTERN PetscErrorCode PetscLayoutReference(PetscLayout,PetscLayout*);
293: PETSC_EXTERN PetscErrorCode PetscLayoutSetLocalSize(PetscLayout,PetscInt);
294: PETSC_EXTERN PetscErrorCode PetscLayoutGetLocalSize(PetscLayout,PetscInt *);
295: PETSC_EXTERN PetscErrorCode PetscLayoutSetSize(PetscLayout,PetscInt);
296: PETSC_EXTERN PetscErrorCode PetscLayoutGetSize(PetscLayout,PetscInt *);
297: PETSC_EXTERN PetscErrorCode PetscLayoutSetBlockSize(PetscLayout,PetscInt);
298: PETSC_EXTERN PetscErrorCode PetscLayoutGetBlockSize(PetscLayout,PetscInt*);
299: PETSC_EXTERN PetscErrorCode PetscLayoutGetRange(PetscLayout,PetscInt *,PetscInt *);
300: PETSC_EXTERN PetscErrorCode PetscLayoutGetRanges(PetscLayout,const PetscInt *[]);
301: PETSC_EXTERN PetscErrorCode PetscLayoutCompare(PetscLayout,PetscLayout,PetscBool*);
302: PETSC_EXTERN PetscErrorCode PetscLayoutSetISLocalToGlobalMapping(PetscLayout,ISLocalToGlobalMapping);
303: PETSC_EXTERN PetscErrorCode PetscSFSetGraphLayout(PetscSF,PetscLayout,PetscInt,const PetscInt*,PetscCopyMode,const PetscInt*);

305: PETSC_EXTERN PetscClassId PETSC_SECTION_CLASSID;

307: PETSC_EXTERN PetscErrorCode PetscSectionCreate(MPI_Comm,PetscSection*);
308: PETSC_EXTERN PetscErrorCode PetscSectionClone(PetscSection, PetscSection*);
309: PETSC_EXTERN PetscErrorCode PetscSectionCopy(PetscSection, PetscSection);
310: PETSC_EXTERN PetscErrorCode PetscSectionGetNumFields(PetscSection, PetscInt *);
311: PETSC_EXTERN PetscErrorCode PetscSectionSetNumFields(PetscSection, PetscInt);
312: PETSC_EXTERN PetscErrorCode PetscSectionGetFieldName(PetscSection, PetscInt, const char *[]);
313: PETSC_EXTERN PetscErrorCode PetscSectionSetFieldName(PetscSection, PetscInt, const char []);
314: PETSC_EXTERN PetscErrorCode PetscSectionGetFieldComponents(PetscSection, PetscInt, PetscInt *);
315: PETSC_EXTERN PetscErrorCode PetscSectionSetFieldComponents(PetscSection, PetscInt, PetscInt);
316: PETSC_EXTERN PetscErrorCode PetscSectionGetChart(PetscSection, PetscInt *, PetscInt *);
317: PETSC_EXTERN PetscErrorCode PetscSectionSetChart(PetscSection, PetscInt, PetscInt);
318: PETSC_EXTERN PetscErrorCode PetscSectionGetPermutation(PetscSection, IS *);
319: PETSC_EXTERN PetscErrorCode PetscSectionSetPermutation(PetscSection, IS);
320: PETSC_EXTERN PetscErrorCode PetscSectionGetDof(PetscSection, PetscInt, PetscInt*);
321: PETSC_EXTERN PetscErrorCode PetscSectionSetDof(PetscSection, PetscInt, PetscInt);
322: PETSC_EXTERN PetscErrorCode PetscSectionAddDof(PetscSection, PetscInt, PetscInt);
323: PETSC_EXTERN PetscErrorCode PetscSectionGetFieldDof(PetscSection, PetscInt, PetscInt, PetscInt*);
324: PETSC_EXTERN PetscErrorCode PetscSectionSetFieldDof(PetscSection, PetscInt, PetscInt, PetscInt);
325: PETSC_EXTERN PetscErrorCode PetscSectionAddFieldDof(PetscSection, PetscInt, PetscInt, PetscInt);
326: PETSC_EXTERN PetscErrorCode PetscSectionHasConstraints(PetscSection, PetscBool *);
327: PETSC_EXTERN PetscErrorCode PetscSectionGetConstraintDof(PetscSection, PetscInt, PetscInt*);
328: PETSC_EXTERN PetscErrorCode PetscSectionSetConstraintDof(PetscSection, PetscInt, PetscInt);
329: PETSC_EXTERN PetscErrorCode PetscSectionAddConstraintDof(PetscSection, PetscInt, PetscInt);
330: PETSC_EXTERN PetscErrorCode PetscSectionGetFieldConstraintDof(PetscSection, PetscInt, PetscInt, PetscInt*);
331: PETSC_EXTERN PetscErrorCode PetscSectionSetFieldConstraintDof(PetscSection, PetscInt, PetscInt, PetscInt);
332: PETSC_EXTERN PetscErrorCode PetscSectionAddFieldConstraintDof(PetscSection, PetscInt, PetscInt, PetscInt);
333: PETSC_EXTERN PetscErrorCode PetscSectionGetConstraintIndices(PetscSection, PetscInt, const PetscInt**);
334: PETSC_EXTERN PetscErrorCode PetscSectionSetConstraintIndices(PetscSection, PetscInt, const PetscInt*);
335: PETSC_EXTERN PetscErrorCode PetscSectionGetFieldConstraintIndices(PetscSection, PetscInt, PetscInt, const PetscInt**);
336: PETSC_EXTERN PetscErrorCode PetscSectionSetFieldConstraintIndices(PetscSection, PetscInt, PetscInt, const PetscInt*);
337: PETSC_EXTERN PetscErrorCode PetscSectionSetUpBC(PetscSection);
338: PETSC_EXTERN PetscErrorCode PetscSectionSetUp(PetscSection);
339: PETSC_EXTERN PetscErrorCode PetscSectionGetMaxDof(PetscSection, PetscInt*);
340: PETSC_EXTERN PetscErrorCode PetscSectionGetStorageSize(PetscSection, PetscInt*);
341: PETSC_EXTERN PetscErrorCode PetscSectionGetConstrainedStorageSize(PetscSection, PetscInt*);
342: PETSC_EXTERN PetscErrorCode PetscSectionGetOffset(PetscSection, PetscInt, PetscInt*);
343: PETSC_EXTERN PetscErrorCode PetscSectionSetOffset(PetscSection, PetscInt, PetscInt);
344: PETSC_EXTERN PetscErrorCode PetscSectionGetFieldOffset(PetscSection, PetscInt, PetscInt, PetscInt*);
345: PETSC_EXTERN PetscErrorCode PetscSectionSetFieldOffset(PetscSection, PetscInt, PetscInt, PetscInt);
346: PETSC_EXTERN PetscErrorCode PetscSectionGetOffsetRange(PetscSection, PetscInt *, PetscInt *);
347: PETSC_EXTERN PetscErrorCode PetscSectionView(PetscSection, PetscViewer);
348: PETSC_STATIC_INLINE PetscErrorCode PetscSectionViewFromOptions(PetscSection A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);}
349: PETSC_EXTERN PetscErrorCode PetscSectionReset(PetscSection);
350: PETSC_EXTERN PetscErrorCode PetscSectionDestroy(PetscSection*);
351: PETSC_EXTERN PetscErrorCode PetscSectionCreateGlobalSection(PetscSection, PetscSF, PetscBool, PetscBool, PetscSection *);
352: PETSC_EXTERN PetscErrorCode PetscSectionCreateGlobalSectionCensored(PetscSection, PetscSF, PetscBool, PetscInt, const PetscInt [], PetscSection *);
353: PETSC_EXTERN PetscErrorCode PetscSectionCreateSubsection(PetscSection, PetscInt, PetscInt [], PetscSection *);
354: PETSC_EXTERN PetscErrorCode PetscSectionCreateSubmeshSection(PetscSection, IS, PetscSection *);
355: PETSC_EXTERN PetscErrorCode PetscSectionGetPointLayout(MPI_Comm, PetscSection, PetscLayout *);
356: PETSC_EXTERN PetscErrorCode PetscSectionGetValueLayout(MPI_Comm, PetscSection, PetscLayout *);
357: PETSC_EXTERN PetscErrorCode PetscSectionPermute(PetscSection, IS, PetscSection *);
358: PETSC_EXTERN PetscErrorCode PetscSectionGetField(PetscSection, PetscInt, PetscSection *);

360: PETSC_EXTERN PetscErrorCode PetscSectionSetClosureIndex(PetscSection, PetscObject, PetscSection, IS);
361: PETSC_EXTERN PetscErrorCode PetscSectionGetClosureIndex(PetscSection, PetscObject, PetscSection *, IS *);
362: PETSC_EXTERN PetscErrorCode PetscSectionSetClosurePermutation(PetscSection, PetscObject, IS);
363: PETSC_EXTERN PetscErrorCode PetscSectionGetClosurePermutation(PetscSection, PetscObject, IS *);
364: PETSC_EXTERN PetscErrorCode PetscSectionGetClosureInversePermutation(PetscSection, PetscObject, IS *);

366: PETSC_EXTERN PetscClassId PETSC_SECTION_SYM_CLASSID;

368: /*J
369:   PetscSectionSymType - String with the name of a PetscSectionSym type.

371:   Level: developer

373:   Notes: PetscSectionSym has no default implementation, but is used by DM in PetscSectionSymCreateLabel().

375: .seealso: PetscSectionSymSetType(), PetscSectionSym, PetscSectionSymCreate(), PetscSectionSymRegister()
376: J*/
377: typedef const char *PetscSectionSymType;

379: PETSC_EXTERN PetscFunctionList PetscSectionSymList;
380: PETSC_EXTERN PetscErrorCode PetscSectionSymSetType(PetscSectionSym, PetscSectionSymType);
381: PETSC_EXTERN PetscErrorCode PetscSectionSymGetType(PetscSectionSym, PetscSectionSymType*);
382: PETSC_EXTERN PetscErrorCode PetscSectionSymRegister(const char[],PetscErrorCode (*)(PetscSectionSym));

384: PETSC_EXTERN PetscErrorCode PetscSectionSymCreate(MPI_Comm, PetscSectionSym*);
385: PETSC_EXTERN PetscErrorCode PetscSectionSymDestroy(PetscSectionSym*);
386: PETSC_EXTERN PetscErrorCode PetscSectionSymView(PetscSectionSym,PetscViewer);

388: PETSC_EXTERN PetscErrorCode PetscSectionSetSym(PetscSection, PetscSectionSym);
389: PETSC_EXTERN PetscErrorCode PetscSectionGetSym(PetscSection, PetscSectionSym*);
390: PETSC_EXTERN PetscErrorCode PetscSectionSetFieldSym(PetscSection, PetscInt, PetscSectionSym);
391: PETSC_EXTERN PetscErrorCode PetscSectionGetFieldSym(PetscSection, PetscInt, PetscSectionSym*);

393: PETSC_EXTERN PetscErrorCode PetscSectionGetPointSyms(PetscSection, PetscInt, const PetscInt *, const PetscInt ***, const PetscScalar ***);
394: PETSC_EXTERN PetscErrorCode PetscSectionRestorePointSyms(PetscSection, PetscInt, const PetscInt *, const PetscInt ***, const PetscScalar ***);
395: PETSC_EXTERN PetscErrorCode PetscSectionGetFieldPointSyms(PetscSection, PetscInt, PetscInt, const PetscInt *, const PetscInt ***, const PetscScalar ***);
396: PETSC_EXTERN PetscErrorCode PetscSectionRestoreFieldPointSyms(PetscSection, PetscInt, PetscInt, const PetscInt *, const PetscInt ***, const PetscScalar ***);


399: /* PetscSF support */
400: PETSC_EXTERN PetscErrorCode PetscSFConvertPartition(PetscSF, PetscSection, IS, ISLocalToGlobalMapping *, PetscSF *);
401: PETSC_EXTERN PetscErrorCode PetscSFCreateRemoteOffsets(PetscSF, PetscSection, PetscSection, PetscInt **);
402: PETSC_EXTERN PetscErrorCode PetscSFDistributeSection(PetscSF, PetscSection, PetscInt **, PetscSection);
403: PETSC_EXTERN PetscErrorCode PetscSFCreateSectionSF(PetscSF, PetscSection, PetscInt [], PetscSection, PetscSF *);

405: #endif