Actual source code: dacreate.c
2: #include <petsc/private/dmdaimpl.h>
4: PetscErrorCode DMSetFromOptions_DA(PetscOptionItems *PetscOptionsObject,DM da)
5: {
7: DM_DA *dd = (DM_DA*)da->data;
8: PetscInt refine = 0,dim = da->dim,maxnlevels = 100,refx[100],refy[100],refz[100],n,i;
9: PetscBool flg;
14: if (dd->M < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"Dimension must be non-negative, call DMSetFromOptions() if you want to change the value at runtime");
15: if (dd->N < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"Dimension must be non-negative, call DMSetFromOptions() if you want to change the value at runtime");
16: if (dd->P < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"Dimension must be non-negative, call DMSetFromOptions() if you want to change the value at runtime");
18: PetscOptionsHead(PetscOptionsObject,"DMDA Options");
19: PetscOptionsBoundedInt("-da_grid_x","Number of grid points in x direction","DMDASetSizes",dd->M,&dd->M,NULL,1);
20: if (dim > 1) {PetscOptionsBoundedInt("-da_grid_y","Number of grid points in y direction","DMDASetSizes",dd->N,&dd->N,NULL,1);}
21: if (dim > 2) {PetscOptionsBoundedInt("-da_grid_z","Number of grid points in z direction","DMDASetSizes",dd->P,&dd->P,NULL,1);}
23: PetscOptionsBoundedInt("-da_overlap","Decomposition overlap in all directions","DMDASetOverlap",dd->xol,&dd->xol,&flg,0);
24: if (flg) {DMDASetOverlap(da,dd->xol,dd->xol,dd->xol);}
25: PetscOptionsBoundedInt("-da_overlap_x","Decomposition overlap in x direction","DMDASetOverlap",dd->xol,&dd->xol,NULL,0);
26: if (dim > 1) {PetscOptionsBoundedInt("-da_overlap_y","Decomposition overlap in y direction","DMDASetOverlap",dd->yol,&dd->yol,NULL,0);}
27: if (dim > 2) {PetscOptionsBoundedInt("-da_overlap_z","Decomposition overlap in z direction","DMDASetOverlap",dd->zol,&dd->zol,NULL,0);}
29: PetscOptionsBoundedInt("-da_local_subdomains","","DMDASetNumLocalSubdomains",dd->Nsub,&dd->Nsub,&flg,PETSC_DECIDE);
30: if (flg) {DMDASetNumLocalSubDomains(da,dd->Nsub);}
32: /* Handle DMDA parallel distribution */
33: PetscOptionsBoundedInt("-da_processors_x","Number of processors in x direction","DMDASetNumProcs",dd->m,&dd->m,NULL,PETSC_DECIDE);
34: if (dim > 1) {PetscOptionsBoundedInt("-da_processors_y","Number of processors in y direction","DMDASetNumProcs",dd->n,&dd->n,NULL,PETSC_DECIDE);}
35: if (dim > 2) {PetscOptionsBoundedInt("-da_processors_z","Number of processors in z direction","DMDASetNumProcs",dd->p,&dd->p,NULL,PETSC_DECIDE);}
36: /* Handle DMDA refinement */
37: PetscOptionsBoundedInt("-da_refine_x","Refinement ratio in x direction","DMDASetRefinementFactor",dd->refine_x,&dd->refine_x,NULL,1);
38: if (dim > 1) {PetscOptionsBoundedInt("-da_refine_y","Refinement ratio in y direction","DMDASetRefinementFactor",dd->refine_y,&dd->refine_y,NULL,1);}
39: if (dim > 2) {PetscOptionsBoundedInt("-da_refine_z","Refinement ratio in z direction","DMDASetRefinementFactor",dd->refine_z,&dd->refine_z,NULL,1);}
40: dd->coarsen_x = dd->refine_x; dd->coarsen_y = dd->refine_y; dd->coarsen_z = dd->refine_z;
42: /* Get refinement factors, defaults taken from the coarse DMDA */
43: DMDAGetRefinementFactor(da,&refx[0],&refy[0],&refz[0]);
44: for (i=1; i<maxnlevels; i++) {
45: refx[i] = refx[0];
46: refy[i] = refy[0];
47: refz[i] = refz[0];
48: }
49: n = maxnlevels;
50: PetscOptionsIntArray("-da_refine_hierarchy_x","Refinement factor for each level","None",refx,&n,&flg);
51: if (flg) {
52: dd->refine_x = refx[0];
53: dd->refine_x_hier_n = n;
54: PetscMalloc1(n,&dd->refine_x_hier);
55: PetscArraycpy(dd->refine_x_hier,refx,n);
56: }
57: if (dim > 1) {
58: n = maxnlevels;
59: PetscOptionsIntArray("-da_refine_hierarchy_y","Refinement factor for each level","None",refy,&n,&flg);
60: if (flg) {
61: dd->refine_y = refy[0];
62: dd->refine_y_hier_n = n;
63: PetscMalloc1(n,&dd->refine_y_hier);
64: PetscArraycpy(dd->refine_y_hier,refy,n);
65: }
66: }
67: if (dim > 2) {
68: n = maxnlevels;
69: PetscOptionsIntArray("-da_refine_hierarchy_z","Refinement factor for each level","None",refz,&n,&flg);
70: if (flg) {
71: dd->refine_z = refz[0];
72: dd->refine_z_hier_n = n;
73: PetscMalloc1(n,&dd->refine_z_hier);
74: PetscArraycpy(dd->refine_z_hier,refz,n);
75: }
76: }
78: PetscOptionsBoundedInt("-da_refine","Uniformly refine DA one or more times","None",refine,&refine,NULL,0);
79: PetscOptionsTail();
81: while (refine--) {
82: if (dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
83: dd->M = dd->refine_x*dd->M;
84: } else {
85: dd->M = 1 + dd->refine_x*(dd->M - 1);
86: }
87: if (dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
88: dd->N = dd->refine_y*dd->N;
89: } else {
90: dd->N = 1 + dd->refine_y*(dd->N - 1);
91: }
92: if (dd->bz == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
93: dd->P = dd->refine_z*dd->P;
94: } else {
95: dd->P = 1 + dd->refine_z*(dd->P - 1);
96: }
97: da->levelup++;
98: if (da->levelup - da->leveldown >= 0) {
99: dd->refine_x = refx[da->levelup - da->leveldown];
100: dd->refine_y = refy[da->levelup - da->leveldown];
101: dd->refine_z = refz[da->levelup - da->leveldown];
102: }
103: if (da->levelup - da->leveldown >= 1) {
104: dd->coarsen_x = refx[da->levelup - da->leveldown - 1];
105: dd->coarsen_y = refy[da->levelup - da->leveldown - 1];
106: dd->coarsen_z = refz[da->levelup - da->leveldown - 1];
107: }
108: }
109: return(0);
110: }
112: extern PetscErrorCode DMCreateGlobalVector_DA(DM,Vec*);
113: extern PetscErrorCode DMCreateLocalVector_DA(DM,Vec*);
114: extern PetscErrorCode DMGlobalToLocalBegin_DA(DM,Vec,InsertMode,Vec);
115: extern PetscErrorCode DMGlobalToLocalEnd_DA(DM,Vec,InsertMode,Vec);
116: extern PetscErrorCode DMLocalToGlobalBegin_DA(DM,Vec,InsertMode,Vec);
117: extern PetscErrorCode DMLocalToGlobalEnd_DA(DM,Vec,InsertMode,Vec);
118: extern PetscErrorCode DMLocalToLocalBegin_DA(DM,Vec,InsertMode,Vec);
119: extern PetscErrorCode DMLocalToLocalEnd_DA(DM,Vec,InsertMode,Vec);
120: extern PetscErrorCode DMCreateInterpolation_DA(DM,DM,Mat*,Vec*);
121: extern PetscErrorCode DMCreateColoring_DA(DM,ISColoringType,ISColoring*);
122: extern PetscErrorCode DMCreateMatrix_DA(DM,Mat*);
123: extern PetscErrorCode DMCreateCoordinateDM_DA(DM,DM*);
124: extern PetscErrorCode DMRefine_DA(DM,MPI_Comm,DM*);
125: extern PetscErrorCode DMCoarsen_DA(DM,MPI_Comm,DM*);
126: extern PetscErrorCode DMRefineHierarchy_DA(DM,PetscInt,DM[]);
127: extern PetscErrorCode DMCoarsenHierarchy_DA(DM,PetscInt,DM[]);
128: extern PetscErrorCode DMCreateInjection_DA(DM,DM,Mat*);
129: extern PetscErrorCode DMView_DA(DM,PetscViewer);
130: extern PetscErrorCode DMSetUp_DA(DM);
131: extern PetscErrorCode DMDestroy_DA(DM);
132: extern PetscErrorCode DMCreateDomainDecomposition_DA(DM,PetscInt*,char***,IS**,IS**,DM**);
133: extern PetscErrorCode DMCreateDomainDecompositionScatters_DA(DM,PetscInt,DM*,VecScatter**,VecScatter**,VecScatter**);
134: PETSC_INTERN PetscErrorCode DMGetCompatibility_DA(DM,DM,PetscBool*,PetscBool*);
136: PetscErrorCode DMLoad_DA(DM da,PetscViewer viewer)
137: {
138: PetscErrorCode ierr;
139: PetscInt dim,m,n,p,dof,swidth;
140: DMDAStencilType stencil;
141: DMBoundaryType bx,by,bz;
142: PetscBool coors;
143: DM dac;
144: Vec c;
147: PetscViewerBinaryRead(viewer,&dim,1,NULL,PETSC_INT);
148: PetscViewerBinaryRead(viewer,&m,1,NULL,PETSC_INT);
149: PetscViewerBinaryRead(viewer,&n,1,NULL,PETSC_INT);
150: PetscViewerBinaryRead(viewer,&p,1,NULL,PETSC_INT);
151: PetscViewerBinaryRead(viewer,&dof,1,NULL,PETSC_INT);
152: PetscViewerBinaryRead(viewer,&swidth,1,NULL,PETSC_INT);
153: PetscViewerBinaryRead(viewer,&bx,1,NULL,PETSC_ENUM);
154: PetscViewerBinaryRead(viewer,&by,1,NULL,PETSC_ENUM);
155: PetscViewerBinaryRead(viewer,&bz,1,NULL,PETSC_ENUM);
156: PetscViewerBinaryRead(viewer,&stencil,1,NULL,PETSC_ENUM);
158: DMSetDimension(da, dim);
159: DMDASetSizes(da, m,n,p);
160: DMDASetBoundaryType(da, bx, by, bz);
161: DMDASetDof(da, dof);
162: DMDASetStencilType(da, stencil);
163: DMDASetStencilWidth(da, swidth);
164: DMSetUp(da);
165: PetscViewerBinaryRead(viewer,&coors,1,NULL,PETSC_ENUM);
166: if (coors) {
167: DMGetCoordinateDM(da,&dac);
168: DMCreateGlobalVector(dac,&c);
169: VecLoad(c,viewer);
170: DMSetCoordinates(da,c);
171: VecDestroy(&c);
172: }
173: return(0);
174: }
176: PetscErrorCode DMCreateSubDM_DA(DM dm, PetscInt numFields, const PetscInt fields[], IS *is, DM *subdm)
177: {
178: DM_DA *da = (DM_DA*) dm->data;
182: if (subdm) {
183: PetscSF sf;
184: Vec coords;
185: void *ctx;
186: /* Cannot use DMClone since the dof stuff is mixed in. Ugh
187: DMClone(dm, subdm); */
188: DMCreate(PetscObjectComm((PetscObject)dm), subdm);
189: DMGetPointSF(dm, &sf);
190: DMSetPointSF(*subdm, sf);
191: DMGetApplicationContext(dm, &ctx);
192: DMSetApplicationContext(*subdm, ctx);
193: DMGetCoordinatesLocal(dm, &coords);
194: if (coords) {
195: DMSetCoordinatesLocal(*subdm, coords);
196: } else {
197: DMGetCoordinates(dm, &coords);
198: if (coords) {DMSetCoordinates(*subdm, coords);}
199: }
201: DMSetType(*subdm, DMDA);
202: DMSetDimension(*subdm, dm->dim);
203: DMDASetSizes(*subdm, da->M, da->N, da->P);
204: DMDASetNumProcs(*subdm, da->m, da->n, da->p);
205: DMDASetBoundaryType(*subdm, da->bx, da->by, da->bz);
206: DMDASetDof(*subdm, numFields);
207: DMDASetStencilType(*subdm, da->stencil_type);
208: DMDASetStencilWidth(*subdm, da->s);
209: DMDASetOwnershipRanges(*subdm, da->lx, da->ly, da->lz);
210: }
211: if (is) {
212: PetscInt *indices, cnt = 0, dof = da->w, i, j;
214: PetscMalloc1(da->Nlocal*numFields/dof, &indices);
215: for (i = da->base/dof; i < (da->base+da->Nlocal)/dof; ++i) {
216: for (j = 0; j < numFields; ++j) {
217: indices[cnt++] = dof*i + fields[j];
218: }
219: }
220: if (cnt != da->Nlocal*numFields/dof) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Count %D does not equal expected value %D", cnt, da->Nlocal*numFields/dof);
221: ISCreateGeneral(PetscObjectComm((PetscObject) dm), cnt, indices, PETSC_OWN_POINTER, is);
222: }
223: return(0);
224: }
226: PetscErrorCode DMCreateFieldDecomposition_DA(DM dm, PetscInt *len,char ***namelist, IS **islist, DM **dmlist)
227: {
228: PetscInt i;
230: DM_DA *dd = (DM_DA*)dm->data;
231: PetscInt dof = dd->w;
234: if (len) *len = dof;
235: if (islist) {
236: Vec v;
237: PetscInt rstart,n;
239: DMGetGlobalVector(dm,&v);
240: VecGetOwnershipRange(v,&rstart,NULL);
241: VecGetLocalSize(v,&n);
242: DMRestoreGlobalVector(dm,&v);
243: PetscMalloc1(dof,islist);
244: for (i=0; i<dof; i++) {
245: ISCreateStride(PetscObjectComm((PetscObject)dm),n/dof,rstart+i,dof,&(*islist)[i]);
246: }
247: }
248: if (namelist) {
249: PetscMalloc1(dof, namelist);
250: if (dd->fieldname) {
251: for (i=0; i<dof; i++) {
252: PetscStrallocpy(dd->fieldname[i],&(*namelist)[i]);
253: }
254: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Currently DMDA must have fieldnames");
255: }
256: if (dmlist) {
257: DM da;
259: DMDACreate(PetscObjectComm((PetscObject)dm), &da);
260: DMSetDimension(da, dm->dim);
261: DMDASetSizes(da, dd->M, dd->N, dd->P);
262: DMDASetNumProcs(da, dd->m, dd->n, dd->p);
263: DMDASetBoundaryType(da, dd->bx, dd->by, dd->bz);
264: DMDASetDof(da, 1);
265: DMDASetStencilType(da, dd->stencil_type);
266: DMDASetStencilWidth(da, dd->s);
267: DMSetUp(da);
268: PetscMalloc1(dof,dmlist);
269: for (i=0; i<dof-1; i++) {PetscObjectReference((PetscObject)da);}
270: for (i=0; i<dof; i++) (*dmlist)[i] = da;
271: }
272: return(0);
273: }
275: PetscErrorCode DMClone_DA(DM dm, DM *newdm)
276: {
277: DM_DA *da = (DM_DA *) dm->data;
281: DMSetType(*newdm, DMDA);
282: DMSetDimension(*newdm, dm->dim);
283: DMDASetSizes(*newdm, da->M, da->N, da->P);
284: DMDASetNumProcs(*newdm, da->m, da->n, da->p);
285: DMDASetBoundaryType(*newdm, da->bx, da->by, da->bz);
286: DMDASetDof(*newdm, da->w);
287: DMDASetStencilType(*newdm, da->stencil_type);
288: DMDASetStencilWidth(*newdm, da->s);
289: DMDASetOwnershipRanges(*newdm, da->lx, da->ly, da->lz);
290: DMSetUp(*newdm);
291: return(0);
292: }
294: static PetscErrorCode DMHasCreateInjection_DA(DM dm, PetscBool *flg)
295: {
296: DM_DA *da = (DM_DA *)dm->data;
301: *flg = da->interptype == DMDA_Q1 ? PETSC_TRUE : PETSC_FALSE;
302: return(0);
303: }
305: static PetscErrorCode DMGetDimPoints_DA(DM dm, PetscInt dim, PetscInt *pStart, PetscInt *pEnd)
306: {
310: DMDAGetDepthStratum(dm, dim, pStart, pEnd);
311: return(0);
312: }
314: static PetscErrorCode DMGetNeighbors_DA(DM dm, PetscInt *nranks, const PetscMPIInt *ranks[])
315: {
317: PetscInt dim;
318: DMDAStencilType st;
321: DMDAGetNeighbors(dm,ranks);
322: DMDAGetInfo(dm,&dim,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,&st);
324: switch (dim) {
325: case 1:
326: *nranks = 3;
327: /* if (st == DMDA_STENCIL_STAR) { *nranks = 3; } */
328: break;
329: case 2:
330: *nranks = 9;
331: /* if (st == DMDA_STENCIL_STAR) { *nranks = 5; } */
332: break;
333: case 3:
334: *nranks = 27;
335: /* if (st == DMDA_STENCIL_STAR) { *nranks = 7; } */
336: break;
337: default:
338: break;
339: }
340: return(0);
341: }
343: /*MC
344: DMDA = "da" - A DM object that is used to manage data for a structured grid in 1, 2, or 3 dimensions.
345: In the global representation of the vector each process stores a non-overlapping rectangular (or slab in 3d) portion of the grid points.
346: In the local representation these rectangular regions (slabs) are extended in all directions by a stencil width.
348: The vectors can be thought of as either cell centered or vertex centered on the mesh. But some variables cannot be cell centered and others
349: vertex centered; see the documentation for DMSTAG, a similar DM implementation which supports these staggered grids.
351: Level: intermediate
353: .seealso: DMType, DMCOMPOSITE, DMSTAG, DMDACreate(), DMCreate(), DMSetType()
354: M*/
356: extern PetscErrorCode DMLocatePoints_DA_Regular(DM,Vec,DMPointLocationType,PetscSF);
357: PETSC_INTERN PetscErrorCode DMSetUpGLVisViewer_DMDA(PetscObject,PetscViewer);
359: PETSC_EXTERN PetscErrorCode DMCreate_DA(DM da)
360: {
362: DM_DA *dd;
366: PetscNewLog(da,&dd);
367: da->data = dd;
369: da->dim = -1;
370: dd->interptype = DMDA_Q1;
371: dd->refine_x = 2;
372: dd->refine_y = 2;
373: dd->refine_z = 2;
374: dd->coarsen_x = 2;
375: dd->coarsen_y = 2;
376: dd->coarsen_z = 2;
377: dd->fieldname = NULL;
378: dd->nlocal = -1;
379: dd->Nlocal = -1;
380: dd->M = -1;
381: dd->N = -1;
382: dd->P = -1;
383: dd->m = -1;
384: dd->n = -1;
385: dd->p = -1;
386: dd->w = -1;
387: dd->s = -1;
389: dd->xs = -1; dd->xe = -1; dd->ys = -1; dd->ye = -1; dd->zs = -1; dd->ze = -1;
390: dd->Xs = -1; dd->Xe = -1; dd->Ys = -1; dd->Ye = -1; dd->Zs = -1; dd->Ze = -1;
392: dd->Nsub = 1;
393: dd->xol = 0;
394: dd->yol = 0;
395: dd->zol = 0;
396: dd->xo = 0;
397: dd->yo = 0;
398: dd->zo = 0;
399: dd->Mo = -1;
400: dd->No = -1;
401: dd->Po = -1;
403: dd->gtol = NULL;
404: dd->ltol = NULL;
405: dd->ao = NULL;
406: PetscStrallocpy(AOBASIC,(char**)&dd->aotype);
407: dd->base = -1;
408: dd->bx = DM_BOUNDARY_NONE;
409: dd->by = DM_BOUNDARY_NONE;
410: dd->bz = DM_BOUNDARY_NONE;
411: dd->stencil_type = DMDA_STENCIL_BOX;
412: dd->interptype = DMDA_Q1;
413: dd->lx = NULL;
414: dd->ly = NULL;
415: dd->lz = NULL;
417: dd->elementtype = DMDA_ELEMENT_Q1;
419: da->ops->globaltolocalbegin = DMGlobalToLocalBegin_DA;
420: da->ops->globaltolocalend = DMGlobalToLocalEnd_DA;
421: da->ops->localtoglobalbegin = DMLocalToGlobalBegin_DA;
422: da->ops->localtoglobalend = DMLocalToGlobalEnd_DA;
423: da->ops->localtolocalbegin = DMLocalToLocalBegin_DA;
424: da->ops->localtolocalend = DMLocalToLocalEnd_DA;
425: da->ops->createglobalvector = DMCreateGlobalVector_DA;
426: da->ops->createlocalvector = DMCreateLocalVector_DA;
427: da->ops->createinterpolation = DMCreateInterpolation_DA;
428: da->ops->getcoloring = DMCreateColoring_DA;
429: da->ops->creatematrix = DMCreateMatrix_DA;
430: da->ops->refine = DMRefine_DA;
431: da->ops->coarsen = DMCoarsen_DA;
432: da->ops->refinehierarchy = DMRefineHierarchy_DA;
433: da->ops->coarsenhierarchy = DMCoarsenHierarchy_DA;
434: da->ops->createinjection = DMCreateInjection_DA;
435: da->ops->hascreateinjection = DMHasCreateInjection_DA;
436: da->ops->destroy = DMDestroy_DA;
437: da->ops->view = NULL;
438: da->ops->setfromoptions = DMSetFromOptions_DA;
439: da->ops->setup = DMSetUp_DA;
440: da->ops->clone = DMClone_DA;
441: da->ops->load = DMLoad_DA;
442: da->ops->createcoordinatedm = DMCreateCoordinateDM_DA;
443: da->ops->createsubdm = DMCreateSubDM_DA;
444: da->ops->createfielddecomposition = DMCreateFieldDecomposition_DA;
445: da->ops->createdomaindecomposition = DMCreateDomainDecomposition_DA;
446: da->ops->createddscatters = DMCreateDomainDecompositionScatters_DA;
447: da->ops->getdimpoints = DMGetDimPoints_DA;
448: da->ops->getneighbors = DMGetNeighbors_DA;
449: da->ops->locatepoints = DMLocatePoints_DA_Regular;
450: da->ops->getcompatibility = DMGetCompatibility_DA;
451: PetscObjectComposeFunction((PetscObject)da,"DMSetUpGLVisViewer_C",DMSetUpGLVisViewer_DMDA);
452: return(0);
453: }
455: /*@
456: DMDACreate - Creates a DMDA object.
458: Collective
460: Input Parameter:
461: . comm - The communicator for the DMDA object
463: Output Parameter:
464: . da - The DMDA object
466: Level: advanced
468: Developers Note:
469: Since there exists DMDACreate1/2/3d() should this routine even exist?
471: .seealso: DMDASetSizes(), DMClone(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d()
472: @*/
473: PetscErrorCode DMDACreate(MPI_Comm comm, DM *da)
474: {
479: DMCreate(comm,da);
480: DMSetType(*da,DMDA);
481: return(0);
482: }