Actual source code: petscdm.h

petsc-3.9.4 2018-09-11
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 DMCreateMassMatrix(DM,DM,Mat*);
 78: PETSC_EXTERN PetscErrorCode DMGetWorkArray(DM,PetscInt,MPI_Datatype,void*);
 79: PETSC_EXTERN PetscErrorCode DMRestoreWorkArray(DM,PetscInt,MPI_Datatype,void*);
 80: PETSC_EXTERN PetscErrorCode DMRefine(DM,MPI_Comm,DM*);
 81: PETSC_EXTERN PetscErrorCode DMCoarsen(DM,MPI_Comm,DM*);
 82: PETSC_EXTERN PetscErrorCode DMGetCoarseDM(DM,DM*);
 83: PETSC_EXTERN PetscErrorCode DMSetCoarseDM(DM,DM);
 84: PETSC_EXTERN PetscErrorCode DMGetFineDM(DM,DM*);
 85: PETSC_EXTERN PetscErrorCode DMSetFineDM(DM,DM);
 86: PETSC_EXTERN PetscErrorCode DMRefineHierarchy(DM,PetscInt,DM[]);
 87: PETSC_EXTERN PetscErrorCode DMCoarsenHierarchy(DM,PetscInt,DM[]);
 88: PETSC_EXTERN PetscErrorCode DMCoarsenHookAdd(DM,PetscErrorCode (*)(DM,DM,void*),PetscErrorCode (*)(DM,Mat,Vec,Mat,DM,void*),void*);
 89: PETSC_EXTERN PetscErrorCode DMCoarsenHookRemove(DM,PetscErrorCode (*)(DM,DM,void*),PetscErrorCode (*)(DM,Mat,Vec,Mat,DM,void*),void*);
 90: PETSC_EXTERN PetscErrorCode DMRefineHookAdd(DM,PetscErrorCode (*)(DM,DM,void*),PetscErrorCode (*)(DM,Mat,DM,void*),void*);
 91: PETSC_EXTERN PetscErrorCode DMRefineHookRemove(DM,PetscErrorCode (*)(DM,DM,void*),PetscErrorCode (*)(DM,Mat,DM,void*),void*);
 92: PETSC_EXTERN PetscErrorCode DMRestrict(DM,Mat,Vec,Mat,DM);
 93: PETSC_EXTERN PetscErrorCode DMInterpolate(DM,Mat,DM);
 94: PETSC_EXTERN PetscErrorCode DMSetFromOptions(DM);
 95: PETSC_STATIC_INLINE PetscErrorCode DMViewFromOptions(DM A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);}

 97: PETSC_EXTERN PetscErrorCode DMAdaptLabel(DM,DMLabel,DM*);
 98: PETSC_EXTERN PetscErrorCode DMAdaptMetric(DM, Vec, DMLabel, DM *);

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

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

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

139: /* block hook interface */
140: PETSC_EXTERN PetscErrorCode DMSubDomainHookAdd(DM,PetscErrorCode (*)(DM,DM,void*),PetscErrorCode (*)(DM,VecScatter,VecScatter,DM,void*),void*);
141: PETSC_EXTERN PetscErrorCode DMSubDomainHookRemove(DM,PetscErrorCode (*)(DM,DM,void*),PetscErrorCode (*)(DM,VecScatter,VecScatter,DM,void*),void*);
142: PETSC_EXTERN PetscErrorCode DMSubDomainRestrict(DM,VecScatter,VecScatter,DM);

