Actual source code: petscdm.h

  1: /*
  2:       Objects to manage the interactions between the mesh data structures and the algebraic objects
  3: */
  4: #if !defined(PETSCDM_H)
  5: #define PETSCDM_H
  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;
 15: PETSC_EXTERN PetscClassId DMLABEL_CLASSID;

 17: #define DMLOCATEPOINT_POINT_NOT_FOUND -367

 19: /*J
 20:     DMType - String with the name of a PETSc DM

 22:    Level: beginner

 24: .seealso: DMSetType(), DM
 25: J*/
 26: typedef const char* DMType;
 27: #define DMDA        "da"
 28: #define DMCOMPOSITE "composite"
 29: #define DMSLICED    "sliced"
 30: #define DMSHELL     "shell"
 31: #define DMPLEX      "plex"
 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"
 40: #define DMPRODUCT   "product"
 41: #define DMSTAG      "stag"

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

 53: PETSC_EXTERN PetscErrorCode DMView(DM,PetscViewer);
 54: PETSC_EXTERN PetscErrorCode DMLoad(DM,PetscViewer);
 55: PETSC_EXTERN PetscErrorCode DMDestroy(DM*);
 56: PETSC_EXTERN PetscErrorCode DMCreateGlobalVector(DM,Vec*);
 57: PETSC_EXTERN PetscErrorCode DMCreateLocalVector(DM,Vec*);
 58: PETSC_EXTERN PetscErrorCode DMGetLocalVector(DM,Vec *);
 59: PETSC_EXTERN PetscErrorCode DMRestoreLocalVector(DM,Vec *);
 60: PETSC_EXTERN PetscErrorCode DMGetGlobalVector(DM,Vec *);
 61: PETSC_EXTERN PetscErrorCode DMRestoreGlobalVector(DM,Vec *);
 62: PETSC_EXTERN PetscErrorCode DMClearGlobalVectors(DM);
 63: PETSC_EXTERN PetscErrorCode DMClearLocalVectors(DM);
 64: PETSC_EXTERN PetscErrorCode DMHasNamedGlobalVector(DM,const char*,PetscBool*);
 65: PETSC_EXTERN PetscErrorCode DMGetNamedGlobalVector(DM,const char*,Vec*);
 66: PETSC_EXTERN PetscErrorCode DMRestoreNamedGlobalVector(DM,const char*,Vec*);
 67: PETSC_EXTERN PetscErrorCode DMHasNamedLocalVector(DM,const char*,PetscBool*);
 68: PETSC_EXTERN PetscErrorCode DMGetNamedLocalVector(DM,const char*,Vec*);
 69: PETSC_EXTERN PetscErrorCode DMRestoreNamedLocalVector(DM,const char*,Vec*);
 70: PETSC_EXTERN PetscErrorCode DMGetLocalToGlobalMapping(DM,ISLocalToGlobalMapping*);
 71: PETSC_EXTERN PetscErrorCode DMCreateFieldIS(DM,PetscInt*,char***,IS**);
 72: PETSC_EXTERN PetscErrorCode DMGetBlockSize(DM,PetscInt*);
 73: PETSC_EXTERN PetscErrorCode DMCreateColoring(DM,ISColoringType,ISColoring*);
 74: PETSC_EXTERN PetscErrorCode DMCreateMatrix(DM,Mat*);
 75: PETSC_EXTERN PetscErrorCode DMSetMatrixPreallocateOnly(DM,PetscBool);
 76: PETSC_EXTERN PetscErrorCode DMSetMatrixStructureOnly(DM,PetscBool);
 77: PETSC_EXTERN PetscErrorCode DMCreateInterpolation(DM,DM,Mat*,Vec*);
 78: PETSC_EXTERN PetscErrorCode DMCreateRestriction(DM,DM,Mat*);
 79: PETSC_EXTERN PetscErrorCode DMCreateInjection(DM,DM,Mat*);
 80: PETSC_EXTERN PetscErrorCode DMCreateMassMatrix(DM,DM,Mat*);
 81: PETSC_EXTERN PetscErrorCode DMGetWorkArray(DM,PetscInt,MPI_Datatype,void*);
 82: PETSC_EXTERN PetscErrorCode DMRestoreWorkArray(DM,PetscInt,MPI_Datatype,void*);
 83: PETSC_EXTERN PetscErrorCode DMRefine(DM,MPI_Comm,DM*);
 84: PETSC_EXTERN PetscErrorCode DMCoarsen(DM,MPI_Comm,DM*);
 85: PETSC_EXTERN PetscErrorCode DMGetCoarseDM(DM,DM*);
 86: PETSC_EXTERN PetscErrorCode DMSetCoarseDM(DM,DM);
 87: PETSC_EXTERN PetscErrorCode DMGetFineDM(DM,DM*);
 88: PETSC_EXTERN PetscErrorCode DMSetFineDM(DM,DM);
 89: PETSC_EXTERN PetscErrorCode DMRefineHierarchy(DM,PetscInt,DM[]);
 90: PETSC_EXTERN PetscErrorCode DMCoarsenHierarchy(DM,PetscInt,DM[]);
 91: PETSC_EXTERN PetscErrorCode DMCoarsenHookAdd(DM,PetscErrorCode (*)(DM,DM,void*),PetscErrorCode (*)(DM,Mat,Vec,Mat,DM,void*),void*);
 92: PETSC_EXTERN PetscErrorCode DMCoarsenHookRemove(DM,PetscErrorCode (*)(DM,DM,void*),PetscErrorCode (*)(DM,Mat,Vec,Mat,DM,void*),void*);
 93: PETSC_EXTERN PetscErrorCode DMRefineHookAdd(DM,PetscErrorCode (*)(DM,DM,void*),PetscErrorCode (*)(DM,Mat,DM,void*),void*);
 94: PETSC_EXTERN PetscErrorCode DMRefineHookRemove(DM,PetscErrorCode (*)(DM,DM,void*),PetscErrorCode (*)(DM,Mat,DM,void*),void*);
 95: PETSC_EXTERN PetscErrorCode DMRestrict(DM,Mat,Vec,Mat,DM);
 96: PETSC_EXTERN PetscErrorCode DMInterpolate(DM,Mat,DM);
 97: PETSC_EXTERN PetscErrorCode DMInterpolateSolution(DM,DM,Mat,Vec,Vec);
 98: PETSC_EXTERN PetscErrorCode DMSetFromOptions(DM);
 99: PETSC_EXTERN PetscErrorCode DMViewFromOptions(DM,PetscObject,const char[]);

