xref: /petsc/src/dm/impls/da/da2.c (revision fcfd50eb5fbe12291db06fc1b33ad004ead64926)
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);
3247c6ae99SBarry 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);
3347c6ae99SBarry 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);
3447c6ae99SBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
353da9ae13SJed Brown     } else {
363da9ae13SJed Brown       ierr = DMView_DA_VTK(da,viewer);CHKERRQ(ierr);
3747c6ae99SBarry Smith     }
3847c6ae99SBarry Smith   } else if (isdraw) {
3947c6ae99SBarry Smith     PetscDraw  draw;
4047c6ae99SBarry Smith     double     ymin = -1*dd->s-1,ymax = dd->N+dd->s;
4147c6ae99SBarry Smith     double     xmin = -1*dd->s-1,xmax = dd->M+dd->s;
4247c6ae99SBarry Smith     double     x,y;
4347c6ae99SBarry Smith     PetscInt   base,*idx;
4447c6ae99SBarry Smith     char       node[10];
4547c6ae99SBarry Smith     PetscBool  isnull;
4647c6ae99SBarry Smith 
4747c6ae99SBarry Smith     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
4847c6ae99SBarry Smith     ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
4947c6ae99SBarry Smith     if (!dd->coordinates) {
5047c6ae99SBarry Smith       ierr = PetscDrawSetCoordinates(draw,xmin,ymin,xmax,ymax);CHKERRQ(ierr);
5147c6ae99SBarry Smith     }
5247c6ae99SBarry Smith     ierr = PetscDrawSynchronizedClear(draw);CHKERRQ(ierr);
5347c6ae99SBarry Smith 
5447c6ae99SBarry Smith     /* first processor draw all node lines */
5547c6ae99SBarry Smith     if (!rank) {
5647c6ae99SBarry Smith       ymin = 0.0; ymax = dd->N - 1;
5747c6ae99SBarry Smith       for (xmin=0; xmin<dd->M; xmin++) {
5847c6ae99SBarry Smith         ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_BLACK);CHKERRQ(ierr);
5947c6ae99SBarry Smith       }
6047c6ae99SBarry Smith       xmin = 0.0; xmax = dd->M - 1;
6147c6ae99SBarry Smith       for (ymin=0; ymin<dd->N; ymin++) {
6247c6ae99SBarry Smith         ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_BLACK);CHKERRQ(ierr);
6347c6ae99SBarry Smith       }
6447c6ae99SBarry Smith     }
6547c6ae99SBarry Smith     ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr);
6647c6ae99SBarry Smith     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
6747c6ae99SBarry Smith 
6847c6ae99SBarry Smith     /* draw my box */
6947c6ae99SBarry Smith     ymin = dd->ys; ymax = dd->ye - 1; xmin = dd->xs/dd->w;
7047c6ae99SBarry Smith     xmax =(dd->xe-1)/dd->w;
7147c6ae99SBarry Smith     ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_RED);CHKERRQ(ierr);
7247c6ae99SBarry Smith     ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_RED);CHKERRQ(ierr);
7347c6ae99SBarry Smith     ierr = PetscDrawLine(draw,xmin,ymax,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr);
7447c6ae99SBarry Smith     ierr = PetscDrawLine(draw,xmax,ymin,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr);
7547c6ae99SBarry Smith 
7647c6ae99SBarry Smith     /* put in numbers */
7747c6ae99SBarry Smith     base = (dd->base)/dd->w;
7847c6ae99SBarry Smith     for (y=ymin; y<=ymax; y++) {
7947c6ae99SBarry Smith       for (x=xmin; x<=xmax; x++) {
8047c6ae99SBarry Smith         sprintf(node,"%d",(int)base++);
8147c6ae99SBarry Smith         ierr = PetscDrawString(draw,x,y,PETSC_DRAW_BLACK,node);CHKERRQ(ierr);
8247c6ae99SBarry Smith       }
8347c6ae99SBarry Smith     }
8447c6ae99SBarry Smith 
8547c6ae99SBarry Smith     ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr);
8647c6ae99SBarry Smith     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
8747c6ae99SBarry Smith     /* overlay ghost numbers, useful for error checking */
8847c6ae99SBarry Smith     /* put in numbers */
8947c6ae99SBarry Smith 
9047c6ae99SBarry Smith     base = 0; idx = dd->idx;
9147c6ae99SBarry Smith     ymin = dd->Ys; ymax = dd->Ye; xmin = dd->Xs; xmax = dd->Xe;
9247c6ae99SBarry Smith     for (y=ymin; y<ymax; y++) {
9347c6ae99SBarry Smith       for (x=xmin; x<xmax; x++) {
9447c6ae99SBarry Smith         if ((base % dd->w) == 0) {
9547c6ae99SBarry Smith           sprintf(node,"%d",(int)(idx[base]/dd->w));
9647c6ae99SBarry Smith           ierr = PetscDrawString(draw,x/dd->w,y,PETSC_DRAW_BLUE,node);CHKERRQ(ierr);
9747c6ae99SBarry Smith         }
9847c6ae99SBarry Smith         base++;
9947c6ae99SBarry Smith       }
10047c6ae99SBarry Smith     }
10147c6ae99SBarry Smith     ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr);
10247c6ae99SBarry Smith     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
1039a42bb27SBarry Smith   } else if (isbinary){
1049a42bb27SBarry Smith     ierr = DMView_DA_Binary(da,viewer);CHKERRQ(ierr);
1059a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
1069a42bb27SBarry Smith   } else if (ismatlab) {
1079a42bb27SBarry Smith     ierr = DMView_DA_Matlab(da,viewer);CHKERRQ(ierr);
1089a42bb27SBarry Smith #endif
109aa219208SBarry Smith   } else SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for DMDA 1d",((PetscObject)viewer)->type_name);
11047c6ae99SBarry Smith   PetscFunctionReturn(0);
11147c6ae99SBarry Smith }
11247c6ae99SBarry Smith 
11347c6ae99SBarry Smith /*
11447c6ae99SBarry Smith       M is number of grid points
11547c6ae99SBarry Smith       m is number of processors
11647c6ae99SBarry Smith 
11747c6ae99SBarry Smith */
11847c6ae99SBarry Smith #undef __FUNCT__
119aa219208SBarry Smith #define __FUNCT__ "DMDASplitComm2d"
1207087cfbeSBarry Smith PetscErrorCode  DMDASplitComm2d(MPI_Comm comm,PetscInt M,PetscInt N,PetscInt sw,MPI_Comm *outcomm)
12147c6ae99SBarry Smith {
12247c6ae99SBarry Smith   PetscErrorCode ierr;
12347c6ae99SBarry Smith   PetscInt       m,n = 0,x = 0,y = 0;
12447c6ae99SBarry Smith   PetscMPIInt    size,csize,rank;
12547c6ae99SBarry Smith 
12647c6ae99SBarry Smith   PetscFunctionBegin;
12747c6ae99SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
12847c6ae99SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
12947c6ae99SBarry Smith 
13047c6ae99SBarry Smith   csize = 4*size;
13147c6ae99SBarry Smith   do {
13247c6ae99SBarry 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);
13347c6ae99SBarry Smith     csize   = csize/4;
13447c6ae99SBarry Smith 
13547c6ae99SBarry Smith     m = (PetscInt)(0.5 + sqrt(((double)M)*((double)csize)/((double)N)));
13647c6ae99SBarry Smith     if (!m) m = 1;
13747c6ae99SBarry Smith     while (m > 0) {
13847c6ae99SBarry Smith       n = csize/m;
13947c6ae99SBarry Smith       if (m*n == csize) break;
14047c6ae99SBarry Smith       m--;
14147c6ae99SBarry Smith     }
14247c6ae99SBarry Smith     if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;}
14347c6ae99SBarry Smith 
14447c6ae99SBarry Smith     x = M/m + ((M % m) > ((csize-1) % m));
14547c6ae99SBarry Smith     y = (N + (csize-1)/m)/n;
14647c6ae99SBarry Smith   } while ((x < 4 || y < 4) && csize > 1);
14747c6ae99SBarry Smith   if (size != csize) {
14847c6ae99SBarry Smith     MPI_Group    entire_group,sub_group;
14947c6ae99SBarry Smith     PetscMPIInt  i,*groupies;
15047c6ae99SBarry Smith 
15147c6ae99SBarry Smith     ierr = MPI_Comm_group(comm,&entire_group);CHKERRQ(ierr);
15247c6ae99SBarry Smith     ierr = PetscMalloc(csize*sizeof(PetscInt),&groupies);CHKERRQ(ierr);
15347c6ae99SBarry Smith     for (i=0; i<csize; i++) {
15447c6ae99SBarry Smith       groupies[i] = (rank/csize)*csize + i;
15547c6ae99SBarry Smith     }
15647c6ae99SBarry Smith     ierr = MPI_Group_incl(entire_group,csize,groupies,&sub_group);CHKERRQ(ierr);
15747c6ae99SBarry Smith     ierr = PetscFree(groupies);CHKERRQ(ierr);
15847c6ae99SBarry Smith     ierr = MPI_Comm_create(comm,sub_group,outcomm);CHKERRQ(ierr);
15947c6ae99SBarry Smith     ierr = MPI_Group_free(&entire_group);CHKERRQ(ierr);
16047c6ae99SBarry Smith     ierr = MPI_Group_free(&sub_group);CHKERRQ(ierr);
161aa219208SBarry Smith     ierr = PetscInfo1(0,"DMDASplitComm2d:Creating redundant coarse problems of size %d\n",csize);CHKERRQ(ierr);
16247c6ae99SBarry Smith   } else {
16347c6ae99SBarry Smith     *outcomm = comm;
16447c6ae99SBarry Smith   }
16547c6ae99SBarry Smith   PetscFunctionReturn(0);
16647c6ae99SBarry Smith }
16747c6ae99SBarry Smith 
16847c6ae99SBarry Smith #undef __FUNCT__
169aa219208SBarry Smith #define __FUNCT__ "DMDASetLocalFunction"
17047c6ae99SBarry Smith /*@C
171aa219208SBarry Smith        DMDASetLocalFunction - Caches in a DM a local function.
17247c6ae99SBarry Smith 
173aa219208SBarry Smith    Logically Collective on DMDA
17447c6ae99SBarry Smith 
17547c6ae99SBarry Smith    Input Parameter:
17647c6ae99SBarry Smith +  da - initial distributed array
17747c6ae99SBarry Smith -  lf - the local function
17847c6ae99SBarry Smith 
17947c6ae99SBarry Smith    Level: intermediate
18047c6ae99SBarry Smith 
18147c6ae99SBarry Smith    Notes: The routine SNESDAFormFunction() uses this the cached function to evaluate the user provided function.
18247c6ae99SBarry Smith 
18347c6ae99SBarry Smith .keywords:  distributed array, refine
18447c6ae99SBarry Smith 
185aa219208SBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalFunctioni()
18647c6ae99SBarry Smith @*/
1877087cfbeSBarry Smith PetscErrorCode  DMDASetLocalFunction(DM da,DMDALocalFunction1 lf)
18847c6ae99SBarry Smith {
18947c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
19047c6ae99SBarry Smith   PetscFunctionBegin;
19147c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
19247c6ae99SBarry Smith   dd->lf    = lf;
19347c6ae99SBarry Smith   PetscFunctionReturn(0);
19447c6ae99SBarry Smith }
19547c6ae99SBarry Smith 
19647c6ae99SBarry Smith #undef __FUNCT__
197aa219208SBarry Smith #define __FUNCT__ "DMDASetLocalFunctioni"
19847c6ae99SBarry Smith /*@C
199aa219208SBarry Smith        DMDASetLocalFunctioni - Caches in a DM a local function that evaluates a single component
20047c6ae99SBarry Smith 
201aa219208SBarry Smith    Logically Collective on DMDA
20247c6ae99SBarry Smith 
20347c6ae99SBarry Smith    Input Parameter:
20447c6ae99SBarry Smith +  da - initial distributed array
20547c6ae99SBarry Smith -  lfi - the local function
20647c6ae99SBarry Smith 
20747c6ae99SBarry Smith    Level: intermediate
20847c6ae99SBarry Smith 
20947c6ae99SBarry Smith .keywords:  distributed array, refine
21047c6ae99SBarry Smith 
211aa219208SBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalFunction()
21247c6ae99SBarry Smith @*/
2137087cfbeSBarry Smith PetscErrorCode  DMDASetLocalFunctioni(DM da,PetscErrorCode (*lfi)(DMDALocalInfo*,MatStencil*,void*,PetscScalar*,void*))
21447c6ae99SBarry Smith {
21547c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
21647c6ae99SBarry Smith   PetscFunctionBegin;
21747c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
21847c6ae99SBarry Smith   dd->lfi = lfi;
21947c6ae99SBarry Smith   PetscFunctionReturn(0);
22047c6ae99SBarry Smith }
22147c6ae99SBarry Smith 
22247c6ae99SBarry Smith #undef __FUNCT__
223aa219208SBarry Smith #define __FUNCT__ "DMDASetLocalFunctionib"
22447c6ae99SBarry Smith /*@C
225aa219208SBarry Smith        DMDASetLocalFunctionib - Caches in a DM a block local function that evaluates a single component
22647c6ae99SBarry Smith 
227aa219208SBarry Smith    Logically Collective on DMDA
22847c6ae99SBarry Smith 
22947c6ae99SBarry Smith    Input Parameter:
23047c6ae99SBarry Smith +  da - initial distributed array
23147c6ae99SBarry Smith -  lfi - the local function
23247c6ae99SBarry Smith 
23347c6ae99SBarry Smith    Level: intermediate
23447c6ae99SBarry Smith 
23547c6ae99SBarry Smith .keywords:  distributed array, refine
23647c6ae99SBarry Smith 
237aa219208SBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalFunction()
23847c6ae99SBarry Smith @*/
2397087cfbeSBarry Smith PetscErrorCode  DMDASetLocalFunctionib(DM da,PetscErrorCode (*lfi)(DMDALocalInfo*,MatStencil*,void*,PetscScalar*,void*))
24047c6ae99SBarry Smith {
24147c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
24247c6ae99SBarry Smith   PetscFunctionBegin;
24347c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
24447c6ae99SBarry Smith   dd->lfib = lfi;
24547c6ae99SBarry Smith   PetscFunctionReturn(0);
24647c6ae99SBarry Smith }
24747c6ae99SBarry Smith 
24847c6ae99SBarry Smith #undef __FUNCT__
249aa219208SBarry Smith #define __FUNCT__ "DMDASetLocalAdicFunction_Private"
250aa219208SBarry Smith PetscErrorCode DMDASetLocalAdicFunction_Private(DM da,DMDALocalFunction1 ad_lf)
25147c6ae99SBarry Smith {
25247c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
25347c6ae99SBarry Smith   PetscFunctionBegin;
25447c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
25547c6ae99SBarry Smith   dd->adic_lf = ad_lf;
25647c6ae99SBarry Smith   PetscFunctionReturn(0);
25747c6ae99SBarry Smith }
25847c6ae99SBarry Smith 
25947c6ae99SBarry Smith /*MC
260aa219208SBarry Smith        DMDASetLocalAdicFunctioni - Caches in a DM a local functioni computed by ADIC/ADIFOR
26147c6ae99SBarry Smith 
26247c6ae99SBarry Smith    Synopsis:
263aa219208SBarry Smith    PetscErrorCode DMDASetLocalAdicFunctioni(DM da,PetscInt (ad_lf*)(DMDALocalInfo*,MatStencil*,void*,void*,void*)
26447c6ae99SBarry Smith 
265aa219208SBarry Smith    Logically Collective on DMDA
26647c6ae99SBarry Smith 
26747c6ae99SBarry Smith    Input Parameter:
26847c6ae99SBarry Smith +  da - initial distributed array
26947c6ae99SBarry Smith -  ad_lfi - the local function as computed by ADIC/ADIFOR
27047c6ae99SBarry Smith 
27147c6ae99SBarry Smith    Level: intermediate
27247c6ae99SBarry Smith 
27347c6ae99SBarry Smith .keywords:  distributed array, refine
27447c6ae99SBarry Smith 
275aa219208SBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalFunction(),
276aa219208SBarry Smith           DMDASetLocalJacobian(), DMDASetLocalFunctioni()
27747c6ae99SBarry Smith M*/
27847c6ae99SBarry Smith 
27947c6ae99SBarry Smith #undef __FUNCT__
280aa219208SBarry Smith #define __FUNCT__ "DMDASetLocalAdicFunctioni_Private"
281aa219208SBarry Smith PetscErrorCode DMDASetLocalAdicFunctioni_Private(DM da,PetscErrorCode (*ad_lfi)(DMDALocalInfo*,MatStencil*,void*,void*,void*))
28247c6ae99SBarry Smith {
28347c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
28447c6ae99SBarry Smith   PetscFunctionBegin;
28547c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
28647c6ae99SBarry Smith   dd->adic_lfi = ad_lfi;
28747c6ae99SBarry Smith   PetscFunctionReturn(0);
28847c6ae99SBarry Smith }
28947c6ae99SBarry Smith 
29047c6ae99SBarry Smith /*MC
291aa219208SBarry Smith        DMDASetLocalAdicMFFunctioni - Caches in a DM a local functioni computed by ADIC/ADIFOR
29247c6ae99SBarry Smith 
29347c6ae99SBarry Smith    Synopsis:
294aa219208SBarry Smith    PetscErrorCode  DMDASetLocalAdicFunctioni(DM da,int (ad_lf*)(DMDALocalInfo*,MatStencil*,void*,void*,void*)
29547c6ae99SBarry Smith 
296aa219208SBarry Smith    Logically Collective on DMDA
29747c6ae99SBarry Smith 
29847c6ae99SBarry Smith    Input Parameter:
29947c6ae99SBarry Smith +  da - initial distributed array
30047c6ae99SBarry Smith -  admf_lfi - the local matrix-free function as computed by ADIC/ADIFOR
30147c6ae99SBarry Smith 
30247c6ae99SBarry Smith    Level: intermediate
30347c6ae99SBarry Smith 
30447c6ae99SBarry Smith .keywords:  distributed array, refine
30547c6ae99SBarry Smith 
306aa219208SBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalFunction(),
307aa219208SBarry Smith           DMDASetLocalJacobian(), DMDASetLocalFunctioni()
30847c6ae99SBarry Smith M*/
30947c6ae99SBarry Smith 
31047c6ae99SBarry Smith #undef __FUNCT__
311aa219208SBarry Smith #define __FUNCT__ "DMDASetLocalAdicMFFunctioni_Private"
312aa219208SBarry Smith PetscErrorCode DMDASetLocalAdicMFFunctioni_Private(DM da,PetscErrorCode (*admf_lfi)(DMDALocalInfo*,MatStencil*,void*,void*,void*))
31347c6ae99SBarry Smith {
31447c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
31547c6ae99SBarry Smith   PetscFunctionBegin;
31647c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
31747c6ae99SBarry Smith   dd->adicmf_lfi = admf_lfi;
31847c6ae99SBarry Smith   PetscFunctionReturn(0);
31947c6ae99SBarry Smith }
32047c6ae99SBarry Smith 
32147c6ae99SBarry Smith /*MC
322aa219208SBarry Smith        DMDASetLocalAdicFunctionib - Caches in a DM a block local functioni computed by ADIC/ADIFOR
32347c6ae99SBarry Smith 
32447c6ae99SBarry Smith    Synopsis:
325aa219208SBarry Smith    PetscErrorCode DMDASetLocalAdicFunctionib(DM da,PetscInt (ad_lf*)(DMDALocalInfo*,MatStencil*,void*,void*,void*)
32647c6ae99SBarry Smith 
327aa219208SBarry Smith    Logically Collective on DMDA
32847c6ae99SBarry Smith 
32947c6ae99SBarry Smith    Input Parameter:
33047c6ae99SBarry Smith +  da - initial distributed array
33147c6ae99SBarry Smith -  ad_lfi - the local function as computed by ADIC/ADIFOR
33247c6ae99SBarry Smith 
33347c6ae99SBarry Smith    Level: intermediate
33447c6ae99SBarry Smith 
33547c6ae99SBarry Smith .keywords:  distributed array, refine
33647c6ae99SBarry Smith 
337aa219208SBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalFunction(),
338aa219208SBarry Smith           DMDASetLocalJacobian(), DMDASetLocalFunctionib()
33947c6ae99SBarry Smith M*/
34047c6ae99SBarry Smith 
34147c6ae99SBarry Smith #undef __FUNCT__
342aa219208SBarry Smith #define __FUNCT__ "DMDASetLocalAdicFunctionib_Private"
343aa219208SBarry Smith PetscErrorCode DMDASetLocalAdicFunctionib_Private(DM da,PetscErrorCode (*ad_lfi)(DMDALocalInfo*,MatStencil*,void*,void*,void*))
34447c6ae99SBarry Smith {
34547c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
34647c6ae99SBarry Smith   PetscFunctionBegin;
34747c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
34847c6ae99SBarry Smith   dd->adic_lfib = ad_lfi;
34947c6ae99SBarry Smith   PetscFunctionReturn(0);
35047c6ae99SBarry Smith }
35147c6ae99SBarry Smith 
35247c6ae99SBarry Smith /*MC
353aa219208SBarry Smith        DMDASetLocalAdicMFFunctionib - Caches in a DM a block local functioni computed by ADIC/ADIFOR
35447c6ae99SBarry Smith 
35547c6ae99SBarry Smith    Synopsis:
356aa219208SBarry Smith    PetscErrorCode  DMDASetLocalAdicFunctionib(DM da,int (ad_lf*)(DMDALocalInfo*,MatStencil*,void*,void*,void*)
35747c6ae99SBarry Smith 
358aa219208SBarry Smith    Logically Collective on DMDA
35947c6ae99SBarry Smith 
36047c6ae99SBarry Smith    Input Parameter:
36147c6ae99SBarry Smith +  da - initial distributed array
36247c6ae99SBarry Smith -  admf_lfi - the local matrix-free function as computed by ADIC/ADIFOR
36347c6ae99SBarry Smith 
36447c6ae99SBarry Smith    Level: intermediate
36547c6ae99SBarry Smith 
36647c6ae99SBarry Smith .keywords:  distributed array, refine
36747c6ae99SBarry Smith 
368aa219208SBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalFunction(),
369aa219208SBarry Smith           DMDASetLocalJacobian(), DMDASetLocalFunctionib()
37047c6ae99SBarry Smith M*/
37147c6ae99SBarry Smith 
37247c6ae99SBarry Smith #undef __FUNCT__
373aa219208SBarry Smith #define __FUNCT__ "DMDASetLocalAdicMFFunctionib_Private"
374aa219208SBarry Smith PetscErrorCode DMDASetLocalAdicMFFunctionib_Private(DM da,PetscErrorCode (*admf_lfi)(DMDALocalInfo*,MatStencil*,void*,void*,void*))
37547c6ae99SBarry Smith {
37647c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
37747c6ae99SBarry Smith   PetscFunctionBegin;
37847c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
37947c6ae99SBarry Smith   dd->adicmf_lfib = admf_lfi;
38047c6ae99SBarry Smith   PetscFunctionReturn(0);
38147c6ae99SBarry Smith }
38247c6ae99SBarry Smith 
38347c6ae99SBarry Smith /*MC
384aa219208SBarry Smith        DMDASetLocalAdicMFFunction - Caches in a DM a local function computed by ADIC/ADIFOR
38547c6ae99SBarry Smith 
38647c6ae99SBarry Smith    Synopsis:
387aa219208SBarry Smith    PetscErrorCode DMDASetLocalAdicMFFunction(DM da,DMDALocalFunction1 ad_lf)
38847c6ae99SBarry Smith 
389aa219208SBarry Smith    Logically Collective on DMDA
39047c6ae99SBarry Smith 
39147c6ae99SBarry Smith    Input Parameter:
39247c6ae99SBarry Smith +  da - initial distributed array
39347c6ae99SBarry Smith -  ad_lf - the local function as computed by ADIC/ADIFOR
39447c6ae99SBarry Smith 
39547c6ae99SBarry Smith    Level: intermediate
39647c6ae99SBarry Smith 
39747c6ae99SBarry Smith .keywords:  distributed array, refine
39847c6ae99SBarry Smith 
399aa219208SBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalFunction(),
400aa219208SBarry Smith           DMDASetLocalJacobian()
40147c6ae99SBarry Smith M*/
40247c6ae99SBarry Smith 
40347c6ae99SBarry Smith #undef __FUNCT__
404aa219208SBarry Smith #define __FUNCT__ "DMDASetLocalAdicMFFunction_Private"
405aa219208SBarry Smith PetscErrorCode DMDASetLocalAdicMFFunction_Private(DM da,DMDALocalFunction1 ad_lf)
40647c6ae99SBarry Smith {
40747c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
40847c6ae99SBarry Smith   PetscFunctionBegin;
40947c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
41047c6ae99SBarry Smith   dd->adicmf_lf = ad_lf;
41147c6ae99SBarry Smith   PetscFunctionReturn(0);
41247c6ae99SBarry Smith }
41347c6ae99SBarry Smith 
41447c6ae99SBarry Smith /*@C
415aa219208SBarry Smith        DMDASetLocalJacobian - Caches in a DM a local Jacobian computation function
41647c6ae99SBarry Smith 
417aa219208SBarry Smith    Logically Collective on DMDA
41847c6ae99SBarry Smith 
41947c6ae99SBarry Smith 
42047c6ae99SBarry Smith    Input Parameter:
42147c6ae99SBarry Smith +  da - initial distributed array
42247c6ae99SBarry Smith -  lj - the local Jacobian
42347c6ae99SBarry Smith 
42447c6ae99SBarry Smith    Level: intermediate
42547c6ae99SBarry Smith 
42647c6ae99SBarry Smith    Notes: The routine SNESDAFormFunction() uses this the cached function to evaluate the user provided function.
42747c6ae99SBarry Smith 
42847c6ae99SBarry Smith .keywords:  distributed array, refine
42947c6ae99SBarry Smith 
430aa219208SBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalFunction()
43147c6ae99SBarry Smith @*/
43247c6ae99SBarry Smith #undef __FUNCT__
433aa219208SBarry Smith #define __FUNCT__ "DMDASetLocalJacobian"
4347087cfbeSBarry Smith PetscErrorCode  DMDASetLocalJacobian(DM da,DMDALocalFunction1 lj)
43547c6ae99SBarry Smith {
43647c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
43747c6ae99SBarry Smith   PetscFunctionBegin;
43847c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
43947c6ae99SBarry Smith   dd->lj    = lj;
44047c6ae99SBarry Smith   PetscFunctionReturn(0);
44147c6ae99SBarry Smith }
44247c6ae99SBarry Smith 
44347c6ae99SBarry Smith #undef __FUNCT__
444aa219208SBarry Smith #define __FUNCT__ "DMDAGetLocalFunction"
44547c6ae99SBarry Smith /*@C
446aa219208SBarry Smith        DMDAGetLocalFunction - Gets from a DM a local function and its ADIC/ADIFOR Jacobian
44747c6ae99SBarry Smith 
44847c6ae99SBarry Smith    Note Collective
44947c6ae99SBarry Smith 
45047c6ae99SBarry Smith    Input Parameter:
45147c6ae99SBarry Smith .  da - initial distributed array
45247c6ae99SBarry Smith 
45347c6ae99SBarry Smith    Output Parameter:
45447c6ae99SBarry Smith .  lf - the local function
45547c6ae99SBarry Smith 
45647c6ae99SBarry Smith    Level: intermediate
45747c6ae99SBarry Smith 
45847c6ae99SBarry Smith .keywords:  distributed array, refine
45947c6ae99SBarry Smith 
460aa219208SBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalJacobian(), DMDASetLocalFunction()
46147c6ae99SBarry Smith @*/
4627087cfbeSBarry Smith PetscErrorCode  DMDAGetLocalFunction(DM da,DMDALocalFunction1 *lf)
46347c6ae99SBarry Smith {
46447c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
46547c6ae99SBarry Smith   PetscFunctionBegin;
46647c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
46747c6ae99SBarry Smith   if (lf) *lf = dd->lf;
46847c6ae99SBarry Smith   PetscFunctionReturn(0);
46947c6ae99SBarry Smith }
47047c6ae99SBarry Smith 
47147c6ae99SBarry Smith #undef __FUNCT__
472aa219208SBarry Smith #define __FUNCT__ "DMDAGetLocalJacobian"
47347c6ae99SBarry Smith /*@C
474aa219208SBarry Smith        DMDAGetLocalJacobian - Gets from a DM a local jacobian
47547c6ae99SBarry Smith 
47647c6ae99SBarry Smith    Not Collective
47747c6ae99SBarry Smith 
47847c6ae99SBarry Smith    Input Parameter:
47947c6ae99SBarry Smith .  da - initial distributed array
48047c6ae99SBarry Smith 
48147c6ae99SBarry Smith    Output Parameter:
48247c6ae99SBarry Smith .  lj - the local jacobian
48347c6ae99SBarry Smith 
48447c6ae99SBarry Smith    Level: intermediate
48547c6ae99SBarry Smith 
48647c6ae99SBarry Smith .keywords:  distributed array, refine
48747c6ae99SBarry Smith 
488aa219208SBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalJacobian()
48947c6ae99SBarry Smith @*/
4907087cfbeSBarry Smith PetscErrorCode  DMDAGetLocalJacobian(DM da,DMDALocalFunction1 *lj)
49147c6ae99SBarry Smith {
49247c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
49347c6ae99SBarry Smith   PetscFunctionBegin;
49447c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
49547c6ae99SBarry Smith   if (lj) *lj = dd->lj;
49647c6ae99SBarry Smith   PetscFunctionReturn(0);
49747c6ae99SBarry Smith }
49847c6ae99SBarry Smith 
49947c6ae99SBarry Smith #undef __FUNCT__
500aa219208SBarry Smith #define __FUNCT__ "DMDAFormFunction"
50147c6ae99SBarry Smith /*@
502aa219208SBarry Smith     DMDAFormFunction - Evaluates a user provided function on each processor that
503aa219208SBarry Smith         share a DMDA
50447c6ae99SBarry Smith 
50547c6ae99SBarry Smith    Input Parameters:
5069a42bb27SBarry Smith +    da - the DM that defines the grid
50747c6ae99SBarry Smith .    vu - input vector
50847c6ae99SBarry Smith .    vfu - output vector
50947c6ae99SBarry Smith -    w - any user data
51047c6ae99SBarry Smith 
51147c6ae99SBarry Smith     Notes: Does NOT do ghost updates on vu upon entry
51247c6ae99SBarry Smith 
513aa219208SBarry Smith            This should eventually replace DMDAFormFunction1
51447c6ae99SBarry Smith 
51547c6ae99SBarry Smith     Level: advanced
51647c6ae99SBarry Smith 
517aa219208SBarry Smith .seealso: DMDAComputeJacobian1WithAdic()
51847c6ae99SBarry Smith 
51947c6ae99SBarry Smith @*/
5207087cfbeSBarry Smith PetscErrorCode  DMDAFormFunction(DM da,PetscErrorCode (*lf)(void),Vec vu,Vec vfu,void *w)
52147c6ae99SBarry Smith {
52247c6ae99SBarry Smith   PetscErrorCode ierr;
52347c6ae99SBarry Smith   void           *u,*fu;
524aa219208SBarry Smith   DMDALocalInfo  info;
525aa219208SBarry Smith   PetscErrorCode (*f)(DMDALocalInfo*,void*,void*,void*) = (PetscErrorCode (*)(DMDALocalInfo*,void*,void*,void*))lf;
52647c6ae99SBarry Smith 
52747c6ae99SBarry Smith   PetscFunctionBegin;
528aa219208SBarry Smith   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
529aa219208SBarry Smith   ierr = DMDAVecGetArray(da,vu,&u);CHKERRQ(ierr);
530aa219208SBarry Smith   ierr = DMDAVecGetArray(da,vfu,&fu);CHKERRQ(ierr);
53147c6ae99SBarry Smith 
53239d508bbSBarry Smith   ierr = (*f)(&info,u,fu,w);CHKERRQ(ierr);
53347c6ae99SBarry Smith 
534aa219208SBarry Smith   ierr = DMDAVecRestoreArray(da,vu,&u);CHKERRQ(ierr);
535aa219208SBarry Smith   ierr = DMDAVecRestoreArray(da,vfu,&fu);CHKERRQ(ierr);
53647c6ae99SBarry Smith   PetscFunctionReturn(0);
53747c6ae99SBarry Smith }
53847c6ae99SBarry Smith 
53947c6ae99SBarry Smith #undef __FUNCT__
540aa219208SBarry Smith #define __FUNCT__ "DMDAFormFunctionLocal"
54147c6ae99SBarry Smith /*@C
542aa219208SBarry Smith    DMDAFormFunctionLocal - This is a universal function evaluation routine for
5439a42bb27SBarry Smith    a local DM function.
54447c6ae99SBarry Smith 
545aa219208SBarry Smith    Collective on DMDA
54647c6ae99SBarry Smith 
54747c6ae99SBarry Smith    Input Parameters:
5489a42bb27SBarry Smith +  da - the DM context
54947c6ae99SBarry Smith .  func - The local function
55047c6ae99SBarry Smith .  X - input vector
55147c6ae99SBarry Smith .  F - function vector
55247c6ae99SBarry Smith -  ctx - A user context
55347c6ae99SBarry Smith 
55447c6ae99SBarry Smith    Level: intermediate
55547c6ae99SBarry Smith 
556aa219208SBarry Smith .seealso: DMDASetLocalFunction(), DMDASetLocalJacobian(), DMDASetLocalAdicFunction(), DMDASetLocalAdicMFFunction(),
55747c6ae99SBarry Smith           SNESSetFunction(), SNESSetJacobian()
55847c6ae99SBarry Smith 
55947c6ae99SBarry Smith @*/
5607087cfbeSBarry Smith PetscErrorCode  DMDAFormFunctionLocal(DM da, DMDALocalFunction1 func, Vec X, Vec F, void *ctx)
56147c6ae99SBarry Smith {
56247c6ae99SBarry Smith   Vec            localX;
563aa219208SBarry Smith   DMDALocalInfo  info;
56447c6ae99SBarry Smith   void           *u;
56547c6ae99SBarry Smith   void           *fu;
56647c6ae99SBarry Smith   PetscErrorCode ierr;
56747c6ae99SBarry Smith 
56847c6ae99SBarry Smith   PetscFunctionBegin;
5699a42bb27SBarry Smith   ierr = DMGetLocalVector(da,&localX);CHKERRQ(ierr);
57047c6ae99SBarry Smith   /*
57147c6ae99SBarry Smith      Scatter ghost points to local vector, using the 2-step process
5729a42bb27SBarry Smith         DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
57347c6ae99SBarry Smith   */
5749a42bb27SBarry Smith   ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
5759a42bb27SBarry Smith   ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
576aa219208SBarry Smith   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
577aa219208SBarry Smith   ierr = DMDAVecGetArray(da,localX,&u);CHKERRQ(ierr);
578aa219208SBarry Smith   ierr = DMDAVecGetArray(da,F,&fu);CHKERRQ(ierr);
57939d508bbSBarry Smith   ierr = (*func)(&info,u,fu,ctx);CHKERRQ(ierr);
580aa219208SBarry Smith   ierr = DMDAVecRestoreArray(da,localX,&u);CHKERRQ(ierr);
581aa219208SBarry Smith   ierr = DMDAVecRestoreArray(da,F,&fu);CHKERRQ(ierr);
5829a42bb27SBarry Smith   ierr = DMRestoreLocalVector(da,&localX);CHKERRQ(ierr);
58347c6ae99SBarry Smith   PetscFunctionReturn(0);
58447c6ae99SBarry Smith }
58547c6ae99SBarry Smith 
58647c6ae99SBarry Smith #undef __FUNCT__
587aa219208SBarry Smith #define __FUNCT__ "DMDAFormFunctionLocalGhost"
58847c6ae99SBarry Smith /*@C
589aa219208SBarry Smith    DMDAFormFunctionLocalGhost - This is a universal function evaluation routine for
5909a42bb27SBarry Smith    a local DM function, but the ghost values of the output are communicated and added.
59147c6ae99SBarry Smith 
592aa219208SBarry Smith    Collective on DMDA
59347c6ae99SBarry Smith 
59447c6ae99SBarry Smith    Input Parameters:
5959a42bb27SBarry Smith +  da - the DM context
59647c6ae99SBarry Smith .  func - The local function
59747c6ae99SBarry Smith .  X - input vector
59847c6ae99SBarry Smith .  F - function vector
59947c6ae99SBarry Smith -  ctx - A user context
60047c6ae99SBarry Smith 
60147c6ae99SBarry Smith    Level: intermediate
60247c6ae99SBarry Smith 
603aa219208SBarry Smith .seealso: DMDASetLocalFunction(), DMDASetLocalJacobian(), DMDASetLocalAdicFunction(), DMDASetLocalAdicMFFunction(),
60447c6ae99SBarry Smith           SNESSetFunction(), SNESSetJacobian()
60547c6ae99SBarry Smith 
60647c6ae99SBarry Smith @*/
6077087cfbeSBarry Smith PetscErrorCode  DMDAFormFunctionLocalGhost(DM da, DMDALocalFunction1 func, Vec X, Vec F, void *ctx)
60847c6ae99SBarry Smith {
60947c6ae99SBarry Smith   Vec            localX, localF;
610aa219208SBarry Smith   DMDALocalInfo  info;
61147c6ae99SBarry Smith   void           *u;
61247c6ae99SBarry Smith   void           *fu;
61347c6ae99SBarry Smith   PetscErrorCode ierr;
61447c6ae99SBarry Smith 
61547c6ae99SBarry Smith   PetscFunctionBegin;
6169a42bb27SBarry Smith   ierr = DMGetLocalVector(da,&localX);CHKERRQ(ierr);
6179a42bb27SBarry Smith   ierr = DMGetLocalVector(da,&localF);CHKERRQ(ierr);
61847c6ae99SBarry Smith   /*
61947c6ae99SBarry Smith      Scatter ghost points to local vector, using the 2-step process
6209a42bb27SBarry Smith         DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
62147c6ae99SBarry Smith   */
6229a42bb27SBarry Smith   ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
6239a42bb27SBarry Smith   ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
62447c6ae99SBarry Smith   ierr = VecSet(F, 0.0);CHKERRQ(ierr);
62547c6ae99SBarry Smith   ierr = VecSet(localF, 0.0);CHKERRQ(ierr);
626aa219208SBarry Smith   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
627aa219208SBarry Smith   ierr = DMDAVecGetArray(da,localX,&u);CHKERRQ(ierr);
628aa219208SBarry Smith   ierr = DMDAVecGetArray(da,localF,&fu);CHKERRQ(ierr);
62939d508bbSBarry Smith   ierr = (*func)(&info,u,fu,ctx);CHKERRQ(ierr);
6309a42bb27SBarry Smith   ierr = DMLocalToGlobalBegin(da,localF,ADD_VALUES,F);CHKERRQ(ierr);
6319a42bb27SBarry Smith   ierr = DMLocalToGlobalEnd(da,localF,ADD_VALUES,F);CHKERRQ(ierr);
632aa219208SBarry Smith   ierr = DMDAVecRestoreArray(da,localX,&u);CHKERRQ(ierr);
633aa219208SBarry Smith   ierr = DMDAVecRestoreArray(da,localF,&fu);CHKERRQ(ierr);
6349a42bb27SBarry Smith   ierr = DMRestoreLocalVector(da,&localX);CHKERRQ(ierr);
6359a42bb27SBarry Smith   ierr = DMRestoreLocalVector(da,&localF);CHKERRQ(ierr);
63647c6ae99SBarry Smith   PetscFunctionReturn(0);
63747c6ae99SBarry Smith }
63847c6ae99SBarry Smith 
63947c6ae99SBarry Smith #undef __FUNCT__
640aa219208SBarry Smith #define __FUNCT__ "DMDAFormFunction1"
64147c6ae99SBarry Smith /*@
642aa219208SBarry Smith     DMDAFormFunction1 - Evaluates a user provided function on each processor that
643aa219208SBarry Smith         share a DMDA
64447c6ae99SBarry Smith 
64547c6ae99SBarry Smith    Input Parameters:
6469a42bb27SBarry Smith +    da - the DM that defines the grid
64747c6ae99SBarry Smith .    vu - input vector
64847c6ae99SBarry Smith .    vfu - output vector
64947c6ae99SBarry Smith -    w - any user data
65047c6ae99SBarry Smith 
65147c6ae99SBarry Smith     Notes: Does NOT do ghost updates on vu upon entry
65247c6ae99SBarry Smith 
65347c6ae99SBarry Smith     Level: advanced
65447c6ae99SBarry Smith 
655aa219208SBarry Smith .seealso: DMDAComputeJacobian1WithAdic()
65647c6ae99SBarry Smith 
65747c6ae99SBarry Smith @*/
6587087cfbeSBarry Smith PetscErrorCode  DMDAFormFunction1(DM da,Vec vu,Vec vfu,void *w)
65947c6ae99SBarry Smith {
66047c6ae99SBarry Smith   PetscErrorCode ierr;
66147c6ae99SBarry Smith   void           *u,*fu;
662aa219208SBarry Smith   DMDALocalInfo  info;
66347c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
66447c6ae99SBarry Smith 
66547c6ae99SBarry Smith   PetscFunctionBegin;
666aa219208SBarry Smith   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
667aa219208SBarry Smith   ierr = DMDAVecGetArray(da,vu,&u);CHKERRQ(ierr);
668aa219208SBarry Smith   ierr = DMDAVecGetArray(da,vfu,&fu);CHKERRQ(ierr);
66947c6ae99SBarry Smith 
67047c6ae99SBarry Smith   CHKMEMQ;
67139d508bbSBarry Smith   ierr = (*dd->lf)(&info,u,fu,w);CHKERRQ(ierr);
67247c6ae99SBarry Smith   CHKMEMQ;
67347c6ae99SBarry Smith 
674aa219208SBarry Smith   ierr = DMDAVecRestoreArray(da,vu,&u);CHKERRQ(ierr);
675aa219208SBarry Smith   ierr = DMDAVecRestoreArray(da,vfu,&fu);CHKERRQ(ierr);
67647c6ae99SBarry Smith   PetscFunctionReturn(0);
67747c6ae99SBarry Smith }
67847c6ae99SBarry Smith 
67947c6ae99SBarry Smith #undef __FUNCT__
680aa219208SBarry Smith #define __FUNCT__ "DMDAFormFunctioniTest1"
6817087cfbeSBarry Smith PetscErrorCode  DMDAFormFunctioniTest1(DM da,void *w)
68247c6ae99SBarry Smith {
68347c6ae99SBarry Smith   Vec            vu,fu,fui;
68447c6ae99SBarry Smith   PetscErrorCode ierr;
68547c6ae99SBarry Smith   PetscInt       i,n;
68647c6ae99SBarry Smith   PetscScalar    *ui;
68747c6ae99SBarry Smith   PetscRandom    rnd;
68847c6ae99SBarry Smith   PetscReal      norm;
68947c6ae99SBarry Smith 
69047c6ae99SBarry Smith   PetscFunctionBegin;
6919a42bb27SBarry Smith   ierr = DMGetLocalVector(da,&vu);CHKERRQ(ierr);
69247c6ae99SBarry Smith   ierr = PetscRandomCreate(PETSC_COMM_SELF,&rnd);CHKERRQ(ierr);
69347c6ae99SBarry Smith   ierr = PetscRandomSetFromOptions(rnd);CHKERRQ(ierr);
69447c6ae99SBarry Smith   ierr = VecSetRandom(vu,rnd);CHKERRQ(ierr);
695*fcfd50ebSBarry Smith   ierr = PetscRandomDestroy(&rnd);CHKERRQ(ierr);
69647c6ae99SBarry Smith 
6979a42bb27SBarry Smith   ierr = DMGetGlobalVector(da,&fu);CHKERRQ(ierr);
6989a42bb27SBarry Smith   ierr = DMGetGlobalVector(da,&fui);CHKERRQ(ierr);
69947c6ae99SBarry Smith 
700aa219208SBarry Smith   ierr = DMDAFormFunction1(da,vu,fu,w);CHKERRQ(ierr);
70147c6ae99SBarry Smith 
70247c6ae99SBarry Smith   ierr = VecGetArray(fui,&ui);CHKERRQ(ierr);
70347c6ae99SBarry Smith   ierr = VecGetLocalSize(fui,&n);CHKERRQ(ierr);
70447c6ae99SBarry Smith   for (i=0; i<n; i++) {
705aa219208SBarry Smith     ierr = DMDAFormFunctioni1(da,i,vu,ui+i,w);CHKERRQ(ierr);
70647c6ae99SBarry Smith   }
70747c6ae99SBarry Smith   ierr = VecRestoreArray(fui,&ui);CHKERRQ(ierr);
70847c6ae99SBarry Smith 
70947c6ae99SBarry Smith   ierr = VecAXPY(fui,-1.0,fu);CHKERRQ(ierr);
71047c6ae99SBarry Smith   ierr = VecNorm(fui,NORM_2,&norm);CHKERRQ(ierr);
71147c6ae99SBarry Smith   ierr = PetscPrintf(((PetscObject)da)->comm,"Norm of difference in vectors %G\n",norm);CHKERRQ(ierr);
71247c6ae99SBarry Smith   ierr = VecView(fu,0);CHKERRQ(ierr);
71347c6ae99SBarry Smith   ierr = VecView(fui,0);CHKERRQ(ierr);
71447c6ae99SBarry Smith 
7159a42bb27SBarry Smith   ierr = DMRestoreLocalVector(da,&vu);CHKERRQ(ierr);
7169a42bb27SBarry Smith   ierr = DMRestoreGlobalVector(da,&fu);CHKERRQ(ierr);
7179a42bb27SBarry Smith   ierr = DMRestoreGlobalVector(da,&fui);CHKERRQ(ierr);
71847c6ae99SBarry Smith   PetscFunctionReturn(0);
71947c6ae99SBarry Smith }
72047c6ae99SBarry Smith 
72147c6ae99SBarry Smith #undef __FUNCT__
722aa219208SBarry Smith #define __FUNCT__ "DMDAFormFunctioni1"
72347c6ae99SBarry Smith /*@
724aa219208SBarry Smith     DMDAFormFunctioni1 - Evaluates a user provided point-wise function
72547c6ae99SBarry Smith 
72647c6ae99SBarry Smith    Input Parameters:
7279a42bb27SBarry Smith +    da - the DM that defines the grid
72847c6ae99SBarry Smith .    i - the component of the function we wish to compute (must be local)
72947c6ae99SBarry Smith .    vu - input vector
73047c6ae99SBarry Smith .    vfu - output value
73147c6ae99SBarry Smith -    w - any user data
73247c6ae99SBarry Smith 
73347c6ae99SBarry Smith     Notes: Does NOT do ghost updates on vu upon entry
73447c6ae99SBarry Smith 
73547c6ae99SBarry Smith     Level: advanced
73647c6ae99SBarry Smith 
737aa219208SBarry Smith .seealso: DMDAComputeJacobian1WithAdic()
73847c6ae99SBarry Smith 
73947c6ae99SBarry Smith @*/
7407087cfbeSBarry Smith PetscErrorCode  DMDAFormFunctioni1(DM da,PetscInt i,Vec vu,PetscScalar *vfu,void *w)
74147c6ae99SBarry Smith {
74247c6ae99SBarry Smith   PetscErrorCode ierr;
74347c6ae99SBarry Smith   void           *u;
744aa219208SBarry Smith   DMDALocalInfo  info;
74547c6ae99SBarry Smith   MatStencil     stencil;
74647c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
74747c6ae99SBarry Smith 
74847c6ae99SBarry Smith   PetscFunctionBegin;
74947c6ae99SBarry Smith 
750aa219208SBarry Smith   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
751aa219208SBarry Smith   ierr = DMDAVecGetArray(da,vu,&u);CHKERRQ(ierr);
75247c6ae99SBarry Smith 
75347c6ae99SBarry Smith   /* figure out stencil value from i */
75447c6ae99SBarry Smith   stencil.c = i % info.dof;
75547c6ae99SBarry Smith   stencil.i = (i % (info.xm*info.dof))/info.dof;
75647c6ae99SBarry Smith   stencil.j = (i % (info.xm*info.ym*info.dof))/(info.xm*info.dof);
75747c6ae99SBarry Smith   stencil.k = i/(info.xm*info.ym*info.dof);
75847c6ae99SBarry Smith 
75947c6ae99SBarry Smith   ierr = (*dd->lfi)(&info,&stencil,u,vfu,w);CHKERRQ(ierr);
76047c6ae99SBarry Smith 
761aa219208SBarry Smith   ierr = DMDAVecRestoreArray(da,vu,&u);CHKERRQ(ierr);
76247c6ae99SBarry Smith   PetscFunctionReturn(0);
76347c6ae99SBarry Smith }
76447c6ae99SBarry Smith 
76547c6ae99SBarry Smith #undef __FUNCT__
766aa219208SBarry Smith #define __FUNCT__ "DMDAFormFunctionib1"
76747c6ae99SBarry Smith /*@
768aa219208SBarry Smith     DMDAFormFunctionib1 - Evaluates a user provided point-block function
76947c6ae99SBarry Smith 
77047c6ae99SBarry Smith    Input Parameters:
7719a42bb27SBarry Smith +    da - the DM that defines the grid
77247c6ae99SBarry Smith .    i - the component of the function we wish to compute (must be local)
77347c6ae99SBarry Smith .    vu - input vector
77447c6ae99SBarry Smith .    vfu - output value
77547c6ae99SBarry Smith -    w - any user data
77647c6ae99SBarry Smith 
77747c6ae99SBarry Smith     Notes: Does NOT do ghost updates on vu upon entry
77847c6ae99SBarry Smith 
77947c6ae99SBarry Smith     Level: advanced
78047c6ae99SBarry Smith 
781aa219208SBarry Smith .seealso: DMDAComputeJacobian1WithAdic()
78247c6ae99SBarry Smith 
78347c6ae99SBarry Smith @*/
7847087cfbeSBarry Smith PetscErrorCode  DMDAFormFunctionib1(DM da,PetscInt i,Vec vu,PetscScalar *vfu,void *w)
78547c6ae99SBarry Smith {
78647c6ae99SBarry Smith   PetscErrorCode ierr;
78747c6ae99SBarry Smith   void           *u;
788aa219208SBarry Smith   DMDALocalInfo  info;
78947c6ae99SBarry Smith   MatStencil     stencil;
79047c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
79147c6ae99SBarry Smith 
79247c6ae99SBarry Smith   PetscFunctionBegin;
793aa219208SBarry Smith   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
794aa219208SBarry Smith   ierr = DMDAVecGetArray(da,vu,&u);CHKERRQ(ierr);
79547c6ae99SBarry Smith 
79647c6ae99SBarry Smith   /* figure out stencil value from i */
79747c6ae99SBarry Smith   stencil.c = i % info.dof;
79847c6ae99SBarry Smith   if (stencil.c) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONG,"Point-block functions can only be called for the entire block");
79947c6ae99SBarry Smith   stencil.i = (i % (info.xm*info.dof))/info.dof;
80047c6ae99SBarry Smith   stencil.j = (i % (info.xm*info.ym*info.dof))/(info.xm*info.dof);
80147c6ae99SBarry Smith   stencil.k = i/(info.xm*info.ym*info.dof);
80247c6ae99SBarry Smith 
80347c6ae99SBarry Smith   ierr = (*dd->lfib)(&info,&stencil,u,vfu,w);CHKERRQ(ierr);
80447c6ae99SBarry Smith 
805aa219208SBarry Smith   ierr = DMDAVecRestoreArray(da,vu,&u);CHKERRQ(ierr);
80647c6ae99SBarry Smith   PetscFunctionReturn(0);
80747c6ae99SBarry Smith }
80847c6ae99SBarry Smith 
80947c6ae99SBarry Smith #if defined(new)
81047c6ae99SBarry Smith #undef __FUNCT__
811aa219208SBarry Smith #define __FUNCT__ "DMDAGetDiagonal_MFFD"
81247c6ae99SBarry Smith /*
813aa219208SBarry Smith   DMDAGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix where local
814aa219208SBarry Smith     function lives on a DMDA
81547c6ae99SBarry Smith 
81647c6ae99SBarry Smith         y ~= (F(u + ha) - F(u))/h,
81747c6ae99SBarry Smith   where F = nonlinear function, as set by SNESSetFunction()
81847c6ae99SBarry Smith         u = current iterate
81947c6ae99SBarry Smith         h = difference interval
82047c6ae99SBarry Smith */
821aa219208SBarry Smith PetscErrorCode DMDAGetDiagonal_MFFD(DM da,Vec U,Vec a)
82247c6ae99SBarry Smith {
82347c6ae99SBarry Smith   PetscScalar    h,*aa,*ww,v;
82447c6ae99SBarry Smith   PetscReal      epsilon = PETSC_SQRT_MACHINE_EPSILON,umin = 100.0*PETSC_SQRT_MACHINE_EPSILON;
82547c6ae99SBarry Smith   PetscErrorCode ierr;
82647c6ae99SBarry Smith   PetscInt       gI,nI;
82747c6ae99SBarry Smith   MatStencil     stencil;
828aa219208SBarry Smith   DMDALocalInfo  info;
82947c6ae99SBarry Smith 
83047c6ae99SBarry Smith   PetscFunctionBegin;
83147c6ae99SBarry Smith   ierr = (*ctx->func)(0,U,a,ctx->funcctx);CHKERRQ(ierr);
83247c6ae99SBarry Smith   ierr = (*ctx->funcisetbase)(U,ctx->funcctx);CHKERRQ(ierr);
83347c6ae99SBarry Smith 
83447c6ae99SBarry Smith   ierr = VecGetArray(U,&ww);CHKERRQ(ierr);
83547c6ae99SBarry Smith   ierr = VecGetArray(a,&aa);CHKERRQ(ierr);
83647c6ae99SBarry Smith 
83747c6ae99SBarry Smith   nI = 0;
83847c6ae99SBarry Smith     h  = ww[gI];
83947c6ae99SBarry Smith     if (h == 0.0) h = 1.0;
84047c6ae99SBarry Smith #if !defined(PETSC_USE_COMPLEX)
84147c6ae99SBarry Smith     if (h < umin && h >= 0.0)      h = umin;
84247c6ae99SBarry Smith     else if (h < 0.0 && h > -umin) h = -umin;
84347c6ae99SBarry Smith #else
84447c6ae99SBarry Smith     if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0)     h = umin;
84547c6ae99SBarry Smith     else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin;
84647c6ae99SBarry Smith #endif
84747c6ae99SBarry Smith     h     *= epsilon;
84847c6ae99SBarry Smith 
84947c6ae99SBarry Smith     ww[gI] += h;
85047c6ae99SBarry Smith     ierr          = (*ctx->funci)(i,w,&v,ctx->funcctx);CHKERRQ(ierr);
85147c6ae99SBarry Smith     aa[nI]  = (v - aa[nI])/h;
85247c6ae99SBarry Smith     ww[gI] -= h;
85347c6ae99SBarry Smith     nI++;
85447c6ae99SBarry Smith   }
85547c6ae99SBarry Smith   ierr = VecRestoreArray(U,&ww);CHKERRQ(ierr);
85647c6ae99SBarry Smith   ierr = VecRestoreArray(a,&aa);CHKERRQ(ierr);
85747c6ae99SBarry Smith   PetscFunctionReturn(0);
85847c6ae99SBarry Smith }
85947c6ae99SBarry Smith #endif
86047c6ae99SBarry Smith 
86147c6ae99SBarry Smith #if defined(PETSC_HAVE_ADIC)
86247c6ae99SBarry Smith EXTERN_C_BEGIN
863c6db04a5SJed Brown #include <adic/ad_utils.h>
86447c6ae99SBarry Smith EXTERN_C_END
86547c6ae99SBarry Smith 
86647c6ae99SBarry Smith #undef __FUNCT__
867aa219208SBarry Smith #define __FUNCT__ "DMDAComputeJacobian1WithAdic"
86847c6ae99SBarry Smith /*@C
869aa219208SBarry Smith     DMDAComputeJacobian1WithAdic - Evaluates a adiC provided Jacobian function on each processor that
870aa219208SBarry Smith         share a DMDA
87147c6ae99SBarry Smith 
87247c6ae99SBarry Smith    Input Parameters:
8739a42bb27SBarry Smith +    da - the DM that defines the grid
87447c6ae99SBarry Smith .    vu - input vector (ghosted)
87547c6ae99SBarry Smith .    J - output matrix
87647c6ae99SBarry Smith -    w - any user data
87747c6ae99SBarry Smith 
87847c6ae99SBarry Smith    Level: advanced
87947c6ae99SBarry Smith 
88047c6ae99SBarry Smith     Notes: Does NOT do ghost updates on vu upon entry
88147c6ae99SBarry Smith 
882aa219208SBarry Smith .seealso: DMDAFormFunction1()
88347c6ae99SBarry Smith 
88447c6ae99SBarry Smith @*/
8857087cfbeSBarry Smith PetscErrorCode  DMDAComputeJacobian1WithAdic(DM da,Vec vu,Mat J,void *w)
88647c6ae99SBarry Smith {
88747c6ae99SBarry Smith   PetscErrorCode ierr;
88847c6ae99SBarry Smith   PetscInt       gtdof,tdof;
88947c6ae99SBarry Smith   PetscScalar    *ustart;
890aa219208SBarry Smith   DMDALocalInfo  info;
89147c6ae99SBarry Smith   void           *ad_u,*ad_f,*ad_ustart,*ad_fstart;
89247c6ae99SBarry Smith   ISColoring     iscoloring;
89347c6ae99SBarry Smith 
89447c6ae99SBarry Smith   PetscFunctionBegin;
895aa219208SBarry Smith   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
89647c6ae99SBarry Smith 
89747c6ae99SBarry Smith   PetscADResetIndep();
89847c6ae99SBarry Smith 
89947c6ae99SBarry Smith   /* get space for derivative objects.  */
900aa219208SBarry Smith   ierr = DMDAGetAdicArray(da,PETSC_TRUE,&ad_u,&ad_ustart,&gtdof);CHKERRQ(ierr);
901aa219208SBarry Smith   ierr = DMDAGetAdicArray(da,PETSC_FALSE,&ad_f,&ad_fstart,&tdof);CHKERRQ(ierr);
90247c6ae99SBarry Smith   ierr = VecGetArray(vu,&ustart);CHKERRQ(ierr);
90394013140SBarry Smith   ierr = DMGetColoring(da,IS_COLORING_GHOSTED,MATAIJ,&iscoloring);CHKERRQ(ierr);
90447c6ae99SBarry Smith 
90547c6ae99SBarry Smith   PetscADSetValueAndColor(ad_ustart,gtdof,iscoloring->colors,ustart);
90647c6ae99SBarry Smith 
90747c6ae99SBarry Smith   ierr = VecRestoreArray(vu,&ustart);CHKERRQ(ierr);
908*fcfd50ebSBarry Smith   ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr);
90947c6ae99SBarry Smith   ierr = PetscADIncrementTotalGradSize(iscoloring->n);CHKERRQ(ierr);
91047c6ae99SBarry Smith   PetscADSetIndepDone();
91147c6ae99SBarry Smith 
912aa219208SBarry Smith   ierr = PetscLogEventBegin(DMDA_LocalADFunction,0,0,0,0);CHKERRQ(ierr);
91347c6ae99SBarry Smith   ierr = (*dd->adic_lf)(&info,ad_u,ad_f,w);CHKERRQ(ierr);
914aa219208SBarry Smith   ierr = PetscLogEventEnd(DMDA_LocalADFunction,0,0,0,0);CHKERRQ(ierr);
91547c6ae99SBarry Smith 
91647c6ae99SBarry Smith   /* stick the values into the matrix */
91747c6ae99SBarry Smith   ierr = MatSetValuesAdic(J,(PetscScalar**)ad_fstart);CHKERRQ(ierr);
91847c6ae99SBarry Smith 
91947c6ae99SBarry Smith   /* return space for derivative objects.  */
920aa219208SBarry Smith   ierr = DMDARestoreAdicArray(da,PETSC_TRUE,&ad_u,&ad_ustart,&gtdof);CHKERRQ(ierr);
921aa219208SBarry Smith   ierr = DMDARestoreAdicArray(da,PETSC_FALSE,&ad_f,&ad_fstart,&tdof);CHKERRQ(ierr);
92247c6ae99SBarry Smith   PetscFunctionReturn(0);
92347c6ae99SBarry Smith }
92447c6ae99SBarry Smith 
92547c6ae99SBarry Smith #undef __FUNCT__
926aa219208SBarry Smith #define __FUNCT__ "DMDAMultiplyByJacobian1WithAdic"
92747c6ae99SBarry Smith /*@C
928aa219208SBarry Smith     DMDAMultiplyByJacobian1WithAdic - Applies an ADIC-provided Jacobian function to a vector on
929aa219208SBarry Smith     each processor that shares a DMDA.
93047c6ae99SBarry Smith 
93147c6ae99SBarry Smith     Input Parameters:
9329a42bb27SBarry Smith +   da - the DM that defines the grid
93347c6ae99SBarry Smith .   vu - Jacobian is computed at this point (ghosted)
93447c6ae99SBarry Smith .   v - product is done on this vector (ghosted)
93547c6ae99SBarry Smith .   fu - output vector = J(vu)*v (not ghosted)
93647c6ae99SBarry Smith -   w - any user data
93747c6ae99SBarry Smith 
93847c6ae99SBarry Smith     Notes:
93947c6ae99SBarry Smith     This routine does NOT do ghost updates on vu upon entry.
94047c6ae99SBarry Smith 
94147c6ae99SBarry Smith    Level: advanced
94247c6ae99SBarry Smith 
943aa219208SBarry Smith .seealso: DMDAFormFunction1()
94447c6ae99SBarry Smith 
94547c6ae99SBarry Smith @*/
9467087cfbeSBarry Smith PetscErrorCode  DMDAMultiplyByJacobian1WithAdic(DM da,Vec vu,Vec v,Vec f,void *w)
94747c6ae99SBarry Smith {
94847c6ae99SBarry Smith   PetscErrorCode ierr;
94947c6ae99SBarry Smith   PetscInt       i,gtdof,tdof;
95047c6ae99SBarry Smith   PetscScalar    *avu,*av,*af,*ad_vustart,*ad_fstart;
951aa219208SBarry Smith   DMDALocalInfo  info;
95247c6ae99SBarry Smith   void           *ad_vu,*ad_f;
95347c6ae99SBarry Smith 
95447c6ae99SBarry Smith   PetscFunctionBegin;
955aa219208SBarry Smith   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
95647c6ae99SBarry Smith 
95747c6ae99SBarry Smith   /* get space for derivative objects.  */
958aa219208SBarry Smith   ierr = DMDAGetAdicMFArray(da,PETSC_TRUE,&ad_vu,&ad_vustart,&gtdof);CHKERRQ(ierr);
959aa219208SBarry Smith   ierr = DMDAGetAdicMFArray(da,PETSC_FALSE,&ad_f,&ad_fstart,&tdof);CHKERRQ(ierr);
96047c6ae99SBarry Smith 
96147c6ae99SBarry Smith   /* copy input vector into derivative object */
96247c6ae99SBarry Smith   ierr = VecGetArray(vu,&avu);CHKERRQ(ierr);
96347c6ae99SBarry Smith   ierr = VecGetArray(v,&av);CHKERRQ(ierr);
96447c6ae99SBarry Smith   for (i=0; i<gtdof; i++) {
96547c6ae99SBarry Smith     ad_vustart[2*i]   = avu[i];
96647c6ae99SBarry Smith     ad_vustart[2*i+1] = av[i];
96747c6ae99SBarry Smith   }
96847c6ae99SBarry Smith   ierr = VecRestoreArray(vu,&avu);CHKERRQ(ierr);
96947c6ae99SBarry Smith   ierr = VecRestoreArray(v,&av);CHKERRQ(ierr);
97047c6ae99SBarry Smith 
97147c6ae99SBarry Smith   PetscADResetIndep();
97247c6ae99SBarry Smith   ierr = PetscADIncrementTotalGradSize(1);CHKERRQ(ierr);
97347c6ae99SBarry Smith   PetscADSetIndepDone();
97447c6ae99SBarry Smith 
97547c6ae99SBarry Smith   ierr = (*dd->adicmf_lf)(&info,ad_vu,ad_f,w);CHKERRQ(ierr);
97647c6ae99SBarry Smith 
97747c6ae99SBarry Smith   /* stick the values into the vector */
97847c6ae99SBarry Smith   ierr = VecGetArray(f,&af);CHKERRQ(ierr);
97947c6ae99SBarry Smith   for (i=0; i<tdof; i++) {
98047c6ae99SBarry Smith     af[i] = ad_fstart[2*i+1];
98147c6ae99SBarry Smith   }
98247c6ae99SBarry Smith   ierr = VecRestoreArray(f,&af);CHKERRQ(ierr);
98347c6ae99SBarry Smith 
98447c6ae99SBarry Smith   /* return space for derivative objects.  */
985aa219208SBarry Smith   ierr = DMDARestoreAdicMFArray(da,PETSC_TRUE,&ad_vu,&ad_vustart,&gtdof);CHKERRQ(ierr);
986aa219208SBarry Smith   ierr = DMDARestoreAdicMFArray(da,PETSC_FALSE,&ad_f,&ad_fstart,&tdof);CHKERRQ(ierr);
98747c6ae99SBarry Smith   PetscFunctionReturn(0);
98847c6ae99SBarry Smith }
98947c6ae99SBarry Smith #endif
99047c6ae99SBarry Smith 
99147c6ae99SBarry Smith #undef __FUNCT__
992aa219208SBarry Smith #define __FUNCT__ "DMDAComputeJacobian1"
99347c6ae99SBarry Smith /*@
994aa219208SBarry Smith     DMDAComputeJacobian1 - Evaluates a local Jacobian function on each processor that
995aa219208SBarry Smith         share a DMDA
99647c6ae99SBarry Smith 
99747c6ae99SBarry Smith    Input Parameters:
9989a42bb27SBarry Smith +    da - the DM that defines the grid
99947c6ae99SBarry Smith .    vu - input vector (ghosted)
100047c6ae99SBarry Smith .    J - output matrix
100147c6ae99SBarry Smith -    w - any user data
100247c6ae99SBarry Smith 
100347c6ae99SBarry Smith     Notes: Does NOT do ghost updates on vu upon entry
100447c6ae99SBarry Smith 
100547c6ae99SBarry Smith     Level: advanced
100647c6ae99SBarry Smith 
1007aa219208SBarry Smith .seealso: DMDAFormFunction1()
100847c6ae99SBarry Smith 
100947c6ae99SBarry Smith @*/
10107087cfbeSBarry Smith PetscErrorCode  DMDAComputeJacobian1(DM da,Vec vu,Mat J,void *w)
101147c6ae99SBarry Smith {
101247c6ae99SBarry Smith   PetscErrorCode ierr;
101347c6ae99SBarry Smith   void           *u;
1014aa219208SBarry Smith   DMDALocalInfo  info;
101547c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
101647c6ae99SBarry Smith 
101747c6ae99SBarry Smith   PetscFunctionBegin;
1018aa219208SBarry Smith   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
1019aa219208SBarry Smith   ierr = DMDAVecGetArray(da,vu,&u);CHKERRQ(ierr);
102047c6ae99SBarry Smith   ierr = (*dd->lj)(&info,u,J,w);CHKERRQ(ierr);
1021aa219208SBarry Smith   ierr = DMDAVecRestoreArray(da,vu,&u);CHKERRQ(ierr);
102247c6ae99SBarry Smith   PetscFunctionReturn(0);
102347c6ae99SBarry Smith }
102447c6ae99SBarry Smith 
102547c6ae99SBarry Smith 
102647c6ae99SBarry Smith #undef __FUNCT__
1027aa219208SBarry Smith #define __FUNCT__ "DMDAComputeJacobian1WithAdifor"
102847c6ae99SBarry Smith /*
1029aa219208SBarry Smith     DMDAComputeJacobian1WithAdifor - Evaluates a ADIFOR provided Jacobian local function on each processor that
1030aa219208SBarry Smith         share a DMDA
103147c6ae99SBarry Smith 
103247c6ae99SBarry Smith    Input Parameters:
10339a42bb27SBarry Smith +    da - the DM that defines the grid
103447c6ae99SBarry Smith .    vu - input vector (ghosted)
103547c6ae99SBarry Smith .    J - output matrix
103647c6ae99SBarry Smith -    w - any user data
103747c6ae99SBarry Smith 
103847c6ae99SBarry Smith     Notes: Does NOT do ghost updates on vu upon entry
103947c6ae99SBarry Smith 
1040aa219208SBarry Smith .seealso: DMDAFormFunction1()
104147c6ae99SBarry Smith 
104247c6ae99SBarry Smith */
10437087cfbeSBarry Smith PetscErrorCode  DMDAComputeJacobian1WithAdifor(DM da,Vec vu,Mat J,void *w)
104447c6ae99SBarry Smith {
104547c6ae99SBarry Smith   PetscErrorCode  ierr;
104647c6ae99SBarry Smith   PetscInt        i,Nc,N;
104747c6ae99SBarry Smith   ISColoringValue *color;
1048aa219208SBarry Smith   DMDALocalInfo   info;
104947c6ae99SBarry Smith   PetscScalar     *u,*g_u,*g_f,*f = 0,*p_u;
105047c6ae99SBarry Smith   ISColoring      iscoloring;
105147c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
1052aa219208SBarry Smith   void            (*lf)(PetscInt*,DMDALocalInfo*,PetscScalar*,PetscScalar*,PetscInt*,PetscScalar*,PetscScalar*,PetscInt*,void*,PetscErrorCode*) =
1053aa219208SBarry Smith                   (void (*)(PetscInt*,DMDALocalInfo*,PetscScalar*,PetscScalar*,PetscInt*,PetscScalar*,PetscScalar*,PetscInt*,void*,PetscErrorCode*))*dd->adifor_lf;
105447c6ae99SBarry Smith 
105547c6ae99SBarry Smith   PetscFunctionBegin;
105694013140SBarry Smith   ierr = DMGetColoring(da,IS_COLORING_GHOSTED,MATAIJ,&iscoloring);CHKERRQ(ierr);
105747c6ae99SBarry Smith   Nc   = iscoloring->n;
1058aa219208SBarry Smith   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
105947c6ae99SBarry Smith   N    = info.gxm*info.gym*info.gzm*info.dof;
106047c6ae99SBarry Smith 
106147c6ae99SBarry Smith   /* get space for derivative objects.  */
106247c6ae99SBarry Smith   ierr  = PetscMalloc(Nc*info.gxm*info.gym*info.gzm*info.dof*sizeof(PetscScalar),&g_u);CHKERRQ(ierr);
106347c6ae99SBarry Smith   ierr  = PetscMemzero(g_u,Nc*info.gxm*info.gym*info.gzm*info.dof*sizeof(PetscScalar));CHKERRQ(ierr);
106447c6ae99SBarry Smith   p_u   = g_u;
106547c6ae99SBarry Smith   color = iscoloring->colors;
106647c6ae99SBarry Smith   for (i=0; i<N; i++) {
106747c6ae99SBarry Smith     p_u[*color++] = 1.0;
106847c6ae99SBarry Smith     p_u          += Nc;
106947c6ae99SBarry Smith   }
1070*fcfd50ebSBarry Smith   ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr);
107147c6ae99SBarry 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);
107247c6ae99SBarry Smith 
107347c6ae99SBarry Smith   /* Seed the input array g_u with coloring information */
107447c6ae99SBarry Smith 
107547c6ae99SBarry Smith   ierr = VecGetArray(vu,&u);CHKERRQ(ierr);
107647c6ae99SBarry Smith   (lf)(&Nc,&info,u,g_u,&Nc,f,g_f,&Nc,w,&ierr);CHKERRQ(ierr);
107747c6ae99SBarry Smith   ierr = VecRestoreArray(vu,&u);CHKERRQ(ierr);
107847c6ae99SBarry Smith 
107947c6ae99SBarry Smith   /* stick the values into the matrix */
108047c6ae99SBarry Smith   /* PetscScalarView(Nc*info.xm*info.ym,g_f,0); */
108147c6ae99SBarry Smith   ierr = MatSetValuesAdifor(J,Nc,g_f);CHKERRQ(ierr);
108247c6ae99SBarry Smith 
108347c6ae99SBarry Smith   /* return space for derivative objects.  */
108447c6ae99SBarry Smith   ierr = PetscFree(g_u);CHKERRQ(ierr);
108547c6ae99SBarry Smith   ierr = PetscFree2(g_f,f);CHKERRQ(ierr);
108647c6ae99SBarry Smith   PetscFunctionReturn(0);
108747c6ae99SBarry Smith }
108847c6ae99SBarry Smith 
108947c6ae99SBarry Smith #undef __FUNCT__
1090aa219208SBarry Smith #define __FUNCT__ "DMDAFormJacobianLocal"
109147c6ae99SBarry Smith /*@C
1092aa219208SBarry Smith    DMDAFormjacobianLocal - This is a universal Jacobian evaluation routine for
10939a42bb27SBarry Smith    a local DM function.
109447c6ae99SBarry Smith 
1095aa219208SBarry Smith    Collective on DMDA
109647c6ae99SBarry Smith 
109747c6ae99SBarry Smith    Input Parameters:
10989a42bb27SBarry Smith +  da - the DM context
109947c6ae99SBarry Smith .  func - The local function
110047c6ae99SBarry Smith .  X - input vector
110147c6ae99SBarry Smith .  J - Jacobian matrix
110247c6ae99SBarry Smith -  ctx - A user context
110347c6ae99SBarry Smith 
110447c6ae99SBarry Smith    Level: intermediate
110547c6ae99SBarry Smith 
1106aa219208SBarry Smith .seealso: DMDASetLocalFunction(), DMDASetLocalJacobian(), DMDASetLocalAdicFunction(), DMDASetLocalAdicMFFunction(),
110747c6ae99SBarry Smith           SNESSetFunction(), SNESSetJacobian()
110847c6ae99SBarry Smith 
110947c6ae99SBarry Smith @*/
11107087cfbeSBarry Smith PetscErrorCode  DMDAFormJacobianLocal(DM da, DMDALocalFunction1 func, Vec X, Mat J, void *ctx)
111147c6ae99SBarry Smith {
111247c6ae99SBarry Smith   Vec            localX;
1113aa219208SBarry Smith   DMDALocalInfo  info;
111447c6ae99SBarry Smith   void           *u;
111547c6ae99SBarry Smith   PetscErrorCode ierr;
111647c6ae99SBarry Smith 
111747c6ae99SBarry Smith   PetscFunctionBegin;
11189a42bb27SBarry Smith   ierr = DMGetLocalVector(da,&localX);CHKERRQ(ierr);
111947c6ae99SBarry Smith   /*
112047c6ae99SBarry Smith      Scatter ghost points to local vector, using the 2-step process
11219a42bb27SBarry Smith         DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
112247c6ae99SBarry Smith   */
11239a42bb27SBarry Smith   ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
11249a42bb27SBarry Smith   ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
1125aa219208SBarry Smith   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
1126aa219208SBarry Smith   ierr = DMDAVecGetArray(da,localX,&u);CHKERRQ(ierr);
112739d508bbSBarry Smith   ierr = (*func)(&info,u,J,ctx);CHKERRQ(ierr);
1128aa219208SBarry Smith   ierr = DMDAVecRestoreArray(da,localX,&u);CHKERRQ(ierr);
11299a42bb27SBarry Smith   ierr = DMRestoreLocalVector(da,&localX);CHKERRQ(ierr);
113047c6ae99SBarry Smith   PetscFunctionReturn(0);
113147c6ae99SBarry Smith }
113247c6ae99SBarry Smith 
113347c6ae99SBarry Smith #undef __FUNCT__
1134aa219208SBarry Smith #define __FUNCT__ "DMDAMultiplyByJacobian1WithAD"
113547c6ae99SBarry Smith /*@C
1136aa219208SBarry Smith     DMDAMultiplyByJacobian1WithAD - Applies a Jacobian function supplied by ADIFOR or ADIC
1137aa219208SBarry Smith     to a vector on each processor that shares a DMDA.
113847c6ae99SBarry Smith 
113947c6ae99SBarry Smith    Input Parameters:
11409a42bb27SBarry Smith +    da - the DM that defines the grid
114147c6ae99SBarry Smith .    vu - Jacobian is computed at this point (ghosted)
114247c6ae99SBarry Smith .    v - product is done on this vector (ghosted)
114347c6ae99SBarry Smith .    fu - output vector = J(vu)*v (not ghosted)
114447c6ae99SBarry Smith -    w - any user data
114547c6ae99SBarry Smith 
114647c6ae99SBarry Smith     Notes:
114747c6ae99SBarry Smith     This routine does NOT do ghost updates on vu and v upon entry.
114847c6ae99SBarry Smith 
1149aa219208SBarry Smith     Automatically calls DMDAMultiplyByJacobian1WithAdifor() or DMDAMultiplyByJacobian1WithAdic()
1150aa219208SBarry Smith     depending on whether DMDASetLocalAdicMFFunction() or DMDASetLocalAdiforMFFunction() was called.
115147c6ae99SBarry Smith 
115247c6ae99SBarry Smith    Level: advanced
115347c6ae99SBarry Smith 
1154aa219208SBarry Smith .seealso: DMDAFormFunction1(), DMDAMultiplyByJacobian1WithAdifor(), DMDAMultiplyByJacobian1WithAdic()
115547c6ae99SBarry Smith 
115647c6ae99SBarry Smith @*/
11577087cfbeSBarry Smith PetscErrorCode  DMDAMultiplyByJacobian1WithAD(DM da,Vec u,Vec v,Vec f,void *w)
115847c6ae99SBarry Smith {
115947c6ae99SBarry Smith   PetscErrorCode ierr;
116047c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
116147c6ae99SBarry Smith 
116247c6ae99SBarry Smith   PetscFunctionBegin;
116347c6ae99SBarry Smith   if (dd->adicmf_lf) {
116447c6ae99SBarry Smith #if defined(PETSC_HAVE_ADIC)
1165aa219208SBarry Smith     ierr = DMDAMultiplyByJacobian1WithAdic(da,u,v,f,w);CHKERRQ(ierr);
116647c6ae99SBarry Smith #else
116747c6ae99SBarry Smith     SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP_SYS,"Requires ADIC to be installed and cannot use complex numbers");
116847c6ae99SBarry Smith #endif
116947c6ae99SBarry Smith   } else if (dd->adiformf_lf) {
1170aa219208SBarry Smith     ierr = DMDAMultiplyByJacobian1WithAdifor(da,u,v,f,w);CHKERRQ(ierr);
117147c6ae99SBarry Smith   } else {
1172aa219208SBarry Smith     SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ORDER,"Must call DMDASetLocalAdiforMFFunction() or DMDASetLocalAdicMFFunction() before using");
117347c6ae99SBarry Smith   }
117447c6ae99SBarry Smith   PetscFunctionReturn(0);
117547c6ae99SBarry Smith }
117647c6ae99SBarry Smith 
117747c6ae99SBarry Smith 
117847c6ae99SBarry Smith #undef __FUNCT__
1179aa219208SBarry Smith #define __FUNCT__ "DMDAMultiplyByJacobian1WithAdifor"
118047c6ae99SBarry Smith /*@C
1181aa219208SBarry Smith     DMDAMultiplyByJacobian1WithAdifor - Applies a ADIFOR provided Jacobian function on each processor that
11829a42bb27SBarry Smith         share a DM to a vector
118347c6ae99SBarry Smith 
118447c6ae99SBarry Smith    Input Parameters:
11859a42bb27SBarry Smith +    da - the DM that defines the grid
118647c6ae99SBarry Smith .    vu - Jacobian is computed at this point (ghosted)
118747c6ae99SBarry Smith .    v - product is done on this vector (ghosted)
118847c6ae99SBarry Smith .    fu - output vector = J(vu)*v (not ghosted)
118947c6ae99SBarry Smith -    w - any user data
119047c6ae99SBarry Smith 
119147c6ae99SBarry Smith     Notes: Does NOT do ghost updates on vu and v upon entry
119247c6ae99SBarry Smith 
119347c6ae99SBarry Smith    Level: advanced
119447c6ae99SBarry Smith 
1195aa219208SBarry Smith .seealso: DMDAFormFunction1()
119647c6ae99SBarry Smith 
119747c6ae99SBarry Smith @*/
11987087cfbeSBarry Smith PetscErrorCode  DMDAMultiplyByJacobian1WithAdifor(DM da,Vec u,Vec v,Vec f,void *w)
119947c6ae99SBarry Smith {
120047c6ae99SBarry Smith   PetscErrorCode ierr;
120147c6ae99SBarry Smith   PetscScalar    *au,*av,*af,*awork;
120247c6ae99SBarry Smith   Vec            work;
1203aa219208SBarry Smith   DMDALocalInfo  info;
120447c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
1205aa219208SBarry Smith   void           (*lf)(DMDALocalInfo*,PetscScalar*,PetscScalar*,PetscScalar*,PetscScalar*,void*,PetscErrorCode*) =
1206aa219208SBarry Smith                  (void (*)(DMDALocalInfo*,PetscScalar*,PetscScalar*,PetscScalar*,PetscScalar*,void*,PetscErrorCode*))*dd->adiformf_lf;
120747c6ae99SBarry Smith 
120847c6ae99SBarry Smith   PetscFunctionBegin;
1209aa219208SBarry Smith   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
121047c6ae99SBarry Smith 
12119a42bb27SBarry Smith   ierr = DMGetGlobalVector(da,&work);CHKERRQ(ierr);
121247c6ae99SBarry Smith   ierr = VecGetArray(u,&au);CHKERRQ(ierr);
121347c6ae99SBarry Smith   ierr = VecGetArray(v,&av);CHKERRQ(ierr);
121447c6ae99SBarry Smith   ierr = VecGetArray(f,&af);CHKERRQ(ierr);
121547c6ae99SBarry Smith   ierr = VecGetArray(work,&awork);CHKERRQ(ierr);
121647c6ae99SBarry Smith   (lf)(&info,au,av,awork,af,w,&ierr);CHKERRQ(ierr);
121747c6ae99SBarry Smith   ierr = VecRestoreArray(u,&au);CHKERRQ(ierr);
121847c6ae99SBarry Smith   ierr = VecRestoreArray(v,&av);CHKERRQ(ierr);
121947c6ae99SBarry Smith   ierr = VecRestoreArray(f,&af);CHKERRQ(ierr);
122047c6ae99SBarry Smith   ierr = VecRestoreArray(work,&awork);CHKERRQ(ierr);
12219a42bb27SBarry Smith   ierr = DMRestoreGlobalVector(da,&work);CHKERRQ(ierr);
122247c6ae99SBarry Smith 
122347c6ae99SBarry Smith   PetscFunctionReturn(0);
122447c6ae99SBarry Smith }
122547c6ae99SBarry Smith 
122647c6ae99SBarry Smith #undef __FUNCT__
12279a42bb27SBarry Smith #define __FUNCT__ "DMSetUp_DA_2D"
12287087cfbeSBarry Smith PetscErrorCode  DMSetUp_DA_2D(DM da)
122947c6ae99SBarry Smith {
123047c6ae99SBarry Smith   DM_DA                  *dd = (DM_DA*)da->data;
123147c6ae99SBarry Smith   const PetscInt         M            = dd->M;
123247c6ae99SBarry Smith   const PetscInt         N            = dd->N;
123347c6ae99SBarry Smith   PetscInt               m            = dd->m;
123447c6ae99SBarry Smith   PetscInt               n            = dd->n;
123547c6ae99SBarry Smith   const PetscInt         dof          = dd->w;
123647c6ae99SBarry Smith   const PetscInt         s            = dd->s;
12371321219cSEthan Coon   const DMDABoundaryType bx         = dd->bx;
12381321219cSEthan Coon   const DMDABoundaryType by         = dd->by;
1239aa219208SBarry Smith   const DMDAStencilType  stencil_type = dd->stencil_type;
124047c6ae99SBarry Smith   PetscInt               *lx           = dd->lx;
124147c6ae99SBarry Smith   PetscInt               *ly           = dd->ly;
124247c6ae99SBarry Smith   MPI_Comm               comm;
124347c6ae99SBarry Smith   PetscMPIInt            rank,size;
1244ce00eea3SSatish Balay   PetscInt               xs,xe,ys,ye,x,y,Xs,Xe,Ys,Ye,start,end,IXs,IXe,IYs,IYe;
1245ce00eea3SSatish Balay   PetscInt               up,down,left,right,i,n0,n1,n2,n3,n5,n6,n7,n8,*idx,nn,*idx_cpy;
1246ce00eea3SSatish Balay   const PetscInt         *idx_full;
1247db87c5ecSEthan Coon   PetscInt               xbase,*bases,*ldims,j,x_t,y_t,s_t,base,count;
124847c6ae99SBarry Smith   PetscInt               s_x,s_y; /* s proportionalized to w */
124947c6ae99SBarry Smith   PetscInt               sn0 = 0,sn2 = 0,sn6 = 0,sn8 = 0;
125047c6ae99SBarry Smith   Vec                    local,global;
125147c6ae99SBarry Smith   VecScatter             ltog,gtol;
1252ce00eea3SSatish Balay   IS                     to,from,ltogis;
125347c6ae99SBarry Smith   PetscErrorCode         ierr;
125447c6ae99SBarry Smith 
125547c6ae99SBarry Smith   PetscFunctionBegin;
125647c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
125747c6ae99SBarry Smith 
125847c6ae99SBarry Smith   if (dof < 1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Must have 1 or more degrees of freedom per node: %D",dof);
125947c6ae99SBarry Smith   if (s < 0) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Stencil width cannot be negative: %D",s);
126047c6ae99SBarry Smith 
126147c6ae99SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
126247c6ae99SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
126347c6ae99SBarry Smith 
126447c6ae99SBarry Smith   dd->dim         = 2;
126547c6ae99SBarry Smith   ierr = PetscMalloc(dof*sizeof(char*),&dd->fieldname);CHKERRQ(ierr);
126647c6ae99SBarry Smith   ierr = PetscMemzero(dd->fieldname,dof*sizeof(char*));CHKERRQ(ierr);
126747c6ae99SBarry Smith 
126847c6ae99SBarry Smith   if (m != PETSC_DECIDE) {
126947c6ae99SBarry Smith     if (m < 1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in X direction: %D",m);
127047c6ae99SBarry Smith     else if (m > size) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in X direction: %D %d",m,size);
127147c6ae99SBarry Smith   }
127247c6ae99SBarry Smith   if (n != PETSC_DECIDE) {
127347c6ae99SBarry Smith     if (n < 1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Y direction: %D",n);
127447c6ae99SBarry Smith     else if (n > size) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Y direction: %D %d",n,size);
127547c6ae99SBarry Smith   }
127647c6ae99SBarry Smith 
127747c6ae99SBarry Smith   if (m == PETSC_DECIDE || n == PETSC_DECIDE) {
127847c6ae99SBarry Smith     if (n != PETSC_DECIDE) {
127947c6ae99SBarry Smith       m = size/n;
128047c6ae99SBarry Smith     } else if (m != PETSC_DECIDE) {
128147c6ae99SBarry Smith       n = size/m;
128247c6ae99SBarry Smith     } else {
128347c6ae99SBarry Smith       /* try for squarish distribution */
128447c6ae99SBarry Smith       m = (PetscInt)(0.5 + sqrt(((double)M)*((double)size)/((double)N)));
128547c6ae99SBarry Smith       if (!m) m = 1;
128647c6ae99SBarry Smith       while (m > 0) {
128747c6ae99SBarry Smith 	n = size/m;
128847c6ae99SBarry Smith 	if (m*n == size) break;
128947c6ae99SBarry Smith 	m--;
129047c6ae99SBarry Smith       }
129147c6ae99SBarry Smith       if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;}
129247c6ae99SBarry Smith     }
129347c6ae99SBarry 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 ");
129447c6ae99SBarry Smith   } else if (m*n != size) SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"Given Bad partition");
129547c6ae99SBarry Smith 
129647c6ae99SBarry Smith   if (M < m) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in x direction is too fine! %D %D",M,m);
129747c6ae99SBarry Smith   if (N < n) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in y direction is too fine! %D %D",N,n);
129847c6ae99SBarry Smith 
129947c6ae99SBarry Smith   /*
130047c6ae99SBarry Smith      Determine locally owned region
130147c6ae99SBarry Smith      xs is the first local node number, x is the number of local nodes
130247c6ae99SBarry Smith   */
130347c6ae99SBarry Smith   if (!lx) {
130447c6ae99SBarry Smith     ierr = PetscMalloc(m*sizeof(PetscInt), &dd->lx);CHKERRQ(ierr);
130547c6ae99SBarry Smith     lx = dd->lx;
130647c6ae99SBarry Smith     for (i=0; i<m; i++) {
130747c6ae99SBarry Smith       lx[i] = M/m + ((M % m) > i);
130847c6ae99SBarry Smith     }
130947c6ae99SBarry Smith   }
131047c6ae99SBarry Smith   x  = lx[rank % m];
131147c6ae99SBarry Smith   xs = 0;
131247c6ae99SBarry Smith   for (i=0; i<(rank % m); i++) {
131347c6ae99SBarry Smith     xs += lx[i];
131447c6ae99SBarry Smith   }
131547c6ae99SBarry Smith #if defined(PETSC_USE_DEBUG)
131647c6ae99SBarry Smith   left = xs;
131747c6ae99SBarry Smith   for (i=(rank % m); i<m; i++) {
131847c6ae99SBarry Smith     left += lx[i];
131947c6ae99SBarry Smith   }
132047c6ae99SBarry 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);
132147c6ae99SBarry Smith #endif
132247c6ae99SBarry Smith 
132347c6ae99SBarry Smith   /*
132447c6ae99SBarry Smith      Determine locally owned region
132547c6ae99SBarry Smith      ys is the first local node number, y is the number of local nodes
132647c6ae99SBarry Smith   */
132747c6ae99SBarry Smith   if (!ly) {
132847c6ae99SBarry Smith     ierr = PetscMalloc(n*sizeof(PetscInt), &dd->ly);CHKERRQ(ierr);
132947c6ae99SBarry Smith     ly = dd->ly;
133047c6ae99SBarry Smith     for (i=0; i<n; i++) {
133147c6ae99SBarry Smith       ly[i] = N/n + ((N % n) > i);
133247c6ae99SBarry Smith     }
133347c6ae99SBarry Smith   }
133447c6ae99SBarry Smith   y  = ly[rank/m];
133547c6ae99SBarry Smith   ys = 0;
133647c6ae99SBarry Smith   for (i=0; i<(rank/m); i++) {
133747c6ae99SBarry Smith     ys += ly[i];
133847c6ae99SBarry Smith   }
133947c6ae99SBarry Smith #if defined(PETSC_USE_DEBUG)
134047c6ae99SBarry Smith   left = ys;
134147c6ae99SBarry Smith   for (i=(rank/m); i<n; i++) {
134247c6ae99SBarry Smith     left += ly[i];
134347c6ae99SBarry Smith   }
134447c6ae99SBarry 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);
134547c6ae99SBarry Smith #endif
134647c6ae99SBarry Smith 
134747c6ae99SBarry Smith   if (x < s) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local x-width of domain x %D is smaller than stencil width s %D",x,s);
134847c6ae99SBarry Smith   if (y < s) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local y-width of domain y %D is smaller than stencil width s %D",y,s);
134947c6ae99SBarry Smith   xe = xs + x;
135047c6ae99SBarry Smith   ye = ys + y;
135147c6ae99SBarry Smith 
1352ce00eea3SSatish Balay   /* determine ghost region (Xs) and region scattered into (IXs)  */
135347c6ae99SBarry Smith   /* Assume No Periodicity */
1354ce00eea3SSatish Balay   if (xs-s > 0) { Xs = xs - s; IXs = xs - s; } else { Xs = 0; IXs = 0; }
1355ce00eea3SSatish Balay   if (xe+s <= M) { Xe = xe + s; IXe = xe + s; } else { Xe = M; IXe = M; }
1356ce00eea3SSatish Balay   if (ys-s > 0) { Ys = ys - s; IYs = ys - s; } else { Ys = 0; IYs = 0; }
1357ce00eea3SSatish Balay   if (ye+s <= N) { Ye = ye + s; IYe = ye + s; } else { Ye = N; IYe = N; }
135847c6ae99SBarry Smith 
1359ce00eea3SSatish Balay   /* fix for periodicity/ghosted */
13601321219cSEthan Coon   if (bx) { Xs = xs - s; Xe = xe + s; }
13611321219cSEthan Coon   if (bx == DMDA_BOUNDARY_PERIODIC) { IXs = xs - s; IXe = xe + s; }
13621321219cSEthan Coon   if (by) { Ys = ys - s; Ye = ye + s; }
13631321219cSEthan Coon   if (by == DMDA_BOUNDARY_PERIODIC) { IYs = ys - s; IYe = ye + s; }
136447c6ae99SBarry Smith 
136547c6ae99SBarry Smith   /* Resize all X parameters to reflect w */
1366ce00eea3SSatish Balay   s_x = s;
136747c6ae99SBarry Smith   s_y = s;
136847c6ae99SBarry Smith 
136947c6ae99SBarry Smith   /* determine starting point of each processor */
137047c6ae99SBarry Smith   nn    = x*y;
137147c6ae99SBarry Smith   ierr  = PetscMalloc2(size+1,PetscInt,&bases,size,PetscInt,&ldims);CHKERRQ(ierr);
137247c6ae99SBarry Smith   ierr  = MPI_Allgather(&nn,1,MPIU_INT,ldims,1,MPIU_INT,comm);CHKERRQ(ierr);
137347c6ae99SBarry Smith   bases[0] = 0;
137447c6ae99SBarry Smith   for (i=1; i<=size; i++) {
137547c6ae99SBarry Smith     bases[i] = ldims[i-1];
137647c6ae99SBarry Smith   }
137747c6ae99SBarry Smith   for (i=1; i<=size; i++) {
137847c6ae99SBarry Smith     bases[i] += bases[i-1];
137947c6ae99SBarry Smith   }
1380ce00eea3SSatish Balay   base = bases[rank]*dof;
138147c6ae99SBarry Smith 
138247c6ae99SBarry Smith   /* allocate the base parallel and sequential vectors */
1383ce00eea3SSatish Balay   dd->Nlocal = x*y*dof;
138447c6ae99SBarry Smith   ierr = VecCreateMPIWithArray(comm,dd->Nlocal,PETSC_DECIDE,0,&global);CHKERRQ(ierr);
138547c6ae99SBarry Smith   ierr = VecSetBlockSize(global,dof);CHKERRQ(ierr);
1386ce00eea3SSatish Balay   dd->nlocal = (Xe-Xs)*(Ye-Ys)*dof;
138747c6ae99SBarry Smith   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,dd->nlocal,0,&local);CHKERRQ(ierr);
138847c6ae99SBarry Smith   ierr = VecSetBlockSize(local,dof);CHKERRQ(ierr);
138947c6ae99SBarry Smith 
139047c6ae99SBarry Smith   /* generate appropriate vector scatters */
139147c6ae99SBarry Smith   /* local to global inserts non-ghost point region into global */
139247c6ae99SBarry Smith   ierr = VecGetOwnershipRange(global,&start,&end);CHKERRQ(ierr);
1393ce00eea3SSatish Balay   ierr = ISCreateStride(comm,x*y*dof,start,1,&to);CHKERRQ(ierr);
139447c6ae99SBarry Smith 
1395db87c5ecSEthan Coon   count = x*y;
1396ce00eea3SSatish Balay   ierr = PetscMalloc(x*y*sizeof(PetscInt),&idx);CHKERRQ(ierr);
1397ce00eea3SSatish Balay   left = xs - Xs; right = left + x;
1398ce00eea3SSatish Balay   down = ys - Ys; up = down + y;
139947c6ae99SBarry Smith   count = 0;
140047c6ae99SBarry Smith   for (i=down; i<up; i++) {
1401ce00eea3SSatish Balay     for (j=left; j<right; j++) {
1402ce00eea3SSatish Balay       idx[count++] = i*(Xe-Xs) + j;
140347c6ae99SBarry Smith     }
140447c6ae99SBarry Smith   }
140547c6ae99SBarry Smith 
1406ce00eea3SSatish Balay   ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&from);CHKERRQ(ierr);
140747c6ae99SBarry Smith   ierr = VecScatterCreate(local,from,global,to,&ltog);CHKERRQ(ierr);
140847c6ae99SBarry Smith   ierr = PetscLogObjectParent(dd,ltog);CHKERRQ(ierr);
1409*fcfd50ebSBarry Smith   ierr = ISDestroy(&from);CHKERRQ(ierr);
1410*fcfd50ebSBarry Smith   ierr = ISDestroy(&to);CHKERRQ(ierr);
141147c6ae99SBarry Smith 
1412ce00eea3SSatish Balay   /* global to local must include ghost points within the domain,
1413ce00eea3SSatish Balay      but not ghost points outside the domain that aren't periodic */
1414aa219208SBarry Smith   if (stencil_type == DMDA_STENCIL_BOX) {
1415db87c5ecSEthan Coon     count = (IXe-IXs)*(IYe-IYs);
1416db87c5ecSEthan Coon     ierr  = PetscMalloc(count*sizeof(PetscInt),&idx);CHKERRQ(ierr);
1417ce00eea3SSatish Balay 
1418ce00eea3SSatish Balay     left = IXs - Xs; right = left + (IXe-IXs);
1419ce00eea3SSatish Balay     down = IYs - Ys; up = down + (IYe-IYs);
1420ce00eea3SSatish Balay     count = 0;
1421ce00eea3SSatish Balay     for (i=down; i<up; i++) {
1422ce00eea3SSatish Balay       for (j=left; j<right; j++) {
1423ce00eea3SSatish Balay         idx[count++] = j + i*(Xe-Xs);
1424ce00eea3SSatish Balay       }
1425ce00eea3SSatish Balay     }
1426ce00eea3SSatish Balay     ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr);
1427ce00eea3SSatish Balay 
142847c6ae99SBarry Smith   } else {
142947c6ae99SBarry Smith     /* must drop into cross shape region */
143047c6ae99SBarry Smith     /*       ---------|
143147c6ae99SBarry Smith             |  top    |
1432ce00eea3SSatish Balay          |---         ---| up
143347c6ae99SBarry Smith          |   middle      |
143447c6ae99SBarry Smith          |               |
1435ce00eea3SSatish Balay          ----         ---- down
143647c6ae99SBarry Smith             | bottom  |
143747c6ae99SBarry Smith             -----------
143847c6ae99SBarry Smith          Xs xs        xe Xe */
1439db87c5ecSEthan Coon     count = (ys-IYs)*x + y*(IXe-IXs) + (IYe-ye)*x;
1440db87c5ecSEthan Coon     ierr  = PetscMalloc(count*sizeof(PetscInt),&idx);CHKERRQ(ierr);
1441ce00eea3SSatish Balay 
1442ce00eea3SSatish Balay     left = xs - Xs; right = left + x;
1443ce00eea3SSatish Balay     down = ys - Ys; up = down + y;
144447c6ae99SBarry Smith     count = 0;
1445ce00eea3SSatish Balay     /* bottom */
1446ce00eea3SSatish Balay     for (i=(IYs-Ys); i<down; i++) {
1447ce00eea3SSatish Balay       for (j=left; j<right; j++) {
1448ce00eea3SSatish Balay         idx[count++] = j + i*(Xe-Xs);
144947c6ae99SBarry Smith       }
145047c6ae99SBarry Smith     }
145147c6ae99SBarry Smith     /* middle */
145247c6ae99SBarry Smith     for (i=down; i<up; i++) {
1453ce00eea3SSatish Balay       for (j=(IXs-Xs); j<(IXe-Xs); j++) {
1454ce00eea3SSatish Balay         idx[count++] = j + i*(Xe-Xs);
145547c6ae99SBarry Smith       }
145647c6ae99SBarry Smith     }
145747c6ae99SBarry Smith     /* top */
1458ce00eea3SSatish Balay     for (i=up; i<up+IYe-ye; i++) {
1459ce00eea3SSatish Balay       for (j=left; j<right; j++) {
1460ce00eea3SSatish Balay         idx[count++] = j + i*(Xe-Xs);
146147c6ae99SBarry Smith       }
146247c6ae99SBarry Smith     }
146347c6ae99SBarry Smith     ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr);
146447c6ae99SBarry Smith   }
146547c6ae99SBarry Smith 
146647c6ae99SBarry Smith 
146747c6ae99SBarry Smith   /* determine who lies on each side of us stored in    n6 n7 n8
146847c6ae99SBarry Smith                                                         n3    n5
146947c6ae99SBarry Smith                                                         n0 n1 n2
147047c6ae99SBarry Smith   */
147147c6ae99SBarry Smith 
147247c6ae99SBarry Smith   /* Assume the Non-Periodic Case */
147347c6ae99SBarry Smith   n1 = rank - m;
147447c6ae99SBarry Smith   if (rank % m) {
147547c6ae99SBarry Smith     n0 = n1 - 1;
147647c6ae99SBarry Smith   } else {
147747c6ae99SBarry Smith     n0 = -1;
147847c6ae99SBarry Smith   }
147947c6ae99SBarry Smith   if ((rank+1) % m) {
148047c6ae99SBarry Smith     n2 = n1 + 1;
148147c6ae99SBarry Smith     n5 = rank + 1;
148247c6ae99SBarry Smith     n8 = rank + m + 1; if (n8 >= m*n) n8 = -1;
148347c6ae99SBarry Smith   } else {
148447c6ae99SBarry Smith     n2 = -1; n5 = -1; n8 = -1;
148547c6ae99SBarry Smith   }
148647c6ae99SBarry Smith   if (rank % m) {
148747c6ae99SBarry Smith     n3 = rank - 1;
148847c6ae99SBarry Smith     n6 = n3 + m; if (n6 >= m*n) n6 = -1;
148947c6ae99SBarry Smith   } else {
149047c6ae99SBarry Smith     n3 = -1; n6 = -1;
149147c6ae99SBarry Smith   }
149247c6ae99SBarry Smith   n7 = rank + m; if (n7 >= m*n) n7 = -1;
149347c6ae99SBarry Smith 
14941321219cSEthan Coon   if (bx == DMDA_BOUNDARY_PERIODIC && by == DMDA_BOUNDARY_PERIODIC) {
149547c6ae99SBarry Smith   /* Modify for Periodic Cases */
149647c6ae99SBarry Smith     /* Handle all four corners */
149747c6ae99SBarry Smith     if ((n6 < 0) && (n7 < 0) && (n3 < 0)) n6 = m-1;
149847c6ae99SBarry Smith     if ((n8 < 0) && (n7 < 0) && (n5 < 0)) n8 = 0;
149947c6ae99SBarry Smith     if ((n2 < 0) && (n5 < 0) && (n1 < 0)) n2 = size-m;
150047c6ae99SBarry Smith     if ((n0 < 0) && (n3 < 0) && (n1 < 0)) n0 = size-1;
150147c6ae99SBarry Smith 
150247c6ae99SBarry Smith     /* Handle Top and Bottom Sides */
150347c6ae99SBarry Smith     if (n1 < 0) n1 = rank + m * (n-1);
150447c6ae99SBarry Smith     if (n7 < 0) n7 = rank - m * (n-1);
150547c6ae99SBarry Smith     if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1;
150647c6ae99SBarry Smith     if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1;
150747c6ae99SBarry Smith     if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1;
150847c6ae99SBarry Smith     if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1;
150947c6ae99SBarry Smith 
151047c6ae99SBarry Smith     /* Handle Left and Right Sides */
151147c6ae99SBarry Smith     if (n3 < 0) n3 = rank + (m-1);
151247c6ae99SBarry Smith     if (n5 < 0) n5 = rank - (m-1);
151347c6ae99SBarry Smith     if ((n1 >= 0) && (n0 < 0)) n0 = rank-1;
151447c6ae99SBarry Smith     if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1;
151547c6ae99SBarry Smith     if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1;
151647c6ae99SBarry Smith     if ((n7 >= 0) && (n8 < 0)) n8 = rank+1;
15171321219cSEthan Coon   } else if (by == DMDA_BOUNDARY_PERIODIC) {  /* Handle Top and Bottom Sides */
1518ce00eea3SSatish Balay     if (n1 < 0) n1 = rank + m * (n-1);
1519ce00eea3SSatish Balay     if (n7 < 0) n7 = rank - m * (n-1);
1520ce00eea3SSatish Balay     if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1;
1521ce00eea3SSatish Balay     if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1;
1522ce00eea3SSatish Balay     if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1;
1523ce00eea3SSatish Balay     if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1;
15241321219cSEthan Coon   } else if (bx == DMDA_BOUNDARY_PERIODIC) { /* Handle Left and Right Sides */
1525ce00eea3SSatish Balay     if (n3 < 0) n3 = rank + (m-1);
1526ce00eea3SSatish Balay     if (n5 < 0) n5 = rank - (m-1);
1527ce00eea3SSatish Balay     if ((n1 >= 0) && (n0 < 0)) n0 = rank-1;
1528ce00eea3SSatish Balay     if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1;
1529ce00eea3SSatish Balay     if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1;
1530ce00eea3SSatish Balay     if ((n7 >= 0) && (n8 < 0)) n8 = rank+1;
153147c6ae99SBarry Smith   }
1532ce00eea3SSatish Balay 
153347c6ae99SBarry Smith   ierr = PetscMalloc(9*sizeof(PetscInt),&dd->neighbors);CHKERRQ(ierr);
153447c6ae99SBarry Smith   dd->neighbors[0] = n0;
153547c6ae99SBarry Smith   dd->neighbors[1] = n1;
153647c6ae99SBarry Smith   dd->neighbors[2] = n2;
153747c6ae99SBarry Smith   dd->neighbors[3] = n3;
153847c6ae99SBarry Smith   dd->neighbors[4] = rank;
153947c6ae99SBarry Smith   dd->neighbors[5] = n5;
154047c6ae99SBarry Smith   dd->neighbors[6] = n6;
154147c6ae99SBarry Smith   dd->neighbors[7] = n7;
154247c6ae99SBarry Smith   dd->neighbors[8] = n8;
154347c6ae99SBarry Smith 
1544aa219208SBarry Smith   if (stencil_type == DMDA_STENCIL_STAR) {
154547c6ae99SBarry Smith     /* save corner processor numbers */
154647c6ae99SBarry Smith     sn0 = n0; sn2 = n2; sn6 = n6; sn8 = n8;
154747c6ae99SBarry Smith     n0 = n2 = n6 = n8 = -1;
154847c6ae99SBarry Smith   }
154947c6ae99SBarry Smith 
1550ce00eea3SSatish Balay   ierr = PetscMalloc((Xe-Xs)*(Ye-Ys)*sizeof(PetscInt),&idx);CHKERRQ(ierr);
1551ce00eea3SSatish Balay   ierr = PetscLogObjectMemory(da,(Xe-Xs)*(Ye-Ys)*sizeof(PetscInt));CHKERRQ(ierr);
155247c6ae99SBarry Smith 
1553ce00eea3SSatish Balay   nn = 0;
155447c6ae99SBarry Smith   xbase = bases[rank];
155547c6ae99SBarry Smith   for (i=1; i<=s_y; i++) {
155647c6ae99SBarry Smith     if (n0 >= 0) { /* left below */
1557ce00eea3SSatish Balay       x_t = lx[n0 % m];
155847c6ae99SBarry Smith       y_t = ly[(n0/m)];
155947c6ae99SBarry Smith       s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x;
156047c6ae99SBarry Smith       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
156147c6ae99SBarry Smith     }
156247c6ae99SBarry Smith     if (n1 >= 0) { /* directly below */
156347c6ae99SBarry Smith       x_t = x;
156447c6ae99SBarry Smith       y_t = ly[(n1/m)];
156547c6ae99SBarry Smith       s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t;
156647c6ae99SBarry Smith       for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
156747c6ae99SBarry Smith     }
156847c6ae99SBarry Smith     if (n2 >= 0) { /* right below */
1569ce00eea3SSatish Balay       x_t = lx[n2 % m];
157047c6ae99SBarry Smith       y_t = ly[(n2/m)];
157147c6ae99SBarry Smith       s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t;
157247c6ae99SBarry Smith       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
157347c6ae99SBarry Smith     }
157447c6ae99SBarry Smith   }
157547c6ae99SBarry Smith 
157647c6ae99SBarry Smith   for (i=0; i<y; i++) {
157747c6ae99SBarry Smith     if (n3 >= 0) { /* directly left */
1578ce00eea3SSatish Balay       x_t = lx[n3 % m];
157947c6ae99SBarry Smith       /* y_t = y; */
158047c6ae99SBarry Smith       s_t = bases[n3] + (i+1)*x_t - s_x;
158147c6ae99SBarry Smith       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
158247c6ae99SBarry Smith     }
158347c6ae99SBarry Smith 
158447c6ae99SBarry Smith     for (j=0; j<x; j++) { idx[nn++] = xbase++; } /* interior */
158547c6ae99SBarry Smith 
158647c6ae99SBarry Smith     if (n5 >= 0) { /* directly right */
1587ce00eea3SSatish Balay       x_t = lx[n5 % m];
158847c6ae99SBarry Smith       /* y_t = y; */
158947c6ae99SBarry Smith       s_t = bases[n5] + (i)*x_t;
159047c6ae99SBarry Smith       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
159147c6ae99SBarry Smith     }
159247c6ae99SBarry Smith   }
159347c6ae99SBarry Smith 
159447c6ae99SBarry Smith   for (i=1; i<=s_y; i++) {
159547c6ae99SBarry Smith     if (n6 >= 0) { /* left above */
1596ce00eea3SSatish Balay       x_t = lx[n6 % m];
159747c6ae99SBarry Smith       /* y_t = ly[(n6/m)]; */
159847c6ae99SBarry Smith       s_t = bases[n6] + (i)*x_t - s_x;
159947c6ae99SBarry Smith       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
160047c6ae99SBarry Smith     }
160147c6ae99SBarry Smith     if (n7 >= 0) { /* directly above */
160247c6ae99SBarry Smith       x_t = x;
160347c6ae99SBarry Smith       /* y_t = ly[(n7/m)]; */
160447c6ae99SBarry Smith       s_t = bases[n7] + (i-1)*x_t;
160547c6ae99SBarry Smith       for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
160647c6ae99SBarry Smith     }
160747c6ae99SBarry Smith     if (n8 >= 0) { /* right above */
1608ce00eea3SSatish Balay       x_t = lx[n8 % m];
160947c6ae99SBarry Smith       /* y_t = ly[(n8/m)]; */
161047c6ae99SBarry Smith       s_t = bases[n8] + (i-1)*x_t;
161147c6ae99SBarry Smith       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
161247c6ae99SBarry Smith     }
161347c6ae99SBarry Smith   }
161447c6ae99SBarry Smith 
1615ce00eea3SSatish Balay   ierr = ISCreateBlock(comm,dof,nn,idx,PETSC_COPY_VALUES,&from);CHKERRQ(ierr);
161647c6ae99SBarry Smith   ierr = VecScatterCreate(global,from,local,to,&gtol);CHKERRQ(ierr);
161747c6ae99SBarry Smith   ierr = PetscLogObjectParent(da,gtol);CHKERRQ(ierr);
1618*fcfd50ebSBarry Smith   ierr = ISDestroy(&to);CHKERRQ(ierr);
1619*fcfd50ebSBarry Smith   ierr = ISDestroy(&from);CHKERRQ(ierr);
162047c6ae99SBarry Smith 
1621aa219208SBarry Smith   if (stencil_type == DMDA_STENCIL_STAR) {
1622ce00eea3SSatish Balay     n0 = sn0; n2 = sn2; n6 = sn6; n8 = sn8;
1623ce00eea3SSatish Balay   }
1624ce00eea3SSatish Balay 
1625ce00eea3SSatish Balay   if ((stencil_type == DMDA_STENCIL_STAR) ||
16261321219cSEthan Coon       (bx && bx != DMDA_BOUNDARY_PERIODIC) ||
16271321219cSEthan Coon       (by && by != DMDA_BOUNDARY_PERIODIC)) {
162847c6ae99SBarry Smith     /*
162947c6ae99SBarry Smith         Recompute the local to global mappings, this time keeping the
1630ce00eea3SSatish Balay       information about the cross corner processor numbers and any ghosted
1631ce00eea3SSatish Balay       but not periodic indices.
163247c6ae99SBarry Smith     */
163347c6ae99SBarry Smith     nn = 0;
163447c6ae99SBarry Smith     xbase = bases[rank];
163547c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
163647c6ae99SBarry Smith       if (n0 >= 0) { /* left below */
1637ce00eea3SSatish Balay         x_t = lx[n0 % m];
163847c6ae99SBarry Smith         y_t = ly[(n0/m)];
163947c6ae99SBarry Smith         s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x;
164047c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1641ce00eea3SSatish Balay       } else if (xs-Xs > 0 && ys-Ys > 0) {
1642ce00eea3SSatish Balay         for (j=0; j<s_x; j++) { idx[nn++] = -1;}
164347c6ae99SBarry Smith       }
164447c6ae99SBarry Smith       if (n1 >= 0) { /* directly below */
164547c6ae99SBarry Smith         x_t = x;
164647c6ae99SBarry Smith         y_t = ly[(n1/m)];
164747c6ae99SBarry Smith         s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t;
164847c6ae99SBarry Smith         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1649ce00eea3SSatish Balay       } else if (ys-Ys > 0) {
1650ce00eea3SSatish Balay         for (j=0; j<x; j++) { idx[nn++] = -1;}
165147c6ae99SBarry Smith       }
165247c6ae99SBarry Smith       if (n2 >= 0) { /* right below */
1653ce00eea3SSatish Balay         x_t = lx[n2 % m];
165447c6ae99SBarry Smith         y_t = ly[(n2/m)];
165547c6ae99SBarry Smith         s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t;
165647c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1657ce00eea3SSatish Balay       } else if (Xe-xe> 0 && ys-Ys > 0) {
1658ce00eea3SSatish Balay         for (j=0; j<s_x; j++) { idx[nn++] = -1;}
165947c6ae99SBarry Smith       }
166047c6ae99SBarry Smith     }
166147c6ae99SBarry Smith 
166247c6ae99SBarry Smith     for (i=0; i<y; i++) {
166347c6ae99SBarry Smith       if (n3 >= 0) { /* directly left */
1664ce00eea3SSatish Balay         x_t = lx[n3 % m];
166547c6ae99SBarry Smith         /* y_t = y; */
166647c6ae99SBarry Smith         s_t = bases[n3] + (i+1)*x_t - s_x;
166747c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1668ce00eea3SSatish Balay       } else if (xs-Xs > 0) {
1669ce00eea3SSatish Balay         for (j=0; j<s_x; j++) { idx[nn++] = -1;}
167047c6ae99SBarry Smith       }
167147c6ae99SBarry Smith 
167247c6ae99SBarry Smith       for (j=0; j<x; j++) { idx[nn++] = xbase++; } /* interior */
167347c6ae99SBarry Smith 
167447c6ae99SBarry Smith       if (n5 >= 0) { /* directly right */
1675ce00eea3SSatish Balay         x_t = lx[n5 % m];
167647c6ae99SBarry Smith         /* y_t = y; */
167747c6ae99SBarry Smith         s_t = bases[n5] + (i)*x_t;
167847c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1679ce00eea3SSatish Balay       } else if (Xe-xe > 0) {
1680ce00eea3SSatish Balay         for (j=0; j<s_x; j++) { idx[nn++] = -1;}
168147c6ae99SBarry Smith       }
168247c6ae99SBarry Smith     }
168347c6ae99SBarry Smith 
168447c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
168547c6ae99SBarry Smith       if (n6 >= 0) { /* left above */
1686ce00eea3SSatish Balay         x_t = lx[n6 % m];
168747c6ae99SBarry Smith         /* y_t = ly[(n6/m)]; */
168847c6ae99SBarry Smith         s_t = bases[n6] + (i)*x_t - s_x;
168947c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1690ce00eea3SSatish Balay       } else if (xs-Xs > 0 && Ye-ye > 0) {
1691ce00eea3SSatish Balay         for (j=0; j<s_x; j++) { idx[nn++] = -1;}
169247c6ae99SBarry Smith       }
169347c6ae99SBarry Smith       if (n7 >= 0) { /* directly above */
169447c6ae99SBarry Smith         x_t = x;
169547c6ae99SBarry Smith         /* y_t = ly[(n7/m)]; */
169647c6ae99SBarry Smith         s_t = bases[n7] + (i-1)*x_t;
169747c6ae99SBarry Smith         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1698ce00eea3SSatish Balay       } else if (Ye-ye > 0) {
1699ce00eea3SSatish Balay         for (j=0; j<x; j++) { idx[nn++] = -1;}
170047c6ae99SBarry Smith       }
170147c6ae99SBarry Smith       if (n8 >= 0) { /* right above */
1702ce00eea3SSatish Balay         x_t = lx[n8 % m];
170347c6ae99SBarry Smith         /* y_t = ly[(n8/m)]; */
170447c6ae99SBarry Smith         s_t = bases[n8] + (i-1)*x_t;
170547c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1706ce00eea3SSatish Balay       } else if (Xe-xe > 0 && Ye-ye > 0) {
1707ce00eea3SSatish Balay         for (j=0; j<s_x; j++) { idx[nn++] = -1;}
170847c6ae99SBarry Smith       }
170947c6ae99SBarry Smith     }
171047c6ae99SBarry Smith   }
1711ce00eea3SSatish Balay   /*
1712ce00eea3SSatish Balay      Set the local to global ordering in the global vector, this allows use
1713ce00eea3SSatish Balay      of VecSetValuesLocal().
1714ce00eea3SSatish Balay   */
1715ce00eea3SSatish Balay   ierr = ISCreateBlock(comm,dof,nn,idx,PETSC_OWN_POINTER,&ltogis);CHKERRQ(ierr);
1716ce00eea3SSatish Balay   ierr = PetscMalloc(nn*dof*sizeof(PetscInt),&idx_cpy);CHKERRQ(ierr);
1717db87c5ecSEthan Coon   ierr = PetscLogObjectMemory(da,nn*dof*sizeof(PetscInt));CHKERRQ(ierr);
1718ce00eea3SSatish Balay   ierr = ISGetIndices(ltogis, &idx_full);
1719ce00eea3SSatish Balay   ierr = PetscMemcpy(idx_cpy,idx_full,nn*dof*sizeof(PetscInt));CHKERRQ(ierr);
1720ce00eea3SSatish Balay   ierr = ISRestoreIndices(ltogis, &idx_full);
1721ce00eea3SSatish Balay   ierr = ISLocalToGlobalMappingCreateIS(ltogis,&da->ltogmap);CHKERRQ(ierr);
1722ce00eea3SSatish Balay   ierr = PetscLogObjectParent(da,da->ltogmap);CHKERRQ(ierr);
1723*fcfd50ebSBarry Smith   ierr = ISDestroy(&ltogis);CHKERRQ(ierr);
1724ce00eea3SSatish Balay   ierr = ISLocalToGlobalMappingBlock(da->ltogmap,dd->w,&da->ltogmapb);CHKERRQ(ierr);
1725ce00eea3SSatish Balay   ierr = PetscLogObjectParent(da,da->ltogmap);CHKERRQ(ierr);
172647c6ae99SBarry Smith 
1727ce00eea3SSatish Balay   ierr = PetscFree2(bases,ldims);CHKERRQ(ierr);
172847c6ae99SBarry Smith   dd->m  = m;  dd->n  = n;
1729ce00eea3SSatish Balay   /* note petsc expects xs/xe/Xs/Xe to be multiplied by #dofs in many places */
1730ce00eea3SSatish Balay   dd->xs = xs*dof; dd->xe = xe*dof; dd->ys = ys; dd->ye = ye; dd->zs = 0; dd->ze = 1;
1731ce00eea3SSatish Balay   dd->Xs = Xs*dof; dd->Xe = Xe*dof; dd->Ys = Ys; dd->Ye = Ye; dd->Zs = 0; dd->Ze = 1;
173247c6ae99SBarry Smith 
1733*fcfd50ebSBarry Smith   ierr = VecDestroy(&local);CHKERRQ(ierr);
1734*fcfd50ebSBarry Smith   ierr = VecDestroy(&global);CHKERRQ(ierr);
173547c6ae99SBarry Smith 
173647c6ae99SBarry Smith   dd->gtol      = gtol;
173747c6ae99SBarry Smith   dd->ltog      = ltog;
1738ce00eea3SSatish Balay   dd->idx       = idx_cpy;
1739ce00eea3SSatish Balay   dd->Nl        = nn*dof;
174047c6ae99SBarry Smith   dd->base      = base;
17419a42bb27SBarry Smith   da->ops->view = DMView_DA_2d;
174247c6ae99SBarry Smith   dd->ltol = PETSC_NULL;
174347c6ae99SBarry Smith   dd->ao   = PETSC_NULL;
174447c6ae99SBarry Smith 
174547c6ae99SBarry Smith   PetscFunctionReturn(0);
174647c6ae99SBarry Smith }
174747c6ae99SBarry Smith 
174847c6ae99SBarry Smith #undef __FUNCT__
1749aa219208SBarry Smith #define __FUNCT__ "DMDACreate2d"
175047c6ae99SBarry Smith /*@C
1751aa219208SBarry Smith    DMDACreate2d -  Creates an object that will manage the communication of  two-dimensional
175247c6ae99SBarry Smith    regular array data that is distributed across some processors.
175347c6ae99SBarry Smith 
175447c6ae99SBarry Smith    Collective on MPI_Comm
175547c6ae99SBarry Smith 
175647c6ae99SBarry Smith    Input Parameters:
175747c6ae99SBarry Smith +  comm - MPI communicator
17581321219cSEthan Coon .  bx,by - type of ghost nodes the array have.
17591321219cSEthan Coon          Use one of DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_GHOSTED, DMDA_BOUNDARY_PERIODIC.
1760aa219208SBarry Smith .  stencil_type - stencil type.  Use either DMDA_STENCIL_BOX or DMDA_STENCIL_STAR.
176147c6ae99SBarry 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
176247c6ae99SBarry Smith             from the command line with -da_grid_x <M> -da_grid_y <N>)
176347c6ae99SBarry Smith .  m,n - corresponding number of processors in each dimension
176447c6ae99SBarry Smith          (or PETSC_DECIDE to have calculated)
176547c6ae99SBarry Smith .  dof - number of degrees of freedom per node
176647c6ae99SBarry Smith .  s - stencil width
176747c6ae99SBarry Smith -  lx, ly - arrays containing the number of nodes in each cell along
176847c6ae99SBarry Smith            the x and y coordinates, or PETSC_NULL. If non-null, these
176947c6ae99SBarry Smith            must be of length as m and n, and the corresponding
177047c6ae99SBarry Smith            m and n cannot be PETSC_DECIDE. The sum of the lx[] entries
177147c6ae99SBarry Smith            must be M, and the sum of the ly[] entries must be N.
177247c6ae99SBarry Smith 
177347c6ae99SBarry Smith    Output Parameter:
177447c6ae99SBarry Smith .  da - the resulting distributed array object
177547c6ae99SBarry Smith 
177647c6ae99SBarry Smith    Options Database Key:
1777aa219208SBarry Smith +  -da_view - Calls DMView() at the conclusion of DMDACreate2d()
177847c6ae99SBarry Smith .  -da_grid_x <nx> - number of grid points in x direction, if M < 0
177947c6ae99SBarry Smith .  -da_grid_y <ny> - number of grid points in y direction, if N < 0
178047c6ae99SBarry Smith .  -da_processors_x <nx> - number of processors in x direction
178147c6ae99SBarry Smith .  -da_processors_y <ny> - number of processors in y direction
178247c6ae99SBarry Smith .  -da_refine_x - refinement ratio in x direction
178347c6ae99SBarry Smith -  -da_refine_y - refinement ratio in y direction
178447c6ae99SBarry Smith 
178547c6ae99SBarry Smith    Level: beginner
178647c6ae99SBarry Smith 
178747c6ae99SBarry Smith    Notes:
1788aa219208SBarry Smith    The stencil type DMDA_STENCIL_STAR with width 1 corresponds to the
1789aa219208SBarry Smith    standard 5-pt stencil, while DMDA_STENCIL_BOX with width 1 denotes
179047c6ae99SBarry Smith    the standard 9-pt stencil.
179147c6ae99SBarry Smith 
1792aa219208SBarry Smith    The array data itself is NOT stored in the DMDA, it is stored in Vec objects;
1793564755cdSBarry Smith    The appropriate vector objects can be obtained with calls to DMCreateGlobalVector()
1794564755cdSBarry Smith    and DMCreateLocalVector() and calls to VecDuplicate() if more are needed.
179547c6ae99SBarry Smith 
179647c6ae99SBarry Smith .keywords: distributed array, create, two-dimensional
179747c6ae99SBarry Smith 
1798aa219208SBarry Smith .seealso: DMDestroy(), DMView(), DMDACreate1d(), DMDACreate3d(), DMGlobalToLocalBegin(), DMDAGetRefinementFactor(),
1799aa219208SBarry Smith           DMGlobalToLocalEnd(), DMLocalToGlobalBegin(), DMDALocalToLocalBegin(), DMDALocalToLocalEnd(), DMDASetRefinementFactor(),
1800aa219208SBarry Smith           DMDAGetInfo(), DMCreateGlobalVector(), DMCreateLocalVector(), DMDACreateNaturalVector(), DMDALoad(), DMDAGetOwnershipRanges()
180147c6ae99SBarry Smith 
180247c6ae99SBarry Smith @*/
18031321219cSEthan Coon PetscErrorCode  DMDACreate2d(MPI_Comm comm,DMDABoundaryType bx,DMDABoundaryType by,DMDAStencilType stencil_type,
18049a42bb27SBarry Smith                           PetscInt M,PetscInt N,PetscInt m,PetscInt n,PetscInt dof,PetscInt s,const PetscInt lx[],const PetscInt ly[],DM *da)
180547c6ae99SBarry Smith {
180647c6ae99SBarry Smith   PetscErrorCode ierr;
180747c6ae99SBarry Smith 
180847c6ae99SBarry Smith   PetscFunctionBegin;
1809aa219208SBarry Smith   ierr = DMDACreate(comm, da);CHKERRQ(ierr);
1810aa219208SBarry Smith   ierr = DMDASetDim(*da, 2);CHKERRQ(ierr);
1811aa219208SBarry Smith   ierr = DMDASetSizes(*da, M, N, 1);CHKERRQ(ierr);
1812aa219208SBarry Smith   ierr = DMDASetNumProcs(*da, m, n, PETSC_DECIDE);CHKERRQ(ierr);
1813755f258dSLisandro Dalcin   ierr = DMDASetBoundaryType(*da, bx, by, DMDA_BOUNDARY_NONE);CHKERRQ(ierr);
1814aa219208SBarry Smith   ierr = DMDASetDof(*da, dof);CHKERRQ(ierr);
1815aa219208SBarry Smith   ierr = DMDASetStencilType(*da, stencil_type);CHKERRQ(ierr);
1816aa219208SBarry Smith   ierr = DMDASetStencilWidth(*da, s);CHKERRQ(ierr);
1817aa219208SBarry Smith   ierr = DMDASetOwnershipRanges(*da, lx, ly, PETSC_NULL);CHKERRQ(ierr);
181847c6ae99SBarry Smith   /* This violates the behavior for other classes, but right now users expect negative dimensions to be handled this way */
18199a42bb27SBarry Smith   ierr = DMSetFromOptions(*da);CHKERRQ(ierr);
18209a42bb27SBarry Smith   ierr = DMSetUp(*da);CHKERRQ(ierr);
18217242296bSJed Brown   ierr = DMView_DA_Private(*da);CHKERRQ(ierr);
182247c6ae99SBarry Smith   PetscFunctionReturn(0);
182347c6ae99SBarry Smith }
1824