xref: /petsc/src/dm/impls/da/da.c (revision 4e278199b78715991f5c71ebbd945c1489263e6c)
1 #include <petsc/private/dmdaimpl.h>    /*I   "petscdmda.h"   I*/
2 
3 /*@
4   DMDASetSizes - Sets the number of grid points in the three dimensional directions
5 
6   Logically Collective on da
7 
8   Input Parameters:
9 + da - the DMDA
10 . M - the global X size
11 . N - the global Y size
12 - P - the global Z size
13 
14   Level: intermediate
15 
16   Developer Notes:
17   Since the dimension may not yet have been set the code cannot error check for non-positive Y and Z number of grid points
18 
19 .seealso: PetscSplitOwnership()
20 @*/
21 PetscErrorCode  DMDASetSizes(DM da, PetscInt M, PetscInt N, PetscInt P)
22 {
23   DM_DA *dd = (DM_DA*)da->data;
24 
25   PetscFunctionBegin;
26   PetscValidHeaderSpecificType(da, DM_CLASSID, 1,DMDA);
27   PetscValidLogicalCollectiveInt(da,M,2);
28   PetscValidLogicalCollectiveInt(da,N,3);
29   PetscValidLogicalCollectiveInt(da,P,4);
30   if (da->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
31   if (M < 1) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_SIZ,"Number of grid points in X direction must be positive");
32   if (N < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_SIZ,"Number of grid points in Y direction must be positive");
33   if (P < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_SIZ,"Number of grid points in Z direction must be positive");
34 
35   dd->M = M;
36   dd->N = N;
37   dd->P = P;
38   PetscFunctionReturn(0);
39 }
40 
41 /*@
42   DMDASetNumProcs - Sets the number of processes in each dimension
43 
44   Logically Collective on da
45 
46   Input Parameters:
47 + da - the DMDA
48 . m - the number of X procs (or PETSC_DECIDE)
49 . n - the number of Y procs (or PETSC_DECIDE)
50 - p - the number of Z procs (or PETSC_DECIDE)
51 
52   Level: intermediate
53 
54 .seealso: DMDASetSizes(), PetscSplitOwnership()
55 @*/
56 PetscErrorCode  DMDASetNumProcs(DM da, PetscInt m, PetscInt n, PetscInt p)
57 {
58   DM_DA          *dd = (DM_DA*)da->data;
59   PetscErrorCode ierr;
60 
61   PetscFunctionBegin;
62   PetscValidHeaderSpecificType(da, DM_CLASSID, 1,DMDA);
63   PetscValidLogicalCollectiveInt(da,m,2);
64   PetscValidLogicalCollectiveInt(da,n,3);
65   PetscValidLogicalCollectiveInt(da,p,4);
66   if (da->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
67   dd->m = m;
68   dd->n = n;
69   dd->p = p;
70   if (da->dim == 2) {
71     PetscMPIInt size;
72     ierr = MPI_Comm_size(PetscObjectComm((PetscObject)da),&size);CHKERRMPI(ierr);
73     if ((dd->m > 0) && (dd->n < 0)) {
74       dd->n = size/dd->m;
75       if (dd->n*dd->m != size) SETERRQ2(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"%D processes in X direction not divisible into comm size %d",m,size);
76     }
77     if ((dd->n > 0) && (dd->m < 0)) {
78       dd->m = size/dd->n;
79       if (dd->n*dd->m != size) SETERRQ2(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"%D processes in Y direction not divisible into comm size %d",n,size);
80     }
81   }
82   PetscFunctionReturn(0);
83 }
84 
85 /*@
86   DMDASetBoundaryType - Sets the type of ghost nodes on domain boundaries.
87 
88   Not collective
89 
90   Input Parameter:
91 + da    - The DMDA
92 - bx,by,bz - One of DM_BOUNDARY_NONE, DM_BOUNDARY_GHOSTED, DM_BOUNDARY_PERIODIC
93 
94   Level: intermediate
95 
96 .seealso: DMDACreate(), DMDestroy(), DMDA, DMBoundaryType
97 @*/
98 PetscErrorCode  DMDASetBoundaryType(DM da,DMBoundaryType bx,DMBoundaryType by,DMBoundaryType bz)
99 {
100   DM_DA *dd = (DM_DA*)da->data;
101 
102   PetscFunctionBegin;
103   PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA);
104   PetscValidLogicalCollectiveEnum(da,bx,2);
105   PetscValidLogicalCollectiveEnum(da,by,3);
106   PetscValidLogicalCollectiveEnum(da,bz,4);
107   if (da->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
108   dd->bx = bx;
109   dd->by = by;
110   dd->bz = bz;
111   PetscFunctionReturn(0);
112 }
113 
114 /*@
115   DMDASetDof - Sets the number of degrees of freedom per vertex
116 
117   Not collective
118 
119   Input Parameters:
120 + da  - The DMDA
121 - dof - Number of degrees of freedom
122 
123   Level: intermediate
124 
125 .seealso: DMDAGetDof(), DMDACreate(), DMDestroy(), DMDA
126 @*/
127 PetscErrorCode  DMDASetDof(DM da, PetscInt dof)
128 {
129   DM_DA *dd = (DM_DA*)da->data;
130 
131   PetscFunctionBegin;
132   PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA);
133   PetscValidLogicalCollectiveInt(da,dof,2);
134   if (da->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
135   dd->w  = dof;
136   da->bs = dof;
137   PetscFunctionReturn(0);
138 }
139 
140 /*@
141   DMDAGetDof - Gets the number of degrees of freedom per vertex
142 
143   Not collective
144 
145   Input Parameter:
146 . da  - The DMDA
147 
148   Output Parameter:
149 . dof - Number of degrees of freedom
150 
151   Level: intermediate
152 
153 .seealso: DMDASetDof(), DMDACreate(), DMDestroy(), DMDA
154 @*/
155 PetscErrorCode DMDAGetDof(DM da, PetscInt *dof)
156 {
157   DM_DA *dd = (DM_DA *) da->data;
158 
159   PetscFunctionBegin;
160   PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA);
161   PetscValidPointer(dof,2);
162   *dof = dd->w;
163   PetscFunctionReturn(0);
164 }
165 
166 /*@
167   DMDAGetOverlap - Gets the size of the per-processor overlap.
168 
169   Not collective
170 
171   Input Parameters:
172 . da  - The DMDA
173 
174   Output Parameters:
175 + x   - Overlap in the x direction
176 . y   - Overlap in the y direction
177 - z   - Overlap in the z direction
178 
179   Level: intermediate
180 
181 .seealso: DMDACreateDomainDecomposition(), DMDASetOverlap(), DMDA
182 @*/
183 PetscErrorCode  DMDAGetOverlap(DM da,PetscInt *x,PetscInt *y,PetscInt *z)
184 {
185   DM_DA *dd = (DM_DA*)da->data;
186 
187   PetscFunctionBegin;
188   PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA);
189   if (x) *x = dd->xol;
190   if (y) *y = dd->yol;
191   if (z) *z = dd->zol;
192   PetscFunctionReturn(0);
193 }
194 
195 /*@
196   DMDASetOverlap - Sets the size of the per-processor overlap.
197 
198   Not collective
199 
200   Input Parameters:
201 + da  - The DMDA
202 . x   - Overlap in the x direction
203 . y   - Overlap in the y direction
204 - z   - Overlap in the z direction
205 
206   Level: intermediate
207 
208 .seealso: DMDACreateDomainDecomposition(), DMDAGetOverlap(), DMDA
209 @*/
210 PetscErrorCode  DMDASetOverlap(DM da,PetscInt x,PetscInt y,PetscInt z)
211 {
212   DM_DA *dd = (DM_DA*)da->data;
213 
214   PetscFunctionBegin;
215   PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA);
216   PetscValidLogicalCollectiveInt(da,x,2);
217   PetscValidLogicalCollectiveInt(da,y,3);
218   PetscValidLogicalCollectiveInt(da,z,4);
219   dd->xol = x;
220   dd->yol = y;
221   dd->zol = z;
222   PetscFunctionReturn(0);
223 }
224 
225 /*@
226   DMDAGetNumLocalSubDomains - Gets the number of local subdomains created upon decomposition.
227 
228   Not collective
229 
230   Input Parameters:
231 . da  - The DMDA
232 
233   Output Parameters:
234 . Nsub   - Number of local subdomains created upon decomposition
235 
236   Level: intermediate
237 
238 .seealso: DMDACreateDomainDecomposition(), DMDASetNumLocalSubDomains(), DMDA
239 @*/
240 PetscErrorCode  DMDAGetNumLocalSubDomains(DM da,PetscInt *Nsub)
241 {
242   DM_DA *dd = (DM_DA*)da->data;
243 
244   PetscFunctionBegin;
245   PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA);
246   if (Nsub) *Nsub = dd->Nsub;
247   PetscFunctionReturn(0);
248 }
249 
250 /*@
251   DMDASetNumLocalSubDomains - Sets the number of local subdomains created upon decomposition.
252 
253   Not collective
254 
255   Input Parameters:
256 + da  - The DMDA
257 - Nsub - The number of local subdomains requested
258 
259   Level: intermediate
260 
261 .seealso: DMDACreateDomainDecomposition(), DMDAGetNumLocalSubDomains(), DMDA
262 @*/
263 PetscErrorCode  DMDASetNumLocalSubDomains(DM da,PetscInt Nsub)
264 {
265   DM_DA *dd = (DM_DA*)da->data;
266 
267   PetscFunctionBegin;
268   PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA);
269   PetscValidLogicalCollectiveInt(da,Nsub,2);
270   dd->Nsub = Nsub;
271   PetscFunctionReturn(0);
272 }
273 
274 /*@
275   DMDASetOffset - Sets the index offset of the DA.
276 
277   Collective on DA
278 
279   Input Parameter:
280 + da  - The DMDA
281 . xo  - The offset in the x direction
282 . yo  - The offset in the y direction
283 - zo  - The offset in the z direction
284 
285   Level: intermediate
286 
287   Notes:
288     This is used primarily to overlap a computation on a local DA with that on a global DA without
289   changing boundary conditions or subdomain features that depend upon the global offsets.
290 
291 .seealso: DMDAGetOffset(), DMDAVecGetArray()
292 @*/
293 PetscErrorCode  DMDASetOffset(DM da, PetscInt xo, PetscInt yo, PetscInt zo, PetscInt Mo, PetscInt No, PetscInt Po)
294 {
295   PetscErrorCode ierr;
296   DM_DA          *dd = (DM_DA*)da->data;
297 
298   PetscFunctionBegin;
299   PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA);
300   PetscValidLogicalCollectiveInt(da,xo,2);
301   PetscValidLogicalCollectiveInt(da,yo,3);
302   PetscValidLogicalCollectiveInt(da,zo,4);
303   PetscValidLogicalCollectiveInt(da,Mo,5);
304   PetscValidLogicalCollectiveInt(da,No,6);
305   PetscValidLogicalCollectiveInt(da,Po,7);
306   dd->xo = xo;
307   dd->yo = yo;
308   dd->zo = zo;
309   dd->Mo = Mo;
310   dd->No = No;
311   dd->Po = Po;
312 
313   if (da->coordinateDM) {
314     ierr = DMDASetOffset(da->coordinateDM,xo,yo,zo,Mo,No,Po);CHKERRQ(ierr);
315   }
316   PetscFunctionReturn(0);
317 }
318 
319 /*@
320   DMDAGetOffset - Gets the index offset of the DA.
321 
322   Not collective
323 
324   Input Parameter:
325 . da  - The DMDA
326 
327   Output Parameters:
328 + xo  - The offset in the x direction
329 . yo  - The offset in the y direction
330 . zo  - The offset in the z direction
331 . Mo  - The global size in the x direction
332 . No  - The global size in the y direction
333 - Po  - The global size in the z direction
334 
335   Level: intermediate
336 
337 .seealso: DMDASetOffset(), DMDAVecGetArray()
338 @*/
339 PetscErrorCode  DMDAGetOffset(DM da,PetscInt *xo,PetscInt *yo,PetscInt *zo,PetscInt *Mo,PetscInt *No,PetscInt *Po)
340 {
341   DM_DA *dd = (DM_DA*)da->data;
342 
343   PetscFunctionBegin;
344   PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA);
345   if (xo) *xo = dd->xo;
346   if (yo) *yo = dd->yo;
347   if (zo) *zo = dd->zo;
348   if (Mo) *Mo = dd->Mo;
349   if (No) *No = dd->No;
350   if (Po) *Po = dd->Po;
351   PetscFunctionReturn(0);
352 }
353 
354 /*@
355   DMDAGetNonOverlappingRegion - Gets the indices of the nonoverlapping region of a subdomain DM.
356 
357   Not collective
358 
359   Input Parameter:
360 . da  - The DMDA
361 
362   Output Parameters:
363 + xs  - The start of the region in x
364 . ys  - The start of the region in y
365 . zs  - The start of the region in z
366 . xs  - The size of the region in x
367 . ys  - The size of the region in y
368 - zs  - The size of the region in z
369 
370   Level: intermediate
371 
372 .seealso: DMDAGetOffset(), DMDAVecGetArray()
373 @*/
374 PetscErrorCode  DMDAGetNonOverlappingRegion(DM da, PetscInt *xs, PetscInt *ys, PetscInt *zs, PetscInt *xm, PetscInt *ym, PetscInt *zm)
375 {
376   DM_DA          *dd = (DM_DA*)da->data;
377 
378   PetscFunctionBegin;
379   PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA);
380   if (xs) *xs = dd->nonxs;
381   if (ys) *ys = dd->nonys;
382   if (zs) *zs = dd->nonzs;
383   if (xm) *xm = dd->nonxm;
384   if (ym) *ym = dd->nonym;
385   if (zm) *zm = dd->nonzm;
386   PetscFunctionReturn(0);
387 }
388 
389 /*@
390   DMDASetNonOverlappingRegion - Sets the indices of the nonoverlapping region of a subdomain DM.
391 
392   Collective on DA
393 
394   Input Parameter:
395 + da  - The DMDA
396 . xs  - The start of the region in x
397 . ys  - The start of the region in y
398 . zs  - The start of the region in z
399 . xs  - The size of the region in x
400 . ys  - The size of the region in y
401 - zs  - The size of the region in z
402 
403   Level: intermediate
404 
405 .seealso: DMDAGetOffset(), DMDAVecGetArray()
406 @*/
407 PetscErrorCode  DMDASetNonOverlappingRegion(DM da, PetscInt xs, PetscInt ys, PetscInt zs, PetscInt xm, PetscInt ym, PetscInt zm)
408 {
409   DM_DA          *dd = (DM_DA*)da->data;
410 
411   PetscFunctionBegin;
412   PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA);
413   PetscValidLogicalCollectiveInt(da,xs,2);
414   PetscValidLogicalCollectiveInt(da,ys,3);
415   PetscValidLogicalCollectiveInt(da,zs,4);
416   PetscValidLogicalCollectiveInt(da,xm,5);
417   PetscValidLogicalCollectiveInt(da,ym,6);
418   PetscValidLogicalCollectiveInt(da,zm,7);
419   dd->nonxs = xs;
420   dd->nonys = ys;
421   dd->nonzs = zs;
422   dd->nonxm = xm;
423   dd->nonym = ym;
424   dd->nonzm = zm;
425 
426   PetscFunctionReturn(0);
427 }
428 
429 /*@
430   DMDASetStencilType - Sets the type of the communication stencil
431 
432   Logically Collective on da
433 
434   Input Parameter:
435 + da    - The DMDA
436 - stype - The stencil type, use either DMDA_STENCIL_BOX or DMDA_STENCIL_STAR.
437 
438   Level: intermediate
439 
440 .seealso: DMDACreate(), DMDestroy(), DMDA
441 @*/
442 PetscErrorCode  DMDASetStencilType(DM da, DMDAStencilType stype)
443 {
444   DM_DA *dd = (DM_DA*)da->data;
445 
446   PetscFunctionBegin;
447   PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA);
448   PetscValidLogicalCollectiveEnum(da,stype,2);
449   if (da->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
450   dd->stencil_type = stype;
451   PetscFunctionReturn(0);
452 }
453 
454 /*@
455   DMDAGetStencilType - Gets the type of the communication stencil
456 
457   Not collective
458 
459   Input Parameter:
460 . da    - The DMDA
461 
462   Output Parameter:
463 . stype - The stencil type, use either DMDA_STENCIL_BOX or DMDA_STENCIL_STAR.
464 
465   Level: intermediate
466 
467 .seealso: DMDACreate(), DMDestroy(), DMDA
468 @*/
469 PetscErrorCode DMDAGetStencilType(DM da, DMDAStencilType *stype)
470 {
471   DM_DA *dd = (DM_DA*)da->data;
472 
473   PetscFunctionBegin;
474   PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA);
475   PetscValidPointer(stype,2);
476   *stype = dd->stencil_type;
477   PetscFunctionReturn(0);
478 }
479 
480 /*@
481   DMDASetStencilWidth - Sets the width of the communication stencil
482 
483   Logically Collective on da
484 
485   Input Parameter:
486 + da    - The DMDA
487 - width - The stencil width
488 
489   Level: intermediate
490 
491 .seealso: DMDACreate(), DMDestroy(), DMDA
492 @*/
493 PetscErrorCode  DMDASetStencilWidth(DM da, PetscInt width)
494 {
495   DM_DA *dd = (DM_DA*)da->data;
496 
497   PetscFunctionBegin;
498   PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA);
499   PetscValidLogicalCollectiveInt(da,width,2);
500   if (da->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
501   dd->s = width;
502   PetscFunctionReturn(0);
503 }
504 
505 /*@
506   DMDAGetStencilWidth - Gets the width of the communication stencil
507 
508   Not collective
509 
510   Input Parameter:
511 . da    - The DMDA
512 
513   Output Parameter:
514 . width - The stencil width
515 
516   Level: intermediate
517 
518 .seealso: DMDACreate(), DMDestroy(), DMDA
519 @*/
520 PetscErrorCode DMDAGetStencilWidth(DM da, PetscInt *width)
521 {
522   DM_DA *dd = (DM_DA *) da->data;
523 
524   PetscFunctionBegin;
525   PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA);
526   PetscValidPointer(width,2);
527   *width = dd->s;
528   PetscFunctionReturn(0);
529 }
530 
531 static PetscErrorCode DMDACheckOwnershipRanges_Private(DM da,PetscInt M,PetscInt m,const PetscInt lx[])
532 {
533   PetscInt i,sum;
534 
535   PetscFunctionBegin;
536   if (M < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"Global dimension not set");
537   for (i=sum=0; i<m; i++) sum += lx[i];
538   if (sum != M) SETERRQ2(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_INCOMP,"Ownership ranges sum to %D but global dimension is %D",sum,M);
539   PetscFunctionReturn(0);
540 }
541 
542 /*@
543   DMDASetOwnershipRanges - Sets the number of nodes in each direction on each process
544 
545   Logically Collective on da
546 
547   Input Parameter:
548 + da - The DMDA
549 . lx - array containing number of nodes in the X direction on each process, or NULL. If non-null, must be of length da->m
550 . ly - array containing number of nodes in the Y direction on each process, or NULL. If non-null, must be of length da->n
551 - lz - array containing number of nodes in the Z direction on each process, or NULL. If non-null, must be of length da->p.
552 
553   Level: intermediate
554 
555   Note: these numbers are NOT multiplied by the number of dof per node.
556 
557 .seealso: DMDACreate(), DMDestroy(), DMDA
558 @*/
559 PetscErrorCode  DMDASetOwnershipRanges(DM da, const PetscInt lx[], const PetscInt ly[], const PetscInt lz[])
560 {
561   PetscErrorCode ierr;
562   DM_DA          *dd = (DM_DA*)da->data;
563 
564   PetscFunctionBegin;
565   PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA);
566   if (da->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"This function must be called before DMSetUp()");
567   if (lx) {
568     if (dd->m < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs");
569     ierr = DMDACheckOwnershipRanges_Private(da,dd->M,dd->m,lx);CHKERRQ(ierr);
570     if (!dd->lx) {
571       ierr = PetscMalloc1(dd->m, &dd->lx);CHKERRQ(ierr);
572     }
573     ierr = PetscArraycpy(dd->lx, lx, dd->m);CHKERRQ(ierr);
574   }
575   if (ly) {
576     if (dd->n < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs");
577     ierr = DMDACheckOwnershipRanges_Private(da,dd->N,dd->n,ly);CHKERRQ(ierr);
578     if (!dd->ly) {
579       ierr = PetscMalloc1(dd->n, &dd->ly);CHKERRQ(ierr);
580     }
581     ierr = PetscArraycpy(dd->ly, ly, dd->n);CHKERRQ(ierr);
582   }
583   if (lz) {
584     if (dd->p < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs");
585     ierr = DMDACheckOwnershipRanges_Private(da,dd->P,dd->p,lz);CHKERRQ(ierr);
586     if (!dd->lz) {
587       ierr = PetscMalloc1(dd->p, &dd->lz);CHKERRQ(ierr);
588     }
589     ierr = PetscArraycpy(dd->lz, lz, dd->p);CHKERRQ(ierr);
590   }
591   PetscFunctionReturn(0);
592 }
593 
594 /*@
595        DMDASetInterpolationType - Sets the type of interpolation that will be
596           returned by DMCreateInterpolation()
597 
598    Logically Collective on da
599 
600    Input Parameter:
601 +  da - initial distributed array
602 -  ctype - DMDA_Q1 and DMDA_Q0 are currently the only supported forms
603 
604    Level: intermediate
605 
606    Notes:
607     you should call this on the coarser of the two DMDAs you pass to DMCreateInterpolation()
608 
609 .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDA, DMDAInterpolationType
610 @*/
611 PetscErrorCode  DMDASetInterpolationType(DM da,DMDAInterpolationType ctype)
612 {
613   DM_DA *dd = (DM_DA*)da->data;
614 
615   PetscFunctionBegin;
616   PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA);
617   PetscValidLogicalCollectiveEnum(da,ctype,2);
618   dd->interptype = ctype;
619   PetscFunctionReturn(0);
620 }
621 
622 /*@
623        DMDAGetInterpolationType - Gets the type of interpolation that will be
624           used by DMCreateInterpolation()
625 
626    Not Collective
627 
628    Input Parameter:
629 .  da - distributed array
630 
631    Output Parameter:
632 .  ctype - interpolation type (DMDA_Q1 and DMDA_Q0 are currently the only supported forms)
633 
634    Level: intermediate
635 
636 .seealso: DMDA, DMDAInterpolationType, DMDASetInterpolationType(), DMCreateInterpolation()
637 @*/
638 PetscErrorCode  DMDAGetInterpolationType(DM da,DMDAInterpolationType *ctype)
639 {
640   DM_DA *dd = (DM_DA*)da->data;
641 
642   PetscFunctionBegin;
643   PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA);
644   PetscValidPointer(ctype,2);
645   *ctype = dd->interptype;
646   PetscFunctionReturn(0);
647 }
648 
649 /*@C
650       DMDAGetNeighbors - Gets an array containing the MPI rank of all the current
651         processes neighbors.
652 
653     Not Collective
654 
655    Input Parameter:
656 .     da - the DMDA object
657 
658    Output Parameters:
659 .     ranks - the neighbors ranks, stored with the x index increasing most rapidly.
660               this process itself is in the list
661 
662    Notes:
663     In 2d the array is of length 9, in 3d of length 27
664           Not supported in 1d
665           Do not free the array, it is freed when the DMDA is destroyed.
666 
667    Fortran Notes:
668     In fortran you must pass in an array of the appropriate length.
669 
670    Level: intermediate
671 
672 @*/
673 PetscErrorCode  DMDAGetNeighbors(DM da,const PetscMPIInt *ranks[])
674 {
675   DM_DA *dd = (DM_DA*)da->data;
676 
677   PetscFunctionBegin;
678   PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA);
679   *ranks = dd->neighbors;
680   PetscFunctionReturn(0);
681 }
682 
683 /*@C
684       DMDAGetOwnershipRanges - Gets the ranges of indices in the x, y and z direction that are owned by each process
685 
686     Not Collective
687 
688    Input Parameter:
689 .     da - the DMDA object
690 
691    Output Parameter:
692 +     lx - ownership along x direction (optional)
693 .     ly - ownership along y direction (optional)
694 -     lz - ownership along z direction (optional)
695 
696    Level: intermediate
697 
698     Note: these correspond to the optional final arguments passed to DMDACreate(), DMDACreate2d(), DMDACreate3d()
699 
700     In Fortran one must pass in arrays lx, ly, and lz that are long enough to hold the values; the sixth, seventh and
701     eighth arguments from DMDAGetInfo()
702 
703      In C you should not free these arrays, nor change the values in them. They will only have valid values while the
704     DMDA they came from still exists (has not been destroyed).
705 
706     These numbers are NOT multiplied by the number of dof per node.
707 
708 .seealso: DMDAGetCorners(), DMDAGetGhostCorners(), DMDACreate(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), VecGetOwnershipRanges()
709 @*/
710 PetscErrorCode  DMDAGetOwnershipRanges(DM da,const PetscInt *lx[],const PetscInt *ly[],const PetscInt *lz[])
711 {
712   DM_DA *dd = (DM_DA*)da->data;
713 
714   PetscFunctionBegin;
715   PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA);
716   if (lx) *lx = dd->lx;
717   if (ly) *ly = dd->ly;
718   if (lz) *lz = dd->lz;
719   PetscFunctionReturn(0);
720 }
721 
722 /*@
723      DMDASetRefinementFactor - Set the ratios that the DMDA grid is refined
724 
725     Logically Collective on da
726 
727   Input Parameters:
728 +    da - the DMDA object
729 .    refine_x - ratio of fine grid to coarse in x direction (2 by default)
730 .    refine_y - ratio of fine grid to coarse in y direction (2 by default)
731 -    refine_z - ratio of fine grid to coarse in z direction (2 by default)
732 
733   Options Database:
734 +  -da_refine_x refine_x - refinement ratio in x direction
735 .  -da_refine_y rafine_y - refinement ratio in y direction
736 .  -da_refine_z refine_z - refinement ratio in z direction
737 -  -da_refine <n> - refine the DMDA object n times when it is created.
738 
739   Level: intermediate
740 
741     Notes:
742     Pass PETSC_IGNORE to leave a value unchanged
743 
744 .seealso: DMRefine(), DMDAGetRefinementFactor()
745 @*/
746 PetscErrorCode  DMDASetRefinementFactor(DM da, PetscInt refine_x, PetscInt refine_y,PetscInt refine_z)
747 {
748   DM_DA *dd = (DM_DA*)da->data;
749 
750   PetscFunctionBegin;
751   PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA);
752   PetscValidLogicalCollectiveInt(da,refine_x,2);
753   PetscValidLogicalCollectiveInt(da,refine_y,3);
754   PetscValidLogicalCollectiveInt(da,refine_z,4);
755 
756   if (refine_x > 0) dd->refine_x = refine_x;
757   if (refine_y > 0) dd->refine_y = refine_y;
758   if (refine_z > 0) dd->refine_z = refine_z;
759   PetscFunctionReturn(0);
760 }
761 
762 /*@C
763      DMDAGetRefinementFactor - Gets the ratios that the DMDA grid is refined
764 
765     Not Collective
766 
767   Input Parameter:
768 .    da - the DMDA object
769 
770   Output Parameters:
771 +    refine_x - ratio of fine grid to coarse in x direction (2 by default)
772 .    refine_y - ratio of fine grid to coarse in y direction (2 by default)
773 -    refine_z - ratio of fine grid to coarse in z direction (2 by default)
774 
775   Level: intermediate
776 
777     Notes:
778     Pass NULL for values you do not need
779 
780 .seealso: DMRefine(), DMDASetRefinementFactor()
781 @*/
782 PetscErrorCode  DMDAGetRefinementFactor(DM da, PetscInt *refine_x, PetscInt *refine_y,PetscInt *refine_z)
783 {
784   DM_DA *dd = (DM_DA*)da->data;
785 
786   PetscFunctionBegin;
787   PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA);
788   if (refine_x) *refine_x = dd->refine_x;
789   if (refine_y) *refine_y = dd->refine_y;
790   if (refine_z) *refine_z = dd->refine_z;
791   PetscFunctionReturn(0);
792 }
793 
794 /*@C
795      DMDASetGetMatrix - Sets the routine used by the DMDA to allocate a matrix.
796 
797     Logically Collective on da
798 
799   Input Parameters:
800 +    da - the DMDA object
801 -    f - the function that allocates the matrix for that specific DMDA
802 
803   Level: developer
804 
805    Notes:
806     See DMDASetBlockFills() that provides a simple way to provide the nonzero structure for
807        the diagonal and off-diagonal blocks of the matrix
808 
809    Not supported from Fortran
810 
811 .seealso: DMCreateMatrix(), DMDASetBlockFills()
812 @*/
813 PetscErrorCode  DMDASetGetMatrix(DM da,PetscErrorCode (*f)(DM, Mat*))
814 {
815   PetscFunctionBegin;
816   PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA);
817   da->ops->creatematrix = f;
818   PetscFunctionReturn(0);
819 }
820 
821 /*
822   Creates "balanced" ownership ranges after refinement, constrained by the need for the
823   fine grid boundaries to fall within one stencil width of the coarse partition.
824 
825   Uses a greedy algorithm to handle non-ideal layouts, could probably do something better.
826 */
827 static PetscErrorCode DMDARefineOwnershipRanges(DM da,PetscBool periodic,PetscInt stencil_width,PetscInt ratio,PetscInt m,const PetscInt lc[],PetscInt lf[])
828 {
829   PetscInt       i,totalc = 0,remaining,startc = 0,startf = 0;
830   PetscErrorCode ierr;
831 
832   PetscFunctionBegin;
833   if (ratio < 1) SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_USER,"Requested refinement ratio %D must be at least 1",ratio);
834   if (ratio == 1) {
835     ierr = PetscArraycpy(lf,lc,m);CHKERRQ(ierr);
836     PetscFunctionReturn(0);
837   }
838   for (i=0; i<m; i++) totalc += lc[i];
839   remaining = (!periodic) + ratio * (totalc - (!periodic));
840   for (i=0; i<m; i++) {
841     PetscInt want = remaining/(m-i) + !!(remaining%(m-i));
842     if (i == m-1) lf[i] = want;
843     else {
844       const PetscInt nextc = startc + lc[i];
845       /* Move the first fine node of the next subdomain to the right until the coarse node on its left is within one
846        * coarse stencil width of the first coarse node in the next subdomain. */
847       while ((startf+want)/ratio < nextc - stencil_width) want++;
848       /* Move the last fine node in the current subdomain to the left until the coarse node on its right is within one
849        * coarse stencil width of the last coarse node in the current subdomain. */
850       while ((startf+want-1+ratio-1)/ratio > nextc-1+stencil_width) want--;
851       /* Make sure all constraints are satisfied */
852       if (want < 0 || want > remaining || ((startf+want)/ratio < nextc - stencil_width)
853           || ((startf+want-1+ratio-1)/ratio > nextc-1+stencil_width)) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_SIZ,"Could not find a compatible refined ownership range");
854     }
855     lf[i]      = want;
856     startc    += lc[i];
857     startf    += lf[i];
858     remaining -= lf[i];
859   }
860   PetscFunctionReturn(0);
861 }
862 
863 /*
864   Creates "balanced" ownership ranges after coarsening, constrained by the need for the
865   fine grid boundaries to fall within one stencil width of the coarse partition.
866 
867   Uses a greedy algorithm to handle non-ideal layouts, could probably do something better.
868 */
869 static PetscErrorCode DMDACoarsenOwnershipRanges(DM da,PetscBool periodic,PetscInt stencil_width,PetscInt ratio,PetscInt m,const PetscInt lf[],PetscInt lc[])
870 {
871   PetscInt       i,totalf,remaining,startc,startf;
872   PetscErrorCode ierr;
873 
874   PetscFunctionBegin;
875   if (ratio < 1) SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_USER,"Requested refinement ratio %D must be at least 1",ratio);
876   if (ratio == 1) {
877     ierr = PetscArraycpy(lc,lf,m);CHKERRQ(ierr);
878     PetscFunctionReturn(0);
879   }
880   for (i=0,totalf=0; i<m; i++) totalf += lf[i];
881   remaining = (!periodic) + (totalf - (!periodic)) / ratio;
882   for (i=0,startc=0,startf=0; i<m; i++) {
883     PetscInt want = remaining/(m-i) + !!(remaining%(m-i));
884     if (i == m-1) lc[i] = want;
885     else {
886       const PetscInt nextf = startf+lf[i];
887       /* Slide first coarse node of next subdomain to the left until the coarse node to the left of the first fine
888        * node is within one stencil width. */
889       while (nextf/ratio < startc+want-stencil_width) want--;
890       /* Slide the last coarse node of the current subdomain to the right until the coarse node to the right of the last
891        * fine node is within one stencil width. */
892       while ((nextf-1+ratio-1)/ratio > startc+want-1+stencil_width) want++;
893       if (want < 0 || want > remaining
894           || (nextf/ratio < startc+want-stencil_width) || ((nextf-1+ratio-1)/ratio > startc+want-1+stencil_width)) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_SIZ,"Could not find a compatible coarsened ownership range");
895     }
896     lc[i]      = want;
897     startc    += lc[i];
898     startf    += lf[i];
899     remaining -= lc[i];
900   }
901   PetscFunctionReturn(0);
902 }
903 
904 PetscErrorCode  DMRefine_DA(DM da,MPI_Comm comm,DM *daref)
905 {
906   PetscErrorCode ierr;
907   PetscInt       M,N,P,i,dim;
908   DM             da2;
909   DM_DA          *dd = (DM_DA*)da->data,*dd2;
910 
911   PetscFunctionBegin;
912   PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA);
913   PetscValidPointer(daref,3);
914 
915   ierr = DMGetDimension(da, &dim);CHKERRQ(ierr);
916   if (dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
917     M = dd->refine_x*dd->M;
918   } else {
919     M = 1 + dd->refine_x*(dd->M - 1);
920   }
921   if (dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
922     if (dim > 1) {
923       N = dd->refine_y*dd->N;
924     } else {
925       N = 1;
926     }
927   } else {
928     N = 1 + dd->refine_y*(dd->N - 1);
929   }
930   if (dd->bz == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
931     if (dim > 2) {
932       P = dd->refine_z*dd->P;
933     } else {
934       P = 1;
935     }
936   } else {
937     P = 1 + dd->refine_z*(dd->P - 1);
938   }
939   ierr = DMDACreate(PetscObjectComm((PetscObject)da),&da2);CHKERRQ(ierr);
940   ierr = DMSetOptionsPrefix(da2,((PetscObject)da)->prefix);CHKERRQ(ierr);
941   ierr = DMSetDimension(da2,dim);CHKERRQ(ierr);
942   ierr = DMDASetSizes(da2,M,N,P);CHKERRQ(ierr);
943   ierr = DMDASetNumProcs(da2,dd->m,dd->n,dd->p);CHKERRQ(ierr);
944   ierr = DMDASetBoundaryType(da2,dd->bx,dd->by,dd->bz);CHKERRQ(ierr);
945   ierr = DMDASetDof(da2,dd->w);CHKERRQ(ierr);
946   ierr = DMDASetStencilType(da2,dd->stencil_type);CHKERRQ(ierr);
947   ierr = DMDASetStencilWidth(da2,dd->s);CHKERRQ(ierr);
948   if (dim == 3) {
949     PetscInt *lx,*ly,*lz;
950     ierr = PetscMalloc3(dd->m,&lx,dd->n,&ly,dd->p,&lz);CHKERRQ(ierr);
951     ierr = DMDARefineOwnershipRanges(da,(PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_x,dd->m,dd->lx,lx);CHKERRQ(ierr);
952     ierr = DMDARefineOwnershipRanges(da,(PetscBool)(dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_y,dd->n,dd->ly,ly);CHKERRQ(ierr);
953     ierr = DMDARefineOwnershipRanges(da,(PetscBool)(dd->bz == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_z,dd->p,dd->lz,lz);CHKERRQ(ierr);
954     ierr = DMDASetOwnershipRanges(da2,lx,ly,lz);CHKERRQ(ierr);
955     ierr = PetscFree3(lx,ly,lz);CHKERRQ(ierr);
956   } else if (dim == 2) {
957     PetscInt *lx,*ly;
958     ierr = PetscMalloc2(dd->m,&lx,dd->n,&ly);CHKERRQ(ierr);
959     ierr = DMDARefineOwnershipRanges(da,(PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_x,dd->m,dd->lx,lx);CHKERRQ(ierr);
960     ierr = DMDARefineOwnershipRanges(da,(PetscBool)(dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_y,dd->n,dd->ly,ly);CHKERRQ(ierr);
961     ierr = DMDASetOwnershipRanges(da2,lx,ly,NULL);CHKERRQ(ierr);
962     ierr = PetscFree2(lx,ly);CHKERRQ(ierr);
963   } else if (dim == 1) {
964     PetscInt *lx;
965     ierr = PetscMalloc1(dd->m,&lx);CHKERRQ(ierr);
966     ierr = DMDARefineOwnershipRanges(da,(PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->refine_x,dd->m,dd->lx,lx);CHKERRQ(ierr);
967     ierr = DMDASetOwnershipRanges(da2,lx,NULL,NULL);CHKERRQ(ierr);
968     ierr = PetscFree(lx);CHKERRQ(ierr);
969   }
970   dd2 = (DM_DA*)da2->data;
971 
972   /* allow overloaded (user replaced) operations to be inherited by refinement clones */
973   da2->ops->creatematrix = da->ops->creatematrix;
974   /* da2->ops->createinterpolation = da->ops->createinterpolation; this causes problem with SNESVI */
975   da2->ops->getcoloring = da->ops->getcoloring;
976   dd2->interptype       = dd->interptype;
977 
978   /* copy fill information if given */
979   if (dd->dfill) {
980     ierr = PetscMalloc1(dd->dfill[dd->w]+dd->w+1,&dd2->dfill);CHKERRQ(ierr);
981     ierr = PetscArraycpy(dd2->dfill,dd->dfill,dd->dfill[dd->w]+dd->w+1);CHKERRQ(ierr);
982   }
983   if (dd->ofill) {
984     ierr = PetscMalloc1(dd->ofill[dd->w]+dd->w+1,&dd2->ofill);CHKERRQ(ierr);
985     ierr = PetscArraycpy(dd2->ofill,dd->ofill,dd->ofill[dd->w]+dd->w+1);CHKERRQ(ierr);
986   }
987   /* copy the refine information */
988   dd2->coarsen_x = dd2->refine_x = dd->refine_x;
989   dd2->coarsen_y = dd2->refine_y = dd->refine_y;
990   dd2->coarsen_z = dd2->refine_z = dd->refine_z;
991 
992   if (dd->refine_z_hier) {
993     if (da->levelup - da->leveldown + 1 > -1 && da->levelup - da->leveldown + 1 < dd->refine_z_hier_n) {
994       dd2->refine_z = dd->refine_z_hier[da->levelup - da->leveldown + 1];
995     }
996     if (da->levelup - da->leveldown > -1 && da->levelup - da->leveldown < dd->refine_z_hier_n) {
997       dd2->coarsen_z = dd->refine_z_hier[da->levelup - da->leveldown];
998     }
999     dd2->refine_z_hier_n = dd->refine_z_hier_n;
1000     ierr = PetscMalloc1(dd2->refine_z_hier_n,&dd2->refine_z_hier);CHKERRQ(ierr);
1001     ierr = PetscArraycpy(dd2->refine_z_hier,dd->refine_z_hier,dd2->refine_z_hier_n);CHKERRQ(ierr);
1002   }
1003   if (dd->refine_y_hier) {
1004     if (da->levelup - da->leveldown + 1 > -1 && da->levelup - da->leveldown + 1 < dd->refine_y_hier_n) {
1005       dd2->refine_y = dd->refine_y_hier[da->levelup - da->leveldown + 1];
1006     }
1007     if (da->levelup - da->leveldown > -1 && da->levelup - da->leveldown < dd->refine_y_hier_n) {
1008       dd2->coarsen_y = dd->refine_y_hier[da->levelup - da->leveldown];
1009     }
1010     dd2->refine_y_hier_n = dd->refine_y_hier_n;
1011     ierr = PetscMalloc1(dd2->refine_y_hier_n,&dd2->refine_y_hier);CHKERRQ(ierr);
1012     ierr = PetscArraycpy(dd2->refine_y_hier,dd->refine_y_hier,dd2->refine_y_hier_n);CHKERRQ(ierr);
1013   }
1014   if (dd->refine_x_hier) {
1015     if (da->levelup - da->leveldown + 1 > -1 && da->levelup - da->leveldown + 1 < dd->refine_x_hier_n) {
1016       dd2->refine_x = dd->refine_x_hier[da->levelup - da->leveldown + 1];
1017     }
1018     if (da->levelup - da->leveldown > -1 && da->levelup - da->leveldown < dd->refine_x_hier_n) {
1019       dd2->coarsen_x = dd->refine_x_hier[da->levelup - da->leveldown];
1020     }
1021     dd2->refine_x_hier_n = dd->refine_x_hier_n;
1022     ierr = PetscMalloc1(dd2->refine_x_hier_n,&dd2->refine_x_hier);CHKERRQ(ierr);
1023     ierr = PetscArraycpy(dd2->refine_x_hier,dd->refine_x_hier,dd2->refine_x_hier_n);CHKERRQ(ierr);
1024   }
1025 
1026   /* copy vector type information */
1027   ierr = DMSetVecType(da2,da->vectype);CHKERRQ(ierr);
1028 
1029   dd2->lf = dd->lf;
1030   dd2->lj = dd->lj;
1031 
1032   da2->leveldown = da->leveldown;
1033   da2->levelup   = da->levelup + 1;
1034 
1035   ierr = DMSetUp(da2);CHKERRQ(ierr);
1036 
1037   /* interpolate coordinates if they are set on the coarse grid */
1038   if (da->coordinates) {
1039     DM  cdaf,cdac;
1040     Vec coordsc,coordsf;
1041     Mat II;
1042 
1043     ierr = DMGetCoordinateDM(da,&cdac);CHKERRQ(ierr);
1044     ierr = DMGetCoordinates(da,&coordsc);CHKERRQ(ierr);
1045     ierr = DMGetCoordinateDM(da2,&cdaf);CHKERRQ(ierr);
1046     /* force creation of the coordinate vector */
1047     ierr = DMDASetUniformCoordinates(da2,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr);
1048     ierr = DMGetCoordinates(da2,&coordsf);CHKERRQ(ierr);
1049     ierr = DMCreateInterpolation(cdac,cdaf,&II,NULL);CHKERRQ(ierr);
1050     ierr = MatInterpolate(II,coordsc,coordsf);CHKERRQ(ierr);
1051     ierr = MatDestroy(&II);CHKERRQ(ierr);
1052   }
1053 
1054   for (i=0; i<da->bs; i++) {
1055     const char *fieldname;
1056     ierr = DMDAGetFieldName(da,i,&fieldname);CHKERRQ(ierr);
1057     ierr = DMDASetFieldName(da2,i,fieldname);CHKERRQ(ierr);
1058   }
1059 
1060   *daref = da2;
1061   PetscFunctionReturn(0);
1062 }
1063 
1064 PetscErrorCode  DMCoarsen_DA(DM dmf, MPI_Comm comm,DM *dmc)
1065 {
1066   PetscErrorCode ierr;
1067   PetscInt       M,N,P,i,dim;
1068   DM             dmc2;
1069   DM_DA          *dd = (DM_DA*)dmf->data,*dd2;
1070 
1071   PetscFunctionBegin;
1072   PetscValidHeaderSpecificType(dmf,DM_CLASSID,1,DMDA);
1073   PetscValidPointer(dmc,3);
1074 
1075   ierr = DMGetDimension(dmf, &dim);CHKERRQ(ierr);
1076   if (dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
1077     M = dd->M / dd->coarsen_x;
1078   } else {
1079     M = 1 + (dd->M - 1) / dd->coarsen_x;
1080   }
1081   if (dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
1082     if (dim > 1) {
1083       N = dd->N / dd->coarsen_y;
1084     } else {
1085       N = 1;
1086     }
1087   } else {
1088     N = 1 + (dd->N - 1) / dd->coarsen_y;
1089   }
1090   if (dd->bz == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
1091     if (dim > 2) {
1092       P = dd->P / dd->coarsen_z;
1093     } else {
1094       P = 1;
1095     }
1096   } else {
1097     P = 1 + (dd->P - 1) / dd->coarsen_z;
1098   }
1099   ierr = DMDACreate(PetscObjectComm((PetscObject)dmf),&dmc2);CHKERRQ(ierr);
1100   ierr = DMSetOptionsPrefix(dmc2,((PetscObject)dmf)->prefix);CHKERRQ(ierr);
1101   ierr = DMSetDimension(dmc2,dim);CHKERRQ(ierr);
1102   ierr = DMDASetSizes(dmc2,M,N,P);CHKERRQ(ierr);
1103   ierr = DMDASetNumProcs(dmc2,dd->m,dd->n,dd->p);CHKERRQ(ierr);
1104   ierr = DMDASetBoundaryType(dmc2,dd->bx,dd->by,dd->bz);CHKERRQ(ierr);
1105   ierr = DMDASetDof(dmc2,dd->w);CHKERRQ(ierr);
1106   ierr = DMDASetStencilType(dmc2,dd->stencil_type);CHKERRQ(ierr);
1107   ierr = DMDASetStencilWidth(dmc2,dd->s);CHKERRQ(ierr);
1108   if (dim == 3) {
1109     PetscInt *lx,*ly,*lz;
1110     ierr = PetscMalloc3(dd->m,&lx,dd->n,&ly,dd->p,&lz);CHKERRQ(ierr);
1111     ierr = DMDACoarsenOwnershipRanges(dmf,(PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_x,dd->m,dd->lx,lx);CHKERRQ(ierr);
1112     ierr = DMDACoarsenOwnershipRanges(dmf,(PetscBool)(dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_y,dd->n,dd->ly,ly);CHKERRQ(ierr);
1113     ierr = DMDACoarsenOwnershipRanges(dmf,(PetscBool)(dd->bz == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_z,dd->p,dd->lz,lz);CHKERRQ(ierr);
1114     ierr = DMDASetOwnershipRanges(dmc2,lx,ly,lz);CHKERRQ(ierr);
1115     ierr = PetscFree3(lx,ly,lz);CHKERRQ(ierr);
1116   } else if (dim == 2) {
1117     PetscInt *lx,*ly;
1118     ierr = PetscMalloc2(dd->m,&lx,dd->n,&ly);CHKERRQ(ierr);
1119     ierr = DMDACoarsenOwnershipRanges(dmf,(PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_x,dd->m,dd->lx,lx);CHKERRQ(ierr);
1120     ierr = DMDACoarsenOwnershipRanges(dmf,(PetscBool)(dd->by == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_y,dd->n,dd->ly,ly);CHKERRQ(ierr);
1121     ierr = DMDASetOwnershipRanges(dmc2,lx,ly,NULL);CHKERRQ(ierr);
1122     ierr = PetscFree2(lx,ly);CHKERRQ(ierr);
1123   } else if (dim == 1) {
1124     PetscInt *lx;
1125     ierr = PetscMalloc1(dd->m,&lx);CHKERRQ(ierr);
1126     ierr = DMDACoarsenOwnershipRanges(dmf,(PetscBool)(dd->bx == DM_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_x,dd->m,dd->lx,lx);CHKERRQ(ierr);
1127     ierr = DMDASetOwnershipRanges(dmc2,lx,NULL,NULL);CHKERRQ(ierr);
1128     ierr = PetscFree(lx);CHKERRQ(ierr);
1129   }
1130   dd2 = (DM_DA*)dmc2->data;
1131 
1132   /* allow overloaded (user replaced) operations to be inherited by refinement clones; why are only some inherited and not all? */
1133   /* dmc2->ops->createinterpolation = dmf->ops->createinterpolation; copying this one causes trouble for DMSetVI */
1134   dmc2->ops->creatematrix = dmf->ops->creatematrix;
1135   dmc2->ops->getcoloring  = dmf->ops->getcoloring;
1136   dd2->interptype         = dd->interptype;
1137 
1138   /* copy fill information if given */
1139   if (dd->dfill) {
1140     ierr = PetscMalloc1(dd->dfill[dd->w]+dd->w+1,&dd2->dfill);CHKERRQ(ierr);
1141     ierr = PetscArraycpy(dd2->dfill,dd->dfill,dd->dfill[dd->w]+dd->w+1);CHKERRQ(ierr);
1142   }
1143   if (dd->ofill) {
1144     ierr = PetscMalloc1(dd->ofill[dd->w]+dd->w+1,&dd2->ofill);CHKERRQ(ierr);
1145     ierr = PetscArraycpy(dd2->ofill,dd->ofill,dd->ofill[dd->w]+dd->w+1);CHKERRQ(ierr);
1146   }
1147   /* copy the refine information */
1148   dd2->coarsen_x = dd2->refine_x = dd->coarsen_x;
1149   dd2->coarsen_y = dd2->refine_y = dd->coarsen_y;
1150   dd2->coarsen_z = dd2->refine_z = dd->coarsen_z;
1151 
1152   if (dd->refine_z_hier) {
1153     if (dmf->levelup - dmf->leveldown -1 > -1 && dmf->levelup - dmf->leveldown - 1< dd->refine_z_hier_n) {
1154       dd2->refine_z = dd->refine_z_hier[dmf->levelup - dmf->leveldown - 1];
1155     }
1156     if (dmf->levelup - dmf->leveldown - 2 > -1 && dmf->levelup - dmf->leveldown - 2 < dd->refine_z_hier_n) {
1157       dd2->coarsen_z = dd->refine_z_hier[dmf->levelup - dmf->leveldown - 2];
1158     }
1159     dd2->refine_z_hier_n = dd->refine_z_hier_n;
1160     ierr = PetscMalloc1(dd2->refine_z_hier_n,&dd2->refine_z_hier);CHKERRQ(ierr);
1161     ierr = PetscArraycpy(dd2->refine_z_hier,dd->refine_z_hier,dd2->refine_z_hier_n);CHKERRQ(ierr);
1162   }
1163   if (dd->refine_y_hier) {
1164     if (dmf->levelup - dmf->leveldown - 1 > -1 && dmf->levelup - dmf->leveldown - 1< dd->refine_y_hier_n) {
1165       dd2->refine_y = dd->refine_y_hier[dmf->levelup - dmf->leveldown - 1];
1166     }
1167     if (dmf->levelup - dmf->leveldown - 2 > -1 && dmf->levelup - dmf->leveldown - 2 < dd->refine_y_hier_n) {
1168       dd2->coarsen_y = dd->refine_y_hier[dmf->levelup - dmf->leveldown - 2];
1169     }
1170     dd2->refine_y_hier_n = dd->refine_y_hier_n;
1171     ierr = PetscMalloc1(dd2->refine_y_hier_n,&dd2->refine_y_hier);CHKERRQ(ierr);
1172     ierr = PetscArraycpy(dd2->refine_y_hier,dd->refine_y_hier,dd2->refine_y_hier_n);CHKERRQ(ierr);
1173   }
1174   if (dd->refine_x_hier) {
1175     if (dmf->levelup - dmf->leveldown - 1 > -1 && dmf->levelup - dmf->leveldown - 1 < dd->refine_x_hier_n) {
1176       dd2->refine_x = dd->refine_x_hier[dmf->levelup - dmf->leveldown - 1];
1177     }
1178     if (dmf->levelup - dmf->leveldown - 2 > -1 && dmf->levelup - dmf->leveldown - 2 < dd->refine_x_hier_n) {
1179       dd2->coarsen_x = dd->refine_x_hier[dmf->levelup - dmf->leveldown - 2];
1180     }
1181     dd2->refine_x_hier_n = dd->refine_x_hier_n;
1182     ierr = PetscMalloc1(dd2->refine_x_hier_n,&dd2->refine_x_hier);CHKERRQ(ierr);
1183     ierr = PetscArraycpy(dd2->refine_x_hier,dd->refine_x_hier,dd2->refine_x_hier_n);CHKERRQ(ierr);
1184   }
1185 
1186   /* copy vector type information */
1187   ierr = DMSetVecType(dmc2,dmf->vectype);CHKERRQ(ierr);
1188 
1189   dd2->lf = dd->lf;
1190   dd2->lj = dd->lj;
1191 
1192   dmc2->leveldown = dmf->leveldown + 1;
1193   dmc2->levelup   = dmf->levelup;
1194 
1195   ierr = DMSetUp(dmc2);CHKERRQ(ierr);
1196 
1197   /* inject coordinates if they are set on the fine grid */
1198   if (dmf->coordinates) {
1199     DM         cdaf,cdac;
1200     Vec        coordsc,coordsf;
1201     Mat        inject;
1202     VecScatter vscat;
1203 
1204     ierr = DMGetCoordinateDM(dmf,&cdaf);CHKERRQ(ierr);
1205     ierr = DMGetCoordinates(dmf,&coordsf);CHKERRQ(ierr);
1206     ierr = DMGetCoordinateDM(dmc2,&cdac);CHKERRQ(ierr);
1207     /* force creation of the coordinate vector */
1208     ierr = DMDASetUniformCoordinates(dmc2,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr);
1209     ierr = DMGetCoordinates(dmc2,&coordsc);CHKERRQ(ierr);
1210 
1211     ierr = DMCreateInjection(cdac,cdaf,&inject);CHKERRQ(ierr);
1212     ierr = MatScatterGetVecScatter(inject,&vscat);CHKERRQ(ierr);
1213     ierr = VecScatterBegin(vscat,coordsf,coordsc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1214     ierr = VecScatterEnd(vscat,coordsf,coordsc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1215     ierr = MatDestroy(&inject);CHKERRQ(ierr);
1216   }
1217 
1218   for (i=0; i<dmf->bs; i++) {
1219     const char *fieldname;
1220     ierr = DMDAGetFieldName(dmf,i,&fieldname);CHKERRQ(ierr);
1221     ierr = DMDASetFieldName(dmc2,i,fieldname);CHKERRQ(ierr);
1222   }
1223 
1224   *dmc = dmc2;
1225   PetscFunctionReturn(0);
1226 }
1227 
1228 PetscErrorCode  DMRefineHierarchy_DA(DM da,PetscInt nlevels,DM daf[])
1229 {
1230   PetscErrorCode ierr;
1231   PetscInt       i,n,*refx,*refy,*refz;
1232 
1233   PetscFunctionBegin;
1234   PetscValidHeaderSpecific(da,DM_CLASSID,1);
1235   if (nlevels < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
1236   if (nlevels == 0) PetscFunctionReturn(0);
1237   PetscValidPointer(daf,3);
1238 
1239   /* Get refinement factors, defaults taken from the coarse DMDA */
1240   ierr = PetscMalloc3(nlevels,&refx,nlevels,&refy,nlevels,&refz);CHKERRQ(ierr);
1241   for (i=0; i<nlevels; i++) {
1242     ierr = DMDAGetRefinementFactor(da,&refx[i],&refy[i],&refz[i]);CHKERRQ(ierr);
1243   }
1244   n    = nlevels;
1245   ierr = PetscOptionsGetIntArray(((PetscObject)da)->options,((PetscObject)da)->prefix,"-da_refine_hierarchy_x",refx,&n,NULL);CHKERRQ(ierr);
1246   n    = nlevels;
1247   ierr = PetscOptionsGetIntArray(((PetscObject)da)->options,((PetscObject)da)->prefix,"-da_refine_hierarchy_y",refy,&n,NULL);CHKERRQ(ierr);
1248   n    = nlevels;
1249   ierr = PetscOptionsGetIntArray(((PetscObject)da)->options,((PetscObject)da)->prefix,"-da_refine_hierarchy_z",refz,&n,NULL);CHKERRQ(ierr);
1250 
1251   ierr = DMDASetRefinementFactor(da,refx[0],refy[0],refz[0]);CHKERRQ(ierr);
1252   ierr = DMRefine(da,PetscObjectComm((PetscObject)da),&daf[0]);CHKERRQ(ierr);
1253   for (i=1; i<nlevels; i++) {
1254     ierr = DMDASetRefinementFactor(daf[i-1],refx[i],refy[i],refz[i]);CHKERRQ(ierr);
1255     ierr = DMRefine(daf[i-1],PetscObjectComm((PetscObject)da),&daf[i]);CHKERRQ(ierr);
1256   }
1257   ierr = PetscFree3(refx,refy,refz);CHKERRQ(ierr);
1258   PetscFunctionReturn(0);
1259 }
1260 
1261 PetscErrorCode  DMCoarsenHierarchy_DA(DM da,PetscInt nlevels,DM dac[])
1262 {
1263   PetscErrorCode ierr;
1264   PetscInt       i;
1265 
1266   PetscFunctionBegin;
1267   PetscValidHeaderSpecific(da,DM_CLASSID,1);
1268   if (nlevels < 0) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
1269   if (nlevels == 0) PetscFunctionReturn(0);
1270   PetscValidPointer(dac,3);
1271   ierr = DMCoarsen(da,PetscObjectComm((PetscObject)da),&dac[0]);CHKERRQ(ierr);
1272   for (i=1; i<nlevels; i++) {
1273     ierr = DMCoarsen(dac[i-1],PetscObjectComm((PetscObject)da),&dac[i]);CHKERRQ(ierr);
1274   }
1275   PetscFunctionReturn(0);
1276 }
1277 
1278 PetscErrorCode DMDASetGLLCoordinates_1d(DM dm,PetscInt n,PetscReal *nodes)
1279 {
1280   PetscErrorCode ierr;
1281   PetscInt       i,j,xs,xn,q;
1282   PetscScalar    *xx;
1283   PetscReal      h;
1284   Vec            x;
1285   DM_DA          *da = (DM_DA*)dm->data;
1286 
1287   PetscFunctionBegin;
1288   if (da->bx != DM_BOUNDARY_PERIODIC) {
1289     ierr = DMDAGetInfo(dm,NULL,&q,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
1290     q    = (q-1)/(n-1);  /* number of spectral elements */
1291     h    = 2.0/q;
1292     ierr = DMDAGetCorners(dm,&xs,NULL,NULL,&xn,NULL,NULL);CHKERRQ(ierr);
1293     xs   = xs/(n-1);
1294     xn   = xn/(n-1);
1295     ierr = DMDASetUniformCoordinates(dm,-1.,1.,0.,0.,0.,0.);CHKERRQ(ierr);
1296     ierr = DMGetCoordinates(dm,&x);CHKERRQ(ierr);
1297     ierr = DMDAVecGetArray(dm,x,&xx);CHKERRQ(ierr);
1298 
1299     /* loop over local spectral elements */
1300     for (j=xs; j<xs+xn; j++) {
1301       /*
1302        Except for the first process, each process starts on the second GLL point of the first element on that process
1303        */
1304       for (i= (j == xs && xs > 0)? 1 : 0; i<n; i++) {
1305         xx[j*(n-1) + i] = -1.0 + h*j + h*(nodes[i]+1.0)/2.;
1306       }
1307     }
1308     ierr = DMDAVecRestoreArray(dm,x,&xx);CHKERRQ(ierr);
1309   } else SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Not yet implemented for periodic");
1310   PetscFunctionReturn(0);
1311 }
1312 
1313 /*@
1314 
1315      DMDASetGLLCoordinates - Sets the global coordinates from -1 to 1 to the GLL points of as many GLL elements that fit the number of grid points
1316 
1317    Collective on da
1318 
1319    Input Parameters:
1320 +   da - the DMDA object
1321 -   n - the number of GLL nodes
1322 -   nodes - the GLL nodes
1323 
1324    Notes:
1325     the parallel decomposition of grid points must correspond to the degree of the GLL. That is, the number of grid points
1326           on each process much be divisible by the number of GLL elements needed per process. This depends on whether the DM is
1327           periodic or not.
1328 
1329    Level: advanced
1330 
1331 .seealso:   DMDACreate(), PetscDTGaussLobattoLegendreQuadrature(), DMGetCoordinates()
1332 @*/
1333 PetscErrorCode DMDASetGLLCoordinates(DM da,PetscInt n,PetscReal *nodes)
1334 {
1335   PetscErrorCode ierr;
1336 
1337   PetscFunctionBegin;
1338   if (da->dim == 1) {
1339     ierr = DMDASetGLLCoordinates_1d(da,n,nodes);CHKERRQ(ierr);
1340   } else SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Not yet implemented for 2 or 3d");
1341   PetscFunctionReturn(0);
1342 }
1343 
1344 PETSC_INTERN PetscErrorCode DMGetCompatibility_DA(DM da1,DM dm2,PetscBool *compatible,PetscBool *set)
1345 {
1346   PetscErrorCode ierr;
1347   DM_DA          *dd1 = (DM_DA*)da1->data,*dd2;
1348   DM             da2;
1349   DMType         dmtype2;
1350   PetscBool      isda,compatibleLocal;
1351   PetscInt       i;
1352 
1353   PetscFunctionBegin;
1354   if (!da1->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da1),PETSC_ERR_ARG_WRONGSTATE,"DMSetUp() must be called on first DM before DMGetCompatibility()");
1355   ierr = DMGetType(dm2,&dmtype2);CHKERRQ(ierr);
1356   ierr = PetscStrcmp(dmtype2,DMDA,&isda);CHKERRQ(ierr);
1357   if (isda) {
1358     da2 = dm2;
1359     dd2 = (DM_DA*)da2->data;
1360     if (!da2->setupcalled) SETERRQ(PetscObjectComm((PetscObject)da2),PETSC_ERR_ARG_WRONGSTATE,"DMSetUp() must be called on second DM before DMGetCompatibility()");
1361     compatibleLocal = (PetscBool)(da1->dim == da2->dim);
1362     if (compatibleLocal) compatibleLocal = (PetscBool)(compatibleLocal && (dd1->s == dd2->s)); /* Stencil width */
1363     /*                                                                           Global size              ranks               Boundary type */
1364     if (compatibleLocal)                 compatibleLocal = (PetscBool)(compatibleLocal && (dd1->M == dd2->M) && (dd1->m == dd2->m) && (dd1->bx == dd2->bx));
1365     if (compatibleLocal && da1->dim > 1) compatibleLocal = (PetscBool)(compatibleLocal && (dd1->N == dd2->N) && (dd1->n == dd2->n) && (dd1->by == dd2->by));
1366     if (compatibleLocal && da1->dim > 2) compatibleLocal = (PetscBool)(compatibleLocal && (dd1->P == dd2->P) && (dd1->p == dd2->p) && (dd1->bz == dd2->bz));
1367     if (compatibleLocal) {
1368       for (i=0; i<dd1->m; ++i) {
1369         compatibleLocal = (PetscBool)(compatibleLocal && (dd1->lx[i] == dd2->lx[i]));           /* Local size     */
1370       }
1371     }
1372     if (compatibleLocal && da1->dim > 1) {
1373       for (i=0; i<dd1->n; ++i) {
1374         compatibleLocal = (PetscBool)(compatibleLocal && (dd1->ly[i] == dd2->ly[i]));
1375       }
1376     }
1377     if (compatibleLocal && da1->dim > 2) {
1378       for (i=0; i<dd1->p; ++i) {
1379         compatibleLocal = (PetscBool)(compatibleLocal && (dd1->lz[i] == dd2->lz[i]));
1380       }
1381     }
1382     *compatible = compatibleLocal;
1383     *set = PETSC_TRUE;
1384   } else {
1385     /* Decline to determine compatibility with other DM types */
1386     *set = PETSC_FALSE;
1387   }
1388   PetscFunctionReturn(0);
1389 }
1390