101: PETSC_EXTERN PetscErrorCode DMAdaptLabel(DM,DMLabel,DM*);
102: PETSC_EXTERN PetscErrorCode DMAdaptMetric(DM, Vec, DMLabel, DM *);

104: PETSC_EXTERN PetscErrorCode DMSetUp(DM);
105: PETSC_EXTERN PetscErrorCode DMCreateInterpolationScale(DM,DM,Mat,Vec*);
106: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use DMDACreateAggregates() or DMCreateRestriction() (since version 3.12)") PetscErrorCode DMCreateAggregates(DM,DM,Mat*);
107: PETSC_EXTERN PetscErrorCode DMGlobalToLocalHookAdd(DM,PetscErrorCode (*)(DM,Vec,InsertMode,Vec,void*),PetscErrorCode (*)(DM,Vec,InsertMode,Vec,void*),void*);
108: PETSC_EXTERN PetscErrorCode DMLocalToGlobalHookAdd(DM,PetscErrorCode (*)(DM,Vec,InsertMode,Vec,void*),PetscErrorCode (*)(DM,Vec,InsertMode,Vec,void*),void*);
109: PETSC_EXTERN PetscErrorCode DMGlobalToLocal(DM,Vec,InsertMode,Vec);
110: PETSC_EXTERN PetscErrorCode DMGlobalToLocalBegin(DM,Vec,InsertMode,Vec);
111: PETSC_EXTERN PetscErrorCode DMGlobalToLocalEnd(DM,Vec,InsertMode,Vec);
112: PETSC_EXTERN PetscErrorCode DMLocalToGlobal(DM,Vec,InsertMode,Vec);
113: PETSC_EXTERN PetscErrorCode DMLocalToGlobalBegin(DM,Vec,InsertMode,Vec);
114: PETSC_EXTERN PetscErrorCode DMLocalToGlobalEnd(DM,Vec,InsertMode,Vec);
115: PETSC_EXTERN PetscErrorCode DMLocalToLocalBegin(DM,Vec,InsertMode,Vec);
116: PETSC_EXTERN PetscErrorCode DMLocalToLocalEnd(DM,Vec,InsertMode,Vec);
117: PETSC_EXTERN PetscErrorCode DMConvert(DM,DMType,DM*);

119: /* Topology support */
120: PETSC_EXTERN PetscErrorCode DMGetDimension(DM,PetscInt*);
121: PETSC_EXTERN PetscErrorCode DMSetDimension(DM,PetscInt);
122: PETSC_EXTERN PetscErrorCode DMGetDimPoints(DM,PetscInt,PetscInt*,PetscInt*);
123: PETSC_EXTERN PetscErrorCode DMGetUseNatural(DM,PetscBool*);
124: PETSC_EXTERN PetscErrorCode DMSetUseNatural(DM,PetscBool);

