xref: /petsc/src/dm/impls/da/da2.c (revision d461ba9780c55a0df6dc85ac9e4e0589fbf99e79)
19a42bb27SBarry Smith 
2c6db04a5SJed Brown #include <private/daimpl.h>    /*I   "petscdmda.h"   I*/
347c6ae99SBarry Smith 
447c6ae99SBarry Smith #undef __FUNCT__
59a42bb27SBarry Smith #define __FUNCT__ "DMView_DA_2d"
69a42bb27SBarry Smith PetscErrorCode DMView_DA_2d(DM da,PetscViewer viewer)
747c6ae99SBarry Smith {
847c6ae99SBarry Smith   PetscErrorCode ierr;
947c6ae99SBarry Smith   PetscMPIInt    rank;
109a42bb27SBarry Smith   PetscBool      iascii,isdraw,isbinary;
1147c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
129a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
139a42bb27SBarry Smith   PetscBool      ismatlab;
149a42bb27SBarry Smith #endif
1547c6ae99SBarry Smith 
1647c6ae99SBarry Smith   PetscFunctionBegin;
1747c6ae99SBarry Smith   ierr = MPI_Comm_rank(((PetscObject)da)->comm,&rank);CHKERRQ(ierr);
1847c6ae99SBarry Smith 
1947c6ae99SBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
2047c6ae99SBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
219a42bb27SBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
229a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
239a42bb27SBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERMATLAB,&ismatlab);CHKERRQ(ierr);
249a42bb27SBarry Smith #endif
2547c6ae99SBarry Smith   if (iascii) {
2647c6ae99SBarry Smith     PetscViewerFormat format;
2747c6ae99SBarry Smith 
2847c6ae99SBarry Smith     ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr);
2947c6ae99SBarry Smith     if (format != PETSC_VIEWER_ASCII_VTK && format != PETSC_VIEWER_ASCII_VTK_CELL) {
30aa219208SBarry Smith       DMDALocalInfo info;
31aa219208SBarry Smith       ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
327b23a99aSBarry Smith       ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
3347c6ae99SBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Processor [%d] M %D N %D m %D n %D w %D s %D\n",rank,dd->M,dd->N,dd->m,dd->n,dd->w,dd->s);CHKERRQ(ierr);
3447c6ae99SBarry Smith       ierr = 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);CHKERRQ(ierr);
3547c6ae99SBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
367b23a99aSBarry Smith       ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr);
373da9ae13SJed Brown     } else {
383da9ae13SJed Brown       ierr = DMView_DA_VTK(da,viewer);CHKERRQ(ierr);
3947c6ae99SBarry Smith     }
4047c6ae99SBarry Smith   } else if (isdraw) {
4147c6ae99SBarry Smith     PetscDraw  draw;
4247c6ae99SBarry Smith     double     ymin = -1*dd->s-1,ymax = dd->N+dd->s;
4347c6ae99SBarry Smith     double     xmin = -1*dd->s-1,xmax = dd->M+dd->s;
4447c6ae99SBarry Smith     double     x,y;
4547c6ae99SBarry Smith     PetscInt   base,*idx;
4647c6ae99SBarry Smith     char       node[10];
4747c6ae99SBarry Smith     PetscBool  isnull;
4847c6ae99SBarry Smith 
4947c6ae99SBarry Smith     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
5047c6ae99SBarry Smith     ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
5147c6ae99SBarry Smith     if (!dd->coordinates) {
5247c6ae99SBarry Smith       ierr = PetscDrawSetCoordinates(draw,xmin,ymin,xmax,ymax);CHKERRQ(ierr);
5347c6ae99SBarry Smith     }
5447c6ae99SBarry Smith     ierr = PetscDrawSynchronizedClear(draw);CHKERRQ(ierr);
5547c6ae99SBarry Smith 
5647c6ae99SBarry Smith     /* first processor draw all node lines */
5747c6ae99SBarry Smith     if (!rank) {
5847c6ae99SBarry Smith       ymin = 0.0; ymax = dd->N - 1;
5947c6ae99SBarry Smith       for (xmin=0; xmin<dd->M; xmin++) {
6047c6ae99SBarry Smith         ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_BLACK);CHKERRQ(ierr);
6147c6ae99SBarry Smith       }
6247c6ae99SBarry Smith       xmin = 0.0; xmax = dd->M - 1;
6347c6ae99SBarry Smith       for (ymin=0; ymin<dd->N; ymin++) {
6447c6ae99SBarry Smith         ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_BLACK);CHKERRQ(ierr);
6547c6ae99SBarry Smith       }
6647c6ae99SBarry Smith     }
6747c6ae99SBarry Smith     ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr);
6847c6ae99SBarry Smith     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
6947c6ae99SBarry Smith 
7047c6ae99SBarry Smith     /* draw my box */
7147c6ae99SBarry Smith     ymin = dd->ys; ymax = dd->ye - 1; xmin = dd->xs/dd->w;
7247c6ae99SBarry Smith     xmax =(dd->xe-1)/dd->w;
7347c6ae99SBarry Smith     ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_RED);CHKERRQ(ierr);
7447c6ae99SBarry Smith     ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_RED);CHKERRQ(ierr);
7547c6ae99SBarry Smith     ierr = PetscDrawLine(draw,xmin,ymax,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr);
7647c6ae99SBarry Smith     ierr = PetscDrawLine(draw,xmax,ymin,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr);
7747c6ae99SBarry Smith 
7847c6ae99SBarry Smith     /* put in numbers */
7947c6ae99SBarry Smith     base = (dd->base)/dd->w;
8047c6ae99SBarry Smith     for (y=ymin; y<=ymax; y++) {
8147c6ae99SBarry Smith       for (x=xmin; x<=xmax; x++) {
8247c6ae99SBarry Smith         sprintf(node,"%d",(int)base++);
8347c6ae99SBarry Smith         ierr = PetscDrawString(draw,x,y,PETSC_DRAW_BLACK,node);CHKERRQ(ierr);
8447c6ae99SBarry Smith       }
8547c6ae99SBarry Smith     }
8647c6ae99SBarry Smith 
8747c6ae99SBarry Smith     ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr);
8847c6ae99SBarry Smith     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
8947c6ae99SBarry Smith     /* overlay ghost numbers, useful for error checking */
9047c6ae99SBarry Smith     /* put in numbers */
9147c6ae99SBarry Smith 
9247c6ae99SBarry Smith     base = 0; idx = dd->idx;
9347c6ae99SBarry Smith     ymin = dd->Ys; ymax = dd->Ye; xmin = dd->Xs; xmax = dd->Xe;
9447c6ae99SBarry Smith     for (y=ymin; y<ymax; y++) {
9547c6ae99SBarry Smith       for (x=xmin; x<xmax; x++) {
9647c6ae99SBarry Smith         if ((base % dd->w) == 0) {
9747c6ae99SBarry Smith           sprintf(node,"%d",(int)(idx[base]/dd->w));
9847c6ae99SBarry Smith           ierr = PetscDrawString(draw,x/dd->w,y,PETSC_DRAW_BLUE,node);CHKERRQ(ierr);
9947c6ae99SBarry Smith         }
10047c6ae99SBarry Smith         base++;
10147c6ae99SBarry Smith       }
10247c6ae99SBarry Smith     }
10347c6ae99SBarry Smith     ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr);
10447c6ae99SBarry Smith     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
1059a42bb27SBarry Smith   } else if (isbinary){
1069a42bb27SBarry Smith     ierr = DMView_DA_Binary(da,viewer);CHKERRQ(ierr);
1079a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
1089a42bb27SBarry Smith   } else if (ismatlab) {
1099a42bb27SBarry Smith     ierr = DMView_DA_Matlab(da,viewer);CHKERRQ(ierr);
1109a42bb27SBarry Smith #endif
111aa219208SBarry Smith   } else SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for DMDA 1d",((PetscObject)viewer)->type_name);
11247c6ae99SBarry Smith   PetscFunctionReturn(0);
11347c6ae99SBarry Smith }
11447c6ae99SBarry Smith 
11547c6ae99SBarry Smith /*
11647c6ae99SBarry Smith       M is number of grid points
11747c6ae99SBarry Smith       m is number of processors
11847c6ae99SBarry Smith 
11947c6ae99SBarry Smith */
12047c6ae99SBarry Smith #undef __FUNCT__
121aa219208SBarry Smith #define __FUNCT__ "DMDASplitComm2d"
1227087cfbeSBarry Smith PetscErrorCode  DMDASplitComm2d(MPI_Comm comm,PetscInt M,PetscInt N,PetscInt sw,MPI_Comm *outcomm)
12347c6ae99SBarry Smith {
12447c6ae99SBarry Smith   PetscErrorCode ierr;
12547c6ae99SBarry Smith   PetscInt       m,n = 0,x = 0,y = 0;
12647c6ae99SBarry Smith   PetscMPIInt    size,csize,rank;
12747c6ae99SBarry Smith 
12847c6ae99SBarry Smith   PetscFunctionBegin;
12947c6ae99SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
13047c6ae99SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
13147c6ae99SBarry Smith 
13247c6ae99SBarry Smith   csize = 4*size;
13347c6ae99SBarry Smith   do {
13447c6ae99SBarry Smith     if (csize % 4) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Cannot split communicator of size %d tried %d %D %D",size,csize,x,y);
13547c6ae99SBarry Smith     csize   = csize/4;
13647c6ae99SBarry Smith 
13747c6ae99SBarry Smith     m = (PetscInt)(0.5 + sqrt(((double)M)*((double)csize)/((double)N)));
13847c6ae99SBarry Smith     if (!m) m = 1;
13947c6ae99SBarry Smith     while (m > 0) {
14047c6ae99SBarry Smith       n = csize/m;
14147c6ae99SBarry Smith       if (m*n == csize) break;
14247c6ae99SBarry Smith       m--;
14347c6ae99SBarry Smith     }
14447c6ae99SBarry Smith     if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;}
14547c6ae99SBarry Smith 
14647c6ae99SBarry Smith     x = M/m + ((M % m) > ((csize-1) % m));
14747c6ae99SBarry Smith     y = (N + (csize-1)/m)/n;
14847c6ae99SBarry Smith   } while ((x < 4 || y < 4) && csize > 1);
14947c6ae99SBarry Smith   if (size != csize) {
15047c6ae99SBarry Smith     MPI_Group    entire_group,sub_group;
15147c6ae99SBarry Smith     PetscMPIInt  i,*groupies;
15247c6ae99SBarry Smith 
15347c6ae99SBarry Smith     ierr = MPI_Comm_group(comm,&entire_group);CHKERRQ(ierr);
15447c6ae99SBarry Smith     ierr = PetscMalloc(csize*sizeof(PetscInt),&groupies);CHKERRQ(ierr);
15547c6ae99SBarry Smith     for (i=0; i<csize; i++) {
15647c6ae99SBarry Smith       groupies[i] = (rank/csize)*csize + i;
15747c6ae99SBarry Smith     }
15847c6ae99SBarry Smith     ierr = MPI_Group_incl(entire_group,csize,groupies,&sub_group);CHKERRQ(ierr);
15947c6ae99SBarry Smith     ierr = PetscFree(groupies);CHKERRQ(ierr);
16047c6ae99SBarry Smith     ierr = MPI_Comm_create(comm,sub_group,outcomm);CHKERRQ(ierr);
16147c6ae99SBarry Smith     ierr = MPI_Group_free(&entire_group);CHKERRQ(ierr);
16247c6ae99SBarry Smith     ierr = MPI_Group_free(&sub_group);CHKERRQ(ierr);
163aa219208SBarry Smith     ierr = PetscInfo1(0,"DMDASplitComm2d:Creating redundant coarse problems of size %d\n",csize);CHKERRQ(ierr);
16447c6ae99SBarry Smith   } else {
16547c6ae99SBarry Smith     *outcomm = comm;
16647c6ae99SBarry Smith   }
16747c6ae99SBarry Smith   PetscFunctionReturn(0);
16847c6ae99SBarry Smith }
16947c6ae99SBarry Smith 
17047c6ae99SBarry Smith #undef __FUNCT__
171644e2e5bSBarry Smith #define __FUNCT__ "DMDAFunction"
172644e2e5bSBarry Smith static PetscErrorCode DMDAFunction(DM dm,Vec x,Vec F)
173644e2e5bSBarry Smith {
174644e2e5bSBarry Smith   PetscErrorCode ierr;
175644e2e5bSBarry Smith   Vec            localX;
176644e2e5bSBarry Smith 
177644e2e5bSBarry Smith   PetscFunctionBegin;
178644e2e5bSBarry Smith   ierr = DMGetLocalVector(dm,&localX);CHKERRQ(ierr);
179644e2e5bSBarry Smith   ierr = DMGlobalToLocalBegin(dm,x,INSERT_VALUES,localX);CHKERRQ(ierr);
180644e2e5bSBarry Smith   ierr = DMGlobalToLocalEnd(dm,x,INSERT_VALUES,localX);CHKERRQ(ierr);
181644e2e5bSBarry Smith   ierr = DMDAFormFunction1(dm,localX,F,dm->ctx);CHKERRQ(ierr);
182644e2e5bSBarry Smith   ierr = DMRestoreLocalVector(dm,&localX);CHKERRQ(ierr);
183644e2e5bSBarry Smith   PetscFunctionReturn(0);
184644e2e5bSBarry Smith }
185644e2e5bSBarry Smith 
186644e2e5bSBarry Smith #undef __FUNCT__
187aa219208SBarry Smith #define __FUNCT__ "DMDASetLocalFunction"
18847c6ae99SBarry Smith /*@C
189aa219208SBarry Smith        DMDASetLocalFunction - Caches in a DM a local function.
19047c6ae99SBarry Smith 
191aa219208SBarry Smith    Logically Collective on DMDA
19247c6ae99SBarry Smith 
19347c6ae99SBarry Smith    Input Parameter:
19447c6ae99SBarry Smith +  da - initial distributed array
19547c6ae99SBarry Smith -  lf - the local function
19647c6ae99SBarry Smith 
19747c6ae99SBarry Smith    Level: intermediate
19847c6ae99SBarry Smith 
19947c6ae99SBarry Smith    Notes: The routine SNESDAFormFunction() uses this the cached function to evaluate the user provided function.
20047c6ae99SBarry Smith 
20147c6ae99SBarry Smith .keywords:  distributed array, refine
20247c6ae99SBarry Smith 
203aa219208SBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalFunctioni()
20447c6ae99SBarry Smith @*/
2057087cfbeSBarry Smith PetscErrorCode  DMDASetLocalFunction(DM da,DMDALocalFunction1 lf)
20647c6ae99SBarry Smith {
207644e2e5bSBarry Smith   PetscErrorCode ierr;
20847c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
209644e2e5bSBarry Smith 
21047c6ae99SBarry Smith   PetscFunctionBegin;
21147c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
212644e2e5bSBarry Smith   ierr = DMSetFunction(da,DMDAFunction);CHKERRQ(ierr);
21347c6ae99SBarry Smith   dd->lf       = lf;
21447c6ae99SBarry Smith   PetscFunctionReturn(0);
21547c6ae99SBarry Smith }
21647c6ae99SBarry Smith 
21747c6ae99SBarry Smith #undef __FUNCT__
218aa219208SBarry Smith #define __FUNCT__ "DMDASetLocalFunctioni"
21947c6ae99SBarry Smith /*@C
220aa219208SBarry Smith        DMDASetLocalFunctioni - Caches in a DM a local function that evaluates a single component
22147c6ae99SBarry Smith 
222aa219208SBarry Smith    Logically Collective on DMDA
22347c6ae99SBarry Smith 
22447c6ae99SBarry Smith    Input Parameter:
22547c6ae99SBarry Smith +  da - initial distributed array
22647c6ae99SBarry Smith -  lfi - the local function
22747c6ae99SBarry Smith 
22847c6ae99SBarry Smith    Level: intermediate
22947c6ae99SBarry Smith 
23047c6ae99SBarry Smith .keywords:  distributed array, refine
23147c6ae99SBarry Smith 
232aa219208SBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalFunction()
23347c6ae99SBarry Smith @*/
2347087cfbeSBarry Smith PetscErrorCode  DMDASetLocalFunctioni(DM da,PetscErrorCode (*lfi)(DMDALocalInfo*,MatStencil*,void*,PetscScalar*,void*))
23547c6ae99SBarry Smith {
23647c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
23747c6ae99SBarry Smith   PetscFunctionBegin;
23847c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
23947c6ae99SBarry Smith   dd->lfi = lfi;
24047c6ae99SBarry Smith   PetscFunctionReturn(0);
24147c6ae99SBarry Smith }
24247c6ae99SBarry Smith 
24347c6ae99SBarry Smith #undef __FUNCT__
244aa219208SBarry Smith #define __FUNCT__ "DMDASetLocalFunctionib"
24547c6ae99SBarry Smith /*@C
246aa219208SBarry Smith        DMDASetLocalFunctionib - Caches in a DM a block local function that evaluates a single component
24747c6ae99SBarry Smith 
248aa219208SBarry Smith    Logically Collective on DMDA
24947c6ae99SBarry Smith 
25047c6ae99SBarry Smith    Input Parameter:
25147c6ae99SBarry Smith +  da - initial distributed array
25247c6ae99SBarry Smith -  lfi - the local function
25347c6ae99SBarry Smith 
25447c6ae99SBarry Smith    Level: intermediate
25547c6ae99SBarry Smith 
25647c6ae99SBarry Smith .keywords:  distributed array, refine
25747c6ae99SBarry Smith 
258aa219208SBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalFunction()
25947c6ae99SBarry Smith @*/
2607087cfbeSBarry Smith PetscErrorCode  DMDASetLocalFunctionib(DM da,PetscErrorCode (*lfi)(DMDALocalInfo*,MatStencil*,void*,PetscScalar*,void*))
26147c6ae99SBarry Smith {
26247c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
26347c6ae99SBarry Smith   PetscFunctionBegin;
26447c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
26547c6ae99SBarry Smith   dd->lfib = lfi;
26647c6ae99SBarry Smith   PetscFunctionReturn(0);
26747c6ae99SBarry Smith }
26847c6ae99SBarry Smith 
26947c6ae99SBarry Smith #undef __FUNCT__
270aa219208SBarry Smith #define __FUNCT__ "DMDASetLocalAdicFunction_Private"
271aa219208SBarry Smith PetscErrorCode DMDASetLocalAdicFunction_Private(DM da,DMDALocalFunction1 ad_lf)
27247c6ae99SBarry Smith {
27347c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
27447c6ae99SBarry Smith   PetscFunctionBegin;
27547c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
27647c6ae99SBarry Smith   dd->adic_lf = ad_lf;
27747c6ae99SBarry Smith   PetscFunctionReturn(0);
27847c6ae99SBarry Smith }
27947c6ae99SBarry Smith 
28047c6ae99SBarry Smith /*MC
281aa219208SBarry Smith        DMDASetLocalAdicFunctioni - Caches in a DM a local functioni computed by ADIC/ADIFOR
28247c6ae99SBarry Smith 
28347c6ae99SBarry Smith    Synopsis:
284aa219208SBarry Smith    PetscErrorCode DMDASetLocalAdicFunctioni(DM da,PetscInt (ad_lf*)(DMDALocalInfo*,MatStencil*,void*,void*,void*)
28547c6ae99SBarry Smith 
286aa219208SBarry Smith    Logically Collective on DMDA
28747c6ae99SBarry Smith 
28847c6ae99SBarry Smith    Input Parameter:
28947c6ae99SBarry Smith +  da - initial distributed array
29047c6ae99SBarry Smith -  ad_lfi - the local function as computed by ADIC/ADIFOR
29147c6ae99SBarry Smith 
29247c6ae99SBarry Smith    Level: intermediate
29347c6ae99SBarry Smith 
29447c6ae99SBarry Smith .keywords:  distributed array, refine
29547c6ae99SBarry Smith 
296aa219208SBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalFunction(),
297aa219208SBarry Smith           DMDASetLocalJacobian(), DMDASetLocalFunctioni()
29847c6ae99SBarry Smith M*/
29947c6ae99SBarry Smith 
30047c6ae99SBarry Smith #undef __FUNCT__
301aa219208SBarry Smith #define __FUNCT__ "DMDASetLocalAdicFunctioni_Private"
302aa219208SBarry Smith PetscErrorCode DMDASetLocalAdicFunctioni_Private(DM da,PetscErrorCode (*ad_lfi)(DMDALocalInfo*,MatStencil*,void*,void*,void*))
30347c6ae99SBarry Smith {
30447c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
30547c6ae99SBarry Smith   PetscFunctionBegin;
30647c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
30747c6ae99SBarry Smith   dd->adic_lfi = ad_lfi;
30847c6ae99SBarry Smith   PetscFunctionReturn(0);
30947c6ae99SBarry Smith }
31047c6ae99SBarry Smith 
31147c6ae99SBarry Smith /*MC
312aa219208SBarry Smith        DMDASetLocalAdicMFFunctioni - Caches in a DM a local functioni computed by ADIC/ADIFOR
31347c6ae99SBarry Smith 
31447c6ae99SBarry Smith    Synopsis:
315aa219208SBarry Smith    PetscErrorCode  DMDASetLocalAdicFunctioni(DM da,int (ad_lf*)(DMDALocalInfo*,MatStencil*,void*,void*,void*)
31647c6ae99SBarry Smith 
317aa219208SBarry Smith    Logically Collective on DMDA
31847c6ae99SBarry Smith 
31947c6ae99SBarry Smith    Input Parameter:
32047c6ae99SBarry Smith +  da - initial distributed array
32147c6ae99SBarry Smith -  admf_lfi - the local matrix-free function as computed by ADIC/ADIFOR
32247c6ae99SBarry Smith 
32347c6ae99SBarry Smith    Level: intermediate
32447c6ae99SBarry Smith 
32547c6ae99SBarry Smith .keywords:  distributed array, refine
32647c6ae99SBarry Smith 
327aa219208SBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalFunction(),
328aa219208SBarry Smith           DMDASetLocalJacobian(), DMDASetLocalFunctioni()
32947c6ae99SBarry Smith M*/
33047c6ae99SBarry Smith 
33147c6ae99SBarry Smith #undef __FUNCT__
332aa219208SBarry Smith #define __FUNCT__ "DMDASetLocalAdicMFFunctioni_Private"
333aa219208SBarry Smith PetscErrorCode DMDASetLocalAdicMFFunctioni_Private(DM da,PetscErrorCode (*admf_lfi)(DMDALocalInfo*,MatStencil*,void*,void*,void*))
33447c6ae99SBarry Smith {
33547c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
33647c6ae99SBarry Smith   PetscFunctionBegin;
33747c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
33847c6ae99SBarry Smith   dd->adicmf_lfi = admf_lfi;
33947c6ae99SBarry Smith   PetscFunctionReturn(0);
34047c6ae99SBarry Smith }
34147c6ae99SBarry Smith 
34247c6ae99SBarry Smith /*MC
343aa219208SBarry Smith        DMDASetLocalAdicFunctionib - Caches in a DM a block local functioni computed by ADIC/ADIFOR
34447c6ae99SBarry Smith 
34547c6ae99SBarry Smith    Synopsis:
346aa219208SBarry Smith    PetscErrorCode DMDASetLocalAdicFunctionib(DM da,PetscInt (ad_lf*)(DMDALocalInfo*,MatStencil*,void*,void*,void*)
34747c6ae99SBarry Smith 
348aa219208SBarry Smith    Logically Collective on DMDA
34947c6ae99SBarry Smith 
35047c6ae99SBarry Smith    Input Parameter:
35147c6ae99SBarry Smith +  da - initial distributed array
35247c6ae99SBarry Smith -  ad_lfi - the local function as computed by ADIC/ADIFOR
35347c6ae99SBarry Smith 
35447c6ae99SBarry Smith    Level: intermediate
35547c6ae99SBarry Smith 
35647c6ae99SBarry Smith .keywords:  distributed array, refine
35747c6ae99SBarry Smith 
358aa219208SBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalFunction(),
359aa219208SBarry Smith           DMDASetLocalJacobian(), DMDASetLocalFunctionib()
36047c6ae99SBarry Smith M*/
36147c6ae99SBarry Smith 
36247c6ae99SBarry Smith #undef __FUNCT__
363aa219208SBarry Smith #define __FUNCT__ "DMDASetLocalAdicFunctionib_Private"
364aa219208SBarry Smith PetscErrorCode DMDASetLocalAdicFunctionib_Private(DM da,PetscErrorCode (*ad_lfi)(DMDALocalInfo*,MatStencil*,void*,void*,void*))
36547c6ae99SBarry Smith {
36647c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
36747c6ae99SBarry Smith   PetscFunctionBegin;
36847c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
36947c6ae99SBarry Smith   dd->adic_lfib = ad_lfi;
37047c6ae99SBarry Smith   PetscFunctionReturn(0);
37147c6ae99SBarry Smith }
37247c6ae99SBarry Smith 
37347c6ae99SBarry Smith /*MC
374aa219208SBarry Smith        DMDASetLocalAdicMFFunctionib - Caches in a DM a block local functioni computed by ADIC/ADIFOR
37547c6ae99SBarry Smith 
37647c6ae99SBarry Smith    Synopsis:
377aa219208SBarry Smith    PetscErrorCode  DMDASetLocalAdicFunctionib(DM da,int (ad_lf*)(DMDALocalInfo*,MatStencil*,void*,void*,void*)
37847c6ae99SBarry Smith 
379aa219208SBarry Smith    Logically Collective on DMDA
38047c6ae99SBarry Smith 
38147c6ae99SBarry Smith    Input Parameter:
38247c6ae99SBarry Smith +  da - initial distributed array
38347c6ae99SBarry Smith -  admf_lfi - the local matrix-free function as computed by ADIC/ADIFOR
38447c6ae99SBarry Smith 
38547c6ae99SBarry Smith    Level: intermediate
38647c6ae99SBarry Smith 
38747c6ae99SBarry Smith .keywords:  distributed array, refine
38847c6ae99SBarry Smith 
389aa219208SBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalFunction(),
390aa219208SBarry Smith           DMDASetLocalJacobian(), DMDASetLocalFunctionib()
39147c6ae99SBarry Smith M*/
39247c6ae99SBarry Smith 
39347c6ae99SBarry Smith #undef __FUNCT__
394aa219208SBarry Smith #define __FUNCT__ "DMDASetLocalAdicMFFunctionib_Private"
395aa219208SBarry Smith PetscErrorCode DMDASetLocalAdicMFFunctionib_Private(DM da,PetscErrorCode (*admf_lfi)(DMDALocalInfo*,MatStencil*,void*,void*,void*))
39647c6ae99SBarry Smith {
39747c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
39847c6ae99SBarry Smith   PetscFunctionBegin;
39947c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
40047c6ae99SBarry Smith   dd->adicmf_lfib = admf_lfi;
40147c6ae99SBarry Smith   PetscFunctionReturn(0);
40247c6ae99SBarry Smith }
40347c6ae99SBarry Smith 
40447c6ae99SBarry Smith /*MC
405aa219208SBarry Smith        DMDASetLocalAdicMFFunction - Caches in a DM a local function computed by ADIC/ADIFOR
40647c6ae99SBarry Smith 
40747c6ae99SBarry Smith    Synopsis:
408aa219208SBarry Smith    PetscErrorCode DMDASetLocalAdicMFFunction(DM da,DMDALocalFunction1 ad_lf)
40947c6ae99SBarry Smith 
410aa219208SBarry Smith    Logically Collective on DMDA
41147c6ae99SBarry Smith 
41247c6ae99SBarry Smith    Input Parameter:
41347c6ae99SBarry Smith +  da - initial distributed array
41447c6ae99SBarry Smith -  ad_lf - the local function as computed by ADIC/ADIFOR
41547c6ae99SBarry Smith 
41647c6ae99SBarry Smith    Level: intermediate
41747c6ae99SBarry Smith 
41847c6ae99SBarry Smith .keywords:  distributed array, refine
41947c6ae99SBarry Smith 
420aa219208SBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalFunction(),
421aa219208SBarry Smith           DMDASetLocalJacobian()
42247c6ae99SBarry Smith M*/
42347c6ae99SBarry Smith 
42447c6ae99SBarry Smith #undef __FUNCT__
425aa219208SBarry Smith #define __FUNCT__ "DMDASetLocalAdicMFFunction_Private"
426aa219208SBarry Smith PetscErrorCode DMDASetLocalAdicMFFunction_Private(DM da,DMDALocalFunction1 ad_lf)
42747c6ae99SBarry Smith {
42847c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
42947c6ae99SBarry Smith   PetscFunctionBegin;
43047c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
43147c6ae99SBarry Smith   dd->adicmf_lf = ad_lf;
43247c6ae99SBarry Smith   PetscFunctionReturn(0);
43347c6ae99SBarry Smith }
43447c6ae99SBarry Smith 
435644e2e5bSBarry Smith #undef __FUNCT__
4362533e041SBarry Smith #define __FUNCT__ "DMDAJacobianDefaultLocal"
4372533e041SBarry Smith PetscErrorCode DMDAJacobianLocal(DM dm,Vec x,Mat A,Mat B, MatStructure *str)
4382533e041SBarry Smith {
4392533e041SBarry Smith   PetscErrorCode ierr;
4402533e041SBarry Smith   Vec            localX;
4412533e041SBarry Smith 
4422533e041SBarry Smith   PetscFunctionBegin;
4432533e041SBarry Smith   ierr = DMGetLocalVector(dm,&localX);CHKERRQ(ierr);
4442533e041SBarry Smith   ierr = DMGlobalToLocalBegin(dm,x,INSERT_VALUES,localX);CHKERRQ(ierr);
4452533e041SBarry Smith   ierr = DMGlobalToLocalEnd(dm,x,INSERT_VALUES,localX);CHKERRQ(ierr);
4462533e041SBarry Smith   ierr = MatFDColoringApply(B,dm->fd,localX,str,dm);CHKERRQ(ierr);
4472533e041SBarry Smith   ierr = DMRestoreLocalVector(dm,&localX);CHKERRQ(ierr);
4482533e041SBarry Smith   /* Assemble true Jacobian; if it is different */
4492533e041SBarry Smith   if (A != B) {
4502533e041SBarry Smith     ierr  = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
4512533e041SBarry Smith     ierr  = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
4522533e041SBarry Smith   }
4532533e041SBarry Smith   ierr  = MatSetOption(B,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
4542533e041SBarry Smith   *str = SAME_NONZERO_PATTERN;
4552533e041SBarry Smith   PetscFunctionReturn(0);
4562533e041SBarry Smith }
4572533e041SBarry Smith 
4582533e041SBarry Smith 
4592533e041SBarry Smith #undef __FUNCT__
460644e2e5bSBarry Smith #define __FUNCT__ "DMDAJacobian"
461644e2e5bSBarry Smith static PetscErrorCode DMDAJacobian(DM dm,Vec x,Mat A,Mat B, MatStructure *str)
462644e2e5bSBarry Smith {
463644e2e5bSBarry Smith   PetscErrorCode ierr;
464644e2e5bSBarry Smith   Vec            localX;
465644e2e5bSBarry Smith 
466644e2e5bSBarry Smith   PetscFunctionBegin;
467644e2e5bSBarry Smith   ierr = DMGetLocalVector(dm,&localX);CHKERRQ(ierr);
468644e2e5bSBarry Smith   ierr = DMGlobalToLocalBegin(dm,x,INSERT_VALUES,localX);CHKERRQ(ierr);
469644e2e5bSBarry Smith   ierr = DMGlobalToLocalEnd(dm,x,INSERT_VALUES,localX);CHKERRQ(ierr);
470644e2e5bSBarry Smith   ierr = DMDAComputeJacobian1(dm,localX,B,dm->ctx);CHKERRQ(ierr);
471644e2e5bSBarry Smith   ierr = DMRestoreLocalVector(dm,&localX);CHKERRQ(ierr);
472644e2e5bSBarry Smith   /* Assemble true Jacobian; if it is different */
473644e2e5bSBarry Smith   if (A != B) {
474644e2e5bSBarry Smith     ierr  = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
475644e2e5bSBarry Smith     ierr  = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
476644e2e5bSBarry Smith   }
477644e2e5bSBarry Smith   ierr  = MatSetOption(B,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
478644e2e5bSBarry Smith   *str = SAME_NONZERO_PATTERN;
479644e2e5bSBarry Smith   PetscFunctionReturn(0);
480644e2e5bSBarry Smith }
481644e2e5bSBarry Smith 
48247c6ae99SBarry Smith /*@C
483aa219208SBarry Smith        DMDASetLocalJacobian - Caches in a DM a local Jacobian computation function
48447c6ae99SBarry Smith 
485aa219208SBarry Smith    Logically Collective on DMDA
48647c6ae99SBarry Smith 
48747c6ae99SBarry Smith 
48847c6ae99SBarry Smith    Input Parameter:
48947c6ae99SBarry Smith +  da - initial distributed array
49047c6ae99SBarry Smith -  lj - the local Jacobian
49147c6ae99SBarry Smith 
49247c6ae99SBarry Smith    Level: intermediate
49347c6ae99SBarry Smith 
49447c6ae99SBarry Smith    Notes: The routine SNESDAFormFunction() uses this the cached function to evaluate the user provided function.
49547c6ae99SBarry Smith 
49647c6ae99SBarry Smith .keywords:  distributed array, refine
49747c6ae99SBarry Smith 
498aa219208SBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalFunction()
49947c6ae99SBarry Smith @*/
50047c6ae99SBarry Smith #undef __FUNCT__
501aa219208SBarry Smith #define __FUNCT__ "DMDASetLocalJacobian"
5027087cfbeSBarry Smith PetscErrorCode  DMDASetLocalJacobian(DM da,DMDALocalFunction1 lj)
50347c6ae99SBarry Smith {
504644e2e5bSBarry Smith   PetscErrorCode ierr;
50547c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
506644e2e5bSBarry Smith 
50747c6ae99SBarry Smith   PetscFunctionBegin;
50847c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
509644e2e5bSBarry Smith   ierr = DMSetJacobian(da,DMDAJacobian);CHKERRQ(ierr);
51047c6ae99SBarry Smith   dd->lj    = lj;
51147c6ae99SBarry Smith   PetscFunctionReturn(0);
51247c6ae99SBarry Smith }
51347c6ae99SBarry Smith 
51447c6ae99SBarry Smith #undef __FUNCT__
515aa219208SBarry Smith #define __FUNCT__ "DMDAGetLocalFunction"
51647c6ae99SBarry Smith /*@C
517aa219208SBarry Smith        DMDAGetLocalFunction - Gets from a DM a local function and its ADIC/ADIFOR Jacobian
51847c6ae99SBarry Smith 
51947c6ae99SBarry Smith    Note Collective
52047c6ae99SBarry Smith 
52147c6ae99SBarry Smith    Input Parameter:
52247c6ae99SBarry Smith .  da - initial distributed array
52347c6ae99SBarry Smith 
52447c6ae99SBarry Smith    Output Parameter:
52547c6ae99SBarry Smith .  lf - the local function
52647c6ae99SBarry Smith 
52747c6ae99SBarry Smith    Level: intermediate
52847c6ae99SBarry Smith 
52947c6ae99SBarry Smith .keywords:  distributed array, refine
53047c6ae99SBarry Smith 
531aa219208SBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalJacobian(), DMDASetLocalFunction()
53247c6ae99SBarry Smith @*/
5337087cfbeSBarry Smith PetscErrorCode  DMDAGetLocalFunction(DM da,DMDALocalFunction1 *lf)
53447c6ae99SBarry Smith {
53547c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
53647c6ae99SBarry Smith   PetscFunctionBegin;
53747c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
53847c6ae99SBarry Smith   if (lf) *lf = dd->lf;
53947c6ae99SBarry Smith   PetscFunctionReturn(0);
54047c6ae99SBarry Smith }
54147c6ae99SBarry Smith 
54247c6ae99SBarry Smith #undef __FUNCT__
543aa219208SBarry Smith #define __FUNCT__ "DMDAGetLocalJacobian"
54447c6ae99SBarry Smith /*@C
545aa219208SBarry Smith        DMDAGetLocalJacobian - Gets from a DM a local jacobian
54647c6ae99SBarry Smith 
54747c6ae99SBarry Smith    Not Collective
54847c6ae99SBarry Smith 
54947c6ae99SBarry Smith    Input Parameter:
55047c6ae99SBarry Smith .  da - initial distributed array
55147c6ae99SBarry Smith 
55247c6ae99SBarry Smith    Output Parameter:
55347c6ae99SBarry Smith .  lj - the local jacobian
55447c6ae99SBarry Smith 
55547c6ae99SBarry Smith    Level: intermediate
55647c6ae99SBarry Smith 
55747c6ae99SBarry Smith .keywords:  distributed array, refine
55847c6ae99SBarry Smith 
559aa219208SBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalJacobian()
56047c6ae99SBarry Smith @*/
5617087cfbeSBarry Smith PetscErrorCode  DMDAGetLocalJacobian(DM da,DMDALocalFunction1 *lj)
56247c6ae99SBarry Smith {
56347c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
56447c6ae99SBarry Smith   PetscFunctionBegin;
56547c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
56647c6ae99SBarry Smith   if (lj) *lj = dd->lj;
56747c6ae99SBarry Smith   PetscFunctionReturn(0);
56847c6ae99SBarry Smith }
56947c6ae99SBarry Smith 
57047c6ae99SBarry Smith #undef __FUNCT__
571aa219208SBarry Smith #define __FUNCT__ "DMDAFormFunction"
57247c6ae99SBarry Smith /*@
573aa219208SBarry Smith     DMDAFormFunction - Evaluates a user provided function on each processor that
574aa219208SBarry Smith         share a DMDA
57547c6ae99SBarry Smith 
57647c6ae99SBarry Smith    Input Parameters:
5779a42bb27SBarry Smith +    da - the DM that defines the grid
57847c6ae99SBarry Smith .    vu - input vector
57947c6ae99SBarry Smith .    vfu - output vector
58047c6ae99SBarry Smith -    w - any user data
58147c6ae99SBarry Smith 
58247c6ae99SBarry Smith     Notes: Does NOT do ghost updates on vu upon entry
58347c6ae99SBarry Smith 
584aa219208SBarry Smith            This should eventually replace DMDAFormFunction1
58547c6ae99SBarry Smith 
58647c6ae99SBarry Smith     Level: advanced
58747c6ae99SBarry Smith 
588aa219208SBarry Smith .seealso: DMDAComputeJacobian1WithAdic()
58947c6ae99SBarry Smith 
59047c6ae99SBarry Smith @*/
5917087cfbeSBarry Smith PetscErrorCode  DMDAFormFunction(DM da,PetscErrorCode (*lf)(void),Vec vu,Vec vfu,void *w)
59247c6ae99SBarry Smith {
59347c6ae99SBarry Smith   PetscErrorCode ierr;
59447c6ae99SBarry Smith   void           *u,*fu;
595aa219208SBarry Smith   DMDALocalInfo  info;
596aa219208SBarry Smith   PetscErrorCode (*f)(DMDALocalInfo*,void*,void*,void*) = (PetscErrorCode (*)(DMDALocalInfo*,void*,void*,void*))lf;
59747c6ae99SBarry Smith 
59847c6ae99SBarry Smith   PetscFunctionBegin;
599aa219208SBarry Smith   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
600aa219208SBarry Smith   ierr = DMDAVecGetArray(da,vu,&u);CHKERRQ(ierr);
601aa219208SBarry Smith   ierr = DMDAVecGetArray(da,vfu,&fu);CHKERRQ(ierr);
60247c6ae99SBarry Smith 
60339d508bbSBarry Smith   ierr = (*f)(&info,u,fu,w);CHKERRQ(ierr);
60447c6ae99SBarry Smith 
605aa219208SBarry Smith   ierr = DMDAVecRestoreArray(da,vu,&u);CHKERRQ(ierr);
606aa219208SBarry Smith   ierr = DMDAVecRestoreArray(da,vfu,&fu);CHKERRQ(ierr);
60747c6ae99SBarry Smith   PetscFunctionReturn(0);
60847c6ae99SBarry Smith }
60947c6ae99SBarry Smith 
61047c6ae99SBarry Smith #undef __FUNCT__
611aa219208SBarry Smith #define __FUNCT__ "DMDAFormFunctionLocal"
61247c6ae99SBarry Smith /*@C
613aa219208SBarry Smith    DMDAFormFunctionLocal - This is a universal function evaluation routine for
6149a42bb27SBarry Smith    a local DM function.
61547c6ae99SBarry Smith 
616aa219208SBarry Smith    Collective on DMDA
61747c6ae99SBarry Smith 
61847c6ae99SBarry Smith    Input Parameters:
6199a42bb27SBarry Smith +  da - the DM context
62047c6ae99SBarry Smith .  func - The local function
62147c6ae99SBarry Smith .  X - input vector
62247c6ae99SBarry Smith .  F - function vector
62347c6ae99SBarry Smith -  ctx - A user context
62447c6ae99SBarry Smith 
62547c6ae99SBarry Smith    Level: intermediate
62647c6ae99SBarry Smith 
627aa219208SBarry Smith .seealso: DMDASetLocalFunction(), DMDASetLocalJacobian(), DMDASetLocalAdicFunction(), DMDASetLocalAdicMFFunction(),
62847c6ae99SBarry Smith           SNESSetFunction(), SNESSetJacobian()
62947c6ae99SBarry Smith 
63047c6ae99SBarry Smith @*/
6317087cfbeSBarry Smith PetscErrorCode  DMDAFormFunctionLocal(DM da, DMDALocalFunction1 func, Vec X, Vec F, void *ctx)
63247c6ae99SBarry Smith {
63347c6ae99SBarry Smith   Vec            localX;
634aa219208SBarry Smith   DMDALocalInfo  info;
63547c6ae99SBarry Smith   void           *u;
63647c6ae99SBarry Smith   void           *fu;
63747c6ae99SBarry Smith   PetscErrorCode ierr;
63847c6ae99SBarry Smith 
63947c6ae99SBarry Smith   PetscFunctionBegin;
6409a42bb27SBarry Smith   ierr = DMGetLocalVector(da,&localX);CHKERRQ(ierr);
64147c6ae99SBarry Smith   /*
64247c6ae99SBarry Smith      Scatter ghost points to local vector, using the 2-step process
6439a42bb27SBarry Smith         DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
64447c6ae99SBarry Smith   */
6459a42bb27SBarry Smith   ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
6469a42bb27SBarry Smith   ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
647aa219208SBarry Smith   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
648aa219208SBarry Smith   ierr = DMDAVecGetArray(da,localX,&u);CHKERRQ(ierr);
649aa219208SBarry Smith   ierr = DMDAVecGetArray(da,F,&fu);CHKERRQ(ierr);
65039d508bbSBarry Smith   ierr = (*func)(&info,u,fu,ctx);CHKERRQ(ierr);
651aa219208SBarry Smith   ierr = DMDAVecRestoreArray(da,localX,&u);CHKERRQ(ierr);
652aa219208SBarry Smith   ierr = DMDAVecRestoreArray(da,F,&fu);CHKERRQ(ierr);
6539a42bb27SBarry Smith   ierr = DMRestoreLocalVector(da,&localX);CHKERRQ(ierr);
65447c6ae99SBarry Smith   PetscFunctionReturn(0);
65547c6ae99SBarry Smith }
65647c6ae99SBarry Smith 
65747c6ae99SBarry Smith #undef __FUNCT__
658aa219208SBarry Smith #define __FUNCT__ "DMDAFormFunctionLocalGhost"
65947c6ae99SBarry Smith /*@C
660aa219208SBarry Smith    DMDAFormFunctionLocalGhost - This is a universal function evaluation routine for
6619a42bb27SBarry Smith    a local DM function, but the ghost values of the output are communicated and added.
66247c6ae99SBarry Smith 
663aa219208SBarry Smith    Collective on DMDA
66447c6ae99SBarry Smith 
66547c6ae99SBarry Smith    Input Parameters:
6669a42bb27SBarry Smith +  da - the DM context
66747c6ae99SBarry Smith .  func - The local function
66847c6ae99SBarry Smith .  X - input vector
66947c6ae99SBarry Smith .  F - function vector
67047c6ae99SBarry Smith -  ctx - A user context
67147c6ae99SBarry Smith 
67247c6ae99SBarry Smith    Level: intermediate
67347c6ae99SBarry Smith 
674aa219208SBarry Smith .seealso: DMDASetLocalFunction(), DMDASetLocalJacobian(), DMDASetLocalAdicFunction(), DMDASetLocalAdicMFFunction(),
67547c6ae99SBarry Smith           SNESSetFunction(), SNESSetJacobian()
67647c6ae99SBarry Smith 
67747c6ae99SBarry Smith @*/
6787087cfbeSBarry Smith PetscErrorCode  DMDAFormFunctionLocalGhost(DM da, DMDALocalFunction1 func, Vec X, Vec F, void *ctx)
67947c6ae99SBarry Smith {
68047c6ae99SBarry Smith   Vec            localX, localF;
681aa219208SBarry Smith   DMDALocalInfo  info;
68247c6ae99SBarry Smith   void           *u;
68347c6ae99SBarry Smith   void           *fu;
68447c6ae99SBarry Smith   PetscErrorCode ierr;
68547c6ae99SBarry Smith 
68647c6ae99SBarry Smith   PetscFunctionBegin;
6879a42bb27SBarry Smith   ierr = DMGetLocalVector(da,&localX);CHKERRQ(ierr);
6889a42bb27SBarry Smith   ierr = DMGetLocalVector(da,&localF);CHKERRQ(ierr);
68947c6ae99SBarry Smith   /*
69047c6ae99SBarry Smith      Scatter ghost points to local vector, using the 2-step process
6919a42bb27SBarry Smith         DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
69247c6ae99SBarry Smith   */
6939a42bb27SBarry Smith   ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
6949a42bb27SBarry Smith   ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
69547c6ae99SBarry Smith   ierr = VecSet(F, 0.0);CHKERRQ(ierr);
69647c6ae99SBarry Smith   ierr = VecSet(localF, 0.0);CHKERRQ(ierr);
697aa219208SBarry Smith   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
698aa219208SBarry Smith   ierr = DMDAVecGetArray(da,localX,&u);CHKERRQ(ierr);
699aa219208SBarry Smith   ierr = DMDAVecGetArray(da,localF,&fu);CHKERRQ(ierr);
70039d508bbSBarry Smith   ierr = (*func)(&info,u,fu,ctx);CHKERRQ(ierr);
7019a42bb27SBarry Smith   ierr = DMLocalToGlobalBegin(da,localF,ADD_VALUES,F);CHKERRQ(ierr);
7029a42bb27SBarry Smith   ierr = DMLocalToGlobalEnd(da,localF,ADD_VALUES,F);CHKERRQ(ierr);
703aa219208SBarry Smith   ierr = DMDAVecRestoreArray(da,localX,&u);CHKERRQ(ierr);
704aa219208SBarry Smith   ierr = DMDAVecRestoreArray(da,localF,&fu);CHKERRQ(ierr);
7059a42bb27SBarry Smith   ierr = DMRestoreLocalVector(da,&localX);CHKERRQ(ierr);
7069a42bb27SBarry Smith   ierr = DMRestoreLocalVector(da,&localF);CHKERRQ(ierr);
70747c6ae99SBarry Smith   PetscFunctionReturn(0);
70847c6ae99SBarry Smith }
70947c6ae99SBarry Smith 
71047c6ae99SBarry Smith #undef __FUNCT__
711aa219208SBarry Smith #define __FUNCT__ "DMDAFormFunction1"
71247c6ae99SBarry Smith /*@
713aa219208SBarry Smith     DMDAFormFunction1 - Evaluates a user provided function on each processor that
714aa219208SBarry Smith         share a DMDA
71547c6ae99SBarry Smith 
71647c6ae99SBarry Smith    Input Parameters:
7179a42bb27SBarry Smith +    da - the DM that defines the grid
71847c6ae99SBarry Smith .    vu - input vector
71947c6ae99SBarry Smith .    vfu - output vector
72047c6ae99SBarry Smith -    w - any user data
72147c6ae99SBarry Smith 
72247c6ae99SBarry Smith     Notes: Does NOT do ghost updates on vu upon entry
72347c6ae99SBarry Smith 
72447c6ae99SBarry Smith     Level: advanced
72547c6ae99SBarry Smith 
726aa219208SBarry Smith .seealso: DMDAComputeJacobian1WithAdic()
72747c6ae99SBarry Smith 
72847c6ae99SBarry Smith @*/
7297087cfbeSBarry Smith PetscErrorCode  DMDAFormFunction1(DM da,Vec vu,Vec vfu,void *w)
73047c6ae99SBarry Smith {
73147c6ae99SBarry Smith   PetscErrorCode ierr;
73247c6ae99SBarry Smith   void           *u,*fu;
733aa219208SBarry Smith   DMDALocalInfo  info;
73447c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
73547c6ae99SBarry Smith 
73647c6ae99SBarry Smith   PetscFunctionBegin;
737aa219208SBarry Smith   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
738aa219208SBarry Smith   ierr = DMDAVecGetArray(da,vu,&u);CHKERRQ(ierr);
739aa219208SBarry Smith   ierr = DMDAVecGetArray(da,vfu,&fu);CHKERRQ(ierr);
74047c6ae99SBarry Smith 
74147c6ae99SBarry Smith   CHKMEMQ;
74239d508bbSBarry Smith   ierr = (*dd->lf)(&info,u,fu,w);CHKERRQ(ierr);
74347c6ae99SBarry Smith   CHKMEMQ;
74447c6ae99SBarry Smith 
745aa219208SBarry Smith   ierr = DMDAVecRestoreArray(da,vu,&u);CHKERRQ(ierr);
746aa219208SBarry Smith   ierr = DMDAVecRestoreArray(da,vfu,&fu);CHKERRQ(ierr);
74747c6ae99SBarry Smith   PetscFunctionReturn(0);
74847c6ae99SBarry Smith }
74947c6ae99SBarry Smith 
75047c6ae99SBarry Smith #undef __FUNCT__
751aa219208SBarry Smith #define __FUNCT__ "DMDAFormFunctioniTest1"
7527087cfbeSBarry Smith PetscErrorCode  DMDAFormFunctioniTest1(DM da,void *w)
75347c6ae99SBarry Smith {
75447c6ae99SBarry Smith   Vec            vu,fu,fui;
75547c6ae99SBarry Smith   PetscErrorCode ierr;
75647c6ae99SBarry Smith   PetscInt       i,n;
75747c6ae99SBarry Smith   PetscScalar    *ui;
75847c6ae99SBarry Smith   PetscRandom    rnd;
75947c6ae99SBarry Smith   PetscReal      norm;
76047c6ae99SBarry Smith 
76147c6ae99SBarry Smith   PetscFunctionBegin;
7629a42bb27SBarry Smith   ierr = DMGetLocalVector(da,&vu);CHKERRQ(ierr);
76347c6ae99SBarry Smith   ierr = PetscRandomCreate(PETSC_COMM_SELF,&rnd);CHKERRQ(ierr);
76447c6ae99SBarry Smith   ierr = PetscRandomSetFromOptions(rnd);CHKERRQ(ierr);
76547c6ae99SBarry Smith   ierr = VecSetRandom(vu,rnd);CHKERRQ(ierr);
766fcfd50ebSBarry Smith   ierr = PetscRandomDestroy(&rnd);CHKERRQ(ierr);
76747c6ae99SBarry Smith 
7689a42bb27SBarry Smith   ierr = DMGetGlobalVector(da,&fu);CHKERRQ(ierr);
7699a42bb27SBarry Smith   ierr = DMGetGlobalVector(da,&fui);CHKERRQ(ierr);
77047c6ae99SBarry Smith 
771aa219208SBarry Smith   ierr = DMDAFormFunction1(da,vu,fu,w);CHKERRQ(ierr);
77247c6ae99SBarry Smith 
77347c6ae99SBarry Smith   ierr = VecGetArray(fui,&ui);CHKERRQ(ierr);
77447c6ae99SBarry Smith   ierr = VecGetLocalSize(fui,&n);CHKERRQ(ierr);
77547c6ae99SBarry Smith   for (i=0; i<n; i++) {
776aa219208SBarry Smith     ierr = DMDAFormFunctioni1(da,i,vu,ui+i,w);CHKERRQ(ierr);
77747c6ae99SBarry Smith   }
77847c6ae99SBarry Smith   ierr = VecRestoreArray(fui,&ui);CHKERRQ(ierr);
77947c6ae99SBarry Smith 
78047c6ae99SBarry Smith   ierr = VecAXPY(fui,-1.0,fu);CHKERRQ(ierr);
78147c6ae99SBarry Smith   ierr = VecNorm(fui,NORM_2,&norm);CHKERRQ(ierr);
78247c6ae99SBarry Smith   ierr = PetscPrintf(((PetscObject)da)->comm,"Norm of difference in vectors %G\n",norm);CHKERRQ(ierr);
78347c6ae99SBarry Smith   ierr = VecView(fu,0);CHKERRQ(ierr);
78447c6ae99SBarry Smith   ierr = VecView(fui,0);CHKERRQ(ierr);
78547c6ae99SBarry Smith 
7869a42bb27SBarry Smith   ierr = DMRestoreLocalVector(da,&vu);CHKERRQ(ierr);
7879a42bb27SBarry Smith   ierr = DMRestoreGlobalVector(da,&fu);CHKERRQ(ierr);
7889a42bb27SBarry Smith   ierr = DMRestoreGlobalVector(da,&fui);CHKERRQ(ierr);
78947c6ae99SBarry Smith   PetscFunctionReturn(0);
79047c6ae99SBarry Smith }
79147c6ae99SBarry Smith 
79247c6ae99SBarry Smith #undef __FUNCT__
793aa219208SBarry Smith #define __FUNCT__ "DMDAFormFunctioni1"
79447c6ae99SBarry Smith /*@
795aa219208SBarry Smith     DMDAFormFunctioni1 - Evaluates a user provided point-wise function
79647c6ae99SBarry Smith 
79747c6ae99SBarry Smith    Input Parameters:
7989a42bb27SBarry Smith +    da - the DM that defines the grid
79947c6ae99SBarry Smith .    i - the component of the function we wish to compute (must be local)
80047c6ae99SBarry Smith .    vu - input vector
80147c6ae99SBarry Smith .    vfu - output value
80247c6ae99SBarry Smith -    w - any user data
80347c6ae99SBarry Smith 
80447c6ae99SBarry Smith     Notes: Does NOT do ghost updates on vu upon entry
80547c6ae99SBarry Smith 
80647c6ae99SBarry Smith     Level: advanced
80747c6ae99SBarry Smith 
808aa219208SBarry Smith .seealso: DMDAComputeJacobian1WithAdic()
80947c6ae99SBarry Smith 
81047c6ae99SBarry Smith @*/
8117087cfbeSBarry Smith PetscErrorCode  DMDAFormFunctioni1(DM da,PetscInt i,Vec vu,PetscScalar *vfu,void *w)
81247c6ae99SBarry Smith {
81347c6ae99SBarry Smith   PetscErrorCode ierr;
81447c6ae99SBarry Smith   void           *u;
815aa219208SBarry Smith   DMDALocalInfo  info;
81647c6ae99SBarry Smith   MatStencil     stencil;
81747c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
81847c6ae99SBarry Smith 
81947c6ae99SBarry Smith   PetscFunctionBegin;
82047c6ae99SBarry Smith 
821aa219208SBarry Smith   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
822aa219208SBarry Smith   ierr = DMDAVecGetArray(da,vu,&u);CHKERRQ(ierr);
82347c6ae99SBarry Smith 
82447c6ae99SBarry Smith   /* figure out stencil value from i */
82547c6ae99SBarry Smith   stencil.c = i % info.dof;
82647c6ae99SBarry Smith   stencil.i = (i % (info.xm*info.dof))/info.dof;
82747c6ae99SBarry Smith   stencil.j = (i % (info.xm*info.ym*info.dof))/(info.xm*info.dof);
82847c6ae99SBarry Smith   stencil.k = i/(info.xm*info.ym*info.dof);
82947c6ae99SBarry Smith 
83047c6ae99SBarry Smith   ierr = (*dd->lfi)(&info,&stencil,u,vfu,w);CHKERRQ(ierr);
83147c6ae99SBarry Smith 
832aa219208SBarry Smith   ierr = DMDAVecRestoreArray(da,vu,&u);CHKERRQ(ierr);
83347c6ae99SBarry Smith   PetscFunctionReturn(0);
83447c6ae99SBarry Smith }
83547c6ae99SBarry Smith 
83647c6ae99SBarry Smith #undef __FUNCT__
837aa219208SBarry Smith #define __FUNCT__ "DMDAFormFunctionib1"
83847c6ae99SBarry Smith /*@
839aa219208SBarry Smith     DMDAFormFunctionib1 - Evaluates a user provided point-block function
84047c6ae99SBarry Smith 
84147c6ae99SBarry Smith    Input Parameters:
8429a42bb27SBarry Smith +    da - the DM that defines the grid
84347c6ae99SBarry Smith .    i - the component of the function we wish to compute (must be local)
84447c6ae99SBarry Smith .    vu - input vector
84547c6ae99SBarry Smith .    vfu - output value
84647c6ae99SBarry Smith -    w - any user data
84747c6ae99SBarry Smith 
84847c6ae99SBarry Smith     Notes: Does NOT do ghost updates on vu upon entry
84947c6ae99SBarry Smith 
85047c6ae99SBarry Smith     Level: advanced
85147c6ae99SBarry Smith 
852aa219208SBarry Smith .seealso: DMDAComputeJacobian1WithAdic()
85347c6ae99SBarry Smith 
85447c6ae99SBarry Smith @*/
8557087cfbeSBarry Smith PetscErrorCode  DMDAFormFunctionib1(DM da,PetscInt i,Vec vu,PetscScalar *vfu,void *w)
85647c6ae99SBarry Smith {
85747c6ae99SBarry Smith   PetscErrorCode ierr;
85847c6ae99SBarry Smith   void           *u;
859aa219208SBarry Smith   DMDALocalInfo  info;
86047c6ae99SBarry Smith   MatStencil     stencil;
86147c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
86247c6ae99SBarry Smith 
86347c6ae99SBarry Smith   PetscFunctionBegin;
864aa219208SBarry Smith   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
865aa219208SBarry Smith   ierr = DMDAVecGetArray(da,vu,&u);CHKERRQ(ierr);
86647c6ae99SBarry Smith 
86747c6ae99SBarry Smith   /* figure out stencil value from i */
86847c6ae99SBarry Smith   stencil.c = i % info.dof;
86947c6ae99SBarry Smith   if (stencil.c) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONG,"Point-block functions can only be called for the entire block");
87047c6ae99SBarry Smith   stencil.i = (i % (info.xm*info.dof))/info.dof;
87147c6ae99SBarry Smith   stencil.j = (i % (info.xm*info.ym*info.dof))/(info.xm*info.dof);
87247c6ae99SBarry Smith   stencil.k = i/(info.xm*info.ym*info.dof);
87347c6ae99SBarry Smith 
87447c6ae99SBarry Smith   ierr = (*dd->lfib)(&info,&stencil,u,vfu,w);CHKERRQ(ierr);
87547c6ae99SBarry Smith 
876aa219208SBarry Smith   ierr = DMDAVecRestoreArray(da,vu,&u);CHKERRQ(ierr);
87747c6ae99SBarry Smith   PetscFunctionReturn(0);
87847c6ae99SBarry Smith }
87947c6ae99SBarry Smith 
88047c6ae99SBarry Smith #if defined(new)
88147c6ae99SBarry Smith #undef __FUNCT__
882aa219208SBarry Smith #define __FUNCT__ "DMDAGetDiagonal_MFFD"
88347c6ae99SBarry Smith /*
884aa219208SBarry Smith   DMDAGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix where local
885aa219208SBarry Smith     function lives on a DMDA
88647c6ae99SBarry Smith 
88747c6ae99SBarry Smith         y ~= (F(u + ha) - F(u))/h,
88847c6ae99SBarry Smith   where F = nonlinear function, as set by SNESSetFunction()
88947c6ae99SBarry Smith         u = current iterate
89047c6ae99SBarry Smith         h = difference interval
89147c6ae99SBarry Smith */
892aa219208SBarry Smith PetscErrorCode DMDAGetDiagonal_MFFD(DM da,Vec U,Vec a)
89347c6ae99SBarry Smith {
89447c6ae99SBarry Smith   PetscScalar    h,*aa,*ww,v;
89547c6ae99SBarry Smith   PetscReal      epsilon = PETSC_SQRT_MACHINE_EPSILON,umin = 100.0*PETSC_SQRT_MACHINE_EPSILON;
89647c6ae99SBarry Smith   PetscErrorCode ierr;
89747c6ae99SBarry Smith   PetscInt       gI,nI;
89847c6ae99SBarry Smith   MatStencil     stencil;
899aa219208SBarry Smith   DMDALocalInfo  info;
90047c6ae99SBarry Smith 
90147c6ae99SBarry Smith   PetscFunctionBegin;
90247c6ae99SBarry Smith   ierr = (*ctx->func)(0,U,a,ctx->funcctx);CHKERRQ(ierr);
90347c6ae99SBarry Smith   ierr = (*ctx->funcisetbase)(U,ctx->funcctx);CHKERRQ(ierr);
90447c6ae99SBarry Smith 
90547c6ae99SBarry Smith   ierr = VecGetArray(U,&ww);CHKERRQ(ierr);
90647c6ae99SBarry Smith   ierr = VecGetArray(a,&aa);CHKERRQ(ierr);
90747c6ae99SBarry Smith 
90847c6ae99SBarry Smith   nI = 0;
90947c6ae99SBarry Smith     h  = ww[gI];
91047c6ae99SBarry Smith     if (h == 0.0) h = 1.0;
91147c6ae99SBarry Smith #if !defined(PETSC_USE_COMPLEX)
91247c6ae99SBarry Smith     if (h < umin && h >= 0.0)      h = umin;
91347c6ae99SBarry Smith     else if (h < 0.0 && h > -umin) h = -umin;
91447c6ae99SBarry Smith #else
91547c6ae99SBarry Smith     if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0)     h = umin;
91647c6ae99SBarry Smith     else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin;
91747c6ae99SBarry Smith #endif
91847c6ae99SBarry Smith     h     *= epsilon;
91947c6ae99SBarry Smith 
92047c6ae99SBarry Smith     ww[gI] += h;
92147c6ae99SBarry Smith     ierr          = (*ctx->funci)(i,w,&v,ctx->funcctx);CHKERRQ(ierr);
92247c6ae99SBarry Smith     aa[nI]  = (v - aa[nI])/h;
92347c6ae99SBarry Smith     ww[gI] -= h;
92447c6ae99SBarry Smith     nI++;
92547c6ae99SBarry Smith   }
92647c6ae99SBarry Smith   ierr = VecRestoreArray(U,&ww);CHKERRQ(ierr);
92747c6ae99SBarry Smith   ierr = VecRestoreArray(a,&aa);CHKERRQ(ierr);
92847c6ae99SBarry Smith   PetscFunctionReturn(0);
92947c6ae99SBarry Smith }
93047c6ae99SBarry Smith #endif
93147c6ae99SBarry Smith 
93247c6ae99SBarry Smith #if defined(PETSC_HAVE_ADIC)
93347c6ae99SBarry Smith EXTERN_C_BEGIN
934c6db04a5SJed Brown #include <adic/ad_utils.h>
93547c6ae99SBarry Smith EXTERN_C_END
93647c6ae99SBarry Smith 
93747c6ae99SBarry Smith #undef __FUNCT__
938aa219208SBarry Smith #define __FUNCT__ "DMDAComputeJacobian1WithAdic"
93947c6ae99SBarry Smith /*@C
940aa219208SBarry Smith     DMDAComputeJacobian1WithAdic - Evaluates a adiC provided Jacobian function on each processor that
941aa219208SBarry Smith         share a DMDA
94247c6ae99SBarry Smith 
94347c6ae99SBarry Smith    Input Parameters:
9449a42bb27SBarry Smith +    da - the DM that defines the grid
94547c6ae99SBarry Smith .    vu - input vector (ghosted)
94647c6ae99SBarry Smith .    J - output matrix
94747c6ae99SBarry Smith -    w - any user data
94847c6ae99SBarry Smith 
94947c6ae99SBarry Smith    Level: advanced
95047c6ae99SBarry Smith 
95147c6ae99SBarry Smith     Notes: Does NOT do ghost updates on vu upon entry
95247c6ae99SBarry Smith 
953aa219208SBarry Smith .seealso: DMDAFormFunction1()
95447c6ae99SBarry Smith 
95547c6ae99SBarry Smith @*/
9567087cfbeSBarry Smith PetscErrorCode  DMDAComputeJacobian1WithAdic(DM da,Vec vu,Mat J,void *w)
95747c6ae99SBarry Smith {
95847c6ae99SBarry Smith   PetscErrorCode ierr;
95947c6ae99SBarry Smith   PetscInt       gtdof,tdof;
96047c6ae99SBarry Smith   PetscScalar    *ustart;
961aa219208SBarry Smith   DMDALocalInfo  info;
96247c6ae99SBarry Smith   void           *ad_u,*ad_f,*ad_ustart,*ad_fstart;
96347c6ae99SBarry Smith   ISColoring     iscoloring;
96447c6ae99SBarry Smith 
96547c6ae99SBarry Smith   PetscFunctionBegin;
966aa219208SBarry Smith   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
96747c6ae99SBarry Smith 
96847c6ae99SBarry Smith   PetscADResetIndep();
96947c6ae99SBarry Smith 
97047c6ae99SBarry Smith   /* get space for derivative objects.  */
971aa219208SBarry Smith   ierr = DMDAGetAdicArray(da,PETSC_TRUE,&ad_u,&ad_ustart,&gtdof);CHKERRQ(ierr);
972aa219208SBarry Smith   ierr = DMDAGetAdicArray(da,PETSC_FALSE,&ad_f,&ad_fstart,&tdof);CHKERRQ(ierr);
97347c6ae99SBarry Smith   ierr = VecGetArray(vu,&ustart);CHKERRQ(ierr);
97494013140SBarry Smith   ierr = DMGetColoring(da,IS_COLORING_GHOSTED,MATAIJ,&iscoloring);CHKERRQ(ierr);
97547c6ae99SBarry Smith 
97647c6ae99SBarry Smith   PetscADSetValueAndColor(ad_ustart,gtdof,iscoloring->colors,ustart);
97747c6ae99SBarry Smith 
97847c6ae99SBarry Smith   ierr = VecRestoreArray(vu,&ustart);CHKERRQ(ierr);
979fcfd50ebSBarry Smith   ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr);
98047c6ae99SBarry Smith   ierr = PetscADIncrementTotalGradSize(iscoloring->n);CHKERRQ(ierr);
98147c6ae99SBarry Smith   PetscADSetIndepDone();
98247c6ae99SBarry Smith 
983aa219208SBarry Smith   ierr = PetscLogEventBegin(DMDA_LocalADFunction,0,0,0,0);CHKERRQ(ierr);
98447c6ae99SBarry Smith   ierr = (*dd->adic_lf)(&info,ad_u,ad_f,w);CHKERRQ(ierr);
985aa219208SBarry Smith   ierr = PetscLogEventEnd(DMDA_LocalADFunction,0,0,0,0);CHKERRQ(ierr);
98647c6ae99SBarry Smith 
98747c6ae99SBarry Smith   /* stick the values into the matrix */
98847c6ae99SBarry Smith   ierr = MatSetValuesAdic(J,(PetscScalar**)ad_fstart);CHKERRQ(ierr);
98947c6ae99SBarry Smith 
99047c6ae99SBarry Smith   /* return space for derivative objects.  */
991aa219208SBarry Smith   ierr = DMDARestoreAdicArray(da,PETSC_TRUE,&ad_u,&ad_ustart,&gtdof);CHKERRQ(ierr);
992aa219208SBarry Smith   ierr = DMDARestoreAdicArray(da,PETSC_FALSE,&ad_f,&ad_fstart,&tdof);CHKERRQ(ierr);
99347c6ae99SBarry Smith   PetscFunctionReturn(0);
99447c6ae99SBarry Smith }
99547c6ae99SBarry Smith 
99647c6ae99SBarry Smith #undef __FUNCT__
997aa219208SBarry Smith #define __FUNCT__ "DMDAMultiplyByJacobian1WithAdic"
99847c6ae99SBarry Smith /*@C
999aa219208SBarry Smith     DMDAMultiplyByJacobian1WithAdic - Applies an ADIC-provided Jacobian function to a vector on
1000aa219208SBarry Smith     each processor that shares a DMDA.
100147c6ae99SBarry Smith 
100247c6ae99SBarry Smith     Input Parameters:
10039a42bb27SBarry Smith +   da - the DM that defines the grid
100447c6ae99SBarry Smith .   vu - Jacobian is computed at this point (ghosted)
100547c6ae99SBarry Smith .   v - product is done on this vector (ghosted)
100647c6ae99SBarry Smith .   fu - output vector = J(vu)*v (not ghosted)
100747c6ae99SBarry Smith -   w - any user data
100847c6ae99SBarry Smith 
100947c6ae99SBarry Smith     Notes:
101047c6ae99SBarry Smith     This routine does NOT do ghost updates on vu upon entry.
101147c6ae99SBarry Smith 
101247c6ae99SBarry Smith    Level: advanced
101347c6ae99SBarry Smith 
1014aa219208SBarry Smith .seealso: DMDAFormFunction1()
101547c6ae99SBarry Smith 
101647c6ae99SBarry Smith @*/
10177087cfbeSBarry Smith PetscErrorCode  DMDAMultiplyByJacobian1WithAdic(DM da,Vec vu,Vec v,Vec f,void *w)
101847c6ae99SBarry Smith {
101947c6ae99SBarry Smith   PetscErrorCode ierr;
102047c6ae99SBarry Smith   PetscInt       i,gtdof,tdof;
102147c6ae99SBarry Smith   PetscScalar    *avu,*av,*af,*ad_vustart,*ad_fstart;
1022aa219208SBarry Smith   DMDALocalInfo  info;
102347c6ae99SBarry Smith   void           *ad_vu,*ad_f;
102447c6ae99SBarry Smith 
102547c6ae99SBarry Smith   PetscFunctionBegin;
1026aa219208SBarry Smith   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
102747c6ae99SBarry Smith 
102847c6ae99SBarry Smith   /* get space for derivative objects.  */
1029aa219208SBarry Smith   ierr = DMDAGetAdicMFArray(da,PETSC_TRUE,&ad_vu,&ad_vustart,&gtdof);CHKERRQ(ierr);
1030aa219208SBarry Smith   ierr = DMDAGetAdicMFArray(da,PETSC_FALSE,&ad_f,&ad_fstart,&tdof);CHKERRQ(ierr);
103147c6ae99SBarry Smith 
103247c6ae99SBarry Smith   /* copy input vector into derivative object */
103347c6ae99SBarry Smith   ierr = VecGetArray(vu,&avu);CHKERRQ(ierr);
103447c6ae99SBarry Smith   ierr = VecGetArray(v,&av);CHKERRQ(ierr);
103547c6ae99SBarry Smith   for (i=0; i<gtdof; i++) {
103647c6ae99SBarry Smith     ad_vustart[2*i]   = avu[i];
103747c6ae99SBarry Smith     ad_vustart[2*i+1] = av[i];
103847c6ae99SBarry Smith   }
103947c6ae99SBarry Smith   ierr = VecRestoreArray(vu,&avu);CHKERRQ(ierr);
104047c6ae99SBarry Smith   ierr = VecRestoreArray(v,&av);CHKERRQ(ierr);
104147c6ae99SBarry Smith 
104247c6ae99SBarry Smith   PetscADResetIndep();
104347c6ae99SBarry Smith   ierr = PetscADIncrementTotalGradSize(1);CHKERRQ(ierr);
104447c6ae99SBarry Smith   PetscADSetIndepDone();
104547c6ae99SBarry Smith 
104647c6ae99SBarry Smith   ierr = (*dd->adicmf_lf)(&info,ad_vu,ad_f,w);CHKERRQ(ierr);
104747c6ae99SBarry Smith 
104847c6ae99SBarry Smith   /* stick the values into the vector */
104947c6ae99SBarry Smith   ierr = VecGetArray(f,&af);CHKERRQ(ierr);
105047c6ae99SBarry Smith   for (i=0; i<tdof; i++) {
105147c6ae99SBarry Smith     af[i] = ad_fstart[2*i+1];
105247c6ae99SBarry Smith   }
105347c6ae99SBarry Smith   ierr = VecRestoreArray(f,&af);CHKERRQ(ierr);
105447c6ae99SBarry Smith 
105547c6ae99SBarry Smith   /* return space for derivative objects.  */
1056aa219208SBarry Smith   ierr = DMDARestoreAdicMFArray(da,PETSC_TRUE,&ad_vu,&ad_vustart,&gtdof);CHKERRQ(ierr);
1057aa219208SBarry Smith   ierr = DMDARestoreAdicMFArray(da,PETSC_FALSE,&ad_f,&ad_fstart,&tdof);CHKERRQ(ierr);
105847c6ae99SBarry Smith   PetscFunctionReturn(0);
105947c6ae99SBarry Smith }
106047c6ae99SBarry Smith #endif
106147c6ae99SBarry Smith 
106247c6ae99SBarry Smith #undef __FUNCT__
1063aa219208SBarry Smith #define __FUNCT__ "DMDAComputeJacobian1"
106447c6ae99SBarry Smith /*@
1065aa219208SBarry Smith     DMDAComputeJacobian1 - Evaluates a local Jacobian function on each processor that
1066aa219208SBarry Smith         share a DMDA
106747c6ae99SBarry Smith 
106847c6ae99SBarry Smith    Input Parameters:
10699a42bb27SBarry Smith +    da - the DM that defines the grid
107047c6ae99SBarry Smith .    vu - input vector (ghosted)
107147c6ae99SBarry Smith .    J - output matrix
107247c6ae99SBarry Smith -    w - any user data
107347c6ae99SBarry Smith 
107447c6ae99SBarry Smith     Notes: Does NOT do ghost updates on vu upon entry
107547c6ae99SBarry Smith 
107647c6ae99SBarry Smith     Level: advanced
107747c6ae99SBarry Smith 
1078aa219208SBarry Smith .seealso: DMDAFormFunction1()
107947c6ae99SBarry Smith 
108047c6ae99SBarry Smith @*/
10817087cfbeSBarry Smith PetscErrorCode  DMDAComputeJacobian1(DM da,Vec vu,Mat J,void *w)
108247c6ae99SBarry Smith {
108347c6ae99SBarry Smith   PetscErrorCode ierr;
108447c6ae99SBarry Smith   void           *u;
1085aa219208SBarry Smith   DMDALocalInfo  info;
108647c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
108747c6ae99SBarry Smith 
108847c6ae99SBarry Smith   PetscFunctionBegin;
1089aa219208SBarry Smith   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
1090aa219208SBarry Smith   ierr = DMDAVecGetArray(da,vu,&u);CHKERRQ(ierr);
109147c6ae99SBarry Smith   ierr = (*dd->lj)(&info,u,J,w);CHKERRQ(ierr);
1092aa219208SBarry Smith   ierr = DMDAVecRestoreArray(da,vu,&u);CHKERRQ(ierr);
109347c6ae99SBarry Smith   PetscFunctionReturn(0);
109447c6ae99SBarry Smith }
109547c6ae99SBarry Smith 
109647c6ae99SBarry Smith 
109747c6ae99SBarry Smith #undef __FUNCT__
1098aa219208SBarry Smith #define __FUNCT__ "DMDAComputeJacobian1WithAdifor"
109947c6ae99SBarry Smith /*
1100aa219208SBarry Smith     DMDAComputeJacobian1WithAdifor - Evaluates a ADIFOR provided Jacobian local function on each processor that
1101aa219208SBarry Smith         share a DMDA
110247c6ae99SBarry Smith 
110347c6ae99SBarry Smith    Input Parameters:
11049a42bb27SBarry Smith +    da - the DM that defines the grid
110547c6ae99SBarry Smith .    vu - input vector (ghosted)
110647c6ae99SBarry Smith .    J - output matrix
110747c6ae99SBarry Smith -    w - any user data
110847c6ae99SBarry Smith 
110947c6ae99SBarry Smith     Notes: Does NOT do ghost updates on vu upon entry
111047c6ae99SBarry Smith 
1111aa219208SBarry Smith .seealso: DMDAFormFunction1()
111247c6ae99SBarry Smith 
111347c6ae99SBarry Smith */
11147087cfbeSBarry Smith PetscErrorCode  DMDAComputeJacobian1WithAdifor(DM da,Vec vu,Mat J,void *w)
111547c6ae99SBarry Smith {
111647c6ae99SBarry Smith   PetscErrorCode  ierr;
111747c6ae99SBarry Smith   PetscInt        i,Nc,N;
111847c6ae99SBarry Smith   ISColoringValue *color;
1119aa219208SBarry Smith   DMDALocalInfo   info;
112047c6ae99SBarry Smith   PetscScalar     *u,*g_u,*g_f,*f = 0,*p_u;
112147c6ae99SBarry Smith   ISColoring      iscoloring;
112247c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
1123aa219208SBarry Smith   void            (*lf)(PetscInt*,DMDALocalInfo*,PetscScalar*,PetscScalar*,PetscInt*,PetscScalar*,PetscScalar*,PetscInt*,void*,PetscErrorCode*) =
1124aa219208SBarry Smith                   (void (*)(PetscInt*,DMDALocalInfo*,PetscScalar*,PetscScalar*,PetscInt*,PetscScalar*,PetscScalar*,PetscInt*,void*,PetscErrorCode*))*dd->adifor_lf;
112547c6ae99SBarry Smith 
112647c6ae99SBarry Smith   PetscFunctionBegin;
112794013140SBarry Smith   ierr = DMGetColoring(da,IS_COLORING_GHOSTED,MATAIJ,&iscoloring);CHKERRQ(ierr);
112847c6ae99SBarry Smith   Nc   = iscoloring->n;
1129aa219208SBarry Smith   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
113047c6ae99SBarry Smith   N    = info.gxm*info.gym*info.gzm*info.dof;
113147c6ae99SBarry Smith 
113247c6ae99SBarry Smith   /* get space for derivative objects.  */
113347c6ae99SBarry Smith   ierr  = PetscMalloc(Nc*info.gxm*info.gym*info.gzm*info.dof*sizeof(PetscScalar),&g_u);CHKERRQ(ierr);
113447c6ae99SBarry Smith   ierr  = PetscMemzero(g_u,Nc*info.gxm*info.gym*info.gzm*info.dof*sizeof(PetscScalar));CHKERRQ(ierr);
113547c6ae99SBarry Smith   p_u   = g_u;
113647c6ae99SBarry Smith   color = iscoloring->colors;
113747c6ae99SBarry Smith   for (i=0; i<N; i++) {
113847c6ae99SBarry Smith     p_u[*color++] = 1.0;
113947c6ae99SBarry Smith     p_u          += Nc;
114047c6ae99SBarry Smith   }
1141fcfd50ebSBarry Smith   ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr);
114247c6ae99SBarry Smith   ierr = PetscMalloc2(Nc*info.xm*info.ym*info.zm*info.dof,PetscScalar,&g_f,info.xm*info.ym*info.zm*info.dof,PetscScalar,&f);CHKERRQ(ierr);
114347c6ae99SBarry Smith 
114447c6ae99SBarry Smith   /* Seed the input array g_u with coloring information */
114547c6ae99SBarry Smith 
114647c6ae99SBarry Smith   ierr = VecGetArray(vu,&u);CHKERRQ(ierr);
114747c6ae99SBarry Smith   (lf)(&Nc,&info,u,g_u,&Nc,f,g_f,&Nc,w,&ierr);CHKERRQ(ierr);
114847c6ae99SBarry Smith   ierr = VecRestoreArray(vu,&u);CHKERRQ(ierr);
114947c6ae99SBarry Smith 
115047c6ae99SBarry Smith   /* stick the values into the matrix */
115147c6ae99SBarry Smith   /* PetscScalarView(Nc*info.xm*info.ym,g_f,0); */
115247c6ae99SBarry Smith   ierr = MatSetValuesAdifor(J,Nc,g_f);CHKERRQ(ierr);
115347c6ae99SBarry Smith 
115447c6ae99SBarry Smith   /* return space for derivative objects.  */
115547c6ae99SBarry Smith   ierr = PetscFree(g_u);CHKERRQ(ierr);
115647c6ae99SBarry Smith   ierr = PetscFree2(g_f,f);CHKERRQ(ierr);
115747c6ae99SBarry Smith   PetscFunctionReturn(0);
115847c6ae99SBarry Smith }
115947c6ae99SBarry Smith 
116047c6ae99SBarry Smith #undef __FUNCT__
1161aa219208SBarry Smith #define __FUNCT__ "DMDAFormJacobianLocal"
116247c6ae99SBarry Smith /*@C
1163aa219208SBarry Smith    DMDAFormjacobianLocal - This is a universal Jacobian evaluation routine for
11649a42bb27SBarry Smith    a local DM function.
116547c6ae99SBarry Smith 
1166aa219208SBarry Smith    Collective on DMDA
116747c6ae99SBarry Smith 
116847c6ae99SBarry Smith    Input Parameters:
11699a42bb27SBarry Smith +  da - the DM context
117047c6ae99SBarry Smith .  func - The local function
117147c6ae99SBarry Smith .  X - input vector
117247c6ae99SBarry Smith .  J - Jacobian matrix
117347c6ae99SBarry Smith -  ctx - A user context
117447c6ae99SBarry Smith 
117547c6ae99SBarry Smith    Level: intermediate
117647c6ae99SBarry Smith 
1177aa219208SBarry Smith .seealso: DMDASetLocalFunction(), DMDASetLocalJacobian(), DMDASetLocalAdicFunction(), DMDASetLocalAdicMFFunction(),
117847c6ae99SBarry Smith           SNESSetFunction(), SNESSetJacobian()
117947c6ae99SBarry Smith 
118047c6ae99SBarry Smith @*/
11817087cfbeSBarry Smith PetscErrorCode  DMDAFormJacobianLocal(DM da, DMDALocalFunction1 func, Vec X, Mat J, void *ctx)
118247c6ae99SBarry Smith {
118347c6ae99SBarry Smith   Vec            localX;
1184aa219208SBarry Smith   DMDALocalInfo  info;
118547c6ae99SBarry Smith   void           *u;
118647c6ae99SBarry Smith   PetscErrorCode ierr;
118747c6ae99SBarry Smith 
118847c6ae99SBarry Smith   PetscFunctionBegin;
11899a42bb27SBarry Smith   ierr = DMGetLocalVector(da,&localX);CHKERRQ(ierr);
119047c6ae99SBarry Smith   /*
119147c6ae99SBarry Smith      Scatter ghost points to local vector, using the 2-step process
11929a42bb27SBarry Smith         DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
119347c6ae99SBarry Smith   */
11949a42bb27SBarry Smith   ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
11959a42bb27SBarry Smith   ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
1196aa219208SBarry Smith   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
1197aa219208SBarry Smith   ierr = DMDAVecGetArray(da,localX,&u);CHKERRQ(ierr);
119839d508bbSBarry Smith   ierr = (*func)(&info,u,J,ctx);CHKERRQ(ierr);
1199aa219208SBarry Smith   ierr = DMDAVecRestoreArray(da,localX,&u);CHKERRQ(ierr);
12009a42bb27SBarry Smith   ierr = DMRestoreLocalVector(da,&localX);CHKERRQ(ierr);
120147c6ae99SBarry Smith   PetscFunctionReturn(0);
120247c6ae99SBarry Smith }
120347c6ae99SBarry Smith 
120447c6ae99SBarry Smith #undef __FUNCT__
1205aa219208SBarry Smith #define __FUNCT__ "DMDAMultiplyByJacobian1WithAD"
120647c6ae99SBarry Smith /*@C
1207aa219208SBarry Smith     DMDAMultiplyByJacobian1WithAD - Applies a Jacobian function supplied by ADIFOR or ADIC
1208aa219208SBarry Smith     to a vector on each processor that shares a DMDA.
120947c6ae99SBarry Smith 
121047c6ae99SBarry Smith    Input Parameters:
12119a42bb27SBarry Smith +    da - the DM that defines the grid
121247c6ae99SBarry Smith .    vu - Jacobian is computed at this point (ghosted)
121347c6ae99SBarry Smith .    v - product is done on this vector (ghosted)
121447c6ae99SBarry Smith .    fu - output vector = J(vu)*v (not ghosted)
121547c6ae99SBarry Smith -    w - any user data
121647c6ae99SBarry Smith 
121747c6ae99SBarry Smith     Notes:
121847c6ae99SBarry Smith     This routine does NOT do ghost updates on vu and v upon entry.
121947c6ae99SBarry Smith 
1220aa219208SBarry Smith     Automatically calls DMDAMultiplyByJacobian1WithAdifor() or DMDAMultiplyByJacobian1WithAdic()
1221aa219208SBarry Smith     depending on whether DMDASetLocalAdicMFFunction() or DMDASetLocalAdiforMFFunction() was called.
122247c6ae99SBarry Smith 
122347c6ae99SBarry Smith    Level: advanced
122447c6ae99SBarry Smith 
1225aa219208SBarry Smith .seealso: DMDAFormFunction1(), DMDAMultiplyByJacobian1WithAdifor(), DMDAMultiplyByJacobian1WithAdic()
122647c6ae99SBarry Smith 
122747c6ae99SBarry Smith @*/
12287087cfbeSBarry Smith PetscErrorCode  DMDAMultiplyByJacobian1WithAD(DM da,Vec u,Vec v,Vec f,void *w)
122947c6ae99SBarry Smith {
123047c6ae99SBarry Smith   PetscErrorCode ierr;
123147c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
123247c6ae99SBarry Smith 
123347c6ae99SBarry Smith   PetscFunctionBegin;
123447c6ae99SBarry Smith   if (dd->adicmf_lf) {
123547c6ae99SBarry Smith #if defined(PETSC_HAVE_ADIC)
1236aa219208SBarry Smith     ierr = DMDAMultiplyByJacobian1WithAdic(da,u,v,f,w);CHKERRQ(ierr);
123747c6ae99SBarry Smith #else
123847c6ae99SBarry Smith     SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP_SYS,"Requires ADIC to be installed and cannot use complex numbers");
123947c6ae99SBarry Smith #endif
124047c6ae99SBarry Smith   } else if (dd->adiformf_lf) {
1241aa219208SBarry Smith     ierr = DMDAMultiplyByJacobian1WithAdifor(da,u,v,f,w);CHKERRQ(ierr);
124247c6ae99SBarry Smith   } else {
1243aa219208SBarry Smith     SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ORDER,"Must call DMDASetLocalAdiforMFFunction() or DMDASetLocalAdicMFFunction() before using");
124447c6ae99SBarry Smith   }
124547c6ae99SBarry Smith   PetscFunctionReturn(0);
124647c6ae99SBarry Smith }
124747c6ae99SBarry Smith 
124847c6ae99SBarry Smith 
124947c6ae99SBarry Smith #undef __FUNCT__
1250aa219208SBarry Smith #define __FUNCT__ "DMDAMultiplyByJacobian1WithAdifor"
125147c6ae99SBarry Smith /*@C
1252aa219208SBarry Smith     DMDAMultiplyByJacobian1WithAdifor - Applies a ADIFOR provided Jacobian function on each processor that
12539a42bb27SBarry Smith         share a DM to a vector
125447c6ae99SBarry Smith 
125547c6ae99SBarry Smith    Input Parameters:
12569a42bb27SBarry Smith +    da - the DM that defines the grid
125747c6ae99SBarry Smith .    vu - Jacobian is computed at this point (ghosted)
125847c6ae99SBarry Smith .    v - product is done on this vector (ghosted)
125947c6ae99SBarry Smith .    fu - output vector = J(vu)*v (not ghosted)
126047c6ae99SBarry Smith -    w - any user data
126147c6ae99SBarry Smith 
126247c6ae99SBarry Smith     Notes: Does NOT do ghost updates on vu and v upon entry
126347c6ae99SBarry Smith 
126447c6ae99SBarry Smith    Level: advanced
126547c6ae99SBarry Smith 
1266aa219208SBarry Smith .seealso: DMDAFormFunction1()
126747c6ae99SBarry Smith 
126847c6ae99SBarry Smith @*/
12697087cfbeSBarry Smith PetscErrorCode  DMDAMultiplyByJacobian1WithAdifor(DM da,Vec u,Vec v,Vec f,void *w)
127047c6ae99SBarry Smith {
127147c6ae99SBarry Smith   PetscErrorCode ierr;
127247c6ae99SBarry Smith   PetscScalar    *au,*av,*af,*awork;
127347c6ae99SBarry Smith   Vec            work;
1274aa219208SBarry Smith   DMDALocalInfo  info;
127547c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
1276aa219208SBarry Smith   void           (*lf)(DMDALocalInfo*,PetscScalar*,PetscScalar*,PetscScalar*,PetscScalar*,void*,PetscErrorCode*) =
1277aa219208SBarry Smith                  (void (*)(DMDALocalInfo*,PetscScalar*,PetscScalar*,PetscScalar*,PetscScalar*,void*,PetscErrorCode*))*dd->adiformf_lf;
127847c6ae99SBarry Smith 
127947c6ae99SBarry Smith   PetscFunctionBegin;
1280aa219208SBarry Smith   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
128147c6ae99SBarry Smith 
12829a42bb27SBarry Smith   ierr = DMGetGlobalVector(da,&work);CHKERRQ(ierr);
128347c6ae99SBarry Smith   ierr = VecGetArray(u,&au);CHKERRQ(ierr);
128447c6ae99SBarry Smith   ierr = VecGetArray(v,&av);CHKERRQ(ierr);
128547c6ae99SBarry Smith   ierr = VecGetArray(f,&af);CHKERRQ(ierr);
128647c6ae99SBarry Smith   ierr = VecGetArray(work,&awork);CHKERRQ(ierr);
128747c6ae99SBarry Smith   (lf)(&info,au,av,awork,af,w,&ierr);CHKERRQ(ierr);
128847c6ae99SBarry Smith   ierr = VecRestoreArray(u,&au);CHKERRQ(ierr);
128947c6ae99SBarry Smith   ierr = VecRestoreArray(v,&av);CHKERRQ(ierr);
129047c6ae99SBarry Smith   ierr = VecRestoreArray(f,&af);CHKERRQ(ierr);
129147c6ae99SBarry Smith   ierr = VecRestoreArray(work,&awork);CHKERRQ(ierr);
12929a42bb27SBarry Smith   ierr = DMRestoreGlobalVector(da,&work);CHKERRQ(ierr);
129347c6ae99SBarry Smith 
129447c6ae99SBarry Smith   PetscFunctionReturn(0);
129547c6ae99SBarry Smith }
129647c6ae99SBarry Smith 
129747c6ae99SBarry Smith #undef __FUNCT__
12989a42bb27SBarry Smith #define __FUNCT__ "DMSetUp_DA_2D"
12997087cfbeSBarry Smith PetscErrorCode  DMSetUp_DA_2D(DM da)
130047c6ae99SBarry Smith {
130147c6ae99SBarry Smith   DM_DA                  *dd = (DM_DA*)da->data;
130247c6ae99SBarry Smith   const PetscInt         M            = dd->M;
130347c6ae99SBarry Smith   const PetscInt         N            = dd->N;
130447c6ae99SBarry Smith   PetscInt               m            = dd->m;
130547c6ae99SBarry Smith   PetscInt               n            = dd->n;
130647c6ae99SBarry Smith   const PetscInt         dof          = dd->w;
130747c6ae99SBarry Smith   const PetscInt         s            = dd->s;
13081321219cSEthan Coon   const DMDABoundaryType bx         = dd->bx;
13091321219cSEthan Coon   const DMDABoundaryType by         = dd->by;
1310aa219208SBarry Smith   const DMDAStencilType  stencil_type = dd->stencil_type;
131147c6ae99SBarry Smith   PetscInt               *lx           = dd->lx;
131247c6ae99SBarry Smith   PetscInt               *ly           = dd->ly;
131347c6ae99SBarry Smith   MPI_Comm               comm;
131447c6ae99SBarry Smith   PetscMPIInt            rank,size;
1315ce00eea3SSatish Balay   PetscInt               xs,xe,ys,ye,x,y,Xs,Xe,Ys,Ye,start,end,IXs,IXe,IYs,IYe;
1316ce00eea3SSatish Balay   PetscInt               up,down,left,right,i,n0,n1,n2,n3,n5,n6,n7,n8,*idx,nn,*idx_cpy;
1317ce00eea3SSatish Balay   const PetscInt         *idx_full;
1318db87c5ecSEthan Coon   PetscInt               xbase,*bases,*ldims,j,x_t,y_t,s_t,base,count;
131947c6ae99SBarry Smith   PetscInt               s_x,s_y; /* s proportionalized to w */
132047c6ae99SBarry Smith   PetscInt               sn0 = 0,sn2 = 0,sn6 = 0,sn8 = 0;
132147c6ae99SBarry Smith   Vec                    local,global;
132247c6ae99SBarry Smith   VecScatter             ltog,gtol;
1323ce00eea3SSatish Balay   IS                     to,from,ltogis;
132447c6ae99SBarry Smith   PetscErrorCode         ierr;
132547c6ae99SBarry Smith 
132647c6ae99SBarry Smith   PetscFunctionBegin;
132747c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
132847c6ae99SBarry Smith 
132947c6ae99SBarry Smith   if (dof < 1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Must have 1 or more degrees of freedom per node: %D",dof);
133047c6ae99SBarry Smith   if (s < 0) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Stencil width cannot be negative: %D",s);
133147c6ae99SBarry Smith 
133247c6ae99SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
133347c6ae99SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
133447c6ae99SBarry Smith 
133547c6ae99SBarry Smith   dd->dim         = 2;
133647c6ae99SBarry Smith   ierr = PetscMalloc(dof*sizeof(char*),&dd->fieldname);CHKERRQ(ierr);
133747c6ae99SBarry Smith   ierr = PetscMemzero(dd->fieldname,dof*sizeof(char*));CHKERRQ(ierr);
133847c6ae99SBarry Smith 
133947c6ae99SBarry Smith   if (m != PETSC_DECIDE) {
134047c6ae99SBarry Smith     if (m < 1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in X direction: %D",m);
134147c6ae99SBarry Smith     else if (m > size) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in X direction: %D %d",m,size);
134247c6ae99SBarry Smith   }
134347c6ae99SBarry Smith   if (n != PETSC_DECIDE) {
134447c6ae99SBarry Smith     if (n < 1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Y direction: %D",n);
134547c6ae99SBarry Smith     else if (n > size) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Y direction: %D %d",n,size);
134647c6ae99SBarry Smith   }
134747c6ae99SBarry Smith 
134847c6ae99SBarry Smith   if (m == PETSC_DECIDE || n == PETSC_DECIDE) {
134947c6ae99SBarry Smith     if (n != PETSC_DECIDE) {
135047c6ae99SBarry Smith       m = size/n;
135147c6ae99SBarry Smith     } else if (m != PETSC_DECIDE) {
135247c6ae99SBarry Smith       n = size/m;
135347c6ae99SBarry Smith     } else {
135447c6ae99SBarry Smith       /* try for squarish distribution */
135547c6ae99SBarry Smith       m = (PetscInt)(0.5 + sqrt(((double)M)*((double)size)/((double)N)));
135647c6ae99SBarry Smith       if (!m) m = 1;
135747c6ae99SBarry Smith       while (m > 0) {
135847c6ae99SBarry Smith 	n = size/m;
135947c6ae99SBarry Smith 	if (m*n == size) break;
136047c6ae99SBarry Smith 	m--;
136147c6ae99SBarry Smith       }
136247c6ae99SBarry Smith       if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;}
136347c6ae99SBarry Smith     }
136447c6ae99SBarry Smith     if (m*n != size) SETERRQ(comm,PETSC_ERR_PLIB,"Unable to create partition, check the size of the communicator and input m and n ");
136547c6ae99SBarry Smith   } else if (m*n != size) SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"Given Bad partition");
136647c6ae99SBarry Smith 
136747c6ae99SBarry Smith   if (M < m) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in x direction is too fine! %D %D",M,m);
136847c6ae99SBarry Smith   if (N < n) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in y direction is too fine! %D %D",N,n);
136947c6ae99SBarry Smith 
137047c6ae99SBarry Smith   /*
137147c6ae99SBarry Smith      Determine locally owned region
137247c6ae99SBarry Smith      xs is the first local node number, x is the number of local nodes
137347c6ae99SBarry Smith   */
137447c6ae99SBarry Smith   if (!lx) {
137547c6ae99SBarry Smith     ierr = PetscMalloc(m*sizeof(PetscInt), &dd->lx);CHKERRQ(ierr);
137647c6ae99SBarry Smith     lx = dd->lx;
137747c6ae99SBarry Smith     for (i=0; i<m; i++) {
137847c6ae99SBarry Smith       lx[i] = M/m + ((M % m) > i);
137947c6ae99SBarry Smith     }
138047c6ae99SBarry Smith   }
138147c6ae99SBarry Smith   x  = lx[rank % m];
138247c6ae99SBarry Smith   xs = 0;
138347c6ae99SBarry Smith   for (i=0; i<(rank % m); i++) {
138447c6ae99SBarry Smith     xs += lx[i];
138547c6ae99SBarry Smith   }
138647c6ae99SBarry Smith #if defined(PETSC_USE_DEBUG)
138747c6ae99SBarry Smith   left = xs;
138847c6ae99SBarry Smith   for (i=(rank % m); i<m; i++) {
138947c6ae99SBarry Smith     left += lx[i];
139047c6ae99SBarry Smith   }
139147c6ae99SBarry Smith   if (left != M) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Sum of lx across processors not equal to M: %D %D",left,M);
139247c6ae99SBarry Smith #endif
139347c6ae99SBarry Smith 
139447c6ae99SBarry Smith   /*
139547c6ae99SBarry Smith      Determine locally owned region
139647c6ae99SBarry Smith      ys is the first local node number, y is the number of local nodes
139747c6ae99SBarry Smith   */
139847c6ae99SBarry Smith   if (!ly) {
139947c6ae99SBarry Smith     ierr = PetscMalloc(n*sizeof(PetscInt), &dd->ly);CHKERRQ(ierr);
140047c6ae99SBarry Smith     ly = dd->ly;
140147c6ae99SBarry Smith     for (i=0; i<n; i++) {
140247c6ae99SBarry Smith       ly[i] = N/n + ((N % n) > i);
140347c6ae99SBarry Smith     }
140447c6ae99SBarry Smith   }
140547c6ae99SBarry Smith   y  = ly[rank/m];
140647c6ae99SBarry Smith   ys = 0;
140747c6ae99SBarry Smith   for (i=0; i<(rank/m); i++) {
140847c6ae99SBarry Smith     ys += ly[i];
140947c6ae99SBarry Smith   }
141047c6ae99SBarry Smith #if defined(PETSC_USE_DEBUG)
141147c6ae99SBarry Smith   left = ys;
141247c6ae99SBarry Smith   for (i=(rank/m); i<n; i++) {
141347c6ae99SBarry Smith     left += ly[i];
141447c6ae99SBarry Smith   }
141547c6ae99SBarry Smith   if (left != N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Sum of ly across processors not equal to N: %D %D",left,N);
141647c6ae99SBarry Smith #endif
141747c6ae99SBarry Smith 
1418bcea557cSEthan Coon   /*
1419bcea557cSEthan Coon    check if the scatter requires more than one process neighbor or wraps around
1420bcea557cSEthan Coon    the domain more than once
1421bcea557cSEthan Coon   */
1422bcea557cSEthan Coon   if ((x < s) && ((m > 1) || (bx == DMDA_BOUNDARY_PERIODIC))) {
1423bcea557cSEthan Coon     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local x-width of domain x %D is smaller than stencil width s %D",x,s);
1424bcea557cSEthan Coon   }
1425bcea557cSEthan Coon   if ((y < s) && ((n > 1) || (by == DMDA_BOUNDARY_PERIODIC))) {
1426bcea557cSEthan Coon     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local y-width of domain y %D is smaller than stencil width s %D",y,s);
1427bcea557cSEthan Coon   }
142847c6ae99SBarry Smith   xe = xs + x;
142947c6ae99SBarry Smith   ye = ys + y;
143047c6ae99SBarry Smith 
1431ce00eea3SSatish Balay   /* determine ghost region (Xs) and region scattered into (IXs)  */
143247c6ae99SBarry Smith   /* Assume No Periodicity */
1433ce00eea3SSatish Balay   if (xs-s > 0) { Xs = xs - s; IXs = xs - s; } else { Xs = 0; IXs = 0; }
1434ce00eea3SSatish Balay   if (xe+s <= M) { Xe = xe + s; IXe = xe + s; } else { Xe = M; IXe = M; }
1435ce00eea3SSatish Balay   if (ys-s > 0) { Ys = ys - s; IYs = ys - s; } else { Ys = 0; IYs = 0; }
1436ce00eea3SSatish Balay   if (ye+s <= N) { Ye = ye + s; IYe = ye + s; } else { Ye = N; IYe = N; }
143747c6ae99SBarry Smith 
1438ce00eea3SSatish Balay   /* fix for periodicity/ghosted */
14391321219cSEthan Coon   if (bx) { Xs = xs - s; Xe = xe + s; }
14401321219cSEthan Coon   if (bx == DMDA_BOUNDARY_PERIODIC) { IXs = xs - s; IXe = xe + s; }
14411321219cSEthan Coon   if (by) { Ys = ys - s; Ye = ye + s; }
14421321219cSEthan Coon   if (by == DMDA_BOUNDARY_PERIODIC) { IYs = ys - s; IYe = ye + s; }
144347c6ae99SBarry Smith 
144447c6ae99SBarry Smith   /* Resize all X parameters to reflect w */
1445ce00eea3SSatish Balay   s_x = s;
144647c6ae99SBarry Smith   s_y = s;
144747c6ae99SBarry Smith 
144847c6ae99SBarry Smith   /* determine starting point of each processor */
144947c6ae99SBarry Smith   nn    = x*y;
145047c6ae99SBarry Smith   ierr  = PetscMalloc2(size+1,PetscInt,&bases,size,PetscInt,&ldims);CHKERRQ(ierr);
145147c6ae99SBarry Smith   ierr  = MPI_Allgather(&nn,1,MPIU_INT,ldims,1,MPIU_INT,comm);CHKERRQ(ierr);
145247c6ae99SBarry Smith   bases[0] = 0;
145347c6ae99SBarry Smith   for (i=1; i<=size; i++) {
145447c6ae99SBarry Smith     bases[i] = ldims[i-1];
145547c6ae99SBarry Smith   }
145647c6ae99SBarry Smith   for (i=1; i<=size; i++) {
145747c6ae99SBarry Smith     bases[i] += bases[i-1];
145847c6ae99SBarry Smith   }
1459ce00eea3SSatish Balay   base = bases[rank]*dof;
146047c6ae99SBarry Smith 
146147c6ae99SBarry Smith   /* allocate the base parallel and sequential vectors */
1462ce00eea3SSatish Balay   dd->Nlocal = x*y*dof;
146347c6ae99SBarry Smith   ierr = VecCreateMPIWithArray(comm,dd->Nlocal,PETSC_DECIDE,0,&global);CHKERRQ(ierr);
146447c6ae99SBarry Smith   ierr = VecSetBlockSize(global,dof);CHKERRQ(ierr);
1465ce00eea3SSatish Balay   dd->nlocal = (Xe-Xs)*(Ye-Ys)*dof;
146647c6ae99SBarry Smith   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,dd->nlocal,0,&local);CHKERRQ(ierr);
146747c6ae99SBarry Smith   ierr = VecSetBlockSize(local,dof);CHKERRQ(ierr);
146847c6ae99SBarry Smith 
146947c6ae99SBarry Smith   /* generate appropriate vector scatters */
147047c6ae99SBarry Smith   /* local to global inserts non-ghost point region into global */
147147c6ae99SBarry Smith   ierr = VecGetOwnershipRange(global,&start,&end);CHKERRQ(ierr);
1472ce00eea3SSatish Balay   ierr = ISCreateStride(comm,x*y*dof,start,1,&to);CHKERRQ(ierr);
147347c6ae99SBarry Smith 
1474db87c5ecSEthan Coon   count = x*y;
1475ce00eea3SSatish Balay   ierr = PetscMalloc(x*y*sizeof(PetscInt),&idx);CHKERRQ(ierr);
1476ce00eea3SSatish Balay   left = xs - Xs; right = left + x;
1477ce00eea3SSatish Balay   down = ys - Ys; up = down + y;
147847c6ae99SBarry Smith   count = 0;
147947c6ae99SBarry Smith   for (i=down; i<up; i++) {
1480ce00eea3SSatish Balay     for (j=left; j<right; j++) {
1481ce00eea3SSatish Balay       idx[count++] = i*(Xe-Xs) + j;
148247c6ae99SBarry Smith     }
148347c6ae99SBarry Smith   }
148447c6ae99SBarry Smith 
1485ce00eea3SSatish Balay   ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&from);CHKERRQ(ierr);
148647c6ae99SBarry Smith   ierr = VecScatterCreate(local,from,global,to,&ltog);CHKERRQ(ierr);
148747c6ae99SBarry Smith   ierr = PetscLogObjectParent(dd,ltog);CHKERRQ(ierr);
1488fcfd50ebSBarry Smith   ierr = ISDestroy(&from);CHKERRQ(ierr);
1489fcfd50ebSBarry Smith   ierr = ISDestroy(&to);CHKERRQ(ierr);
149047c6ae99SBarry Smith 
1491ce00eea3SSatish Balay   /* global to local must include ghost points within the domain,
1492ce00eea3SSatish Balay      but not ghost points outside the domain that aren't periodic */
1493aa219208SBarry Smith   if (stencil_type == DMDA_STENCIL_BOX) {
1494db87c5ecSEthan Coon     count = (IXe-IXs)*(IYe-IYs);
1495db87c5ecSEthan Coon     ierr  = PetscMalloc(count*sizeof(PetscInt),&idx);CHKERRQ(ierr);
1496ce00eea3SSatish Balay 
1497ce00eea3SSatish Balay     left = IXs - Xs; right = left + (IXe-IXs);
1498ce00eea3SSatish Balay     down = IYs - Ys; up = down + (IYe-IYs);
1499ce00eea3SSatish Balay     count = 0;
1500ce00eea3SSatish Balay     for (i=down; i<up; i++) {
1501ce00eea3SSatish Balay       for (j=left; j<right; j++) {
1502ce00eea3SSatish Balay         idx[count++] = j + i*(Xe-Xs);
1503ce00eea3SSatish Balay       }
1504ce00eea3SSatish Balay     }
1505ce00eea3SSatish Balay     ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr);
1506ce00eea3SSatish Balay 
150747c6ae99SBarry Smith   } else {
150847c6ae99SBarry Smith     /* must drop into cross shape region */
150947c6ae99SBarry Smith     /*       ---------|
151047c6ae99SBarry Smith             |  top    |
1511ce00eea3SSatish Balay          |---         ---| up
151247c6ae99SBarry Smith          |   middle      |
151347c6ae99SBarry Smith          |               |
1514ce00eea3SSatish Balay          ----         ---- down
151547c6ae99SBarry Smith             | bottom  |
151647c6ae99SBarry Smith             -----------
151747c6ae99SBarry Smith          Xs xs        xe Xe */
1518db87c5ecSEthan Coon     count = (ys-IYs)*x + y*(IXe-IXs) + (IYe-ye)*x;
1519db87c5ecSEthan Coon     ierr  = PetscMalloc(count*sizeof(PetscInt),&idx);CHKERRQ(ierr);
1520ce00eea3SSatish Balay 
1521ce00eea3SSatish Balay     left = xs - Xs; right = left + x;
1522ce00eea3SSatish Balay     down = ys - Ys; up = down + y;
152347c6ae99SBarry Smith     count = 0;
1524ce00eea3SSatish Balay     /* bottom */
1525ce00eea3SSatish Balay     for (i=(IYs-Ys); i<down; i++) {
1526ce00eea3SSatish Balay       for (j=left; j<right; j++) {
1527ce00eea3SSatish Balay         idx[count++] = j + i*(Xe-Xs);
152847c6ae99SBarry Smith       }
152947c6ae99SBarry Smith     }
153047c6ae99SBarry Smith     /* middle */
153147c6ae99SBarry Smith     for (i=down; i<up; i++) {
1532ce00eea3SSatish Balay       for (j=(IXs-Xs); j<(IXe-Xs); j++) {
1533ce00eea3SSatish Balay         idx[count++] = j + i*(Xe-Xs);
153447c6ae99SBarry Smith       }
153547c6ae99SBarry Smith     }
153647c6ae99SBarry Smith     /* top */
1537ce00eea3SSatish Balay     for (i=up; i<up+IYe-ye; i++) {
1538ce00eea3SSatish Balay       for (j=left; j<right; j++) {
1539ce00eea3SSatish Balay         idx[count++] = j + i*(Xe-Xs);
154047c6ae99SBarry Smith       }
154147c6ae99SBarry Smith     }
154247c6ae99SBarry Smith     ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr);
154347c6ae99SBarry Smith   }
154447c6ae99SBarry Smith 
154547c6ae99SBarry Smith 
154647c6ae99SBarry Smith   /* determine who lies on each side of us stored in    n6 n7 n8
154747c6ae99SBarry Smith                                                         n3    n5
154847c6ae99SBarry Smith                                                         n0 n1 n2
154947c6ae99SBarry Smith   */
155047c6ae99SBarry Smith 
155147c6ae99SBarry Smith   /* Assume the Non-Periodic Case */
155247c6ae99SBarry Smith   n1 = rank - m;
155347c6ae99SBarry Smith   if (rank % m) {
155447c6ae99SBarry Smith     n0 = n1 - 1;
155547c6ae99SBarry Smith   } else {
155647c6ae99SBarry Smith     n0 = -1;
155747c6ae99SBarry Smith   }
155847c6ae99SBarry Smith   if ((rank+1) % m) {
155947c6ae99SBarry Smith     n2 = n1 + 1;
156047c6ae99SBarry Smith     n5 = rank + 1;
156147c6ae99SBarry Smith     n8 = rank + m + 1; if (n8 >= m*n) n8 = -1;
156247c6ae99SBarry Smith   } else {
156347c6ae99SBarry Smith     n2 = -1; n5 = -1; n8 = -1;
156447c6ae99SBarry Smith   }
156547c6ae99SBarry Smith   if (rank % m) {
156647c6ae99SBarry Smith     n3 = rank - 1;
156747c6ae99SBarry Smith     n6 = n3 + m; if (n6 >= m*n) n6 = -1;
156847c6ae99SBarry Smith   } else {
156947c6ae99SBarry Smith     n3 = -1; n6 = -1;
157047c6ae99SBarry Smith   }
157147c6ae99SBarry Smith   n7 = rank + m; if (n7 >= m*n) n7 = -1;
157247c6ae99SBarry Smith 
15731321219cSEthan Coon   if (bx == DMDA_BOUNDARY_PERIODIC && by == DMDA_BOUNDARY_PERIODIC) {
157447c6ae99SBarry Smith   /* Modify for Periodic Cases */
157547c6ae99SBarry Smith     /* Handle all four corners */
157647c6ae99SBarry Smith     if ((n6 < 0) && (n7 < 0) && (n3 < 0)) n6 = m-1;
157747c6ae99SBarry Smith     if ((n8 < 0) && (n7 < 0) && (n5 < 0)) n8 = 0;
157847c6ae99SBarry Smith     if ((n2 < 0) && (n5 < 0) && (n1 < 0)) n2 = size-m;
157947c6ae99SBarry Smith     if ((n0 < 0) && (n3 < 0) && (n1 < 0)) n0 = size-1;
158047c6ae99SBarry Smith 
158147c6ae99SBarry Smith     /* Handle Top and Bottom Sides */
158247c6ae99SBarry Smith     if (n1 < 0) n1 = rank + m * (n-1);
158347c6ae99SBarry Smith     if (n7 < 0) n7 = rank - m * (n-1);
158447c6ae99SBarry Smith     if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1;
158547c6ae99SBarry Smith     if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1;
158647c6ae99SBarry Smith     if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1;
158747c6ae99SBarry Smith     if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1;
158847c6ae99SBarry Smith 
158947c6ae99SBarry Smith     /* Handle Left and Right Sides */
159047c6ae99SBarry Smith     if (n3 < 0) n3 = rank + (m-1);
159147c6ae99SBarry Smith     if (n5 < 0) n5 = rank - (m-1);
159247c6ae99SBarry Smith     if ((n1 >= 0) && (n0 < 0)) n0 = rank-1;
159347c6ae99SBarry Smith     if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1;
159447c6ae99SBarry Smith     if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1;
159547c6ae99SBarry Smith     if ((n7 >= 0) && (n8 < 0)) n8 = rank+1;
15961321219cSEthan Coon   } else if (by == DMDA_BOUNDARY_PERIODIC) {  /* Handle Top and Bottom Sides */
1597ce00eea3SSatish Balay     if (n1 < 0) n1 = rank + m * (n-1);
1598ce00eea3SSatish Balay     if (n7 < 0) n7 = rank - m * (n-1);
1599ce00eea3SSatish Balay     if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1;
1600ce00eea3SSatish Balay     if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1;
1601ce00eea3SSatish Balay     if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1;
1602ce00eea3SSatish Balay     if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1;
16031321219cSEthan Coon   } else if (bx == DMDA_BOUNDARY_PERIODIC) { /* Handle Left and Right Sides */
1604ce00eea3SSatish Balay     if (n3 < 0) n3 = rank + (m-1);
1605ce00eea3SSatish Balay     if (n5 < 0) n5 = rank - (m-1);
1606ce00eea3SSatish Balay     if ((n1 >= 0) && (n0 < 0)) n0 = rank-1;
1607ce00eea3SSatish Balay     if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1;
1608ce00eea3SSatish Balay     if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1;
1609ce00eea3SSatish Balay     if ((n7 >= 0) && (n8 < 0)) n8 = rank+1;
161047c6ae99SBarry Smith   }
1611ce00eea3SSatish Balay 
161247c6ae99SBarry Smith   ierr = PetscMalloc(9*sizeof(PetscInt),&dd->neighbors);CHKERRQ(ierr);
161347c6ae99SBarry Smith   dd->neighbors[0] = n0;
161447c6ae99SBarry Smith   dd->neighbors[1] = n1;
161547c6ae99SBarry Smith   dd->neighbors[2] = n2;
161647c6ae99SBarry Smith   dd->neighbors[3] = n3;
161747c6ae99SBarry Smith   dd->neighbors[4] = rank;
161847c6ae99SBarry Smith   dd->neighbors[5] = n5;
161947c6ae99SBarry Smith   dd->neighbors[6] = n6;
162047c6ae99SBarry Smith   dd->neighbors[7] = n7;
162147c6ae99SBarry Smith   dd->neighbors[8] = n8;
162247c6ae99SBarry Smith 
1623aa219208SBarry Smith   if (stencil_type == DMDA_STENCIL_STAR) {
162447c6ae99SBarry Smith     /* save corner processor numbers */
162547c6ae99SBarry Smith     sn0 = n0; sn2 = n2; sn6 = n6; sn8 = n8;
162647c6ae99SBarry Smith     n0 = n2 = n6 = n8 = -1;
162747c6ae99SBarry Smith   }
162847c6ae99SBarry Smith 
1629ce00eea3SSatish Balay   ierr = PetscMalloc((Xe-Xs)*(Ye-Ys)*sizeof(PetscInt),&idx);CHKERRQ(ierr);
1630ce00eea3SSatish Balay   ierr = PetscLogObjectMemory(da,(Xe-Xs)*(Ye-Ys)*sizeof(PetscInt));CHKERRQ(ierr);
163147c6ae99SBarry Smith 
1632ce00eea3SSatish Balay   nn = 0;
163347c6ae99SBarry Smith   xbase = bases[rank];
163447c6ae99SBarry Smith   for (i=1; i<=s_y; i++) {
163547c6ae99SBarry Smith     if (n0 >= 0) { /* left below */
1636ce00eea3SSatish Balay       x_t = lx[n0 % m];
163747c6ae99SBarry Smith       y_t = ly[(n0/m)];
163847c6ae99SBarry Smith       s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x;
163947c6ae99SBarry Smith       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
164047c6ae99SBarry Smith     }
164147c6ae99SBarry Smith     if (n1 >= 0) { /* directly below */
164247c6ae99SBarry Smith       x_t = x;
164347c6ae99SBarry Smith       y_t = ly[(n1/m)];
164447c6ae99SBarry Smith       s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t;
164547c6ae99SBarry Smith       for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
164647c6ae99SBarry Smith     }
164747c6ae99SBarry Smith     if (n2 >= 0) { /* right below */
1648ce00eea3SSatish Balay       x_t = lx[n2 % m];
164947c6ae99SBarry Smith       y_t = ly[(n2/m)];
165047c6ae99SBarry Smith       s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t;
165147c6ae99SBarry Smith       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
165247c6ae99SBarry Smith     }
165347c6ae99SBarry Smith   }
165447c6ae99SBarry Smith 
165547c6ae99SBarry Smith   for (i=0; i<y; i++) {
165647c6ae99SBarry Smith     if (n3 >= 0) { /* directly left */
1657ce00eea3SSatish Balay       x_t = lx[n3 % m];
165847c6ae99SBarry Smith       /* y_t = y; */
165947c6ae99SBarry Smith       s_t = bases[n3] + (i+1)*x_t - s_x;
166047c6ae99SBarry Smith       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
166147c6ae99SBarry Smith     }
166247c6ae99SBarry Smith 
166347c6ae99SBarry Smith     for (j=0; j<x; j++) { idx[nn++] = xbase++; } /* interior */
166447c6ae99SBarry Smith 
166547c6ae99SBarry Smith     if (n5 >= 0) { /* directly right */
1666ce00eea3SSatish Balay       x_t = lx[n5 % m];
166747c6ae99SBarry Smith       /* y_t = y; */
166847c6ae99SBarry Smith       s_t = bases[n5] + (i)*x_t;
166947c6ae99SBarry Smith       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
167047c6ae99SBarry Smith     }
167147c6ae99SBarry Smith   }
167247c6ae99SBarry Smith 
167347c6ae99SBarry Smith   for (i=1; i<=s_y; i++) {
167447c6ae99SBarry Smith     if (n6 >= 0) { /* left above */
1675ce00eea3SSatish Balay       x_t = lx[n6 % m];
167647c6ae99SBarry Smith       /* y_t = ly[(n6/m)]; */
167747c6ae99SBarry Smith       s_t = bases[n6] + (i)*x_t - s_x;
167847c6ae99SBarry Smith       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
167947c6ae99SBarry Smith     }
168047c6ae99SBarry Smith     if (n7 >= 0) { /* directly above */
168147c6ae99SBarry Smith       x_t = x;
168247c6ae99SBarry Smith       /* y_t = ly[(n7/m)]; */
168347c6ae99SBarry Smith       s_t = bases[n7] + (i-1)*x_t;
168447c6ae99SBarry Smith       for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
168547c6ae99SBarry Smith     }
168647c6ae99SBarry Smith     if (n8 >= 0) { /* right above */
1687ce00eea3SSatish Balay       x_t = lx[n8 % m];
168847c6ae99SBarry Smith       /* y_t = ly[(n8/m)]; */
168947c6ae99SBarry Smith       s_t = bases[n8] + (i-1)*x_t;
169047c6ae99SBarry Smith       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
169147c6ae99SBarry Smith     }
169247c6ae99SBarry Smith   }
169347c6ae99SBarry Smith 
1694ce00eea3SSatish Balay   ierr = ISCreateBlock(comm,dof,nn,idx,PETSC_COPY_VALUES,&from);CHKERRQ(ierr);
169547c6ae99SBarry Smith   ierr = VecScatterCreate(global,from,local,to,&gtol);CHKERRQ(ierr);
169647c6ae99SBarry Smith   ierr = PetscLogObjectParent(da,gtol);CHKERRQ(ierr);
1697fcfd50ebSBarry Smith   ierr = ISDestroy(&to);CHKERRQ(ierr);
1698fcfd50ebSBarry Smith   ierr = ISDestroy(&from);CHKERRQ(ierr);
169947c6ae99SBarry Smith 
1700aa219208SBarry Smith   if (stencil_type == DMDA_STENCIL_STAR) {
1701ce00eea3SSatish Balay     n0 = sn0; n2 = sn2; n6 = sn6; n8 = sn8;
1702ce00eea3SSatish Balay   }
1703ce00eea3SSatish Balay 
1704ce00eea3SSatish Balay   if ((stencil_type == DMDA_STENCIL_STAR) ||
17051321219cSEthan Coon       (bx && bx != DMDA_BOUNDARY_PERIODIC) ||
17061321219cSEthan Coon       (by && by != DMDA_BOUNDARY_PERIODIC)) {
170747c6ae99SBarry Smith     /*
170847c6ae99SBarry Smith         Recompute the local to global mappings, this time keeping the
1709ce00eea3SSatish Balay       information about the cross corner processor numbers and any ghosted
1710ce00eea3SSatish Balay       but not periodic indices.
171147c6ae99SBarry Smith     */
171247c6ae99SBarry Smith     nn = 0;
171347c6ae99SBarry Smith     xbase = bases[rank];
171447c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
171547c6ae99SBarry Smith       if (n0 >= 0) { /* left below */
1716ce00eea3SSatish Balay         x_t = lx[n0 % m];
171747c6ae99SBarry Smith         y_t = ly[(n0/m)];
171847c6ae99SBarry Smith         s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x;
171947c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1720ce00eea3SSatish Balay       } else if (xs-Xs > 0 && ys-Ys > 0) {
1721ce00eea3SSatish Balay         for (j=0; j<s_x; j++) { idx[nn++] = -1;}
172247c6ae99SBarry Smith       }
172347c6ae99SBarry Smith       if (n1 >= 0) { /* directly below */
172447c6ae99SBarry Smith         x_t = x;
172547c6ae99SBarry Smith         y_t = ly[(n1/m)];
172647c6ae99SBarry Smith         s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t;
172747c6ae99SBarry Smith         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1728ce00eea3SSatish Balay       } else if (ys-Ys > 0) {
1729ce00eea3SSatish Balay         for (j=0; j<x; j++) { idx[nn++] = -1;}
173047c6ae99SBarry Smith       }
173147c6ae99SBarry Smith       if (n2 >= 0) { /* right below */
1732ce00eea3SSatish Balay         x_t = lx[n2 % m];
173347c6ae99SBarry Smith         y_t = ly[(n2/m)];
173447c6ae99SBarry Smith         s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t;
173547c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1736ce00eea3SSatish Balay       } else if (Xe-xe> 0 && ys-Ys > 0) {
1737ce00eea3SSatish Balay         for (j=0; j<s_x; j++) { idx[nn++] = -1;}
173847c6ae99SBarry Smith       }
173947c6ae99SBarry Smith     }
174047c6ae99SBarry Smith 
174147c6ae99SBarry Smith     for (i=0; i<y; i++) {
174247c6ae99SBarry Smith       if (n3 >= 0) { /* directly left */
1743ce00eea3SSatish Balay         x_t = lx[n3 % m];
174447c6ae99SBarry Smith         /* y_t = y; */
174547c6ae99SBarry Smith         s_t = bases[n3] + (i+1)*x_t - s_x;
174647c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1747ce00eea3SSatish Balay       } else if (xs-Xs > 0) {
1748ce00eea3SSatish Balay         for (j=0; j<s_x; j++) { idx[nn++] = -1;}
174947c6ae99SBarry Smith       }
175047c6ae99SBarry Smith 
175147c6ae99SBarry Smith       for (j=0; j<x; j++) { idx[nn++] = xbase++; } /* interior */
175247c6ae99SBarry Smith 
175347c6ae99SBarry Smith       if (n5 >= 0) { /* directly right */
1754ce00eea3SSatish Balay         x_t = lx[n5 % m];
175547c6ae99SBarry Smith         /* y_t = y; */
175647c6ae99SBarry Smith         s_t = bases[n5] + (i)*x_t;
175747c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1758ce00eea3SSatish Balay       } else if (Xe-xe > 0) {
1759ce00eea3SSatish Balay         for (j=0; j<s_x; j++) { idx[nn++] = -1;}
176047c6ae99SBarry Smith       }
176147c6ae99SBarry Smith     }
176247c6ae99SBarry Smith 
176347c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
176447c6ae99SBarry Smith       if (n6 >= 0) { /* left above */
1765ce00eea3SSatish Balay         x_t = lx[n6 % m];
176647c6ae99SBarry Smith         /* y_t = ly[(n6/m)]; */
176747c6ae99SBarry Smith         s_t = bases[n6] + (i)*x_t - s_x;
176847c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1769ce00eea3SSatish Balay       } else if (xs-Xs > 0 && Ye-ye > 0) {
1770ce00eea3SSatish Balay         for (j=0; j<s_x; j++) { idx[nn++] = -1;}
177147c6ae99SBarry Smith       }
177247c6ae99SBarry Smith       if (n7 >= 0) { /* directly above */
177347c6ae99SBarry Smith         x_t = x;
177447c6ae99SBarry Smith         /* y_t = ly[(n7/m)]; */
177547c6ae99SBarry Smith         s_t = bases[n7] + (i-1)*x_t;
177647c6ae99SBarry Smith         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1777ce00eea3SSatish Balay       } else if (Ye-ye > 0) {
1778ce00eea3SSatish Balay         for (j=0; j<x; j++) { idx[nn++] = -1;}
177947c6ae99SBarry Smith       }
178047c6ae99SBarry Smith       if (n8 >= 0) { /* right above */
1781ce00eea3SSatish Balay         x_t = lx[n8 % m];
178247c6ae99SBarry Smith         /* y_t = ly[(n8/m)]; */
178347c6ae99SBarry Smith         s_t = bases[n8] + (i-1)*x_t;
178447c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1785ce00eea3SSatish Balay       } else if (Xe-xe > 0 && Ye-ye > 0) {
1786ce00eea3SSatish Balay         for (j=0; j<s_x; j++) { idx[nn++] = -1;}
178747c6ae99SBarry Smith       }
178847c6ae99SBarry Smith     }
178947c6ae99SBarry Smith   }
1790ce00eea3SSatish Balay   /*
1791ce00eea3SSatish Balay      Set the local to global ordering in the global vector, this allows use
1792ce00eea3SSatish Balay      of VecSetValuesLocal().
1793ce00eea3SSatish Balay   */
1794ce00eea3SSatish Balay   ierr = ISCreateBlock(comm,dof,nn,idx,PETSC_OWN_POINTER,&ltogis);CHKERRQ(ierr);
1795ce00eea3SSatish Balay   ierr = PetscMalloc(nn*dof*sizeof(PetscInt),&idx_cpy);CHKERRQ(ierr);
1796db87c5ecSEthan Coon   ierr = PetscLogObjectMemory(da,nn*dof*sizeof(PetscInt));CHKERRQ(ierr);
1797ce00eea3SSatish Balay   ierr = ISGetIndices(ltogis, &idx_full);
1798ce00eea3SSatish Balay   ierr = PetscMemcpy(idx_cpy,idx_full,nn*dof*sizeof(PetscInt));CHKERRQ(ierr);
1799ce00eea3SSatish Balay   ierr = ISRestoreIndices(ltogis, &idx_full);
1800ce00eea3SSatish Balay   ierr = ISLocalToGlobalMappingCreateIS(ltogis,&da->ltogmap);CHKERRQ(ierr);
1801ce00eea3SSatish Balay   ierr = PetscLogObjectParent(da,da->ltogmap);CHKERRQ(ierr);
1802fcfd50ebSBarry Smith   ierr = ISDestroy(&ltogis);CHKERRQ(ierr);
1803ce00eea3SSatish Balay   ierr = ISLocalToGlobalMappingBlock(da->ltogmap,dd->w,&da->ltogmapb);CHKERRQ(ierr);
1804ce00eea3SSatish Balay   ierr = PetscLogObjectParent(da,da->ltogmap);CHKERRQ(ierr);
180547c6ae99SBarry Smith 
1806ce00eea3SSatish Balay   ierr = PetscFree2(bases,ldims);CHKERRQ(ierr);
180747c6ae99SBarry Smith   dd->m  = m;  dd->n  = n;
1808ce00eea3SSatish Balay   /* note petsc expects xs/xe/Xs/Xe to be multiplied by #dofs in many places */
1809ce00eea3SSatish Balay   dd->xs = xs*dof; dd->xe = xe*dof; dd->ys = ys; dd->ye = ye; dd->zs = 0; dd->ze = 1;
1810ce00eea3SSatish Balay   dd->Xs = Xs*dof; dd->Xe = Xe*dof; dd->Ys = Ys; dd->Ye = Ye; dd->Zs = 0; dd->Ze = 1;
181147c6ae99SBarry Smith 
1812fcfd50ebSBarry Smith   ierr = VecDestroy(&local);CHKERRQ(ierr);
1813fcfd50ebSBarry Smith   ierr = VecDestroy(&global);CHKERRQ(ierr);
181447c6ae99SBarry Smith 
181547c6ae99SBarry Smith   dd->gtol      = gtol;
181647c6ae99SBarry Smith   dd->ltog      = ltog;
1817ce00eea3SSatish Balay   dd->idx       = idx_cpy;
1818ce00eea3SSatish Balay   dd->Nl        = nn*dof;
181947c6ae99SBarry Smith   dd->base      = base;
18209a42bb27SBarry Smith   da->ops->view = DMView_DA_2d;
182147c6ae99SBarry Smith   dd->ltol = PETSC_NULL;
182247c6ae99SBarry Smith   dd->ao   = PETSC_NULL;
182347c6ae99SBarry Smith 
182447c6ae99SBarry Smith   PetscFunctionReturn(0);
182547c6ae99SBarry Smith }
182647c6ae99SBarry Smith 
182747c6ae99SBarry Smith #undef __FUNCT__
1828aa219208SBarry Smith #define __FUNCT__ "DMDACreate2d"
182947c6ae99SBarry Smith /*@C
1830aa219208SBarry Smith    DMDACreate2d -  Creates an object that will manage the communication of  two-dimensional
183147c6ae99SBarry Smith    regular array data that is distributed across some processors.
183247c6ae99SBarry Smith 
183347c6ae99SBarry Smith    Collective on MPI_Comm
183447c6ae99SBarry Smith 
183547c6ae99SBarry Smith    Input Parameters:
183647c6ae99SBarry Smith +  comm - MPI communicator
18371321219cSEthan Coon .  bx,by - type of ghost nodes the array have.
18381321219cSEthan Coon          Use one of DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_GHOSTED, DMDA_BOUNDARY_PERIODIC.
1839aa219208SBarry Smith .  stencil_type - stencil type.  Use either DMDA_STENCIL_BOX or DMDA_STENCIL_STAR.
184047c6ae99SBarry Smith .  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
184147c6ae99SBarry Smith             from the command line with -da_grid_x <M> -da_grid_y <N>)
184247c6ae99SBarry Smith .  m,n - corresponding number of processors in each dimension
184347c6ae99SBarry Smith          (or PETSC_DECIDE to have calculated)
184447c6ae99SBarry Smith .  dof - number of degrees of freedom per node
184547c6ae99SBarry Smith .  s - stencil width
184647c6ae99SBarry Smith -  lx, ly - arrays containing the number of nodes in each cell along
184747c6ae99SBarry Smith            the x and y coordinates, or PETSC_NULL. If non-null, these
184847c6ae99SBarry Smith            must be of length as m and n, and the corresponding
184947c6ae99SBarry Smith            m and n cannot be PETSC_DECIDE. The sum of the lx[] entries
185047c6ae99SBarry Smith            must be M, and the sum of the ly[] entries must be N.
185147c6ae99SBarry Smith 
185247c6ae99SBarry Smith    Output Parameter:
185347c6ae99SBarry Smith .  da - the resulting distributed array object
185447c6ae99SBarry Smith 
185547c6ae99SBarry Smith    Options Database Key:
1856aa219208SBarry Smith +  -da_view - Calls DMView() at the conclusion of DMDACreate2d()
185747c6ae99SBarry Smith .  -da_grid_x <nx> - number of grid points in x direction, if M < 0
185847c6ae99SBarry Smith .  -da_grid_y <ny> - number of grid points in y direction, if N < 0
185947c6ae99SBarry Smith .  -da_processors_x <nx> - number of processors in x direction
186047c6ae99SBarry Smith .  -da_processors_y <ny> - number of processors in y direction
1861e0f5d30fSBarry Smith .  -da_refine_x <rx> - refinement ratio in x direction
1862e0f5d30fSBarry Smith .  -da_refine_y <ry> - refinement ratio in y direction
1863e0f5d30fSBarry Smith -  -da_refine <n> - refine the DMDA n times before creating, if M or N < 0
1864e0f5d30fSBarry Smith 
186547c6ae99SBarry Smith 
186647c6ae99SBarry Smith    Level: beginner
186747c6ae99SBarry Smith 
186847c6ae99SBarry Smith    Notes:
1869aa219208SBarry Smith    The stencil type DMDA_STENCIL_STAR with width 1 corresponds to the
1870aa219208SBarry Smith    standard 5-pt stencil, while DMDA_STENCIL_BOX with width 1 denotes
187147c6ae99SBarry Smith    the standard 9-pt stencil.
187247c6ae99SBarry Smith 
1873aa219208SBarry Smith    The array data itself is NOT stored in the DMDA, it is stored in Vec objects;
1874564755cdSBarry Smith    The appropriate vector objects can be obtained with calls to DMCreateGlobalVector()
1875564755cdSBarry Smith    and DMCreateLocalVector() and calls to VecDuplicate() if more are needed.
187647c6ae99SBarry Smith 
187747c6ae99SBarry Smith .keywords: distributed array, create, two-dimensional
187847c6ae99SBarry Smith 
1879aa219208SBarry Smith .seealso: DMDestroy(), DMView(), DMDACreate1d(), DMDACreate3d(), DMGlobalToLocalBegin(), DMDAGetRefinementFactor(),
1880aa219208SBarry Smith           DMGlobalToLocalEnd(), DMLocalToGlobalBegin(), DMDALocalToLocalBegin(), DMDALocalToLocalEnd(), DMDASetRefinementFactor(),
1881*d461ba97SBarry Smith           DMDAGetInfo(), DMCreateGlobalVector(), DMCreateLocalVector(), DMDACreateNaturalVector(), DMLoad(), DMDAGetOwnershipRanges()
188247c6ae99SBarry Smith 
188347c6ae99SBarry Smith @*/
18841321219cSEthan Coon PetscErrorCode  DMDACreate2d(MPI_Comm comm,DMDABoundaryType bx,DMDABoundaryType by,DMDAStencilType stencil_type,
18859a42bb27SBarry Smith                           PetscInt M,PetscInt N,PetscInt m,PetscInt n,PetscInt dof,PetscInt s,const PetscInt lx[],const PetscInt ly[],DM *da)
188647c6ae99SBarry Smith {
188747c6ae99SBarry Smith   PetscErrorCode ierr;
188847c6ae99SBarry Smith 
188947c6ae99SBarry Smith   PetscFunctionBegin;
1890aa219208SBarry Smith   ierr = DMDACreate(comm, da);CHKERRQ(ierr);
1891aa219208SBarry Smith   ierr = DMDASetDim(*da, 2);CHKERRQ(ierr);
1892aa219208SBarry Smith   ierr = DMDASetSizes(*da, M, N, 1);CHKERRQ(ierr);
1893aa219208SBarry Smith   ierr = DMDASetNumProcs(*da, m, n, PETSC_DECIDE);CHKERRQ(ierr);
1894755f258dSLisandro Dalcin   ierr = DMDASetBoundaryType(*da, bx, by, DMDA_BOUNDARY_NONE);CHKERRQ(ierr);
1895aa219208SBarry Smith   ierr = DMDASetDof(*da, dof);CHKERRQ(ierr);
1896aa219208SBarry Smith   ierr = DMDASetStencilType(*da, stencil_type);CHKERRQ(ierr);
1897aa219208SBarry Smith   ierr = DMDASetStencilWidth(*da, s);CHKERRQ(ierr);
1898aa219208SBarry Smith   ierr = DMDASetOwnershipRanges(*da, lx, ly, PETSC_NULL);CHKERRQ(ierr);
189947c6ae99SBarry Smith   /* This violates the behavior for other classes, but right now users expect negative dimensions to be handled this way */
19009a42bb27SBarry Smith   ierr = DMSetFromOptions(*da);CHKERRQ(ierr);
19019a42bb27SBarry Smith   ierr = DMSetUp(*da);CHKERRQ(ierr);
19027242296bSJed Brown   ierr = DMView_DA_Private(*da);CHKERRQ(ierr);
190347c6ae99SBarry Smith   PetscFunctionReturn(0);
190447c6ae99SBarry Smith }
1905