Actual source code: petscdm.h

petsc-3.8.1 2017-11-04
Report Typos and Errors
  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 DMCARTESIAN "cartesian"
 32: #define DMREDUNDANT "redundant"
 33: #define DMPATCH     "patch"
 34: #define DMMOAB      "moab"
 35: #define DMNETWORK   "network"
 36: #define DMFOREST    "forest"
 37: #define DMP4EST     "p4est"
 38: #define DMP8EST     "p8est"
 39: #define DMSWARM     "swarm"

 41: PETSC_EXTERN const char *const DMBoundaryTypes[];
 42: PETSC_EXTERN PetscFunctionList DMList;
 43: PETSC_EXTERN PetscErrorCode DMCreate(MPI_Comm,DM*);
 44: PETSC_EXTERN PetscErrorCode DMClone(DM,DM*);
 45: PETSC_EXTERN PetscErrorCode DMSetType(DM, DMType);
 46: PETSC_EXTERN PetscErrorCode DMGetType(DM, DMType *);
 47: PETSC_EXTERN PetscErrorCode DMRegister(const char[],PetscErrorCode (*)(DM));
 48: PETSC_EXTERN PetscErrorCode DMRegisterDestroy(void);

 50: PETSC_EXTERN PetscErrorCode DMView(DM,PetscViewer);
 51: PETSC_EXTERN PetscErrorCode DMLoad(DM,PetscViewer);
 52: PETSC_EXTERN PetscErrorCode DMDestroy(DM*);
 53: PETSC_EXTERN PetscErrorCode DMCreateGlobalVector(DM,Vec*);
 54: PETSC_EXTERN PetscErrorCode DMCreateLocalVector(DM,Vec*);
 55: PETSC_EXTERN PetscErrorCode DMGetLocalVector(DM,Vec *);
 56: PETSC_EXTERN PetscErrorCode DMRestoreLocalVector(DM,Vec *);
 57: PETSC_EXTERN PetscErrorCode DMGetGlobalVector(DM,Vec *);
 58: PETSC_EXTERN PetscErrorCode DMRestoreGlobalVector(DM,Vec *);
 59: PETSC_EXTERN PetscErrorCode DMClearGlobalVectors(DM);
 60: PETSC_EXTERN PetscErrorCode DMClearLocalVectors(DM);
 61: PETSC_EXTERN PetscErrorCode DMHasNamedGlobalVector(DM,const char*,PetscBool*);
 62: PETSC_EXTERN PetscErrorCode DMGetNamedGlobalVector(DM,const char*,Vec*);
 63: PETSC_EXTERN PetscErrorCode DMRestoreNamedGlobalVector(DM,const char*,Vec*);
 64: PETSC_EXTERN PetscErrorCode DMHasNamedLocalVector(DM,const char*,PetscBool*);
 65: PETSC_EXTERN PetscErrorCode DMGetNamedLocalVector(DM,const char*,Vec*);
 66: PETSC_EXTERN PetscErrorCode DMRestoreNamedLocalVector(DM,const char*,Vec*);
 67: PETSC_EXTERN PetscErrorCode DMGetLocalToGlobalMapping(DM,ISLocalToGlobalMapping*);
 68: PETSC_EXTERN PetscErrorCode DMCreateFieldIS(DM,PetscInt*,char***,IS**);
 69: PETSC_EXTERN PetscErrorCode DMGetBlockSize(DM,PetscInt*);
 70: PETSC_EXTERN PetscErrorCode DMCreateColoring(DM,ISColoringType,ISColoring*);
 71: PETSC_EXTERN PetscErrorCode DMCreateMatrix(DM,Mat*);
 72: PETSC_EXTERN PetscErrorCode DMSetMatrixPreallocateOnly(DM,PetscBool);
 73: PETSC_EXTERN PetscErrorCode DMSetMatrixStructureOnly(DM,PetscBool);
 74: PETSC_EXTERN PetscErrorCode DMCreateInterpolation(DM,DM,Mat*,Vec*);
 75: PETSC_EXTERN PetscErrorCode DMCreateRestriction(DM,DM,Mat*);
 76: PETSC_EXTERN PetscErrorCode DMCreateInjection(DM,DM,Mat*);
 77: PETSC_EXTERN PetscErrorCode DMGetWorkArray(DM,PetscInt,PetscDataType,void*);
 78: PETSC_EXTERN PetscErrorCode DMRestoreWorkArray(DM,PetscInt,PetscDataType,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 DMRefineHookAdd(DM,PetscErrorCode (*)(DM,DM,void*),PetscErrorCode (*)(DM,Mat,DM,void*),void*);
 89: PETSC_EXTERN PetscErrorCode DMRestrict(DM,Mat,Vec,Mat,DM);
 90: PETSC_EXTERN PetscErrorCode DMInterpolate(DM,Mat,DM);
 91: PETSC_EXTERN PetscErrorCode DMSetFromOptions(DM);
 92: PETSC_STATIC_INLINE PetscErrorCode DMViewFromOptions(DM A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);}

 94: typedef enum {DM_ADAPT_DETERMINE = PETSC_DETERMINE, DM_ADAPT_KEEP = 0, DM_ADAPT_REFINE, DM_ADAPT_COARSEN, DM_ADAPT_RESERVED_COUNT} DMAdaptFlag;
 95: PETSC_EXTERN PetscErrorCode DMAdaptLabel(DM,DMLabel,DM*);
 96: PETSC_EXTERN PetscErrorCode DMAdaptMetric(DM, Vec, DMLabel, DM *);

 98: PETSC_EXTERN PetscErrorCode DMSetUp(DM);
 99: PETSC_EXTERN PetscErrorCode DMCreateInterpolationScale(DM,DM,Mat,Vec*);