126: /* Coordinate support */
127: PETSC_EXTERN PetscErrorCode DMGetCoordinateDM(DM,DM*);
128: PETSC_EXTERN PetscErrorCode DMSetCoordinateDM(DM,DM);
129: PETSC_EXTERN PetscErrorCode DMGetCoordinateDim(DM,PetscInt*);
130: PETSC_EXTERN PetscErrorCode DMSetCoordinateDim(DM,PetscInt);
131: PETSC_EXTERN PetscErrorCode DMGetCoordinateSection(DM,PetscSection*);
132: PETSC_EXTERN PetscErrorCode DMSetCoordinateSection(DM,PetscInt,PetscSection);
133: PETSC_EXTERN PetscErrorCode DMGetCoordinates(DM,Vec*);
134: PETSC_EXTERN PetscErrorCode DMSetCoordinates(DM,Vec);
135: PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocal(DM,Vec*);
136: PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocalSetUp(DM);
137: PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocalNoncollective(DM,Vec*);
138: PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocalTuple(DM,IS,PetscSection*,Vec*);
139: PETSC_EXTERN PetscErrorCode DMSetCoordinatesLocal(DM,Vec);
140: PETSC_EXTERN PetscErrorCode DMLocatePoints(DM,Vec,DMPointLocationType,PetscSF*);
141: PETSC_EXTERN PetscErrorCode DMGetPeriodicity(DM,PetscBool*,const PetscReal**,const PetscReal**,const DMBoundaryType**);
142: PETSC_EXTERN PetscErrorCode DMSetPeriodicity(DM,PetscBool,const PetscReal[],const PetscReal[],const DMBoundaryType[]);
143: PETSC_EXTERN PetscErrorCode DMLocalizeCoordinate(DM, const PetscScalar[], PetscBool, PetscScalar[]);
144: PETSC_EXTERN PetscErrorCode DMLocalizeCoordinates(DM);
145: PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocalized(DM,PetscBool*);
146: PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocalizedLocal(DM,PetscBool*);
147: PETSC_EXTERN PetscErrorCode DMGetNeighbors(DM,PetscInt*,const PetscMPIInt**);
148: PETSC_EXTERN PetscErrorCode DMGetCoordinateField(DM,DMField*);
149: PETSC_EXTERN PetscErrorCode DMSetCoordinateField(DM,DMField);
150: PETSC_EXTERN PetscErrorCode DMGetBoundingBox(DM,PetscReal[],PetscReal[]);
151: PETSC_EXTERN PetscErrorCode DMGetLocalBoundingBox(DM,PetscReal[],PetscReal[]);
152: PETSC_EXTERN PetscErrorCode DMProjectCoordinates(DM,PetscFE);

154: /* block hook interface */
155: PETSC_EXTERN PetscErrorCode DMSubDomainHookAdd(DM,PetscErrorCode (*)(DM,DM,void*),PetscErrorCode (*)(DM,VecScatter,VecScatter,DM,void*),void*);
156: PETSC_EXTERN PetscErrorCode DMSubDomainHookRemove(DM,PetscErrorCode (*)(DM,DM,void*),PetscErrorCode (*)(DM,VecScatter,VecScatter,DM,void*),void*);
157: PETSC_EXTERN PetscErrorCode DMSubDomainRestrict(DM,VecScatter,VecScatter,DM);

159: PETSC_EXTERN PetscErrorCode DMSetOptionsPrefix(DM,const char []);
160: PETSC_EXTERN PetscErrorCode DMAppendOptionsPrefix(DM,const char []);
161: PETSC_EXTERN PetscErrorCode DMGetOptionsPrefix(DM,const char*[]);
162: PETSC_EXTERN PetscErrorCode DMSetVecType(DM,VecType);
163: PETSC_EXTERN PetscErrorCode DMGetVecType(DM,VecType*);
164: PETSC_EXTERN PetscErrorCode DMSetMatType(DM,MatType);
165: PETSC_EXTERN PetscErrorCode DMGetMatType(DM,MatType*);
166: PETSC_EXTERN PetscErrorCode DMSetISColoringType(DM,ISColoringType);
167: PETSC_EXTERN PetscErrorCode DMGetISColoringType(DM,ISColoringType*);
168: PETSC_EXTERN PetscErrorCode DMSetApplicationContext(DM,void*);
169: PETSC_EXTERN PetscErrorCode DMSetApplicationContextDestroy(DM,PetscErrorCode (*)(void**));
170: PETSC_EXTERN PetscErrorCode DMGetApplicationContext(DM,void*);
171: PETSC_EXTERN PetscErrorCode DMSetVariableBounds(DM,PetscErrorCode (*)(DM,Vec,Vec));
172: PETSC_EXTERN PetscErrorCode DMHasVariableBounds(DM,PetscBool *);
173: PETSC_EXTERN PetscErrorCode DMHasColoring(DM,PetscBool *);
174: PETSC_EXTERN PetscErrorCode DMHasCreateRestriction(DM,PetscBool *);
175: PETSC_EXTERN PetscErrorCode DMHasCreateInjection(DM,PetscBool *);
176: PETSC_EXTERN PetscErrorCode DMComputeVariableBounds(DM,Vec,Vec);

