xref: /petsc/src/dm/impls/da/da.c (revision efe48dd8cf8b935accbbb9f4bcb20bc83865fa4d)
1 #define PETSCDM_DLL
2 #include "private/daimpl.h"    /*I   "petscdm.h"   I*/
3 
4 
5 #undef __FUNCT__
6 #define __FUNCT__ "DMDASetDim"
7 /*@
8   DMDASetDim - Sets the dimension
9 
10   Collective on DMDA
11 
12   Input Parameters:
13 + da - the DMDA
14 - dim - the dimension (or PETSC_DECIDE)
15 
16   Level: intermediate
17 
18 .seealso: DaGetDim(), DMDASetSizes()
19 @*/
20 PetscErrorCode PETSCDM_DLLEXPORT DMDASetDim(DM da, PetscInt dim)
21 {
22   DM_DA *dd = (DM_DA*)da->data;
23 
24   PetscFunctionBegin;
25   PetscValidHeaderSpecific(da, DM_CLASSID, 1);
26   if (dd->dim > 0 && dim != dd->dim) SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot change DMDA dim from %D after it was set to %D",dd->dim,dim);
27   dd->dim = dim;
28   PetscFunctionReturn(0);
29 }
30 
31 #undef __FUNCT__
32 #define __FUNCT__ "DMDASetSizes"
33 /*@
34   DMDASetSizes - Sets the global sizes
35 
36   Logically Collective on DMDA
37 
38   Input Parameters:
39 + da - the DMDA
40 . M - the global X size (or PETSC_DECIDE)
41 . N - the global Y size (or PETSC_DECIDE)
42 - P - the global Z size (or PETSC_DECIDE)
43 
44   Level: intermediate
45 
46 .seealso: DMDAGetSize(), PetscSplitOwnership()
47 @*/
48 PetscErrorCode PETSCDM_DLLEXPORT DMDASetSizes(DM da, PetscInt M, PetscInt N, PetscInt P)
49 {
50   DM_DA *dd = (DM_DA*)da->data;
51 
52   PetscFunctionBegin;
53   PetscValidHeaderSpecific(da, DM_CLASSID, 1);
54   PetscValidLogicalCollectiveInt(da,M,2);
55   PetscValidLogicalCollectiveInt(da,N,3);
56   PetscValidLogicalCollectiveInt(da,P,4);
57 
58   dd->M = M;
59   dd->N = N;
60   dd->P = P;
61   PetscFunctionReturn(0);
62 }
63 
64 #undef __FUNCT__
65 #define __FUNCT__ "DMDASetNumProcs"
66 /*@
67   DMDASetNumProcs - Sets the number of processes in each dimension
68 
69   Logically Collective on DMDA
70 
71   Input Parameters:
72 + da - the DMDA
73 . m - the number of X procs (or PETSC_DECIDE)
74 . n - the number of Y procs (or PETSC_DECIDE)
75 - p - the number of Z procs (or PETSC_DECIDE)
76 
77   Level: intermediate
78 
79 .seealso: DMDASetSizes(), DMDAGetSize(), PetscSplitOwnership()
80 @*/
81 PetscErrorCode PETSCDM_DLLEXPORT DMDASetNumProcs(DM da, PetscInt m, PetscInt n, PetscInt p)
82 {
83   DM_DA *dd = (DM_DA*)da->data;
84 
85   PetscFunctionBegin;
86   PetscValidHeaderSpecific(da, DM_CLASSID, 1);
87   PetscValidLogicalCollectiveInt(da,m,2);
88   PetscValidLogicalCollectiveInt(da,n,3);
89   PetscValidLogicalCollectiveInt(da,p,4);
90   dd->m = m;
91   dd->n = n;
92   dd->p = p;
93   PetscFunctionReturn(0);
94 }
95 
96 #undef __FUNCT__
97 #define __FUNCT__ "DMDASetPeriodicity"
98 /*@
99   DMDASetPeriodicity - Sets the type of periodicity
100 
101   Not collective
102 
103   Input Parameter:
104 + da    - The DMDA
105 - ptype - One of DMDA_NONPERIODIC, DMDA_XPERIODIC, DMDA_YPERIODIC, DMDA_ZPERIODIC, DMDA_XYPERIODIC, DMDA_XZPERIODIC, DMDA_YZPERIODIC, or DMDA_XYZPERIODIC
106 
107   Level: intermediate
108 
109 .keywords:  distributed array, periodicity
110 .seealso: DMDACreate(), DMDestroy(), DMDA, DMDAPeriodicType
111 @*/
112 PetscErrorCode PETSCDM_DLLEXPORT DMDASetPeriodicity(DM da, DMDAPeriodicType ptype)
113 {
114   DM_DA *dd = (DM_DA*)da->data;
115 
116   PetscFunctionBegin;
117   PetscValidHeaderSpecific(da,DM_CLASSID,1);
118   dd->wrap = ptype;
119   PetscFunctionReturn(0);
120 }
121 
122 #undef __FUNCT__
123 #define __FUNCT__ "DMDASetDof"
124 /*@
125   DMDASetDof - Sets the number of degrees of freedom per vertex
126 
127   Not collective
128 
129   Input Parameter:
130 + da  - The DMDA
131 - dof - Number of degrees of freedom
132 
133   Level: intermediate
134 
135 .keywords:  distributed array, degrees of freedom
136 .seealso: DMDACreate(), DMDestroy(), DMDA
137 @*/
138 PetscErrorCode PETSCDM_DLLEXPORT DMDASetDof(DM da, int dof)
139 {
140   DM_DA *dd = (DM_DA*)da->data;
141 
142   PetscFunctionBegin;
143   PetscValidHeaderSpecific(da,DM_CLASSID,1);
144   dd->w = dof;
145   PetscFunctionReturn(0);
146 }
147 
148 #undef __FUNCT__
149 #define __FUNCT__ "DMDASetStencilType"
150 /*@
151   DMDASetStencilType - Sets the type of the communication stencil
152 
153   Logically Collective on DMDA
154 
155   Input Parameter:
156 + da    - The DMDA
157 - stype - The stencil type, use either DMDA_STENCIL_BOX or DMDA_STENCIL_STAR.
158 
159   Level: intermediate
160 
161 .keywords:  distributed array, stencil
162 .seealso: DMDACreate(), DMDestroy(), DMDA
163 @*/
164 PetscErrorCode PETSCDM_DLLEXPORT DMDASetStencilType(DM da, DMDAStencilType stype)
165 {
166   DM_DA *dd = (DM_DA*)da->data;
167 
168   PetscFunctionBegin;
169   PetscValidHeaderSpecific(da,DM_CLASSID,1);
170   PetscValidLogicalCollectiveEnum(da,stype,2);
171   dd->stencil_type = stype;
172   PetscFunctionReturn(0);
173 }
174 
175 #undef __FUNCT__
176 #define __FUNCT__ "DMDASetStencilWidth"
177 /*@
178   DMDASetStencilWidth - Sets the width of the communication stencil
179 
180   Logically Collective on DMDA
181 
182   Input Parameter:
183 + da    - The DMDA
184 - width - The stencil width
185 
186   Level: intermediate
187 
188 .keywords:  distributed array, stencil
189 .seealso: DMDACreate(), DMDestroy(), DMDA
190 @*/
191 PetscErrorCode PETSCDM_DLLEXPORT DMDASetStencilWidth(DM da, PetscInt width)
192 {
193   DM_DA *dd = (DM_DA*)da->data;
194 
195   PetscFunctionBegin;
196   PetscValidHeaderSpecific(da,DM_CLASSID,1);
197   PetscValidLogicalCollectiveInt(da,width,2);
198   dd->s = width;
199   PetscFunctionReturn(0);
200 }
201 
202 #undef __FUNCT__
203 #define __FUNCT__ "DMDACheckOwnershipRanges_Private"
204 static PetscErrorCode DMDACheckOwnershipRanges_Private(DM da,PetscInt M,PetscInt m,const PetscInt lx[])
205 {
206   PetscInt i,sum;
207 
208   PetscFunctionBegin;
209   if (M < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Global dimension not set");
210   for (i=sum=0; i<m; i++) sum += lx[i];
211   if (sum != M) SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_ARG_INCOMP,"Ownership ranges sum to %D but global dimension is %D",sum,M);
212   PetscFunctionReturn(0);
213 }
214 
215 #undef __FUNCT__
216 #define __FUNCT__ "DMDASetOwnershipRanges"
217 /*@
218   DMDASetOwnershipRanges - Sets the number of nodes in each direction on each process
219 
220   Logically Collective on DMDA
221 
222   Input Parameter:
223 + da - The DMDA
224 . lx - array containing number of nodes in the X direction on each process, or PETSC_NULL. If non-null, must be of length da->m
225 . ly - array containing number of nodes in the Y direction on each process, or PETSC_NULL. If non-null, must be of length da->n
226 - lz - array containing number of nodes in the Z direction on each process, or PETSC_NULL. If non-null, must be of length da->p.
227 
228   Level: intermediate
229 
230 .keywords:  distributed array
231 .seealso: DMDACreate(), DMDestroy(), DMDA
232 @*/
233 PetscErrorCode PETSCDM_DLLEXPORT DMDASetOwnershipRanges(DM da, const PetscInt lx[], const PetscInt ly[], const PetscInt lz[])
234 {
235   PetscErrorCode ierr;
236   DM_DA          *dd = (DM_DA*)da->data;
237 
238   PetscFunctionBegin;
239   PetscValidHeaderSpecific(da,DM_CLASSID,1);
240   if (lx) {
241     if (dd->m < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs");
242     ierr = DMDACheckOwnershipRanges_Private(da,dd->M,dd->m,lx);CHKERRQ(ierr);
243     if (!dd->lx) {
244       ierr = PetscMalloc(dd->m*sizeof(PetscInt), &dd->lx);CHKERRQ(ierr);
245     }
246     ierr = PetscMemcpy(dd->lx, lx, dd->m*sizeof(PetscInt));CHKERRQ(ierr);
247   }
248   if (ly) {
249     if (dd->n < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs");
250     ierr = DMDACheckOwnershipRanges_Private(da,dd->N,dd->n,ly);CHKERRQ(ierr);
251     if (!dd->ly) {
252       ierr = PetscMalloc(dd->n*sizeof(PetscInt), &dd->ly);CHKERRQ(ierr);
253     }
254     ierr = PetscMemcpy(dd->ly, ly, dd->n*sizeof(PetscInt));CHKERRQ(ierr);
255   }
256   if (lz) {
257     if (dd->p < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot set ownership ranges before setting number of procs");
258     ierr = DMDACheckOwnershipRanges_Private(da,dd->P,dd->p,lz);CHKERRQ(ierr);
259     if (!dd->lz) {
260       ierr = PetscMalloc(dd->p*sizeof(PetscInt), &dd->lz);CHKERRQ(ierr);
261     }
262     ierr = PetscMemcpy(dd->lz, lz, dd->p*sizeof(PetscInt));CHKERRQ(ierr);
263   }
264   PetscFunctionReturn(0);
265 }
266 
267 #undef __FUNCT__
268 #define __FUNCT__ "DMDACreateOwnershipRanges"
269 /*
270  Ensure that da->lx, ly, and lz exist.  Collective on DMDA.
271 */
272 PetscErrorCode PETSCDM_DLLEXPORT DMDACreateOwnershipRanges(DM da)
273 {
274   DM_DA          *dd = (DM_DA*)da->data;
275   PetscErrorCode ierr;
276   PetscInt       n;
277   MPI_Comm       comm;
278   PetscMPIInt    size;
279 
280   PetscFunctionBegin;
281   if (!dd->lx) {
282     ierr = PetscMalloc(dd->m*sizeof(PetscInt),&dd->lx);CHKERRQ(ierr);
283     ierr = DMDAGetProcessorSubset(da,DMDA_X,dd->xs,&comm);CHKERRQ(ierr);
284     ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
285     n = (dd->xe-dd->xs)/dd->w;
286     ierr = MPI_Allgather(&n,1,MPIU_INT,dd->lx,1,MPIU_INT,comm);CHKERRQ(ierr);
287   }
288   if (dd->dim > 1 && !dd->ly) {
289     ierr = PetscMalloc(dd->n*sizeof(PetscInt),&dd->ly);CHKERRQ(ierr);
290     ierr = DMDAGetProcessorSubset(da,DMDA_Y,dd->ys,&comm);CHKERRQ(ierr);
291     n = dd->ye-dd->ys;
292     ierr = MPI_Allgather(&n,1,MPIU_INT,dd->ly,1,MPIU_INT,comm);CHKERRQ(ierr);
293   }
294   if (dd->dim > 2 && !dd->lz) {
295     ierr = PetscMalloc(dd->p*sizeof(PetscInt),&dd->lz);CHKERRQ(ierr);
296     ierr = DMDAGetProcessorSubset(da,DMDA_Z,dd->zs,&comm);CHKERRQ(ierr);
297     n = dd->ze-dd->zs;
298     ierr = MPI_Allgather(&n,1,MPIU_INT,dd->lz,1,MPIU_INT,comm);CHKERRQ(ierr);
299   }
300   PetscFunctionReturn(0);
301 }
302 
303 #undef __FUNCT__
304 #define __FUNCT__ "DMDASetInterpolationType"
305 /*@
306        DMDASetInterpolationType - Sets the type of interpolation that will be
307           returned by DMGetInterpolation()
308 
309    Logically Collective on DMDA
310 
311    Input Parameter:
312 +  da - initial distributed array
313 .  ctype - DMDA_Q1 and DMDA_Q0 are currently the only supported forms
314 
315    Level: intermediate
316 
317    Notes: you should call this on the coarser of the two DMDAs you pass to DMGetInterpolation()
318 
319 .keywords:  distributed array, interpolation
320 
321 .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDA, DMDAInterpolationType
322 @*/
323 PetscErrorCode PETSCDM_DLLEXPORT DMDASetInterpolationType(DM da,DMDAInterpolationType ctype)
324 {
325   DM_DA *dd = (DM_DA*)da->data;
326 
327   PetscFunctionBegin;
328   PetscValidHeaderSpecific(da,DM_CLASSID,1);
329   PetscValidLogicalCollectiveEnum(da,ctype,2);
330   dd->interptype = ctype;
331   PetscFunctionReturn(0);
332 }
333 
334 
335 #undef __FUNCT__
336 #define __FUNCT__ "DMDAGetNeighbors"
337 /*@C
338       DMDAGetNeighbors - Gets an array containing the MPI rank of all the current
339         processes neighbors.
340 
341     Not Collective
342 
343    Input Parameter:
344 .     da - the DMDA object
345 
346    Output Parameters:
347 .     ranks - the neighbors ranks, stored with the x index increasing most rapidly.
348               this process itself is in the list
349 
350    Notes: In 2d the array is of length 9, in 3d of length 27
351           Not supported in 1d
352           Do not free the array, it is freed when the DMDA is destroyed.
353 
354    Fortran Notes: In fortran you must pass in an array of the appropriate length.
355 
356    Level: intermediate
357 
358 @*/
359 PetscErrorCode PETSCDM_DLLEXPORT DMDAGetNeighbors(DM da,const PetscMPIInt *ranks[])
360 {
361   DM_DA *dd = (DM_DA*)da->data;
362   PetscFunctionBegin;
363   PetscValidHeaderSpecific(da,DM_CLASSID,1);
364   *ranks = dd->neighbors;
365   PetscFunctionReturn(0);
366 }
367 
368 /*@C
369       DMDASetElementType - Sets the element type to be returned by DMGetElements()
370 
371     Not Collective
372 
373    Input Parameter:
374 .     da - the DMDA object
375 
376    Output Parameters:
377 .     etype - the element type, currently either DMDA_ELEMENT_P1 or ELEMENT_Q1
378 
379    Level: intermediate
380 
381 .seealso: DMDAElementType, DMDAGetElementType(), DMGetElements(), DMRestoreElements()
382 @*/
383 #undef __FUNCT__
384 #define __FUNCT__ "DMDASetElementType"
385 PetscErrorCode PETSCDM_DLLEXPORT DMDASetElementType(DM da, DMDAElementType etype)
386 {
387   DM_DA          *dd = (DM_DA*)da->data;
388   PetscErrorCode ierr;
389 
390   PetscFunctionBegin;
391   PetscValidHeaderSpecific(da,DM_CLASSID,1);
392   PetscValidLogicalCollectiveEnum(da,etype,2);
393   if (dd->elementtype != etype) {
394     ierr = PetscFree(dd->e);CHKERRQ(ierr);
395     dd->elementtype = etype;
396     dd->ne          = 0;
397     dd->e           = PETSC_NULL;
398   }
399   PetscFunctionReturn(0);
400 }
401 
402 /*@C
403       DMDAGetElementType - Gets the element type to be returned by DMGetElements()
404 
405     Not Collective
406 
407    Input Parameter:
408 .     da - the DMDA object
409 
410    Output Parameters:
411 .     etype - the element type, currently either DMDA_ELEMENT_P1 or ELEMENT_Q1
412 
413    Level: intermediate
414 
415 .seealso: DMDAElementType, DMDASetElementType(), DMGetElements(), DMRestoreElements()
416 @*/
417 #undef __FUNCT__
418 #define __FUNCT__ "DMDAGetElementType"
419 PetscErrorCode PETSCDM_DLLEXPORT DMDAGetElementType(DM da, DMDAElementType *etype)
420 {
421   DM_DA *dd = (DM_DA*)da->data;
422   PetscFunctionBegin;
423   PetscValidHeaderSpecific(da,DM_CLASSID,1);
424   PetscValidPointer(etype,2);
425   *etype = dd->elementtype;
426   PetscFunctionReturn(0);
427 }
428 
429 #undef __FUNCT__
430 #define __FUNCT__ "DMGetElements"
431 /*@C
432       DMGetElements - Gets an array containing the indices (in local coordinates)
433                  of all the local elements
434 
435     Not Collective
436 
437    Input Parameter:
438 .     dm - the DM object
439 
440    Output Parameters:
441 +     nel - number of local elements
442 .     nen - number of element nodes
443 -     e - the indices of the elements vertices
444 
445    Level: intermediate
446 
447 .seealso: DMElementType, DMSetElementType(), DMRestoreElements()
448 @*/
449 PetscErrorCode PETSCDM_DLLEXPORT DMGetElements(DM dm,PetscInt *nel,PetscInt *nen,const PetscInt *e[])
450 {
451   PetscErrorCode ierr;
452   PetscFunctionBegin;
453   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
454   PetscValidIntPointer(nel,2);
455   PetscValidIntPointer(nen,3);
456   PetscValidPointer(e,4);
457   ierr = (dm->ops->getelements)(dm,nel,nen,e);CHKERRQ(ierr);
458   PetscFunctionReturn(0);
459 }
460 
461 #undef __FUNCT__
462 #define __FUNCT__ "DMRestoreElements"
463 /*@C
464       DMRestoreElements - Returns an array containing the indices (in local coordinates)
465                  of all the local elements obtained with DMGetElements()
466 
467     Not Collective
468 
469    Input Parameter:
470 +     dm - the DM object
471 .     nel - number of local elements
472 .     nen - number of element nodes
473 -     e - the indices of the elements vertices
474 
475    Level: intermediate
476 
477 .seealso: DMElementType, DMSetElementType(), DMGetElements()
478 @*/
479 PetscErrorCode PETSCDM_DLLEXPORT DMRestoreElements(DM dm,PetscInt *nel,PetscInt *nen,const PetscInt *e[])
480 {
481   PetscErrorCode ierr;
482   PetscFunctionBegin;
483   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
484   PetscValidIntPointer(nel,2);
485   PetscValidIntPointer(nen,3);
486   PetscValidPointer(e,4);
487   if (dm->ops->restoreelements) {
488     ierr = (dm->ops->restoreelements)(dm,nel,nen,e);CHKERRQ(ierr);
489   }
490   PetscFunctionReturn(0);
491 }
492 
493 #undef __FUNCT__
494 #define __FUNCT__ "DMDAGetOwnershipRanges"
495 /*@C
496       DMDAGetOwnershipRanges - Gets the ranges of indices in the x, y and z direction that are owned by each process
497 
498     Not Collective
499 
500    Input Parameter:
501 .     da - the DMDA object
502 
503    Output Parameter:
504 +     lx - ownership along x direction (optional)
505 .     ly - ownership along y direction (optional)
506 -     lz - ownership along z direction (optional)
507 
508    Level: intermediate
509 
510     Note: these correspond to the optional final arguments passed to DMDACreate(), DMDACreate2d(), DMDACreate3d()
511 
512     In Fortran one must pass in arrays lx, ly, and lz that are long enough to hold the values; the sixth, seventh and
513     eighth arguments from DMDAGetInfo()
514 
515      In C you should not free these arrays, nor change the values in them. They will only have valid values while the
516     DMDA they came from still exists (has not been destroyed).
517 
518 .seealso: DMDAGetCorners(), DMDAGetGhostCorners(), DMDACreate(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), VecGetOwnershipRanges()
519 @*/
520 PetscErrorCode PETSCDM_DLLEXPORT DMDAGetOwnershipRanges(DM da,const PetscInt *lx[],const PetscInt *ly[],const PetscInt *lz[])
521 {
522   DM_DA *dd = (DM_DA*)da->data;
523 
524   PetscFunctionBegin;
525   PetscValidHeaderSpecific(da,DM_CLASSID,1);
526   if (lx) *lx = dd->lx;
527   if (ly) *ly = dd->ly;
528   if (lz) *lz = dd->lz;
529   PetscFunctionReturn(0);
530 }
531 
532 #undef __FUNCT__
533 #define __FUNCT__ "DMDASetRefinementFactor"
534 /*@
535      DMDASetRefinementFactor - Set the ratios that the DMDA grid is refined
536 
537     Logically Collective on DMDA
538 
539   Input Parameters:
540 +    da - the DMDA object
541 .    refine_x - ratio of fine grid to coarse in x direction (2 by default)
542 .    refine_y - ratio of fine grid to coarse in y direction (2 by default)
543 -    refine_z - ratio of fine grid to coarse in z direction (2 by default)
544 
545   Options Database:
546 +  -da_refine_x - refinement ratio in x direction
547 .  -da_refine_y - refinement ratio in y direction
548 -  -da_refine_z - refinement ratio in z direction
549 
550   Level: intermediate
551 
552     Notes: Pass PETSC_IGNORE to leave a value unchanged
553 
554 .seealso: DMRefine(), DMDAGetRefinementFactor()
555 @*/
556 PetscErrorCode PETSCDM_DLLEXPORT DMDASetRefinementFactor(DM da, PetscInt refine_x, PetscInt refine_y,PetscInt refine_z)
557 {
558   DM_DA *dd = (DM_DA*)da->data;
559 
560   PetscFunctionBegin;
561   PetscValidHeaderSpecific(da,DM_CLASSID,1);
562   PetscValidLogicalCollectiveInt(da,refine_x,2);
563   PetscValidLogicalCollectiveInt(da,refine_y,3);
564   PetscValidLogicalCollectiveInt(da,refine_z,4);
565 
566   if (refine_x > 0) dd->refine_x = refine_x;
567   if (refine_y > 0) dd->refine_y = refine_y;
568   if (refine_z > 0) dd->refine_z = refine_z;
569   PetscFunctionReturn(0);
570 }
571 
572 #undef __FUNCT__
573 #define __FUNCT__ "DMDAGetRefinementFactor"
574 /*@C
575      DMDAGetRefinementFactor - Gets the ratios that the DMDA grid is refined
576 
577     Not Collective
578 
579   Input Parameter:
580 .    da - the DMDA object
581 
582   Output Parameters:
583 +    refine_x - ratio of fine grid to coarse in x direction (2 by default)
584 .    refine_y - ratio of fine grid to coarse in y direction (2 by default)
585 -    refine_z - ratio of fine grid to coarse in z direction (2 by default)
586 
587   Level: intermediate
588 
589     Notes: Pass PETSC_NULL for values you do not need
590 
591 .seealso: DMRefine(), DMDASetRefinementFactor()
592 @*/
593 PetscErrorCode PETSCDM_DLLEXPORT DMDAGetRefinementFactor(DM da, PetscInt *refine_x, PetscInt *refine_y,PetscInt *refine_z)
594 {
595   DM_DA *dd = (DM_DA*)da->data;
596 
597   PetscFunctionBegin;
598   PetscValidHeaderSpecific(da,DM_CLASSID,1);
599   if (refine_x) *refine_x = dd->refine_x;
600   if (refine_y) *refine_y = dd->refine_y;
601   if (refine_z) *refine_z = dd->refine_z;
602   PetscFunctionReturn(0);
603 }
604 
605 #undef __FUNCT__
606 #define __FUNCT__ "DMDASetGetMatrix"
607 /*@C
608      DMDASetGetMatrix - Sets the routine used by the DMDA to allocate a matrix.
609 
610     Logically Collective on DMDA
611 
612   Input Parameters:
613 +    da - the DMDA object
614 -    f - the function that allocates the matrix for that specific DMDA
615 
616   Level: developer
617 
618    Notes: See DMDASetBlockFills() that provides a simple way to provide the nonzero structure for
619        the diagonal and off-diagonal blocks of the matrix
620 
621 .seealso: DMGetMatrix(), DMDASetBlockFills()
622 @*/
623 PetscErrorCode PETSCDM_DLLEXPORT DMDASetGetMatrix(DM da,PetscErrorCode (*f)(DM, const MatType,Mat*))
624 {
625   PetscFunctionBegin;
626   PetscValidHeaderSpecific(da,DM_CLASSID,1);
627   da->ops->getmatrix = f;
628   PetscFunctionReturn(0);
629 }
630 
631 #undef __FUNCT__
632 #define __FUNCT__ "DMDARefineOwnershipRanges"
633 /*
634   Creates "balanced" ownership ranges after refinement, constrained by the need for the
635   fine grid boundaries to fall within one stencil width of the coarse partition.
636 
637   Uses a greedy algorithm to handle non-ideal layouts, could probably do something better.
638 */
639 static PetscErrorCode DMDARefineOwnershipRanges(DM da,PetscBool periodic,PetscInt stencil_width,PetscInt ratio,PetscInt m,const PetscInt lc[],PetscInt lf[])
640 {
641   PetscInt       i,totalc = 0,remaining,startc = 0,startf = 0;
642   PetscErrorCode ierr;
643 
644   PetscFunctionBegin;
645   if (ratio < 1) SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_USER,"Requested refinement ratio %D must be at least 1",ratio);
646   if (ratio == 1) {
647     ierr = PetscMemcpy(lf,lc,m*sizeof(lc[0]));CHKERRQ(ierr);
648     PetscFunctionReturn(0);
649   }
650   for (i=0; i<m; i++) totalc += lc[i];
651   remaining = (!periodic) + ratio * (totalc - (!periodic));
652   for (i=0; i<m; i++) {
653     PetscInt want = remaining/(m-i) + !!(remaining%(m-i));
654     if (i == m-1) lf[i] = want;
655     else {
656       const PetscInt nextc = startc + lc[i];
657       /* Move the first fine node of the next subdomain to the right until the coarse node on its left is within one
658        * coarse stencil width of the first coarse node in the next subdomain. */
659       while ((startf+want)/ratio < nextc - stencil_width) want++;
660       /* Move the last fine node in the current subdomain to the left until the coarse node on its right is within one
661        * coarse stencil width of the last coarse node in the current subdomain. */
662       while ((startf+want-1+ratio-1)/ratio > nextc-1+stencil_width) want--;
663       /* Make sure all constraints are satisfied */
664       if (want < 0 || want > remaining
665           || ((startf+want)/ratio < nextc - stencil_width)
666           || ((startf+want-1+ratio-1)/ratio > nextc-1+stencil_width))
667         SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_SIZ,"Could not find a compatible refined ownership range");
668     }
669     lf[i] = want;
670     startc += lc[i];
671     startf += lf[i];
672     remaining -= lf[i];
673   }
674   PetscFunctionReturn(0);
675 }
676 
677 #undef __FUNCT__
678 #define __FUNCT__ "DMDACoarsenOwnershipRanges"
679 /*
680   Creates "balanced" ownership ranges after coarsening, constrained by the need for the
681   fine grid boundaries to fall within one stencil width of the coarse partition.
682 
683   Uses a greedy algorithm to handle non-ideal layouts, could probably do something better.
684 */
685 static PetscErrorCode DMDACoarsenOwnershipRanges(DM da,PetscBool periodic,PetscInt stencil_width,PetscInt ratio,PetscInt m,const PetscInt lf[],PetscInt lc[])
686 {
687   PetscInt       i,totalf,remaining,startc,startf;
688   PetscErrorCode ierr;
689 
690   PetscFunctionBegin;
691   if (ratio < 1) SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_USER,"Requested refinement ratio %D must be at least 1",ratio);
692   if (ratio == 1) {
693     ierr = PetscMemcpy(lc,lf,m*sizeof(lf[0]));CHKERRQ(ierr);
694     PetscFunctionReturn(0);
695   }
696   for (i=0,totalf=0; i<m; i++) totalf += lf[i];
697   remaining = (!periodic) + (totalf - (!periodic)) / ratio;
698   for (i=0,startc=0,startf=0; i<m; i++) {
699     PetscInt want = remaining/(m-i) + !!(remaining%(m-i));
700     if (i == m-1) lc[i] = want;
701     else {
702       const PetscInt nextf = startf+lf[i];
703       /* Slide first coarse node of next subdomain to the left until the coarse node to the left of the first fine
704        * node is within one stencil width. */
705       while (nextf/ratio < startc+want-stencil_width) want--;
706       /* Slide the last coarse node of the current subdomain to the right until the coarse node to the right of the last
707        * fine node is within one stencil width. */
708       while ((nextf-1+ratio-1)/ratio > startc+want-1+stencil_width) want++;
709       if (want < 0 || want > remaining
710           || (nextf/ratio < startc+want-stencil_width)
711           || ((nextf-1+ratio-1)/ratio > startc+want-1+stencil_width))
712         SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_SIZ,"Could not find a compatible coarsened ownership range");
713     }
714     lc[i] = want;
715     startc += lc[i];
716     startf += lf[i];
717     remaining -= lc[i];
718   }
719   PetscFunctionReturn(0);
720 }
721 
722 #undef __FUNCT__
723 #define __FUNCT__ "DMRefine_DA"
724 PetscErrorCode PETSCDM_DLLEXPORT DMRefine_DA(DM da,MPI_Comm comm,DM *daref)
725 {
726   PetscErrorCode ierr;
727   PetscInt       M,N,P;
728   DM             da2;
729   DM_DA          *dd = (DM_DA*)da->data,*dd2;
730 
731   PetscFunctionBegin;
732   PetscValidHeaderSpecific(da,DM_CLASSID,1);
733   PetscValidPointer(daref,3);
734 
735   if (DMDAXPeriodic(dd->wrap) || dd->interptype == DMDA_Q0){
736     M = dd->refine_x*dd->M;
737   } else {
738     M = 1 + dd->refine_x*(dd->M - 1);
739   }
740   if (DMDAYPeriodic(dd->wrap) || dd->interptype == DMDA_Q0){
741     N = dd->refine_y*dd->N;
742   } else {
743     N = 1 + dd->refine_y*(dd->N - 1);
744   }
745   if (DMDAZPeriodic(dd->wrap) || dd->interptype == DMDA_Q0){
746     P = dd->refine_z*dd->P;
747   } else {
748     P = 1 + dd->refine_z*(dd->P - 1);
749   }
750   ierr = DMDACreateOwnershipRanges(da);CHKERRQ(ierr);
751   if (dd->dim == 3) {
752     PetscInt *lx,*ly,*lz;
753     ierr = PetscMalloc3(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly,dd->p,PetscInt,&lz);CHKERRQ(ierr);
754     ierr = DMDARefineOwnershipRanges(da,(PetscBool)(DMDAXPeriodic(dd->wrap) || dd->interptype == DMDA_Q0),dd->s,dd->refine_x,dd->m,dd->lx,lx);CHKERRQ(ierr);
755     ierr = DMDARefineOwnershipRanges(da,(PetscBool)(DMDAYPeriodic(dd->wrap) || dd->interptype == DMDA_Q0),dd->s,dd->refine_y,dd->n,dd->ly,ly);CHKERRQ(ierr);
756     ierr = DMDARefineOwnershipRanges(da,(PetscBool)(DMDAZPeriodic(dd->wrap) || dd->interptype == DMDA_Q0),dd->s,dd->refine_z,dd->p,dd->lz,lz);CHKERRQ(ierr);
757     ierr = DMDACreate3d(((PetscObject)da)->comm,dd->wrap,dd->stencil_type,M,N,P,dd->m,dd->n,dd->p,dd->w,dd->s,lx,ly,lz,&da2);CHKERRQ(ierr);
758     ierr = PetscFree3(lx,ly,lz);CHKERRQ(ierr);
759   } else if (dd->dim == 2) {
760     PetscInt *lx,*ly;
761     ierr = PetscMalloc2(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly);CHKERRQ(ierr);
762     ierr = DMDARefineOwnershipRanges(da,(PetscBool)(DMDAXPeriodic(dd->wrap) || dd->interptype == DMDA_Q0),dd->s,dd->refine_x,dd->m,dd->lx,lx);CHKERRQ(ierr);
763     ierr = DMDARefineOwnershipRanges(da,(PetscBool)(DMDAYPeriodic(dd->wrap) || dd->interptype == DMDA_Q0),dd->s,dd->refine_y,dd->n,dd->ly,ly);CHKERRQ(ierr);
764     ierr = DMDACreate2d(((PetscObject)da)->comm,dd->wrap,dd->stencil_type,M,N,dd->m,dd->n,dd->w,dd->s,lx,ly,&da2);CHKERRQ(ierr);
765     ierr = PetscFree2(lx,ly);CHKERRQ(ierr);
766   } else if (dd->dim == 1) {
767     PetscInt *lx;
768     ierr = PetscMalloc(dd->m*sizeof(PetscInt),&lx);CHKERRQ(ierr);
769     ierr = DMDARefineOwnershipRanges(da,(PetscBool)(DMDAXPeriodic(dd->wrap) || dd->interptype == DMDA_Q0),dd->s,dd->refine_x,dd->m,dd->lx,lx);CHKERRQ(ierr);
770     ierr = DMDACreate1d(((PetscObject)da)->comm,dd->wrap,M,dd->w,dd->s,lx,&da2);CHKERRQ(ierr);
771     ierr = PetscFree(lx);CHKERRQ(ierr);
772   }
773   dd2 = (DM_DA*)da2->data;
774 
775   /* allow overloaded (user replaced) operations to be inherited by refinement clones */
776   da2->ops->getmatrix        = da->ops->getmatrix;
777   da2->ops->getinterpolation = da->ops->getinterpolation;
778   da2->ops->getcoloring      = da->ops->getcoloring;
779   dd2->interptype            = dd->interptype;
780 
781   /* copy fill information if given */
782   if (dd->dfill) {
783     ierr = PetscMalloc((dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->dfill);CHKERRQ(ierr);
784     ierr = PetscMemcpy(dd2->dfill,dd->dfill,(dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr);
785   }
786   if (dd->ofill) {
787     ierr = PetscMalloc((dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->ofill);CHKERRQ(ierr);
788     ierr = PetscMemcpy(dd2->ofill,dd->ofill,(dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr);
789   }
790   /* copy the refine information */
791   dd2->refine_x = dd->refine_x;
792   dd2->refine_y = dd->refine_y;
793   dd2->refine_z = dd->refine_z;
794 
795   /* copy vector type information */
796   ierr = PetscFree(da2->vectype);CHKERRQ(ierr);
797   ierr = PetscStrallocpy(da->vectype,&da2->vectype);CHKERRQ(ierr);
798   *daref = da2;
799   PetscFunctionReturn(0);
800 }
801 
802 #undef __FUNCT__
803 #define __FUNCT__ "DMCoarsen_DA"
804 PetscErrorCode PETSCDM_DLLEXPORT DMCoarsen_DA(DM da, MPI_Comm comm,DM *daref)
805 {
806   PetscErrorCode ierr;
807   PetscInt       M,N,P;
808   DM             da2;
809   DM_DA          *dd = (DM_DA*)da->data,*dd2;
810 
811   PetscFunctionBegin;
812   PetscValidHeaderSpecific(da,DM_CLASSID,1);
813   PetscValidPointer(daref,3);
814 
815   if (DMDAXPeriodic(dd->wrap) || dd->interptype == DMDA_Q0){
816     M = dd->M / dd->refine_x;
817   } else {
818     M = 1 + (dd->M - 1) / dd->refine_x;
819   }
820   if (DMDAYPeriodic(dd->wrap) || dd->interptype == DMDA_Q0){
821     N = dd->N / dd->refine_y;
822   } else {
823     N = 1 + (dd->N - 1) / dd->refine_y;
824   }
825   if (DMDAZPeriodic(dd->wrap) || dd->interptype == DMDA_Q0){
826     P = dd->P / dd->refine_z;
827   } else {
828     P = 1 + (dd->P - 1) / dd->refine_z;
829   }
830   ierr = DMDACreateOwnershipRanges(da);CHKERRQ(ierr);
831   if (dd->dim == 3) {
832     PetscInt *lx,*ly,*lz;
833     ierr = PetscMalloc3(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly,dd->p,PetscInt,&lz);CHKERRQ(ierr);
834     ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(DMDAXPeriodic(dd->wrap) || dd->interptype == DMDA_Q0),dd->s,dd->refine_x,dd->m,dd->lx,lx);CHKERRQ(ierr);
835     ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(DMDAYPeriodic(dd->wrap) || dd->interptype == DMDA_Q0),dd->s,dd->refine_y,dd->n,dd->ly,ly);CHKERRQ(ierr);
836     ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(DMDAZPeriodic(dd->wrap) || dd->interptype == DMDA_Q0),dd->s,dd->refine_z,dd->p,dd->lz,lz);CHKERRQ(ierr);
837     ierr = DMDACreate3d(((PetscObject)da)->comm,dd->wrap,dd->stencil_type,M,N,P,dd->m,dd->n,dd->p,dd->w,dd->s,lx,ly,lz,&da2);CHKERRQ(ierr);
838     ierr = PetscFree3(lx,ly,lz);CHKERRQ(ierr);
839   } else if (dd->dim == 2) {
840     PetscInt *lx,*ly;
841     ierr = PetscMalloc2(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly);CHKERRQ(ierr);
842     ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(DMDAXPeriodic(dd->wrap) || dd->interptype == DMDA_Q0),dd->s,dd->refine_x,dd->m,dd->lx,lx);CHKERRQ(ierr);
843     ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(DMDAYPeriodic(dd->wrap) || dd->interptype == DMDA_Q0),dd->s,dd->refine_y,dd->n,dd->ly,ly);CHKERRQ(ierr);
844     ierr = DMDACreate2d(((PetscObject)da)->comm,dd->wrap,dd->stencil_type,M,N,dd->m,dd->n,dd->w,dd->s,lx,ly,&da2);CHKERRQ(ierr);
845     ierr = PetscFree2(lx,ly);CHKERRQ(ierr);
846   } else if (dd->dim == 1) {
847     PetscInt *lx;
848     ierr = PetscMalloc(dd->m*sizeof(PetscInt),&lx);CHKERRQ(ierr);
849     ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(DMDAXPeriodic(dd->wrap) || dd->interptype == DMDA_Q0),dd->s,dd->refine_x,dd->m,dd->lx,lx);CHKERRQ(ierr);
850     ierr = DMDACreate1d(((PetscObject)da)->comm,dd->wrap,M,dd->w,dd->s,lx,&da2);CHKERRQ(ierr);
851     ierr = PetscFree(lx);CHKERRQ(ierr);
852   }
853   dd2 = (DM_DA*)da2->data;
854 
855   /* allow overloaded (user replaced) operations to be inherited by refinement clones */
856   da2->ops->getmatrix        = da->ops->getmatrix;
857   da2->ops->getinterpolation = da->ops->getinterpolation;
858   da2->ops->getcoloring      = da->ops->getcoloring;
859   dd2->interptype            = dd->interptype;
860 
861   /* copy fill information if given */
862   if (dd->dfill) {
863     ierr = PetscMalloc((dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->dfill);CHKERRQ(ierr);
864     ierr = PetscMemcpy(dd2->dfill,dd->dfill,(dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr);
865   }
866   if (dd->ofill) {
867     ierr = PetscMalloc((dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->ofill);CHKERRQ(ierr);
868     ierr = PetscMemcpy(dd2->ofill,dd->ofill,(dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr);
869   }
870   /* copy the refine information */
871   dd2->refine_x = dd->refine_x;
872   dd2->refine_y = dd->refine_y;
873   dd2->refine_z = dd->refine_z;
874 
875   /* copy vector type information */
876   ierr = PetscFree(da2->vectype);CHKERRQ(ierr);
877   ierr = PetscStrallocpy(da->vectype,&da2->vectype);CHKERRQ(ierr);
878 
879   *daref = da2;
880   PetscFunctionReturn(0);
881 }
882 
883 #undef __FUNCT__
884 #define __FUNCT__ "DMRefineHierarchy_DA"
885 PetscErrorCode PETSCDM_DLLEXPORT DMRefineHierarchy_DA(DM da,PetscInt nlevels,DM daf[])
886 {
887   PetscErrorCode ierr;
888   PetscInt       i,n,*refx,*refy,*refz;
889 
890   PetscFunctionBegin;
891   PetscValidHeaderSpecific(da,DM_CLASSID,1);
892   if (nlevels < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
893   if (nlevels == 0) PetscFunctionReturn(0);
894   PetscValidPointer(daf,3);
895 
896   /* Get refinement factors, defaults taken from the coarse DMDA */
897   ierr = PetscMalloc3(nlevels,PetscInt,&refx,nlevels,PetscInt,&refy,nlevels,PetscInt,&refz);CHKERRQ(ierr);
898   for (i=0; i<nlevels; i++) {
899     ierr = DMDAGetRefinementFactor(da,&refx[i],&refy[i],&refz[i]);CHKERRQ(ierr);
900   }
901   n = nlevels;
902   ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_x",refx,&n,PETSC_NULL);CHKERRQ(ierr);
903   n = nlevels;
904   ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_y",refy,&n,PETSC_NULL);CHKERRQ(ierr);
905   n = nlevels;
906   ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_z",refz,&n,PETSC_NULL);CHKERRQ(ierr);
907 
908   ierr = DMDASetRefinementFactor(da,refx[0],refy[0],refz[0]);CHKERRQ(ierr);
909   ierr = DMRefine(da,((PetscObject)da)->comm,&daf[0]);CHKERRQ(ierr);
910   for (i=1; i<nlevels; i++) {
911     ierr = DMDASetRefinementFactor(daf[i-1],refx[i],refy[i],refz[i]);CHKERRQ(ierr);
912     ierr = DMRefine(daf[i-1],((PetscObject)da)->comm,&daf[i]);CHKERRQ(ierr);
913   }
914   ierr = PetscFree3(refx,refy,refz);CHKERRQ(ierr);
915   PetscFunctionReturn(0);
916 }
917 
918 #undef __FUNCT__
919 #define __FUNCT__ "DMCoarsenHierarchy_DA"
920 PetscErrorCode PETSCDM_DLLEXPORT DMCoarsenHierarchy_DA(DM da,PetscInt nlevels,DM dac[])
921 {
922   PetscErrorCode ierr;
923   PetscInt       i;
924 
925   PetscFunctionBegin;
926   PetscValidHeaderSpecific(da,DM_CLASSID,1);
927   if (nlevels < 0) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative");
928   if (nlevels == 0) PetscFunctionReturn(0);
929   PetscValidPointer(dac,3);
930   ierr = DMCoarsen(da,((PetscObject)da)->comm,&dac[0]);CHKERRQ(ierr);
931   for (i=1; i<nlevels; i++) {
932     ierr = DMCoarsen(dac[i-1],((PetscObject)da)->comm,&dac[i]);CHKERRQ(ierr);
933   }
934   PetscFunctionReturn(0);
935 }
936