100: PETSC_EXTERN PetscErrorCode DMCreateAggregates(DM,DM,Mat*);
101: PETSC_EXTERN PetscErrorCode DMGlobalToLocalHookAdd(DM,PetscErrorCode (*)(DM,Vec,InsertMode,Vec,void*),PetscErrorCode (*)(DM,Vec,InsertMode,Vec,void*),void*);
102: PETSC_EXTERN PetscErrorCode DMLocalToGlobalHookAdd(DM,PetscErrorCode (*)(DM,Vec,InsertMode,Vec,void*),PetscErrorCode (*)(DM,Vec,InsertMode,Vec,void*),void*);
103: PETSC_EXTERN PetscErrorCode DMGlobalToLocalBegin(DM,Vec,InsertMode,Vec);
104: PETSC_EXTERN PetscErrorCode DMGlobalToLocalEnd(DM,Vec,InsertMode,Vec);
105: PETSC_EXTERN PetscErrorCode DMLocalToGlobalBegin(DM,Vec,InsertMode,Vec);
106: PETSC_EXTERN PetscErrorCode DMLocalToGlobalEnd(DM,Vec,InsertMode,Vec);
107: PETSC_EXTERN PetscErrorCode DMLocalToLocalBegin(DM,Vec,InsertMode,Vec);
108: PETSC_EXTERN PetscErrorCode DMLocalToLocalEnd(DM,Vec,InsertMode,Vec);
109: PETSC_EXTERN PetscErrorCode DMConvert(DM,DMType,DM*);

111: /* Topology support */
112: PETSC_EXTERN PetscErrorCode DMGetDimension(DM,PetscInt*);
113: PETSC_EXTERN PetscErrorCode DMSetDimension(DM,PetscInt);
114: PETSC_EXTERN PetscErrorCode DMGetDimPoints(DM,PetscInt,PetscInt*,PetscInt*);
115: PETSC_EXTERN PetscErrorCode DMGetUseNatural(DM,PetscBool*);
116: PETSC_EXTERN PetscErrorCode DMSetUseNatural(DM,PetscBool);