178: PETSC_EXTERN PetscErrorCode DMCreateSubDM(DM, PetscInt, const PetscInt[], IS *, DM *);
179: PETSC_EXTERN PetscErrorCode DMCreateSuperDM(DM[], PetscInt, IS **, DM *);
180: PETSC_EXTERN PetscErrorCode DMCreateSectionSubDM(DM,PetscInt,const PetscInt[],IS*,DM*);
181: PETSC_EXTERN PetscErrorCode DMCreateSectionSuperDM(DM[],PetscInt,IS**,DM*);
182: PETSC_EXTERN PetscErrorCode DMCreateFieldDecomposition(DM,PetscInt*,char***,IS**,DM**);
183: PETSC_EXTERN PetscErrorCode DMCreateDomainDecomposition(DM,PetscInt*,char***,IS**,IS**,DM**);
184: PETSC_EXTERN PetscErrorCode DMCreateDomainDecompositionScatters(DM,PetscInt,DM*,VecScatter**,VecScatter**,VecScatter**);

186: PETSC_EXTERN PetscErrorCode DMGetRefineLevel(DM,PetscInt*);
187: PETSC_EXTERN PetscErrorCode DMSetRefineLevel(DM,PetscInt);
188: PETSC_EXTERN PetscErrorCode DMGetCoarsenLevel(DM,PetscInt*);
189: PETSC_EXTERN PetscErrorCode DMSetCoarsenLevel(DM,PetscInt);
190: PETSC_EXTERN PetscErrorCode DMFinalizePackage(void);

192: PETSC_EXTERN PetscErrorCode VecGetDM(Vec, DM*);
193: PETSC_EXTERN PetscErrorCode VecSetDM(Vec, DM);
194: PETSC_EXTERN PetscErrorCode MatGetDM(Mat, DM*);
195: PETSC_EXTERN PetscErrorCode MatSetDM(Mat, DM);
196: PETSC_EXTERN PetscErrorCode MatFDColoringUseDM(Mat,MatFDColoring);

198: typedef struct NLF_DAAD* NLF;

200: #define DM_FILE_CLASSID 1211221

202: /* FEM support */
203: PETSC_EXTERN PetscErrorCode DMPrintCellVector(PetscInt, const char [], PetscInt, const PetscScalar []);
204: PETSC_EXTERN PetscErrorCode DMPrintCellMatrix(PetscInt, const char [], PetscInt, PetscInt, const PetscScalar []);
205: PETSC_EXTERN PetscErrorCode DMPrintLocalVec(DM, const char [], PetscReal, Vec);

207: PETSC_EXTERN PetscErrorCode DMSetNullSpaceConstructor(DM, PetscInt, PetscErrorCode (*)(DM, PetscInt, PetscInt, MatNullSpace *));
208: PETSC_EXTERN PetscErrorCode DMGetNullSpaceConstructor(DM, PetscInt, PetscErrorCode (**)(DM, PetscInt, PetscInt, MatNullSpace *));
209: PETSC_EXTERN PetscErrorCode DMSetNearNullSpaceConstructor(DM, PetscInt, PetscErrorCode (*)(DM, PetscInt, PetscInt, MatNullSpace *));
210: PETSC_EXTERN PetscErrorCode DMGetNearNullSpaceConstructor(DM, PetscInt, PetscErrorCode (**)(DM, PetscInt, PetscInt, MatNullSpace *));

