Actual source code: da2.c

  1: /*$Id: da2.c,v 1.180 2001/09/07 20:12:17 bsmith Exp $*/
  2: 
 3:  #include src/dm/da/daimpl.h

  7: int DAGetOwnershipRange(DA da,int **lx,int **ly,int **lz)
  8: {
 11:   if (lx) *lx = da->lx;
 12:   if (ly) *ly = da->ly;
 13:   if (lz) *lz = da->lz;
 14:   return(0);
 15: }

 19: int DAView_2d(DA da,PetscViewer viewer)
 20: {
 21:   int        rank,ierr;
 22:   PetscTruth isascii,isdraw,isbinary;

 25:   MPI_Comm_rank(da->comm,&rank);

 27:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
 28:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);
 29:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);
 30:   if (isascii) {
 31:     PetscViewerASCIISynchronizedPrintf(viewer,"Processor [%d] M %d N %d m %d n %d w %d s %d\n",rank,da->M,
 32:                              da->N,da->m,da->n,da->w,da->s);
 33:     PetscViewerASCIISynchronizedPrintf(viewer,"X range of indices: %d %d, Y range of indices: %d %d\n",da->xs,da->xe,da->ys,da->ye);
 34:     PetscViewerFlush(viewer);
 35:   } else if (isdraw) {
 36:     PetscDraw       draw;
 37:     double     ymin = -1*da->s-1,ymax = da->N+da->s;
 38:     double     xmin = -1*da->s-1,xmax = da->M+da->s;
 39:     double     x,y;
 40:     int        base,*idx;
 41:     char       node[10];
 42:     PetscTruth isnull;
 43: 
 44:     PetscViewerDrawGetDraw(viewer,0,&draw);
 45:     PetscDrawIsNull(draw,&isnull); if (isnull) return(0);
 46:     PetscDrawSetCoordinates(draw,xmin,ymin,xmax,ymax);
 47:     PetscDrawSynchronizedClear(draw);

 49:     /* first processor draw all node lines */
 50:     if (!rank) {
 51:       ymin = 0.0; ymax = da->N - 1;
 52:       for (xmin=0; xmin<da->M; xmin++) {
 53:         PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_BLACK);
 54:       }
 55:       xmin = 0.0; xmax = da->M - 1;
 56:       for (ymin=0; ymin<da->N; ymin++) {
 57:         PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_BLACK);
 58:       }
 59:     }
 60:     PetscDrawSynchronizedFlush(draw);
 61:     PetscDrawPause(draw);

 63:     /* draw my box */
 64:     ymin = da->ys; ymax = da->ye - 1; xmin = da->xs/da->w;
 65:     xmax =(da->xe-1)/da->w;
 66:     PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_RED);
 67:     PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_RED);
 68:     PetscDrawLine(draw,xmin,ymax,xmax,ymax,PETSC_DRAW_RED);
 69:     PetscDrawLine(draw,xmax,ymin,xmax,ymax,PETSC_DRAW_RED);

 71:     /* put in numbers */
 72:     base = (da->base)/da->w;
 73:     for (y=ymin; y<=ymax; y++) {
 74:       for (x=xmin; x<=xmax; x++) {
 75:         sprintf(node,"%d",base++);
 76:         PetscDrawString(draw,x,y,PETSC_DRAW_BLACK,node);
 77:       }
 78:     }

 80:     PetscDrawSynchronizedFlush(draw);
 81:     PetscDrawPause(draw);
 82:     /* overlay ghost numbers, useful for error checking */
 83:     /* put in numbers */

 85:     base = 0; idx = da->idx;
 86:     ymin = da->Ys; ymax = da->Ye; xmin = da->Xs; xmax = da->Xe;
 87:     for (y=ymin; y<ymax; y++) {
 88:       for (x=xmin; x<xmax; x++) {
 89:         if ((base % da->w) == 0) {
 90:           sprintf(node,"%d",idx[base]/da->w);
 91:           PetscDrawString(draw,x/da->w,y,PETSC_DRAW_BLUE,node);
 92:         }
 93:         base++;
 94:       }
 95:     }
 96:     PetscDrawSynchronizedFlush(draw);
 97:     PetscDrawPause(draw);
 98:   } else if (isbinary) {
 99:     DAView_Binary(da,viewer);
100:   } else {
101:     SETERRQ1(1,"Viewer type %s not supported for DA2d",((PetscObject)viewer)->type_name);
102:   }
103:   return(0);
104: }

106: #if defined(PETSC_HAVE_AMS)
107: /*
108:       This function tells the AMS the layout of the vectors, it is called
109:    in the VecPublish_xx routines.
110: */
111: EXTERN_C_BEGIN
114: int AMSSetFieldBlock_DA(AMS_Memory amem,char *name,Vec vec)
115: {
116:   int        ierr,dof,dim,ends[4],shift = 0,starts[] = {0,0,0,0};
117:   DA         da = 0;
118:   PetscTruth isseq,ismpi;

121:   if (((PetscObject)vec)->amem < 0) return(0); /* return if not published */

123:   PetscObjectQuery((PetscObject)vec,"DA",(PetscObject*)&da);
124:   if (!da) return(0);
125:   DAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0);
126:   if (dof > 1) {dim++; shift = 1; ends[0] = dof;}

128:   PetscTypeCompare((PetscObject)vec,VECSEQ,&isseq);
129:   PetscTypeCompare((PetscObject)vec,VECMPI,&ismpi);
130:   if (isseq) {
131:     DAGetGhostCorners(da,0,0,0,ends+shift,ends+shift+1,ends+shift+2);
132:     ends[shift]   += starts[shift]-1;
133:     ends[shift+1] += starts[shift+1]-1;
134:     ends[shift+2] += starts[shift+2]-1;
135:     AMS_Memory_set_field_block(amem,name,dim,starts,ends);
136:     if (ierr) {
137:       char *message;
138:       AMS_Explain_error(ierr,&message);
139:       SETERRQ(ierr,message);
140:     }
141:   } else if (ismpi) {
142:     DAGetCorners(da,starts+shift,starts+shift+1,starts+shift+2,
143:                            ends+shift,ends+shift+1,ends+shift+2);
144:     ends[shift]   += starts[shift]-1;
145:     ends[shift+1] += starts[shift+1]-1;
146:     ends[shift+2] += starts[shift+2]-1;
147:     AMS_Memory_set_field_block(amem,name,dim,starts,ends);
148:     if (ierr) {
149:       char *message;
150:       AMS_Explain_error(ierr,&message);
151:       SETERRQ(ierr,message);
152:     }
153:   } else {
154:     SETERRQ1(1,"Wrong vector type %s for this call",((PetscObject)vec)->type_name);
155:   }

157:   return(0);
158: }
159: EXTERN_C_END
160: #endif

164: int DAPublish_Petsc(PetscObject obj)
165: {
166: #if defined(PETSC_HAVE_AMS)
167:   DA          v = (DA) obj;
168:   int         ierr;
169: #endif


173: #if defined(PETSC_HAVE_AMS)
174:   /* if it is already published then return */
175:   if (v->amem >=0) return(0);

177:   PetscObjectPublishBaseBegin(obj);
178:   PetscObjectPublishBaseEnd(obj);
179: #endif

181:   return(0);
182: }


