xref: /petsc/src/dm/impls/da/da1.c (revision 8cbdbec6f8317ddf7886f91eb9c6bd083b543c50)
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,j,*idx,nn,left,xs,xe,x,Xs,Xe,start,end,m,IXs,IXe;
125   PetscErrorCode   ierr;
126 
127   PetscFunctionBegin;
128   ierr = PetscObjectGetComm((PetscObject) da, &comm);CHKERRQ(ierr);
129   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
130   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
131 
132   dd->m = size;
133   m     = dd->m;
134 
135   if (s > 0) {
136     /* if not communicating data then should be ok to have nothing on some processes */
137     if (M < m)     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"More processes than data points! %D %D",m,M);
138     if ((M-1) < s) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Array is too small for stencil! %D %D",M-1,s);
139   }
140 
141   /*
142      Determine locally owned region
143      xs is the first local node number, x is the number of local nodes
144   */
145   if (!lx) {
146     ierr = PetscMalloc(m*sizeof(PetscInt), &dd->lx);CHKERRQ(ierr);
147     ierr = PetscOptionsGetBool(PETSC_NULL,"-da_partition_blockcomm",&flg1,PETSC_NULL);CHKERRQ(ierr);
148     ierr = PetscOptionsGetBool(PETSC_NULL,"-da_partition_nodes_at_end",&flg2,PETSC_NULL);CHKERRQ(ierr);
149     if (flg1) {      /* Block Comm type Distribution */
150       xs = rank*M/m;
151       x  = (rank + 1)*M/m - xs;
152     } else if (flg2) { /* The odd nodes are evenly distributed across last nodes */
153       x = (M + rank)/m;
154       if (M/m == x) { xs = rank*x; }
155       else          { xs = rank*(x-1) + (M+rank)%(x*m); }
156     } else { /* The odd nodes are evenly distributed across the first k nodes */
157       /* Regular PETSc Distribution */
158       x = M/m + ((M % m) > rank);
159       if (rank >= (M % m)) {xs = (rank * (PetscInt)(M/m) + M % m);}
160       else                 {xs = rank * (PetscInt)(M/m) + rank;}
161     }
162     ierr = MPI_Allgather(&xs,1,MPIU_INT,dd->lx,1,MPIU_INT,comm);CHKERRQ(ierr);
163     for (i=0; i<m-1; i++) dd->lx[i] = dd->lx[i+1] - dd->lx[i];
164     dd->lx[m-1] = M - dd->lx[m-1];
165   } else {
166     x  = lx[rank];
167     xs = 0;
168     for (i=0; i<rank; i++) {
169       xs += lx[i];
170     }
171     /* verify that data user provided is consistent */
172     left = xs;
173     for (i=rank; i<size; i++) {
174       left += lx[i];
175     }
176     if (left != M) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Sum of lx across processors not equal to M %D %D",left,M);
177   }
178 
179   /*
180    check if the scatter requires more than one process neighbor or wraps around
181    the domain more than once
182   */
183   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);
184 
185   /* From now on x,xs,xe,Xs,Xe are the exact location in the array */
186   x  *= dof;
187   xs *= dof;
188   xe  = xs + x;
189 
190   /* determine ghost region (Xs) and region scattered into (IXs)  */
191   if (xs-sDist-oDist > 0) {
192     Xs = xs - sDist - oDist;
193     IXs = xs - sDist - oDist;
194   } else {
195     if (bx) {
196       Xs = xs - sDist;
197     } else {
198       Xs = 0;
199     }
200     IXs = 0;
201   }
202   if (xe+sDist+oDist <= M*dof) {
203     Xe = xe + sDist + oDist;
204     IXe = xe + sDist + oDist;
205   } else {
206     if (bx) {
207       Xe = xe + sDist;
208     } else {
209       Xe = M*dof;
210     }
211     IXe = M*dof;
212   }
213 
214   if (bx == DMDA_BOUNDARY_PERIODIC || bx == DMDA_BOUNDARY_MIRROR) {
215     Xs = xs - sDist - oDist;
216     Xe = xe + sDist + oDist;
217     IXs = xs - sDist - oDist;
218     IXe = xe + sDist + oDist;
219   }
220 
221   /* allocate the base parallel and sequential vectors */
222   dd->Nlocal = x;
223   ierr = VecCreateMPIWithArray(comm,dof,dd->Nlocal,PETSC_DECIDE,0,&global);CHKERRQ(ierr);
224   dd->nlocal = (Xe-Xs);
225   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,dof,dd->nlocal,0,&local);CHKERRQ(ierr);
226 
227   /* Create Local to Global Vector Scatter Context */
228   /* local to global inserts non-ghost point region into global */
229   ierr = VecGetOwnershipRange(global,&start,&end);CHKERRQ(ierr);
230   ierr = ISCreateStride(comm,x,start,1,&to);CHKERRQ(ierr);
231   ierr = ISCreateStride(comm,x,xs-Xs,1,&from);CHKERRQ(ierr);
232   ierr = VecScatterCreate(local,from,global,to,&ltog);CHKERRQ(ierr);
233   ierr = PetscLogObjectParent(da,ltog);CHKERRQ(ierr);
234   ierr = ISDestroy(&from);CHKERRQ(ierr);
235   ierr = ISDestroy(&to);CHKERRQ(ierr);
236 
237   /* Create Global to Local Vector Scatter Context */
238   /* global to local must retrieve ghost points */
239   ierr = ISCreateStride(comm,(IXe-IXs),IXs-Xs,1,&to);CHKERRQ(ierr);
240 
241   ierr = PetscMalloc((x+2*(sDist+oDist))*sizeof(PetscInt),&idx);CHKERRQ(ierr);
242   ierr = PetscLogObjectMemory(da,(x+2*(sDist+oDist))*sizeof(PetscInt));CHKERRQ(ierr);
243 
244   for (i=0; i<IXs-Xs; i++) {idx[i] = -1; } /* prepend with -1s if needed for ghosted case*/
245 
246   nn = IXs-Xs;
247   if (bx == DMDA_BOUNDARY_PERIODIC) { /* Handle all cases with periodic first */
248     for (i=0; i<sDist+oDist; i++) {  /* Left ghost points */
249       if ((xs-sDist-oDist+i)>=0) {idx[nn++] = xs-sDist-oDist+i;}
250       else                       {idx[nn++] = M*dof+(xs-sDist-oDist+i);}
251     }
252 
253     for (i=0; i<x; i++) { idx [nn++] = xs + i;}  /* Non-ghost points */
254 
255     for (i=0; i<sDist+oDist; i++) { /* Right ghost points */
256       if ((xe+i)<M*dof) {idx [nn++] =  xe+i; }
257       else              {idx [nn++] = (xe+i) - M*dof;}
258     }
259   } else if (bx == DMDA_BOUNDARY_MIRROR) { /* Handle all cases with periodic first */
260     for (i=0; i<(sDist+oDist)/dof; i++) {  /* Left ghost points */
261       for (j=0; j<dof; j++) {
262         if ((xs-sDist-oDist+i*dof + j)>=0) {idx[nn++] = xs-sDist-oDist+i*dof +j;}
263         else                               {idx[nn++] = sDist+oDist - dof*(i) + j;}
264       }
265     }
266 
267     for (i=0; i<x; i++) { idx [nn++] = xs + i;}  /* Non-ghost points */
268 
269     for (i=0; i<(sDist+oDist)/dof; i++) { /* Right ghost points */
270       for (j=0; j<dof; j++) {
271         if ((xe+i)<M*dof) {idx[nn++] =  xe+i*dof+j; }
272         else              {idx[nn++] = M*dof - dof*(i + 2) + j ;}
273       }
274     }
275   } else {      /* Now do all cases with no periodicity */
276     if (0 <= xs-sDist-oDist) {for (i=0; i<sDist+oDist; i++) {idx[nn++] = xs - sDist - oDist + i;}}
277     else                     {for (i=0; i<xs;    i++) {idx[nn++] = i;}}
278 
279     for (i=0; i<x; i++) { idx [nn++] = xs + i;}
280 
281     if ((xe+sDist+oDist)<=M*dof) {for (i=0;  i<sDist+oDist;   i++) {idx[nn++]=xe+i;}}
282     else                         {for (i=xe; i<(M*dof); i++) {idx[nn++]=i;}}
283   }
284 
285   ierr = ISCreateGeneral(comm,nn-IXs+Xs,&idx[IXs-Xs],PETSC_COPY_VALUES,&from);CHKERRQ(ierr);
286   ierr = VecScatterCreate(global,from,local,to,&gtol);CHKERRQ(ierr);
287   ierr = PetscLogObjectParent(da,to);CHKERRQ(ierr);
288   ierr = PetscLogObjectParent(da,from);CHKERRQ(ierr);
289   ierr = PetscLogObjectParent(da,gtol);CHKERRQ(ierr);
290   ierr = ISDestroy(&to);CHKERRQ(ierr);
291   ierr = ISDestroy(&from);CHKERRQ(ierr);
292   ierr = VecDestroy(&local);CHKERRQ(ierr);
293   ierr = VecDestroy(&global);CHKERRQ(ierr);
294 
295   dd->xs = xs; dd->xe = xe; dd->ys = 0; dd->ye = 1; dd->zs = 0; dd->ze = 1;
296   dd->Xs = Xs; dd->Xe = Xe; dd->Ys = 0; dd->Ye = 1; dd->Zs = 0; dd->Ze = 1;
297 
298   dd->gtol      = gtol;
299   dd->ltog      = ltog;
300   dd->base      = xs;
301   da->ops->view = DMView_DA_1d;
302 
303   /*
304      Set the local to global ordering in the global vector, this allows use
305      of VecSetValuesLocal().
306   */
307   for (i=0; i<Xe-IXe; i++) {idx[nn++] = -1; } /* pad with -1s if needed for ghosted case*/
308 
309   ierr = ISLocalToGlobalMappingCreate(comm,nn,idx,PETSC_COPY_VALUES,&da->ltogmap);CHKERRQ(ierr);
310   ierr = ISLocalToGlobalMappingBlock(da->ltogmap,dd->w,&da->ltogmapb);CHKERRQ(ierr);
311   ierr = PetscLogObjectParent(da,da->ltogmap);CHKERRQ(ierr);
312 
313   dd->idx = idx;
314   dd->Nl  = nn;
315 
316   PetscFunctionReturn(0);
317 }
318 
319 
320 #undef __FUNCT__
321 #define __FUNCT__ "DMDACreate1d"
322 /*@C
323    DMDACreate1d - Creates an object that will manage the communication of  one-dimensional
324    regular array data that is distributed across some processors.
325 
326    Collective on MPI_Comm
327 
328    Input Parameters:
329 +  comm - MPI communicator
330 .  bx - type of ghost cells at the boundary the array should have, if any. Use
331           DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_GHOSTED, or DMDA_BOUNDARY_PERIODIC.
332 .  M - global dimension of the array (use -M to indicate that it may be set to a different value
333             from the command line with -da_grid_x <M>)
334 .  dof - number of degrees of freedom per node
335 .  s - stencil width
336 -  lx - array containing number of nodes in the X direction on each processor,
337         or PETSC_NULL. If non-null, must be of length as the number of processes in the MPI_Comm.
338 
339    Output Parameter:
340 .  da - the resulting distributed array object
341 
342    Options Database Key:
343 +  -dm_view - Calls DMView() at the conclusion of DMDACreate1d()
344 .  -da_grid_x <nx> - number of grid points in x direction; can set if M < 0
345 .  -da_refine_x <rx> - refinement factor
346 -  -da_refine <n> - refine the DMDA n times before creating it, if M < 0
347 
348    Level: beginner
349 
350    Notes:
351    The array data itself is NOT stored in the DMDA, it is stored in Vec objects;
352    The appropriate vector objects can be obtained with calls to DMCreateGlobalVector()
353    and DMCreateLocalVector() and calls to VecDuplicate() if more are needed.
354 
355 .keywords: distributed array, create, one-dimensional
356 
357 .seealso: DMDestroy(), DMView(), DMDACreate2d(), DMDACreate3d(), DMGlobalToLocalBegin(), DMDASetRefinementFactor(),
358           DMGlobalToLocalEnd(), DMLocalToGlobalBegin(), DMDALocalToLocalBegin(), DMDALocalToLocalEnd(), DMDAGetRefinementFactor(),
359           DMDAGetInfo(), DMCreateGlobalVector(), DMCreateLocalVector(), DMDACreateNaturalVector(), DMLoad(), DMDAGetOwnershipRanges()
360 
361 @*/
362 PetscErrorCode  DMDACreate1d(MPI_Comm comm, DMDABoundaryType bx, PetscInt M, PetscInt dof, PetscInt s, const PetscInt lx[], DM *da)
363 {
364   PetscErrorCode ierr;
365   PetscMPIInt    size;
366 
367   PetscFunctionBegin;
368   ierr = DMDACreate(comm, da);CHKERRQ(ierr);
369   ierr = DMDASetDim(*da, 1);CHKERRQ(ierr);
370   ierr = DMDASetSizes(*da, M, 1, 1);CHKERRQ(ierr);
371   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
372   ierr = DMDASetNumProcs(*da, size, PETSC_DECIDE, PETSC_DECIDE);CHKERRQ(ierr);
373   ierr = DMDASetBoundaryType(*da, bx, DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_NONE);CHKERRQ(ierr);
374   ierr = DMDASetDof(*da, dof);CHKERRQ(ierr);
375   ierr = DMDASetStencilWidth(*da, s);CHKERRQ(ierr);
376   ierr = DMDASetOwnershipRanges(*da, lx, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
377   /* This violates the behavior for other classes, but right now users expect negative dimensions to be handled this way */
378   ierr = DMSetFromOptions(*da);CHKERRQ(ierr);
379   ierr = DMSetUp(*da);CHKERRQ(ierr);
380   ierr = DMViewFromOptions(*da,"-dm_view");CHKERRQ(ierr);
381   PetscFunctionReturn(0);
382 }
383