Actual source code: da2.c

  1: #define PETSCDM_DLL
  2: 
 3:  #include private/daimpl.h

  7: PetscErrorCode DAView_2d(DA da,PetscViewer viewer)
  8: {
 10:   PetscMPIInt    rank;
 11:   PetscTruth     iascii,isdraw;

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

 16:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
 17:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);
 18:   if (iascii) {
 19:     PetscViewerFormat format;

 21:     PetscViewerGetFormat(viewer, &format);
 22:     if (format != PETSC_VIEWER_ASCII_VTK && format != PETSC_VIEWER_ASCII_VTK_CELL) {
 23:       DALocalInfo info;
 24:       DAGetLocalInfo(da,&info);
 25:       PetscViewerASCIISynchronizedPrintf(viewer,"Processor [%d] M %D N %D m %D n %D w %D s %D\n",rank,da->M,
 26:                                                 da->N,da->m,da->n,da->w,da->s);
 27:       PetscViewerASCIISynchronizedPrintf(viewer,"X range of indices: %D %D, Y range of indices: %D %D\n",info.xs,info.xs+info.xm,info.ys,info.ys+info.ym);
 28:       PetscViewerFlush(viewer);
 29:     }
 30:   } else if (isdraw) {
 31:     PetscDraw       draw;
 32:     double     ymin = -1*da->s-1,ymax = da->N+da->s;
 33:     double     xmin = -1*da->s-1,xmax = da->M+da->s;
 34:     double     x,y;
 35:     PetscInt   base,*idx;
 36:     char       node[10];
 37:     PetscTruth isnull;
 38: 
 39:     PetscViewerDrawGetDraw(viewer,0,&draw);
 40:     PetscDrawIsNull(draw,&isnull); if (isnull) return(0);
 41:     if (!da->coordinates) {
 42:       PetscDrawSetCoordinates(draw,xmin,ymin,xmax,ymax);
 43:     }
 44:     PetscDrawSynchronizedClear(draw);

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

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

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

 77:     PetscDrawSynchronizedFlush(draw);
 78:     PetscDrawPause(draw);
 79:     /* overlay ghost numbers, useful for error checking */
 80:     /* put in numbers */

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

101: /*
102:       M is number of grid points 
103:       m is number of processors

105: */
108: PetscErrorCode  DASplitComm2d(MPI_Comm comm,PetscInt M,PetscInt N,PetscInt sw,MPI_Comm *outcomm)
109: {
111:   PetscInt       m,n = 0,x = 0,y = 0;
112:   PetscMPIInt    size,csize,rank;

115:   MPI_Comm_size(comm,&size);
116:   MPI_Comm_rank(comm,&rank);

118:   csize = 4*size;
119:   do {
120:     if (csize % 4) SETERRQ4(PETSC_ERR_ARG_INCOMP,"Cannot split communicator of size %d tried %d %D %D",size,csize,x,y);
121:     csize   = csize/4;
122: 
123:     m = (PetscInt)(0.5 + sqrt(((double)M)*((double)csize)/((double)N)));
124:     if (!m) m = 1;
125:     while (m > 0) {
126:       n = csize/m;
127:       if (m*n == csize) break;
128:       m--;
129:     }
130:     if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;}

132:     x = M/m + ((M % m) > ((csize-1) % m));
133:     y = (N + (csize-1)/m)/n;
134:   } while ((x < 4 || y < 4) && csize > 1);
135:   if (size != csize) {
136:     MPI_Group    entire_group,sub_group;
137:     PetscMPIInt  i,*groupies;

139:     MPI_Comm_group(comm,&entire_group);
140:     PetscMalloc(csize*sizeof(PetscInt),&groupies);
141:     for (i=0; i<csize; i++) {
142:       groupies[i] = (rank/csize)*csize + i;
143:     }
144:     MPI_Group_incl(entire_group,csize,groupies,&sub_group);
145:     PetscFree(groupies);
146:     MPI_Comm_create(comm,sub_group,outcomm);
147:     MPI_Group_free(&entire_group);
148:     MPI_Group_free(&sub_group);
149:     PetscInfo1(0,"DASplitComm2d:Creating redundant coarse problems of size %d\n",csize);
150:   } else {
151:     *outcomm = comm;
152:   }
153:   return(0);
154: }

158: PetscErrorCode DAGetElements_2d_P1(DA da,PetscInt *n,const PetscInt *e[])
159: {
161:   PetscInt       i,j,cnt,xs,xe = da->xe,ys,ye = da->ye,Xs = da->Xs, Xe = da->Xe, Ys = da->Ys;

164:   if (!da->e) {
165:     if (da->xs == Xs) xs = da->xs; else xs = da->xs - 1;
166:     if (da->ys == Ys) ys = da->ys; else ys = da->ys - 1;
167:     da->ne = 2*(xe - xs - 1)*(ye - ys - 1);
168:     PetscMalloc((1 + 3*da->ne)*sizeof(PetscInt),&da->e);
169:     cnt    = 0;
170:     for (j=ys; j<ye-1; j++) {
171:       for (i=xs; i<xe-1; i++) {
172:         da->e[cnt]   = i - Xs + (j - Ys)*(Xe - Xs);
173:         da->e[cnt+1] = i - Xs + 1 + (j - Ys)*(Xe - Xs);
174:         da->e[cnt+2] = i - Xs + (j - Ys + 1)*(Xe - Xs);

176:         da->e[cnt+3] = i - Xs + 1 + (j - Ys + 1)*(Xe - Xs);
177:         da->e[cnt+4] = i - Xs + (j - Ys + 1)*(Xe - Xs);
178:         da->e[cnt+5] = i - Xs + 1 + (j - Ys)*(Xe - Xs);
179:         cnt += 6;
180:       }
181:     }
182:   }
183:   *n = da->ne;
184:   *e = da->e;
185:   return(0);
186: }

190: /*@C
191:        DASetLocalFunction - Caches in a DA a local function. 

193:    Collective on DA

195:    Input Parameter:
196: +  da - initial distributed array
197: -  lf - the local function

199:    Level: intermediate

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

203: .keywords:  distributed array, refine

205: .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DAGetLocalFunction(), DASetLocalFunctioni()
206: @*/
207: PetscErrorCode  DASetLocalFunction(DA da,DALocalFunction1 lf)
208: {
211:   da->lf    = lf;
212:   return(0);
213: }

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

220:    Collective on DA

222:    Input Parameter:
223: +  da - initial distributed array
224: -  lfi - the local function

226:    Level: intermediate

228: .keywords:  distributed array, refine

230: .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DAGetLocalFunction(), DASetLocalFunction()
231: @*/
232: PetscErrorCode  DASetLocalFunctioni(DA da,PetscErrorCode (*lfi)(DALocalInfo*,MatStencil*,void*,PetscScalar*,void*))
233: {
236:   da->lfi = lfi;
237:   return(0);
238: }

242: /*@C
243:        DASetLocalFunctionib - Caches in a DA a block local function that evaluates a single component

245:    Collective on DA

247:    Input Parameter:
248: +  da - initial distributed array
249: -  lfi - the local function

251:    Level: intermediate

253: .keywords:  distributed array, refine

255: .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DAGetLocalFunction(), DASetLocalFunction()
256: @*/
257: PetscErrorCode  DASetLocalFunctionib(DA da,PetscErrorCode (*lfi)(DALocalInfo*,MatStencil*,void*,PetscScalar*,void*))
258: {
261:   da->lfib = lfi;
262:   return(0);
263: }

267: PetscErrorCode DASetLocalAdicFunction_Private(DA da,DALocalFunction1 ad_lf)
268: {
271:   da->adic_lf = ad_lf;
272:   return(0);
273: }

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

278:    Collective on DA