187: /*@C
188:    DACreate2d -  Creates an object that will manage the communication of  two-dimensional 
189:    regular array data that is distributed across some processors.

191:    Collective on MPI_Comm

193:    Input Parameters:
194: +  comm - MPI communicator
195: .  wrap - type of periodicity should the array have. 
196:          Use one of DA_NONPERIODIC, DA_XPERIODIC, DA_YPERIODIC, or DA_XYPERIODIC.
197: .  stencil_type - stencil type.  Use either DA_STENCIL_BOX or DA_STENCIL_STAR.
198: .  M,N - global dimension in each direction of the array
199: .  m,n - corresponding number of processors in each dimension 
200:          (or PETSC_DECIDE to have calculated)
201: .  dof - number of degrees of freedom per node
202: .  s - stencil width
203: -  lx, ly - arrays containing the number of nodes in each cell along
204:            the x and y coordinates, or PETSC_NULL. If non-null, these
205:            must be of length as m and n, and the corresponding
206:            m and n cannot be PETSC_DECIDE. The sum of the lx[] entries
207:            must be M, and the sum of the ly[] entries must be N.

209:    Output Parameter:
210: .  inra - the resulting distributed array object

212:    Options Database Key:
213: +  -da_view - Calls DAView() at the conclusion of DACreate2d()
214: .  -da_grid_x <nx> - number of grid points in x direction, if M < 0
215: .  -da_grid_y <ny> - number of grid points in y direction, if N < 0
216: .  -da_processors_x <nx> - number of processors in x direction
217: -  -da_processors_y <ny> - number of processors in y direction

219:    Level: beginner

221:    Notes:
222:    The stencil type DA_STENCIL_STAR with width 1 corresponds to the 
223:    standard 5-pt stencil, while DA_STENCIL_BOX with width 1 denotes
224:    the standard 9-pt stencil.

226:    The array data itself is NOT stored in the DA, it is stored in Vec objects;
227:    The appropriate vector objects can be obtained with calls to DACreateGlobalVector()
228:    and DACreateLocalVector() and calls to VecDuplicate() if more are needed.

230: .keywords: distributed array, create, two-dimensional

232: .seealso: DADestroy(), DAView(), DACreate1d(), DACreate3d(), DAGlobalToLocalBegin(),
233:           DAGlobalToLocalEnd(), DALocalToGlobal(), DALocalToLocalBegin(), DALocalToLocalEnd(),
234:           DAGetInfo(), DACreateGlobalVector(), DACreateLocalVector(), DACreateNaturalVector(), DALoad(), DAView()

236: @*/
237: int DACreate2d(MPI_Comm comm,DAPeriodicType wrap,DAStencilType stencil_type,
238:                 int M,int N,int m,int n,int dof,int s,int *lx,int *ly,DA *inra)
239: {
240:   int           rank,size,xs,xe,ys,ye,x,y,Xs,Xe,Ys,Ye,ierr,start,end;
241:   int           up,down,left,i,n0,n1,n2,n3,n5,n6,n7,n8,*idx,nn;
242:   int           xbase,*bases,*ldims,j,x_t,y_t,s_t,base,count;
243:   int           s_x,s_y; /* s proportionalized to w */
244:   int           *flx = 0,*fly = 0;
245:   int           sn0 = 0,sn2 = 0,sn6 = 0,sn8 = 0,refine_x = 2, refine_y = 2,tM = M,tN = N;
246:   PetscTruth    flg1,flg2;
247:   DA            da;
248:   Vec           local,global;
249:   VecScatter    ltog,gtol;
250:   IS            to,from;

254:   *inra = 0;
255: #ifndef PETSC_USE_DYNAMIC_LIBRARIES
256:   DMInitializePackage(PETSC_NULL);
257: #endif

259:   if (dof < 1) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Must have 1 or more degrees of freedom per node: %d",dof);
260:   if (s < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Stencil width cannot be negative: %d",s);

262:   PetscOptionsBegin(comm,PETSC_NULL,"2d DA Options","DA");
263:     if (M < 0){
264:       tM = -M;
265:       PetscOptionsInt("-da_grid_x","Number of grid points in x direction","DACreate2d",tM,&tM,PETSC_NULL);
266:     }
267:     if (N < 0){
268:       tN = -N;
269:       PetscOptionsInt("-da_grid_y","Number of grid points in y direction","DACreate2d",tN,&tN,PETSC_NULL);
270:     }
271:     PetscOptionsInt("-da_processors_x","Number of processors in x direction","DACreate2d",m,&m,PETSC_NULL);
272:     PetscOptionsInt("-da_processors_y","Number of processors in y direction","DACreate2d",n,&n,PETSC_NULL);
273:     PetscOptionsInt("-da_refine_x","Refinement ratio in x direction","DACreate2d",refine_x,&refine_x,PETSC_NULL);
274:     PetscOptionsInt("-da_refine_y","Refinement ratio in y direction","DACreate2d",refine_y,&refine_y,PETSC_NULL);
275:   PetscOptionsEnd();
276:   M = tM; N = tN;

278:   PetscHeaderCreate(da,_p_DA,struct _DAOps,DA_COOKIE,0,"DA",comm,DADestroy,DAView);
279:   PetscLogObjectCreate(da);
280:   da->bops->publish           = DAPublish_Petsc;
281:   da->ops->createglobalvector = DACreateGlobalVector;
282:   da->ops->getinterpolation   = DAGetInterpolation;
283:   da->ops->getcoloring        = DAGetColoring;
284:   da->ops->getmatrix          = DAGetMatrix;
285:   da->ops->refine             = DARefine;
286:   da->ops->getinjection       = DAGetInjection;
287:   PetscLogObjectMemory(da,sizeof(struct _p_DA));
288:   da->dim        = 2;
289:   da->interptype = DA_Q1;
290:   da->refine_x   = refine_x;
291:   da->refine_y   = refine_y;
292:   PetscMalloc(dof*sizeof(char*),&da->fieldname);
293:   PetscMemzero(da->fieldname,dof*sizeof(char*));

295:   MPI_Comm_size(comm,&size);
296:   MPI_Comm_rank(comm,&rank);

298:   if (m != PETSC_DECIDE) {
299:     if (m < 1) {SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in X direction: %d",m);}
300:     else if (m > size) {SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in X direction: %d %d",m,size);}
301:   }
302:   if (n != PETSC_DECIDE) {
303:     if (n < 1) {SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Y direction: %d",n);}
304:     else if (n > size) {SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Y direction: %d %d",n,size);}
305:   }

307:   if (m == PETSC_DECIDE || n == PETSC_DECIDE) {
308:     /* try for squarish distribution */
309:     /* This should use MPI_Dims_create instead */
310:     m = (int)(0.5 + sqrt(((double)M)*((double)size)/((double)N)));
311:     if (!m) m = 1;
312:     while (m > 0) {
313:       n = size/m;
314:       if (m*n == size) break;
315:       m--;
316:     }
317:     if (M > N && m < n) {int _m = m; m = n; n = _m;}
318:     if (m*n != size) SETERRQ(PETSC_ERR_PLIB,"Internally Created Bad Partition");
319:   } else if (m*n != size) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Given Bad partition");

321:   if (M < m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Partition in x direction is too fine! %d %d",M,m);
322:   if (N < n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Partition in y direction is too fine! %d %d",N,n);

324:   /*
325:      We should create an MPI Cartesian topology here, with reorder
326:      set to true.  That would create a NEW communicator that we would
327:      need to use for operations on this distributed array 
328:   */
329:   PetscOptionsHasName(PETSC_NULL,"-da_partition_nodes_at_end",&flg2);

331:   /* 
332:      Determine locally owned region 
333:      xs is the first local node number, x is the number of local nodes 
334:   */
335:   if (lx) { /* user sets distribution */
336:     x  = lx[rank % m];
337:     xs = 0;
338:     for (i=0; i<(rank % m); i++) {
339:       xs += lx[i];
340:     }
341:     left = xs;
342:     for (i=(rank % m); i<m; i++) {
343:       left += lx[i];
344:     }
345:     if (left != M) {
346:       SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Sum of lx across processors not equal to M: %d %d",left,M);
347:     }
348:   } else if (flg2) {
349:     x = (M + rank%m)/m;
350:     if (m > 1 && x < s) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column width is too thin for stencil! %d %d",x,s);
351:     if (M/m == x) { xs = (rank % m)*x; }
352:     else          { xs = (rank % m)*(x-1) + (M+(rank % m))%(x*m); }
353:     SETERRQ(PETSC_ERR_SUP,"-da_partition_nodes_at_end not supported");
354:   } else { /* Normal PETSc distribution */
355:     x = M/m + ((M % m) > (rank % m));
356:     if (m > 1 && x < s) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column width is too thin for stencil! %d %d",x,s);
357:     if ((M % m) > (rank % m)) { xs = (rank % m)*x; }
358:     else                      { xs = (M % m)*(x+1) + ((rank % m)-(M % m))*x; }
359:     PetscMalloc(m*sizeof(int),&lx);
360:     flx = lx;
361:     for (i=0; i<m; i++) {
362:       lx[i] = M/m + ((M % m) > i);
363:     }
364:   }

366:   /* 
367:      Determine locally owned region 
368:      ys is the first local node number, y is the number of local nodes 
369:   */
370:   if (ly) { /* user sets distribution */
371:     y  = ly[rank/m];
372:     ys = 0;
373:     for (i=0; i<(rank/m); i++) {
374:       ys += ly[i];
375:     }
376:     left = ys;
377:     for (i=(rank/m); i<n; i++) {
378:       left += ly[i];
379:     }
380:     if (left != N) {
381:       SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Sum of ly across processors not equal to N: %d %d",left,N);
382:     }
383:   } else if (flg2) {
384:     y = (N + rank/m)/n;
385:     if (n > 1 && y < s) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row width is too thin for stencil! %d %d",y,s);
386:     if (N/n == y) { ys = (rank/m)*y;  }
387:     else          { ys = (rank/m)*(y-1) + (N+(rank/m))%(y*n); }
388:     SETERRQ(PETSC_ERR_SUP,"-da_partition_nodes_at_end not supported");
389:   } else { /* Normal PETSc distribution */
390:     y = N/n + ((N % n) > (rank/m));
391:     if (n > 1 && y < s) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row width is too thin for stencil! %d %d",y,s);
392:     if ((N % n) > (rank/m)) { ys = (rank/m)*y; }
393:     else                    { ys = (N % n)*(y+1) + ((rank/m)-(N % n))*y; }
394:     PetscMalloc(n*sizeof(int),&ly);
395:     fly  = ly;
396:     for (i=0; i<n; i++) {
397:       ly[i] = N/n + ((N % n) > i);
398:     }
399:   }

401:   xe = xs + x;
402:   ye = ys + y;

404:   /* determine ghost region */
405:   /* Assume No Periodicity */
406:   if (xs-s > 0) Xs = xs - s; else Xs = 0;
407:   if (ys-s > 0) Ys = ys - s; else Ys = 0;
408:   if (xe+s <= M) Xe = xe + s; else Xe = M;
409:   if (ye+s <= N) Ye = ye + s; else Ye = N;

411:   /* X Periodic */
412:   if (DAXPeriodic(wrap)){
413:     Xs = xs - s;
414:     Xe = xe + s;
415:   }

417:   /* Y Periodic */
418:   if (DAYPeriodic(wrap)){
419:     Ys = ys - s;
420:     Ye = ye + s;
421:   }

423:   /* Resize all X parameters to reflect w */
424:   x   *= dof;
425:   xs  *= dof;
426:   xe  *= dof;
427:   Xs  *= dof;
428:   Xe  *= dof;
429:   s_x = s*dof;
430:   s_y = s;

432:   /* determine starting point of each processor */
433:   nn = x*y;
434:   PetscMalloc((2*size+1)*sizeof(int),&bases);
435:   ldims = (int*)(bases+size+1);
436:   MPI_Allgather(&nn,1,MPI_INT,ldims,1,MPI_INT,comm);
437:   bases[0] = 0;
438:   for (i=1; i<=size; i++) {
439:     bases[i] = ldims[i-1];
440:   }
441:   for (i=1; i<=size; i++) {
442:     bases[i] += bases[i-1];
443:   }

445:   /* allocate the base parallel and sequential vectors */
446:   da->Nlocal = x*y;
447:   VecCreateMPIWithArray(comm,da->Nlocal,PETSC_DECIDE,0,&global);
448:   VecSetBlockSize(global,dof);
449:   da->nlocal = (Xe-Xs)*(Ye-Ys) ;
450:   VecCreateSeqWithArray(PETSC_COMM_SELF,da->nlocal,0,&local);
451:   VecSetBlockSize(local,dof);


454:   /* generate appropriate vector scatters */
455:   /* local to global inserts non-ghost point region into global */
456:   VecGetOwnershipRange(global,&start,&end);
457:   ISCreateStride(comm,x*y,start,1,&to);

459:   left  = xs - Xs; down  = ys - Ys; up    = down + y;
460:   PetscMalloc(x*(up - down)*sizeof(int),&idx);
461:   count = 0;
462:   for (i=down; i<up; i++) {
463:     for (j=0; j<x; j++) {
464:       idx[count++] = left + i*(Xe-Xs) + j;
465:     }
466:   }
467:   ISCreateGeneral(comm,count,idx,&from);
468:   PetscFree(idx);

470:   VecScatterCreate(local,from,global,to,&ltog);
471:   PetscLogObjectParent(da,to);
472:   PetscLogObjectParent(da,from);
473:   PetscLogObjectParent(da,ltog);
474:   ISDestroy(from);
475:   ISDestroy(to);

477:   /* global to local must include ghost points */
478:   if (stencil_type == DA_STENCIL_BOX) {
479:     ISCreateStride(comm,(Xe-Xs)*(Ye-Ys),0,1,&to);
480:   } else {
481:     /* must drop into cross shape region */
482:     /*       ---------|
483:             |  top    |
484:          |---         ---|
485:          |   middle      |
486:          |               |
487:          ----         ----
488:             | bottom  |
489:             -----------
490:         Xs xs        xe  Xe */
491:     /* bottom */
492:     left  = xs - Xs; down = ys - Ys; up    = down + y;
493:     count = down*(xe-xs) + (up-down)*(Xe-Xs) + (Ye-Ys-up)*(xe-xs);
494:     PetscMalloc(count*sizeof(int),&idx);
495:     count = 0;
496:     for (i=0; i<down; i++) {
497:       for (j=0; j<xe-xs; j++) {
498:         idx[count++] = left + i*(Xe-Xs) + j;
499:       }
500:     }
501:     /* middle */
502:     for (i=down; i<up; i++) {
503:       for (j=0; j<Xe-Xs; j++) {
504:         idx[count++] = i*(Xe-Xs) + j;
505:       }
506:     }
507:     /* top */
508:     for (i=up; i<Ye-Ys; i++) {
509:       for (j=0; j<xe-xs; j++) {
510:         idx[count++] = left + i*(Xe-Xs) + j;
511:       }
512:     }
513:     ISCreateGeneral(comm,count,idx,&to);
514:     PetscFree(idx);
515:   }


518:   /* determine who lies on each side of us stored in    n6 n7 n8
519:                                                         n3    n5
520:                                                         n0 n1 n2
521:   */

523:   /* Assume the Non-Periodic Case */
524:   n1 = rank - m;
525:   if (rank % m) {
526:     n0 = n1 - 1;
527:   } else {
528:     n0 = -1;
529:   }
530:   if ((rank+1) % m) {
531:     n2 = n1 + 1;
532:     n5 = rank + 1;
533:     n8 = rank + m + 1; if (n8 >= m*n) n8 = -1;
534:   } else {
535:     n2 = -1; n5 = -1; n8 = -1;
536:   }
537:   if (rank % m) {
538:     n3 = rank - 1;
539:     n6 = n3 + m; if (n6 >= m*n) n6 = -1;
540:   } else {
541:     n3 = -1; n6 = -1;
542:   }
543:   n7 = rank + m; if (n7 >= m*n) n7 = -1;


546:   /* Modify for Periodic Cases */
547:   if (wrap == DA_YPERIODIC) {  /* Handle Top and Bottom Sides */
548:     if (n1 < 0) n1 = rank + m * (n-1);
549:     if (n7 < 0) n7 = rank - m * (n-1);
550:     if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1;
551:     if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1;
552:     if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1;
553:     if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1;
554:   } else if (wrap == DA_XPERIODIC) { /* Handle Left and Right Sides */
555:     if (n3 < 0) n3 = rank + (m-1);
556:     if (n5 < 0) n5 = rank - (m-1);
557:     if ((n1 >= 0) && (n0 < 0)) n0 = rank-1;
558:     if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1;
559:     if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1;
560:     if ((n7 >= 0) && (n8 < 0)) n8 = rank+1;
561:   } else if (wrap == DA_XYPERIODIC) {

563:     /* Handle all four corners */
564:     if ((n6 < 0) && (n7 < 0) && (n3 < 0)) n6 = m-1;
565:     if ((n8 < 0) && (n7 < 0) && (n5 < 0)) n8 = 0;
566:     if ((n2 < 0) && (n5 < 0) && (n1 < 0)) n2 = size-m;
567:     if ((n0 < 0) && (n3 < 0) && (n1 < 0)) n0 = size-1;

569:     /* Handle Top and Bottom Sides */
570:     if (n1 < 0) n1 = rank + m * (n-1);
571:     if (n7 < 0) n7 = rank - m * (n-1);
572:     if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1;
573:     if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1;
574:     if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1;
575:     if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1;

577:     /* Handle Left and Right Sides */
578:     if (n3 < 0) n3 = rank + (m-1);
579:     if (n5 < 0) n5 = rank - (m-1);
580:     if ((n1 >= 0) && (n0 < 0)) n0 = rank-1;
581:     if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1;
582:     if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1;
583:     if ((n7 >= 0) && (n8 < 0)) n8 = rank+1;
584:   }

586:   if (stencil_type == DA_STENCIL_STAR) {
587:     /* save corner processor numbers */
588:     sn0 = n0; sn2 = n2; sn6 = n6; sn8 = n8;
589:     n0 = n2 = n6 = n8 = -1;
590:   }

592:   PetscMalloc((x+2*s_x)*(y+2*s_y)*sizeof(int),&idx);
593:   PetscLogObjectMemory(da,(x+2*s_x)*(y+2*s_y)*sizeof(int));
594:   nn = 0;

596:   xbase = bases[rank];
597:   for (i=1; i<=s_y; i++) {
598:     if (n0 >= 0) { /* left below */
599:       x_t = lx[n0 % m]*dof;
600:       y_t = ly[(n0/m)];
601:       s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x;
602:       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
603:     }
604:     if (n1 >= 0) { /* directly below */
605:       x_t = x;
606:       y_t = ly[(n1/m)];
607:       s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t;
608:       for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
609:     }
610:     if (n2 >= 0) { /* right below */
611:       x_t = lx[n2 % m]*dof;
612:       y_t = ly[(n2/m)];
613:       s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t;
614:       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
615:     }
616:   }

618:   for (i=0; i<y; i++) {
619:     if (n3 >= 0) { /* directly left */
620:       x_t = lx[n3 % m]*dof;
621:       /* y_t = y; */
622:       s_t = bases[n3] + (i+1)*x_t - s_x;
623:       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
624:     }

626:     for (j=0; j<x; j++) { idx[nn++] = xbase++; } /* interior */

628:     if (n5 >= 0) { /* directly right */
629:       x_t = lx[n5 % m]*dof;
630:       /* y_t = y; */
631:       s_t = bases[n5] + (i)*x_t;
632:       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
633:     }
634:   }

636:   for (i=1; i<=s_y; i++) {
637:     if (n6 >= 0) { /* left above */
638:       x_t = lx[n6 % m]*dof;
639:       /* y_t = ly[(n6/m)]; */
640:       s_t = bases[n6] + (i)*x_t - s_x;
641:       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
642:     }
643:     if (n7 >= 0) { /* directly above */
644:       x_t = x;
645:       /* y_t = ly[(n7/m)]; */
646:       s_t = bases[n7] + (i-1)*x_t;
647:       for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
648:     }
649:     if (n8 >= 0) { /* right above */
650:       x_t = lx[n8 % m]*dof;
651:       /* y_t = ly[(n8/m)]; */
652:       s_t = bases[n8] + (i-1)*x_t;
653:       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
654:     }
655:   }

657:   base = bases[rank];
658:   ISCreateGeneral(comm,nn,idx,&from);
659:   VecScatterCreate(global,from,local,to,&gtol);
660:   PetscLogObjectParent(da,to);
661:   PetscLogObjectParent(da,from);
662:   PetscLogObjectParent(da,gtol);
663:   ISDestroy(to);
664:   ISDestroy(from);

666:   if (stencil_type == DA_STENCIL_STAR) {
667:     /*
668:         Recompute the local to global mappings, this time keeping the 
669:       information about the cross corner processor numbers.
670:     */
671:     n0 = sn0; n2 = sn2; n6 = sn6; n8 = sn8;
672:     nn = 0;
673:     xbase = bases[rank];
674:     for (i=1; i<=s_y; i++) {
675:       if (n0 >= 0) { /* left below */
676:         x_t = lx[n0 % m]*dof;
677:         y_t = ly[(n0/m)];
678:         s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x;
679:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
680:       }
681:       if (n1 >= 0) { /* directly below */
682:         x_t = x;
683:         y_t = ly[(n1/m)];
684:         s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t;
685:         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
686:       }
687:       if (n2 >= 0) { /* right below */
688:         x_t = lx[n2 % m]*dof;
689:         y_t = ly[(n2/m)];
690:         s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t;
691:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
692:       }
693:     }

695:     for (i=0; i<y; i++) {
696:       if (n3 >= 0) { /* directly left */
697:         x_t = lx[n3 % m]*dof;
698:         /* y_t = y; */
699:         s_t = bases[n3] + (i+1)*x_t - s_x;
700:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
701:       }

703:       for (j=0; j<x; j++) { idx[nn++] = xbase++; } /* interior */

705:       if (n5 >= 0) { /* directly right */
706:         x_t = lx[n5 % m]*dof;
707:         /* y_t = y; */
708:         s_t = bases[n5] + (i)*x_t;
709:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
710:       }
711:     }

713:     for (i=1; i<=s_y; i++) {
714:       if (n6 >= 0) { /* left above */
715:         x_t = lx[n6 % m]*dof;
716:         /* y_t = ly[(n6/m)]; */
717:         s_t = bases[n6] + (i)*x_t - s_x;
718:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
719:       }
720:       if (n7 >= 0) { /* directly above */
721:         x_t = x;
722:         /* y_t = ly[(n7/m)]; */
723:         s_t = bases[n7] + (i-1)*x_t;
724:         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
725:       }
726:       if (n8 >= 0) { /* right above */
727:         x_t = lx[n8 % m]*dof;
728:         /* y_t = ly[(n8/m)]; */
729:         s_t = bases[n8] + (i-1)*x_t;
730:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
731:       }
732:     }
733:   }
734:   PetscFree(bases);

736:   da->M  = M;  da->N  = N;  da->m  = m;  da->n  = n;  da->w = dof;  da->s = s;
737:   da->xs = xs; da->xe = xe; da->ys = ys; da->ye = ye; da->zs = 0; da->ze = 1;
738:   da->Xs = Xs; da->Xe = Xe; da->Ys = Ys; da->Ye = Ye; da->Zs = 0; da->Ze = 1;
739:   da->P  = 1;  da->p  = 1;

741:   VecDestroy(local);
742:   VecDestroy(global);

744:   da->gtol         = gtol;
745:   da->ltog         = ltog;
746:   da->idx          = idx;
747:   da->Nl           = nn;
748:   da->base         = base;
749:   da->wrap         = wrap;
750:   da->ops->view    = DAView_2d;
751:   da->stencil_type = stencil_type;

753:   /* 
754:      Set the local to global ordering in the global vector, this allows use
755:      of VecSetValuesLocal().
756:   */
757:   ISLocalToGlobalMappingCreateNC(comm,nn,idx,&da->ltogmap);
758:   ISLocalToGlobalMappingBlock(da->ltogmap,da->w,&da->ltogmapb);
759:   PetscLogObjectParent(da,da->ltogmap);

761:   *inra = da;

763:   da->ltol = PETSC_NULL;
764:   da->ao   = PETSC_NULL;


767:   if (!flx) {
768:     PetscMalloc(m*sizeof(int),&flx);
769:     PetscMemcpy(flx,lx,m*sizeof(int));
770:   }
771:   if (!fly) {
772:     PetscMalloc(n*sizeof(int),&fly);
773:     PetscMemcpy(fly,ly,n*sizeof(int));
774:   }
775:   da->lx = flx;
776:   da->ly = fly;

778:   PetscOptionsHasName(PETSC_NULL,"-da_view",&flg1);
779:   if (flg1) {DAView(da,PETSC_VIEWER_STDOUT_(da->comm));}
780:   PetscOptionsHasName(PETSC_NULL,"-da_view_draw",&flg1);
781:   if (flg1) {DAView(da,PETSC_VIEWER_DRAW_(da->comm));}
782:   PetscOptionsHasName(PETSC_NULL,"-help",&flg1);
783:   if (flg1) {DAPrintHelp(da);}

785:   PetscPublishAll(da);
786:   return(0);
787: }

791: /*@
792:    DAPrintHelp - Prints command line options for DA.

794:    Collective on DA

796:    Input Parameters:
797: .  da - the distributed array

799:    Level: intermediate

801: .seealso: DACreate1d(), DACreate2d(), DACreate3d()

803: .keywords: DA, help

805: @*/
806: int DAPrintHelp(DA da)
807: {
808:   static PetscTruth called = PETSC_FALSE;
809:   MPI_Comm          comm;
810:   int               ierr;


815:   comm = da->comm;
816:   if (!called) {
817:     (*PetscHelpPrintf)(comm,"General Distributed Array (DA) options:\n");
818:     (*PetscHelpPrintf)(comm,"  -da_view: print DA distribution to screen\n");
819:     (*PetscHelpPrintf)(comm,"  -da_view_draw: display DA in window\n");
820:     called = PETSC_TRUE;
821:   }
822:   return(0);
823: }

827: /*@C
828:    DARefine - Creates a new distributed array that is a refinement of a given
829:    distributed array.

831:    Collective on DA

833:    Input Parameter:
834: +  da - initial distributed array
835: -  comm - communicator to contain refined DA, must be either same as the da communicator or include the 
836:           da communicator and be 2, 4, or 8 times larger. Currently ignored

838:    Output Parameter:
839: .  daref - refined distributed array

841:    Level: advanced

843:    Note:
844:    Currently, refinement consists of just doubling the number of grid spaces
845:    in each dimension of the DA.

847: .keywords:  distributed array, refine

849: .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy()
850: @*/
851: int DARefine(DA da,MPI_Comm comm,DA *daref)
852: {
853:   int M,N,P,ierr;
854:   DA  da2;


859:   if (DAXPeriodic(da->wrap) || da->interptype == DA_Q0){
860:     M = da->refine_x*da->M;
861:   } else {
862:     M = 1 + da->refine_x*(da->M - 1);
863:   }
864:   if (DAYPeriodic(da->wrap) || da->interptype == DA_Q0){
865:     N = da->refine_y*da->N;
866:   } else {
867:     N = 1 + da->refine_y*(da->N - 1);
868:   }
869:   if (DAZPeriodic(da->wrap) || da->interptype == DA_Q0){
870:     P = da->refine_z*da->P;
871:   } else {
872:     P = 1 + da->refine_z*(da->P - 1);
873:   }
874:   if (da->dim == 1) {
875:     DACreate1d(da->comm,da->wrap,M,da->w,da->s,PETSC_NULL,&da2);
876:   } else if (da->dim == 2) {
877:     DACreate2d(da->comm,da->wrap,da->stencil_type,M,N,da->m,da->n,da->w,da->s,PETSC_NULL,PETSC_NULL,&da2);
878:   } else if (da->dim == 3) {
879:     DACreate3d(da->comm,da->wrap,da->stencil_type,M,N,P,da->m,da->n,da->p,da->w,da->s,0,0,0,&da2);
880:   }
881:   /* allow overloaded (user replaced) operations to be inherited by refinement clones */
882:   da2->ops->getmatrix        = da->ops->getmatrix;
883:   da2->ops->getinterpolation = da->ops->getinterpolation;
884:   da2->ops->getcoloring      = da->ops->getcoloring;
885: 
886:   /* copy fill information if given */
887:   if (da->dfill) {
888:     PetscMalloc((da->dfill[da->w]+da->w+1)*sizeof(int),&da2->dfill);
889:     PetscMemcpy(da2->dfill,da->dfill,(da->dfill[da->w]+da->w+1)*sizeof(int));
890:   }
891:   if (da->ofill) {
892:     PetscMalloc((da->ofill[da->w]+da->w+1)*sizeof(int),&da2->ofill);
893:     PetscMemcpy(da2->ofill,da->ofill,(da->ofill[da->w]+da->w+1)*sizeof(int));
894:   }
895:   *daref = da2;
896:   return(0);
897: }

899: /*@C
900:      DASetGetMatrix - Sets the routine used by the DA to allocate a matrix.

902:     Collective on DA

904:   Input Parameters:
905: +    da - the DA object
906: -    f - the function that allocates the matrix for that specific DA

908:   Level: developer

910:    Notes: See DASetBlockFills() that provides a simple way to provide the nonzero structure for 
911:        the diagonal and off-diagonal blocks of the matrix

913: .seealso: DAGetMatrix(), DASetBlockFills()
914: @*/
915: int DASetGetMatrix(DA da,int (*f)(DA,MatType,Mat*))
916: {
918:   da->ops->getmatrix = f;
919:   return(0);
920: }

922: /*
923:       M is number of grid points 
924:       m is number of processors

926: */
929: int DASplitComm2d(MPI_Comm comm,int M,int N,int sw,MPI_Comm *outcomm)
930: {
931:   int ierr,m,n = 0,csize,size,rank,x = 0,y = 0;

934:   MPI_Comm_size(comm,&size);
935:   MPI_Comm_rank(comm,&rank);

937:   csize = 4*size;
938:   do {
939:     if (csize % 4) SETERRQ4(1,"Cannot split communicator of size %d tried %d %d %d",size,csize,x,y);
940:     csize   = csize/4;
941: 
942:     m = (int)(0.5 + sqrt(((double)M)*((double)csize)/((double)N)));
943:     if (!m) m = 1;
944:     while (m > 0) {
945:       n = csize/m;
946:       if (m*n == csize) break;
947:       m--;
948:     }
949:     if (M > N && m < n) {int _m = m; m = n; n = _m;}

951:     x = M/m + ((M % m) > ((csize-1) % m));
952:     y = (N + (csize-1)/m)/n;
953:   } while ((x < 4 || y < 4) && csize > 1);
954:   if (size != csize) {
955:     MPI_Group entire_group,sub_group;
956:     int       i,*groupies;

958:     MPI_Comm_group(comm,&entire_group);
959:     PetscMalloc(csize*sizeof(int),&groupies);
960:     for (i=0; i<csize; i++) {
961:       groupies[i] = (rank/csize)*csize + i;
962:     }
963:     MPI_Group_incl(entire_group,csize,groupies,&sub_group);
964:     PetscFree(groupies);
965:     MPI_Comm_create(comm,sub_group,outcomm);
966:     MPI_Group_free(&entire_group);
967:     MPI_Group_free(&sub_group);
968:     PetscLogInfo(0,"Creating redundant coarse problems of size %d\n",csize);
969:   } else {
970:     *outcomm = comm;
971:   }
972:   return(0);
973: }

977: /*@C
978:        DASetLocalFunction - Caches in a DA a local function. 

980:    Collective on DA

982:    Input Parameter:
983: +  da - initial distributed array
984: -  lf - the local function

986:    Level: intermediate

988:    Notes: The routine SNESDAFormFunction() uses this the cached function to evaluate the user provided function.

990: .keywords:  distributed array, refine

992: .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DAGetLocalFunction(), DASetLocalFunctioni()
993: @*/
994: int DASetLocalFunction(DA da,DALocalFunction1 lf)
995: {
998:   da->lf    = lf;
999:   return(0);
1000: }

1004: /*@C
1005:        DASetLocalFunctioni - Caches in a DA a local function that evaluates a single component

1007:    Collective on DA

1009:    Input Parameter:
1010: +  da - initial distributed array
1011: -  lfi - the local function

1013:    Level: intermediate

1015: .keywords:  distributed array, refine

1017: .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DAGetLocalFunction(), DASetLocalFunction()
1018: @*/
1019: int DASetLocalFunctioni(DA da,int (*lfi)(DALocalInfo*,MatStencil*,void*,PetscScalar*,void*))
1020: {
1023:   da->lfi = lfi;
1024:   return(0);
1025: }


1030: int DASetLocalAdicFunction_Private(DA da,DALocalFunction1 ad_lf)
1031: {
1034:   da->adic_lf = ad_lf;
1035:   return(0);
1036: }

1038: /*MC
1039:        DASetLocalAdicFunctioni - Caches in a DA a local functioni computed by ADIC/ADIFOR

1041:    Collective on DA

1043:    Synopsis:
1044:    int int DASetLocalAdicFunctioni(DA da,int (ad_lf*)(DALocalInfo*,MatStencil*,void*,void*,void*)
1045:    
1046:    Input Parameter:
1047: +  da - initial distributed array
1048: -  ad_lfi - the local function as computed by ADIC/ADIFOR

1050:    Level: intermediate

1052: .keywords:  distributed array, refine

1054: .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DAGetLocalFunction(), DASetLocalFunction(),
1055:           DASetLocalJacobian(), DASetLocalFunctioni()
1056: M*/

1060: int DASetLocalAdicFunctioni_Private(DA da,int (*ad_lfi)(DALocalInfo*,MatStencil*,void*,void*,void*))
1061: {
1064:   da->adic_lfi = ad_lfi;
1065:   return(0);
1066: }

1068: /*MC
1069:        DASetLocalAdicMFFunctioni - Caches in a DA a local functioni computed by ADIC/ADIFOR

1071:    Collective on DA

1073:    Synopsis:
1074:    int int DASetLocalAdicFunctioni(DA da,int (ad_lf*)(DALocalInfo*,MatStencil*,void*,void*,void*)
1075:    
1076:    Input Parameter:
1077: +  da - initial distributed array
1078: -  admf_lfi - the local matrix-free function as computed by ADIC/ADIFOR

1080:    Level: intermediate

1082: .keywords:  distributed array, refine

1084: .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DAGetLocalFunction(), DASetLocalFunction(),
1085:           DASetLocalJacobian(), DASetLocalFunctioni()
1086: M*/

1090: int DASetLocalAdicMFFunctioni_Private(DA da,int (*admf_lfi)(DALocalInfo*,MatStencil*,void*,void*,void*))
1091: {
1094:   da->adicmf_lfi = admf_lfi;
1095:   return(0);
1096: }

1098: /*MC
1099:        DASetLocalAdicMFFunction - Caches in a DA a local function computed by ADIC/ADIFOR

1101:    Collective on DA

1103:    Synopsis:
1104:    int int DASetLocalAdicMFFunction(DA da,DALocalFunction1 ad_lf)
1105:    
1106:    Input Parameter:
1107: +  da - initial distributed array
1108: -  ad_lf - the local function as computed by ADIC/ADIFOR

1110:    Level: intermediate

1112: .keywords:  distributed array, refine

1114: .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DAGetLocalFunction(), DASetLocalFunction(),
1115:           DASetLocalJacobian()
1116: M*/

1120: int DASetLocalAdicMFFunction_Private(DA da,DALocalFunction1 ad_lf)
1121: {
1124:   da->adicmf_lf = ad_lf;
1125:   return(0);
1126: }

1128: /*@C
1129:        DASetLocalJacobian - Caches in a DA a local Jacobian

1131:    Collective on DA

1133:    
1134:    Input Parameter:
1135: +  da - initial distributed array
1136: -  lj - the local Jacobian

1138:    Level: intermediate

1140:    Notes: The routine SNESDAFormFunction() uses this the cached function to evaluate the user provided function.

1142: .keywords:  distributed array, refine

1144: .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DAGetLocalFunction(), DASetLocalFunction()
1145: @*/
1148: int DASetLocalJacobian(DA da,DALocalFunction1 lj)
1149: {
1152:   da->lj    = lj;
1153:   return(0);
1154: }

1158: /*@C
1159:        DAGetLocalFunction - Gets from a DA a local function and its ADIC/ADIFOR Jacobian

1161:    Collective on DA

1163:    Input Parameter:
1164: .  da - initial distributed array

1166:    Output Parameters:
1167: .  lf - the local function

1169:    Level: intermediate

1171: .keywords:  distributed array, refine

1173: .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DASetLocalFunction()
1174: @*/
1175: int DAGetLocalFunction(DA da,DALocalFunction1 *lf)
1176: {
1179:   if (lf)       *lf = da->lf;
1180:   return(0);
1181: }

1185: /*@
1186:     DAFormFunction1 - Evaluates a user provided function on each processor that 
1187:         share a DA

1189:    Input Parameters:
1190: +    da - the DA that defines the grid
1191: .    vu - input vector
1192: .    vfu - output vector 
1193: -    w - any user data

1195:     Notes: Does NOT do ghost updates on vu upon entry

1197:     Level: advanced

1199: .seealso: DAComputeJacobian1WithAdic()

1201: @*/
1202: int DAFormFunction1(DA da,Vec vu,Vec vfu,void *w)
1203: {
1204:   int         ierr;
1205:   void        *u,*fu;
1206:   DALocalInfo info;
1207: 

1210:   DAGetLocalInfo(da,&info);
1211:   DAVecGetArray(da,vu,&u);
1212:   DAVecGetArray(da,vfu,&fu);

1214:   (*da->lf)(&info,u,fu,w);

1216:   DAVecRestoreArray(da,vu,&u);
1217:   DAVecRestoreArray(da,vfu,&fu);
1218:   return(0);
1219: }

1223: int DAFormFunctioniTest1(DA da,void *w)
1224: {
1225:   Vec         vu,fu,fui;
1226:   int         ierr,i,n;
1227:   PetscScalar *ui,mone = -1.0;
1228:   PetscRandom rnd;
1229:   PetscReal   norm;

1232:   DAGetLocalVector(da,&vu);
1233:   PetscRandomCreate(PETSC_COMM_SELF,RANDOM_DEFAULT,&rnd);
1234:   VecSetRandom(rnd,vu);
1235:   PetscRandomDestroy(rnd);

1237:   DAGetGlobalVector(da,&fu);
1238:   DAGetGlobalVector(da,&fui);
1239: 
1240:   DAFormFunction1(da,vu,fu,w);

1242:   VecGetArray(fui,&ui);
1243:   VecGetLocalSize(fui,&n);
1244:   for (i=0; i<n; i++) {
1245:     DAFormFunctioni1(da,i,vu,ui+i,w);
1246:   }
1247:   VecRestoreArray(fui,&ui);

1249:   VecAXPY(&mone,fu,fui);
1250:   VecNorm(fui,NORM_2,&norm);
1251:   PetscPrintf(da->comm,"Norm of difference in vectors %g\n",norm);
1252:   VecView(fu,0);
1253:   VecView(fui,0);

1255:   DARestoreLocalVector(da,&vu);
1256:   DARestoreGlobalVector(da,&fu);
1257:   DARestoreGlobalVector(da,&fui);
1258:   return(0);
1259: }

1263: /*@
1264:     DAFormFunctioni1 - Evaluates a user provided function

1266:    Input Parameters:
1267: +    da - the DA that defines the grid
1268: .    i - the component of the function we wish to compute (must be local)
1269: .    vu - input vector
1270: .    vfu - output value
1271: -    w - any user data

1273:     Notes: Does NOT do ghost updates on vu upon entry

1275:     Level: advanced

1277: .seealso: DAComputeJacobian1WithAdic()

1279: @*/
1280: int DAFormFunctioni1(DA da,int i,Vec vu,PetscScalar *vfu,void *w)
1281: {
1282:   int         ierr;
1283:   void        *u;
1284:   DALocalInfo info;
1285:   MatStencil  stencil;
1286: 

1289:   DAGetLocalInfo(da,&info);
1290:   DAVecGetArray(da,vu,&u);

1292:   /* figure out stencil value from i */
1293:   stencil.c = i % info.dof;
1294:   stencil.i = (i % (info.xm*info.dof))/info.dof;
1295:   stencil.j = (i % (info.xm*info.ym*info.dof))/(info.xm*info.dof);
1296:   stencil.k = i/(info.xm*info.ym*info.dof);

1298:   (*da->lfi)(&info,&stencil,u,vfu,w);

1300:   DAVecRestoreArray(da,vu,&u);
1301:   return(0);
1302: }

1304: #if defined(new)
1307: /*
1308:   DAGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix where local
1309:     function lives on a DA

1311:         y ~= (F(u + ha) - F(u))/h, 
1312:   where F = nonlinear function, as set by SNESSetFunction()
1313:         u = current iterate
1314:         h = difference interval
1315: */
1316: int DAGetDiagonal_MFFD(DA da,Vec U,Vec a)
1317: {
1318:   PetscScalar  h,*aa,*ww,v;
1319:   PetscReal    epsilon = PETSC_SQRT_MACHINE_EPSILON,umin = 100.0*PETSC_SQRT_MACHINE_EPSILON;
1320:   int          ierr,gI,nI;
1321:   MatStencil   stencil;
1322:   DALocalInfo  info;
1323: 
1325:   (*ctx->func)(0,U,a,ctx->funcctx);
1326:   (*ctx->funcisetbase)(U,ctx->funcctx);

1328:   VecGetArray(U,&ww);
1329:   VecGetArray(a,&aa);
1330: 
1331:   nI = 0;
1332:     h  = ww[gI];
1333:     if (h == 0.0) h = 1.0;
1334: #if !defined(PETSC_USE_COMPLEX)
1335:     if (h < umin && h >= 0.0)      h = umin;
1336:     else if (h < 0.0 && h > -umin) h = -umin;
1337: #else
1338:     if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0)     h = umin;
1339:     else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin;
1340: #endif
1341:     h     *= epsilon;
1342: 
1343:     ww[gI += h;
1344:     (*ctx->funci)(i,w,&v,ctx->funcctx);
1345:     aa[nI]  = (v - aa[nI])/h;
1346:     ww[gI] -= h;
1347:     nI++;
1348:   }
1349:   VecRestoreArray(U,&ww);
1350:   VecRestoreArray(a,&aa);
1351:   return(0);
1352: }
1353: #endif

1355: #if defined(PETSC_HAVE_ADIC) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE)
1356: EXTERN_C_BEGIN
1357: #include "adic/ad_utils.h"
1358: EXTERN_C_END

1362: /*@
1363:     DAComputeJacobian1WithAdic - Evaluates a adiC provided Jacobian function on each processor that 
1364:         share a DA

1366:    Input Parameters:
1367: +    da - the DA that defines the grid
1368: .    vu - input vector (ghosted)
1369: .    J - output matrix
1370: -    w - any user data

1372:    Level: advanced

1374:     Notes: Does NOT do ghost updates on vu upon entry

1376: .seealso: DAFormFunction1()

1378: @*/
1379: int DAComputeJacobian1WithAdic(DA da,Vec vu,Mat J,void *w)
1380: {
1381:   int         ierr,gtdof,tdof;
1382:   PetscScalar *ustart;
1383:   DALocalInfo info;
1384:   void        *ad_u,*ad_f,*ad_ustart,*ad_fstart;
1385:   ISColoring  iscoloring;

1388:   DAGetLocalInfo(da,&info);

1390:   PetscADResetIndep();

1392:   /* get space for derivative objects.  */
1393:   DAGetAdicArray(da,PETSC_TRUE,(void **)&ad_u,&ad_ustart,&gtdof);
1394:   DAGetAdicArray(da,PETSC_FALSE,(void **)&ad_f,&ad_fstart,&tdof);
1395:   VecGetArray(vu,&ustart);
1396:   DAGetColoring(da,IS_COLORING_GHOSTED,&iscoloring);

1398:   PetscADSetValueAndColor(ad_ustart,gtdof,iscoloring->colors,ustart);

1400:   VecRestoreArray(vu,&ustart);
1401:   ISColoringDestroy(iscoloring);
1402:   PetscADIncrementTotalGradSize(iscoloring->n);
1403:   PetscADSetIndepDone();

1405:   DALogEventBegin(DA_LocalADFunction,0,0,0,0);
1406:   (*da->adic_lf)(&info,ad_u,ad_f,w);
1407:   DALogEventEnd(DA_LocalADFunction,0,0,0,0);

1409:   /* stick the values into the matrix */
1410:   MatSetValuesAdic(J,(PetscScalar**)ad_fstart);

1412:   /* return space for derivative objects.  */
1413:   DARestoreAdicArray(da,PETSC_TRUE,(void **)&ad_u,&ad_ustart,&gtdof);
1414:   DARestoreAdicArray(da,PETSC_FALSE,(void **)&ad_f,&ad_fstart,&tdof);
1415:   return(0);
1416: }

1420: /*@C
1421:     DAMultiplyByJacobian1WithAdic - Applies an ADIC-provided Jacobian function to a vector on 
1422:     each processor that shares a DA.

1424:     Input Parameters:
1425: +   da - the DA that defines the grid
1426: .   vu - Jacobian is computed at this point (ghosted)
1427: .   v - product is done on this vector (ghosted)
1428: .   fu - output vector = J(vu)*v (not ghosted)
1429: -   w - any user data

1431:     Notes: 
1432:     This routine does NOT do ghost updates on vu upon entry.

1434:    Level: advanced

1436: .seealso: DAFormFunction1()

1438: @*/
1439: int DAMultiplyByJacobian1WithAdic(DA da,Vec vu,Vec v,Vec f,void *w)
1440: {
1441:   int         ierr,i,gtdof,tdof;
1442:   PetscScalar *avu,*av,*af,*ad_vustart,*ad_fstart;
1443:   DALocalInfo info;
1444:   void        *ad_vu,*ad_f;

1447:   DAGetLocalInfo(da,&info);

1449:   /* get space for derivative objects.  */
1450:   DAGetAdicMFArray(da,PETSC_TRUE,(void **)&ad_vu,(void**)&ad_vustart,&gtdof);
1451:   DAGetAdicMFArray(da,PETSC_FALSE,(void **)&ad_f,(void**)&ad_fstart,&tdof);

1453:   /* copy input vector into derivative object */
1454:   VecGetArray(vu,&avu);
1455:   VecGetArray(v,&av);
1456:   for (i=0; i<gtdof; i++) {
1457:     ad_vustart[2*i]   = avu[i];
1458:     ad_vustart[2*i+1] = av[i];
1459:   }
1460:   VecRestoreArray(vu,&avu);
1461:   VecRestoreArray(v,&av);

1463:   PetscADResetIndep();
1464:   PetscADIncrementTotalGradSize(1);
1465:   PetscADSetIndepDone();

1467:   (*da->adicmf_lf)(&info,ad_vu,ad_f,w);

1469:   /* stick the values into the vector */
1470:   VecGetArray(f,&af);
1471:   for (i=0; i<tdof; i++) {
1472:     af[i] = ad_fstart[2*i+1];
1473:   }
1474:   VecRestoreArray(f,&af);

1476:   /* return space for derivative objects.  */
1477:   DARestoreAdicMFArray(da,PETSC_TRUE,(void **)&ad_vu,(void**)&ad_vustart,&gtdof);
1478:   DARestoreAdicMFArray(da,PETSC_FALSE,(void **)&ad_f,(void**)&ad_fstart,&tdof);
1479:   return(0);
1480: }


1483: #else

1487: int DAComputeJacobian1WithAdic(DA da,Vec vu,Mat J,void *w)
1488: {
1490:   SETERRQ(1,"Must compile with bmake/PETSC_ARCH/packages flag PETSC_HAVE_ADIC for this routine");
1491: }

1495: int DAMultiplyByJacobian1WithAdic(DA da,Vec vu,Vec v,Vec f,void *w)
1496: {
1498:   SETERRQ(1,"Must compile with bmake/PETSC_ARCH/packages flag PETSC_HAVE_ADIC for this routine");
1499: }

1501: #endif

1505: /*@
1506:     DAComputeJacobian1 - Evaluates a local Jacobian function on each processor that 
1507:         share a DA

1509:    Input Parameters:
1510: +    da - the DA that defines the grid
1511: .    vu - input vector (ghosted)
1512: .    J - output matrix
1513: -    w - any user data

1515:     Notes: Does NOT do ghost updates on vu upon entry

1517:     Level: advanced

1519: .seealso: DAFormFunction1()

1521: @*/
1522: int DAComputeJacobian1(DA da,Vec vu,Mat J,void *w)
1523: {
1524:   int         ierr;
1525:   void        *u;
1526:   DALocalInfo info;

1529:   DAGetLocalInfo(da,&info);
1530:   DAVecGetArray(da,vu,&u);
1531:   (*da->lj)(&info,u,J,w);
1532:   DAVecRestoreArray(da,vu,&u);
1533:   return(0);
1534: }


1539: /*
1540:     DAComputeJacobian1WithAdifor - Evaluates a ADIFOR provided Jacobian local function on each processor that 
1541:         share a DA

1543:    Input Parameters:
1544: +    da - the DA that defines the grid
1545: .    vu - input vector (ghosted)
1546: .    J - output matrix
1547: -    w - any user data

1549:     Notes: Does NOT do ghost updates on vu upon entry

1551: .seealso: DAFormFunction1()

1553: */
1554: int DAComputeJacobian1WithAdifor(DA da,Vec vu,Mat J,void *w)
1555: {
1556:   int             i,ierr,Nc,N;
1557:   ISColoringValue *color;
1558:   DALocalInfo     info;
1559:   PetscScalar     *u,*g_u,*g_f,*f,*p_u;
1560:   ISColoring      iscoloring;
1561:   void            (*lf)(int *,DALocalInfo*,PetscScalar*,PetscScalar*,int*,PetscScalar*,PetscScalar*,int*,void*,int*) =
1562:                   (void (*)(int *,DALocalInfo*,PetscScalar*,PetscScalar*,int*,PetscScalar*,PetscScalar*,int*,void*,int*))*da->adifor_lf;

1565:   DAGetColoring(da,IS_COLORING_GHOSTED,&iscoloring);
1566:   Nc   = iscoloring->n;
1567:   DAGetLocalInfo(da,&info);
1568:   N    = info.gxm*info.gym*info.gzm*info.dof;

1570:   /* get space for derivative objects.  */
1571:   PetscMalloc(Nc*info.gxm*info.gym*info.gzm*info.dof*sizeof(PetscScalar),&g_u);
1572:   PetscMemzero(g_u,Nc*info.gxm*info.gym*info.gzm*info.dof*sizeof(PetscScalar));
1573:   p_u   = g_u;
1574:   color = iscoloring->colors;
1575:   for (i=0; i<N; i++) {
1576:     p_u[*color++] = 1.0;
1577:     p_u          += Nc;
1578:   }
1579:   ISColoringDestroy(iscoloring);
1580:   PetscMalloc(Nc*info.xm*info.ym*info.zm*info.dof*sizeof(PetscScalar),&g_f);
1581:   PetscMalloc(info.xm*info.ym*info.zm*info.dof*sizeof(PetscScalar),&f);

1583:   /* Seed the input array g_u with coloring information */
1584: 
1585:   VecGetArray(vu,&u);
1586:   (lf)(&Nc,&info,u,g_u,&Nc,f,g_f,&Nc,w,&ierr);
1587:   VecRestoreArray(vu,&u);

1589:   /* stick the values into the matrix */
1590:   /* PetscScalarView(Nc*info.xm*info.ym,g_f,0); */
1591:   MatSetValuesAdifor(J,Nc,g_f);

1593:   /* return space for derivative objects.  */
1594:   PetscFree(g_u);
1595:   PetscFree(g_f);
1596:   PetscFree(f);
1597:   return(0);
1598: }

1602: /*@C
1603:     DAMultiplyByJacobian1WithAD - Applies a Jacobian function supplied by ADIFOR or ADIC
1604:     to a vector on each processor that shares a DA.

1606:    Input Parameters:
1607: +    da - the DA that defines the grid
1608: .    vu - Jacobian is computed at this point (ghosted)
1609: .    v - product is done on this vector (ghosted)
1610: .    fu - output vector = J(vu)*v (not ghosted)
1611: -    w - any user data

1613:     Notes: 
1614:     This routine does NOT do ghost updates on vu and v upon entry.
1615:            
1616:     Automatically calls DAMultiplyByJacobian1WithAdifor() or DAMultiplyByJacobian1WithAdic()
1617:     depending on whether DASetLocalAdicMFFunction() or DASetLocalAdiforMFFunction() was called.

1619:    Level: advanced

1621: .seealso: DAFormFunction1(), DAMultiplyByJacobian1WithAdifor(), DAMultiplyByJacobian1WithAdic()

1623: @*/
1624: int DAMultiplyByJacobian1WithAD(DA da,Vec u,Vec v,Vec f,void *w)
1625: {
1626:   int         ierr;

1629:   if (da->adicmf_lf) {
1630: #if defined(PETSC_HAVE_ADIC) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE)
1631:     DAMultiplyByJacobian1WithAdic(da,u,v,f,w);
1632: #else
1633:     SETERRQ(1,"Requires ADIC to be installed and cannot use complex numbers");
1634: #endif
1635:   } else if (da->adiformf_lf) {
1636:     DAMultiplyByJacobian1WithAdifor(da,u,v,f,w);
1637:   } else {
1638:     SETERRQ(1,"Must call DASetLocalAdiforMFFunction() or DASetLocalAdicMFFunction() before using");
1639:   }
1640:   return(0);
1641: }


1646: /*@C
1647:     DAMultiplyByJacobian1WithAdifor - Applies a ADIFOR provided Jacobian function on each processor that 
1648:         share a DA to a vector

1650:    Input Parameters:
1651: +    da - the DA that defines the grid
1652: .    vu - Jacobian is computed at this point (ghosted)
1653: .    v - product is done on this vector (ghosted)
1654: .    fu - output vector = J(vu)*v (not ghosted)
1655: -    w - any user data

1657:     Notes: Does NOT do ghost updates on vu and v upon entry

1659:    Level: advanced

1661: .seealso: DAFormFunction1()

1663: @*/
1664: int DAMultiplyByJacobian1WithAdifor(DA da,Vec u,Vec v,Vec f,void *w)
1665: {
1666:   int         ierr;
1667:   PetscScalar *au,*av,*af,*awork;
1668:   Vec         work;
1669:   DALocalInfo info;
1670:   void        (*lf)(DALocalInfo*,PetscScalar*,PetscScalar*,PetscScalar*,PetscScalar*,void*,int*) =
1671:               (void (*)(DALocalInfo*,PetscScalar*,PetscScalar*,PetscScalar*,PetscScalar*,void*,int*))*da->adiformf_lf;

1674:   DAGetLocalInfo(da,&info);

1676:   DAGetGlobalVector(da,&work);
1677:   VecGetArray(u,&au);
1678:   VecGetArray(v,&av);
1679:   VecGetArray(f,&af);
1680:   VecGetArray(work,&awork);
1681:   (lf)(&info,au,av,awork,af,w,&ierr);
1682:   VecRestoreArray(u,&au);
1683:   VecRestoreArray(v,&av);
1684:   VecRestoreArray(f,&af);
1685:   VecRestoreArray(work,&awork);
1686:   DARestoreGlobalVector(da,&work);

1688:   return(0);
1689: }

1693: /*@C
1694:        DASetInterpolationType - Sets the type of interpolation that will be 
1695:           returned by DAGetInterpolation()

1697:    Collective on DA

1699:    Input Parameter:
1700: +  da - initial distributed array
1701: .  ctype - DA_Q1 and DA_Q0 are currently the only supported forms

1703:    Level: intermediate

1705:    Notes: you should call this on the coarser of the two DAs you pass to DAGetInterpolation()

1707: .keywords:  distributed array, interpolation

1709: .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DA, DAInterpolationType
1710: @*/
1711: int DASetInterpolationType(DA da,DAInterpolationType ctype)
1712: {
1715:   da->interptype = ctype;
1716:   return(0);
1717: }