212: PETSC_EXTERN PetscErrorCode DMGetSection(DM, PetscSection *); /* Use DMGetLocalSection() in new code (since v3.12) */
213: PETSC_EXTERN PetscErrorCode DMSetSection(DM, PetscSection);   /* Use DMSetLocalSection() in new code (since v3.12) */
214: PETSC_EXTERN PetscErrorCode DMGetLocalSection(DM, PetscSection *);
215: PETSC_EXTERN PetscErrorCode DMSetLocalSection(DM, PetscSection);
216: PETSC_EXTERN PetscErrorCode DMGetGlobalSection(DM, PetscSection *);
217: PETSC_EXTERN PetscErrorCode DMSetGlobalSection(DM, PetscSection);
218: PETSC_STATIC_INLINE PETSC_DEPRECATED_FUNCTION("Use DMGetSection() (since v3.9)") PetscErrorCode DMGetDefaultSection(DM dm, PetscSection *s) {return DMGetSection(dm,s);}
219: PETSC_STATIC_INLINE PETSC_DEPRECATED_FUNCTION("Use DMSetSection() (since v3.9)") PetscErrorCode DMSetDefaultSection(DM dm, PetscSection s) {return DMSetSection(dm,s);}
220: PETSC_STATIC_INLINE PETSC_DEPRECATED_FUNCTION("Use DMGetGlobalSection() (since v3.9)") PetscErrorCode DMGetDefaultGlobalSection(DM dm, PetscSection *s) {return DMGetGlobalSection(dm,s);}
221: PETSC_STATIC_INLINE PETSC_DEPRECATED_FUNCTION("Use DMSetGlobalSection() (since v3.9)") PetscErrorCode DMSetDefaultGlobalSection(DM dm, PetscSection s) {return DMSetGlobalSection(dm,s);}

223: PETSC_EXTERN PetscErrorCode DMGetSectionSF(DM, PetscSF*);
224: PETSC_EXTERN PetscErrorCode DMSetSectionSF(DM, PetscSF);
225: PETSC_EXTERN PetscErrorCode DMCreateSectionSF(DM, PetscSection, PetscSection);
226: PETSC_STATIC_INLINE PETSC_DEPRECATED_FUNCTION("Use DMGetSectionSF() (since v3.12)") PetscErrorCode DMGetDefaultSF(DM dm, PetscSF *s) {return DMGetSectionSF(dm,s);}
227: PETSC_STATIC_INLINE PETSC_DEPRECATED_FUNCTION("Use DMSetSectionSF() (since v3.12)") PetscErrorCode DMSetDefaultSF(DM dm, PetscSF s) {return DMSetSectionSF(dm,s);}
228: PETSC_STATIC_INLINE PETSC_DEPRECATED_FUNCTION("Use DMCreateSectionSF() (since v3.12)") PetscErrorCode DMCreateDefaultSF(DM dm, PetscSection l, PetscSection g) {return DMCreateSectionSF(dm,l,g);}
229: PETSC_EXTERN PetscErrorCode DMGetPointSF(DM, PetscSF *);
230: PETSC_EXTERN PetscErrorCode DMSetPointSF(DM, PetscSF);


233: PETSC_EXTERN PetscErrorCode DMGetDefaultConstraints(DM, PetscSection *, Mat *);
234: PETSC_EXTERN PetscErrorCode DMSetDefaultConstraints(DM, PetscSection, Mat);

236: PETSC_EXTERN PetscErrorCode DMGetOutputDM(DM, DM *);
237: PETSC_EXTERN PetscErrorCode DMGetOutputSequenceNumber(DM, PetscInt *, PetscReal *);
238: PETSC_EXTERN PetscErrorCode DMSetOutputSequenceNumber(DM, PetscInt, PetscReal);
239: PETSC_EXTERN PetscErrorCode DMOutputSequenceLoad(DM, PetscViewer, const char *, PetscInt, PetscReal *);

241: PETSC_EXTERN PetscErrorCode DMGetNumFields(DM, PetscInt *);
242: PETSC_EXTERN PetscErrorCode DMSetNumFields(DM, PetscInt);
243: PETSC_EXTERN PetscErrorCode DMGetField(DM, PetscInt, DMLabel *, PetscObject *);
244: PETSC_EXTERN PetscErrorCode DMSetField(DM, PetscInt, DMLabel, PetscObject);
245: PETSC_EXTERN PetscErrorCode DMAddField(DM, DMLabel, PetscObject);
246: PETSC_EXTERN PetscErrorCode DMSetFieldAvoidTensor(DM, PetscInt, PetscBool);
247: PETSC_EXTERN PetscErrorCode DMGetFieldAvoidTensor(DM, PetscInt, PetscBool *);
248: PETSC_EXTERN PetscErrorCode DMClearFields(DM);
249: PETSC_EXTERN PetscErrorCode DMCopyFields(DM, DM);
250: PETSC_EXTERN PetscErrorCode DMGetAdjacency(DM, PetscInt, PetscBool *, PetscBool *);
251: PETSC_EXTERN PetscErrorCode DMSetAdjacency(DM, PetscInt, PetscBool, PetscBool);
252: PETSC_EXTERN PetscErrorCode DMGetBasicAdjacency(DM, PetscBool *, PetscBool *);
253: PETSC_EXTERN PetscErrorCode DMSetBasicAdjacency(DM, PetscBool, PetscBool);