280:    Synopsis:
281:    PetscErrorCode DASetLocalAdicFunctioni(DA da,PetscInt (ad_lf*)(DALocalInfo*,MatStencil*,void*,void*,void*)
282:    
283:    Input Parameter:
284: +  da - initial distributed array
285: -  ad_lfi - the local function as computed by ADIC/ADIFOR

287:    Level: intermediate

289: .keywords:  distributed array, refine

291: .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DAGetLocalFunction(), DASetLocalFunction(),
292:           DASetLocalJacobian(), DASetLocalFunctioni()
293: M*/

297: PetscErrorCode DASetLocalAdicFunctioni_Private(DA da,PetscErrorCode (*ad_lfi)(DALocalInfo*,MatStencil*,void*,void*,void*))
298: {
301:   da->adic_lfi = ad_lfi;
302:   return(0);
303: }

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

308:    Collective on DA

310:    Synopsis:
311:    PetscErrorCode  DASetLocalAdicFunctioni(DA da,int (ad_lf*)(DALocalInfo*,MatStencil*,void*,void*,void*)
312:    
313:    Input Parameter:
314: +  da - initial distributed array
315: -  admf_lfi - the local matrix-free function as computed by ADIC/ADIFOR

317:    Level: intermediate

319: .keywords:  distributed array, refine

321: .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DAGetLocalFunction(), DASetLocalFunction(),
322:           DASetLocalJacobian(), DASetLocalFunctioni()
323: M*/

327: PetscErrorCode DASetLocalAdicMFFunctioni_Private(DA da,PetscErrorCode (*admf_lfi)(DALocalInfo*,MatStencil*,void*,void*,void*))
328: {
331:   da->adicmf_lfi = admf_lfi;
332:   return(0);
333: }

335: /*MC
336:        DASetLocalAdicFunctionib - Caches in a DA a block local functioni computed by ADIC/ADIFOR

338:    Collective on DA

340:    Synopsis:
341:    PetscErrorCode DASetLocalAdicFunctionib(DA da,PetscInt (ad_lf*)(DALocalInfo*,MatStencil*,void*,void*,void*)
342:    
343:    Input Parameter:
344: +  da - initial distributed array
345: -  ad_lfi - the local function as computed by ADIC/ADIFOR

347:    Level: intermediate

349: .keywords:  distributed array, refine

351: .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DAGetLocalFunction(), DASetLocalFunction(),
352:           DASetLocalJacobian(), DASetLocalFunctionib()
353: M*/

357: PetscErrorCode DASetLocalAdicFunctionib_Private(DA da,PetscErrorCode (*ad_lfi)(DALocalInfo*,MatStencil*,void*,void*,void*))
358: {
361:   da->adic_lfib = ad_lfi;
362:   return(0);
363: }

365: /*MC
366:        DASetLocalAdicMFFunctionib - Caches in a DA a block local functioni computed by ADIC/ADIFOR

368:    Collective on DA

370:    Synopsis:
371:    PetscErrorCode  DASetLocalAdicFunctionib(DA da,int (ad_lf*)(DALocalInfo*,MatStencil*,void*,void*,void*)
372:    
373:    Input Parameter:
374: +  da - initial distributed array
375: -  admf_lfi - the local matrix-free function as computed by ADIC/ADIFOR

377:    Level: intermediate

379: .keywords:  distributed array, refine

381: .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DAGetLocalFunction(), DASetLocalFunction(),
382:           DASetLocalJacobian(), DASetLocalFunctionib()
383: M*/

387: PetscErrorCode DASetLocalAdicMFFunctionib_Private(DA da,PetscErrorCode (*admf_lfi)(DALocalInfo*,MatStencil*,void*,void*,void*))
388: {
391:   da->adicmf_lfib = admf_lfi;
392:   return(0);
393: }

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

398:    Collective on DA

400:    Synopsis:
401:    PetscErrorCode DASetLocalAdicMFFunction(DA da,DALocalFunction1 ad_lf)
402:    
403:    Input Parameter:
404: +  da - initial distributed array
405: -  ad_lf - the local function as computed by ADIC/ADIFOR

407:    Level: intermediate

409: .keywords:  distributed array, refine

411: .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DAGetLocalFunction(), DASetLocalFunction(),
412:           DASetLocalJacobian()
413: M*/

417: PetscErrorCode DASetLocalAdicMFFunction_Private(DA da,DALocalFunction1 ad_lf)
418: {
421:   da->adicmf_lf = ad_lf;
422:   return(0);
423: }

425: /*@C
426:        DASetLocalJacobian - Caches in a DA a local Jacobian

428:    Collective on DA

430:    
431:    Input Parameter:
432: +  da - initial distributed array
433: -  lj - the local Jacobian

435:    Level: intermediate

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

439: .keywords:  distributed array, refine

441: .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DAGetLocalFunction(), DASetLocalFunction()
442: @*/
445: PetscErrorCode  DASetLocalJacobian(DA da,DALocalFunction1 lj)
446: {
449:   da->lj    = lj;
450:   return(0);
451: }

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

458:    Collective on DA

460:    Input Parameter:
461: .  da - initial distributed array

463:    Output Parameter:
464: .  lf - the local function

466:    Level: intermediate

468: .keywords:  distributed array, refine

470: .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DAGetLocalJacobian(), DASetLocalFunction()
471: @*/
472: PetscErrorCode  DAGetLocalFunction(DA da,DALocalFunction1 *lf)
473: {
476:   if (lf)       *lf = da->lf;
477:   return(0);
478: }

482: /*@C
483:        DAGetLocalJacobian - Gets from a DA a local jacobian

485:    Collective on DA

487:    Input Parameter:
488: .  da - initial distributed array

490:    Output Parameter:
491: .  lj - the local jacobian

493:    Level: intermediate

495: .keywords:  distributed array, refine

497: .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DAGetLocalFunction(), DASetLocalJacobian()
498: @*/
499: PetscErrorCode  DAGetLocalJacobian(DA da,DALocalFunction1 *lj)
500: {
503:   if (lj) *lj = da->lj;
504:   return(0);
505: }

509: /*@
510:     DAFormFunction - Evaluates a user provided function on each processor that 
511:         share a DA

513:    Input Parameters:
514: +    da - the DA that defines the grid
515: .    vu - input vector
516: .    vfu - output vector 
517: -    w - any user data

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

521:            This should eventually replace DAFormFunction1

523:     Level: advanced

525: .seealso: DAComputeJacobian1WithAdic()

527: @*/
528: PetscErrorCode  DAFormFunction(DA da,PetscErrorCode (*lf)(void),Vec vu,Vec vfu,void *w)
529: {
531:   void           *u,*fu;
532:   DALocalInfo    info;
533:   PetscErrorCode (*f)(DALocalInfo*,void*,void*,void*) = (PetscErrorCode (*)(DALocalInfo*,void*,void*,void*))lf;
534: 
536:   DAGetLocalInfo(da,&info);
537:   DAVecGetArray(da,vu,&u);
538:   DAVecGetArray(da,vfu,&fu);

540:   (*f)(&info,u,fu,w);
541:   if (PetscExceptionValue(ierr)) {
542:     PetscErrorCode pDAVecRestoreArray(da,vu,&u);CHKERRQ(pierr);
543:     pDAVecRestoreArray(da,vfu,&fu);CHKERRQ(pierr);
544:   }
545: 

547:   DAVecRestoreArray(da,vu,&u);
548:   DAVecRestoreArray(da,vfu,&fu);
549:   return(0);
550: }

554: /*@C 
555:    DAFormFunctionLocal - This is a universal function evaluation routine for
556:    a local DA function.

558:    Collective on DA

560:    Input Parameters:
561: +  da - the DA context
562: .  func - The local function
563: .  X - input vector
564: .  F - function vector
565: -  ctx - A user context

567:    Level: intermediate

569: .seealso: DASetLocalFunction(), DASetLocalJacobian(), DASetLocalAdicFunction(), DASetLocalAdicMFFunction(),
570:           SNESSetFunction(), SNESSetJacobian()

572: @*/
573: PetscErrorCode  DAFormFunctionLocal(DA da, DALocalFunction1 func, Vec X, Vec F, void *ctx)
574: {
575:   Vec            localX;
576:   DALocalInfo    info;
577:   void          *u;
578:   void          *fu;

582:   DAGetLocalVector(da,&localX);
583:   /*
584:      Scatter ghost points to local vector, using the 2-step process
585:         DAGlobalToLocalBegin(), DAGlobalToLocalEnd().
586:   */
587:   DAGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
588:   DAGlobalToLocalEnd(da,X,INSERT_VALUES,localX);
589:   DAGetLocalInfo(da,&info);
590:   DAVecGetArray(da,localX,&u);
591:   DAVecGetArray(da,F,&fu);
592:   (*func)(&info,u,fu,ctx);
593:   if (PetscExceptionValue(ierr)) {
594:     PetscErrorCode pDAVecRestoreArray(da,localX,&u);CHKERRQ(pierr);
595:     pDAVecRestoreArray(da,F,&fu);CHKERRQ(pierr);
596:   }
597: 
598:   DAVecRestoreArray(da,localX,&u);
599:   DAVecRestoreArray(da,F,&fu);
600:   if (PetscExceptionValue(ierr)) {
601:     PetscErrorCode pDARestoreLocalVector(da,&localX);CHKERRQ(pierr);
602:   }
603: 
604:   DARestoreLocalVector(da,&localX);
605:   return(0);
606: }

610: /*@C 
611:    DAFormFunctionLocalGhost - This is a universal function evaluation routine for
612:    a local DA function, but the ghost values of the output are communicated and added.

614:    Collective on DA

616:    Input Parameters:
617: +  da - the DA context
618: .  func - The local function
619: .  X - input vector
620: .  F - function vector
621: -  ctx - A user context

623:    Level: intermediate

625: .seealso: DASetLocalFunction(), DASetLocalJacobian(), DASetLocalAdicFunction(), DASetLocalAdicMFFunction(),
626:           SNESSetFunction(), SNESSetJacobian()

628: @*/
629: PetscErrorCode  DAFormFunctionLocalGhost(DA da, DALocalFunction1 func, Vec X, Vec F, void *ctx)
630: {
631:   Vec            localX, localF;
632:   DALocalInfo    info;
633:   void          *u;
634:   void          *fu;

638:   DAGetLocalVector(da,&localX);
639:   DAGetLocalVector(da,&localF);
640:   /*
641:      Scatter ghost points to local vector, using the 2-step process
642:         DAGlobalToLocalBegin(), DAGlobalToLocalEnd().
643:   */
644:   DAGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
645:   DAGlobalToLocalEnd(da,X,INSERT_VALUES,localX);
646:   VecSet(F, 0.0);
647:   VecSet(localF, 0.0);
648:   DAGetLocalInfo(da,&info);
649:   DAVecGetArray(da,localX,&u);
650:   DAVecGetArray(da,localF,&fu);
651:   (*func)(&info,u,fu,ctx);
652:   if (PetscExceptionValue(ierr)) {
653:     PetscErrorCode pDAVecRestoreArray(da,localX,&u);CHKERRQ(pierr);
654:     pDAVecRestoreArray(da,localF,&fu);CHKERRQ(pierr);
655:   }
656: 
657:   DALocalToGlobalBegin(da,localF,F);
658:   DALocalToGlobalEnd(da,localF,F);
659:   DAVecRestoreArray(da,localX,&u);
660:   DAVecRestoreArray(da,localF,&fu);
661:   if (PetscExceptionValue(ierr)) {
662:     PetscErrorCode pDARestoreLocalVector(da,&localX);CHKERRQ(pierr);
663:   DARestoreLocalVector(da,&localF);
664:   }
665: 
666:   DARestoreLocalVector(da,&localX);
667:   DARestoreLocalVector(da,&localF);
668:   return(0);
669: }

673: /*@
674:     DAFormFunction1 - Evaluates a user provided function on each processor that 
675:         share a DA

677:    Input Parameters:
678: +    da - the DA that defines the grid
679: .    vu - input vector
680: .    vfu - output vector 
681: -    w - any user data

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

685:     Level: advanced

687: .seealso: DAComputeJacobian1WithAdic()

689: @*/
690: PetscErrorCode  DAFormFunction1(DA da,Vec vu,Vec vfu,void *w)
691: {
693:   void           *u,*fu;
694:   DALocalInfo    info;
695: 

698:   DAGetLocalInfo(da,&info);
699:   DAVecGetArray(da,vu,&u);
700:   DAVecGetArray(da,vfu,&fu);

702:   CHKMEMQ;
703:   (*da->lf)(&info,u,fu,w);
704:   if (PetscExceptionValue(ierr)) {
705:     PetscErrorCode pDAVecRestoreArray(da,vu,&u);CHKERRQ(pierr);
706:     pDAVecRestoreArray(da,vfu,&fu);CHKERRQ(pierr);
707:   }
708: 
709:   CHKMEMQ;

711:   DAVecRestoreArray(da,vu,&u);
712:   DAVecRestoreArray(da,vfu,&fu);
713:   return(0);
714: }

718: PetscErrorCode  DAFormFunctioniTest1(DA da,void *w)
719: {
720:   Vec            vu,fu,fui;
722:   PetscInt       i,n;
723:   PetscScalar    *ui;
724:   PetscRandom    rnd;
725:   PetscReal      norm;

728:   DAGetLocalVector(da,&vu);
729:   PetscRandomCreate(PETSC_COMM_SELF,&rnd);
730:   PetscRandomSetFromOptions(rnd);
731:   VecSetRandom(vu,rnd);
732:   PetscRandomDestroy(rnd);

734:   DAGetGlobalVector(da,&fu);
735:   DAGetGlobalVector(da,&fui);
736: 
737:   DAFormFunction1(da,vu,fu,w);

739:   VecGetArray(fui,&ui);
740:   VecGetLocalSize(fui,&n);
741:   for (i=0; i<n; i++) {
742:     DAFormFunctioni1(da,i,vu,ui+i,w);
743:   }
744:   VecRestoreArray(fui,&ui);

746:   VecAXPY(fui,-1.0,fu);
747:   VecNorm(fui,NORM_2,&norm);
748:   PetscPrintf(((PetscObject)da)->comm,"Norm of difference in vectors %G\n",norm);
749:   VecView(fu,0);
750:   VecView(fui,0);

752:   DARestoreLocalVector(da,&vu);
753:   DARestoreGlobalVector(da,&fu);
754:   DARestoreGlobalVector(da,&fui);
755:   return(0);
756: }

760: /*@
761:     DAFormFunctioni1 - Evaluates a user provided point-wise function

763:    Input Parameters:
764: +    da - the DA that defines the grid
765: .    i - the component of the function we wish to compute (must be local)
766: .    vu - input vector
767: .    vfu - output value
768: -    w - any user data

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

772:     Level: advanced

774: .seealso: DAComputeJacobian1WithAdic()

776: @*/
777: PetscErrorCode  DAFormFunctioni1(DA da,PetscInt i,Vec vu,PetscScalar *vfu,void *w)
778: {
780:   void           *u;
781:   DALocalInfo    info;
782:   MatStencil     stencil;
783: 

786:   DAGetLocalInfo(da,&info);
787:   DAVecGetArray(da,vu,&u);

789:   /* figure out stencil value from i */
790:   stencil.c = i % info.dof;
791:   stencil.i = (i % (info.xm*info.dof))/info.dof;
792:   stencil.j = (i % (info.xm*info.ym*info.dof))/(info.xm*info.dof);
793:   stencil.k = i/(info.xm*info.ym*info.dof);

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

797:   DAVecRestoreArray(da,vu,&u);
798:   return(0);
799: }

803: /*@
804:     DAFormFunctionib1 - Evaluates a user provided point-block function

806:    Input Parameters:
807: +    da - the DA that defines the grid
808: .    i - the component of the function we wish to compute (must be local)
809: .    vu - input vector
810: .    vfu - output value
811: -    w - any user data

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

815:     Level: advanced

817: .seealso: DAComputeJacobian1WithAdic()

819: @*/
820: PetscErrorCode  DAFormFunctionib1(DA da,PetscInt i,Vec vu,PetscScalar *vfu,void *w)
821: {
823:   void           *u;
824:   DALocalInfo    info;
825:   MatStencil     stencil;
826: 
828:   DAGetLocalInfo(da,&info);
829:   DAVecGetArray(da,vu,&u);

831:   /* figure out stencil value from i */
832:   stencil.c = i % info.dof;
833:   if (stencil.c) SETERRQ(PETSC_ERR_ARG_WRONG,"Point-block functions can only be called for the entire block");
834:   stencil.i = (i % (info.xm*info.dof))/info.dof;
835:   stencil.j = (i % (info.xm*info.ym*info.dof))/(info.xm*info.dof);
836:   stencil.k = i/(info.xm*info.ym*info.dof);

838:   (*da->lfib)(&info,&stencil,u,vfu,w);

840:   DAVecRestoreArray(da,vu,&u);
841:   return(0);
842: }

844: #if defined(new)
847: /*
848:   DAGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix where local
849:     function lives on a DA

851:         y ~= (F(u + ha) - F(u))/h, 
852:   where F = nonlinear function, as set by SNESSetFunction()
853:         u = current iterate
854:         h = difference interval
855: */
856: PetscErrorCode DAGetDiagonal_MFFD(DA da,Vec U,Vec a)
857: {
858:   PetscScalar    h,*aa,*ww,v;
859:   PetscReal      epsilon = PETSC_SQRT_MACHINE_EPSILON,umin = 100.0*PETSC_SQRT_MACHINE_EPSILON;
861:   PetscInt       gI,nI;
862:   MatStencil     stencil;
863:   DALocalInfo    info;
864: 
866:   (*ctx->func)(0,U,a,ctx->funcctx);
867:   (*ctx->funcisetbase)(U,ctx->funcctx);

869:   VecGetArray(U,&ww);
870:   VecGetArray(a,&aa);
871: 
872:   nI = 0;
873:     h  = ww[gI];
874:     if (h == 0.0) h = 1.0;
875: #if !defined(PETSC_USE_COMPLEX)
876:     if (h < umin && h >= 0.0)      h = umin;
877:     else if (h < 0.0 && h > -umin) h = -umin;
878: #else
879:     if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0)     h = umin;
880:     else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin;
881: #endif
882:     h     *= epsilon;
883: 
884:     ww[gI += h;
885:     (*ctx->funci)(i,w,&v,ctx->funcctx);
886:     aa[nI]  = (v - aa[nI])/h;
887:     ww[gI] -= h;
888:     nI++;
889:   }
890:   VecRestoreArray(U,&ww);
891:   VecRestoreArray(a,&aa);
892:   return(0);
893: }
894: #endif

896: #if defined(PETSC_HAVE_ADIC)
898: #include "adic/ad_utils.h"

903: /*@C
904:     DAComputeJacobian1WithAdic - Evaluates a adiC provided Jacobian function on each processor that 
905:         share a DA

907:    Input Parameters:
908: +    da - the DA that defines the grid
909: .    vu - input vector (ghosted)
910: .    J - output matrix
911: -    w - any user data

913:    Level: advanced

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

917: .seealso: DAFormFunction1()

919: @*/
920: PetscErrorCode  DAComputeJacobian1WithAdic(DA da,Vec vu,Mat J,void *w)
921: {
923:   PetscInt       gtdof,tdof;
924:   PetscScalar    *ustart;
925:   DALocalInfo    info;
926:   void           *ad_u,*ad_f,*ad_ustart,*ad_fstart;
927:   ISColoring     iscoloring;

930:   DAGetLocalInfo(da,&info);

932:   PetscADResetIndep();

934:   /* get space for derivative objects.  */
935:   DAGetAdicArray(da,PETSC_TRUE,(void **)&ad_u,&ad_ustart,&gtdof);
936:   DAGetAdicArray(da,PETSC_FALSE,(void **)&ad_f,&ad_fstart,&tdof);
937:   VecGetArray(vu,&ustart);
938:   DAGetColoring(da,IS_COLORING_GHOSTED,MATAIJ,&iscoloring);

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

942:   VecRestoreArray(vu,&ustart);
943:   ISColoringDestroy(iscoloring);
944:   PetscADIncrementTotalGradSize(iscoloring->n);
945:   PetscADSetIndepDone();

947:   PetscLogEventBegin(DA_LocalADFunction,0,0,0,0);
948:   (*da->adic_lf)(&info,ad_u,ad_f,w);
949:   PetscLogEventEnd(DA_LocalADFunction,0,0,0,0);

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

954:   /* return space for derivative objects.  */
955:   DARestoreAdicArray(da,PETSC_TRUE,(void **)&ad_u,&ad_ustart,&gtdof);
956:   DARestoreAdicArray(da,PETSC_FALSE,(void **)&ad_f,&ad_fstart,&tdof);
957:   return(0);
958: }

962: /*@C
963:     DAMultiplyByJacobian1WithAdic - Applies an ADIC-provided Jacobian function to a vector on 
964:     each processor that shares a DA.

966:     Input Parameters:
967: +   da - the DA that defines the grid
968: .   vu - Jacobian is computed at this point (ghosted)
969: .   v - product is done on this vector (ghosted)
970: .   fu - output vector = J(vu)*v (not ghosted)
971: -   w - any user data

973:     Notes: 
974:     This routine does NOT do ghost updates on vu upon entry.

976:    Level: advanced

978: .seealso: DAFormFunction1()

980: @*/
981: PetscErrorCode  DAMultiplyByJacobian1WithAdic(DA da,Vec vu,Vec v,Vec f,void *w)
982: {
984:   PetscInt       i,gtdof,tdof;
985:   PetscScalar    *avu,*av,*af,*ad_vustart,*ad_fstart;
986:   DALocalInfo    info;
987:   void           *ad_vu,*ad_f;

990:   DAGetLocalInfo(da,&info);

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

996:   /* copy input vector into derivative object */
997:   VecGetArray(vu,&avu);
998:   VecGetArray(v,&av);
999:   for (i=0; i<gtdof; i++) {
1000:     ad_vustart[2*i]   = avu[i];
1001:     ad_vustart[2*i+1] = av[i];
1002:   }
1003:   VecRestoreArray(vu,&avu);
1004:   VecRestoreArray(v,&av);

1006:   PetscADResetIndep();
1007:   PetscADIncrementTotalGradSize(1);
1008:   PetscADSetIndepDone();

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

1012:   /* stick the values into the vector */
1013:   VecGetArray(f,&af);
1014:   for (i=0; i<tdof; i++) {
1015:     af[i] = ad_fstart[2*i+1];
1016:   }
1017:   VecRestoreArray(f,&af);

1019:   /* return space for derivative objects.  */
1020:   DARestoreAdicMFArray(da,PETSC_TRUE,(void **)&ad_vu,(void**)&ad_vustart,&gtdof);
1021:   DARestoreAdicMFArray(da,PETSC_FALSE,(void **)&ad_f,(void**)&ad_fstart,&tdof);
1022:   return(0);
1023: }
1024: #endif

1028: /*@
1029:     DAComputeJacobian1 - Evaluates a local Jacobian function on each processor that 
1030:         share a DA

1032:    Input Parameters:
1033: +    da - the DA that defines the grid
1034: .    vu - input vector (ghosted)
1035: .    J - output matrix
1036: -    w - any user data

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

1040:     Level: advanced

1042: .seealso: DAFormFunction1()

1044: @*/
1045: PetscErrorCode  DAComputeJacobian1(DA da,Vec vu,Mat J,void *w)
1046: {
1048:   void           *u;
1049:   DALocalInfo    info;

1052:   DAGetLocalInfo(da,&info);
1053:   DAVecGetArray(da,vu,&u);
1054:   (*da->lj)(&info,u,J,w);
1055:   DAVecRestoreArray(da,vu,&u);
1056:   return(0);
1057: }


1062: /*
1063:     DAComputeJacobian1WithAdifor - Evaluates a ADIFOR provided Jacobian local function on each processor that 
1064:         share a DA

1066:    Input Parameters:
1067: +    da - the DA that defines the grid
1068: .    vu - input vector (ghosted)
1069: .    J - output matrix
1070: -    w - any user data

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

1074: .seealso: DAFormFunction1()

1076: */
1077: PetscErrorCode  DAComputeJacobian1WithAdifor(DA da,Vec vu,Mat J,void *w)
1078: {
1079:   PetscErrorCode  ierr;
1080:   PetscInt        i,Nc,N;
1081:   ISColoringValue *color;
1082:   DALocalInfo     info;
1083:   PetscScalar     *u,*g_u,*g_f,*f = 0,*p_u;
1084:   ISColoring      iscoloring;
1085:   void            (*lf)(PetscInt*,DALocalInfo*,PetscScalar*,PetscScalar*,PetscInt*,PetscScalar*,PetscScalar*,PetscInt*,void*,PetscErrorCode*) =
1086:                   (void (*)(PetscInt*,DALocalInfo*,PetscScalar*,PetscScalar*,PetscInt*,PetscScalar*,PetscScalar*,PetscInt*,void*,PetscErrorCode*))*da->adifor_lf;

1089:   DAGetColoring(da,IS_COLORING_GHOSTED,MATAIJ,&iscoloring);
1090:   Nc   = iscoloring->n;
1091:   DAGetLocalInfo(da,&info);
1092:   N    = info.gxm*info.gym*info.gzm*info.dof;

1094:   /* get space for derivative objects.  */
1095:   PetscMalloc(Nc*info.gxm*info.gym*info.gzm*info.dof*sizeof(PetscScalar),&g_u);
1096:   PetscMemzero(g_u,Nc*info.gxm*info.gym*info.gzm*info.dof*sizeof(PetscScalar));
1097:   p_u   = g_u;
1098:   color = iscoloring->colors;
1099:   for (i=0; i<N; i++) {
1100:     p_u[*color++] = 1.0;
1101:     p_u          += Nc;
1102:   }
1103:   ISColoringDestroy(iscoloring);
1104:   PetscMalloc2(Nc*info.xm*info.ym*info.zm*info.dof,PetscScalar,&g_f,info.xm*info.ym*info.zm*info.dof,PetscScalar,&f);

1106:   /* Seed the input array g_u with coloring information */
1107: 
1108:   VecGetArray(vu,&u);
1109:   (lf)(&Nc,&info,u,g_u,&Nc,f,g_f,&Nc,w,&ierr);
1110:   VecRestoreArray(vu,&u);

1112:   /* stick the values into the matrix */
1113:   /* PetscScalarView(Nc*info.xm*info.ym,g_f,0); */
1114:   MatSetValuesAdifor(J,Nc,g_f);

1116:   /* return space for derivative objects.  */
1117:   PetscFree(g_u);
1118:   PetscFree2(g_f,f);
1119:   return(0);
1120: }

1124: /*@C 
1125:    DAFormjacobianLocal - This is a universal Jacobian evaluation routine for
1126:    a local DA function.

1128:    Collective on DA

1130:    Input Parameters:
1131: +  da - the DA context
1132: .  func - The local function
1133: .  X - input vector
1134: .  J - Jacobian matrix
1135: -  ctx - A user context

1137:    Level: intermediate

1139: .seealso: DASetLocalFunction(), DASetLocalJacobian(), DASetLocalAdicFunction(), DASetLocalAdicMFFunction(),
1140:           SNESSetFunction(), SNESSetJacobian()

1142: @*/
1143: PetscErrorCode  DAFormJacobianLocal(DA da, DALocalFunction1 func, Vec X, Mat J, void *ctx)
1144: {
1145:   Vec            localX;
1146:   DALocalInfo    info;
1147:   void          *u;

1151:   DAGetLocalVector(da,&localX);
1152:   /*
1153:      Scatter ghost points to local vector, using the 2-step process
1154:         DAGlobalToLocalBegin(), DAGlobalToLocalEnd().
1155:   */
1156:   DAGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
1157:   DAGlobalToLocalEnd(da,X,INSERT_VALUES,localX);
1158:   DAGetLocalInfo(da,&info);
1159:   DAVecGetArray(da,localX,&u);
1160:   (*func)(&info,u,J,ctx);
1161:   if (PetscExceptionValue(ierr)) {
1162:     PetscErrorCode pDAVecRestoreArray(da,localX,&u);CHKERRQ(pierr);
1163:   }
1164: 
1165:   DAVecRestoreArray(da,localX,&u);
1166:   if (PetscExceptionValue(ierr)) {
1167:     PetscErrorCode pDARestoreLocalVector(da,&localX);CHKERRQ(pierr);
1168:   }
1169: 
1170:   DARestoreLocalVector(da,&localX);
1171:   return(0);
1172: }

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

1180:    Input Parameters:
1181: +    da - the DA that defines the grid
1182: .    vu - Jacobian is computed at this point (ghosted)
1183: .    v - product is done on this vector (ghosted)
1184: .    fu - output vector = J(vu)*v (not ghosted)
1185: -    w - any user data

1187:     Notes: 
1188:     This routine does NOT do ghost updates on vu and v upon entry.
1189:            
1190:     Automatically calls DAMultiplyByJacobian1WithAdifor() or DAMultiplyByJacobian1WithAdic()
1191:     depending on whether DASetLocalAdicMFFunction() or DASetLocalAdiforMFFunction() was called.

1193:    Level: advanced

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

1197: @*/
1198: PetscErrorCode  DAMultiplyByJacobian1WithAD(DA da,Vec u,Vec v,Vec f,void *w)
1199: {

1203:   if (da->adicmf_lf) {
1204: #if defined(PETSC_HAVE_ADIC)
1205:     DAMultiplyByJacobian1WithAdic(da,u,v,f,w);
1206: #else
1207:     SETERRQ(PETSC_ERR_SUP_SYS,"Requires ADIC to be installed and cannot use complex numbers");
1208: #endif
1209:   } else if (da->adiformf_lf) {
1210:     DAMultiplyByJacobian1WithAdifor(da,u,v,f,w);
1211:   } else {
1212:     SETERRQ(PETSC_ERR_ORDER,"Must call DASetLocalAdiforMFFunction() or DASetLocalAdicMFFunction() before using");
1213:   }
1214:   return(0);
1215: }


1220: /*@C
1221:     DAMultiplyByJacobian1WithAdifor - Applies a ADIFOR provided Jacobian function on each processor that 
1222:         share a DA to a vector

1224:    Input Parameters:
1225: +    da - the DA that defines the grid
1226: .    vu - Jacobian is computed at this point (ghosted)
1227: .    v - product is done on this vector (ghosted)
1228: .    fu - output vector = J(vu)*v (not ghosted)
1229: -    w - any user data

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

1233:    Level: advanced

1235: .seealso: DAFormFunction1()

1237: @*/
1238: PetscErrorCode  DAMultiplyByJacobian1WithAdifor(DA da,Vec u,Vec v,Vec f,void *w)
1239: {
1241:   PetscScalar    *au,*av,*af,*awork;
1242:   Vec            work;
1243:   DALocalInfo    info;
1244:   void           (*lf)(DALocalInfo*,PetscScalar*,PetscScalar*,PetscScalar*,PetscScalar*,void*,PetscErrorCode*) =
1245:                  (void (*)(DALocalInfo*,PetscScalar*,PetscScalar*,PetscScalar*,PetscScalar*,void*,PetscErrorCode*))*da->adiformf_lf;

1248:   DAGetLocalInfo(da,&info);

1250:   DAGetGlobalVector(da,&work);
1251:   VecGetArray(u,&au);
1252:   VecGetArray(v,&av);
1253:   VecGetArray(f,&af);
1254:   VecGetArray(work,&awork);
1255:   (lf)(&info,au,av,awork,af,w,&ierr);
1256:   VecRestoreArray(u,&au);
1257:   VecRestoreArray(v,&av);
1258:   VecRestoreArray(f,&af);
1259:   VecRestoreArray(work,&awork);
1260:   DARestoreGlobalVector(da,&work);

1262:   return(0);
1263: }

1268: PetscErrorCode  DACreate_2D(DA da)
1269: {
1270:   const PetscInt       dim          = da->dim;
1271:   const PetscInt       M            = da->M;
1272:   const PetscInt       N            = da->N;
1273:   PetscInt             m            = da->m;
1274:   PetscInt             n            = da->n;
1275:   const PetscInt       dof          = da->w;
1276:   const PetscInt       s            = da->s;
1277:   const DAPeriodicType wrap         = da->wrap;
1278:   const DAStencilType  stencil_type = da->stencil_type;
1279:   PetscInt            *lx           = da->lx;
1280:   PetscInt            *ly           = da->ly;
1281:   MPI_Comm             comm;
1282:   PetscMPIInt    rank,size;
1283:   PetscInt       xs,xe,ys,ye,x,y,Xs,Xe,Ys,Ye,start,end;
1284:   PetscInt       up,down,left,i,n0,n1,n2,n3,n5,n6,n7,n8,*idx,nn;
1285:   PetscInt       xbase,*bases,*ldims,j,x_t,y_t,s_t,base,count;
1286:   PetscInt       s_x,s_y; /* s proportionalized to w */
1287:   PetscInt       sn0 = 0,sn2 = 0,sn6 = 0,sn8 = 0;
1288:   Vec            local,global;
1289:   VecScatter     ltog,gtol;
1290:   IS             to,from;

1294: #ifndef PETSC_USE_DYNAMIC_LIBRARIES
1295:   DMInitializePackage(PETSC_NULL);
1296: #endif

1298:   if (dim != PETSC_DECIDE && dim != 2) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Dimension should be 2: %D",dim);
1299:   if (dof < 1) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Must have 1 or more degrees of freedom per node: %D",dof);
1300:   if (s < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Stencil width cannot be negative: %D",s);

1302:   PetscObjectGetComm((PetscObject) da, &comm);
1303:   MPI_Comm_size(comm,&size);
1304:   MPI_Comm_rank(comm,&rank);

1306:   da->ops->getelements = DAGetElements_2d_P1;

1308:   da->dim         = 2;
1309:   da->elementtype = DA_ELEMENT_P1;
1310:   PetscMalloc(dof*sizeof(char*),&da->fieldname);
1311:   PetscMemzero(da->fieldname,dof*sizeof(char*));

1313:   if (m != PETSC_DECIDE) {
1314:     if (m < 1) {SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in X direction: %D",m);}
1315:     else if (m > size) {SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in X direction: %D %d",m,size);}
1316:   }
1317:   if (n != PETSC_DECIDE) {
1318:     if (n < 1) {SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Y direction: %D",n);}
1319:     else if (n > size) {SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Y direction: %D %d",n,size);}
1320:   }

1322:   if (m == PETSC_DECIDE || n == PETSC_DECIDE) {
1323:     if (n != PETSC_DECIDE) {
1324:       m = size/n;
1325:     } else if (m != PETSC_DECIDE) {
1326:       n = size/m;
1327:     } else {
1328:       /* try for squarish distribution */
1329:       m = (PetscInt)(0.5 + sqrt(((double)M)*((double)size)/((double)N)));
1330:       if (!m) m = 1;
1331:       while (m > 0) {
1332:         n = size/m;
1333:         if (m*n == size) break;
1334:         m--;
1335:       }
1336:       if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;}
1337:     }
1338:     if (m*n != size) SETERRQ(PETSC_ERR_PLIB,"Unable to create partition, check the size of the communicator and input m and n ");
1339:   } else if (m*n != size) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Given Bad partition");

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

1344:   /* 
1345:      Determine locally owned region 
1346:      xs is the first local node number, x is the number of local nodes 
1347:   */
1348:   if (!lx) {
1349:     PetscMalloc(m*sizeof(PetscInt), &da->lx);
1350:     lx = da->lx;
1351:     for (i=0; i<m; i++) {
1352:       lx[i] = M/m + ((M % m) > i);
1353:     }
1354:   }
1355:   x  = lx[rank % m];
1356:   xs = 0;
1357:   for (i=0; i<(rank % m); i++) {
1358:     xs += lx[i];
1359:   }
1360: #if defined(PETSC_USE_DEBUG)
1361:   left = xs;
1362:   for (i=(rank % m); i<m; i++) {
1363:     left += lx[i];
1364:   }
1365:   if (left != M) {
1366:     SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Sum of lx across processors not equal to M: %D %D",left,M);
1367:   }
1368: #endif

1370:   /* 
1371:      Determine locally owned region 
1372:      ys is the first local node number, y is the number of local nodes 
1373:   */
1374:   if (!ly) {
1375:     PetscMalloc(n*sizeof(PetscInt), &da->ly);
1376:     ly = da->ly;
1377:     for (i=0; i<n; i++) {
1378:       ly[i] = N/n + ((N % n) > i);
1379:     }
1380:   }
1381:   y  = ly[rank/m];
1382:   ys = 0;
1383:   for (i=0; i<(rank/m); i++) {
1384:     ys += ly[i];
1385:   }
1386: #if defined(PETSC_USE_DEBUG)
1387:   left = ys;
1388:   for (i=(rank/m); i<n; i++) {
1389:     left += ly[i];
1390:   }
1391:   if (left != N) {
1392:     SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Sum of ly across processors not equal to N: %D %D",left,N);
1393:   }
1394: #endif

1396:   if (x < s) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Local x-width of domain x %D is smaller than stencil width s %D",x,s);
1397:   if (y < s) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Local y-width of domain y %D is smaller than stencil width s %D",y,s);
1398:   xe = xs + x;
1399:   ye = ys + y;

1401:   /* determine ghost region */
1402:   /* Assume No Periodicity */
1403:   if (xs-s > 0) Xs = xs - s; else Xs = 0;
1404:   if (ys-s > 0) Ys = ys - s; else Ys = 0;
1405:   if (xe+s <= M) Xe = xe + s; else Xe = M;
1406:   if (ye+s <= N) Ye = ye + s; else Ye = N;

1408:   /* X Periodic */
1409:   if (DAXPeriodic(wrap)){
1410:     Xs = xs - s;
1411:     Xe = xe + s;
1412:   }

1414:   /* Y Periodic */
1415:   if (DAYPeriodic(wrap)){
1416:     Ys = ys - s;
1417:     Ye = ye + s;
1418:   }

1420:   /* Resize all X parameters to reflect w */
1421:   x   *= dof;
1422:   xs  *= dof;
1423:   xe  *= dof;
1424:   Xs  *= dof;
1425:   Xe  *= dof;
1426:   s_x = s*dof;
1427:   s_y = s;

1429:   /* determine starting point of each processor */
1430:   nn    = x*y;
1431:   PetscMalloc2(size+1,PetscInt,&bases,size,PetscInt,&ldims);
1432:   MPI_Allgather(&nn,1,MPIU_INT,ldims,1,MPIU_INT,comm);
1433:   bases[0] = 0;
1434:   for (i=1; i<=size; i++) {
1435:     bases[i] = ldims[i-1];
1436:   }
1437:   for (i=1; i<=size; i++) {
1438:     bases[i] += bases[i-1];
1439:   }

1441:   /* allocate the base parallel and sequential vectors */
1442:   da->Nlocal = x*y;
1443:   VecCreateMPIWithArray(comm,da->Nlocal,PETSC_DECIDE,0,&global);
1444:   VecSetBlockSize(global,dof);
1445:   da->nlocal = (Xe-Xs)*(Ye-Ys);
1446:   VecCreateSeqWithArray(PETSC_COMM_SELF,da->nlocal,0,&local);
1447:   VecSetBlockSize(local,dof);


1450:   /* generate appropriate vector scatters */
1451:   /* local to global inserts non-ghost point region into global */
1452:   VecGetOwnershipRange(global,&start,&end);
1453:   ISCreateStride(comm,x*y,start,1,&to);

1455:   left  = xs - Xs; down  = ys - Ys; up    = down + y;
1456:   PetscMalloc(x*(up - down)*sizeof(PetscInt),&idx);
1457:   count = 0;
1458:   for (i=down; i<up; i++) {
1459:     for (j=0; j<x/dof; j++) {
1460:       idx[count++] = left + i*(Xe-Xs) + j*dof;
1461:     }
1462:   }
1463:   ISCreateBlock(comm,dof,count,idx,&from);
1464:   PetscFree(idx);

1466:   VecScatterCreate(local,from,global,to,&ltog);
1467:   PetscLogObjectParent(da,ltog);
1468:   ISDestroy(from);
1469:   ISDestroy(to);

1471:   /* global to local must include ghost points */
1472:   if (stencil_type == DA_STENCIL_BOX) {
1473:     ISCreateStride(comm,(Xe-Xs)*(Ye-Ys),0,1,&to);
1474:   } else {
1475:     /* must drop into cross shape region */
1476:     /*       ---------|
1477:             |  top    |
1478:          |---         ---|
1479:          |   middle      |
1480:          |               |
1481:          ----         ----
1482:             | bottom  |
1483:             -----------
1484:         Xs xs        xe  Xe */
1485:     /* bottom */
1486:     left  = xs - Xs; down = ys - Ys; up    = down + y;
1487:     count = down*(xe-xs) + (up-down)*(Xe-Xs) + (Ye-Ys-up)*(xe-xs);
1488:     PetscMalloc(count*sizeof(PetscInt)/dof,&idx);
1489:     count = 0;
1490:     for (i=0; i<down; i++) {
1491:       for (j=0; j<xe-xs; j += dof) {
1492:         idx[count++] = left + i*(Xe-Xs) + j;
1493:       }
1494:     }
1495:     /* middle */
1496:     for (i=down; i<up; i++) {
1497:       for (j=0; j<Xe-Xs; j += dof) {
1498:         idx[count++] = i*(Xe-Xs) + j;
1499:       }
1500:     }
1501:     /* top */
1502:     for (i=up; i<Ye-Ys; i++) {
1503:       for (j=0; j<xe-xs; j += dof) {
1504:         idx[count++] = left + i*(Xe-Xs) + j;
1505:       }
1506:     }
1507:     ISCreateBlock(comm,dof,count,idx,&to);
1508:     PetscFree(idx);
1509:   }


1512:   /* determine who lies on each side of us stored in    n6 n7 n8
1513:                                                         n3    n5
1514:                                                         n0 n1 n2
1515:   */

1517:   /* Assume the Non-Periodic Case */
1518:   n1 = rank - m;
1519:   if (rank % m) {
1520:     n0 = n1 - 1;
1521:   } else {
1522:     n0 = -1;
1523:   }
1524:   if ((rank+1) % m) {
1525:     n2 = n1 + 1;
1526:     n5 = rank + 1;
1527:     n8 = rank + m + 1; if (n8 >= m*n) n8 = -1;
1528:   } else {
1529:     n2 = -1; n5 = -1; n8 = -1;
1530:   }
1531:   if (rank % m) {
1532:     n3 = rank - 1;
1533:     n6 = n3 + m; if (n6 >= m*n) n6 = -1;
1534:   } else {
1535:     n3 = -1; n6 = -1;
1536:   }
1537:   n7 = rank + m; if (n7 >= m*n) n7 = -1;


1540:   /* Modify for Periodic Cases */
1541:   if (wrap == DA_YPERIODIC) {  /* Handle Top and Bottom Sides */
1542:     if (n1 < 0) n1 = rank + m * (n-1);
1543:     if (n7 < 0) n7 = rank - m * (n-1);
1544:     if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1;
1545:     if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1;
1546:     if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1;
1547:     if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1;
1548:   } else if (wrap == DA_XPERIODIC) { /* Handle Left and Right Sides */
1549:     if (n3 < 0) n3 = rank + (m-1);
1550:     if (n5 < 0) n5 = rank - (m-1);
1551:     if ((n1 >= 0) && (n0 < 0)) n0 = rank-1;
1552:     if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1;
1553:     if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1;
1554:     if ((n7 >= 0) && (n8 < 0)) n8 = rank+1;
1555:   } else if (wrap == DA_XYPERIODIC) {

1557:     /* Handle all four corners */
1558:     if ((n6 < 0) && (n7 < 0) && (n3 < 0)) n6 = m-1;
1559:     if ((n8 < 0) && (n7 < 0) && (n5 < 0)) n8 = 0;
1560:     if ((n2 < 0) && (n5 < 0) && (n1 < 0)) n2 = size-m;
1561:     if ((n0 < 0) && (n3 < 0) && (n1 < 0)) n0 = size-1;

1563:     /* Handle Top and Bottom Sides */
1564:     if (n1 < 0) n1 = rank + m * (n-1);
1565:     if (n7 < 0) n7 = rank - m * (n-1);
1566:     if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1;
1567:     if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1;
1568:     if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1;
1569:     if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1;

1571:     /* Handle Left and Right Sides */
1572:     if (n3 < 0) n3 = rank + (m-1);
1573:     if (n5 < 0) n5 = rank - (m-1);
1574:     if ((n1 >= 0) && (n0 < 0)) n0 = rank-1;
1575:     if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1;
1576:     if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1;
1577:     if ((n7 >= 0) && (n8 < 0)) n8 = rank+1;
1578:   }
1579:   PetscMalloc(9*sizeof(PetscInt),&da->neighbors);
1580:   da->neighbors[0] = n0;
1581:   da->neighbors[1] = n1;
1582:   da->neighbors[2] = n2;
1583:   da->neighbors[3] = n3;
1584:   da->neighbors[4] = rank;
1585:   da->neighbors[5] = n5;
1586:   da->neighbors[6] = n6;
1587:   da->neighbors[7] = n7;
1588:   da->neighbors[8] = n8;

1590:   if (stencil_type == DA_STENCIL_STAR) {
1591:     /* save corner processor numbers */
1592:     sn0 = n0; sn2 = n2; sn6 = n6; sn8 = n8;
1593:     n0 = n2 = n6 = n8 = -1;
1594:   }

1596:   PetscMalloc((x+2*s_x)*(y+2*s_y)*sizeof(PetscInt),&idx);
1597:   PetscLogObjectMemory(da,(x+2*s_x)*(y+2*s_y)*sizeof(PetscInt));
1598:   nn = 0;

1600:   xbase = bases[rank];
1601:   for (i=1; i<=s_y; i++) {
1602:     if (n0 >= 0) { /* left below */
1603:       x_t = lx[n0 % m]*dof;
1604:       y_t = ly[(n0/m)];
1605:       s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x;
1606:       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1607:     }
1608:     if (n1 >= 0) { /* directly below */
1609:       x_t = x;
1610:       y_t = ly[(n1/m)];
1611:       s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t;
1612:       for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1613:     }
1614:     if (n2 >= 0) { /* right below */
1615:       x_t = lx[n2 % m]*dof;
1616:       y_t = ly[(n2/m)];
1617:       s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t;
1618:       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1619:     }
1620:   }

1622:   for (i=0; i<y; i++) {
1623:     if (n3 >= 0) { /* directly left */
1624:       x_t = lx[n3 % m]*dof;
1625:       /* y_t = y; */
1626:       s_t = bases[n3] + (i+1)*x_t - s_x;
1627:       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1628:     }

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

1632:     if (n5 >= 0) { /* directly right */
1633:       x_t = lx[n5 % m]*dof;
1634:       /* y_t = y; */
1635:       s_t = bases[n5] + (i)*x_t;
1636:       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1637:     }
1638:   }

1640:   for (i=1; i<=s_y; i++) {
1641:     if (n6 >= 0) { /* left above */
1642:       x_t = lx[n6 % m]*dof;
1643:       /* y_t = ly[(n6/m)]; */
1644:       s_t = bases[n6] + (i)*x_t - s_x;
1645:       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1646:     }
1647:     if (n7 >= 0) { /* directly above */
1648:       x_t = x;
1649:       /* y_t = ly[(n7/m)]; */
1650:       s_t = bases[n7] + (i-1)*x_t;
1651:       for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1652:     }
1653:     if (n8 >= 0) { /* right above */
1654:       x_t = lx[n8 % m]*dof;
1655:       /* y_t = ly[(n8/m)]; */
1656:       s_t = bases[n8] + (i-1)*x_t;
1657:       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1658:     }
1659:   }

1661:   base = bases[rank];
1662:   {
1663:     PetscInt nnn = nn/dof,*iidx;
1664:     PetscMalloc(nnn*sizeof(PetscInt),&iidx);
1665:     for (i=0; i<nnn; i++) {
1666:       iidx[i] = idx[dof*i];
1667:     }
1668:     ISCreateBlock(comm,dof,nnn,iidx,&from);
1669:     PetscFree(iidx);
1670:   }
1671:   VecScatterCreate(global,from,local,to,&gtol);
1672:   PetscLogObjectParent(da,gtol);
1673:   ISDestroy(to);
1674:   ISDestroy(from);

1676:   if (stencil_type == DA_STENCIL_STAR) {
1677:     /*
1678:         Recompute the local to global mappings, this time keeping the 
1679:       information about the cross corner processor numbers.
1680:     */
1681:     n0 = sn0; n2 = sn2; n6 = sn6; n8 = sn8;
1682:     nn = 0;
1683:     xbase = bases[rank];
1684:     for (i=1; i<=s_y; i++) {
1685:       if (n0 >= 0) { /* left below */
1686:         x_t = lx[n0 % m]*dof;
1687:         y_t = ly[(n0/m)];
1688:         s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x;
1689:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1690:       }
1691:       if (n1 >= 0) { /* directly below */
1692:         x_t = x;
1693:         y_t = ly[(n1/m)];
1694:         s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t;
1695:         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1696:       }
1697:       if (n2 >= 0) { /* right below */
1698:         x_t = lx[n2 % m]*dof;
1699:         y_t = ly[(n2/m)];
1700:         s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t;
1701:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1702:       }
1703:     }

1705:     for (i=0; i<y; i++) {
1706:       if (n3 >= 0) { /* directly left */
1707:         x_t = lx[n3 % m]*dof;
1708:         /* y_t = y; */
1709:         s_t = bases[n3] + (i+1)*x_t - s_x;
1710:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1711:       }

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

1715:       if (n5 >= 0) { /* directly right */
1716:         x_t = lx[n5 % m]*dof;
1717:         /* y_t = y; */
1718:         s_t = bases[n5] + (i)*x_t;
1719:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1720:       }
1721:     }

1723:     for (i=1; i<=s_y; i++) {
1724:       if (n6 >= 0) { /* left above */
1725:         x_t = lx[n6 % m]*dof;
1726:         /* y_t = ly[(n6/m)]; */
1727:         s_t = bases[n6] + (i)*x_t - s_x;
1728:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1729:       }
1730:       if (n7 >= 0) { /* directly above */
1731:         x_t = x;
1732:         /* y_t = ly[(n7/m)]; */
1733:         s_t = bases[n7] + (i-1)*x_t;
1734:         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1735:       }
1736:       if (n8 >= 0) { /* right above */
1737:         x_t = lx[n8 % m]*dof;
1738:         /* y_t = ly[(n8/m)]; */
1739:         s_t = bases[n8] + (i-1)*x_t;
1740:         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1741:       }
1742:     }
1743:   }
1744:   PetscFree2(bases,ldims);

1746:   da->m  = m;  da->n  = n;
1747:   da->xs = xs; da->xe = xe; da->ys = ys; da->ye = ye; da->zs = 0; da->ze = 1;
1748:   da->Xs = Xs; da->Xe = Xe; da->Ys = Ys; da->Ye = Ye; da->Zs = 0; da->Ze = 1;

1750:   VecDestroy(local);
1751:   VecDestroy(global);

1753:   da->gtol      = gtol;
1754:   da->ltog      = ltog;
1755:   da->idx       = idx;
1756:   da->Nl        = nn;
1757:   da->base      = base;
1758:   da->ops->view = DAView_2d;

1760:   /* 
1761:      Set the local to global ordering in the global vector, this allows use
1762:      of VecSetValuesLocal().
1763:   */
1764:   ISLocalToGlobalMappingCreateNC(comm,nn,idx,&da->ltogmap);
1765:   ISLocalToGlobalMappingBlock(da->ltogmap,da->w,&da->ltogmapb);
1766:   PetscLogObjectParent(da,da->ltogmap);

1768:   da->ltol = PETSC_NULL;
1769:   da->ao   = PETSC_NULL;

1771:   PetscPublishAll(da);
1772:   return(0);
1773: }

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

1782:    Collective on MPI_Comm

1784:    Input Parameters:
1785: +  comm - MPI communicator
1786: .  wrap - type of periodicity should the array have. 
1787:          Use one of DA_NONPERIODIC, DA_XPERIODIC, DA_YPERIODIC, or DA_XYPERIODIC.
1788: .  stencil_type - stencil type.  Use either DA_STENCIL_BOX or DA_STENCIL_STAR.
1789: .  M,N - global dimension in each direction of the array (use -M and or -N to indicate that it may be set to a different value 
1790:             from the command line with -da_grid_x <M> -da_grid_y <N>)
1791: .  m,n - corresponding number of processors in each dimension 
1792:          (or PETSC_DECIDE to have calculated)
1793: .  dof - number of degrees of freedom per node
1794: .  s - stencil width
1795: -  lx, ly - arrays containing the number of nodes in each cell along
1796:            the x and y coordinates, or PETSC_NULL. If non-null, these
1797:            must be of length as m and n, and the corresponding
1798:            m and n cannot be PETSC_DECIDE. The sum of the lx[] entries
1799:            must be M, and the sum of the ly[] entries must be N.

1801:    Output Parameter:
1802: .  da - the resulting distributed array object

1804:    Options Database Key:
1805: +  -da_view - Calls DAView() at the conclusion of DACreate2d()
1806: .  -da_grid_x <nx> - number of grid points in x direction, if M < 0
1807: .  -da_grid_y <ny> - number of grid points in y direction, if N < 0
1808: .  -da_processors_x <nx> - number of processors in x direction
1809: .  -da_processors_y <ny> - number of processors in y direction
1810: .  -da_refine_x - refinement ratio in x direction
1811: -  -da_refine_y - refinement ratio in y direction

1813:    Level: beginner

1815:    Notes:
1816:    The stencil type DA_STENCIL_STAR with width 1 corresponds to the 
1817:    standard 5-pt stencil, while DA_STENCIL_BOX with width 1 denotes
1818:    the standard 9-pt stencil.

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

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

1826: .seealso: DADestroy(), DAView(), DACreate1d(), DACreate3d(), DAGlobalToLocalBegin(), DAGetRefinementFactor(),
1827:           DAGlobalToLocalEnd(), DALocalToGlobal(), DALocalToLocalBegin(), DALocalToLocalEnd(), DASetRefinementFactor(),
1828:           DAGetInfo(), DACreateGlobalVector(), DACreateLocalVector(), DACreateNaturalVector(), DALoad(), DAView(), DAGetOwnershipRanges()

1830: @*/
1831: PetscErrorCode  DACreate2d(MPI_Comm comm,DAPeriodicType wrap,DAStencilType stencil_type,
1832:                           PetscInt M,PetscInt N,PetscInt m,PetscInt n,PetscInt dof,PetscInt s,const PetscInt lx[],const PetscInt ly[],DA *da)
1833: {

1837:   DACreate(comm, da);
1838:   DASetDim(*da, 2);
1839:   DASetSizes(*da, M, N, PETSC_DECIDE);
1840:   DASetNumProcs(*da, m, n, PETSC_DECIDE);
1841:   DASetPeriodicity(*da, wrap);
1842:   DASetDof(*da, dof);
1843:   DASetStencilType(*da, stencil_type);
1844:   DASetStencilWidth(*da, s);
1845:   DASetVertexDivision(*da, lx, ly, PETSC_NULL);
1846:   /* This violates the behavior for other classes, but right now users expect negative dimensions to be handled this way */
1847:   DASetFromOptions(*da);
1848:   DASetType(*da, DA2D);
1849:   return(0);
1850: }