Actual source code: petscdm.h

petsc-3.7.5 2017-01-01
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: /*J
 17:     DMType - String with the name of a PETSc DM

 19:    Level: beginner

 21: .seealso: DMSetType(), DM
 22: J*/
 23: typedef const char* DMType;
 24: #define DMDA        "da"
 25: #define DMCOMPOSITE "composite"
 26: #define DMSLICED    "sliced"
 27: #define DMSHELL     "shell"
 28: #define DMPLEX      "plex"
 29: #define DMCARTESIAN "cartesian"
 30: #define DMREDUNDANT "redundant"
 31: #define DMPATCH     "patch"
 32: #define DMMOAB      "moab"
 33: #define DMNETWORK   "network"
 34: #define DMFOREST    "forest"
 35: #define DMP4EST     "p4est"
 36: #define DMP8EST     "p8est"

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

 47: PETSC_EXTERN PetscErrorCode DMView(DM,PetscViewer);
 48: PETSC_EXTERN PetscErrorCode DMLoad(DM,PetscViewer);
 49: PETSC_EXTERN PetscErrorCode DMDestroy(DM*);
 50: PETSC_EXTERN PetscErrorCode DMCreateGlobalVector(DM,Vec*);
 51: PETSC_EXTERN PetscErrorCode DMCreateLocalVector(DM,Vec*);
 52: PETSC_EXTERN PetscErrorCode DMGetLocalVector(DM,Vec *);
 53: PETSC_EXTERN PetscErrorCode DMRestoreLocalVector(DM,Vec *);
 54: PETSC_EXTERN PetscErrorCode DMGetGlobalVector(DM,Vec *);
 55: PETSC_EXTERN PetscErrorCode DMRestoreGlobalVector(DM,Vec *);
 56: PETSC_EXTERN PetscErrorCode DMClearGlobalVectors(DM);
 57: PETSC_EXTERN PetscErrorCode DMClearLocalVectors(DM);
 58: PETSC_EXTERN PetscErrorCode DMHasNamedGlobalVector(DM,const char*,PetscBool*);
 59: PETSC_EXTERN PetscErrorCode DMGetNamedGlobalVector(DM,const char*,Vec*);
 60: PETSC_EXTERN PetscErrorCode DMRestoreNamedGlobalVector(DM,const char*,Vec*);
 61: PETSC_EXTERN PetscErrorCode DMHasNamedLocalVector(DM,const char*,PetscBool*);
 62: PETSC_EXTERN PetscErrorCode DMGetNamedLocalVector(DM,const char*,Vec*);
 63: PETSC_EXTERN PetscErrorCode DMRestoreNamedLocalVector(DM,const char*,Vec*);
 64: PETSC_EXTERN PetscErrorCode DMGetLocalToGlobalMapping(DM,ISLocalToGlobalMapping*);
 65: PETSC_EXTERN PetscErrorCode DMCreateFieldIS(DM,PetscInt*,char***,IS**);
 66: PETSC_EXTERN PetscErrorCode DMGetBlockSize(DM,PetscInt*);
 67: PETSC_EXTERN PetscErrorCode DMCreateColoring(DM,ISColoringType,ISColoring*);
 68: PETSC_EXTERN PetscErrorCode DMCreateMatrix(DM,Mat*);
 69: PETSC_EXTERN PetscErrorCode DMSetMatrixPreallocateOnly(DM,PetscBool);
 70: PETSC_EXTERN PetscErrorCode DMCreateInterpolation(DM,DM,Mat*,Vec*);
 71: PETSC_EXTERN PetscErrorCode DMCreateRestriction(DM,DM,Mat*);
 72: PETSC_EXTERN PetscErrorCode DMCreateInjection(DM,DM,Mat*);
 73: PETSC_EXTERN PetscErrorCode DMGetWorkArray(DM,PetscInt,PetscDataType,void*);
 74: PETSC_EXTERN PetscErrorCode DMRestoreWorkArray(DM,PetscInt,PetscDataType,void*);
 75: PETSC_EXTERN PetscErrorCode DMRefine(DM,MPI_Comm,DM*);
 76: PETSC_EXTERN PetscErrorCode DMCoarsen(DM,MPI_Comm,DM*);
 77: PETSC_EXTERN PetscErrorCode DMGetCoarseDM(DM,DM*);
 78: PETSC_EXTERN PetscErrorCode DMSetCoarseDM(DM,DM);
 79: PETSC_EXTERN PetscErrorCode DMGetFineDM(DM,DM*);
 80: PETSC_EXTERN PetscErrorCode DMSetFineDM(DM,DM);
 81: PETSC_EXTERN PetscErrorCode DMRefineHierarchy(DM,PetscInt,DM[]);
 82: PETSC_EXTERN PetscErrorCode DMCoarsenHierarchy(DM,PetscInt,DM[]);
 83: PETSC_EXTERN PetscErrorCode DMCoarsenHookAdd(DM,PetscErrorCode (*)(DM,DM,void*),PetscErrorCode (*)(DM,Mat,Vec,Mat,DM,void*),void*);
 84: PETSC_EXTERN PetscErrorCode DMRefineHookAdd(DM,PetscErrorCode (*)(DM,DM,void*),PetscErrorCode (*)(DM,Mat,DM,void*),void*);
 85: PETSC_EXTERN PetscErrorCode DMRestrict(DM,Mat,Vec,Mat,DM);
 86: PETSC_EXTERN PetscErrorCode DMInterpolate(DM,Mat,DM);
 87: PETSC_EXTERN PetscErrorCode DMSetFromOptions(DM);
 88: PETSC_STATIC_INLINE PetscErrorCode DMViewFromOptions(DM A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);}

 90: PETSC_EXTERN PetscErrorCode DMSetUp(DM);
 91: PETSC_EXTERN PetscErrorCode DMCreateInterpolationScale(DM,DM,Mat,Vec*);
 92: PETSC_EXTERN PetscErrorCode DMCreateAggregates(DM,DM,Mat*);
 93: PETSC_EXTERN PetscErrorCode DMGlobalToLocalHookAdd(DM,PetscErrorCode (*)(DM,Vec,InsertMode,Vec,void*),PetscErrorCode (*)(DM,Vec,InsertMode,Vec,void*),void*);
 94: PETSC_EXTERN PetscErrorCode DMLocalToGlobalHookAdd(DM,PetscErrorCode (*)(DM,Vec,InsertMode,Vec,void*),PetscErrorCode (*)(DM,Vec,InsertMode,Vec,void*),void*);
 95: PETSC_EXTERN PetscErrorCode DMGlobalToLocalBegin(DM,Vec,InsertMode,Vec);
 96: PETSC_EXTERN PetscErrorCode DMGlobalToLocalEnd(DM,Vec,InsertMode,Vec);
 97: PETSC_EXTERN PetscErrorCode DMLocalToGlobalBegin(DM,Vec,InsertMode,Vec);
 98: PETSC_EXTERN PetscErrorCode DMLocalToGlobalEnd(DM,Vec,InsertMode,Vec);
 99: PETSC_EXTERN PetscErrorCode DMLocalToLocalBegin(DM,Vec,InsertMode,Vec);