255: PETSC_EXTERN PetscErrorCode DMGetNumDS(DM, PetscInt *);
256: PETSC_EXTERN PetscErrorCode DMGetDS(DM, PetscDS *);
257: PETSC_EXTERN PetscErrorCode DMGetCellDS(DM, PetscInt, PetscDS *);
258: PETSC_EXTERN PetscErrorCode DMGetRegionDS(DM, DMLabel, IS *, PetscDS *);
259: PETSC_EXTERN PetscErrorCode DMSetRegionDS(DM, DMLabel, IS, PetscDS);
260: PETSC_EXTERN PetscErrorCode DMGetRegionNumDS(DM, PetscInt, DMLabel *, IS *, PetscDS *);
261: PETSC_EXTERN PetscErrorCode DMSetRegionNumDS(DM, PetscInt, DMLabel, IS, PetscDS);
262: PETSC_EXTERN PetscErrorCode DMFindRegionNum(DM, PetscDS, PetscInt *);
263: PETSC_EXTERN PetscErrorCode DMCreateDS(DM);
264: PETSC_EXTERN PetscErrorCode DMClearDS(DM);
265: PETSC_EXTERN PetscErrorCode DMCopyDS(DM, DM);
266: PETSC_EXTERN PetscErrorCode DMCopyDisc(DM, DM);
267: PETSC_EXTERN PetscErrorCode DMComputeExactSolution(DM, PetscReal, Vec, Vec);

269: /*MC
270:   DMInterpolationInfo - Structure for holding information about interpolation on a mesh

272:   Level: intermediate

274:   Synopsis:
275:     comm   - The communicator
276:     dim    - The spatial dimension of points
277:     nInput - The number of input points
278:     points - The input point coordinates
279:     cells  - The cell containing each point
280:     n      - The number of local points
281:     coords - The point coordinates
282:     dof    - The number of components to interpolate

284: .seealso: DMInterpolationCreate(), DMInterpolationEvaluate(), DMInterpolationAddPoints()
285: M*/
286: struct _DMInterpolationInfo {
287:   MPI_Comm   comm;
288:   PetscInt   dim;    /*1 The spatial dimension of points */
289:   PetscInt   nInput; /* The number of input points */
290:   PetscReal *points; /* The input point coordinates */
291:   PetscInt  *cells;  /* The cell containing each point */
292:   PetscInt   n;      /* The number of local points */
293:   Vec        coords; /* The point coordinates */
294:   PetscInt   dof;    /* The number of components to interpolate */
295: };
296: typedef struct _DMInterpolationInfo *DMInterpolationInfo;

298: PETSC_EXTERN PetscErrorCode DMInterpolationCreate(MPI_Comm, DMInterpolationInfo *);
299: PETSC_EXTERN PetscErrorCode DMInterpolationSetDim(DMInterpolationInfo, PetscInt);
300: PETSC_EXTERN PetscErrorCode DMInterpolationGetDim(DMInterpolationInfo, PetscInt *);
301: PETSC_EXTERN PetscErrorCode DMInterpolationSetDof(DMInterpolationInfo, PetscInt);
302: PETSC_EXTERN PetscErrorCode DMInterpolationGetDof(DMInterpolationInfo, PetscInt *);
303: PETSC_EXTERN PetscErrorCode DMInterpolationAddPoints(DMInterpolationInfo, PetscInt, PetscReal[]);
304: PETSC_EXTERN PetscErrorCode DMInterpolationSetUp(DMInterpolationInfo, DM, PetscBool, PetscBool);
305: PETSC_EXTERN PetscErrorCode DMInterpolationGetCoordinates(DMInterpolationInfo, Vec *);
306: PETSC_EXTERN PetscErrorCode DMInterpolationGetVector(DMInterpolationInfo, Vec *);
307: PETSC_EXTERN PetscErrorCode DMInterpolationRestoreVector(DMInterpolationInfo, Vec *);
308: PETSC_EXTERN PetscErrorCode DMInterpolationEvaluate(DMInterpolationInfo, DM, Vec, Vec);
309: PETSC_EXTERN PetscErrorCode DMInterpolationDestroy(DMInterpolationInfo *);

