Actual source code: petscdmda.h
1: #pragma once
3: #include <petscdm.h>
4: #include <petscdmdatypes.h>
5: #include <petscpf.h>
6: #include <petscao.h>
7: #include <petscfe.h>
9: /* SUBMANSEC = DMDA */
11: /*MC
12: DMDA_STENCIL_STAR - "Star"-type stencil. In logical grid coordinates, only (i,j,k), (i+s,j,k), (i,j+s,k),
13: (i,j,k+s) are in the stencil NOT, for example, (i+s,j+s,k)
15: Level: beginner
17: Note:
18: Determines what ghost point values are brought over to each process in `DMGlobalToLocalBegin()`/ `DMGlobalToLocalEnd()`; in this case the "corner" values are not
19: brought over and hence should not be accessed locally
21: .seealso: [](ch_dmbase), `DMDA`, `DMDA_STENCIL_BOX`, `DMDAStencilType`, `DMDASetStencilType()`
22: M*/
24: /*MC
25: DMDA_STENCIL_BOX - "Box"-type stencil. In logical grid coordinates, any of (i,j,k), (i+s,j+r,k+t) may
26: be in the stencil.
28: Level: beginner
30: Note:
31: Determines what ghost point values are brought over to each process in `DMGlobalToLocalBegin()`/ `DMGlobalToLocalEnd()`
33: .seealso: [](ch_dmbase), `DMDA`, `DMDA_STENCIL_STAR`, `DMDAStencilType`, `DMDASetStencilType()`
34: M*/
36: PETSC_EXTERN PetscErrorCode DMDASetInterpolationType(DM, DMDAInterpolationType);
37: PETSC_EXTERN PetscErrorCode DMDAGetInterpolationType(DM, DMDAInterpolationType *);
38: PETSC_EXTERN PetscErrorCode DMDACreateAggregates(DM, DM, Mat *);
40: /* FEM */
41: PETSC_EXTERN PetscErrorCode DMDASetElementType(DM, DMDAElementType);
42: PETSC_EXTERN PetscErrorCode DMDAGetElementType(DM, DMDAElementType *);
43: PETSC_EXTERN PetscErrorCode DMDAGetElements(DM, PetscInt *, PetscInt *, const PetscInt *[]);
44: PETSC_EXTERN PetscErrorCode DMDARestoreElements(DM, PetscInt *, PetscInt *, const PetscInt *[]);
45: PETSC_EXTERN PetscErrorCode DMDAGetElementsSizes(DM, PetscInt *, PetscInt *, PetscInt *);
46: PETSC_EXTERN PetscErrorCode DMDAGetElementsCorners(DM, PetscInt *, PetscInt *, PetscInt *);
47: PETSC_EXTERN PetscErrorCode DMDAGetSubdomainCornersIS(DM, IS *);
48: PETSC_EXTERN PetscErrorCode DMDARestoreSubdomainCornersIS(DM, IS *);
50: #define MATSEQUSFFT "sequsfft"
52: PETSC_EXTERN PetscErrorCode DMDACreate(MPI_Comm, DM *);
53: PETSC_EXTERN PetscErrorCode DMDASetSizes(DM, PetscInt, PetscInt, PetscInt);
54: PETSC_EXTERN PetscErrorCode DMDACreate1d(MPI_Comm, DMBoundaryType, PetscInt, PetscInt, PetscInt, const PetscInt[], DM *);
55: PETSC_EXTERN PetscErrorCode DMDACreate2d(MPI_Comm, DMBoundaryType, DMBoundaryType, DMDAStencilType, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], DM *);
56: PETSC_EXTERN PetscErrorCode DMDACreate3d(MPI_Comm, DMBoundaryType, DMBoundaryType, DMBoundaryType, DMDAStencilType, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscInt[], DM *);
58: PETSC_EXTERN PetscErrorCode DMDAGlobalToNaturalBegin(DM, Vec, InsertMode, Vec);
59: PETSC_EXTERN PetscErrorCode DMDAGlobalToNaturalEnd(DM, Vec, InsertMode, Vec);
60: PETSC_EXTERN PetscErrorCode DMDANaturalToGlobalBegin(DM, Vec, InsertMode, Vec);
61: PETSC_EXTERN PetscErrorCode DMDANaturalToGlobalEnd(DM, Vec, InsertMode, Vec);
62: PETSC_DEPRECATED_FUNCTION(3, 5, 0, "DMLocalToLocalBegin()", ) static inline PetscErrorCode DMDALocalToLocalBegin(DM dm, Vec g, InsertMode mode, Vec l)
63: {
64: return DMLocalToLocalBegin(dm, g, mode, l);
65: }
66: PETSC_DEPRECATED_FUNCTION(3, 5, 0, "DMLocalToLocalEnd()", ) static inline PetscErrorCode DMDALocalToLocalEnd(DM dm, Vec g, InsertMode mode, Vec l)
67: {
68: return DMLocalToLocalEnd(dm, g, mode, l);
69: }
70: PETSC_EXTERN PetscErrorCode DMDACreateNaturalVector(DM, Vec *);
72: PETSC_EXTERN PetscErrorCode DMDAGetCorners(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *);
73: PETSC_EXTERN PetscErrorCode DMDAGetGhostCorners(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *);
74: PETSC_EXTERN PetscErrorCode DMDAGetInfo(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *, DMBoundaryType *, DMBoundaryType *, DMBoundaryType *, DMDAStencilType *);
75: PETSC_EXTERN PetscErrorCode DMDAGetProcessorSubset(DM, DMDirection, PetscInt, MPI_Comm *);
76: PETSC_EXTERN PetscErrorCode DMDAGetProcessorSubsets(DM, DMDirection, MPI_Comm *);
77: PETSC_EXTERN PetscErrorCode DMDAGetRay(DM, DMDirection, PetscInt, Vec *, VecScatter *);
79: PETSC_EXTERN PetscErrorCode DMDAGlobalToNaturalAllCreate(DM, VecScatter *);
80: PETSC_EXTERN PetscErrorCode DMDANaturalAllToGlobalCreate(DM, VecScatter *);
82: PETSC_EXTERN PetscErrorCode DMDAGetScatter(DM, VecScatter *, VecScatter *);
83: PETSC_EXTERN PetscErrorCode DMDAGetNeighbors(DM, const PetscMPIInt **);
85: PETSC_EXTERN PetscErrorCode DMDASetAOType(DM, AOType);
86: PETSC_EXTERN PetscErrorCode DMDAGetAO(DM, AO *);
87: PETSC_EXTERN PetscErrorCode DMDASetUniformCoordinates(DM, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal);
88: PETSC_EXTERN PetscErrorCode DMDASetGLLCoordinates(DM, PetscInt, PetscReal *);
89: PETSC_EXTERN PetscErrorCode DMDAGetCoordinateArray(DM, void *);
90: PETSC_EXTERN PetscErrorCode DMDARestoreCoordinateArray(DM, void *);
91: PETSC_EXTERN PetscErrorCode DMDAGetLogicalCoordinate(DM, PetscScalar, PetscScalar, PetscScalar, PetscInt *, PetscInt *, PetscInt *, PetscScalar *, PetscScalar *, PetscScalar *);
92: /* function to wrap coordinates around boundary */
93: PETSC_EXTERN PetscErrorCode DMDAMapCoordsToPeriodicDomain(DM, PetscScalar *, PetscScalar *);
95: PETSC_EXTERN PetscErrorCode DMDACreateCompatibleDMDA(DM, PetscInt, DM *);
96: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 10, 0, "DMDACreateCompatibleDMDA()", ) PetscErrorCode DMDAGetReducedDMDA(DM, PetscInt, DM *);
98: PETSC_EXTERN PetscErrorCode DMDASetFieldName(DM, PetscInt, const char[]);
99: PETSC_EXTERN PetscErrorCode DMDAGetFieldName(DM, PetscInt, const char **);
100: PETSC_EXTERN PetscErrorCode DMDASetFieldNames(DM, const char *const *);
101: PETSC_EXTERN PetscErrorCode DMDAGetFieldNames(DM, const char *const **);
102: PETSC_EXTERN PetscErrorCode DMDASetCoordinateName(DM, PetscInt, const char[]);
103: PETSC_EXTERN PetscErrorCode DMDAGetCoordinateName(DM, PetscInt, const char **);
105: PETSC_EXTERN PetscErrorCode DMDASetBoundaryType(DM, DMBoundaryType, DMBoundaryType, DMBoundaryType);
106: PETSC_EXTERN PetscErrorCode DMDAGetBoundaryType(DM, DMBoundaryType *, DMBoundaryType *, DMBoundaryType *);
107: PETSC_EXTERN PetscErrorCode DMDASetDof(DM, PetscInt);
108: PETSC_EXTERN PetscErrorCode DMDAGetDof(DM, PetscInt *);
109: PETSC_EXTERN PetscErrorCode DMDASetOverlap(DM, PetscInt, PetscInt, PetscInt);
110: PETSC_EXTERN PetscErrorCode DMDAGetOverlap(DM, PetscInt *, PetscInt *, PetscInt *);
111: PETSC_EXTERN PetscErrorCode DMDASetNumLocalSubDomains(DM, PetscInt);
112: PETSC_EXTERN PetscErrorCode DMDAGetNumLocalSubDomains(DM, PetscInt *);
113: PETSC_EXTERN PetscErrorCode DMDAGetOffset(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *);
114: PETSC_EXTERN PetscErrorCode DMDASetOffset(DM, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt);
115: PETSC_EXTERN PetscErrorCode DMDAGetNonOverlappingRegion(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *);
116: PETSC_EXTERN PetscErrorCode DMDASetNonOverlappingRegion(DM, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt);
117: PETSC_EXTERN PetscErrorCode DMDASetStencilWidth(DM, PetscInt);
118: PETSC_EXTERN PetscErrorCode DMDAGetStencilWidth(DM, PetscInt *);
119: PETSC_EXTERN PetscErrorCode DMDAMapMatStencilToGlobal(DM, PetscInt, const MatStencil[], PetscInt[]);
120: PETSC_EXTERN PetscErrorCode DMDASetOwnershipRanges(DM, const PetscInt[], const PetscInt[], const PetscInt[]);
121: PETSC_EXTERN PetscErrorCode DMDAGetOwnershipRanges(DM, const PetscInt **, const PetscInt **, const PetscInt **);
122: PETSC_EXTERN PetscErrorCode DMDASetNumProcs(DM, PetscInt, PetscInt, PetscInt);
123: PETSC_EXTERN PetscErrorCode DMDASetStencilType(DM, DMDAStencilType);
124: PETSC_EXTERN PetscErrorCode DMDAGetStencilType(DM, DMDAStencilType *);
126: PETSC_EXTERN PetscErrorCode DMDAVecGetArray(DM, Vec, void *);
127: PETSC_EXTERN PetscErrorCode DMDAVecRestoreArray(DM, Vec, void *);
128: PETSC_EXTERN PetscErrorCode DMDAVecGetArrayWrite(DM, Vec, void *);
129: PETSC_EXTERN PetscErrorCode DMDAVecRestoreArrayWrite(DM, Vec, void *);
131: PETSC_EXTERN PetscErrorCode DMDAVecGetArrayDOF(DM, Vec, void *);
132: PETSC_EXTERN PetscErrorCode DMDAVecRestoreArrayDOF(DM, Vec, void *);
134: PETSC_EXTERN PetscErrorCode DMDAVecGetArrayRead(DM, Vec, void *);
135: PETSC_EXTERN PetscErrorCode DMDAVecRestoreArrayRead(DM, Vec, void *);
137: PETSC_EXTERN PetscErrorCode DMDAVecGetArrayDOFRead(DM, Vec, void *);
138: PETSC_EXTERN PetscErrorCode DMDAVecRestoreArrayDOFRead(DM, Vec, void *);
140: PETSC_EXTERN PetscErrorCode DMDAVecGetArrayDOFWrite(DM, Vec, void *);
141: PETSC_EXTERN PetscErrorCode DMDAVecRestoreArrayDOFWrite(DM, Vec, void *);
143: PETSC_EXTERN PetscErrorCode DMDACreatePatchIS(DM, MatStencil *, MatStencil *, IS *, PetscBool);
145: /*MC
146: DMDACoor2d - Structure for holding 2d (x and y) coordinates when working with `DMDA`
148: Synopsis:
149: .vb
150: DMDACoor2d **coors;
151: Vec vcoors;
152: DM cda;
153: DMGetCoordinates(da,&vcoors);
154: DMGetCoordinateDM(da,&cda);
155: DMDAVecGetArray(cda,vcoors,&coors);
156: DMDAGetCorners(cda,&mstart,&nstart,0,&m,&n,0)
157: for (i=mstart; i<mstart+m; i++) {
158: for (j=nstart; j<nstart+n; j++) {
159: x = coors[j][i].x;
160: y = coors[j][i].y;
161: ......
162: }
163: }
164: DMDAVecRestoreArray(dac,vcoors,&coors);
165: .ve
167: Level: intermediate
169: .seealso: [](ch_dmbase), `DMDA`, `DMDACoor3d`, `DMDAVecRestoreArray()`, `DMDAVecGetArray()`, `DMGetCoordinateDM()`, `DMGetCoordinates()`
170: M*/
171: typedef struct {
172: PetscScalar x, y;
173: } DMDACoor2d;
175: /*MC
176: DMDACoor3d - Structure for holding 3d (x, y and z) coordinates coordinates when working with `DMDA`
178: Synopsis:
179: .vb
180: DMDACoor3d ***coors;
181: Vec vcoors;
182: DM cda;
183: DMGetCoordinates(da,&vcoors);
184: DMGetCoordinateDM(da,&cda);
185: DMDAVecGetArray(cda,vcoors,&coors);
186: DMDAGetCorners(cda,&mstart,&nstart,&pstart,&m,&n,&p)
187: for (i=mstart; i<mstart+m; i++) {
188: for (j=nstart; j<nstart+n; j++) {
189: for (k=pstart; k<pstart+p; k++) {
190: x = coors[k][j][i].x;
191: y = coors[k][j][i].y;
192: z = coors[k][j][i].z;
193: ......
194: }
195: }
196: DMDAVecRestoreArray(dac,vcoors,&coors);
197: .ve
199: Level: intermediate
201: .seealso: [](ch_dmbase), `DMDA`, `DMDACoor2d`, `DMDAVecRestoreArray()`, `DMDAVecGetArray()`, `DMGetCoordinateDM()`, `DMGetCoordinates()`
202: M*/
203: typedef struct {
204: PetscScalar x, y, z;
205: } DMDACoor3d;
207: PETSC_EXTERN PetscErrorCode DMDAGetLocalInfo(DM, DMDALocalInfo *);
209: PETSC_EXTERN PetscErrorCode MatRegisterDAAD(void);
210: PETSC_EXTERN PetscErrorCode MatCreateSeqUSFFT(Vec, DM, Mat *);
212: PETSC_EXTERN PetscErrorCode DMDASetGetMatrix(DM, PetscErrorCode (*)(DM, Mat *));
213: PETSC_EXTERN PetscErrorCode DMDASetBlockFills(DM, const PetscInt *, const PetscInt *);
214: PETSC_EXTERN PetscErrorCode DMDASetBlockFillsSparse(DM, const PetscInt *, const PetscInt *);
215: PETSC_EXTERN PetscErrorCode DMDASetRefinementFactor(DM, PetscInt, PetscInt, PetscInt);
216: PETSC_EXTERN PetscErrorCode DMDAGetRefinementFactor(DM, PetscInt *, PetscInt *, PetscInt *);
218: PETSC_EXTERN PetscErrorCode DMDAGetArray(DM, PetscBool, void *);
219: PETSC_EXTERN PetscErrorCode DMDARestoreArray(DM, PetscBool, void *);
221: PETSC_EXTERN PetscErrorCode DMDACreatePF(DM, PF *);
223: /* Emulation of DMPlex */
224: PETSC_EXTERN PetscErrorCode DMDAGetNumCells(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *);
225: PETSC_EXTERN PetscErrorCode DMDAGetCellPoint(DM, PetscInt, PetscInt, PetscInt, PetscInt *);
226: PETSC_EXTERN PetscErrorCode DMDAGetNumVertices(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *);
227: PETSC_EXTERN PetscErrorCode DMDAGetNumFaces(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *);
228: PETSC_EXTERN PetscErrorCode DMDAGetHeightStratum(DM, PetscInt, PetscInt *, PetscInt *);
229: PETSC_EXTERN PetscErrorCode DMDAGetDepthStratum(DM, PetscInt, PetscInt *, PetscInt *);
230: PETSC_EXTERN PetscErrorCode DMDAConvertToCell(DM, MatStencil, PetscInt *);
231: PETSC_EXTERN PetscErrorCode DMDASetVertexCoordinates(DM, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal);
232: PETSC_EXTERN PetscErrorCode DMDASetPreallocationCenterDimension(DM, PetscInt);
233: PETSC_EXTERN PetscErrorCode DMDAGetPreallocationCenterDimension(DM, PetscInt *);
235: PETSC_EXTERN PetscErrorCode DMDAVTKWriteAll(PetscObject, PetscViewer);