100: PETSC_EXTERN PetscErrorCode DMLocalToLocalEnd(DM,Vec,InsertMode,Vec);
101: PETSC_EXTERN PetscErrorCode DMConvert(DM,DMType,DM*);

103: /* Topology support */
104: PETSC_EXTERN PetscErrorCode DMGetDimension(DM,PetscInt*);
105: PETSC_EXTERN PetscErrorCode DMSetDimension(DM,PetscInt);
106: PETSC_EXTERN PetscErrorCode DMGetDimPoints(DM,PetscInt,PetscInt*,PetscInt*);
107: PETSC_EXTERN PetscErrorCode DMGetUseNatural(DM,PetscBool*);
108: PETSC_EXTERN PetscErrorCode DMSetUseNatural(DM,PetscBool);

110: /* Coordinate support */
111: PETSC_EXTERN PetscErrorCode DMGetCoordinateDM(DM,DM*);
112: PETSC_EXTERN PetscErrorCode DMSetCoordinateDM(DM,DM);
113: PETSC_EXTERN PetscErrorCode DMGetCoordinateDim(DM,PetscInt*);
114: PETSC_EXTERN PetscErrorCode DMSetCoordinateDim(DM,PetscInt);
115: PETSC_EXTERN PetscErrorCode DMGetCoordinateSection(DM,PetscSection*);
116: PETSC_EXTERN PetscErrorCode DMSetCoordinateSection(DM,PetscInt,PetscSection);
117: PETSC_EXTERN PetscErrorCode DMGetCoordinates(DM,Vec*);
118: PETSC_EXTERN PetscErrorCode DMSetCoordinates(DM,Vec);
119: PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocal(DM,Vec*);
120: PETSC_EXTERN PetscErrorCode DMSetCoordinatesLocal(DM,Vec);
121: PETSC_EXTERN PetscErrorCode DMLocatePoints(DM,Vec,PetscSF*);
122: PETSC_EXTERN PetscErrorCode DMGetPeriodicity(DM,const PetscReal**,const PetscReal**,const DMBoundaryType**);
123: PETSC_EXTERN PetscErrorCode DMSetPeriodicity(DM,const PetscReal[],const PetscReal[],const DMBoundaryType[]);
124: PETSC_EXTERN PetscErrorCode DMLocalizeCoordinate(DM, const PetscScalar[], PetscScalar[]);
125: PETSC_EXTERN PetscErrorCode DMLocalizeCoordinates(DM);
126: PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocalized(DM,PetscBool*);

