xref: /petsc/src/dm/impls/da/da1.c (revision cffb1e400d6c3ea2f6cb522ae2432dc42cf29e73)
1 
2 /*
3    Code for manipulating distributed regular 1d arrays in parallel.
4    This file was created by Peter Mell   6/30/95
5 */
6 
7 #include <petsc-private/daimpl.h>     /*I  "petscdmda.h"   I*/
8 
9 const char *const DMDABoundaryTypes[] = {"BOUNDARY_NONE","BOUNDARY_GHOSTED","BOUNDARY_PERIODIC","DMDA_",0};
10 
11 #undef __FUNCT__
12 #define __FUNCT__ "DMView_DA_1d"
13 PetscErrorCode DMView_DA_1d(DM da,PetscViewer viewer)
14 {
15   PetscErrorCode ierr;
16   PetscMPIInt    rank;
17   PetscBool      iascii,isdraw,isbinary;
18   DM_DA          *dd = (DM_DA*)da->data;
19 #if defined(PETSC_HAVE_MATLAB_ENGINE)
20   PetscBool      ismatlab;
21 #endif
22 
23   PetscFunctionBegin;
24   ierr = MPI_Comm_rank(((PetscObject)da)->comm,&rank);CHKERRQ(ierr);
25 
26   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
27   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
28   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
29 #if defined(PETSC_HAVE_MATLAB_ENGINE)
30   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERMATLAB,&ismatlab);CHKERRQ(ierr);
31 #endif
32   if (iascii) {
33     PetscViewerFormat format;
34 
35     ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr);
36     if (format != PETSC_VIEWER_ASCII_VTK && format != PETSC_VIEWER_ASCII_VTK_CELL) {
37       DMDALocalInfo info;
38       ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
39       ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
40       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Processor [%d] M %D m %D w %D s %D\n",rank,dd->M,dd->m,dd->w,dd->s);CHKERRQ(ierr);
41       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"X range of indices: %D %D\n",info.xs,info.xs+info.xm);CHKERRQ(ierr);
42       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
43       ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr);
44     } else {
45       ierr = DMView_DA_VTK(da, viewer);CHKERRQ(ierr);
46     }
47   } else if (isdraw) {
48     PetscDraw  draw;
49     double     ymin = -1,ymax = 1,xmin = -1,xmax = dd->M,x;
50     PetscInt   base;
51     char       node[10];
52     PetscBool  isnull;
53 
54     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
55     ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
56 
57     ierr = PetscDrawSetCoordinates(draw,xmin,ymin,xmax,ymax);CHKERRQ(ierr);
58     ierr = PetscDrawSynchronizedClear(draw);CHKERRQ(ierr);
59 
60     /* first processor draws all node lines */
61     if (!rank) {
62       PetscInt xmin_tmp;
63       ymin = 0.0; ymax = 0.3;
64 
65       /* ADIC doesn't like doubles in a for loop */
66       for (xmin_tmp =0; xmin_tmp < dd->M; xmin_tmp++) {
67          ierr = PetscDrawLine(draw,(double)xmin_tmp,ymin,(double)xmin_tmp,ymax,PETSC_DRAW_BLACK);CHKERRQ(ierr);
68       }
69 
70       xmin = 0.0; xmax = dd->M - 1;
71       ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_BLACK);CHKERRQ(ierr);
72       ierr = PetscDrawLine(draw,xmin,ymax,xmax,ymax,PETSC_DRAW_BLACK);CHKERRQ(ierr);
73     }
74 
75     ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr);
76     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
77 
78     /* draw my box */
79     ymin = 0; ymax = 0.3; xmin = dd->xs / dd->w; xmax = (dd->xe / dd->w)  - 1;
80     ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_RED);CHKERRQ(ierr);
81     ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_RED);CHKERRQ(ierr);
82     ierr = PetscDrawLine(draw,xmin,ymax,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr);
83     ierr = PetscDrawLine(draw,xmax,ymin,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr);
84 
85     /* Put in index numbers */
86     base = dd->base / dd->w;
87     for (x=xmin; x<=xmax; x++) {
88       sprintf(node,"%d",(int)base++);
89       ierr = PetscDrawString(draw,x,ymin,PETSC_DRAW_RED,node);CHKERRQ(ierr);
90     }
91 
92     ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr);
93     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
94   } else if (isbinary){
95     ierr = DMView_DA_Binary(da,viewer);CHKERRQ(ierr);
96 #if defined(PETSC_HAVE_MATLAB_ENGINE)
97   } else if (ismatlab) {
98     ierr = DMView_DA_Matlab(da,viewer);CHKERRQ(ierr);
99 #endif
100   } else SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for DMDA 1d",((PetscObject)viewer)->type_name);
101   PetscFunctionReturn(0);
102 }
103 
104 
105 #undef __FUNCT__
106 #define __FUNCT__ "DMSetUp_DA_1D"
107 PetscErrorCode  DMSetUp_DA_1D(DM da)
108 {
109   DM_DA            *dd = (DM_DA*)da->data;
110   const PetscInt   M     = dd->M;
111   const PetscInt   dof   = dd->w;
112   const PetscInt   s     = dd->s;
113   const PetscInt   o     = dd->overlap;
114   const PetscInt   sDist = s*dof;  /* absolute stencil distance */
115   const PetscInt   oDist = o*dof;
116   const PetscInt   *lx    = dd->lx;
117   DMDABoundaryType bx  = dd->bx;
118   MPI_Comm         comm;
119   Vec              local, global;
120   VecScatter       ltog, gtol;
121   IS               to, from;
122   PetscBool        flg1 = PETSC_FALSE, flg2 = PETSC_FALSE;
123   PetscMPIInt      rank, size;
124   PetscInt         i,*idx,nn,left,xs,xe,x,Xs,Xe,start,end,m,IXs,IXe;
125   PetscErrorCode   ierr;
126 
127   PetscFunctionBegin;
128   if (dof < 1) SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Must have 1 or more degrees of freedom per node: %D",dof);
129   if (s < 0) SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Stencil width cannot be negative: %D",s);
130 
131   dd->dim = 1;
132   ierr = PetscMalloc(dof*sizeof(char*),&dd->fieldname);CHKERRQ(ierr);
133   ierr = PetscMemzero(dd->fieldname,dof*sizeof(char*));CHKERRQ(ierr);
134   ierr = PetscObjectGetComm((PetscObject) da, &comm);CHKERRQ(ierr);
135   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
136   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
137 
138   dd->m = size;
139   m     = dd->m;
140 
141   if (s > 0) {
142     /* if not communicating data then should be ok to have nothing on some processes */
143     if (M < m)     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"More processes than data points! %D %D",m,M);
144     if ((M-1) < s) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Array is too small for stencil! %D %D",M-1,s);
145   }
146 
147   /*
148      Determine locally owned region
149      xs is the first local node number, x is the number of local nodes
150   */
151   if (!lx) {
152     ierr = PetscMalloc(m*sizeof(PetscInt), &dd->lx);CHKERRQ(ierr);
153     ierr = PetscOptionsGetBool(PETSC_NULL,"-da_partition_blockcomm",&flg1,PETSC_NULL);CHKERRQ(ierr);
154     ierr = PetscOptionsGetBool(PETSC_NULL,"-da_partition_nodes_at_end",&flg2,PETSC_NULL);CHKERRQ(ierr);
155     if (flg1) {      /* Block Comm type Distribution */
156       xs = rank*M/m;
157       x  = (rank + 1)*M/m - xs;
158     } else if (flg2) { /* The odd nodes are evenly distributed across last nodes */
159       x = (M + rank)/m;
160       if (M/m == x) { xs = rank*x; }
161       else          { xs = rank*(x-1) + (M+rank)%(x*m); }
162     } else { /* The odd nodes are evenly distributed across the first k nodes */
163       /* Regular PETSc Distribution */
164       x = M/m + ((M % m) > rank);
165       if (rank >= (M % m)) {xs = (rank * (PetscInt)(M/m) + M % m);}
166       else                 {xs = rank * (PetscInt)(M/m) + rank;}
167     }
168     ierr = MPI_Allgather(&xs,1,MPIU_INT,dd->lx,1,MPIU_INT,comm);CHKERRQ(ierr);
169     for (i=0; i<m-1; i++) dd->lx[i] = dd->lx[i+1] - dd->lx[i];
170     dd->lx[m-1] = M - dd->lx[m-1];
171   } else {
172     x  = lx[rank];
173     xs = 0;
174     for (i=0; i<rank; i++) {
175       xs += lx[i];
176     }
177     /* verify that data user provided is consistent */
178     left = xs;
179     for (i=rank; i<size; i++) {
180       left += lx[i];
181     }
182     if (left != M) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Sum of lx across processors not equal to M %D %D",left,M);
183   }
184 
185   /*
186    check if the scatter requires more than one process neighbor or wraps around
187    the domain more than once
188   */
189   if ((x < s+o) & ((M > 1) | (bx == DMDA_BOUNDARY_PERIODIC))) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local x-width of domain x %D is smaller than stencil width s %D",x,s+o);
190 
191   /* From now on x,xs,xe,Xs,Xe are the exact location in the array */
192   x  *= dof;
193   xs *= dof;
194   xe  = xs + x;
195 
196   /* determine ghost region (Xs) and region scattered into (IXs)  */
197   if (xs-sDist-oDist > 0) {
198     Xs = xs - sDist - oDist;
199     IXs = xs - sDist - oDist;
200   } else {
201     if (bx) {
202       Xs = xs - sDist;
203     } else {
204       Xs = 0;
205     }
206     IXs = 0;
207   }
208   if (xe+sDist+oDist <= M*dof) {
209     Xe = xe + sDist + oDist;
210     IXe = xe + sDist + oDist;
211   } else {
212     if (bx) {
213       Xe = xe + sDist;
214     } else {
215       Xe = M*dof;
216     }
217     IXe = M*dof;
218   }
219 
220   if (bx == DMDA_BOUNDARY_PERIODIC) {
221     Xs = xs - sDist - oDist;
222     Xe = xe + sDist + oDist;
223     IXs = xs - sDist - oDist;
224     IXe = xe + sDist + oDist;
225   }
226 
227   /* allocate the base parallel and sequential vectors */
228   dd->Nlocal = x;
229   ierr = VecCreateMPIWithArray(comm,dof,dd->Nlocal,PETSC_DECIDE,0,&global);CHKERRQ(ierr);
230   dd->nlocal = (Xe-Xs);
231   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,dof,dd->nlocal,0,&local);CHKERRQ(ierr);
232 
233   /* Create Local to Global Vector Scatter Context */
234   /* local to global inserts non-ghost point region into global */
235   ierr = VecGetOwnershipRange(global,&start,&end);CHKERRQ(ierr);
236   ierr = ISCreateStride(comm,x,start,1,&to);CHKERRQ(ierr);
237   ierr = ISCreateStride(comm,x,xs-Xs,1,&from);CHKERRQ(ierr);
238   ierr = VecScatterCreate(local,from,global,to,&ltog);CHKERRQ(ierr);
239   ierr = PetscLogObjectParent(da,ltog);CHKERRQ(ierr);
240   ierr = ISDestroy(&from);CHKERRQ(ierr);
241   ierr = ISDestroy(&to);CHKERRQ(ierr);
242 
243   /* Create Global to Local Vector Scatter Context */
244   /* global to local must retrieve ghost points */
245   ierr = ISCreateStride(comm,(IXe-IXs),IXs-Xs,1,&to);CHKERRQ(ierr);
246 
247   ierr = PetscMalloc((x+2*(sDist+oDist))*sizeof(PetscInt),&idx);CHKERRQ(ierr);
248   ierr = PetscLogObjectMemory(da,(x+2*(sDist+oDist))*sizeof(PetscInt));CHKERRQ(ierr);
249 
250   for (i=0; i<IXs-Xs; i++) {idx[i] = -1; } /* prepend with -1s if needed for ghosted case*/
251 
252   nn = IXs-Xs;
253   if (bx == DMDA_BOUNDARY_PERIODIC) { /* Handle all cases with wrap first */
254     for (i=0; i<sDist+oDist; i++) {  /* Left ghost points */
255       if ((xs-sDist-oDist+i)>=0) { idx[nn++] = xs-sDist-oDist+i;}
256       else                 { idx[nn++] = M*dof+(xs-sDist-oDist+i);}
257     }
258 
259     for (i=0; i<x; i++) { idx [nn++] = xs + i;}  /* Non-ghost points */
260 
261     for (i=0; i<sDist+oDist; i++) { /* Right ghost points */
262       if ((xe+i)<M*dof) { idx [nn++] =  xe+i; }
263       else              { idx [nn++] = (xe+i) - M*dof;}
264     }
265   } else {      /* Now do all cases with no wrapping */
266     if (0 <= xs-sDist-oDist) {for (i=0; i<sDist+oDist; i++) {idx[nn++] = xs - sDist - oDist + i;}}
267     else               {for (i=0; i<xs;    i++) {idx[nn++] = i;}}
268 
269     for (i=0; i<x; i++) { idx [nn++] = xs + i;}
270 
271     if ((xe+sDist+oDist)<=M*dof) {for (i=0;  i<sDist+oDist;   i++) {idx[nn++]=xe+i;}}
272     else                   {for (i=xe; i<(M*dof); i++) {idx[nn++]=i;}}
273   }
274 
275   ierr = ISCreateGeneral(comm,nn-IXs+Xs,&idx[IXs-Xs],PETSC_COPY_VALUES,&from);CHKERRQ(ierr);
276   ierr = VecScatterCreate(global,from,local,to,&gtol);CHKERRQ(ierr);
277   ierr = PetscLogObjectParent(da,to);CHKERRQ(ierr);
278   ierr = PetscLogObjectParent(da,from);CHKERRQ(ierr);
279   ierr = PetscLogObjectParent(da,gtol);CHKERRQ(ierr);
280   ierr = ISDestroy(&to);CHKERRQ(ierr);
281   ierr = ISDestroy(&from);CHKERRQ(ierr);
282   ierr = VecDestroy(&local);CHKERRQ(ierr);
283   ierr = VecDestroy(&global);CHKERRQ(ierr);
284 
285   dd->xs = xs; dd->xe = xe; dd->ys = 0; dd->ye = 1; dd->zs = 0; dd->ze = 1;
286   dd->Xs = Xs; dd->Xe = Xe; dd->Ys = 0; dd->Ye = 1; dd->Zs = 0; dd->Ze = 1;
287 
288   dd->gtol      = gtol;
289   dd->ltog      = ltog;
290   dd->base      = xs;
291   da->ops->view = DMView_DA_1d;
292 
293   /*
294      Set the local to global ordering in the global vector, this allows use
295      of VecSetValuesLocal().
296   */
297   for (i=0; i<Xe-IXe; i++) {idx[nn++] = -1; } /* pad with -1s if needed for ghosted case*/
298 
299   ierr = ISLocalToGlobalMappingCreate(comm,nn,idx,PETSC_COPY_VALUES,&da->ltogmap);CHKERRQ(ierr);
300   ierr = ISLocalToGlobalMappingBlock(da->ltogmap,dd->w,&da->ltogmapb);CHKERRQ(ierr);
301   ierr = PetscLogObjectParent(da,da->ltogmap);CHKERRQ(ierr);
302 
303   dd->idx = idx;
304   dd->Nl  = nn;
305 
306   PetscFunctionReturn(0);
307 }
308 
309 
310 #undef __FUNCT__
311 #define __FUNCT__ "DMDACreate1d"
312 /*@C
313    DMDACreate1d - Creates an object that will manage the communication of  one-dimensional
314    regular array data that is distributed across some processors.
315 
316    Collective on MPI_Comm
317 
318    Input Parameters:
319 +  comm - MPI communicator
320 .  bx - type of ghost cells at the boundary the array should have, if any. Use
321           DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_GHOSTED, or DMDA_BOUNDARY_PERIODIC.
322 .  M - global dimension of the array (use -M to indicate that it may be set to a different value
323             from the command line with -da_grid_x <M>)
324 .  dof - number of degrees of freedom per node
325 .  s - stencil width
326 -  lx - array containing number of nodes in the X direction on each processor,
327         or PETSC_NULL. If non-null, must be of length as the number of processes in the MPI_Comm.
328 
329    Output Parameter:
330 .  da - the resulting distributed array object
331 
332    Options Database Key:
333 +  -dm_view - Calls DMView() at the conclusion of DMDACreate1d()
334 .  -da_grid_x <nx> - number of grid points in x direction; can set if M < 0
335 .  -da_refine_x <rx> - refinement factor
336 -  -da_refine <n> - refine the DMDA n times before creating it, if M < 0
337 
338    Level: beginner
339 
340    Notes:
341    The array data itself is NOT stored in the DMDA, it is stored in Vec objects;
342    The appropriate vector objects can be obtained with calls to DMCreateGlobalVector()
343    and DMCreateLocalVector() and calls to VecDuplicate() if more are needed.
344 
345 .keywords: distributed array, create, one-dimensional
346 
347 .seealso: DMDestroy(), DMView(), DMDACreate2d(), DMDACreate3d(), DMGlobalToLocalBegin(), DMDASetRefinementFactor(),
348           DMGlobalToLocalEnd(), DMLocalToGlobalBegin(), DMDALocalToLocalBegin(), DMDALocalToLocalEnd(), DMDAGetRefinementFactor(),
349           DMDAGetInfo(), DMCreateGlobalVector(), DMCreateLocalVector(), DMDACreateNaturalVector(), DMLoad(), DMDAGetOwnershipRanges()
350 
351 @*/
352 PetscErrorCode  DMDACreate1d(MPI_Comm comm, DMDABoundaryType bx, PetscInt M, PetscInt dof, PetscInt s, const PetscInt lx[], DM *da)
353 {
354   PetscErrorCode ierr;
355   PetscMPIInt    size;
356 
357   PetscFunctionBegin;
358   ierr = DMDACreate(comm, da);CHKERRQ(ierr);
359   ierr = DMDASetDim(*da, 1);CHKERRQ(ierr);
360   ierr = DMDASetSizes(*da, M, 1, 1);CHKERRQ(ierr);
361   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
362   ierr = DMDASetNumProcs(*da, size, PETSC_DECIDE, PETSC_DECIDE);CHKERRQ(ierr);
363   ierr = DMDASetBoundaryType(*da, bx, DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_NONE);CHKERRQ(ierr);
364   ierr = DMDASetDof(*da, dof);CHKERRQ(ierr);
365   ierr = DMDASetStencilWidth(*da, s);CHKERRQ(ierr);
366   ierr = DMDASetOwnershipRanges(*da, lx, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
367   /* This violates the behavior for other classes, but right now users expect negative dimensions to be handled this way */
368   ierr = DMSetFromOptions(*da);CHKERRQ(ierr);
369   ierr = DMSetUp(*da);CHKERRQ(ierr);
370   ierr = DMViewFromOptions(*da,"-dm_view");CHKERRQ(ierr);
371   PetscFunctionReturn(0);
372 }
373