Actual source code: petscdmda.h
1: #if !defined(PETSCDMDA_H)
2: #define PETSCDMDA_H
4: #include <petscdm.h>
5: #include <petscdmdatypes.h>
6: #include <petscpf.h>
7: #include <petscao.h>
8: #include <petscfe.h>
10: /*MC
11: DMDA_STENCIL_STAR - "Star"-type stencil. In logical grid coordinates, only (i,j,k), (i+s,j,k), (i,j+s,k),
12: (i,j,k+s) are in the stencil NOT, for example, (i+s,j+s,k)
14: Level: beginner
16: Determines what ghost point values are brought over to each process; in this case the "corner" values are not
17: brought over and hence should not be accessed locally
19: .seealso: DMDA_STENCIL_BOX, DMDAStencilType, DMDASetStencilType()
20: M*/
22: /*MC
23: DMDA_STENCIL_BOX - "Box"-type stencil. In logical grid coordinates, any of (i,j,k), (i+s,j+r,k+t) may
24: be in the stencil.
26: Level: beginner
28: .seealso: DMDA_STENCIL_STAR, DMDAStencilType, DMDASetStencilType()
29: M*/
31: PETSC_EXTERN PetscErrorCode DMDASetInterpolationType(DM,DMDAInterpolationType);
32: PETSC_EXTERN PetscErrorCode DMDAGetInterpolationType(DM,DMDAInterpolationType*);
33: PETSC_EXTERN PetscErrorCode DMDACreateAggregates(DM,DM,Mat*);
35: /* FEM */
36: PETSC_EXTERN PetscErrorCode DMDASetElementType(DM,DMDAElementType);
37: PETSC_EXTERN PetscErrorCode DMDAGetElementType(DM,DMDAElementType*);
38: PETSC_EXTERN PetscErrorCode DMDAGetElements(DM,PetscInt*,PetscInt*,const PetscInt*[]);
39: PETSC_EXTERN PetscErrorCode DMDARestoreElements(DM,PetscInt*,PetscInt*,const PetscInt*[]);
40: PETSC_EXTERN PetscErrorCode DMDAGetElementsSizes(DM,PetscInt*,PetscInt*,PetscInt*);
41: PETSC_EXTERN PetscErrorCode DMDAGetElementsCorners(DM,PetscInt*,PetscInt*,PetscInt*);
42: PETSC_EXTERN PetscErrorCode DMDAGetSubdomainCornersIS(DM,IS*);
43: PETSC_EXTERN PetscErrorCode DMDARestoreSubdomainCornersIS(DM,IS*);
45: #define MATSEQUSFFT "sequsfft"
47: PETSC_EXTERN PetscErrorCode DMDACreate(MPI_Comm,DM*);
48: PETSC_EXTERN PetscErrorCode DMDASetSizes(DM,PetscInt,PetscInt,PetscInt);
49: PETSC_EXTERN PetscErrorCode DMDACreate1d(MPI_Comm,DMBoundaryType,PetscInt,PetscInt,PetscInt,const PetscInt[],DM *);
50: PETSC_EXTERN PetscErrorCode DMDACreate2d(MPI_Comm,DMBoundaryType,DMBoundaryType,DMDAStencilType,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],DM*);
51: 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*);
53: PETSC_EXTERN PetscErrorCode DMDAGlobalToNaturalBegin(DM,Vec,InsertMode,Vec);
54: PETSC_EXTERN PetscErrorCode DMDAGlobalToNaturalEnd(DM,Vec,InsertMode,Vec);
55: PETSC_EXTERN PetscErrorCode DMDANaturalToGlobalBegin(DM,Vec,InsertMode,Vec);
56: PETSC_EXTERN PetscErrorCode DMDANaturalToGlobalEnd(DM,Vec,InsertMode,Vec);
57: PETSC_DEPRECATED_FUNCTION("Use DMLocalToLocalBegin() (since version 3.5)") PETSC_STATIC_INLINE PetscErrorCode DMDALocalToLocalBegin(DM dm,Vec g,InsertMode mode,Vec l) {return DMLocalToLocalBegin(dm,g,mode,l);}
58: PETSC_DEPRECATED_FUNCTION("Use DMLocalToLocalEnd() (since version 3.5)") PETSC_STATIC_INLINE PetscErrorCode DMDALocalToLocalEnd(DM dm,Vec g,InsertMode mode,Vec l) {return DMLocalToLocalEnd(dm,g,mode,l);}
59: PETSC_EXTERN PetscErrorCode DMDACreateNaturalVector(DM,Vec *);
61: PETSC_EXTERN PetscErrorCode DMDAGetCorners(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
62: PETSC_EXTERN PetscErrorCode DMDAGetGhostCorners(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
63: PETSC_EXTERN PetscErrorCode DMDAGetInfo(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,DMBoundaryType*,DMBoundaryType*,DMBoundaryType*,DMDAStencilType*);
64: PETSC_EXTERN PetscErrorCode DMDAGetProcessorSubset(DM,DMDirection,PetscInt,MPI_Comm*);
65: PETSC_EXTERN PetscErrorCode DMDAGetProcessorSubsets(DM,DMDirection,MPI_Comm*);
66: PETSC_EXTERN PetscErrorCode DMDAGetRay(DM,DMDirection,PetscInt,Vec*,VecScatter*);
68: PETSC_EXTERN PetscErrorCode DMDAGlobalToNaturalAllCreate(DM,VecScatter*);
69: PETSC_EXTERN PetscErrorCode DMDANaturalAllToGlobalCreate(DM,VecScatter*);
71: PETSC_EXTERN PetscErrorCode DMDAGetScatter(DM,VecScatter*,VecScatter*);
72: PETSC_EXTERN PetscErrorCode DMDAGetNeighbors(DM,const PetscMPIInt**);
74: PETSC_EXTERN PetscErrorCode DMDASetAOType(DM,AOType);
75: PETSC_EXTERN PetscErrorCode DMDAGetAO(DM,AO*);
76: PETSC_EXTERN PetscErrorCode DMDASetUniformCoordinates(DM,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal);
77: PETSC_EXTERN PetscErrorCode DMDASetGLLCoordinates(DM,PetscInt,PetscReal*);
78: PETSC_EXTERN PetscErrorCode DMDAGetCoordinateArray(DM,void*);
79: PETSC_EXTERN PetscErrorCode DMDARestoreCoordinateArray(DM,void*);
80: PETSC_EXTERN PetscErrorCode DMDAGetLogicalCoordinate(DM,PetscScalar,PetscScalar,PetscScalar,PetscInt*,PetscInt*,PetscInt*,PetscScalar*,PetscScalar*,PetscScalar*);
81: /* function to wrap coordinates around boundary */
82: PETSC_EXTERN PetscErrorCode DMDAMapCoordsToPeriodicDomain(DM,PetscScalar*,PetscScalar*);
84: PETSC_EXTERN PetscErrorCode DMDACreateCompatibleDMDA(DM,PetscInt,DM*);
85: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use DMDACreateCompatibleDMDA() (since version 3.10)") PetscErrorCode DMDAGetReducedDMDA(DM,PetscInt,DM*);
87: PETSC_EXTERN PetscErrorCode DMDASetFieldName(DM,PetscInt,const char[]);
88: PETSC_EXTERN PetscErrorCode DMDAGetFieldName(DM,PetscInt,const char**);
89: PETSC_EXTERN PetscErrorCode DMDASetFieldNames(DM,const char * const *);
90: PETSC_EXTERN PetscErrorCode DMDAGetFieldNames(DM,const char * const **);
91: PETSC_EXTERN PetscErrorCode DMDASetCoordinateName(DM,PetscInt,const char[]);
92: PETSC_EXTERN PetscErrorCode DMDAGetCoordinateName(DM,PetscInt,const char**);
94: PETSC_EXTERN PetscErrorCode DMDASetBoundaryType(DM,DMBoundaryType,DMBoundaryType,DMBoundaryType);
95: PETSC_EXTERN PetscErrorCode DMDASetDof(DM, PetscInt);
96: PETSC_EXTERN PetscErrorCode DMDAGetDof(DM, PetscInt*);
97: PETSC_EXTERN PetscErrorCode DMDASetOverlap(DM,PetscInt,PetscInt,PetscInt);
98: PETSC_EXTERN PetscErrorCode DMDAGetOverlap(DM,PetscInt*,PetscInt*,PetscInt*);
99: PETSC_EXTERN PetscErrorCode DMDASetNumLocalSubDomains(DM,PetscInt);
100: PETSC_EXTERN PetscErrorCode DMDAGetNumLocalSubDomains(DM,PetscInt*);
101: PETSC_EXTERN PetscErrorCode DMDAGetOffset(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
102: PETSC_EXTERN PetscErrorCode DMDASetOffset(DM,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt);
103: PETSC_EXTERN PetscErrorCode DMDAGetNonOverlappingRegion(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
104: PETSC_EXTERN PetscErrorCode DMDASetNonOverlappingRegion(DM,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt);
105: PETSC_EXTERN PetscErrorCode DMDASetStencilWidth(DM, PetscInt);
106: PETSC_EXTERN PetscErrorCode DMDAGetStencilWidth(DM, PetscInt*);
107: PETSC_EXTERN PetscErrorCode DMDASetOwnershipRanges(DM,const PetscInt[],const PetscInt[],const PetscInt[]);
108: PETSC_EXTERN PetscErrorCode DMDAGetOwnershipRanges(DM,const PetscInt**,const PetscInt**,const PetscInt**);
109: PETSC_EXTERN PetscErrorCode DMDASetNumProcs(DM, PetscInt, PetscInt, PetscInt);
110: PETSC_EXTERN PetscErrorCode DMDASetStencilType(DM, DMDAStencilType);
111: PETSC_EXTERN PetscErrorCode DMDAGetStencilType(DM, DMDAStencilType*);
113: PETSC_EXTERN PetscErrorCode DMDAVecGetArray(DM,Vec,void *);
114: PETSC_EXTERN PetscErrorCode DMDAVecRestoreArray(DM,Vec,void *);
115: PETSC_EXTERN PetscErrorCode DMDAVecGetArrayWrite(DM,Vec,void *);
116: PETSC_EXTERN PetscErrorCode DMDAVecRestoreArrayWrite(DM,Vec,void *);
118: PETSC_EXTERN PetscErrorCode DMDAVecGetArrayDOF(DM,Vec,void *);
119: PETSC_EXTERN PetscErrorCode DMDAVecRestoreArrayDOF(DM,Vec,void *);
121: PETSC_EXTERN PetscErrorCode DMDAVecGetArrayRead(DM,Vec,void *);
122: PETSC_EXTERN PetscErrorCode DMDAVecRestoreArrayRead(DM,Vec,void *);
124: PETSC_EXTERN PetscErrorCode DMDAVecGetArrayDOFRead(DM,Vec,void *);
125: PETSC_EXTERN PetscErrorCode DMDAVecRestoreArrayDOFRead(DM,Vec,void *);
127: PETSC_EXTERN PetscErrorCode DMDAVecGetArrayDOFWrite(DM,Vec,void *);
128: PETSC_EXTERN PetscErrorCode DMDAVecRestoreArrayDOFWrite(DM,Vec,void *);
130: PETSC_EXTERN PetscErrorCode DMDACreatePatchIS(DM,MatStencil*,MatStencil*,IS*,PetscBool);
132: /*MC
133: DMDACoor2d - Structure for holding 2d (x and y) coordinates.
135: Level: intermediate
137: Synopsis:
138: DMDACoor2d **coors;
139: Vec vcoors;
140: DM cda;
141: DMGetCoordinates(da,&vcoors);
142: DMGetCoordinateDM(da,&cda);
143: DMDAVecGetArray(cda,vcoors,&coors);
144: DMDAGetCorners(cda,&mstart,&nstart,0,&m,&n,0)
145: for (i=mstart; i<mstart+m; i++) {
146: for (j=nstart; j<nstart+n; j++) {
147: x = coors[j][i].x;
148: y = coors[j][i].y;
149: ......
150: }
151: }
152: DMDAVecRestoreArray(dac,vcoors,&coors);
154: .seealso: DMDACoor3d, DMGetCoordinateDM(), DMGetCoordinates()
155: M*/
156: typedef struct {PetscScalar x,y;} DMDACoor2d;
158: /*MC
159: DMDACoor3d - Structure for holding 3d (x, y and z) coordinates.
161: Level: intermediate
163: Synopsis:
164: DMDACoor3d ***coors;
165: Vec vcoors;
166: DM cda;
167: DMGetCoordinates(da,&vcoors);
168: DMGetCoordinateDM(da,&cda);
169: DMDAVecGetArray(cda,vcoors,&coors);
170: DMDAGetCorners(cda,&mstart,&nstart,&pstart,&m,&n,&p)
171: for (i=mstart; i<mstart+m; i++) {
172: for (j=nstart; j<nstart+n; j++) {
173: for (k=pstart; k<pstart+p; k++) {
174: x = coors[k][j][i].x;
175: y = coors[k][j][i].y;
176: z = coors[k][j][i].z;
177: ......
178: }
179: }
180: DMDAVecRestoreArray(dac,vcoors,&coors);
182: .seealso: DMDACoor2d, DMGetCoordinateDM(), DMGetCoordinates()
183: M*/
184: typedef struct {PetscScalar x,y,z;} DMDACoor3d;
186: PETSC_EXTERN PetscErrorCode DMDAGetLocalInfo(DM,DMDALocalInfo*);
188: PETSC_EXTERN PetscErrorCode MatRegisterDAAD(void);
189: PETSC_EXTERN PetscErrorCode MatCreateDAAD(DM,Mat*);
190: PETSC_EXTERN PetscErrorCode MatCreateSeqUSFFT(Vec,DM,Mat*);
192: PETSC_EXTERN PetscErrorCode DMDASetGetMatrix(DM,PetscErrorCode (*)(DM, Mat *));
193: PETSC_EXTERN PetscErrorCode DMDASetBlockFills(DM,const PetscInt*,const PetscInt*);
194: PETSC_EXTERN PetscErrorCode DMDASetBlockFillsSparse(DM,const PetscInt*,const PetscInt*);
195: PETSC_EXTERN PetscErrorCode DMDASetRefinementFactor(DM,PetscInt,PetscInt,PetscInt);
196: PETSC_EXTERN PetscErrorCode DMDAGetRefinementFactor(DM,PetscInt*,PetscInt*,PetscInt*);
198: PETSC_EXTERN PetscErrorCode DMDAGetArray(DM,PetscBool ,void*);
199: PETSC_EXTERN PetscErrorCode DMDARestoreArray(DM,PetscBool ,void*);
201: PETSC_EXTERN PetscErrorCode DMDACreatePF(DM,PF*);
203: /* Emulation of DMPlex */
204: PETSC_EXTERN PetscErrorCode DMDAGetNumCells(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *);
205: PETSC_EXTERN PetscErrorCode DMDAGetCellPoint(DM, PetscInt, PetscInt, PetscInt, PetscInt *);
206: PETSC_EXTERN PetscErrorCode DMDAGetNumVertices(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *);
207: PETSC_EXTERN PetscErrorCode DMDAGetNumFaces(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *);
208: PETSC_EXTERN PetscErrorCode DMDAGetHeightStratum(DM, PetscInt, PetscInt *, PetscInt *);
209: PETSC_EXTERN PetscErrorCode DMDAGetDepthStratum(DM, PetscInt, PetscInt *, PetscInt *);
210: PETSC_EXTERN PetscErrorCode DMDAComputeCellGeometryFEM(DM, PetscInt, PetscQuadrature, PetscReal [], PetscReal [], PetscReal [], PetscReal []);
211: PETSC_EXTERN PetscErrorCode DMDAGetTransitiveClosure(DM, PetscInt, PetscBool, PetscInt *, PetscInt **);
212: PETSC_EXTERN PetscErrorCode DMDARestoreTransitiveClosure(DM, PetscInt, PetscBool, PetscInt *, PetscInt **);
213: PETSC_EXTERN PetscErrorCode DMDAConvertToCell(DM, MatStencil, PetscInt *);
214: PETSC_EXTERN PetscErrorCode DMDASetVertexCoordinates(DM,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal);
215: PETSC_EXTERN PetscErrorCode DMDASetPreallocationCenterDimension(DM, PetscInt);
216: PETSC_EXTERN PetscErrorCode DMDAGetPreallocationCenterDimension(DM, PetscInt*);
218: #endif