128: /* block hook interface */
129: PETSC_EXTERN PetscErrorCode DMSubDomainHookAdd(DM,PetscErrorCode (*)(DM,DM,void*),PetscErrorCode (*)(DM,VecScatter,VecScatter,DM,void*),void*);
130: PETSC_EXTERN PetscErrorCode DMSubDomainRestrict(DM,VecScatter,VecScatter,DM);

132: PETSC_EXTERN PetscErrorCode DMSetOptionsPrefix(DM,const char []);
133: PETSC_EXTERN PetscErrorCode DMAppendOptionsPrefix(DM,const char []);
134: PETSC_EXTERN PetscErrorCode DMGetOptionsPrefix(DM,const char*[]);
135: PETSC_EXTERN PetscErrorCode DMSetVecType(DM,VecType);
136: PETSC_EXTERN PetscErrorCode DMGetVecType(DM,VecType*);
137: PETSC_EXTERN PetscErrorCode DMSetMatType(DM,MatType);
138: PETSC_EXTERN PetscErrorCode DMGetMatType(DM,MatType*);
139: PETSC_EXTERN PetscErrorCode DMSetApplicationContext(DM,void*);
140: PETSC_EXTERN PetscErrorCode DMSetApplicationContextDestroy(DM,PetscErrorCode (*)(void**));
141: PETSC_EXTERN PetscErrorCode DMGetApplicationContext(DM,void*);
142: PETSC_EXTERN PetscErrorCode DMSetVariableBounds(DM,PetscErrorCode (*)(DM,Vec,Vec));
143: PETSC_EXTERN PetscErrorCode DMHasVariableBounds(DM,PetscBool *);
144: PETSC_EXTERN PetscErrorCode DMHasColoring(DM,PetscBool *);
145: PETSC_EXTERN PetscErrorCode DMHasCreateRestriction(DM,PetscBool *);
146: PETSC_EXTERN PetscErrorCode DMComputeVariableBounds(DM,Vec,Vec);

148: PETSC_EXTERN PetscErrorCode DMCreateSubDM(DM, PetscInt, PetscInt[], IS *, DM *);
149: PETSC_EXTERN PetscErrorCode DMCreateFieldDecomposition(DM,PetscInt*,char***,IS**,DM**);
150: PETSC_EXTERN PetscErrorCode DMCreateDomainDecomposition(DM,PetscInt*,char***,IS**,IS**,DM**);
151: PETSC_EXTERN PetscErrorCode DMCreateDomainDecompositionScatters(DM,PetscInt,DM*,VecScatter**,VecScatter**,VecScatter**);

153: PETSC_EXTERN PetscErrorCode DMGetRefineLevel(DM,PetscInt*);
154: PETSC_EXTERN PetscErrorCode DMSetRefineLevel(DM,PetscInt);
155: PETSC_EXTERN PetscErrorCode DMGetCoarsenLevel(DM,PetscInt*);
156: PETSC_EXTERN PetscErrorCode DMFinalizePackage(void);

158: PETSC_EXTERN PetscErrorCode VecGetDM(Vec, DM*);
159: PETSC_EXTERN PetscErrorCode VecSetDM(Vec, DM);
160: PETSC_EXTERN PetscErrorCode MatGetDM(Mat, DM*);
161: PETSC_EXTERN PetscErrorCode MatSetDM(Mat, DM);

163: typedef struct NLF_DAAD* NLF;

165: #define DM_FILE_CLASSID 1211221

