Actual source code: petscdmda.h
petsc-3.4.5 2014-06-29
4: #include <petscdm.h>
5: #include <petscdmdatypes.h>
6: #include <petscpf.h>
7: #include <petscao.h>
9: /*MC
10: DMDA_STENCIL_STAR - "Star"-type stencil. In logical grid coordinates, only (i,j,k), (i+s,j,k), (i,j+s,k),
11: (i,j,k+s) are in the stencil NOT, for example, (i+s,j+s,k)
13: Level: beginner
15: .seealso: DMDA_STENCIL_BOX, DMDAStencilType, DMDASetStencilType()
16: M*/
18: /*MC
19: DMDA_STENCIL_BOX - "Box"-type stencil. In logical grid coordinates, any of (i,j,k), (i+s,j+r,k+t) may
20: be in the stencil.
22: Level: beginner
24: .seealso: DMDA_STENCIL_STAR, DMDAStencilType, DMDASetStencilType()
25: M*/
27: PETSC_EXTERN const char *const DMDABoundaryTypes[];
29: PETSC_EXTERN PetscErrorCode DMDASetInterpolationType(DM,DMDAInterpolationType);
30: PETSC_EXTERN PetscErrorCode DMDAGetInterpolationType(DM,DMDAInterpolationType*);
32: PETSC_EXTERN PetscErrorCode DMDASetElementType(DM,DMDAElementType);
33: PETSC_EXTERN PetscErrorCode DMDAGetElementType(DM,DMDAElementType*);
34: PETSC_EXTERN PetscErrorCode DMDAGetElements(DM,PetscInt *,PetscInt *,const PetscInt*[]);
35: PETSC_EXTERN PetscErrorCode DMDARestoreElements(DM,PetscInt *,PetscInt *,const PetscInt*[]);
37: typedef enum { DMDA_X,DMDA_Y,DMDA_Z } DMDADirection;
39: #define MATSEQUSFFT "sequsfft"
41: PETSC_EXTERN PetscErrorCode DMDACreate(MPI_Comm,DM*);
42: PETSC_EXTERN PetscErrorCode DMDASetDim(DM,PetscInt);
43: PETSC_EXTERN PetscErrorCode DMDASetSizes(DM,PetscInt,PetscInt,PetscInt);
44: PETSC_EXTERN PetscErrorCode DMDACreate1d(MPI_Comm,DMDABoundaryType,PetscInt,PetscInt,PetscInt,const PetscInt[],DM *);
45: PETSC_EXTERN PetscErrorCode DMDACreate2d(MPI_Comm,DMDABoundaryType,DMDABoundaryType,DMDAStencilType,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],DM*);
46: PETSC_EXTERN PetscErrorCode DMDACreate3d(MPI_Comm,DMDABoundaryType,DMDABoundaryType,DMDABoundaryType,DMDAStencilType,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscInt[],DM*);
48: PETSC_EXTERN PetscErrorCode DMDAGlobalToNaturalBegin(DM,Vec,InsertMode,Vec);
49: PETSC_EXTERN PetscErrorCode DMDAGlobalToNaturalEnd(DM,Vec,InsertMode,Vec);
50: PETSC_EXTERN PetscErrorCode DMDANaturalToGlobalBegin(DM,Vec,InsertMode,Vec);
51: PETSC_EXTERN PetscErrorCode DMDANaturalToGlobalEnd(DM,Vec,InsertMode,Vec);
52: PETSC_EXTERN PetscErrorCode DMDALocalToLocalBegin(DM,Vec,InsertMode,Vec);
53: PETSC_EXTERN PetscErrorCode DMDALocalToLocalEnd(DM,Vec,InsertMode,Vec);
54: PETSC_EXTERN PetscErrorCode DMDACreateNaturalVector(DM,Vec *);
56: PETSC_EXTERN PetscErrorCode DMDAGetCorners(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
57: PETSC_EXTERN PetscErrorCode DMDAGetGhostCorners(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
58: PETSC_EXTERN PetscErrorCode DMDAGetInfo(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,DMDABoundaryType*,DMDABoundaryType*,DMDABoundaryType*,DMDAStencilType*);
59: PETSC_EXTERN PetscErrorCode DMDAGetProcessorSubset(DM,DMDADirection,PetscInt,MPI_Comm*);
60: PETSC_EXTERN PetscErrorCode DMDAGetProcessorSubsets(DM,DMDADirection,MPI_Comm*);
61: PETSC_EXTERN PetscErrorCode DMDAGetRay(DM,DMDADirection,PetscInt,Vec*,VecScatter*);
63: PETSC_EXTERN PetscErrorCode DMDAGlobalToNaturalAllCreate(DM,VecScatter*);
64: PETSC_EXTERN PetscErrorCode DMDANaturalAllToGlobalCreate(DM,VecScatter*);
66: PETSC_EXTERN PetscErrorCode DMDAGetGlobalIndices(DM,PetscInt*,PetscInt**);
68: PETSC_EXTERN PetscErrorCode DMDAGetScatter(DM,VecScatter*,VecScatter*,VecScatter*);
69: PETSC_EXTERN PetscErrorCode DMDAGetNeighbors(DM,const PetscMPIInt**);
71: PETSC_EXTERN PetscErrorCode DMDAGetAO(DM,AO*);
72: PETSC_EXTERN PetscErrorCode DMDASetUniformCoordinates(DM,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal);
73: PETSC_EXTERN PetscErrorCode DMDAGetBoundingBox(DM,PetscReal[],PetscReal[]);
74: PETSC_EXTERN PetscErrorCode DMDAGetLocalBoundingBox(DM,PetscReal[],PetscReal[]);
75: PETSC_EXTERN PetscErrorCode DMDAGetLogicalCoordinate(DM,PetscScalar,PetscScalar,PetscScalar,PetscInt*,PetscInt*,PetscInt*,PetscScalar*,PetscScalar*,PetscScalar*);
76: /* function to wrap coordinates around boundary */
77: PETSC_EXTERN PetscErrorCode DMDAMapCoordsToPeriodicDomain(DM,PetscScalar*,PetscScalar*);
79: PETSC_EXTERN PetscErrorCode DMDAGetReducedDMDA(DM,PetscInt,DM*);
81: PETSC_EXTERN PetscErrorCode DMDASetFieldName(DM,PetscInt,const char[]);
82: PETSC_EXTERN PetscErrorCode DMDAGetFieldName(DM,PetscInt,const char**);
83: PETSC_EXTERN PetscErrorCode DMDASetCoordinateName(DM,PetscInt,const char[]);
84: PETSC_EXTERN PetscErrorCode DMDAGetCoordinateName(DM,PetscInt,const char**);
86: PETSC_EXTERN PetscErrorCode DMDASetBoundaryType(DM,DMDABoundaryType,DMDABoundaryType,DMDABoundaryType);
87: PETSC_EXTERN PetscErrorCode DMDASetDof(DM, PetscInt);
88: PETSC_EXTERN PetscErrorCode DMDASetOverlap(DM,PetscInt,PetscInt,PetscInt);
89: PETSC_EXTERN PetscErrorCode DMDAGetOverlap(DM,PetscInt*,PetscInt*,PetscInt*);
90: PETSC_EXTERN PetscErrorCode DMDASetNumLocalSubDomains(DM,PetscInt);
91: PETSC_EXTERN PetscErrorCode DMDAGetNumLocalSubDomains(DM,PetscInt*);
92: PETSC_EXTERN PetscErrorCode DMDAGetOffset(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
93: PETSC_EXTERN PetscErrorCode DMDASetOffset(DM,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt);
94: PETSC_EXTERN PetscErrorCode DMDAGetNonOverlappingRegion(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
95: PETSC_EXTERN PetscErrorCode DMDASetNonOverlappingRegion(DM,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt);
96: PETSC_EXTERN PetscErrorCode DMDASetStencilWidth(DM, PetscInt);
97: PETSC_EXTERN PetscErrorCode DMDASetOwnershipRanges(DM,const PetscInt[],const PetscInt[],const PetscInt[]);
98: PETSC_EXTERN PetscErrorCode DMDAGetOwnershipRanges(DM,const PetscInt**,const PetscInt**,const PetscInt**);
99: PETSC_EXTERN PetscErrorCode DMDASetNumProcs(DM, PetscInt, PetscInt, PetscInt);
100: PETSC_EXTERN PetscErrorCode DMDASetStencilType(DM, DMDAStencilType);
102: PETSC_EXTERN PetscErrorCode DMDAVecGetArray(DM,Vec,void *);
103: PETSC_EXTERN PetscErrorCode DMDAVecRestoreArray(DM,Vec,void *);
105: PETSC_EXTERN PetscErrorCode DMDAVecGetArrayDOF(DM,Vec,void *);
106: PETSC_EXTERN PetscErrorCode DMDAVecRestoreArrayDOF(DM,Vec,void *);
108: PETSC_EXTERN PetscErrorCode DMDASplitComm2d(MPI_Comm,PetscInt,PetscInt,PetscInt,MPI_Comm*);
110: PETSC_EXTERN PetscErrorCode DMDACreatePatchIS(DM,MatStencil*,MatStencil*,IS*);
113: /*MC
114: DMDACoor2d - Structure for holding 2d (x and y) coordinates.
116: Level: intermediate
118: Sample Usage:
119: DMDACoor2d **coors;
120: Vec vcoors;
121: DM cda;
123: DMGetCoordinates(da,&vcoors);
124: DMDAGetCoordinateDA(da,&cda);
125: DMDAVecGetArray(cda,vcoors,&coors);
126: DMDAGetCorners(cda,&mstart,&nstart,0,&m,&n,0)
127: for (i=mstart; i<mstart+m; i++) {
128: for (j=nstart; j<nstart+n; j++) {
129: x = coors[j][i].x;
130: y = coors[j][i].y;
131: ......
132: }
133: }
134: DMDAVecRestoreArray(dac,vcoors,&coors);
136: .seealso: DMDACoor3d, DMDAForEachPointBegin(), DMDAGetCoordinateDA(), DMGetCoordinates(), DMDAGetGhostCoordinates()
137: M*/
138: typedef struct {PetscScalar x,y;} DMDACoor2d;
140: /*MC
141: DMDACoor3d - Structure for holding 3d (x, y and z) coordinates.
143: Level: intermediate
145: Sample Usage:
146: DMDACoor3d ***coors;
147: Vec vcoors;
148: DM cda;
150: DMGetCoordinates(da,&vcoors);
151: DMDAGetCoordinateDA(da,&cda);
152: DMDAVecGetArray(cda,vcoors,&coors);
153: DMDAGetCorners(cda,&mstart,&nstart,&pstart,&m,&n,&p)
154: for (i=mstart; i<mstart+m; i++) {
155: for (j=nstart; j<nstart+n; j++) {
156: for (k=pstart; k<pstart+p; k++) {
157: x = coors[k][j][i].x;
158: y = coors[k][j][i].y;
159: z = coors[k][j][i].z;
160: ......
161: }
162: }
163: DMDAVecRestoreArray(dac,vcoors,&coors);
165: .seealso: DMDACoor2d, DMDAForEachPointBegin(), DMDAGetCoordinateDA(), DMGetCoordinates(), DMDAGetGhostCoordinates()
166: M*/
167: typedef struct {PetscScalar x,y,z;} DMDACoor3d;
169: PETSC_EXTERN PetscErrorCode DMDAGetLocalInfo(DM,DMDALocalInfo*);
171: PETSC_EXTERN PetscErrorCode MatRegisterDAAD(void);
172: PETSC_EXTERN PetscErrorCode MatCreateDAAD(DM,Mat*);
173: PETSC_EXTERN PetscErrorCode MatCreateSeqUSFFT(Vec,DM,Mat*);
175: PETSC_EXTERN PetscErrorCode DMDASetGetMatrix(DM,PetscErrorCode (*)(DM, MatType,Mat *));
176: PETSC_EXTERN PetscErrorCode DMDASetBlockFills(DM,const PetscInt*,const PetscInt*);
177: PETSC_EXTERN PetscErrorCode DMDASetRefinementFactor(DM,PetscInt,PetscInt,PetscInt);
178: PETSC_EXTERN PetscErrorCode DMDAGetRefinementFactor(DM,PetscInt*,PetscInt*,PetscInt*);
180: PETSC_EXTERN PetscErrorCode DMDAGetArray(DM,PetscBool ,void*);
181: PETSC_EXTERN PetscErrorCode DMDARestoreArray(DM,PetscBool ,void*);
183: PETSC_EXTERN PetscErrorCode DMDACreatePF(DM,PF*);
185: PETSC_EXTERN PetscErrorCode DMDAGetNumCells(DM, PetscInt *);
186: PETSC_EXTERN PetscErrorCode DMDAGetNumVertices(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *);
187: PETSC_EXTERN PetscErrorCode DMDAGetNumFaces(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *);
188: PETSC_EXTERN PetscErrorCode DMDAGetHeightStratum(DM, PetscInt, PetscInt *, PetscInt *);
189: PETSC_EXTERN PetscErrorCode DMDACreateSection(DM, PetscInt[], PetscInt[], PetscInt[], PetscInt[]);
190: PETSC_EXTERN PetscErrorCode DMDAComputeCellGeometry(DM, PetscInt, PetscQuadrature *, PetscReal [], PetscReal [], PetscReal [], PetscReal []);
191: PETSC_EXTERN PetscErrorCode DMDAVecGetClosure(DM, PetscSection, Vec, PetscInt, const PetscScalar **);
192: PETSC_EXTERN PetscErrorCode DMDAVecSetClosure(DM, PetscSection, Vec, PetscInt, const PetscScalar *, InsertMode);
193: PETSC_EXTERN PetscErrorCode DMDAGetClosure(DM,PetscSection,PetscInt,PetscInt*,const PetscInt**);
194: PETSC_EXTERN PetscErrorCode DMDARestoreClosure(DM,PetscSection,PetscInt,PetscInt*,const PetscInt**);
195: PETSC_EXTERN PetscErrorCode DMDAGetClosureScalar(DM,PetscSection,PetscInt,PetscScalar*,const PetscScalar**);
196: PETSC_EXTERN PetscErrorCode DMDARestoreClosureScalar(DM,PetscSection,PetscInt,PetscScalar*,const PetscScalar**);
197: PETSC_EXTERN PetscErrorCode DMDASetClosureScalar(DM,PetscSection,PetscInt,PetscScalar*,const PetscScalar*,InsertMode);
198: PETSC_EXTERN PetscErrorCode DMDAConvertToCell(DM, MatStencil, PetscInt *);
200: #endif