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