167: /* FEM support */
168: PETSC_EXTERN PetscErrorCode DMPrintCellVector(PetscInt, const char [], PetscInt, const PetscScalar []);
169: PETSC_EXTERN PetscErrorCode DMPrintCellMatrix(PetscInt, const char [], PetscInt, PetscInt, const PetscScalar []);
170: PETSC_EXTERN PetscErrorCode DMPrintLocalVec(DM, const char [], PetscReal, Vec);

172: PETSC_EXTERN PetscErrorCode DMGetDefaultSection(DM, PetscSection *);
173: PETSC_EXTERN PetscErrorCode DMSetDefaultSection(DM, PetscSection);
174: PETSC_EXTERN PetscErrorCode DMGetDefaultConstraints(DM, PetscSection *, Mat *);
175: PETSC_EXTERN PetscErrorCode DMSetDefaultConstraints(DM, PetscSection, Mat);
176: PETSC_EXTERN PetscErrorCode DMGetDefaultGlobalSection(DM, PetscSection *);
177: PETSC_EXTERN PetscErrorCode DMSetDefaultGlobalSection(DM, PetscSection);
178: PETSC_EXTERN PetscErrorCode DMGetDefaultSF(DM, PetscSF *);
179: PETSC_EXTERN PetscErrorCode DMSetDefaultSF(DM, PetscSF);
180: PETSC_EXTERN PetscErrorCode DMCreateDefaultSF(DM, PetscSection, PetscSection);
181: PETSC_EXTERN PetscErrorCode DMGetPointSF(DM, PetscSF *);
182: PETSC_EXTERN PetscErrorCode DMSetPointSF(DM, PetscSF);

184: PETSC_EXTERN PetscErrorCode DMGetOutputDM(DM, DM *);
185: PETSC_EXTERN PetscErrorCode DMGetOutputSequenceNumber(DM, PetscInt *, PetscReal *);
186: PETSC_EXTERN PetscErrorCode DMSetOutputSequenceNumber(DM, PetscInt, PetscReal);
187: PETSC_EXTERN PetscErrorCode DMOutputSequenceLoad(DM, PetscViewer, const char *, PetscInt, PetscReal *);

189: PETSC_EXTERN PetscErrorCode DMGetDS(DM, PetscDS *);
190: PETSC_EXTERN PetscErrorCode DMSetDS(DM, PetscDS);
191: PETSC_EXTERN PetscErrorCode DMGetNumFields(DM, PetscInt *);
192: PETSC_EXTERN PetscErrorCode DMSetNumFields(DM, PetscInt);
193: PETSC_EXTERN PetscErrorCode DMGetField(DM, PetscInt, PetscObject *);
194: PETSC_EXTERN PetscErrorCode DMSetField(DM, PetscInt, PetscObject);

196: 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;

198: struct _DMInterpolationInfo {
199:   MPI_Comm   comm;
200:   PetscInt   dim;    /*1 The spatial dimension of points */
201:   PetscInt   nInput; /* The number of input points */
202:   PetscReal *points; /* The input point coordinates */
203:   PetscInt  *cells;  /* The cell containing each point */
204:   PetscInt   n;      /* The number of local points */
205:   Vec        coords; /* The point coordinates */
206:   PetscInt   dof;    /* The number of components to interpolate */
207: };
208: typedef struct _DMInterpolationInfo *DMInterpolationInfo;

210: PETSC_EXTERN PetscErrorCode DMInterpolationCreate(MPI_Comm, DMInterpolationInfo *);
211: PETSC_EXTERN PetscErrorCode DMInterpolationSetDim(DMInterpolationInfo, PetscInt);
212: PETSC_EXTERN PetscErrorCode DMInterpolationGetDim(DMInterpolationInfo, PetscInt *);
213: PETSC_EXTERN PetscErrorCode DMInterpolationSetDof(DMInterpolationInfo, PetscInt);
214: PETSC_EXTERN PetscErrorCode DMInterpolationGetDof(DMInterpolationInfo, PetscInt *);
215: PETSC_EXTERN PetscErrorCode DMInterpolationAddPoints(DMInterpolationInfo, PetscInt, PetscReal[]);
216: PETSC_EXTERN PetscErrorCode DMInterpolationSetUp(DMInterpolationInfo, DM, PetscBool);
217: PETSC_EXTERN PetscErrorCode DMInterpolationGetCoordinates(DMInterpolationInfo, Vec *);
218: PETSC_EXTERN PetscErrorCode DMInterpolationGetVector(DMInterpolationInfo, Vec *);
219: PETSC_EXTERN PetscErrorCode DMInterpolationRestoreVector(DMInterpolationInfo, Vec *);
220: PETSC_EXTERN PetscErrorCode DMInterpolationEvaluate(DMInterpolationInfo, DM, Vec, Vec);
221: PETSC_EXTERN PetscErrorCode DMInterpolationDestroy(DMInterpolationInfo *);