144: PETSC_EXTERN PetscErrorCode DMSetOptionsPrefix(DM,const char []);
145: PETSC_EXTERN PetscErrorCode DMAppendOptionsPrefix(DM,const char []);
146: PETSC_EXTERN PetscErrorCode DMGetOptionsPrefix(DM,const char*[]);
147: PETSC_EXTERN PetscErrorCode DMSetVecType(DM,VecType);
148: PETSC_EXTERN PetscErrorCode DMGetVecType(DM,VecType*);
149: PETSC_EXTERN PetscErrorCode DMSetMatType(DM,MatType);
150: PETSC_EXTERN PetscErrorCode DMGetMatType(DM,MatType*);
151: PETSC_EXTERN PetscErrorCode DMSetISColoringType(DM,ISColoringType);
152: PETSC_EXTERN PetscErrorCode DMGetISColoringType(DM,ISColoringType*);
153: PETSC_EXTERN PetscErrorCode DMSetApplicationContext(DM,void*);
154: PETSC_EXTERN PetscErrorCode DMSetApplicationContextDestroy(DM,PetscErrorCode (*)(void**));
155: PETSC_EXTERN PetscErrorCode DMGetApplicationContext(DM,void*);
156: PETSC_EXTERN PetscErrorCode DMSetVariableBounds(DM,PetscErrorCode (*)(DM,Vec,Vec));
157: PETSC_EXTERN PetscErrorCode DMHasVariableBounds(DM,PetscBool *);
158: PETSC_EXTERN PetscErrorCode DMHasColoring(DM,PetscBool *);
159: PETSC_EXTERN PetscErrorCode DMHasCreateRestriction(DM,PetscBool *);
160: PETSC_EXTERN PetscErrorCode DMHasCreateInjection(DM,PetscBool *);
161: PETSC_EXTERN PetscErrorCode DMComputeVariableBounds(DM,Vec,Vec);

163: PETSC_EXTERN PetscErrorCode DMCreateSubDM(DM, PetscInt, const PetscInt[], IS *, DM *);
164: PETSC_EXTERN PetscErrorCode DMCreateSuperDM(DM[], PetscInt, IS **, DM *);
165: PETSC_EXTERN PetscErrorCode DMCreateFieldDecomposition(DM,PetscInt*,char***,IS**,DM**);
166: PETSC_EXTERN PetscErrorCode DMCreateDomainDecomposition(DM,PetscInt*,char***,IS**,IS**,DM**);
167: PETSC_EXTERN PetscErrorCode DMCreateDomainDecompositionScatters(DM,PetscInt,DM*,VecScatter**,VecScatter**,VecScatter**);

169: PETSC_EXTERN PetscErrorCode DMGetRefineLevel(DM,PetscInt*);
170: PETSC_EXTERN PetscErrorCode DMSetRefineLevel(DM,PetscInt);
171: PETSC_EXTERN PetscErrorCode DMGetCoarsenLevel(DM,PetscInt*);
172: PETSC_EXTERN PetscErrorCode DMFinalizePackage(void);

174: PETSC_EXTERN PetscErrorCode VecGetDM(Vec, DM*);
175: PETSC_EXTERN PetscErrorCode VecSetDM(Vec, DM);
176: PETSC_EXTERN PetscErrorCode MatGetDM(Mat, DM*);
177: PETSC_EXTERN PetscErrorCode MatSetDM(Mat, DM);
178: PETSC_EXTERN PetscErrorCode MatFDColoringUseDM(Mat,MatFDColoring);

180: typedef struct NLF_DAAD* NLF;

182: #define DM_FILE_CLASSID 1211221

184: /* FEM support */
185: PETSC_EXTERN PetscErrorCode DMPrintCellVector(PetscInt, const char [], PetscInt, const PetscScalar []);
186: PETSC_EXTERN PetscErrorCode DMPrintCellMatrix(PetscInt, const char [], PetscInt, PetscInt, const PetscScalar []);
187: PETSC_EXTERN PetscErrorCode DMPrintLocalVec(DM, const char [], PetscReal, Vec);

