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