311: PETSC_EXTERN PetscErrorCode DMCreateLabel(DM, const char []);
312: PETSC_EXTERN PetscErrorCode DMGetLabelValue(DM, const char[], PetscInt, PetscInt *);
313: PETSC_EXTERN PetscErrorCode DMSetLabelValue(DM, const char[], PetscInt, PetscInt);
314: PETSC_EXTERN PetscErrorCode DMClearLabelValue(DM, const char[], PetscInt, PetscInt);
315: PETSC_EXTERN PetscErrorCode DMGetLabelSize(DM, const char[], PetscInt *);
316: PETSC_EXTERN PetscErrorCode DMGetLabelIdIS(DM, const char[], IS *);
317: PETSC_EXTERN PetscErrorCode DMGetStratumSize(DM, const char [], PetscInt, PetscInt *);
318: PETSC_EXTERN PetscErrorCode DMGetStratumIS(DM, const char [], PetscInt, IS *);
319: PETSC_EXTERN PetscErrorCode DMSetStratumIS(DM, const char [], PetscInt, IS);
320: PETSC_EXTERN PetscErrorCode DMClearLabelStratum(DM, const char[], PetscInt);
321: PETSC_EXTERN PetscErrorCode DMGetLabelOutput(DM, const char[], PetscBool *);
322: PETSC_EXTERN PetscErrorCode DMSetLabelOutput(DM, const char[], PetscBool);

324: PETSC_EXTERN PetscErrorCode DMGetNumLabels(DM, PetscInt *);
325: PETSC_EXTERN PetscErrorCode DMGetLabelName(DM, PetscInt, const char **);
326: PETSC_EXTERN PetscErrorCode DMHasLabel(DM, const char [], PetscBool *);
327: PETSC_EXTERN PetscErrorCode DMGetLabel(DM, const char *, DMLabel *);
328: PETSC_EXTERN PetscErrorCode DMGetLabelByNum(DM, PetscInt, DMLabel *);
329: PETSC_EXTERN PetscErrorCode DMAddLabel(DM, DMLabel);
330: PETSC_EXTERN PetscErrorCode DMRemoveLabel(DM, const char [], DMLabel *);
331: PETSC_EXTERN PetscErrorCode DMRemoveLabelBySelf(DM, DMLabel *, PetscBool);
332: PETSC_EXTERN PetscErrorCode DMCopyLabels(DM, DM, PetscCopyMode, PetscBool);

334: PETSC_EXTERN PetscErrorCode DMAddBoundary(DM, DMBoundaryConditionType, const char[], const char[], PetscInt, PetscInt, const PetscInt *, void (*)(void), void (*)(void), PetscInt, const PetscInt *, void *);
335: PETSC_EXTERN PetscErrorCode DMGetNumBoundary(DM, PetscInt *);
336: PETSC_EXTERN PetscErrorCode DMGetBoundary(DM, PetscInt, DMBoundaryConditionType *, const char **, const char **, PetscInt *, PetscInt *, const PetscInt **, void (**)(void), void (**)(void), PetscInt *, const PetscInt **, void **);
337: PETSC_EXTERN PetscErrorCode DMIsBoundaryPoint(DM, PetscInt, PetscBool *);
338: PETSC_EXTERN PetscErrorCode DMCopyBoundary(DM, DM);

340: PETSC_EXTERN PetscErrorCode DMProjectFunction(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void**,InsertMode,Vec);
341: PETSC_EXTERN PetscErrorCode DMProjectFunctionLocal(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void**,InsertMode,Vec);
342: PETSC_EXTERN PetscErrorCode DMProjectFunctionLabel(DM, PetscReal, DMLabel, PetscInt, const PetscInt[], PetscInt, const PetscInt[], PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **, InsertMode, Vec);
343: PETSC_EXTERN PetscErrorCode DMProjectFunctionLabelLocal(DM,PetscReal,DMLabel,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,InsertMode,Vec);
344: 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);
345: 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);
346: PETSC_EXTERN PetscErrorCode DMProjectBdFieldLabelLocal(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[],const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),InsertMode,Vec);
347: PETSC_EXTERN PetscErrorCode DMComputeL2Diff(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,Vec,PetscReal *);
348: PETSC_EXTERN PetscErrorCode DMComputeL2GradientDiff(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal [], const PetscReal [], PetscInt, PetscScalar *, void *), void **, Vec, const PetscReal [], PetscReal *);
349: PETSC_EXTERN PetscErrorCode DMComputeL2FieldDiff(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,Vec,PetscReal *);
350: PETSC_EXTERN PetscErrorCode DMComputeError(DM, Vec, PetscReal[], Vec *);
351: PETSC_EXTERN PetscErrorCode DMHasBasisTransform(DM,PetscBool*);
352: PETSC_EXTERN PetscErrorCode DMCopyTransform(DM, DM);

