Actual source code: petscdm.h
petsc-3.10.5 2019-03-28
1: /*
2: Objects to manage the interactions between the mesh data structures and the algebraic objects
3: */
6: #include <petscmat.h>
7: #include <petscdmtypes.h>
8: #include <petscfetypes.h>
9: #include <petscdstypes.h>
10: #include <petscdmlabel.h>
12: PETSC_EXTERN PetscErrorCode DMInitializePackage(void);
14: PETSC_EXTERN PetscClassId DM_CLASSID;
16: #define DMLOCATEPOINT_POINT_NOT_FOUND -367
18: /*J
19: DMType - String with the name of a PETSc DM
21: Level: beginner
23: .seealso: DMSetType(), DM
24: J*/
25: typedef const char* DMType;
26: #define DMDA "da"
27: #define DMCOMPOSITE "composite"
28: #define DMSLICED "sliced"
29: #define DMSHELL "shell"
30: #define DMPLEX "plex"
31: #define DMREDUNDANT "redundant"
32: #define DMPATCH "patch"
33: #define DMMOAB "moab"
34: #define DMNETWORK "network"
35: #define DMFOREST "forest"
36: #define DMP4EST "p4est"
37: #define DMP8EST "p8est"
38: #define DMSWARM "swarm"
40: PETSC_EXTERN const char *const DMBoundaryTypes[];
41: PETSC_EXTERN PetscFunctionList DMList;
42: PETSC_EXTERN PetscErrorCode DMCreate(MPI_Comm,DM*);
43: PETSC_EXTERN PetscErrorCode DMClone(DM,DM*);
44: PETSC_EXTERN PetscErrorCode DMSetType(DM, DMType);
45: PETSC_EXTERN PetscErrorCode DMGetType(DM, DMType *);
46: PETSC_EXTERN PetscErrorCode DMRegister(const char[],PetscErrorCode (*)(DM));
47: PETSC_EXTERN PetscErrorCode DMRegisterDestroy(void);
49: PETSC_EXTERN PetscErrorCode DMView(DM,PetscViewer);
50: PETSC_EXTERN PetscErrorCode DMLoad(DM,PetscViewer);
51: PETSC_EXTERN PetscErrorCode DMDestroy(DM*);
52: PETSC_EXTERN PetscErrorCode DMCreateGlobalVector(DM,Vec*);
53: PETSC_EXTERN PetscErrorCode DMCreateLocalVector(DM,Vec*);
54: PETSC_EXTERN PetscErrorCode DMGetLocalVector(DM,Vec *);
55: PETSC_EXTERN PetscErrorCode DMRestoreLocalVector(DM,Vec *);
56: PETSC_EXTERN PetscErrorCode DMGetGlobalVector(DM,Vec *);
57: PETSC_EXTERN PetscErrorCode DMRestoreGlobalVector(DM,Vec *);
58: PETSC_EXTERN PetscErrorCode DMClearGlobalVectors(DM);
59: PETSC_EXTERN PetscErrorCode DMClearLocalVectors(DM);
60: PETSC_EXTERN PetscErrorCode DMHasNamedGlobalVector(DM,const char*,PetscBool*);
61: PETSC_EXTERN PetscErrorCode DMGetNamedGlobalVector(DM,const char*,Vec*);
62: PETSC_EXTERN PetscErrorCode DMRestoreNamedGlobalVector(DM,const char*,Vec*);
63: PETSC_EXTERN PetscErrorCode DMHasNamedLocalVector(DM,const char*,PetscBool*);
64: PETSC_EXTERN PetscErrorCode DMGetNamedLocalVector(DM,const char*,Vec*);
65: PETSC_EXTERN PetscErrorCode DMRestoreNamedLocalVector(DM,const char*,Vec*);
66: PETSC_EXTERN PetscErrorCode DMGetLocalToGlobalMapping(DM,ISLocalToGlobalMapping*);
67: PETSC_EXTERN PetscErrorCode DMCreateFieldIS(DM,PetscInt*,char***,IS**);
68: PETSC_EXTERN PetscErrorCode DMGetBlockSize(DM,PetscInt*);
69: PETSC_EXTERN PetscErrorCode DMCreateColoring(DM,ISColoringType,ISColoring*);
70: PETSC_EXTERN PetscErrorCode DMCreateMatrix(DM,Mat*);
71: PETSC_EXTERN PetscErrorCode DMSetMatrixPreallocateOnly(DM,PetscBool);
72: PETSC_EXTERN PetscErrorCode DMSetMatrixStructureOnly(DM,PetscBool);
73: PETSC_EXTERN PetscErrorCode DMCreateInterpolation(DM,DM,Mat*,Vec*);
74: PETSC_EXTERN PetscErrorCode DMCreateRestriction(DM,DM,Mat*);
75: PETSC_EXTERN PetscErrorCode DMCreateInjection(DM,DM,Mat*);
76: PETSC_EXTERN PetscErrorCode DMCreateMassMatrix(DM,DM,Mat*);
77: PETSC_EXTERN PetscErrorCode DMGetWorkArray(DM,PetscInt,MPI_Datatype,void*);
78: PETSC_EXTERN PetscErrorCode DMRestoreWorkArray(DM,PetscInt,MPI_Datatype,void*);
79: PETSC_EXTERN PetscErrorCode DMRefine(DM,MPI_Comm,DM*);
80: PETSC_EXTERN PetscErrorCode DMCoarsen(DM,MPI_Comm,DM*);
81: PETSC_EXTERN PetscErrorCode DMGetCoarseDM(DM,DM*);
82: PETSC_EXTERN PetscErrorCode DMSetCoarseDM(DM,DM);
83: PETSC_EXTERN PetscErrorCode DMGetFineDM(DM,DM*);
84: PETSC_EXTERN PetscErrorCode DMSetFineDM(DM,DM);
85: PETSC_EXTERN PetscErrorCode DMRefineHierarchy(DM,PetscInt,DM[]);
86: PETSC_EXTERN PetscErrorCode DMCoarsenHierarchy(DM,PetscInt,DM[]);
87: PETSC_EXTERN PetscErrorCode DMCoarsenHookAdd(DM,PetscErrorCode (*)(DM,DM,void*),PetscErrorCode (*)(DM,Mat,Vec,Mat,DM,void*),void*);
88: PETSC_EXTERN PetscErrorCode DMCoarsenHookRemove(DM,PetscErrorCode (*)(DM,DM,void*),PetscErrorCode (*)(DM,Mat,Vec,Mat,DM,void*),void*);
89: PETSC_EXTERN PetscErrorCode DMRefineHookAdd(DM,PetscErrorCode (*)(DM,DM,void*),PetscErrorCode (*)(DM,Mat,DM,void*),void*);
90: PETSC_EXTERN PetscErrorCode DMRefineHookRemove(DM,PetscErrorCode (*)(DM,DM,void*),PetscErrorCode (*)(DM,Mat,DM,void*),void*);
91: PETSC_EXTERN PetscErrorCode DMRestrict(DM,Mat,Vec,Mat,DM);
92: PETSC_EXTERN PetscErrorCode DMInterpolate(DM,Mat,DM);
93: PETSC_EXTERN PetscErrorCode DMSetFromOptions(DM);
94: PETSC_STATIC_INLINE PetscErrorCode DMViewFromOptions(DM A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);}
96: PETSC_EXTERN PetscErrorCode DMAdaptLabel(DM,DMLabel,DM*);
97: PETSC_EXTERN PetscErrorCode DMAdaptMetric(DM, Vec, DMLabel, DM *);
99: PETSC_EXTERN PetscErrorCode DMSetUp(DM);
100: PETSC_EXTERN PetscErrorCode DMCreateInterpolationScale(DM,DM,Mat,Vec*);
101: PETSC_EXTERN PetscErrorCode DMCreateAggregates(DM,DM,Mat*);
102: PETSC_EXTERN PetscErrorCode DMGlobalToLocalHookAdd(DM,PetscErrorCode (*)(DM,Vec,InsertMode,Vec,void*),PetscErrorCode (*)(DM,Vec,InsertMode,Vec,void*),void*);
103: PETSC_EXTERN PetscErrorCode DMLocalToGlobalHookAdd(DM,PetscErrorCode (*)(DM,Vec,InsertMode,Vec,void*),PetscErrorCode (*)(DM,Vec,InsertMode,Vec,void*),void*);
104: PETSC_EXTERN PetscErrorCode DMGlobalToLocalBegin(DM,Vec,InsertMode,Vec);
105: PETSC_EXTERN PetscErrorCode DMGlobalToLocalEnd(DM,Vec,InsertMode,Vec);
106: PETSC_EXTERN PetscErrorCode DMLocalToGlobalBegin(DM,Vec,InsertMode,Vec);
107: PETSC_EXTERN PetscErrorCode DMLocalToGlobalEnd(DM,Vec,InsertMode,Vec);
108: PETSC_EXTERN PetscErrorCode DMLocalToLocalBegin(DM,Vec,InsertMode,Vec);
109: PETSC_EXTERN PetscErrorCode DMLocalToLocalEnd(DM,Vec,InsertMode,Vec);
110: PETSC_EXTERN PetscErrorCode DMConvert(DM,DMType,DM*);
112: /* Topology support */
113: PETSC_EXTERN PetscErrorCode DMGetDimension(DM,PetscInt*);
114: PETSC_EXTERN PetscErrorCode DMSetDimension(DM,PetscInt);
115: PETSC_EXTERN PetscErrorCode DMGetDimPoints(DM,PetscInt,PetscInt*,PetscInt*);
116: PETSC_EXTERN PetscErrorCode DMGetUseNatural(DM,PetscBool*);
117: PETSC_EXTERN PetscErrorCode DMSetUseNatural(DM,PetscBool);
119: /* Coordinate support */
120: PETSC_EXTERN PetscErrorCode DMGetCoordinateDM(DM,DM*);
121: PETSC_EXTERN PetscErrorCode DMSetCoordinateDM(DM,DM);
122: PETSC_EXTERN PetscErrorCode DMGetCoordinateDim(DM,PetscInt*);
123: PETSC_EXTERN PetscErrorCode DMSetCoordinateDim(DM,PetscInt);
124: PETSC_EXTERN PetscErrorCode DMGetCoordinateSection(DM,PetscSection*);
125: PETSC_EXTERN PetscErrorCode DMSetCoordinateSection(DM,PetscInt,PetscSection);
126: PETSC_EXTERN PetscErrorCode DMGetCoordinates(DM,Vec*);
127: PETSC_EXTERN PetscErrorCode DMSetCoordinates(DM,Vec);
128: PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocal(DM,Vec*);
129: PETSC_EXTERN PetscErrorCode DMSetCoordinatesLocal(DM,Vec);
130: PETSC_EXTERN PetscErrorCode DMLocatePoints(DM,Vec,DMPointLocationType,PetscSF*);
131: PETSC_EXTERN PetscErrorCode DMGetPeriodicity(DM,PetscBool*,const PetscReal**,const PetscReal**,const DMBoundaryType**);
132: PETSC_EXTERN PetscErrorCode DMSetPeriodicity(DM,PetscBool,const PetscReal[],const PetscReal[],const DMBoundaryType[]);
133: PETSC_EXTERN PetscErrorCode DMLocalizeCoordinate(DM, const PetscScalar[], PetscBool, PetscScalar[]);
134: PETSC_EXTERN PetscErrorCode DMLocalizeCoordinates(DM);
135: PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocalized(DM,PetscBool*);
136: PETSC_EXTERN PetscErrorCode DMGetNeighbors(DM,PetscInt*,const PetscMPIInt**);
137: PETSC_EXTERN PetscErrorCode DMGetCoordinateField(DM,DMField*);
138: PETSC_EXTERN PetscErrorCode DMSetCoordinateField(DM,DMField);
140: /* block hook interface */
141: PETSC_EXTERN PetscErrorCode DMSubDomainHookAdd(DM,PetscErrorCode (*)(DM,DM,void*),PetscErrorCode (*)(DM,VecScatter,VecScatter,DM,void*),void*);
142: PETSC_EXTERN PetscErrorCode DMSubDomainHookRemove(DM,PetscErrorCode (*)(DM,DM,void*),PetscErrorCode (*)(DM,VecScatter,VecScatter,DM,void*),void*);
143: PETSC_EXTERN PetscErrorCode DMSubDomainRestrict(DM,VecScatter,VecScatter,DM);
145: PETSC_EXTERN PetscErrorCode DMSetOptionsPrefix(DM,const char []);
146: PETSC_EXTERN PetscErrorCode DMAppendOptionsPrefix(DM,const char []);
147: PETSC_EXTERN PetscErrorCode DMGetOptionsPrefix(DM,const char*[]);
148: PETSC_EXTERN PetscErrorCode DMSetVecType(DM,VecType);
149: PETSC_EXTERN PetscErrorCode DMGetVecType(DM,VecType*);
150: PETSC_EXTERN PetscErrorCode DMSetMatType(DM,MatType);
151: PETSC_EXTERN PetscErrorCode DMGetMatType(DM,MatType*);
152: PETSC_EXTERN PetscErrorCode DMSetISColoringType(DM,ISColoringType);
153: PETSC_EXTERN PetscErrorCode DMGetISColoringType(DM,ISColoringType*);
154: PETSC_EXTERN PetscErrorCode DMSetApplicationContext(DM,void*);
155: PETSC_EXTERN PetscErrorCode DMSetApplicationContextDestroy(DM,PetscErrorCode (*)(void**));
156: PETSC_EXTERN PetscErrorCode DMGetApplicationContext(DM,void*);
157: PETSC_EXTERN PetscErrorCode DMSetVariableBounds(DM,PetscErrorCode (*)(DM,Vec,Vec));
158: PETSC_EXTERN PetscErrorCode DMHasVariableBounds(DM,PetscBool *);
159: PETSC_EXTERN PetscErrorCode DMHasColoring(DM,PetscBool *);
160: PETSC_EXTERN PetscErrorCode DMHasCreateRestriction(DM,PetscBool *);
161: PETSC_EXTERN PetscErrorCode DMHasCreateInjection(DM,PetscBool *);
162: PETSC_EXTERN PetscErrorCode DMComputeVariableBounds(DM,Vec,Vec);
164: PETSC_EXTERN PetscErrorCode DMCreateSubDM(DM, PetscInt, const PetscInt[], IS *, DM *);
165: PETSC_EXTERN PetscErrorCode DMCreateSuperDM(DM[], PetscInt, IS **, DM *);
166: PETSC_EXTERN PetscErrorCode DMCreateFieldDecomposition(DM,PetscInt*,char***,IS**,DM**);
167: PETSC_EXTERN PetscErrorCode DMCreateDomainDecomposition(DM,PetscInt*,char***,IS**,IS**,DM**);
168: PETSC_EXTERN PetscErrorCode DMCreateDomainDecompositionScatters(DM,PetscInt,DM*,VecScatter**,VecScatter**,VecScatter**);
170: PETSC_EXTERN PetscErrorCode DMGetRefineLevel(DM,PetscInt*);
171: PETSC_EXTERN PetscErrorCode DMSetRefineLevel(DM,PetscInt);
172: PETSC_EXTERN PetscErrorCode DMGetCoarsenLevel(DM,PetscInt*);
173: PETSC_EXTERN PetscErrorCode DMFinalizePackage(void);
175: PETSC_EXTERN PetscErrorCode VecGetDM(Vec, DM*);
176: PETSC_EXTERN PetscErrorCode VecSetDM(Vec, DM);
177: PETSC_EXTERN PetscErrorCode MatGetDM(Mat, DM*);
178: PETSC_EXTERN PetscErrorCode MatSetDM(Mat, DM);
179: PETSC_EXTERN PetscErrorCode MatFDColoringUseDM(Mat,MatFDColoring);
181: typedef struct NLF_DAAD* NLF;
183: #define DM_FILE_CLASSID 1211221
185: /* FEM support */
186: PETSC_EXTERN PetscErrorCode DMPrintCellVector(PetscInt, const char [], PetscInt, const PetscScalar []);
187: PETSC_EXTERN PetscErrorCode DMPrintCellMatrix(PetscInt, const char [], PetscInt, PetscInt, const PetscScalar []);
188: PETSC_EXTERN PetscErrorCode DMPrintLocalVec(DM, const char [], PetscReal, Vec);
190: PETSC_EXTERN PetscErrorCode DMGetSection(DM, PetscSection *);
191: PETSC_EXTERN PetscErrorCode DMSetSection(DM, PetscSection);
192: PETSC_EXTERN PetscErrorCode DMGetDefaultConstraints(DM, PetscSection *, Mat *);
193: PETSC_EXTERN PetscErrorCode DMSetDefaultConstraints(DM, PetscSection, Mat);
194: PETSC_EXTERN PetscErrorCode DMGetGlobalSection(DM, PetscSection *);
195: PETSC_EXTERN PetscErrorCode DMSetGlobalSection(DM, PetscSection);
196: PETSC_EXTERN PetscErrorCode DMGetDefaultSF(DM, PetscSF *);
197: PETSC_EXTERN PetscErrorCode DMSetDefaultSF(DM, PetscSF);
198: PETSC_EXTERN PetscErrorCode DMCreateDefaultSF(DM, PetscSection, PetscSection);
199: PETSC_EXTERN PetscErrorCode DMGetPointSF(DM, PetscSF *);
200: PETSC_EXTERN PetscErrorCode DMSetPointSF(DM, PetscSF);
202: PETSC_STATIC_INLINE PetscErrorCode DMGetDefaultSection(DM dm, PetscSection *s) {return DMGetSection(dm,s);}
203: PETSC_STATIC_INLINE PetscErrorCode DMSetDefaultSection(DM dm, PetscSection s) {return DMSetSection(dm,s);}
204: PETSC_STATIC_INLINE PetscErrorCode DMGetDefaultGlobalSection(DM dm, PetscSection *s) {return DMGetGlobalSection(dm,s);}
205: PETSC_STATIC_INLINE PetscErrorCode DMSetDefaultGlobalSection(DM dm, PetscSection s) {return DMSetGlobalSection(dm,s);}
208: PETSC_EXTERN PetscErrorCode DMGetOutputDM(DM, DM *);
209: PETSC_EXTERN PetscErrorCode DMGetOutputSequenceNumber(DM, PetscInt *, PetscReal *);
210: PETSC_EXTERN PetscErrorCode DMSetOutputSequenceNumber(DM, PetscInt, PetscReal);
211: PETSC_EXTERN PetscErrorCode DMOutputSequenceLoad(DM, PetscViewer, const char *, PetscInt, PetscReal *);
213: PETSC_EXTERN PetscErrorCode DMGetDS(DM, PetscDS *);
214: PETSC_EXTERN PetscErrorCode DMSetDS(DM, PetscDS);
215: PETSC_EXTERN PetscErrorCode DMGetNumFields(DM, PetscInt *);
216: PETSC_EXTERN PetscErrorCode DMSetNumFields(DM, PetscInt);
217: PETSC_EXTERN PetscErrorCode DMGetField(DM, PetscInt, PetscObject *);
218: PETSC_EXTERN PetscErrorCode DMSetField(DM, PetscInt, PetscObject);
220: struct _DMInterpolationInfo {
221: MPI_Comm comm;
222: PetscInt dim; /*1 The spatial dimension of points */
223: PetscInt nInput; /* The number of input points */
224: PetscReal *points; /* The input point coordinates */
225: PetscInt *cells; /* The cell containing each point */
226: PetscInt n; /* The number of local points */
227: Vec coords; /* The point coordinates */
228: PetscInt dof; /* The number of components to interpolate */
229: };
230: typedef struct _DMInterpolationInfo *DMInterpolationInfo;
232: PETSC_EXTERN PetscErrorCode DMInterpolationCreate(MPI_Comm, DMInterpolationInfo *);
233: PETSC_EXTERN PetscErrorCode DMInterpolationSetDim(DMInterpolationInfo, PetscInt);
234: PETSC_EXTERN PetscErrorCode DMInterpolationGetDim(DMInterpolationInfo, PetscInt *);
235: PETSC_EXTERN PetscErrorCode DMInterpolationSetDof(DMInterpolationInfo, PetscInt);
236: PETSC_EXTERN PetscErrorCode DMInterpolationGetDof(DMInterpolationInfo, PetscInt *);
237: PETSC_EXTERN PetscErrorCode DMInterpolationAddPoints(DMInterpolationInfo, PetscInt, PetscReal[]);
238: PETSC_EXTERN PetscErrorCode DMInterpolationSetUp(DMInterpolationInfo, DM, PetscBool);
239: PETSC_EXTERN PetscErrorCode DMInterpolationGetCoordinates(DMInterpolationInfo, Vec *);
240: PETSC_EXTERN PetscErrorCode DMInterpolationGetVector(DMInterpolationInfo, Vec *);
241: PETSC_EXTERN PetscErrorCode DMInterpolationRestoreVector(DMInterpolationInfo, Vec *);
242: PETSC_EXTERN PetscErrorCode DMInterpolationEvaluate(DMInterpolationInfo, DM, Vec, Vec);
243: PETSC_EXTERN PetscErrorCode DMInterpolationDestroy(DMInterpolationInfo *);
245: PETSC_EXTERN PetscErrorCode DMCreateLabel(DM, const char []);
246: PETSC_EXTERN PetscErrorCode DMGetLabelValue(DM, const char[], PetscInt, PetscInt *);
247: PETSC_EXTERN PetscErrorCode DMSetLabelValue(DM, const char[], PetscInt, PetscInt);
248: PETSC_EXTERN PetscErrorCode DMClearLabelValue(DM, const char[], PetscInt, PetscInt);
249: PETSC_EXTERN PetscErrorCode DMGetLabelSize(DM, const char[], PetscInt *);
250: PETSC_EXTERN PetscErrorCode DMGetLabelIdIS(DM, const char[], IS *);
251: PETSC_EXTERN PetscErrorCode DMGetStratumSize(DM, const char [], PetscInt, PetscInt *);
252: PETSC_EXTERN PetscErrorCode DMGetStratumIS(DM, const char [], PetscInt, IS *);
253: PETSC_EXTERN PetscErrorCode DMSetStratumIS(DM, const char [], PetscInt, IS);
254: PETSC_EXTERN PetscErrorCode DMClearLabelStratum(DM, const char[], PetscInt);
255: PETSC_EXTERN PetscErrorCode DMGetLabelOutput(DM, const char[], PetscBool *);
256: PETSC_EXTERN PetscErrorCode DMSetLabelOutput(DM, const char[], PetscBool);
258: PETSC_EXTERN PetscErrorCode DMGetNumLabels(DM, PetscInt *);
259: PETSC_EXTERN PetscErrorCode DMGetLabelName(DM, PetscInt, const char **);
260: PETSC_EXTERN PetscErrorCode DMHasLabel(DM, const char [], PetscBool *);
261: PETSC_EXTERN PetscErrorCode DMGetLabel(DM, const char *, DMLabel *);
262: PETSC_EXTERN PetscErrorCode DMGetLabelByNum(DM, PetscInt, DMLabel *);
263: PETSC_EXTERN PetscErrorCode DMAddLabel(DM, DMLabel);
264: PETSC_EXTERN PetscErrorCode DMRemoveLabel(DM, const char [], DMLabel *);
265: PETSC_EXTERN PetscErrorCode DMCopyLabels(DM, DM);
267: PETSC_EXTERN PetscErrorCode DMAddBoundary(DM, DMBoundaryConditionType, const char[], const char[], PetscInt, PetscInt, const PetscInt *, void (*)(void), PetscInt, const PetscInt *, void *);
268: PETSC_EXTERN PetscErrorCode DMGetNumBoundary(DM, PetscInt *);
269: PETSC_EXTERN PetscErrorCode DMGetBoundary(DM, PetscInt, DMBoundaryConditionType *, const char **, const char **, PetscInt *, PetscInt *, const PetscInt **, void (**)(void), PetscInt *, const PetscInt **, void **);
270: PETSC_EXTERN PetscErrorCode DMIsBoundaryPoint(DM, PetscInt, PetscBool *);
271: PETSC_EXTERN PetscErrorCode DMCopyBoundary(DM, DM);
273: PETSC_EXTERN PetscErrorCode DMProjectFunction(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void**,InsertMode,Vec);
274: PETSC_EXTERN PetscErrorCode DMProjectFunctionLocal(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void**,InsertMode,Vec);
275: PETSC_EXTERN PetscErrorCode DMProjectFunctionLabel(DM, PetscReal, DMLabel, PetscInt, const PetscInt[], PetscInt, const PetscInt[], PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **, InsertMode, Vec);
276: PETSC_EXTERN PetscErrorCode DMProjectFunctionLabelLocal(DM,PetscReal,DMLabel,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,InsertMode,Vec);
277: PETSC_EXTERN PetscErrorCode DMProjectFieldLocal(DM,PetscReal,Vec,void (**)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),InsertMode,Vec);
278: PETSC_EXTERN PetscErrorCode DMProjectFieldLabelLocal(DM,PetscReal,DMLabel,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Vec,void (**)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),InsertMode,Vec);
279: PETSC_EXTERN PetscErrorCode DMComputeL2Diff(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,Vec,PetscReal *);
280: PETSC_EXTERN PetscErrorCode DMComputeL2GradientDiff(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal [], const PetscReal [], PetscInt, PetscScalar *, void *), void **, Vec, const PetscReal [], PetscReal *);
281: PETSC_EXTERN PetscErrorCode DMComputeL2FieldDiff(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,Vec,PetscReal *);
283: PETSC_EXTERN PetscErrorCode DMGetNullSpaceConstructor(DM, PetscInt, PetscErrorCode (**)(DM, PetscInt, MatNullSpace *));
284: PETSC_EXTERN PetscErrorCode DMSetNullSpaceConstructor(DM, PetscInt, PetscErrorCode (*)(DM, PetscInt, MatNullSpace *));
286: PETSC_EXTERN PetscErrorCode DMGetCompatibility(DM,DM,PetscBool*,PetscBool*);
287: #endif