118: /* Coordinate support */
119: PETSC_EXTERN PetscErrorCode DMGetCoordinateDM(DM,DM*);
120: PETSC_EXTERN PetscErrorCode DMSetCoordinateDM(DM,DM);
121: PETSC_EXTERN PetscErrorCode DMGetCoordinateDim(DM,PetscInt*);
122: PETSC_EXTERN PetscErrorCode DMSetCoordinateDim(DM,PetscInt);
123: PETSC_EXTERN PetscErrorCode DMGetCoordinateSection(DM,PetscSection*);
124: PETSC_EXTERN PetscErrorCode DMSetCoordinateSection(DM,PetscInt,PetscSection);
125: PETSC_EXTERN PetscErrorCode DMGetCoordinates(DM,Vec*);
126: PETSC_EXTERN PetscErrorCode DMSetCoordinates(DM,Vec);
127: PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocal(DM,Vec*);
128: PETSC_EXTERN PetscErrorCode DMSetCoordinatesLocal(DM,Vec);
129: PETSC_EXTERN PetscErrorCode DMLocatePoints(DM,Vec,DMPointLocationType,PetscSF*);
130: PETSC_EXTERN PetscErrorCode DMGetPeriodicity(DM,PetscBool*,const PetscReal**,const PetscReal**,const DMBoundaryType**);
131: PETSC_EXTERN PetscErrorCode DMSetPeriodicity(DM,PetscBool,const PetscReal[],const PetscReal[],const DMBoundaryType[]);
132: PETSC_EXTERN PetscErrorCode DMLocalizeCoordinate(DM, const PetscScalar[], PetscBool, PetscScalar[]);
133: PETSC_EXTERN PetscErrorCode DMLocalizeCoordinates(DM);
134: PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocalized(DM,PetscBool*);
135: PETSC_EXTERN PetscErrorCode DMGetNeighbors(DM,PetscInt*,const PetscMPIInt**);

137: /* block hook interface */
138: PETSC_EXTERN PetscErrorCode DMSubDomainHookAdd(DM,PetscErrorCode (*)(DM,DM,void*),PetscErrorCode (*)(DM,VecScatter,VecScatter,DM,void*),void*);
139: PETSC_EXTERN PetscErrorCode DMSubDomainRestrict(DM,VecScatter,VecScatter,DM);

141: PETSC_EXTERN PetscErrorCode DMSetOptionsPrefix(DM,const char []);
142: PETSC_EXTERN PetscErrorCode DMAppendOptionsPrefix(DM,const char []);
143: PETSC_EXTERN PetscErrorCode DMGetOptionsPrefix(DM,const char*[]);
144: PETSC_EXTERN PetscErrorCode DMSetVecType(DM,VecType);
145: PETSC_EXTERN PetscErrorCode DMGetVecType(DM,VecType*);
146: PETSC_EXTERN PetscErrorCode DMSetMatType(DM,MatType);
147: PETSC_EXTERN PetscErrorCode DMGetMatType(DM,MatType*);
148: PETSC_EXTERN PetscErrorCode DMSetApplicationContext(DM,void*);
149: PETSC_EXTERN PetscErrorCode DMSetApplicationContextDestroy(DM,PetscErrorCode (*)(void**));
150: PETSC_EXTERN PetscErrorCode DMGetApplicationContext(DM,void*);
151: PETSC_EXTERN PetscErrorCode DMSetVariableBounds(DM,PetscErrorCode (*)(DM,Vec,Vec));
152: PETSC_EXTERN PetscErrorCode DMHasVariableBounds(DM,PetscBool *);
153: PETSC_EXTERN PetscErrorCode DMHasColoring(DM,PetscBool *);
154: PETSC_EXTERN PetscErrorCode DMHasCreateRestriction(DM,PetscBool *);
155: PETSC_EXTERN PetscErrorCode DMComputeVariableBounds(DM,Vec,Vec);