354: PETSC_EXTERN PetscErrorCode DMGetCompatibility(DM,DM,PetscBool*,PetscBool*);

356: PETSC_EXTERN PetscErrorCode DMMonitorSet(DM, PetscErrorCode (*)(DM, void *), void *, PetscErrorCode (*)(void**));
357: PETSC_EXTERN PetscErrorCode DMMonitorCancel(DM);
358: PETSC_EXTERN PetscErrorCode DMMonitorSetFromOptions(DM, const char[], const char[], const char[], PetscErrorCode (*)(DM, void *), PetscErrorCode (*)(DM, PetscViewerAndFormat *), PetscBool *);
359: PETSC_EXTERN PetscErrorCode DMMonitor(DM);

361: PETSC_STATIC_INLINE PetscInt DMPolytopeTypeGetDim(DMPolytopeType ct) {
362:   switch (ct) {
363:     case DM_POLYTOPE_POINT:
364:       return 0;
365:     case DM_POLYTOPE_SEGMENT:
366:     case DM_POLYTOPE_POINT_PRISM_TENSOR:
367:       return 1;
368:     case DM_POLYTOPE_TRIANGLE:
369:     case DM_POLYTOPE_QUADRILATERAL:
370:     case DM_POLYTOPE_SEG_PRISM_TENSOR:
371:       return 2;
372:     case DM_POLYTOPE_TETRAHEDRON:
373:     case DM_POLYTOPE_HEXAHEDRON:
374:     case DM_POLYTOPE_TRI_PRISM:
375:     case DM_POLYTOPE_TRI_PRISM_TENSOR:
376:     case DM_POLYTOPE_QUAD_PRISM_TENSOR:
377:     case DM_POLYTOPE_PYRAMID:
378:       return 3;
379:     default: return -1;
380:   }
381: }

383: PETSC_STATIC_INLINE PetscInt DMPolytopeTypeGetConeSize(DMPolytopeType ct)
384: {
385:   switch (ct) {
386:     case DM_POLYTOPE_POINT:              return 0;
387:     case DM_POLYTOPE_SEGMENT:            return 2;
388:     case DM_POLYTOPE_POINT_PRISM_TENSOR: return 2;
389:     case DM_POLYTOPE_TRIANGLE:           return 3;
390:     case DM_POLYTOPE_QUADRILATERAL:      return 4;
391:     case DM_POLYTOPE_SEG_PRISM_TENSOR:   return 4;
392:     case DM_POLYTOPE_TETRAHEDRON:        return 4;
393:     case DM_POLYTOPE_HEXAHEDRON:         return 6;
394:     case DM_POLYTOPE_TRI_PRISM:          return 5;
395:     case DM_POLYTOPE_TRI_PRISM_TENSOR:   return 5;
396:     case DM_POLYTOPE_QUAD_PRISM_TENSOR:  return 6;
397:     case DM_POLYTOPE_PYRAMID:            return 5;
398:     default: return -1;
399:   }
400: }

402: PETSC_STATIC_INLINE PetscInt DMPolytopeTypeGetNumVertices(DMPolytopeType ct)
403: {
404:   switch (ct) {
405:     case DM_POLYTOPE_POINT:              return 1;
406:     case DM_POLYTOPE_SEGMENT:            return 2;
407:     case DM_POLYTOPE_POINT_PRISM_TENSOR: return 2;
408:     case DM_POLYTOPE_TRIANGLE:           return 3;
409:     case DM_POLYTOPE_QUADRILATERAL:      return 4;
410:     case DM_POLYTOPE_SEG_PRISM_TENSOR:   return 4;
411:     case DM_POLYTOPE_TETRAHEDRON:        return 4;
412:     case DM_POLYTOPE_HEXAHEDRON:         return 8;
413:     case DM_POLYTOPE_TRI_PRISM:          return 6;
414:     case DM_POLYTOPE_TRI_PRISM_TENSOR:   return 6;
415:     case DM_POLYTOPE_QUAD_PRISM_TENSOR:  return 8;
416:     case DM_POLYTOPE_PYRAMID:            return 5;
417:     default: return -1;
418:   }
419: }

421: #endif