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