157: PETSC_EXTERN PetscErrorCode DMCreateSubDM(DM, PetscInt, PetscInt[], IS *, DM *);
158: PETSC_EXTERN PetscErrorCode DMCreateFieldDecomposition(DM,PetscInt*,char***,IS**,DM**);
159: PETSC_EXTERN PetscErrorCode DMCreateDomainDecomposition(DM,PetscInt*,char***,IS**,IS**,DM**);
160: PETSC_EXTERN PetscErrorCode DMCreateDomainDecompositionScatters(DM,PetscInt,DM*,VecScatter**,VecScatter**,VecScatter**);

162: PETSC_EXTERN PetscErrorCode DMGetRefineLevel(DM,PetscInt*);
163: PETSC_EXTERN PetscErrorCode DMSetRefineLevel(DM,PetscInt);
164: PETSC_EXTERN PetscErrorCode DMGetCoarsenLevel(DM,PetscInt*);
165: PETSC_EXTERN PetscErrorCode DMFinalizePackage(void);

167: PETSC_EXTERN PetscErrorCode VecGetDM(Vec, DM*);
168: PETSC_EXTERN PetscErrorCode VecSetDM(Vec, DM);
169: PETSC_EXTERN PetscErrorCode MatGetDM(Mat, DM*);
170: PETSC_EXTERN PetscErrorCode MatSetDM(Mat, DM);
171: PETSC_EXTERN PetscErrorCode MatFDColoringUseDM(Mat,MatFDColoring);

173: typedef struct NLF_DAAD* NLF;

175: #define DM_FILE_CLASSID 1211221

177: /* FEM support */
178: PETSC_EXTERN PetscErrorCode DMPrintCellVector(PetscInt, const char [], PetscInt, const PetscScalar []);
179: PETSC_EXTERN PetscErrorCode DMPrintCellMatrix(PetscInt, const char [], PetscInt, PetscInt, const PetscScalar []);
180: PETSC_EXTERN PetscErrorCode DMPrintLocalVec(DM, const char [], PetscReal, Vec);

182: PETSC_EXTERN PetscErrorCode DMGetDefaultSection(DM, PetscSection *);
183: PETSC_EXTERN PetscErrorCode DMSetDefaultSection(DM, PetscSection);
184: PETSC_EXTERN PetscErrorCode DMGetDefaultConstraints(DM, PetscSection *, Mat *);
185: PETSC_EXTERN PetscErrorCode DMSetDefaultConstraints(DM, PetscSection, Mat);
186: PETSC_EXTERN PetscErrorCode DMGetDefaultGlobalSection(DM, PetscSection *);
187: PETSC_EXTERN PetscErrorCode DMSetDefaultGlobalSection(DM, PetscSection);
188: PETSC_EXTERN PetscErrorCode DMGetDefaultSF(DM, PetscSF *);
189: PETSC_EXTERN PetscErrorCode DMSetDefaultSF(DM, PetscSF);
190: PETSC_EXTERN PetscErrorCode DMCreateDefaultSF(DM, PetscSection, PetscSection);
191: PETSC_EXTERN PetscErrorCode DMGetPointSF(DM, PetscSF *);
192: PETSC_EXTERN PetscErrorCode DMSetPointSF(DM, PetscSF);

194: PETSC_EXTERN PetscErrorCode DMGetOutputDM(DM, DM *);
195: PETSC_EXTERN PetscErrorCode DMGetOutputSequenceNumber(DM, PetscInt *, PetscReal *);
196: PETSC_EXTERN PetscErrorCode DMSetOutputSequenceNumber(DM, PetscInt, PetscReal);
197: PETSC_EXTERN PetscErrorCode DMOutputSequenceLoad(DM, PetscViewer, const char *, PetscInt, PetscReal *);

199: PETSC_EXTERN PetscErrorCode DMGetDS(DM, PetscDS *);
200: PETSC_EXTERN PetscErrorCode DMSetDS(DM, PetscDS);
201: PETSC_EXTERN PetscErrorCode DMGetNumFields(DM, PetscInt *);
202: PETSC_EXTERN PetscErrorCode DMSetNumFields(DM, PetscInt);
203: PETSC_EXTERN PetscErrorCode DMGetField(DM, PetscInt, PetscObject *);
204: PETSC_EXTERN PetscErrorCode DMSetField(DM, PetscInt, PetscObject);

