Actual source code: petscdm.h
petsc-3.6.4 2016-04-12
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>
11: PETSC_EXTERN PetscErrorCode DMInitializePackage(void);
13: PETSC_EXTERN PetscClassId DM_CLASSID;
15: /*J
16: DMType - String with the name of a PETSc DM
18: Level: beginner
20: .seealso: DMSetType(), DM
21: J*/
22: typedef const char* DMType;
23: #define DMDA "da"
24: #define DMCOMPOSITE "composite"
25: #define DMSLICED "sliced"
26: #define DMSHELL "shell"
27: #define DMPLEX "plex"
28: #define DMCARTESIAN "cartesian"
29: #define DMREDUNDANT "redundant"
30: #define DMPATCH "patch"
31: #define DMMOAB "moab"
32: #define DMNETWORK "network"
34: PETSC_EXTERN const char *const DMBoundaryTypes[];
35: PETSC_EXTERN PetscFunctionList DMList;
36: PETSC_EXTERN PetscErrorCode DMCreate(MPI_Comm,DM*);
37: PETSC_EXTERN PetscErrorCode DMClone(DM,DM*);
38: PETSC_EXTERN PetscErrorCode DMSetType(DM, DMType);
39: PETSC_EXTERN PetscErrorCode DMGetType(DM, DMType *);
40: PETSC_EXTERN PetscErrorCode DMRegister(const char[],PetscErrorCode (*)(DM));
41: PETSC_EXTERN PetscErrorCode DMRegisterDestroy(void);
43: PETSC_EXTERN PetscErrorCode DMView(DM,PetscViewer);
44: PETSC_EXTERN PetscErrorCode DMLoad(DM,PetscViewer);
45: PETSC_EXTERN PetscErrorCode DMDestroy(DM*);
46: PETSC_EXTERN PetscErrorCode DMCreateGlobalVector(DM,Vec*);
47: PETSC_EXTERN PetscErrorCode DMCreateLocalVector(DM,Vec*);
48: PETSC_EXTERN PetscErrorCode DMGetLocalVector(DM,Vec *);
49: PETSC_EXTERN PetscErrorCode DMRestoreLocalVector(DM,Vec *);
50: PETSC_EXTERN PetscErrorCode DMGetGlobalVector(DM,Vec *);
51: PETSC_EXTERN PetscErrorCode DMRestoreGlobalVector(DM,Vec *);
52: PETSC_EXTERN PetscErrorCode DMClearGlobalVectors(DM);
53: PETSC_EXTERN PetscErrorCode DMGetNamedGlobalVector(DM,const char*,Vec*);
54: PETSC_EXTERN PetscErrorCode DMRestoreNamedGlobalVector(DM,const char*,Vec*);
55: PETSC_EXTERN PetscErrorCode DMGetNamedLocalVector(DM,const char*,Vec*);
56: PETSC_EXTERN PetscErrorCode DMRestoreNamedLocalVector(DM,const char*,Vec*);
57: PETSC_EXTERN PetscErrorCode DMGetLocalToGlobalMapping(DM,ISLocalToGlobalMapping*);
58: PETSC_EXTERN PetscErrorCode DMCreateFieldIS(DM,PetscInt*,char***,IS**);
59: PETSC_EXTERN PetscErrorCode DMGetBlockSize(DM,PetscInt*);
60: PETSC_EXTERN PetscErrorCode DMCreateColoring(DM,ISColoringType,ISColoring*);
61: PETSC_EXTERN PetscErrorCode DMCreateMatrix(DM,Mat*);
62: PETSC_EXTERN PetscErrorCode DMSetMatrixPreallocateOnly(DM,PetscBool);
63: PETSC_EXTERN PetscErrorCode DMCreateInterpolation(DM,DM,Mat*,Vec*);
64: PETSC_EXTERN PetscErrorCode DMCreateInjection(DM,DM,Mat*);
65: PETSC_EXTERN PetscErrorCode DMGetWorkArray(DM,PetscInt,PetscDataType,void*);
66: PETSC_EXTERN PetscErrorCode DMRestoreWorkArray(DM,PetscInt,PetscDataType,void*);
67: PETSC_EXTERN PetscErrorCode DMRefine(DM,MPI_Comm,DM*);
68: PETSC_EXTERN PetscErrorCode DMCoarsen(DM,MPI_Comm,DM*);
69: PETSC_EXTERN PetscErrorCode DMRefineHierarchy(DM,PetscInt,DM[]);
70: PETSC_EXTERN PetscErrorCode DMCoarsenHierarchy(DM,PetscInt,DM[]);
71: PETSC_EXTERN PetscErrorCode DMCoarsenHookAdd(DM,PetscErrorCode (*)(DM,DM,void*),PetscErrorCode (*)(DM,Mat,Vec,Mat,DM,void*),void*);
72: PETSC_EXTERN PetscErrorCode DMRefineHookAdd(DM,PetscErrorCode (*)(DM,DM,void*),PetscErrorCode (*)(DM,Mat,DM,void*),void*);
73: PETSC_EXTERN PetscErrorCode DMRestrict(DM,Mat,Vec,Mat,DM);
74: PETSC_EXTERN PetscErrorCode DMInterpolate(DM,Mat,DM);
75: PETSC_EXTERN PetscErrorCode DMSetFromOptions(DM);
76: PETSC_STATIC_INLINE PetscErrorCode DMViewFromOptions(DM A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);}
78: PETSC_EXTERN PetscErrorCode DMSetUp(DM);
79: PETSC_EXTERN PetscErrorCode DMCreateInterpolationScale(DM,DM,Mat,Vec*);
80: PETSC_EXTERN PetscErrorCode DMCreateAggregates(DM,DM,Mat*);
81: PETSC_EXTERN PetscErrorCode DMGlobalToLocalHookAdd(DM,PetscErrorCode (*)(DM,Vec,InsertMode,Vec,void*),PetscErrorCode (*)(DM,Vec,InsertMode,Vec,void*),void*);
82: PETSC_EXTERN PetscErrorCode DMLocalToGlobalHookAdd(DM,PetscErrorCode (*)(DM,Vec,InsertMode,Vec,void*),PetscErrorCode (*)(DM,Vec,InsertMode,Vec,void*),void*);
83: PETSC_EXTERN PetscErrorCode DMGlobalToLocalBegin(DM,Vec,InsertMode,Vec);
84: PETSC_EXTERN PetscErrorCode DMGlobalToLocalEnd(DM,Vec,InsertMode,Vec);
85: PETSC_EXTERN PetscErrorCode DMLocalToGlobalBegin(DM,Vec,InsertMode,Vec);
86: PETSC_EXTERN PetscErrorCode DMLocalToGlobalEnd(DM,Vec,InsertMode,Vec);
87: PETSC_EXTERN PetscErrorCode DMLocalToLocalBegin(DM,Vec,InsertMode,Vec);
88: PETSC_EXTERN PetscErrorCode DMLocalToLocalEnd(DM,Vec,InsertMode,Vec);
89: PETSC_EXTERN PetscErrorCode DMConvert(DM,DMType,DM*);
91: /* Topology support */
92: PETSC_EXTERN PetscErrorCode DMGetDimension(DM,PetscInt*);
93: PETSC_EXTERN PetscErrorCode DMSetDimension(DM,PetscInt);
94: PETSC_EXTERN PetscErrorCode DMGetDimPoints(DM,PetscInt,PetscInt*,PetscInt*);
96: /* Coordinate support */
97: PETSC_EXTERN PetscErrorCode DMGetCoordinateDM(DM,DM*);
98: PETSC_EXTERN PetscErrorCode DMSetCoordinateDM(DM,DM);
99: PETSC_EXTERN PetscErrorCode DMGetCoordinateDim(DM,PetscInt*);
100: PETSC_EXTERN PetscErrorCode DMSetCoordinateDim(DM,PetscInt);
101: PETSC_EXTERN PetscErrorCode DMGetCoordinateSection(DM,PetscSection*);
102: PETSC_EXTERN PetscErrorCode DMSetCoordinateSection(DM,PetscInt,PetscSection);
103: PETSC_EXTERN PetscErrorCode DMGetCoordinates(DM,Vec*);
104: PETSC_EXTERN PetscErrorCode DMSetCoordinates(DM,Vec);
105: PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocal(DM,Vec*);
106: PETSC_EXTERN PetscErrorCode DMSetCoordinatesLocal(DM,Vec);
107: PETSC_EXTERN PetscErrorCode DMLocatePoints(DM,Vec,IS*);
108: PETSC_EXTERN PetscErrorCode DMGetPeriodicity(DM,const PetscReal**,const PetscReal**,const DMBoundaryType**);
109: PETSC_EXTERN PetscErrorCode DMSetPeriodicity(DM,const PetscReal[],const PetscReal[],const DMBoundaryType[]);
111: /* block hook interface */
112: PETSC_EXTERN PetscErrorCode DMSubDomainHookAdd(DM,PetscErrorCode (*)(DM,DM,void*),PetscErrorCode (*)(DM,VecScatter,VecScatter,DM,void*),void*);
113: PETSC_EXTERN PetscErrorCode DMSubDomainRestrict(DM,VecScatter,VecScatter,DM);
115: PETSC_EXTERN PetscErrorCode DMSetOptionsPrefix(DM,const char []);
116: PETSC_EXTERN PetscErrorCode DMSetVecType(DM,VecType);
117: PETSC_EXTERN PetscErrorCode DMGetVecType(DM,VecType*);
118: PETSC_EXTERN PetscErrorCode DMSetMatType(DM,MatType);
119: PETSC_EXTERN PetscErrorCode DMGetMatType(DM,MatType*);
120: PETSC_EXTERN PetscErrorCode DMSetApplicationContext(DM,void*);
121: PETSC_EXTERN PetscErrorCode DMSetApplicationContextDestroy(DM,PetscErrorCode (*)(void**));
122: PETSC_EXTERN PetscErrorCode DMGetApplicationContext(DM,void*);
123: PETSC_EXTERN PetscErrorCode DMSetVariableBounds(DM,PetscErrorCode (*)(DM,Vec,Vec));
124: PETSC_EXTERN PetscErrorCode DMHasVariableBounds(DM,PetscBool *);
125: PETSC_EXTERN PetscErrorCode DMHasColoring(DM,PetscBool *);
126: PETSC_EXTERN PetscErrorCode DMComputeVariableBounds(DM,Vec,Vec);
128: PETSC_EXTERN PetscErrorCode DMCreateSubDM(DM, PetscInt, PetscInt[], IS *, DM *);
129: PETSC_EXTERN PetscErrorCode DMCreateFieldDecomposition(DM,PetscInt*,char***,IS**,DM**);
130: PETSC_EXTERN PetscErrorCode DMCreateDomainDecomposition(DM,PetscInt*,char***,IS**,IS**,DM**);
131: PETSC_EXTERN PetscErrorCode DMCreateDomainDecompositionScatters(DM,PetscInt,DM*,VecScatter**,VecScatter**,VecScatter**);
133: PETSC_EXTERN PetscErrorCode DMGetRefineLevel(DM,PetscInt*);
134: PETSC_EXTERN PetscErrorCode DMGetCoarsenLevel(DM,PetscInt*);
135: PETSC_EXTERN PetscErrorCode DMFinalizePackage(void);
137: PETSC_EXTERN PetscErrorCode VecGetDM(Vec, DM*);
138: PETSC_EXTERN PetscErrorCode VecSetDM(Vec, DM);
139: PETSC_EXTERN PetscErrorCode MatGetDM(Mat, DM*);
140: PETSC_EXTERN PetscErrorCode MatSetDM(Mat, DM);
142: typedef struct NLF_DAAD* NLF;
144: #define DM_FILE_CLASSID 1211221
146: /* FEM support */
147: PETSC_EXTERN PetscErrorCode DMPrintCellVector(PetscInt, const char [], PetscInt, const PetscScalar []);
148: PETSC_EXTERN PetscErrorCode DMPrintCellMatrix(PetscInt, const char [], PetscInt, PetscInt, const PetscScalar []);
149: PETSC_EXTERN PetscErrorCode DMPrintLocalVec(DM, const char [], PetscReal, Vec);
151: PETSC_EXTERN PetscErrorCode DMGetDefaultSection(DM, PetscSection *);
152: PETSC_EXTERN PetscErrorCode DMSetDefaultSection(DM, PetscSection);
153: PETSC_EXTERN PetscErrorCode DMGetDefaultConstraints(DM, PetscSection *, Mat *);
154: PETSC_EXTERN PetscErrorCode DMSetDefaultConstraints(DM, PetscSection, Mat);
155: PETSC_EXTERN PetscErrorCode DMGetDefaultGlobalSection(DM, PetscSection *);
156: PETSC_EXTERN PetscErrorCode DMSetDefaultGlobalSection(DM, PetscSection);
157: PETSC_EXTERN PetscErrorCode DMGetDefaultSF(DM, PetscSF *);
158: PETSC_EXTERN PetscErrorCode DMSetDefaultSF(DM, PetscSF);
159: PETSC_EXTERN PetscErrorCode DMCreateDefaultSF(DM, PetscSection, PetscSection);
160: PETSC_EXTERN PetscErrorCode DMGetPointSF(DM, PetscSF *);
161: PETSC_EXTERN PetscErrorCode DMSetPointSF(DM, PetscSF);
163: PETSC_EXTERN PetscErrorCode DMGetOutputDM(DM, DM *);
164: PETSC_EXTERN PetscErrorCode DMGetOutputSequenceNumber(DM, PetscInt *, PetscReal *);
165: PETSC_EXTERN PetscErrorCode DMSetOutputSequenceNumber(DM, PetscInt, PetscReal);
166: PETSC_EXTERN PetscErrorCode DMOutputSequenceLoad(DM, PetscViewer, const char *, PetscInt, PetscReal *);
168: PETSC_EXTERN PetscErrorCode DMGetDS(DM, PetscDS *);
169: PETSC_EXTERN PetscErrorCode DMSetDS(DM, PetscDS);
170: PETSC_EXTERN PetscErrorCode DMGetNumFields(DM, PetscInt *);
171: PETSC_EXTERN PetscErrorCode DMSetNumFields(DM, PetscInt);
172: PETSC_EXTERN PetscErrorCode DMGetField(DM, PetscInt, PetscObject *);
173: PETSC_EXTERN PetscErrorCode DMSetField(DM, PetscInt, PetscObject);
175: 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;
177: struct _DMInterpolationInfo {
178: MPI_Comm comm;
179: PetscInt dim; /*1 The spatial dimension of points */
180: PetscInt nInput; /* The number of input points */
181: PetscReal *points; /* The input point coordinates */
182: PetscInt *cells; /* The cell containing each point */
183: PetscInt n; /* The number of local points */
184: Vec coords; /* The point coordinates */
185: PetscInt dof; /* The number of components to interpolate */
186: };
187: typedef struct _DMInterpolationInfo *DMInterpolationInfo;
189: PETSC_EXTERN PetscErrorCode DMInterpolationCreate(MPI_Comm, DMInterpolationInfo *);
190: PETSC_EXTERN PetscErrorCode DMInterpolationSetDim(DMInterpolationInfo, PetscInt);
191: PETSC_EXTERN PetscErrorCode DMInterpolationGetDim(DMInterpolationInfo, PetscInt *);
192: PETSC_EXTERN PetscErrorCode DMInterpolationSetDof(DMInterpolationInfo, PetscInt);
193: PETSC_EXTERN PetscErrorCode DMInterpolationGetDof(DMInterpolationInfo, PetscInt *);
194: PETSC_EXTERN PetscErrorCode DMInterpolationAddPoints(DMInterpolationInfo, PetscInt, PetscReal[]);
195: PETSC_EXTERN PetscErrorCode DMInterpolationSetUp(DMInterpolationInfo, DM, PetscBool);
196: PETSC_EXTERN PetscErrorCode DMInterpolationGetCoordinates(DMInterpolationInfo, Vec *);
197: PETSC_EXTERN PetscErrorCode DMInterpolationGetVector(DMInterpolationInfo, Vec *);
198: PETSC_EXTERN PetscErrorCode DMInterpolationRestoreVector(DMInterpolationInfo, Vec *);
199: PETSC_EXTERN PetscErrorCode DMInterpolationEvaluate(DMInterpolationInfo, DM, Vec, Vec);
200: PETSC_EXTERN PetscErrorCode DMInterpolationDestroy(DMInterpolationInfo *);
201: #endif