189: PETSC_EXTERN PetscErrorCode DMGetDefaultSection(DM, PetscSection *);
190: PETSC_EXTERN PetscErrorCode DMSetDefaultSection(DM, PetscSection);
191: PETSC_EXTERN PetscErrorCode DMGetDefaultConstraints(DM, PetscSection *, Mat *);
192: PETSC_EXTERN PetscErrorCode DMSetDefaultConstraints(DM, PetscSection, Mat);
193: PETSC_EXTERN PetscErrorCode DMGetDefaultGlobalSection(DM, PetscSection *);
194: PETSC_EXTERN PetscErrorCode DMSetDefaultGlobalSection(DM, PetscSection);
195: PETSC_EXTERN PetscErrorCode DMGetDefaultSF(DM, PetscSF *);
196: PETSC_EXTERN PetscErrorCode DMSetDefaultSF(DM, PetscSF);
197: PETSC_EXTERN PetscErrorCode DMCreateDefaultSF(DM, PetscSection, PetscSection);
198: PETSC_EXTERN PetscErrorCode DMGetPointSF(DM, PetscSF *);
199: PETSC_EXTERN PetscErrorCode DMSetPointSF(DM, PetscSF);

201: PETSC_EXTERN PetscErrorCode DMGetOutputDM(DM, DM *);
202: PETSC_EXTERN PetscErrorCode DMGetOutputSequenceNumber(DM, PetscInt *, PetscReal *);
203: PETSC_EXTERN PetscErrorCode DMSetOutputSequenceNumber(DM, PetscInt, PetscReal);
204: PETSC_EXTERN PetscErrorCode DMOutputSequenceLoad(DM, PetscViewer, const char *, PetscInt, PetscReal *);

206: PETSC_EXTERN PetscErrorCode DMGetDS(DM, PetscDS *);
207: PETSC_EXTERN PetscErrorCode DMSetDS(DM, PetscDS);
208: PETSC_EXTERN PetscErrorCode DMGetNumFields(DM, PetscInt *);
209: PETSC_EXTERN PetscErrorCode DMSetNumFields(DM, PetscInt);
210: PETSC_EXTERN PetscErrorCode DMGetField(DM, PetscInt, PetscObject *);
211: PETSC_EXTERN PetscErrorCode DMSetField(DM, PetscInt, PetscObject);

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

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

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

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

260: PETSC_EXTERN PetscErrorCode DMAddBoundary(DM, DMBoundaryConditionType, const char[], const char[], PetscInt, PetscInt, const PetscInt *, void (*)(void), PetscInt, const PetscInt *, void *);
261: PETSC_EXTERN PetscErrorCode DMGetNumBoundary(DM, PetscInt *);
262: PETSC_EXTERN PetscErrorCode DMGetBoundary(DM, PetscInt, DMBoundaryConditionType *, const char **, const char **, PetscInt *, PetscInt *, const PetscInt **, void (**)(void), PetscInt *, const PetscInt **, void **);
263: PETSC_EXTERN PetscErrorCode DMIsBoundaryPoint(DM, PetscInt, PetscBool *);
264: PETSC_EXTERN PetscErrorCode DMCopyBoundary(DM, DM);

266: PETSC_EXTERN PetscErrorCode DMProjectFunction(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void**,InsertMode,Vec);
267: PETSC_EXTERN PetscErrorCode DMProjectFunctionLocal(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void**,InsertMode,Vec);
268: PETSC_EXTERN PetscErrorCode DMProjectFunctionLabel(DM, PetscReal, DMLabel, PetscInt, const PetscInt[], PetscInt, const PetscInt[], PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **, InsertMode, Vec);
269: PETSC_EXTERN PetscErrorCode DMProjectFunctionLabelLocal(DM,PetscReal,DMLabel,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,InsertMode,Vec);
270: 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);
271: 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);
272: PETSC_EXTERN PetscErrorCode DMComputeL2Diff(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,Vec,PetscReal *);
273: PETSC_EXTERN PetscErrorCode DMComputeL2GradientDiff(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal [], const PetscReal [], PetscInt, PetscScalar *, void *), void **, Vec, const PetscReal [], PetscReal *);
274: PETSC_EXTERN PetscErrorCode DMComputeL2FieldDiff(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,Vec,PetscReal *);
275: #endif