206: /*E
207:   PetscUnit - The seven fundamental SI units

209:   Level: beginner

211: .seealso: DMPlexGetScale(), DMPlexSetScale()
212: E*/
213: typedef enum {PETSC_UNIT_LENGTH, PETSC_UNIT_MASS, PETSC_UNIT_TIME, PETSC_UNIT_CURRENT, PETSC_UNIT_TEMPERATURE, PETSC_UNIT_AMOUNT, PETSC_UNIT_LUMINOSITY, NUM_PETSC_UNITS} PetscUnit;

215: struct _DMInterpolationInfo {
216:   MPI_Comm   comm;
217:   PetscInt   dim;    /*1 The spatial dimension of points */
218:   PetscInt   nInput; /* The number of input points */
219:   PetscReal *points; /* The input point coordinates */
220:   PetscInt  *cells;  /* The cell containing each point */
221:   PetscInt   n;      /* The number of local points */
222:   Vec        coords; /* The point coordinates */
223:   PetscInt   dof;    /* The number of components to interpolate */
224: };
225: typedef struct _DMInterpolationInfo *DMInterpolationInfo;

227: PETSC_EXTERN PetscErrorCode DMInterpolationCreate(MPI_Comm, DMInterpolationInfo *);
228: PETSC_EXTERN PetscErrorCode DMInterpolationSetDim(DMInterpolationInfo, PetscInt);
229: PETSC_EXTERN PetscErrorCode DMInterpolationGetDim(DMInterpolationInfo, PetscInt *);
230: PETSC_EXTERN PetscErrorCode DMInterpolationSetDof(DMInterpolationInfo, PetscInt);
231: PETSC_EXTERN PetscErrorCode DMInterpolationGetDof(DMInterpolationInfo, PetscInt *);
232: PETSC_EXTERN PetscErrorCode DMInterpolationAddPoints(DMInterpolationInfo, PetscInt, PetscReal[]);
233: PETSC_EXTERN PetscErrorCode DMInterpolationSetUp(DMInterpolationInfo, DM, PetscBool);
234: PETSC_EXTERN PetscErrorCode DMInterpolationGetCoordinates(DMInterpolationInfo, Vec *);
235: PETSC_EXTERN PetscErrorCode DMInterpolationGetVector(DMInterpolationInfo, Vec *);
236: PETSC_EXTERN PetscErrorCode DMInterpolationRestoreVector(DMInterpolationInfo, Vec *);
237: PETSC_EXTERN PetscErrorCode DMInterpolationEvaluate(DMInterpolationInfo, DM, Vec, Vec);
238: PETSC_EXTERN PetscErrorCode DMInterpolationDestroy(DMInterpolationInfo *);

240: PETSC_EXTERN PetscErrorCode DMCreateLabel(DM, const char []);
241: PETSC_EXTERN PetscErrorCode DMGetLabelValue(DM, const char[], PetscInt, PetscInt *);
242: PETSC_EXTERN PetscErrorCode DMSetLabelValue(DM, const char[], PetscInt, PetscInt);
243: PETSC_EXTERN PetscErrorCode DMClearLabelValue(DM, const char[], PetscInt, PetscInt);
244: PETSC_EXTERN PetscErrorCode DMGetLabelSize(DM, const char[], PetscInt *);
245: PETSC_EXTERN PetscErrorCode DMGetLabelIdIS(DM, const char[], IS *);
246: PETSC_EXTERN PetscErrorCode DMGetStratumSize(DM, const char [], PetscInt, PetscInt *);
247: PETSC_EXTERN PetscErrorCode DMGetStratumIS(DM, const char [], PetscInt, IS *);
248: PETSC_EXTERN PetscErrorCode DMSetStratumIS(DM, const char [], PetscInt, IS);
249: PETSC_EXTERN PetscErrorCode DMClearLabelStratum(DM, const char[], PetscInt);
250: PETSC_EXTERN PetscErrorCode DMGetLabelOutput(DM, const char[], PetscBool *);
251: PETSC_EXTERN PetscErrorCode DMSetLabelOutput(DM, const char[], PetscBool);