223: PETSC_EXTERN PetscErrorCode DMCreateLabel(DM, const char []);
224: PETSC_EXTERN PetscErrorCode DMGetLabelValue(DM, const char[], PetscInt, PetscInt *);
225: PETSC_EXTERN PetscErrorCode DMSetLabelValue(DM, const char[], PetscInt, PetscInt);
226: PETSC_EXTERN PetscErrorCode DMClearLabelValue(DM, const char[], PetscInt, PetscInt);
227: PETSC_EXTERN PetscErrorCode DMGetLabelSize(DM, const char[], PetscInt *);
228: PETSC_EXTERN PetscErrorCode DMGetLabelIdIS(DM, const char[], IS *);
229: PETSC_EXTERN PetscErrorCode DMGetStratumSize(DM, const char [], PetscInt, PetscInt *);
230: PETSC_EXTERN PetscErrorCode DMGetStratumIS(DM, const char [], PetscInt, IS *);
231: PETSC_EXTERN PetscErrorCode DMClearLabelStratum(DM, const char[], PetscInt);
232: PETSC_EXTERN PetscErrorCode DMGetLabelOutput(DM, const char[], PetscBool *);
233: PETSC_EXTERN PetscErrorCode DMSetLabelOutput(DM, const char[], PetscBool);

235: PETSC_EXTERN PetscErrorCode DMGetNumLabels(DM, PetscInt *);
236: PETSC_EXTERN PetscErrorCode DMGetLabelName(DM, PetscInt, const char **);
237: PETSC_EXTERN PetscErrorCode DMHasLabel(DM, const char [], PetscBool *);
238: PETSC_EXTERN PetscErrorCode DMGetLabel(DM, const char *, DMLabel *);
239: PETSC_EXTERN PetscErrorCode DMGetLabelByNum(DM, PetscInt, DMLabel *);
240: PETSC_EXTERN PetscErrorCode DMAddLabel(DM, DMLabel);
241: PETSC_EXTERN PetscErrorCode DMRemoveLabel(DM, const char [], DMLabel *);
242: PETSC_EXTERN PetscErrorCode DMCopyLabels(DM, DM);

244: PETSC_EXTERN PetscErrorCode DMAddBoundary(DM, PetscBool, const char[], const char[], PetscInt, PetscInt, const PetscInt *, void (*)(), PetscInt, const PetscInt *, void *);
245: PETSC_EXTERN PetscErrorCode DMGetNumBoundary(DM, PetscInt *);
246: PETSC_EXTERN PetscErrorCode DMGetBoundary(DM, PetscInt, PetscBool *, const char **, const char **, PetscInt *, PetscInt *, const PetscInt **, void (**)(), PetscInt *, const PetscInt **, void **);
247: PETSC_EXTERN PetscErrorCode DMIsBoundaryPoint(DM, PetscInt, PetscBool *);
248: PETSC_EXTERN PetscErrorCode DMCopyBoundary(DM, DM);

250: PETSC_EXTERN PetscErrorCode DMProjectFunction(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void**,InsertMode,Vec);
251: PETSC_EXTERN PetscErrorCode DMProjectFunctionLocal(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void**,InsertMode,Vec);
252: PETSC_EXTERN PetscErrorCode DMProjectFunctionLabelLocal(DM,PetscReal,DMLabel,PetscInt,const PetscInt[],PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,InsertMode,Vec);
253: PETSC_EXTERN PetscErrorCode DMProjectFieldLocal(DM,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[],PetscScalar[]),InsertMode,Vec);
254: PETSC_EXTERN PetscErrorCode DMComputeL2Diff(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,Vec,PetscReal *);
255: PETSC_EXTERN PetscErrorCode DMComputeL2GradientDiff(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal [], const PetscReal [], PetscInt, PetscScalar *, void *), void **, Vec, const PetscReal [], PetscReal *);
256: PETSC_EXTERN PetscErrorCode DMComputeL2FieldDiff(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,Vec,PetscReal *);
257: #endif