253: PETSC_EXTERN PetscErrorCode DMGetNumLabels(DM, PetscInt *);
254: PETSC_EXTERN PetscErrorCode DMGetLabelName(DM, PetscInt, const char **);
255: PETSC_EXTERN PetscErrorCode DMHasLabel(DM, const char [], PetscBool *);
256: PETSC_EXTERN PetscErrorCode DMGetLabel(DM, const char *, DMLabel *);
257: PETSC_EXTERN PetscErrorCode DMGetLabelByNum(DM, PetscInt, DMLabel *);
258: PETSC_EXTERN PetscErrorCode DMAddLabel(DM, DMLabel);
259: PETSC_EXTERN PetscErrorCode DMRemoveLabel(DM, const char [], DMLabel *);
260: PETSC_EXTERN PetscErrorCode DMCopyLabels(DM, DM);

262: /*E
263:   DMBoundaryConditionType - indicates what type of boundary condition is to be imposed

265:   Note: This flag indicates the type of function which will define the condition:
266: $ DM_BC_ESSENTIAL       - A Dirichlet condition using a function of the coordinates
267: $ DM_BC_ESSENTIAL_FIELD - A Dirichlet condition using a function of the coordinates and auxiliary field data
268: $ DM_BC_NATURAL         - A Neumann condition using a function of the coordinates
269: $ DM_BC_NATURAL_FIELD   - A Dirichlet condition using a function of the coordinates and auxiliary field data
270: $ DM_BC_NATURAL_RIEMANN - A flux condition which determines the state in ghost cells
271: The user can check whether a boundary condition is essential using (type & DM_BC_ESSENTIAL), and similarly for
272: natural conditions (type & DM_BC_NATURAL)

274:   Level: beginner

276: .seealso: DMAddBoundary(), DMGetBoundary()
277: E*/
278: typedef enum {DM_BC_ESSENTIAL = 1, DM_BC_ESSENTIAL_FIELD = 5, DM_BC_NATURAL = 2, DM_BC_NATURAL_FIELD = 6, DM_BC_NATURAL_RIEMANN = 10} DMBoundaryConditionType;
279: PETSC_EXTERN PetscErrorCode DMAddBoundary(DM, DMBoundaryConditionType, const char[], const char[], PetscInt, PetscInt, const PetscInt *, void (*)(void), PetscInt, const PetscInt *, void *);
280: PETSC_EXTERN PetscErrorCode DMGetNumBoundary(DM, PetscInt *);
281: PETSC_EXTERN PetscErrorCode DMGetBoundary(DM, PetscInt, DMBoundaryConditionType *, const char **, const char **, PetscInt *, PetscInt *, const PetscInt **, void (**)(void), PetscInt *, const PetscInt **, void **);
282: PETSC_EXTERN PetscErrorCode DMIsBoundaryPoint(DM, PetscInt, PetscBool *);
283: PETSC_EXTERN PetscErrorCode DMCopyBoundary(DM, DM);

285: PETSC_EXTERN PetscErrorCode DMProjectFunction(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void**,InsertMode,Vec);
286: PETSC_EXTERN PetscErrorCode DMProjectFunctionLocal(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void**,InsertMode,Vec);
287: PETSC_EXTERN PetscErrorCode DMProjectFunctionLabelLocal(DM,PetscReal,DMLabel,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,InsertMode,Vec);
288: 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);
289: 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);
290: PETSC_EXTERN PetscErrorCode DMComputeL2Diff(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,Vec,PetscReal *);
291: PETSC_EXTERN PetscErrorCode DMComputeL2GradientDiff(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal [], const PetscReal [], PetscInt, PetscScalar *, void *), void **, Vec, const PetscReal [], PetscReal *);
292: PETSC_EXTERN PetscErrorCode DMComputeL2FieldDiff(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,Vec,PetscReal *);
293: #endif