xref: /petsc/src/dm/impls/da/da2.c (revision 1411c6ee19d0933bd2bf6d23611b5e480183c563)
19a42bb27SBarry Smith 
247c6ae99SBarry Smith #define PETSCDM_DLL
347c6ae99SBarry Smith 
4e1589f56SBarry Smith #include "private/daimpl.h"    /*I   "petscdm.h"   I*/
547c6ae99SBarry Smith 
647c6ae99SBarry Smith #undef __FUNCT__
79a42bb27SBarry Smith #define __FUNCT__ "DMView_DA_2d"
89a42bb27SBarry Smith PetscErrorCode DMView_DA_2d(DM da,PetscViewer viewer)
947c6ae99SBarry Smith {
1047c6ae99SBarry Smith   PetscErrorCode ierr;
1147c6ae99SBarry Smith   PetscMPIInt    rank;
129a42bb27SBarry Smith   PetscBool      iascii,isdraw,isbinary;
1347c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
149a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
159a42bb27SBarry Smith   PetscBool      ismatlab;
169a42bb27SBarry Smith #endif
1747c6ae99SBarry Smith 
1847c6ae99SBarry Smith   PetscFunctionBegin;
1947c6ae99SBarry Smith   ierr = MPI_Comm_rank(((PetscObject)da)->comm,&rank);CHKERRQ(ierr);
2047c6ae99SBarry Smith 
2147c6ae99SBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
2247c6ae99SBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
239a42bb27SBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
249a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
259a42bb27SBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERMATLAB,&ismatlab);CHKERRQ(ierr);
269a42bb27SBarry Smith #endif
2747c6ae99SBarry Smith   if (iascii) {
2847c6ae99SBarry Smith     PetscViewerFormat format;
2947c6ae99SBarry Smith 
3047c6ae99SBarry Smith     ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr);
3147c6ae99SBarry Smith     if (format != PETSC_VIEWER_ASCII_VTK && format != PETSC_VIEWER_ASCII_VTK_CELL) {
32aa219208SBarry Smith       DMDALocalInfo info;
33aa219208SBarry Smith       ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
3447c6ae99SBarry 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);
3547c6ae99SBarry 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);
3647c6ae99SBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
373da9ae13SJed Brown     } else {
383da9ae13SJed Brown       ierr = DMView_DA_VTK(da,viewer);CHKERRQ(ierr);
3947c6ae99SBarry Smith     }
4047c6ae99SBarry Smith   } else if (isdraw) {
4147c6ae99SBarry Smith     PetscDraw  draw;
4247c6ae99SBarry Smith     double     ymin = -1*dd->s-1,ymax = dd->N+dd->s;
4347c6ae99SBarry Smith     double     xmin = -1*dd->s-1,xmax = dd->M+dd->s;
4447c6ae99SBarry Smith     double     x,y;
4547c6ae99SBarry Smith     PetscInt   base,*idx;
4647c6ae99SBarry Smith     char       node[10];
4747c6ae99SBarry Smith     PetscBool  isnull;
4847c6ae99SBarry Smith 
4947c6ae99SBarry Smith     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
5047c6ae99SBarry Smith     ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
5147c6ae99SBarry Smith     if (!dd->coordinates) {
5247c6ae99SBarry Smith       ierr = PetscDrawSetCoordinates(draw,xmin,ymin,xmax,ymax);CHKERRQ(ierr);
5347c6ae99SBarry Smith     }
5447c6ae99SBarry Smith     ierr = PetscDrawSynchronizedClear(draw);CHKERRQ(ierr);
5547c6ae99SBarry Smith 
5647c6ae99SBarry Smith     /* first processor draw all node lines */
5747c6ae99SBarry Smith     if (!rank) {
5847c6ae99SBarry Smith       ymin = 0.0; ymax = dd->N - 1;
5947c6ae99SBarry Smith       for (xmin=0; xmin<dd->M; xmin++) {
6047c6ae99SBarry Smith         ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_BLACK);CHKERRQ(ierr);
6147c6ae99SBarry Smith       }
6247c6ae99SBarry Smith       xmin = 0.0; xmax = dd->M - 1;
6347c6ae99SBarry Smith       for (ymin=0; ymin<dd->N; ymin++) {
6447c6ae99SBarry Smith         ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_BLACK);CHKERRQ(ierr);
6547c6ae99SBarry Smith       }
6647c6ae99SBarry Smith     }
6747c6ae99SBarry Smith     ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr);
6847c6ae99SBarry Smith     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
6947c6ae99SBarry Smith 
7047c6ae99SBarry Smith     /* draw my box */
7147c6ae99SBarry Smith     ymin = dd->ys; ymax = dd->ye - 1; xmin = dd->xs/dd->w;
7247c6ae99SBarry Smith     xmax =(dd->xe-1)/dd->w;
7347c6ae99SBarry Smith     ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_RED);CHKERRQ(ierr);
7447c6ae99SBarry Smith     ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_RED);CHKERRQ(ierr);
7547c6ae99SBarry Smith     ierr = PetscDrawLine(draw,xmin,ymax,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr);
7647c6ae99SBarry Smith     ierr = PetscDrawLine(draw,xmax,ymin,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr);
7747c6ae99SBarry Smith 
7847c6ae99SBarry Smith     /* put in numbers */
7947c6ae99SBarry Smith     base = (dd->base)/dd->w;
8047c6ae99SBarry Smith     for (y=ymin; y<=ymax; y++) {
8147c6ae99SBarry Smith       for (x=xmin; x<=xmax; x++) {
8247c6ae99SBarry Smith         sprintf(node,"%d",(int)base++);
8347c6ae99SBarry Smith         ierr = PetscDrawString(draw,x,y,PETSC_DRAW_BLACK,node);CHKERRQ(ierr);
8447c6ae99SBarry Smith       }
8547c6ae99SBarry Smith     }
8647c6ae99SBarry Smith 
8747c6ae99SBarry Smith     ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr);
8847c6ae99SBarry Smith     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
8947c6ae99SBarry Smith     /* overlay ghost numbers, useful for error checking */
9047c6ae99SBarry Smith     /* put in numbers */
9147c6ae99SBarry Smith 
9247c6ae99SBarry Smith     base = 0; idx = dd->idx;
9347c6ae99SBarry Smith     ymin = dd->Ys; ymax = dd->Ye; xmin = dd->Xs; xmax = dd->Xe;
9447c6ae99SBarry Smith     for (y=ymin; y<ymax; y++) {
9547c6ae99SBarry Smith       for (x=xmin; x<xmax; x++) {
9647c6ae99SBarry Smith         if ((base % dd->w) == 0) {
9747c6ae99SBarry Smith           sprintf(node,"%d",(int)(idx[base]/dd->w));
9847c6ae99SBarry Smith           ierr = PetscDrawString(draw,x/dd->w,y,PETSC_DRAW_BLUE,node);CHKERRQ(ierr);
9947c6ae99SBarry Smith         }
10047c6ae99SBarry Smith         base++;
10147c6ae99SBarry Smith       }
10247c6ae99SBarry Smith     }
10347c6ae99SBarry Smith     ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr);
10447c6ae99SBarry Smith     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
1059a42bb27SBarry Smith   } else if (isbinary){
1069a42bb27SBarry Smith     ierr = DMView_DA_Binary(da,viewer);CHKERRQ(ierr);
1079a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
1089a42bb27SBarry Smith   } else if (ismatlab) {
1099a42bb27SBarry Smith     ierr = DMView_DA_Matlab(da,viewer);CHKERRQ(ierr);
1109a42bb27SBarry Smith #endif
111aa219208SBarry Smith   } else SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for DMDA 1d",((PetscObject)viewer)->type_name);
11247c6ae99SBarry Smith   PetscFunctionReturn(0);
11347c6ae99SBarry Smith }
11447c6ae99SBarry Smith 
11547c6ae99SBarry Smith /*
11647c6ae99SBarry Smith       M is number of grid points
11747c6ae99SBarry Smith       m is number of processors
11847c6ae99SBarry Smith 
11947c6ae99SBarry Smith */
12047c6ae99SBarry Smith #undef __FUNCT__
121aa219208SBarry Smith #define __FUNCT__ "DMDASplitComm2d"
122aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDASplitComm2d(MPI_Comm comm,PetscInt M,PetscInt N,PetscInt sw,MPI_Comm *outcomm)
12347c6ae99SBarry Smith {
12447c6ae99SBarry Smith   PetscErrorCode ierr;
12547c6ae99SBarry Smith   PetscInt       m,n = 0,x = 0,y = 0;
12647c6ae99SBarry Smith   PetscMPIInt    size,csize,rank;
12747c6ae99SBarry Smith 
12847c6ae99SBarry Smith   PetscFunctionBegin;
12947c6ae99SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
13047c6ae99SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
13147c6ae99SBarry Smith 
13247c6ae99SBarry Smith   csize = 4*size;
13347c6ae99SBarry Smith   do {
13447c6ae99SBarry Smith     if (csize % 4) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Cannot split communicator of size %d tried %d %D %D",size,csize,x,y);
13547c6ae99SBarry Smith     csize   = csize/4;
13647c6ae99SBarry Smith 
13747c6ae99SBarry Smith     m = (PetscInt)(0.5 + sqrt(((double)M)*((double)csize)/((double)N)));
13847c6ae99SBarry Smith     if (!m) m = 1;
13947c6ae99SBarry Smith     while (m > 0) {
14047c6ae99SBarry Smith       n = csize/m;
14147c6ae99SBarry Smith       if (m*n == csize) break;
14247c6ae99SBarry Smith       m--;
14347c6ae99SBarry Smith     }
14447c6ae99SBarry Smith     if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;}
14547c6ae99SBarry Smith 
14647c6ae99SBarry Smith     x = M/m + ((M % m) > ((csize-1) % m));
14747c6ae99SBarry Smith     y = (N + (csize-1)/m)/n;
14847c6ae99SBarry Smith   } while ((x < 4 || y < 4) && csize > 1);
14947c6ae99SBarry Smith   if (size != csize) {
15047c6ae99SBarry Smith     MPI_Group    entire_group,sub_group;
15147c6ae99SBarry Smith     PetscMPIInt  i,*groupies;
15247c6ae99SBarry Smith 
15347c6ae99SBarry Smith     ierr = MPI_Comm_group(comm,&entire_group);CHKERRQ(ierr);
15447c6ae99SBarry Smith     ierr = PetscMalloc(csize*sizeof(PetscInt),&groupies);CHKERRQ(ierr);
15547c6ae99SBarry Smith     for (i=0; i<csize; i++) {
15647c6ae99SBarry Smith       groupies[i] = (rank/csize)*csize + i;
15747c6ae99SBarry Smith     }
15847c6ae99SBarry Smith     ierr = MPI_Group_incl(entire_group,csize,groupies,&sub_group);CHKERRQ(ierr);
15947c6ae99SBarry Smith     ierr = PetscFree(groupies);CHKERRQ(ierr);
16047c6ae99SBarry Smith     ierr = MPI_Comm_create(comm,sub_group,outcomm);CHKERRQ(ierr);
16147c6ae99SBarry Smith     ierr = MPI_Group_free(&entire_group);CHKERRQ(ierr);
16247c6ae99SBarry Smith     ierr = MPI_Group_free(&sub_group);CHKERRQ(ierr);
163aa219208SBarry Smith     ierr = PetscInfo1(0,"DMDASplitComm2d:Creating redundant coarse problems of size %d\n",csize);CHKERRQ(ierr);
16447c6ae99SBarry Smith   } else {
16547c6ae99SBarry Smith     *outcomm = comm;
16647c6ae99SBarry Smith   }
16747c6ae99SBarry Smith   PetscFunctionReturn(0);
16847c6ae99SBarry Smith }
16947c6ae99SBarry Smith 
17047c6ae99SBarry Smith #undef __FUNCT__
171aa219208SBarry Smith #define __FUNCT__ "DMDASetLocalFunction"
17247c6ae99SBarry Smith /*@C
173aa219208SBarry Smith        DMDASetLocalFunction - Caches in a DM a local function.
17447c6ae99SBarry Smith 
175aa219208SBarry Smith    Logically Collective on DMDA
17647c6ae99SBarry Smith 
17747c6ae99SBarry Smith    Input Parameter:
17847c6ae99SBarry Smith +  da - initial distributed array
17947c6ae99SBarry Smith -  lf - the local function
18047c6ae99SBarry Smith 
18147c6ae99SBarry Smith    Level: intermediate
18247c6ae99SBarry Smith 
18347c6ae99SBarry Smith    Notes: The routine SNESDAFormFunction() uses this the cached function to evaluate the user provided function.
18447c6ae99SBarry Smith 
18547c6ae99SBarry Smith .keywords:  distributed array, refine
18647c6ae99SBarry Smith 
187aa219208SBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalFunctioni()
18847c6ae99SBarry Smith @*/
189aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDASetLocalFunction(DM da,DMDALocalFunction1 lf)
19047c6ae99SBarry Smith {
19147c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
19247c6ae99SBarry Smith   PetscFunctionBegin;
19347c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
19447c6ae99SBarry Smith   dd->lf    = lf;
19547c6ae99SBarry Smith   PetscFunctionReturn(0);
19647c6ae99SBarry Smith }
19747c6ae99SBarry Smith 
19847c6ae99SBarry Smith #undef __FUNCT__
199aa219208SBarry Smith #define __FUNCT__ "DMDASetLocalFunctioni"
20047c6ae99SBarry Smith /*@C
201aa219208SBarry Smith        DMDASetLocalFunctioni - Caches in a DM a local function that evaluates a single component
20247c6ae99SBarry Smith 
203aa219208SBarry Smith    Logically Collective on DMDA
20447c6ae99SBarry Smith 
20547c6ae99SBarry Smith    Input Parameter:
20647c6ae99SBarry Smith +  da - initial distributed array
20747c6ae99SBarry Smith -  lfi - the local function
20847c6ae99SBarry Smith 
20947c6ae99SBarry Smith    Level: intermediate
21047c6ae99SBarry Smith 
21147c6ae99SBarry Smith .keywords:  distributed array, refine
21247c6ae99SBarry Smith 
213aa219208SBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalFunction()
21447c6ae99SBarry Smith @*/
215aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDASetLocalFunctioni(DM da,PetscErrorCode (*lfi)(DMDALocalInfo*,MatStencil*,void*,PetscScalar*,void*))
21647c6ae99SBarry Smith {
21747c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
21847c6ae99SBarry Smith   PetscFunctionBegin;
21947c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
22047c6ae99SBarry Smith   dd->lfi = lfi;
22147c6ae99SBarry Smith   PetscFunctionReturn(0);
22247c6ae99SBarry Smith }
22347c6ae99SBarry Smith 
22447c6ae99SBarry Smith #undef __FUNCT__
225aa219208SBarry Smith #define __FUNCT__ "DMDASetLocalFunctionib"
22647c6ae99SBarry Smith /*@C
227aa219208SBarry Smith        DMDASetLocalFunctionib - Caches in a DM a block local function that evaluates a single component
22847c6ae99SBarry Smith 
229aa219208SBarry Smith    Logically Collective on DMDA
23047c6ae99SBarry Smith 
23147c6ae99SBarry Smith    Input Parameter:
23247c6ae99SBarry Smith +  da - initial distributed array
23347c6ae99SBarry Smith -  lfi - the local function
23447c6ae99SBarry Smith 
23547c6ae99SBarry Smith    Level: intermediate
23647c6ae99SBarry Smith 
23747c6ae99SBarry Smith .keywords:  distributed array, refine
23847c6ae99SBarry Smith 
239aa219208SBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalFunction()
24047c6ae99SBarry Smith @*/
241aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDASetLocalFunctionib(DM da,PetscErrorCode (*lfi)(DMDALocalInfo*,MatStencil*,void*,PetscScalar*,void*))
24247c6ae99SBarry Smith {
24347c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
24447c6ae99SBarry Smith   PetscFunctionBegin;
24547c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
24647c6ae99SBarry Smith   dd->lfib = lfi;
24747c6ae99SBarry Smith   PetscFunctionReturn(0);
24847c6ae99SBarry Smith }
24947c6ae99SBarry Smith 
25047c6ae99SBarry Smith #undef __FUNCT__
251aa219208SBarry Smith #define __FUNCT__ "DMDASetLocalAdicFunction_Private"
252aa219208SBarry Smith PetscErrorCode DMDASetLocalAdicFunction_Private(DM da,DMDALocalFunction1 ad_lf)
25347c6ae99SBarry Smith {
25447c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
25547c6ae99SBarry Smith   PetscFunctionBegin;
25647c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
25747c6ae99SBarry Smith   dd->adic_lf = ad_lf;
25847c6ae99SBarry Smith   PetscFunctionReturn(0);
25947c6ae99SBarry Smith }
26047c6ae99SBarry Smith 
26147c6ae99SBarry Smith /*MC
262aa219208SBarry Smith        DMDASetLocalAdicFunctioni - Caches in a DM a local functioni computed by ADIC/ADIFOR
26347c6ae99SBarry Smith 
26447c6ae99SBarry Smith    Synopsis:
265aa219208SBarry Smith    PetscErrorCode DMDASetLocalAdicFunctioni(DM da,PetscInt (ad_lf*)(DMDALocalInfo*,MatStencil*,void*,void*,void*)
26647c6ae99SBarry Smith 
267aa219208SBarry Smith    Logically Collective on DMDA
26847c6ae99SBarry Smith 
26947c6ae99SBarry Smith    Input Parameter:
27047c6ae99SBarry Smith +  da - initial distributed array
27147c6ae99SBarry Smith -  ad_lfi - the local function as computed by ADIC/ADIFOR
27247c6ae99SBarry Smith 
27347c6ae99SBarry Smith    Level: intermediate
27447c6ae99SBarry Smith 
27547c6ae99SBarry Smith .keywords:  distributed array, refine
27647c6ae99SBarry Smith 
277aa219208SBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalFunction(),
278aa219208SBarry Smith           DMDASetLocalJacobian(), DMDASetLocalFunctioni()
27947c6ae99SBarry Smith M*/
28047c6ae99SBarry Smith 
28147c6ae99SBarry Smith #undef __FUNCT__
282aa219208SBarry Smith #define __FUNCT__ "DMDASetLocalAdicFunctioni_Private"
283aa219208SBarry Smith PetscErrorCode DMDASetLocalAdicFunctioni_Private(DM da,PetscErrorCode (*ad_lfi)(DMDALocalInfo*,MatStencil*,void*,void*,void*))
28447c6ae99SBarry Smith {
28547c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
28647c6ae99SBarry Smith   PetscFunctionBegin;
28747c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
28847c6ae99SBarry Smith   dd->adic_lfi = ad_lfi;
28947c6ae99SBarry Smith   PetscFunctionReturn(0);
29047c6ae99SBarry Smith }
29147c6ae99SBarry Smith 
29247c6ae99SBarry Smith /*MC
293aa219208SBarry Smith        DMDASetLocalAdicMFFunctioni - Caches in a DM a local functioni computed by ADIC/ADIFOR
29447c6ae99SBarry Smith 
29547c6ae99SBarry Smith    Synopsis:
296aa219208SBarry Smith    PetscErrorCode  DMDASetLocalAdicFunctioni(DM da,int (ad_lf*)(DMDALocalInfo*,MatStencil*,void*,void*,void*)
29747c6ae99SBarry Smith 
298aa219208SBarry Smith    Logically Collective on DMDA
29947c6ae99SBarry Smith 
30047c6ae99SBarry Smith    Input Parameter:
30147c6ae99SBarry Smith +  da - initial distributed array
30247c6ae99SBarry Smith -  admf_lfi - the local matrix-free function as computed by ADIC/ADIFOR
30347c6ae99SBarry Smith 
30447c6ae99SBarry Smith    Level: intermediate
30547c6ae99SBarry Smith 
30647c6ae99SBarry Smith .keywords:  distributed array, refine
30747c6ae99SBarry Smith 
308aa219208SBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalFunction(),
309aa219208SBarry Smith           DMDASetLocalJacobian(), DMDASetLocalFunctioni()
31047c6ae99SBarry Smith M*/
31147c6ae99SBarry Smith 
31247c6ae99SBarry Smith #undef __FUNCT__
313aa219208SBarry Smith #define __FUNCT__ "DMDASetLocalAdicMFFunctioni_Private"
314aa219208SBarry Smith PetscErrorCode DMDASetLocalAdicMFFunctioni_Private(DM da,PetscErrorCode (*admf_lfi)(DMDALocalInfo*,MatStencil*,void*,void*,void*))
31547c6ae99SBarry Smith {
31647c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
31747c6ae99SBarry Smith   PetscFunctionBegin;
31847c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
31947c6ae99SBarry Smith   dd->adicmf_lfi = admf_lfi;
32047c6ae99SBarry Smith   PetscFunctionReturn(0);
32147c6ae99SBarry Smith }
32247c6ae99SBarry Smith 
32347c6ae99SBarry Smith /*MC
324aa219208SBarry Smith        DMDASetLocalAdicFunctionib - Caches in a DM a block local functioni computed by ADIC/ADIFOR
32547c6ae99SBarry Smith 
32647c6ae99SBarry Smith    Synopsis:
327aa219208SBarry Smith    PetscErrorCode DMDASetLocalAdicFunctionib(DM da,PetscInt (ad_lf*)(DMDALocalInfo*,MatStencil*,void*,void*,void*)
32847c6ae99SBarry Smith 
329aa219208SBarry Smith    Logically Collective on DMDA
33047c6ae99SBarry Smith 
33147c6ae99SBarry Smith    Input Parameter:
33247c6ae99SBarry Smith +  da - initial distributed array
33347c6ae99SBarry Smith -  ad_lfi - the local function as computed by ADIC/ADIFOR
33447c6ae99SBarry Smith 
33547c6ae99SBarry Smith    Level: intermediate
33647c6ae99SBarry Smith 
33747c6ae99SBarry Smith .keywords:  distributed array, refine
33847c6ae99SBarry Smith 
339aa219208SBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalFunction(),
340aa219208SBarry Smith           DMDASetLocalJacobian(), DMDASetLocalFunctionib()
34147c6ae99SBarry Smith M*/
34247c6ae99SBarry Smith 
34347c6ae99SBarry Smith #undef __FUNCT__
344aa219208SBarry Smith #define __FUNCT__ "DMDASetLocalAdicFunctionib_Private"
345aa219208SBarry Smith PetscErrorCode DMDASetLocalAdicFunctionib_Private(DM da,PetscErrorCode (*ad_lfi)(DMDALocalInfo*,MatStencil*,void*,void*,void*))
34647c6ae99SBarry Smith {
34747c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
34847c6ae99SBarry Smith   PetscFunctionBegin;
34947c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
35047c6ae99SBarry Smith   dd->adic_lfib = ad_lfi;
35147c6ae99SBarry Smith   PetscFunctionReturn(0);
35247c6ae99SBarry Smith }
35347c6ae99SBarry Smith 
35447c6ae99SBarry Smith /*MC
355aa219208SBarry Smith        DMDASetLocalAdicMFFunctionib - Caches in a DM a block local functioni computed by ADIC/ADIFOR
35647c6ae99SBarry Smith 
35747c6ae99SBarry Smith    Synopsis:
358aa219208SBarry Smith    PetscErrorCode  DMDASetLocalAdicFunctionib(DM da,int (ad_lf*)(DMDALocalInfo*,MatStencil*,void*,void*,void*)
35947c6ae99SBarry Smith 
360aa219208SBarry Smith    Logically Collective on DMDA
36147c6ae99SBarry Smith 
36247c6ae99SBarry Smith    Input Parameter:
36347c6ae99SBarry Smith +  da - initial distributed array
36447c6ae99SBarry Smith -  admf_lfi - the local matrix-free function as computed by ADIC/ADIFOR
36547c6ae99SBarry Smith 
36647c6ae99SBarry Smith    Level: intermediate
36747c6ae99SBarry Smith 
36847c6ae99SBarry Smith .keywords:  distributed array, refine
36947c6ae99SBarry Smith 
370aa219208SBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalFunction(),
371aa219208SBarry Smith           DMDASetLocalJacobian(), DMDASetLocalFunctionib()
37247c6ae99SBarry Smith M*/
37347c6ae99SBarry Smith 
37447c6ae99SBarry Smith #undef __FUNCT__
375aa219208SBarry Smith #define __FUNCT__ "DMDASetLocalAdicMFFunctionib_Private"
376aa219208SBarry Smith PetscErrorCode DMDASetLocalAdicMFFunctionib_Private(DM da,PetscErrorCode (*admf_lfi)(DMDALocalInfo*,MatStencil*,void*,void*,void*))
37747c6ae99SBarry Smith {
37847c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
37947c6ae99SBarry Smith   PetscFunctionBegin;
38047c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
38147c6ae99SBarry Smith   dd->adicmf_lfib = admf_lfi;
38247c6ae99SBarry Smith   PetscFunctionReturn(0);
38347c6ae99SBarry Smith }
38447c6ae99SBarry Smith 
38547c6ae99SBarry Smith /*MC
386aa219208SBarry Smith        DMDASetLocalAdicMFFunction - Caches in a DM a local function computed by ADIC/ADIFOR
38747c6ae99SBarry Smith 
38847c6ae99SBarry Smith    Synopsis:
389aa219208SBarry Smith    PetscErrorCode DMDASetLocalAdicMFFunction(DM da,DMDALocalFunction1 ad_lf)
39047c6ae99SBarry Smith 
391aa219208SBarry Smith    Logically Collective on DMDA
39247c6ae99SBarry Smith 
39347c6ae99SBarry Smith    Input Parameter:
39447c6ae99SBarry Smith +  da - initial distributed array
39547c6ae99SBarry Smith -  ad_lf - the local function as computed by ADIC/ADIFOR
39647c6ae99SBarry Smith 
39747c6ae99SBarry Smith    Level: intermediate
39847c6ae99SBarry Smith 
39947c6ae99SBarry Smith .keywords:  distributed array, refine
40047c6ae99SBarry Smith 
401aa219208SBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalFunction(),
402aa219208SBarry Smith           DMDASetLocalJacobian()
40347c6ae99SBarry Smith M*/
40447c6ae99SBarry Smith 
40547c6ae99SBarry Smith #undef __FUNCT__
406aa219208SBarry Smith #define __FUNCT__ "DMDASetLocalAdicMFFunction_Private"
407aa219208SBarry Smith PetscErrorCode DMDASetLocalAdicMFFunction_Private(DM da,DMDALocalFunction1 ad_lf)
40847c6ae99SBarry Smith {
40947c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
41047c6ae99SBarry Smith   PetscFunctionBegin;
41147c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
41247c6ae99SBarry Smith   dd->adicmf_lf = ad_lf;
41347c6ae99SBarry Smith   PetscFunctionReturn(0);
41447c6ae99SBarry Smith }
41547c6ae99SBarry Smith 
41647c6ae99SBarry Smith /*@C
417aa219208SBarry Smith        DMDASetLocalJacobian - Caches in a DM a local Jacobian computation function
41847c6ae99SBarry Smith 
419aa219208SBarry Smith    Logically Collective on DMDA
42047c6ae99SBarry Smith 
42147c6ae99SBarry Smith 
42247c6ae99SBarry Smith    Input Parameter:
42347c6ae99SBarry Smith +  da - initial distributed array
42447c6ae99SBarry Smith -  lj - the local Jacobian
42547c6ae99SBarry Smith 
42647c6ae99SBarry Smith    Level: intermediate
42747c6ae99SBarry Smith 
42847c6ae99SBarry Smith    Notes: The routine SNESDAFormFunction() uses this the cached function to evaluate the user provided function.
42947c6ae99SBarry Smith 
43047c6ae99SBarry Smith .keywords:  distributed array, refine
43147c6ae99SBarry Smith 
432aa219208SBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalFunction()
43347c6ae99SBarry Smith @*/
43447c6ae99SBarry Smith #undef __FUNCT__
435aa219208SBarry Smith #define __FUNCT__ "DMDASetLocalJacobian"
436aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDASetLocalJacobian(DM da,DMDALocalFunction1 lj)
43747c6ae99SBarry Smith {
43847c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
43947c6ae99SBarry Smith   PetscFunctionBegin;
44047c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
44147c6ae99SBarry Smith   dd->lj    = lj;
44247c6ae99SBarry Smith   PetscFunctionReturn(0);
44347c6ae99SBarry Smith }
44447c6ae99SBarry Smith 
44547c6ae99SBarry Smith #undef __FUNCT__
446aa219208SBarry Smith #define __FUNCT__ "DMDAGetLocalFunction"
44747c6ae99SBarry Smith /*@C
448aa219208SBarry Smith        DMDAGetLocalFunction - Gets from a DM a local function and its ADIC/ADIFOR Jacobian
44947c6ae99SBarry Smith 
45047c6ae99SBarry Smith    Note Collective
45147c6ae99SBarry Smith 
45247c6ae99SBarry Smith    Input Parameter:
45347c6ae99SBarry Smith .  da - initial distributed array
45447c6ae99SBarry Smith 
45547c6ae99SBarry Smith    Output Parameter:
45647c6ae99SBarry Smith .  lf - the local function
45747c6ae99SBarry Smith 
45847c6ae99SBarry Smith    Level: intermediate
45947c6ae99SBarry Smith 
46047c6ae99SBarry Smith .keywords:  distributed array, refine
46147c6ae99SBarry Smith 
462aa219208SBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalJacobian(), DMDASetLocalFunction()
46347c6ae99SBarry Smith @*/
464aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDAGetLocalFunction(DM da,DMDALocalFunction1 *lf)
46547c6ae99SBarry Smith {
46647c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
46747c6ae99SBarry Smith   PetscFunctionBegin;
46847c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
46947c6ae99SBarry Smith   if (lf) *lf = dd->lf;
47047c6ae99SBarry Smith   PetscFunctionReturn(0);
47147c6ae99SBarry Smith }
47247c6ae99SBarry Smith 
47347c6ae99SBarry Smith #undef __FUNCT__
474aa219208SBarry Smith #define __FUNCT__ "DMDAGetLocalJacobian"
47547c6ae99SBarry Smith /*@C
476aa219208SBarry Smith        DMDAGetLocalJacobian - Gets from a DM a local jacobian
47747c6ae99SBarry Smith 
47847c6ae99SBarry Smith    Not Collective
47947c6ae99SBarry Smith 
48047c6ae99SBarry Smith    Input Parameter:
48147c6ae99SBarry Smith .  da - initial distributed array
48247c6ae99SBarry Smith 
48347c6ae99SBarry Smith    Output Parameter:
48447c6ae99SBarry Smith .  lj - the local jacobian
48547c6ae99SBarry Smith 
48647c6ae99SBarry Smith    Level: intermediate
48747c6ae99SBarry Smith 
48847c6ae99SBarry Smith .keywords:  distributed array, refine
48947c6ae99SBarry Smith 
490aa219208SBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalJacobian()
49147c6ae99SBarry Smith @*/
492aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDAGetLocalJacobian(DM da,DMDALocalFunction1 *lj)
49347c6ae99SBarry Smith {
49447c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
49547c6ae99SBarry Smith   PetscFunctionBegin;
49647c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
49747c6ae99SBarry Smith   if (lj) *lj = dd->lj;
49847c6ae99SBarry Smith   PetscFunctionReturn(0);
49947c6ae99SBarry Smith }
50047c6ae99SBarry Smith 
50147c6ae99SBarry Smith #undef __FUNCT__
502aa219208SBarry Smith #define __FUNCT__ "DMDAFormFunction"
50347c6ae99SBarry Smith /*@
504aa219208SBarry Smith     DMDAFormFunction - Evaluates a user provided function on each processor that
505aa219208SBarry Smith         share a DMDA
50647c6ae99SBarry Smith 
50747c6ae99SBarry Smith    Input Parameters:
5089a42bb27SBarry Smith +    da - the DM that defines the grid
50947c6ae99SBarry Smith .    vu - input vector
51047c6ae99SBarry Smith .    vfu - output vector
51147c6ae99SBarry Smith -    w - any user data
51247c6ae99SBarry Smith 
51347c6ae99SBarry Smith     Notes: Does NOT do ghost updates on vu upon entry
51447c6ae99SBarry Smith 
515aa219208SBarry Smith            This should eventually replace DMDAFormFunction1
51647c6ae99SBarry Smith 
51747c6ae99SBarry Smith     Level: advanced
51847c6ae99SBarry Smith 
519aa219208SBarry Smith .seealso: DMDAComputeJacobian1WithAdic()
52047c6ae99SBarry Smith 
52147c6ae99SBarry Smith @*/
522aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDAFormFunction(DM da,PetscErrorCode (*lf)(void),Vec vu,Vec vfu,void *w)
52347c6ae99SBarry Smith {
52447c6ae99SBarry Smith   PetscErrorCode ierr;
52547c6ae99SBarry Smith   void           *u,*fu;
526aa219208SBarry Smith   DMDALocalInfo  info;
527aa219208SBarry Smith   PetscErrorCode (*f)(DMDALocalInfo*,void*,void*,void*) = (PetscErrorCode (*)(DMDALocalInfo*,void*,void*,void*))lf;
52847c6ae99SBarry Smith 
52947c6ae99SBarry Smith   PetscFunctionBegin;
530aa219208SBarry Smith   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
531aa219208SBarry Smith   ierr = DMDAVecGetArray(da,vu,&u);CHKERRQ(ierr);
532aa219208SBarry Smith   ierr = DMDAVecGetArray(da,vfu,&fu);CHKERRQ(ierr);
53347c6ae99SBarry Smith 
53447c6ae99SBarry Smith   ierr = (*f)(&info,u,fu,w);
53547c6ae99SBarry Smith   if (PetscExceptionValue(ierr)) {
536aa219208SBarry Smith     PetscErrorCode pierr = DMDAVecRestoreArray(da,vu,&u);CHKERRQ(pierr);
537aa219208SBarry Smith     pierr = DMDAVecRestoreArray(da,vfu,&fu);CHKERRQ(pierr);
53847c6ae99SBarry Smith   }
53947c6ae99SBarry Smith  CHKERRQ(ierr);
54047c6ae99SBarry Smith 
541aa219208SBarry Smith   ierr = DMDAVecRestoreArray(da,vu,&u);CHKERRQ(ierr);
542aa219208SBarry Smith   ierr = DMDAVecRestoreArray(da,vfu,&fu);CHKERRQ(ierr);
54347c6ae99SBarry Smith   PetscFunctionReturn(0);
54447c6ae99SBarry Smith }
54547c6ae99SBarry Smith 
54647c6ae99SBarry Smith #undef __FUNCT__
547aa219208SBarry Smith #define __FUNCT__ "DMDAFormFunctionLocal"
54847c6ae99SBarry Smith /*@C
549aa219208SBarry Smith    DMDAFormFunctionLocal - This is a universal function evaluation routine for
5509a42bb27SBarry Smith    a local DM function.
55147c6ae99SBarry Smith 
552aa219208SBarry Smith    Collective on DMDA
55347c6ae99SBarry Smith 
55447c6ae99SBarry Smith    Input Parameters:
5559a42bb27SBarry Smith +  da - the DM context
55647c6ae99SBarry Smith .  func - The local function
55747c6ae99SBarry Smith .  X - input vector
55847c6ae99SBarry Smith .  F - function vector
55947c6ae99SBarry Smith -  ctx - A user context
56047c6ae99SBarry Smith 
56147c6ae99SBarry Smith    Level: intermediate
56247c6ae99SBarry Smith 
563aa219208SBarry Smith .seealso: DMDASetLocalFunction(), DMDASetLocalJacobian(), DMDASetLocalAdicFunction(), DMDASetLocalAdicMFFunction(),
56447c6ae99SBarry Smith           SNESSetFunction(), SNESSetJacobian()
56547c6ae99SBarry Smith 
56647c6ae99SBarry Smith @*/
567aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDAFormFunctionLocal(DM da, DMDALocalFunction1 func, Vec X, Vec F, void *ctx)
56847c6ae99SBarry Smith {
56947c6ae99SBarry Smith   Vec            localX;
570aa219208SBarry Smith   DMDALocalInfo  info;
57147c6ae99SBarry Smith   void           *u;
57247c6ae99SBarry Smith   void           *fu;
57347c6ae99SBarry Smith   PetscErrorCode ierr;
57447c6ae99SBarry Smith 
57547c6ae99SBarry Smith   PetscFunctionBegin;
5769a42bb27SBarry Smith   ierr = DMGetLocalVector(da,&localX);CHKERRQ(ierr);
57747c6ae99SBarry Smith   /*
57847c6ae99SBarry Smith      Scatter ghost points to local vector, using the 2-step process
5799a42bb27SBarry Smith         DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
58047c6ae99SBarry Smith   */
5819a42bb27SBarry Smith   ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
5829a42bb27SBarry Smith   ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
583aa219208SBarry Smith   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
584aa219208SBarry Smith   ierr = DMDAVecGetArray(da,localX,&u);CHKERRQ(ierr);
585aa219208SBarry Smith   ierr = DMDAVecGetArray(da,F,&fu);CHKERRQ(ierr);
58647c6ae99SBarry Smith   ierr = (*func)(&info,u,fu,ctx);
58747c6ae99SBarry Smith   if (PetscExceptionValue(ierr)) {
588aa219208SBarry Smith     PetscErrorCode pierr = DMDAVecRestoreArray(da,localX,&u);CHKERRQ(pierr);
589aa219208SBarry Smith     pierr = DMDAVecRestoreArray(da,F,&fu);CHKERRQ(pierr);
59047c6ae99SBarry Smith   }
59147c6ae99SBarry Smith   CHKERRQ(ierr);
592aa219208SBarry Smith   ierr = DMDAVecRestoreArray(da,localX,&u);CHKERRQ(ierr);
593aa219208SBarry Smith   ierr = DMDAVecRestoreArray(da,F,&fu);CHKERRQ(ierr);
59447c6ae99SBarry Smith   if (PetscExceptionValue(ierr)) {
5959a42bb27SBarry Smith     PetscErrorCode pierr = DMRestoreLocalVector(da,&localX);CHKERRQ(pierr);
59647c6ae99SBarry Smith   }
59747c6ae99SBarry Smith   CHKERRQ(ierr);
5989a42bb27SBarry Smith   ierr = DMRestoreLocalVector(da,&localX);CHKERRQ(ierr);
59947c6ae99SBarry Smith   PetscFunctionReturn(0);
60047c6ae99SBarry Smith }
60147c6ae99SBarry Smith 
60247c6ae99SBarry Smith #undef __FUNCT__
603aa219208SBarry Smith #define __FUNCT__ "DMDAFormFunctionLocalGhost"
60447c6ae99SBarry Smith /*@C
605aa219208SBarry Smith    DMDAFormFunctionLocalGhost - This is a universal function evaluation routine for
6069a42bb27SBarry Smith    a local DM function, but the ghost values of the output are communicated and added.
60747c6ae99SBarry Smith 
608aa219208SBarry Smith    Collective on DMDA
60947c6ae99SBarry Smith 
61047c6ae99SBarry Smith    Input Parameters:
6119a42bb27SBarry Smith +  da - the DM context
61247c6ae99SBarry Smith .  func - The local function
61347c6ae99SBarry Smith .  X - input vector
61447c6ae99SBarry Smith .  F - function vector
61547c6ae99SBarry Smith -  ctx - A user context
61647c6ae99SBarry Smith 
61747c6ae99SBarry Smith    Level: intermediate
61847c6ae99SBarry Smith 
619aa219208SBarry Smith .seealso: DMDASetLocalFunction(), DMDASetLocalJacobian(), DMDASetLocalAdicFunction(), DMDASetLocalAdicMFFunction(),
62047c6ae99SBarry Smith           SNESSetFunction(), SNESSetJacobian()
62147c6ae99SBarry Smith 
62247c6ae99SBarry Smith @*/
623aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDAFormFunctionLocalGhost(DM da, DMDALocalFunction1 func, Vec X, Vec F, void *ctx)
62447c6ae99SBarry Smith {
62547c6ae99SBarry Smith   Vec            localX, localF;
626aa219208SBarry Smith   DMDALocalInfo  info;
62747c6ae99SBarry Smith   void           *u;
62847c6ae99SBarry Smith   void           *fu;
62947c6ae99SBarry Smith   PetscErrorCode ierr;
63047c6ae99SBarry Smith 
63147c6ae99SBarry Smith   PetscFunctionBegin;
6329a42bb27SBarry Smith   ierr = DMGetLocalVector(da,&localX);CHKERRQ(ierr);
6339a42bb27SBarry Smith   ierr = DMGetLocalVector(da,&localF);CHKERRQ(ierr);
63447c6ae99SBarry Smith   /*
63547c6ae99SBarry Smith      Scatter ghost points to local vector, using the 2-step process
6369a42bb27SBarry Smith         DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
63747c6ae99SBarry Smith   */
6389a42bb27SBarry Smith   ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
6399a42bb27SBarry Smith   ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
64047c6ae99SBarry Smith   ierr = VecSet(F, 0.0);CHKERRQ(ierr);
64147c6ae99SBarry Smith   ierr = VecSet(localF, 0.0);CHKERRQ(ierr);
642aa219208SBarry Smith   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
643aa219208SBarry Smith   ierr = DMDAVecGetArray(da,localX,&u);CHKERRQ(ierr);
644aa219208SBarry Smith   ierr = DMDAVecGetArray(da,localF,&fu);CHKERRQ(ierr);
64547c6ae99SBarry Smith   ierr = (*func)(&info,u,fu,ctx);
64647c6ae99SBarry Smith   if (PetscExceptionValue(ierr)) {
647aa219208SBarry Smith     PetscErrorCode pierr = DMDAVecRestoreArray(da,localX,&u);CHKERRQ(pierr);
648aa219208SBarry Smith     pierr = DMDAVecRestoreArray(da,localF,&fu);CHKERRQ(pierr);
64947c6ae99SBarry Smith   }
65047c6ae99SBarry Smith   CHKERRQ(ierr);
6519a42bb27SBarry Smith   ierr = DMLocalToGlobalBegin(da,localF,ADD_VALUES,F);CHKERRQ(ierr);
6529a42bb27SBarry Smith   ierr = DMLocalToGlobalEnd(da,localF,ADD_VALUES,F);CHKERRQ(ierr);
653aa219208SBarry Smith   ierr = DMDAVecRestoreArray(da,localX,&u);CHKERRQ(ierr);
654aa219208SBarry Smith   ierr = DMDAVecRestoreArray(da,localF,&fu);CHKERRQ(ierr);
65547c6ae99SBarry Smith   if (PetscExceptionValue(ierr)) {
6569a42bb27SBarry Smith     PetscErrorCode pierr = DMRestoreLocalVector(da,&localX);CHKERRQ(pierr);
6579a42bb27SBarry Smith   ierr = DMRestoreLocalVector(da,&localF);CHKERRQ(ierr);
65847c6ae99SBarry Smith   }
65947c6ae99SBarry Smith   CHKERRQ(ierr);
6609a42bb27SBarry Smith   ierr = DMRestoreLocalVector(da,&localX);CHKERRQ(ierr);
6619a42bb27SBarry Smith   ierr = DMRestoreLocalVector(da,&localF);CHKERRQ(ierr);
66247c6ae99SBarry Smith   PetscFunctionReturn(0);
66347c6ae99SBarry Smith }
66447c6ae99SBarry Smith 
66547c6ae99SBarry Smith #undef __FUNCT__
666aa219208SBarry Smith #define __FUNCT__ "DMDAFormFunction1"
66747c6ae99SBarry Smith /*@
668aa219208SBarry Smith     DMDAFormFunction1 - Evaluates a user provided function on each processor that
669aa219208SBarry Smith         share a DMDA
67047c6ae99SBarry Smith 
67147c6ae99SBarry Smith    Input Parameters:
6729a42bb27SBarry Smith +    da - the DM that defines the grid
67347c6ae99SBarry Smith .    vu - input vector
67447c6ae99SBarry Smith .    vfu - output vector
67547c6ae99SBarry Smith -    w - any user data
67647c6ae99SBarry Smith 
67747c6ae99SBarry Smith     Notes: Does NOT do ghost updates on vu upon entry
67847c6ae99SBarry Smith 
67947c6ae99SBarry Smith     Level: advanced
68047c6ae99SBarry Smith 
681aa219208SBarry Smith .seealso: DMDAComputeJacobian1WithAdic()
68247c6ae99SBarry Smith 
68347c6ae99SBarry Smith @*/
684aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDAFormFunction1(DM da,Vec vu,Vec vfu,void *w)
68547c6ae99SBarry Smith {
68647c6ae99SBarry Smith   PetscErrorCode ierr;
68747c6ae99SBarry Smith   void           *u,*fu;
688aa219208SBarry Smith   DMDALocalInfo  info;
68947c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
69047c6ae99SBarry Smith 
69147c6ae99SBarry Smith   PetscFunctionBegin;
69247c6ae99SBarry Smith 
693aa219208SBarry Smith   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
694aa219208SBarry Smith   ierr = DMDAVecGetArray(da,vu,&u);CHKERRQ(ierr);
695aa219208SBarry Smith   ierr = DMDAVecGetArray(da,vfu,&fu);CHKERRQ(ierr);
69647c6ae99SBarry Smith 
69747c6ae99SBarry Smith   CHKMEMQ;
69847c6ae99SBarry Smith   ierr = (*dd->lf)(&info,u,fu,w);
69947c6ae99SBarry Smith   if (PetscExceptionValue(ierr)) {
700aa219208SBarry Smith     PetscErrorCode pierr = DMDAVecRestoreArray(da,vu,&u);CHKERRQ(pierr);
701aa219208SBarry Smith     pierr = DMDAVecRestoreArray(da,vfu,&fu);CHKERRQ(pierr);
70247c6ae99SBarry Smith   }
70347c6ae99SBarry Smith   CHKERRQ(ierr);
70447c6ae99SBarry Smith   CHKMEMQ;
70547c6ae99SBarry Smith 
706aa219208SBarry Smith   ierr = DMDAVecRestoreArray(da,vu,&u);CHKERRQ(ierr);
707aa219208SBarry Smith   ierr = DMDAVecRestoreArray(da,vfu,&fu);CHKERRQ(ierr);
70847c6ae99SBarry Smith   PetscFunctionReturn(0);
70947c6ae99SBarry Smith }
71047c6ae99SBarry Smith 
71147c6ae99SBarry Smith #undef __FUNCT__
712aa219208SBarry Smith #define __FUNCT__ "DMDAFormFunctioniTest1"
713aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDAFormFunctioniTest1(DM da,void *w)
71447c6ae99SBarry Smith {
71547c6ae99SBarry Smith   Vec            vu,fu,fui;
71647c6ae99SBarry Smith   PetscErrorCode ierr;
71747c6ae99SBarry Smith   PetscInt       i,n;
71847c6ae99SBarry Smith   PetscScalar    *ui;
71947c6ae99SBarry Smith   PetscRandom    rnd;
72047c6ae99SBarry Smith   PetscReal      norm;
72147c6ae99SBarry Smith 
72247c6ae99SBarry Smith   PetscFunctionBegin;
7239a42bb27SBarry Smith   ierr = DMGetLocalVector(da,&vu);CHKERRQ(ierr);
72447c6ae99SBarry Smith   ierr = PetscRandomCreate(PETSC_COMM_SELF,&rnd);CHKERRQ(ierr);
72547c6ae99SBarry Smith   ierr = PetscRandomSetFromOptions(rnd);CHKERRQ(ierr);
72647c6ae99SBarry Smith   ierr = VecSetRandom(vu,rnd);CHKERRQ(ierr);
72747c6ae99SBarry Smith   ierr = PetscRandomDestroy(rnd);CHKERRQ(ierr);
72847c6ae99SBarry Smith 
7299a42bb27SBarry Smith   ierr = DMGetGlobalVector(da,&fu);CHKERRQ(ierr);
7309a42bb27SBarry Smith   ierr = DMGetGlobalVector(da,&fui);CHKERRQ(ierr);
73147c6ae99SBarry Smith 
732aa219208SBarry Smith   ierr = DMDAFormFunction1(da,vu,fu,w);CHKERRQ(ierr);
73347c6ae99SBarry Smith 
73447c6ae99SBarry Smith   ierr = VecGetArray(fui,&ui);CHKERRQ(ierr);
73547c6ae99SBarry Smith   ierr = VecGetLocalSize(fui,&n);CHKERRQ(ierr);
73647c6ae99SBarry Smith   for (i=0; i<n; i++) {
737aa219208SBarry Smith     ierr = DMDAFormFunctioni1(da,i,vu,ui+i,w);CHKERRQ(ierr);
73847c6ae99SBarry Smith   }
73947c6ae99SBarry Smith   ierr = VecRestoreArray(fui,&ui);CHKERRQ(ierr);
74047c6ae99SBarry Smith 
74147c6ae99SBarry Smith   ierr = VecAXPY(fui,-1.0,fu);CHKERRQ(ierr);
74247c6ae99SBarry Smith   ierr = VecNorm(fui,NORM_2,&norm);CHKERRQ(ierr);
74347c6ae99SBarry Smith   ierr = PetscPrintf(((PetscObject)da)->comm,"Norm of difference in vectors %G\n",norm);CHKERRQ(ierr);
74447c6ae99SBarry Smith   ierr = VecView(fu,0);CHKERRQ(ierr);
74547c6ae99SBarry Smith   ierr = VecView(fui,0);CHKERRQ(ierr);
74647c6ae99SBarry Smith 
7479a42bb27SBarry Smith   ierr = DMRestoreLocalVector(da,&vu);CHKERRQ(ierr);
7489a42bb27SBarry Smith   ierr = DMRestoreGlobalVector(da,&fu);CHKERRQ(ierr);
7499a42bb27SBarry Smith   ierr = DMRestoreGlobalVector(da,&fui);CHKERRQ(ierr);
75047c6ae99SBarry Smith   PetscFunctionReturn(0);
75147c6ae99SBarry Smith }
75247c6ae99SBarry Smith 
75347c6ae99SBarry Smith #undef __FUNCT__
754aa219208SBarry Smith #define __FUNCT__ "DMDAFormFunctioni1"
75547c6ae99SBarry Smith /*@
756aa219208SBarry Smith     DMDAFormFunctioni1 - Evaluates a user provided point-wise function
75747c6ae99SBarry Smith 
75847c6ae99SBarry Smith    Input Parameters:
7599a42bb27SBarry Smith +    da - the DM that defines the grid
76047c6ae99SBarry Smith .    i - the component of the function we wish to compute (must be local)
76147c6ae99SBarry Smith .    vu - input vector
76247c6ae99SBarry Smith .    vfu - output value
76347c6ae99SBarry Smith -    w - any user data
76447c6ae99SBarry Smith 
76547c6ae99SBarry Smith     Notes: Does NOT do ghost updates on vu upon entry
76647c6ae99SBarry Smith 
76747c6ae99SBarry Smith     Level: advanced
76847c6ae99SBarry Smith 
769aa219208SBarry Smith .seealso: DMDAComputeJacobian1WithAdic()
77047c6ae99SBarry Smith 
77147c6ae99SBarry Smith @*/
772aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDAFormFunctioni1(DM da,PetscInt i,Vec vu,PetscScalar *vfu,void *w)
77347c6ae99SBarry Smith {
77447c6ae99SBarry Smith   PetscErrorCode ierr;
77547c6ae99SBarry Smith   void           *u;
776aa219208SBarry Smith   DMDALocalInfo  info;
77747c6ae99SBarry Smith   MatStencil     stencil;
77847c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
77947c6ae99SBarry Smith 
78047c6ae99SBarry Smith   PetscFunctionBegin;
78147c6ae99SBarry Smith 
782aa219208SBarry Smith   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
783aa219208SBarry Smith   ierr = DMDAVecGetArray(da,vu,&u);CHKERRQ(ierr);
78447c6ae99SBarry Smith 
78547c6ae99SBarry Smith   /* figure out stencil value from i */
78647c6ae99SBarry Smith   stencil.c = i % info.dof;
78747c6ae99SBarry Smith   stencil.i = (i % (info.xm*info.dof))/info.dof;
78847c6ae99SBarry Smith   stencil.j = (i % (info.xm*info.ym*info.dof))/(info.xm*info.dof);
78947c6ae99SBarry Smith   stencil.k = i/(info.xm*info.ym*info.dof);
79047c6ae99SBarry Smith 
79147c6ae99SBarry Smith   ierr = (*dd->lfi)(&info,&stencil,u,vfu,w);CHKERRQ(ierr);
79247c6ae99SBarry Smith 
793aa219208SBarry Smith   ierr = DMDAVecRestoreArray(da,vu,&u);CHKERRQ(ierr);
79447c6ae99SBarry Smith   PetscFunctionReturn(0);
79547c6ae99SBarry Smith }
79647c6ae99SBarry Smith 
79747c6ae99SBarry Smith #undef __FUNCT__
798aa219208SBarry Smith #define __FUNCT__ "DMDAFormFunctionib1"
79947c6ae99SBarry Smith /*@
800aa219208SBarry Smith     DMDAFormFunctionib1 - Evaluates a user provided point-block function
80147c6ae99SBarry Smith 
80247c6ae99SBarry Smith    Input Parameters:
8039a42bb27SBarry Smith +    da - the DM that defines the grid
80447c6ae99SBarry Smith .    i - the component of the function we wish to compute (must be local)
80547c6ae99SBarry Smith .    vu - input vector
80647c6ae99SBarry Smith .    vfu - output value
80747c6ae99SBarry Smith -    w - any user data
80847c6ae99SBarry Smith 
80947c6ae99SBarry Smith     Notes: Does NOT do ghost updates on vu upon entry
81047c6ae99SBarry Smith 
81147c6ae99SBarry Smith     Level: advanced
81247c6ae99SBarry Smith 
813aa219208SBarry Smith .seealso: DMDAComputeJacobian1WithAdic()
81447c6ae99SBarry Smith 
81547c6ae99SBarry Smith @*/
816aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDAFormFunctionib1(DM da,PetscInt i,Vec vu,PetscScalar *vfu,void *w)
81747c6ae99SBarry Smith {
81847c6ae99SBarry Smith   PetscErrorCode ierr;
81947c6ae99SBarry Smith   void           *u;
820aa219208SBarry Smith   DMDALocalInfo  info;
82147c6ae99SBarry Smith   MatStencil     stencil;
82247c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
82347c6ae99SBarry Smith 
82447c6ae99SBarry Smith   PetscFunctionBegin;
825aa219208SBarry Smith   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
826aa219208SBarry Smith   ierr = DMDAVecGetArray(da,vu,&u);CHKERRQ(ierr);
82747c6ae99SBarry Smith 
82847c6ae99SBarry Smith   /* figure out stencil value from i */
82947c6ae99SBarry Smith   stencil.c = i % info.dof;
83047c6ae99SBarry Smith   if (stencil.c) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONG,"Point-block functions can only be called for the entire block");
83147c6ae99SBarry Smith   stencil.i = (i % (info.xm*info.dof))/info.dof;
83247c6ae99SBarry Smith   stencil.j = (i % (info.xm*info.ym*info.dof))/(info.xm*info.dof);
83347c6ae99SBarry Smith   stencil.k = i/(info.xm*info.ym*info.dof);
83447c6ae99SBarry Smith 
83547c6ae99SBarry Smith   ierr = (*dd->lfib)(&info,&stencil,u,vfu,w);CHKERRQ(ierr);
83647c6ae99SBarry Smith 
837aa219208SBarry Smith   ierr = DMDAVecRestoreArray(da,vu,&u);CHKERRQ(ierr);
83847c6ae99SBarry Smith   PetscFunctionReturn(0);
83947c6ae99SBarry Smith }
84047c6ae99SBarry Smith 
84147c6ae99SBarry Smith #if defined(new)
84247c6ae99SBarry Smith #undef __FUNCT__
843aa219208SBarry Smith #define __FUNCT__ "DMDAGetDiagonal_MFFD"
84447c6ae99SBarry Smith /*
845aa219208SBarry Smith   DMDAGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix where local
846aa219208SBarry Smith     function lives on a DMDA
84747c6ae99SBarry Smith 
84847c6ae99SBarry Smith         y ~= (F(u + ha) - F(u))/h,
84947c6ae99SBarry Smith   where F = nonlinear function, as set by SNESSetFunction()
85047c6ae99SBarry Smith         u = current iterate
85147c6ae99SBarry Smith         h = difference interval
85247c6ae99SBarry Smith */
853aa219208SBarry Smith PetscErrorCode DMDAGetDiagonal_MFFD(DM da,Vec U,Vec a)
85447c6ae99SBarry Smith {
85547c6ae99SBarry Smith   PetscScalar    h,*aa,*ww,v;
85647c6ae99SBarry Smith   PetscReal      epsilon = PETSC_SQRT_MACHINE_EPSILON,umin = 100.0*PETSC_SQRT_MACHINE_EPSILON;
85747c6ae99SBarry Smith   PetscErrorCode ierr;
85847c6ae99SBarry Smith   PetscInt       gI,nI;
85947c6ae99SBarry Smith   MatStencil     stencil;
860aa219208SBarry Smith   DMDALocalInfo  info;
86147c6ae99SBarry Smith 
86247c6ae99SBarry Smith   PetscFunctionBegin;
86347c6ae99SBarry Smith   ierr = (*ctx->func)(0,U,a,ctx->funcctx);CHKERRQ(ierr);
86447c6ae99SBarry Smith   ierr = (*ctx->funcisetbase)(U,ctx->funcctx);CHKERRQ(ierr);
86547c6ae99SBarry Smith 
86647c6ae99SBarry Smith   ierr = VecGetArray(U,&ww);CHKERRQ(ierr);
86747c6ae99SBarry Smith   ierr = VecGetArray(a,&aa);CHKERRQ(ierr);
86847c6ae99SBarry Smith 
86947c6ae99SBarry Smith   nI = 0;
87047c6ae99SBarry Smith     h  = ww[gI];
87147c6ae99SBarry Smith     if (h == 0.0) h = 1.0;
87247c6ae99SBarry Smith #if !defined(PETSC_USE_COMPLEX)
87347c6ae99SBarry Smith     if (h < umin && h >= 0.0)      h = umin;
87447c6ae99SBarry Smith     else if (h < 0.0 && h > -umin) h = -umin;
87547c6ae99SBarry Smith #else
87647c6ae99SBarry Smith     if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0)     h = umin;
87747c6ae99SBarry Smith     else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin;
87847c6ae99SBarry Smith #endif
87947c6ae99SBarry Smith     h     *= epsilon;
88047c6ae99SBarry Smith 
88147c6ae99SBarry Smith     ww[gI] += h;
88247c6ae99SBarry Smith     ierr          = (*ctx->funci)(i,w,&v,ctx->funcctx);CHKERRQ(ierr);
88347c6ae99SBarry Smith     aa[nI]  = (v - aa[nI])/h;
88447c6ae99SBarry Smith     ww[gI] -= h;
88547c6ae99SBarry Smith     nI++;
88647c6ae99SBarry Smith   }
88747c6ae99SBarry Smith   ierr = VecRestoreArray(U,&ww);CHKERRQ(ierr);
88847c6ae99SBarry Smith   ierr = VecRestoreArray(a,&aa);CHKERRQ(ierr);
88947c6ae99SBarry Smith   PetscFunctionReturn(0);
89047c6ae99SBarry Smith }
89147c6ae99SBarry Smith #endif
89247c6ae99SBarry Smith 
89347c6ae99SBarry Smith #if defined(PETSC_HAVE_ADIC)
89447c6ae99SBarry Smith EXTERN_C_BEGIN
89547c6ae99SBarry Smith #include "adic/ad_utils.h"
89647c6ae99SBarry Smith EXTERN_C_END
89747c6ae99SBarry Smith 
89847c6ae99SBarry Smith #undef __FUNCT__
899aa219208SBarry Smith #define __FUNCT__ "DMDAComputeJacobian1WithAdic"
90047c6ae99SBarry Smith /*@C
901aa219208SBarry Smith     DMDAComputeJacobian1WithAdic - Evaluates a adiC provided Jacobian function on each processor that
902aa219208SBarry Smith         share a DMDA
90347c6ae99SBarry Smith 
90447c6ae99SBarry Smith    Input Parameters:
9059a42bb27SBarry Smith +    da - the DM that defines the grid
90647c6ae99SBarry Smith .    vu - input vector (ghosted)
90747c6ae99SBarry Smith .    J - output matrix
90847c6ae99SBarry Smith -    w - any user data
90947c6ae99SBarry Smith 
91047c6ae99SBarry Smith    Level: advanced
91147c6ae99SBarry Smith 
91247c6ae99SBarry Smith     Notes: Does NOT do ghost updates on vu upon entry
91347c6ae99SBarry Smith 
914aa219208SBarry Smith .seealso: DMDAFormFunction1()
91547c6ae99SBarry Smith 
91647c6ae99SBarry Smith @*/
917aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDAComputeJacobian1WithAdic(DM da,Vec vu,Mat J,void *w)
91847c6ae99SBarry Smith {
91947c6ae99SBarry Smith   PetscErrorCode ierr;
92047c6ae99SBarry Smith   PetscInt       gtdof,tdof;
92147c6ae99SBarry Smith   PetscScalar    *ustart;
922aa219208SBarry Smith   DMDALocalInfo  info;
92347c6ae99SBarry Smith   void           *ad_u,*ad_f,*ad_ustart,*ad_fstart;
92447c6ae99SBarry Smith   ISColoring     iscoloring;
92547c6ae99SBarry Smith 
92647c6ae99SBarry Smith   PetscFunctionBegin;
927aa219208SBarry Smith   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
92847c6ae99SBarry Smith 
92947c6ae99SBarry Smith   PetscADResetIndep();
93047c6ae99SBarry Smith 
93147c6ae99SBarry Smith   /* get space for derivative objects.  */
932aa219208SBarry Smith   ierr = DMDAGetAdicArray(da,PETSC_TRUE,&ad_u,&ad_ustart,&gtdof);CHKERRQ(ierr);
933aa219208SBarry Smith   ierr = DMDAGetAdicArray(da,PETSC_FALSE,&ad_f,&ad_fstart,&tdof);CHKERRQ(ierr);
93447c6ae99SBarry Smith   ierr = VecGetArray(vu,&ustart);CHKERRQ(ierr);
93594013140SBarry Smith   ierr = DMGetColoring(da,IS_COLORING_GHOSTED,MATAIJ,&iscoloring);CHKERRQ(ierr);
93647c6ae99SBarry Smith 
93747c6ae99SBarry Smith   PetscADSetValueAndColor(ad_ustart,gtdof,iscoloring->colors,ustart);
93847c6ae99SBarry Smith 
93947c6ae99SBarry Smith   ierr = VecRestoreArray(vu,&ustart);CHKERRQ(ierr);
94047c6ae99SBarry Smith   ierr = ISColoringDestroy(iscoloring);CHKERRQ(ierr);
94147c6ae99SBarry Smith   ierr = PetscADIncrementTotalGradSize(iscoloring->n);CHKERRQ(ierr);
94247c6ae99SBarry Smith   PetscADSetIndepDone();
94347c6ae99SBarry Smith 
944aa219208SBarry Smith   ierr = PetscLogEventBegin(DMDA_LocalADFunction,0,0,0,0);CHKERRQ(ierr);
94547c6ae99SBarry Smith   ierr = (*dd->adic_lf)(&info,ad_u,ad_f,w);CHKERRQ(ierr);
946aa219208SBarry Smith   ierr = PetscLogEventEnd(DMDA_LocalADFunction,0,0,0,0);CHKERRQ(ierr);
94747c6ae99SBarry Smith 
94847c6ae99SBarry Smith   /* stick the values into the matrix */
94947c6ae99SBarry Smith   ierr = MatSetValuesAdic(J,(PetscScalar**)ad_fstart);CHKERRQ(ierr);
95047c6ae99SBarry Smith 
95147c6ae99SBarry Smith   /* return space for derivative objects.  */
952aa219208SBarry Smith   ierr = DMDARestoreAdicArray(da,PETSC_TRUE,&ad_u,&ad_ustart,&gtdof);CHKERRQ(ierr);
953aa219208SBarry Smith   ierr = DMDARestoreAdicArray(da,PETSC_FALSE,&ad_f,&ad_fstart,&tdof);CHKERRQ(ierr);
95447c6ae99SBarry Smith   PetscFunctionReturn(0);
95547c6ae99SBarry Smith }
95647c6ae99SBarry Smith 
95747c6ae99SBarry Smith #undef __FUNCT__
958aa219208SBarry Smith #define __FUNCT__ "DMDAMultiplyByJacobian1WithAdic"
95947c6ae99SBarry Smith /*@C
960aa219208SBarry Smith     DMDAMultiplyByJacobian1WithAdic - Applies an ADIC-provided Jacobian function to a vector on
961aa219208SBarry Smith     each processor that shares a DMDA.
96247c6ae99SBarry Smith 
96347c6ae99SBarry Smith     Input Parameters:
9649a42bb27SBarry Smith +   da - the DM that defines the grid
96547c6ae99SBarry Smith .   vu - Jacobian is computed at this point (ghosted)
96647c6ae99SBarry Smith .   v - product is done on this vector (ghosted)
96747c6ae99SBarry Smith .   fu - output vector = J(vu)*v (not ghosted)
96847c6ae99SBarry Smith -   w - any user data
96947c6ae99SBarry Smith 
97047c6ae99SBarry Smith     Notes:
97147c6ae99SBarry Smith     This routine does NOT do ghost updates on vu upon entry.
97247c6ae99SBarry Smith 
97347c6ae99SBarry Smith    Level: advanced
97447c6ae99SBarry Smith 
975aa219208SBarry Smith .seealso: DMDAFormFunction1()
97647c6ae99SBarry Smith 
97747c6ae99SBarry Smith @*/
978aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDAMultiplyByJacobian1WithAdic(DM da,Vec vu,Vec v,Vec f,void *w)
97947c6ae99SBarry Smith {
98047c6ae99SBarry Smith   PetscErrorCode ierr;
98147c6ae99SBarry Smith   PetscInt       i,gtdof,tdof;
98247c6ae99SBarry Smith   PetscScalar    *avu,*av,*af,*ad_vustart,*ad_fstart;
983aa219208SBarry Smith   DMDALocalInfo  info;
98447c6ae99SBarry Smith   void           *ad_vu,*ad_f;
98547c6ae99SBarry Smith 
98647c6ae99SBarry Smith   PetscFunctionBegin;
987aa219208SBarry Smith   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
98847c6ae99SBarry Smith 
98947c6ae99SBarry Smith   /* get space for derivative objects.  */
990aa219208SBarry Smith   ierr = DMDAGetAdicMFArray(da,PETSC_TRUE,&ad_vu,&ad_vustart,&gtdof);CHKERRQ(ierr);
991aa219208SBarry Smith   ierr = DMDAGetAdicMFArray(da,PETSC_FALSE,&ad_f,&ad_fstart,&tdof);CHKERRQ(ierr);
99247c6ae99SBarry Smith 
99347c6ae99SBarry Smith   /* copy input vector into derivative object */
99447c6ae99SBarry Smith   ierr = VecGetArray(vu,&avu);CHKERRQ(ierr);
99547c6ae99SBarry Smith   ierr = VecGetArray(v,&av);CHKERRQ(ierr);
99647c6ae99SBarry Smith   for (i=0; i<gtdof; i++) {
99747c6ae99SBarry Smith     ad_vustart[2*i]   = avu[i];
99847c6ae99SBarry Smith     ad_vustart[2*i+1] = av[i];
99947c6ae99SBarry Smith   }
100047c6ae99SBarry Smith   ierr = VecRestoreArray(vu,&avu);CHKERRQ(ierr);
100147c6ae99SBarry Smith   ierr = VecRestoreArray(v,&av);CHKERRQ(ierr);
100247c6ae99SBarry Smith 
100347c6ae99SBarry Smith   PetscADResetIndep();
100447c6ae99SBarry Smith   ierr = PetscADIncrementTotalGradSize(1);CHKERRQ(ierr);
100547c6ae99SBarry Smith   PetscADSetIndepDone();
100647c6ae99SBarry Smith 
100747c6ae99SBarry Smith   ierr = (*dd->adicmf_lf)(&info,ad_vu,ad_f,w);CHKERRQ(ierr);
100847c6ae99SBarry Smith 
100947c6ae99SBarry Smith   /* stick the values into the vector */
101047c6ae99SBarry Smith   ierr = VecGetArray(f,&af);CHKERRQ(ierr);
101147c6ae99SBarry Smith   for (i=0; i<tdof; i++) {
101247c6ae99SBarry Smith     af[i] = ad_fstart[2*i+1];
101347c6ae99SBarry Smith   }
101447c6ae99SBarry Smith   ierr = VecRestoreArray(f,&af);CHKERRQ(ierr);
101547c6ae99SBarry Smith 
101647c6ae99SBarry Smith   /* return space for derivative objects.  */
1017aa219208SBarry Smith   ierr = DMDARestoreAdicMFArray(da,PETSC_TRUE,&ad_vu,&ad_vustart,&gtdof);CHKERRQ(ierr);
1018aa219208SBarry Smith   ierr = DMDARestoreAdicMFArray(da,PETSC_FALSE,&ad_f,&ad_fstart,&tdof);CHKERRQ(ierr);
101947c6ae99SBarry Smith   PetscFunctionReturn(0);
102047c6ae99SBarry Smith }
102147c6ae99SBarry Smith #endif
102247c6ae99SBarry Smith 
102347c6ae99SBarry Smith #undef __FUNCT__
1024aa219208SBarry Smith #define __FUNCT__ "DMDAComputeJacobian1"
102547c6ae99SBarry Smith /*@
1026aa219208SBarry Smith     DMDAComputeJacobian1 - Evaluates a local Jacobian function on each processor that
1027aa219208SBarry Smith         share a DMDA
102847c6ae99SBarry Smith 
102947c6ae99SBarry Smith    Input Parameters:
10309a42bb27SBarry Smith +    da - the DM that defines the grid
103147c6ae99SBarry Smith .    vu - input vector (ghosted)
103247c6ae99SBarry Smith .    J - output matrix
103347c6ae99SBarry Smith -    w - any user data
103447c6ae99SBarry Smith 
103547c6ae99SBarry Smith     Notes: Does NOT do ghost updates on vu upon entry
103647c6ae99SBarry Smith 
103747c6ae99SBarry Smith     Level: advanced
103847c6ae99SBarry Smith 
1039aa219208SBarry Smith .seealso: DMDAFormFunction1()
104047c6ae99SBarry Smith 
104147c6ae99SBarry Smith @*/
1042aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDAComputeJacobian1(DM da,Vec vu,Mat J,void *w)
104347c6ae99SBarry Smith {
104447c6ae99SBarry Smith   PetscErrorCode ierr;
104547c6ae99SBarry Smith   void           *u;
1046aa219208SBarry Smith   DMDALocalInfo  info;
104747c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
104847c6ae99SBarry Smith 
104947c6ae99SBarry Smith   PetscFunctionBegin;
1050aa219208SBarry Smith   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
1051aa219208SBarry Smith   ierr = DMDAVecGetArray(da,vu,&u);CHKERRQ(ierr);
105247c6ae99SBarry Smith   ierr = (*dd->lj)(&info,u,J,w);CHKERRQ(ierr);
1053aa219208SBarry Smith   ierr = DMDAVecRestoreArray(da,vu,&u);CHKERRQ(ierr);
105447c6ae99SBarry Smith   PetscFunctionReturn(0);
105547c6ae99SBarry Smith }
105647c6ae99SBarry Smith 
105747c6ae99SBarry Smith 
105847c6ae99SBarry Smith #undef __FUNCT__
1059aa219208SBarry Smith #define __FUNCT__ "DMDAComputeJacobian1WithAdifor"
106047c6ae99SBarry Smith /*
1061aa219208SBarry Smith     DMDAComputeJacobian1WithAdifor - Evaluates a ADIFOR provided Jacobian local function on each processor that
1062aa219208SBarry Smith         share a DMDA
106347c6ae99SBarry Smith 
106447c6ae99SBarry Smith    Input Parameters:
10659a42bb27SBarry Smith +    da - the DM that defines the grid
106647c6ae99SBarry Smith .    vu - input vector (ghosted)
106747c6ae99SBarry Smith .    J - output matrix
106847c6ae99SBarry Smith -    w - any user data
106947c6ae99SBarry Smith 
107047c6ae99SBarry Smith     Notes: Does NOT do ghost updates on vu upon entry
107147c6ae99SBarry Smith 
1072aa219208SBarry Smith .seealso: DMDAFormFunction1()
107347c6ae99SBarry Smith 
107447c6ae99SBarry Smith */
1075aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDAComputeJacobian1WithAdifor(DM da,Vec vu,Mat J,void *w)
107647c6ae99SBarry Smith {
107747c6ae99SBarry Smith   PetscErrorCode  ierr;
107847c6ae99SBarry Smith   PetscInt        i,Nc,N;
107947c6ae99SBarry Smith   ISColoringValue *color;
1080aa219208SBarry Smith   DMDALocalInfo   info;
108147c6ae99SBarry Smith   PetscScalar     *u,*g_u,*g_f,*f = 0,*p_u;
108247c6ae99SBarry Smith   ISColoring      iscoloring;
108347c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
1084aa219208SBarry Smith   void            (*lf)(PetscInt*,DMDALocalInfo*,PetscScalar*,PetscScalar*,PetscInt*,PetscScalar*,PetscScalar*,PetscInt*,void*,PetscErrorCode*) =
1085aa219208SBarry Smith                   (void (*)(PetscInt*,DMDALocalInfo*,PetscScalar*,PetscScalar*,PetscInt*,PetscScalar*,PetscScalar*,PetscInt*,void*,PetscErrorCode*))*dd->adifor_lf;
108647c6ae99SBarry Smith 
108747c6ae99SBarry Smith   PetscFunctionBegin;
108894013140SBarry Smith   ierr = DMGetColoring(da,IS_COLORING_GHOSTED,MATAIJ,&iscoloring);CHKERRQ(ierr);
108947c6ae99SBarry Smith   Nc   = iscoloring->n;
1090aa219208SBarry Smith   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
109147c6ae99SBarry Smith   N    = info.gxm*info.gym*info.gzm*info.dof;
109247c6ae99SBarry Smith 
109347c6ae99SBarry Smith   /* get space for derivative objects.  */
109447c6ae99SBarry Smith   ierr  = PetscMalloc(Nc*info.gxm*info.gym*info.gzm*info.dof*sizeof(PetscScalar),&g_u);CHKERRQ(ierr);
109547c6ae99SBarry Smith   ierr  = PetscMemzero(g_u,Nc*info.gxm*info.gym*info.gzm*info.dof*sizeof(PetscScalar));CHKERRQ(ierr);
109647c6ae99SBarry Smith   p_u   = g_u;
109747c6ae99SBarry Smith   color = iscoloring->colors;
109847c6ae99SBarry Smith   for (i=0; i<N; i++) {
109947c6ae99SBarry Smith     p_u[*color++] = 1.0;
110047c6ae99SBarry Smith     p_u          += Nc;
110147c6ae99SBarry Smith   }
110247c6ae99SBarry Smith   ierr = ISColoringDestroy(iscoloring);CHKERRQ(ierr);
110347c6ae99SBarry 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);
110447c6ae99SBarry Smith 
110547c6ae99SBarry Smith   /* Seed the input array g_u with coloring information */
110647c6ae99SBarry Smith 
110747c6ae99SBarry Smith   ierr = VecGetArray(vu,&u);CHKERRQ(ierr);
110847c6ae99SBarry Smith   (lf)(&Nc,&info,u,g_u,&Nc,f,g_f,&Nc,w,&ierr);CHKERRQ(ierr);
110947c6ae99SBarry Smith   ierr = VecRestoreArray(vu,&u);CHKERRQ(ierr);
111047c6ae99SBarry Smith 
111147c6ae99SBarry Smith   /* stick the values into the matrix */
111247c6ae99SBarry Smith   /* PetscScalarView(Nc*info.xm*info.ym,g_f,0); */
111347c6ae99SBarry Smith   ierr = MatSetValuesAdifor(J,Nc,g_f);CHKERRQ(ierr);
111447c6ae99SBarry Smith 
111547c6ae99SBarry Smith   /* return space for derivative objects.  */
111647c6ae99SBarry Smith   ierr = PetscFree(g_u);CHKERRQ(ierr);
111747c6ae99SBarry Smith   ierr = PetscFree2(g_f,f);CHKERRQ(ierr);
111847c6ae99SBarry Smith   PetscFunctionReturn(0);
111947c6ae99SBarry Smith }
112047c6ae99SBarry Smith 
112147c6ae99SBarry Smith #undef __FUNCT__
1122aa219208SBarry Smith #define __FUNCT__ "DMDAFormJacobianLocal"
112347c6ae99SBarry Smith /*@C
1124aa219208SBarry Smith    DMDAFormjacobianLocal - This is a universal Jacobian evaluation routine for
11259a42bb27SBarry Smith    a local DM function.
112647c6ae99SBarry Smith 
1127aa219208SBarry Smith    Collective on DMDA
112847c6ae99SBarry Smith 
112947c6ae99SBarry Smith    Input Parameters:
11309a42bb27SBarry Smith +  da - the DM context
113147c6ae99SBarry Smith .  func - The local function
113247c6ae99SBarry Smith .  X - input vector
113347c6ae99SBarry Smith .  J - Jacobian matrix
113447c6ae99SBarry Smith -  ctx - A user context
113547c6ae99SBarry Smith 
113647c6ae99SBarry Smith    Level: intermediate
113747c6ae99SBarry Smith 
1138aa219208SBarry Smith .seealso: DMDASetLocalFunction(), DMDASetLocalJacobian(), DMDASetLocalAdicFunction(), DMDASetLocalAdicMFFunction(),
113947c6ae99SBarry Smith           SNESSetFunction(), SNESSetJacobian()
114047c6ae99SBarry Smith 
114147c6ae99SBarry Smith @*/
1142aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDAFormJacobianLocal(DM da, DMDALocalFunction1 func, Vec X, Mat J, void *ctx)
114347c6ae99SBarry Smith {
114447c6ae99SBarry Smith   Vec            localX;
1145aa219208SBarry Smith   DMDALocalInfo  info;
114647c6ae99SBarry Smith   void           *u;
114747c6ae99SBarry Smith   PetscErrorCode ierr;
114847c6ae99SBarry Smith 
114947c6ae99SBarry Smith   PetscFunctionBegin;
11509a42bb27SBarry Smith   ierr = DMGetLocalVector(da,&localX);CHKERRQ(ierr);
115147c6ae99SBarry Smith   /*
115247c6ae99SBarry Smith      Scatter ghost points to local vector, using the 2-step process
11539a42bb27SBarry Smith         DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
115447c6ae99SBarry Smith   */
11559a42bb27SBarry Smith   ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
11569a42bb27SBarry Smith   ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
1157aa219208SBarry Smith   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
1158aa219208SBarry Smith   ierr = DMDAVecGetArray(da,localX,&u);CHKERRQ(ierr);
115947c6ae99SBarry Smith   ierr = (*func)(&info,u,J,ctx);
116047c6ae99SBarry Smith   if (PetscExceptionValue(ierr)) {
1161aa219208SBarry Smith     PetscErrorCode pierr = DMDAVecRestoreArray(da,localX,&u);CHKERRQ(pierr);
116247c6ae99SBarry Smith   }
116347c6ae99SBarry Smith   CHKERRQ(ierr);
1164aa219208SBarry Smith   ierr = DMDAVecRestoreArray(da,localX,&u);CHKERRQ(ierr);
116547c6ae99SBarry Smith   if (PetscExceptionValue(ierr)) {
11669a42bb27SBarry Smith     PetscErrorCode pierr = DMRestoreLocalVector(da,&localX);CHKERRQ(pierr);
116747c6ae99SBarry Smith   }
116847c6ae99SBarry Smith   CHKERRQ(ierr);
11699a42bb27SBarry Smith   ierr = DMRestoreLocalVector(da,&localX);CHKERRQ(ierr);
117047c6ae99SBarry Smith   PetscFunctionReturn(0);
117147c6ae99SBarry Smith }
117247c6ae99SBarry Smith 
117347c6ae99SBarry Smith #undef __FUNCT__
1174aa219208SBarry Smith #define __FUNCT__ "DMDAMultiplyByJacobian1WithAD"
117547c6ae99SBarry Smith /*@C
1176aa219208SBarry Smith     DMDAMultiplyByJacobian1WithAD - Applies a Jacobian function supplied by ADIFOR or ADIC
1177aa219208SBarry Smith     to a vector on each processor that shares a DMDA.
117847c6ae99SBarry Smith 
117947c6ae99SBarry Smith    Input Parameters:
11809a42bb27SBarry Smith +    da - the DM that defines the grid
118147c6ae99SBarry Smith .    vu - Jacobian is computed at this point (ghosted)
118247c6ae99SBarry Smith .    v - product is done on this vector (ghosted)
118347c6ae99SBarry Smith .    fu - output vector = J(vu)*v (not ghosted)
118447c6ae99SBarry Smith -    w - any user data
118547c6ae99SBarry Smith 
118647c6ae99SBarry Smith     Notes:
118747c6ae99SBarry Smith     This routine does NOT do ghost updates on vu and v upon entry.
118847c6ae99SBarry Smith 
1189aa219208SBarry Smith     Automatically calls DMDAMultiplyByJacobian1WithAdifor() or DMDAMultiplyByJacobian1WithAdic()
1190aa219208SBarry Smith     depending on whether DMDASetLocalAdicMFFunction() or DMDASetLocalAdiforMFFunction() was called.
119147c6ae99SBarry Smith 
119247c6ae99SBarry Smith    Level: advanced
119347c6ae99SBarry Smith 
1194aa219208SBarry Smith .seealso: DMDAFormFunction1(), DMDAMultiplyByJacobian1WithAdifor(), DMDAMultiplyByJacobian1WithAdic()
119547c6ae99SBarry Smith 
119647c6ae99SBarry Smith @*/
1197aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDAMultiplyByJacobian1WithAD(DM da,Vec u,Vec v,Vec f,void *w)
119847c6ae99SBarry Smith {
119947c6ae99SBarry Smith   PetscErrorCode ierr;
120047c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
120147c6ae99SBarry Smith 
120247c6ae99SBarry Smith   PetscFunctionBegin;
120347c6ae99SBarry Smith   if (dd->adicmf_lf) {
120447c6ae99SBarry Smith #if defined(PETSC_HAVE_ADIC)
1205aa219208SBarry Smith     ierr = DMDAMultiplyByJacobian1WithAdic(da,u,v,f,w);CHKERRQ(ierr);
120647c6ae99SBarry Smith #else
120747c6ae99SBarry Smith     SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP_SYS,"Requires ADIC to be installed and cannot use complex numbers");
120847c6ae99SBarry Smith #endif
120947c6ae99SBarry Smith   } else if (dd->adiformf_lf) {
1210aa219208SBarry Smith     ierr = DMDAMultiplyByJacobian1WithAdifor(da,u,v,f,w);CHKERRQ(ierr);
121147c6ae99SBarry Smith   } else {
1212aa219208SBarry Smith     SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ORDER,"Must call DMDASetLocalAdiforMFFunction() or DMDASetLocalAdicMFFunction() before using");
121347c6ae99SBarry Smith   }
121447c6ae99SBarry Smith   PetscFunctionReturn(0);
121547c6ae99SBarry Smith }
121647c6ae99SBarry Smith 
121747c6ae99SBarry Smith 
121847c6ae99SBarry Smith #undef __FUNCT__
1219aa219208SBarry Smith #define __FUNCT__ "DMDAMultiplyByJacobian1WithAdifor"
122047c6ae99SBarry Smith /*@C
1221aa219208SBarry Smith     DMDAMultiplyByJacobian1WithAdifor - Applies a ADIFOR provided Jacobian function on each processor that
12229a42bb27SBarry Smith         share a DM to a vector
122347c6ae99SBarry Smith 
122447c6ae99SBarry Smith    Input Parameters:
12259a42bb27SBarry Smith +    da - the DM that defines the grid
122647c6ae99SBarry Smith .    vu - Jacobian is computed at this point (ghosted)
122747c6ae99SBarry Smith .    v - product is done on this vector (ghosted)
122847c6ae99SBarry Smith .    fu - output vector = J(vu)*v (not ghosted)
122947c6ae99SBarry Smith -    w - any user data
123047c6ae99SBarry Smith 
123147c6ae99SBarry Smith     Notes: Does NOT do ghost updates on vu and v upon entry
123247c6ae99SBarry Smith 
123347c6ae99SBarry Smith    Level: advanced
123447c6ae99SBarry Smith 
1235aa219208SBarry Smith .seealso: DMDAFormFunction1()
123647c6ae99SBarry Smith 
123747c6ae99SBarry Smith @*/
1238aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDAMultiplyByJacobian1WithAdifor(DM da,Vec u,Vec v,Vec f,void *w)
123947c6ae99SBarry Smith {
124047c6ae99SBarry Smith   PetscErrorCode ierr;
124147c6ae99SBarry Smith   PetscScalar    *au,*av,*af,*awork;
124247c6ae99SBarry Smith   Vec            work;
1243aa219208SBarry Smith   DMDALocalInfo  info;
124447c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
1245aa219208SBarry Smith   void           (*lf)(DMDALocalInfo*,PetscScalar*,PetscScalar*,PetscScalar*,PetscScalar*,void*,PetscErrorCode*) =
1246aa219208SBarry Smith                  (void (*)(DMDALocalInfo*,PetscScalar*,PetscScalar*,PetscScalar*,PetscScalar*,void*,PetscErrorCode*))*dd->adiformf_lf;
124747c6ae99SBarry Smith 
124847c6ae99SBarry Smith   PetscFunctionBegin;
1249aa219208SBarry Smith   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
125047c6ae99SBarry Smith 
12519a42bb27SBarry Smith   ierr = DMGetGlobalVector(da,&work);CHKERRQ(ierr);
125247c6ae99SBarry Smith   ierr = VecGetArray(u,&au);CHKERRQ(ierr);
125347c6ae99SBarry Smith   ierr = VecGetArray(v,&av);CHKERRQ(ierr);
125447c6ae99SBarry Smith   ierr = VecGetArray(f,&af);CHKERRQ(ierr);
125547c6ae99SBarry Smith   ierr = VecGetArray(work,&awork);CHKERRQ(ierr);
125647c6ae99SBarry Smith   (lf)(&info,au,av,awork,af,w,&ierr);CHKERRQ(ierr);
125747c6ae99SBarry Smith   ierr = VecRestoreArray(u,&au);CHKERRQ(ierr);
125847c6ae99SBarry Smith   ierr = VecRestoreArray(v,&av);CHKERRQ(ierr);
125947c6ae99SBarry Smith   ierr = VecRestoreArray(f,&af);CHKERRQ(ierr);
126047c6ae99SBarry Smith   ierr = VecRestoreArray(work,&awork);CHKERRQ(ierr);
12619a42bb27SBarry Smith   ierr = DMRestoreGlobalVector(da,&work);CHKERRQ(ierr);
126247c6ae99SBarry Smith 
126347c6ae99SBarry Smith   PetscFunctionReturn(0);
126447c6ae99SBarry Smith }
126547c6ae99SBarry Smith 
126647c6ae99SBarry Smith #undef __FUNCT__
12679a42bb27SBarry Smith #define __FUNCT__ "DMSetUp_DA_2D"
12689a42bb27SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMSetUp_DA_2D(DM da)
126947c6ae99SBarry Smith {
127047c6ae99SBarry Smith   DM_DA                  *dd = (DM_DA*)da->data;
127147c6ae99SBarry Smith   const PetscInt         M            = dd->M;
127247c6ae99SBarry Smith   const PetscInt         N            = dd->N;
127347c6ae99SBarry Smith   PetscInt               m            = dd->m;
127447c6ae99SBarry Smith   PetscInt               n            = dd->n;
127547c6ae99SBarry Smith   const PetscInt         dof          = dd->w;
127647c6ae99SBarry Smith   const PetscInt         s            = dd->s;
1277aa219208SBarry Smith   const DMDAPeriodicType wrap         = dd->wrap;
1278aa219208SBarry Smith   const DMDAStencilType  stencil_type = dd->stencil_type;
127947c6ae99SBarry Smith   PetscInt               *lx           = dd->lx;
128047c6ae99SBarry Smith   PetscInt               *ly           = dd->ly;
128147c6ae99SBarry Smith   MPI_Comm               comm;
128247c6ae99SBarry Smith   PetscMPIInt            rank,size;
128347c6ae99SBarry Smith   PetscInt               xs,xe,ys,ye,x,y,Xs,Xe,Ys,Ye,start,end;
128447c6ae99SBarry Smith   PetscInt               up,down,left,i,n0,n1,n2,n3,n5,n6,n7,n8,*idx,nn;
128547c6ae99SBarry Smith   PetscInt               xbase,*bases,*ldims,j,x_t,y_t,s_t,base,count;
128647c6ae99SBarry Smith   PetscInt               s_x,s_y; /* s proportionalized to w */
128747c6ae99SBarry Smith   PetscInt               sn0 = 0,sn2 = 0,sn6 = 0,sn8 = 0;
128847c6ae99SBarry Smith   Vec                    local,global;
128947c6ae99SBarry Smith   VecScatter             ltog,gtol;
129047c6ae99SBarry Smith   IS                     to,from;
129147c6ae99SBarry Smith   PetscErrorCode         ierr;
129247c6ae99SBarry Smith 
129347c6ae99SBarry Smith   PetscFunctionBegin;
129447c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
129547c6ae99SBarry Smith 
129647c6ae99SBarry Smith   if (dof < 1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Must have 1 or more degrees of freedom per node: %D",dof);
129747c6ae99SBarry Smith   if (s < 0) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Stencil width cannot be negative: %D",s);
129847c6ae99SBarry Smith 
129947c6ae99SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
130047c6ae99SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
130147c6ae99SBarry Smith 
130247c6ae99SBarry Smith   dd->dim         = 2;
130347c6ae99SBarry Smith   ierr = PetscMalloc(dof*sizeof(char*),&dd->fieldname);CHKERRQ(ierr);
130447c6ae99SBarry Smith   ierr = PetscMemzero(dd->fieldname,dof*sizeof(char*));CHKERRQ(ierr);
130547c6ae99SBarry Smith 
130647c6ae99SBarry Smith   if (m != PETSC_DECIDE) {
130747c6ae99SBarry Smith     if (m < 1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in X direction: %D",m);
130847c6ae99SBarry Smith     else if (m > size) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in X direction: %D %d",m,size);
130947c6ae99SBarry Smith   }
131047c6ae99SBarry Smith   if (n != PETSC_DECIDE) {
131147c6ae99SBarry Smith     if (n < 1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Y direction: %D",n);
131247c6ae99SBarry Smith     else if (n > size) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Y direction: %D %d",n,size);
131347c6ae99SBarry Smith   }
131447c6ae99SBarry Smith 
131547c6ae99SBarry Smith   if (m == PETSC_DECIDE || n == PETSC_DECIDE) {
131647c6ae99SBarry Smith     if (n != PETSC_DECIDE) {
131747c6ae99SBarry Smith       m = size/n;
131847c6ae99SBarry Smith     } else if (m != PETSC_DECIDE) {
131947c6ae99SBarry Smith       n = size/m;
132047c6ae99SBarry Smith     } else {
132147c6ae99SBarry Smith       /* try for squarish distribution */
132247c6ae99SBarry Smith       m = (PetscInt)(0.5 + sqrt(((double)M)*((double)size)/((double)N)));
132347c6ae99SBarry Smith       if (!m) m = 1;
132447c6ae99SBarry Smith       while (m > 0) {
132547c6ae99SBarry Smith 	n = size/m;
132647c6ae99SBarry Smith 	if (m*n == size) break;
132747c6ae99SBarry Smith 	m--;
132847c6ae99SBarry Smith       }
132947c6ae99SBarry Smith       if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;}
133047c6ae99SBarry Smith     }
133147c6ae99SBarry 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 ");
133247c6ae99SBarry Smith   } else if (m*n != size) SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"Given Bad partition");
133347c6ae99SBarry Smith 
133447c6ae99SBarry Smith   if (M < m) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in x direction is too fine! %D %D",M,m);
133547c6ae99SBarry Smith   if (N < n) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in y direction is too fine! %D %D",N,n);
133647c6ae99SBarry Smith 
133747c6ae99SBarry Smith   /*
133847c6ae99SBarry Smith      Determine locally owned region
133947c6ae99SBarry Smith      xs is the first local node number, x is the number of local nodes
134047c6ae99SBarry Smith   */
134147c6ae99SBarry Smith   if (!lx) {
134247c6ae99SBarry Smith     ierr = PetscMalloc(m*sizeof(PetscInt), &dd->lx);CHKERRQ(ierr);
134347c6ae99SBarry Smith     lx = dd->lx;
134447c6ae99SBarry Smith     for (i=0; i<m; i++) {
134547c6ae99SBarry Smith       lx[i] = M/m + ((M % m) > i);
134647c6ae99SBarry Smith     }
134747c6ae99SBarry Smith   }
134847c6ae99SBarry Smith   x  = lx[rank % m];
134947c6ae99SBarry Smith   xs = 0;
135047c6ae99SBarry Smith   for (i=0; i<(rank % m); i++) {
135147c6ae99SBarry Smith     xs += lx[i];
135247c6ae99SBarry Smith   }
135347c6ae99SBarry Smith #if defined(PETSC_USE_DEBUG)
135447c6ae99SBarry Smith   left = xs;
135547c6ae99SBarry Smith   for (i=(rank % m); i<m; i++) {
135647c6ae99SBarry Smith     left += lx[i];
135747c6ae99SBarry Smith   }
135847c6ae99SBarry 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);
135947c6ae99SBarry Smith #endif
136047c6ae99SBarry Smith 
136147c6ae99SBarry Smith   /*
136247c6ae99SBarry Smith      Determine locally owned region
136347c6ae99SBarry Smith      ys is the first local node number, y is the number of local nodes
136447c6ae99SBarry Smith   */
136547c6ae99SBarry Smith   if (!ly) {
136647c6ae99SBarry Smith     ierr = PetscMalloc(n*sizeof(PetscInt), &dd->ly);CHKERRQ(ierr);
136747c6ae99SBarry Smith     ly = dd->ly;
136847c6ae99SBarry Smith     for (i=0; i<n; i++) {
136947c6ae99SBarry Smith       ly[i] = N/n + ((N % n) > i);
137047c6ae99SBarry Smith     }
137147c6ae99SBarry Smith   }
137247c6ae99SBarry Smith   y  = ly[rank/m];
137347c6ae99SBarry Smith   ys = 0;
137447c6ae99SBarry Smith   for (i=0; i<(rank/m); i++) {
137547c6ae99SBarry Smith     ys += ly[i];
137647c6ae99SBarry Smith   }
137747c6ae99SBarry Smith #if defined(PETSC_USE_DEBUG)
137847c6ae99SBarry Smith   left = ys;
137947c6ae99SBarry Smith   for (i=(rank/m); i<n; i++) {
138047c6ae99SBarry Smith     left += ly[i];
138147c6ae99SBarry Smith   }
138247c6ae99SBarry 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);
138347c6ae99SBarry Smith #endif
138447c6ae99SBarry Smith 
138547c6ae99SBarry 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);
138647c6ae99SBarry 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);
138747c6ae99SBarry Smith   xe = xs + x;
138847c6ae99SBarry Smith   ye = ys + y;
138947c6ae99SBarry Smith 
139047c6ae99SBarry Smith   /* determine ghost region */
139147c6ae99SBarry Smith   /* Assume No Periodicity */
139247c6ae99SBarry Smith   if (xs-s > 0) Xs = xs - s; else Xs = 0;
139347c6ae99SBarry Smith   if (ys-s > 0) Ys = ys - s; else Ys = 0;
139447c6ae99SBarry Smith   if (xe+s <= M) Xe = xe + s; else Xe = M;
139547c6ae99SBarry Smith   if (ye+s <= N) Ye = ye + s; else Ye = N;
139647c6ae99SBarry Smith 
139747c6ae99SBarry Smith   /* X Periodic */
1398aa219208SBarry Smith   if (DMDAXPeriodic(wrap)){
139947c6ae99SBarry Smith     Xs = xs - s;
140047c6ae99SBarry Smith     Xe = xe + s;
140147c6ae99SBarry Smith   }
140247c6ae99SBarry Smith 
140347c6ae99SBarry Smith   /* Y Periodic */
1404aa219208SBarry Smith   if (DMDAYPeriodic(wrap)){
140547c6ae99SBarry Smith     Ys = ys - s;
140647c6ae99SBarry Smith     Ye = ye + s;
140747c6ae99SBarry Smith   }
140847c6ae99SBarry Smith 
140947c6ae99SBarry Smith   /* Resize all X parameters to reflect w */
141047c6ae99SBarry Smith   x   *= dof;
141147c6ae99SBarry Smith   xs  *= dof;
141247c6ae99SBarry Smith   xe  *= dof;
141347c6ae99SBarry Smith   Xs  *= dof;
141447c6ae99SBarry Smith   Xe  *= dof;
141547c6ae99SBarry Smith   s_x = s*dof;
141647c6ae99SBarry Smith   s_y = s;
141747c6ae99SBarry Smith 
141847c6ae99SBarry Smith   /* determine starting point of each processor */
141947c6ae99SBarry Smith   nn    = x*y;
142047c6ae99SBarry Smith   ierr  = PetscMalloc2(size+1,PetscInt,&bases,size,PetscInt,&ldims);CHKERRQ(ierr);
142147c6ae99SBarry Smith   ierr  = MPI_Allgather(&nn,1,MPIU_INT,ldims,1,MPIU_INT,comm);CHKERRQ(ierr);
142247c6ae99SBarry Smith   bases[0] = 0;
142347c6ae99SBarry Smith   for (i=1; i<=size; i++) {
142447c6ae99SBarry Smith     bases[i] = ldims[i-1];
142547c6ae99SBarry Smith   }
142647c6ae99SBarry Smith   for (i=1; i<=size; i++) {
142747c6ae99SBarry Smith     bases[i] += bases[i-1];
142847c6ae99SBarry Smith   }
142947c6ae99SBarry Smith 
143047c6ae99SBarry Smith   /* allocate the base parallel and sequential vectors */
143147c6ae99SBarry Smith   dd->Nlocal = x*y;
143247c6ae99SBarry Smith   ierr = VecCreateMPIWithArray(comm,dd->Nlocal,PETSC_DECIDE,0,&global);CHKERRQ(ierr);
143347c6ae99SBarry Smith   ierr = VecSetBlockSize(global,dof);CHKERRQ(ierr);
143447c6ae99SBarry Smith   dd->nlocal = (Xe-Xs)*(Ye-Ys);
143547c6ae99SBarry Smith   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,dd->nlocal,0,&local);CHKERRQ(ierr);
143647c6ae99SBarry Smith   ierr = VecSetBlockSize(local,dof);CHKERRQ(ierr);
143747c6ae99SBarry Smith 
143847c6ae99SBarry Smith 
143947c6ae99SBarry Smith   /* generate appropriate vector scatters */
144047c6ae99SBarry Smith   /* local to global inserts non-ghost point region into global */
144147c6ae99SBarry Smith   ierr = VecGetOwnershipRange(global,&start,&end);CHKERRQ(ierr);
144247c6ae99SBarry Smith   ierr = ISCreateStride(comm,x*y,start,1,&to);CHKERRQ(ierr);
144347c6ae99SBarry Smith 
144447c6ae99SBarry Smith   left  = xs - Xs; down  = ys - Ys; up    = down + y;
144547c6ae99SBarry Smith   ierr = PetscMalloc(x*(up - down)*sizeof(PetscInt),&idx);CHKERRQ(ierr);
144647c6ae99SBarry Smith   count = 0;
144747c6ae99SBarry Smith   for (i=down; i<up; i++) {
144847c6ae99SBarry Smith     for (j=0; j<x/dof; j++) {
144947c6ae99SBarry Smith       idx[count++] = (left + i*(Xe-Xs) + j*dof)/dof;
145047c6ae99SBarry Smith     }
145147c6ae99SBarry Smith   }
145247c6ae99SBarry Smith   ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&from);CHKERRQ(ierr);
145347c6ae99SBarry Smith 
145447c6ae99SBarry Smith   ierr = VecScatterCreate(local,from,global,to,&ltog);CHKERRQ(ierr);
145547c6ae99SBarry Smith   ierr = PetscLogObjectParent(dd,ltog);CHKERRQ(ierr);
145647c6ae99SBarry Smith   ierr = ISDestroy(from);CHKERRQ(ierr);
145747c6ae99SBarry Smith   ierr = ISDestroy(to);CHKERRQ(ierr);
145847c6ae99SBarry Smith 
145947c6ae99SBarry Smith   /* global to local must include ghost points */
1460aa219208SBarry Smith   if (stencil_type == DMDA_STENCIL_BOX) {
146147c6ae99SBarry Smith     ierr = ISCreateStride(comm,(Xe-Xs)*(Ye-Ys),0,1,&to);CHKERRQ(ierr);
146247c6ae99SBarry Smith   } else {
146347c6ae99SBarry Smith     /* must drop into cross shape region */
146447c6ae99SBarry Smith     /*       ---------|
146547c6ae99SBarry Smith             |  top    |
146647c6ae99SBarry Smith          |---         ---|
146747c6ae99SBarry Smith          |   middle      |
146847c6ae99SBarry Smith          |               |
146947c6ae99SBarry Smith          ----         ----
147047c6ae99SBarry Smith             | bottom  |
147147c6ae99SBarry Smith             -----------
147247c6ae99SBarry Smith         Xs xs        xe  Xe */
147347c6ae99SBarry Smith     /* bottom */
147447c6ae99SBarry Smith     left  = xs - Xs; down = ys - Ys; up    = down + y;
147547c6ae99SBarry Smith     count = down*(xe-xs) + (up-down)*(Xe-Xs) + (Ye-Ys-up)*(xe-xs);
147647c6ae99SBarry Smith     ierr  = PetscMalloc(count*sizeof(PetscInt)/dof,&idx);CHKERRQ(ierr);
147747c6ae99SBarry Smith     count = 0;
147847c6ae99SBarry Smith     for (i=0; i<down; i++) {
147947c6ae99SBarry Smith       for (j=0; j<xe-xs; j += dof) {
148047c6ae99SBarry Smith         idx[count++] = left + i*(Xe-Xs) + j;
148147c6ae99SBarry Smith       }
148247c6ae99SBarry Smith     }
148347c6ae99SBarry Smith     /* middle */
148447c6ae99SBarry Smith     for (i=down; i<up; i++) {
148547c6ae99SBarry Smith       for (j=0; j<Xe-Xs; j += dof) {
148647c6ae99SBarry Smith         idx[count++] = i*(Xe-Xs) + j;
148747c6ae99SBarry Smith       }
148847c6ae99SBarry Smith     }
148947c6ae99SBarry Smith     /* top */
149047c6ae99SBarry Smith     for (i=up; i<Ye-Ys; i++) {
149147c6ae99SBarry Smith       for (j=0; j<xe-xs; j += dof) {
149247c6ae99SBarry Smith         idx[count++] = (left + i*(Xe-Xs) + j);
149347c6ae99SBarry Smith       }
149447c6ae99SBarry Smith     }
149547c6ae99SBarry Smith     for (i=0; i<count; i++) idx[i] = idx[i]/dof;
149647c6ae99SBarry Smith     ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr);
149747c6ae99SBarry Smith   }
149847c6ae99SBarry Smith 
149947c6ae99SBarry Smith 
150047c6ae99SBarry Smith   /* determine who lies on each side of us stored in    n6 n7 n8
150147c6ae99SBarry Smith                                                         n3    n5
150247c6ae99SBarry Smith                                                         n0 n1 n2
150347c6ae99SBarry Smith   */
150447c6ae99SBarry Smith 
150547c6ae99SBarry Smith   /* Assume the Non-Periodic Case */
150647c6ae99SBarry Smith   n1 = rank - m;
150747c6ae99SBarry Smith   if (rank % m) {
150847c6ae99SBarry Smith     n0 = n1 - 1;
150947c6ae99SBarry Smith   } else {
151047c6ae99SBarry Smith     n0 = -1;
151147c6ae99SBarry Smith   }
151247c6ae99SBarry Smith   if ((rank+1) % m) {
151347c6ae99SBarry Smith     n2 = n1 + 1;
151447c6ae99SBarry Smith     n5 = rank + 1;
151547c6ae99SBarry Smith     n8 = rank + m + 1; if (n8 >= m*n) n8 = -1;
151647c6ae99SBarry Smith   } else {
151747c6ae99SBarry Smith     n2 = -1; n5 = -1; n8 = -1;
151847c6ae99SBarry Smith   }
151947c6ae99SBarry Smith   if (rank % m) {
152047c6ae99SBarry Smith     n3 = rank - 1;
152147c6ae99SBarry Smith     n6 = n3 + m; if (n6 >= m*n) n6 = -1;
152247c6ae99SBarry Smith   } else {
152347c6ae99SBarry Smith     n3 = -1; n6 = -1;
152447c6ae99SBarry Smith   }
152547c6ae99SBarry Smith   n7 = rank + m; if (n7 >= m*n) n7 = -1;
152647c6ae99SBarry Smith 
152747c6ae99SBarry Smith 
152847c6ae99SBarry Smith   /* Modify for Periodic Cases */
1529aa219208SBarry Smith   if (wrap == DMDA_YPERIODIC) {  /* Handle Top and Bottom Sides */
153047c6ae99SBarry Smith     if (n1 < 0) n1 = rank + m * (n-1);
153147c6ae99SBarry Smith     if (n7 < 0) n7 = rank - m * (n-1);
153247c6ae99SBarry Smith     if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1;
153347c6ae99SBarry Smith     if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1;
153447c6ae99SBarry Smith     if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1;
153547c6ae99SBarry Smith     if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1;
1536aa219208SBarry Smith   } else if (wrap == DMDA_XPERIODIC) { /* Handle Left and Right Sides */
153747c6ae99SBarry Smith     if (n3 < 0) n3 = rank + (m-1);
153847c6ae99SBarry Smith     if (n5 < 0) n5 = rank - (m-1);
153947c6ae99SBarry Smith     if ((n1 >= 0) && (n0 < 0)) n0 = rank-1;
154047c6ae99SBarry Smith     if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1;
154147c6ae99SBarry Smith     if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1;
154247c6ae99SBarry Smith     if ((n7 >= 0) && (n8 < 0)) n8 = rank+1;
1543aa219208SBarry Smith   } else if (wrap == DMDA_XYPERIODIC) {
154447c6ae99SBarry Smith 
154547c6ae99SBarry Smith     /* Handle all four corners */
154647c6ae99SBarry Smith     if ((n6 < 0) && (n7 < 0) && (n3 < 0)) n6 = m-1;
154747c6ae99SBarry Smith     if ((n8 < 0) && (n7 < 0) && (n5 < 0)) n8 = 0;
154847c6ae99SBarry Smith     if ((n2 < 0) && (n5 < 0) && (n1 < 0)) n2 = size-m;
154947c6ae99SBarry Smith     if ((n0 < 0) && (n3 < 0) && (n1 < 0)) n0 = size-1;
155047c6ae99SBarry Smith 
155147c6ae99SBarry Smith     /* Handle Top and Bottom Sides */
155247c6ae99SBarry Smith     if (n1 < 0) n1 = rank + m * (n-1);
155347c6ae99SBarry Smith     if (n7 < 0) n7 = rank - m * (n-1);
155447c6ae99SBarry Smith     if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1;
155547c6ae99SBarry Smith     if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1;
155647c6ae99SBarry Smith     if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1;
155747c6ae99SBarry Smith     if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1;
155847c6ae99SBarry Smith 
155947c6ae99SBarry Smith     /* Handle Left and Right Sides */
156047c6ae99SBarry Smith     if (n3 < 0) n3 = rank + (m-1);
156147c6ae99SBarry Smith     if (n5 < 0) n5 = rank - (m-1);
156247c6ae99SBarry Smith     if ((n1 >= 0) && (n0 < 0)) n0 = rank-1;
156347c6ae99SBarry Smith     if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1;
156447c6ae99SBarry Smith     if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1;
156547c6ae99SBarry Smith     if ((n7 >= 0) && (n8 < 0)) n8 = rank+1;
156647c6ae99SBarry Smith   }
156747c6ae99SBarry Smith   ierr = PetscMalloc(9*sizeof(PetscInt),&dd->neighbors);CHKERRQ(ierr);
156847c6ae99SBarry Smith   dd->neighbors[0] = n0;
156947c6ae99SBarry Smith   dd->neighbors[1] = n1;
157047c6ae99SBarry Smith   dd->neighbors[2] = n2;
157147c6ae99SBarry Smith   dd->neighbors[3] = n3;
157247c6ae99SBarry Smith   dd->neighbors[4] = rank;
157347c6ae99SBarry Smith   dd->neighbors[5] = n5;
157447c6ae99SBarry Smith   dd->neighbors[6] = n6;
157547c6ae99SBarry Smith   dd->neighbors[7] = n7;
157647c6ae99SBarry Smith   dd->neighbors[8] = n8;
157747c6ae99SBarry Smith 
1578aa219208SBarry Smith   if (stencil_type == DMDA_STENCIL_STAR) {
157947c6ae99SBarry Smith     /* save corner processor numbers */
158047c6ae99SBarry Smith     sn0 = n0; sn2 = n2; sn6 = n6; sn8 = n8;
158147c6ae99SBarry Smith     n0 = n2 = n6 = n8 = -1;
158247c6ae99SBarry Smith   }
158347c6ae99SBarry Smith 
158447c6ae99SBarry Smith   ierr = PetscMalloc((x+2*s_x)*(y+2*s_y)*sizeof(PetscInt),&idx);CHKERRQ(ierr);
158547c6ae99SBarry Smith   ierr = PetscLogObjectMemory(da,(x+2*s_x)*(y+2*s_y)*sizeof(PetscInt));CHKERRQ(ierr);
158647c6ae99SBarry Smith   nn = 0;
158747c6ae99SBarry Smith 
158847c6ae99SBarry Smith   xbase = bases[rank];
158947c6ae99SBarry Smith   for (i=1; i<=s_y; i++) {
159047c6ae99SBarry Smith     if (n0 >= 0) { /* left below */
159147c6ae99SBarry Smith       x_t = lx[n0 % m]*dof;
159247c6ae99SBarry Smith       y_t = ly[(n0/m)];
159347c6ae99SBarry Smith       s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x;
159447c6ae99SBarry Smith       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
159547c6ae99SBarry Smith     }
159647c6ae99SBarry Smith     if (n1 >= 0) { /* directly below */
159747c6ae99SBarry Smith       x_t = x;
159847c6ae99SBarry Smith       y_t = ly[(n1/m)];
159947c6ae99SBarry Smith       s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t;
160047c6ae99SBarry Smith       for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
160147c6ae99SBarry Smith     }
160247c6ae99SBarry Smith     if (n2 >= 0) { /* right below */
160347c6ae99SBarry Smith       x_t = lx[n2 % m]*dof;
160447c6ae99SBarry Smith       y_t = ly[(n2/m)];
160547c6ae99SBarry Smith       s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t;
160647c6ae99SBarry Smith       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
160747c6ae99SBarry Smith     }
160847c6ae99SBarry Smith   }
160947c6ae99SBarry Smith 
161047c6ae99SBarry Smith   for (i=0; i<y; i++) {
161147c6ae99SBarry Smith     if (n3 >= 0) { /* directly left */
161247c6ae99SBarry Smith       x_t = lx[n3 % m]*dof;
161347c6ae99SBarry Smith       /* y_t = y; */
161447c6ae99SBarry Smith       s_t = bases[n3] + (i+1)*x_t - s_x;
161547c6ae99SBarry Smith       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
161647c6ae99SBarry Smith     }
161747c6ae99SBarry Smith 
161847c6ae99SBarry Smith     for (j=0; j<x; j++) { idx[nn++] = xbase++; } /* interior */
161947c6ae99SBarry Smith 
162047c6ae99SBarry Smith     if (n5 >= 0) { /* directly right */
162147c6ae99SBarry Smith       x_t = lx[n5 % m]*dof;
162247c6ae99SBarry Smith       /* y_t = y; */
162347c6ae99SBarry Smith       s_t = bases[n5] + (i)*x_t;
162447c6ae99SBarry Smith       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
162547c6ae99SBarry Smith     }
162647c6ae99SBarry Smith   }
162747c6ae99SBarry Smith 
162847c6ae99SBarry Smith   for (i=1; i<=s_y; i++) {
162947c6ae99SBarry Smith     if (n6 >= 0) { /* left above */
163047c6ae99SBarry Smith       x_t = lx[n6 % m]*dof;
163147c6ae99SBarry Smith       /* y_t = ly[(n6/m)]; */
163247c6ae99SBarry Smith       s_t = bases[n6] + (i)*x_t - s_x;
163347c6ae99SBarry Smith       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
163447c6ae99SBarry Smith     }
163547c6ae99SBarry Smith     if (n7 >= 0) { /* directly above */
163647c6ae99SBarry Smith       x_t = x;
163747c6ae99SBarry Smith       /* y_t = ly[(n7/m)]; */
163847c6ae99SBarry Smith       s_t = bases[n7] + (i-1)*x_t;
163947c6ae99SBarry Smith       for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
164047c6ae99SBarry Smith     }
164147c6ae99SBarry Smith     if (n8 >= 0) { /* right above */
164247c6ae99SBarry Smith       x_t = lx[n8 % m]*dof;
164347c6ae99SBarry Smith       /* y_t = ly[(n8/m)]; */
164447c6ae99SBarry Smith       s_t = bases[n8] + (i-1)*x_t;
164547c6ae99SBarry Smith       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
164647c6ae99SBarry Smith     }
164747c6ae99SBarry Smith   }
164847c6ae99SBarry Smith 
164947c6ae99SBarry Smith   base = bases[rank];
165047c6ae99SBarry Smith   {
165147c6ae99SBarry Smith     PetscInt nnn = nn/dof,*iidx;
165247c6ae99SBarry Smith     ierr = PetscMalloc(nnn*sizeof(PetscInt),&iidx);CHKERRQ(ierr);
165347c6ae99SBarry Smith     for (i=0; i<nnn; i++) {
165447c6ae99SBarry Smith       iidx[i] = idx[dof*i]/dof;
165547c6ae99SBarry Smith     }
165647c6ae99SBarry Smith     ierr = ISCreateBlock(comm,dof,nnn,iidx,PETSC_OWN_POINTER,&from);CHKERRQ(ierr);
165747c6ae99SBarry Smith   }
165847c6ae99SBarry Smith   ierr = VecScatterCreate(global,from,local,to,&gtol);CHKERRQ(ierr);
165947c6ae99SBarry Smith   ierr = PetscLogObjectParent(da,gtol);CHKERRQ(ierr);
166047c6ae99SBarry Smith   ierr = ISDestroy(to);CHKERRQ(ierr);
166147c6ae99SBarry Smith   ierr = ISDestroy(from);CHKERRQ(ierr);
166247c6ae99SBarry Smith 
1663aa219208SBarry Smith   if (stencil_type == DMDA_STENCIL_STAR) {
166447c6ae99SBarry Smith     /*
166547c6ae99SBarry Smith         Recompute the local to global mappings, this time keeping the
166647c6ae99SBarry Smith       information about the cross corner processor numbers.
166747c6ae99SBarry Smith     */
166847c6ae99SBarry Smith     n0 = sn0; n2 = sn2; n6 = sn6; n8 = sn8;
166947c6ae99SBarry Smith     nn = 0;
167047c6ae99SBarry Smith     xbase = bases[rank];
167147c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
167247c6ae99SBarry Smith       if (n0 >= 0) { /* left below */
167347c6ae99SBarry Smith         x_t = lx[n0 % m]*dof;
167447c6ae99SBarry Smith         y_t = ly[(n0/m)];
167547c6ae99SBarry Smith         s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x;
167647c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
167747c6ae99SBarry Smith       }
167847c6ae99SBarry Smith       if (n1 >= 0) { /* directly below */
167947c6ae99SBarry Smith         x_t = x;
168047c6ae99SBarry Smith         y_t = ly[(n1/m)];
168147c6ae99SBarry Smith         s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t;
168247c6ae99SBarry Smith         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
168347c6ae99SBarry Smith       }
168447c6ae99SBarry Smith       if (n2 >= 0) { /* right below */
168547c6ae99SBarry Smith         x_t = lx[n2 % m]*dof;
168647c6ae99SBarry Smith         y_t = ly[(n2/m)];
168747c6ae99SBarry Smith         s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t;
168847c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
168947c6ae99SBarry Smith       }
169047c6ae99SBarry Smith     }
169147c6ae99SBarry Smith 
169247c6ae99SBarry Smith     for (i=0; i<y; i++) {
169347c6ae99SBarry Smith       if (n3 >= 0) { /* directly left */
169447c6ae99SBarry Smith         x_t = lx[n3 % m]*dof;
169547c6ae99SBarry Smith         /* y_t = y; */
169647c6ae99SBarry Smith         s_t = bases[n3] + (i+1)*x_t - s_x;
169747c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
169847c6ae99SBarry Smith       }
169947c6ae99SBarry Smith 
170047c6ae99SBarry Smith       for (j=0; j<x; j++) { idx[nn++] = xbase++; } /* interior */
170147c6ae99SBarry Smith 
170247c6ae99SBarry Smith       if (n5 >= 0) { /* directly right */
170347c6ae99SBarry Smith         x_t = lx[n5 % m]*dof;
170447c6ae99SBarry Smith         /* y_t = y; */
170547c6ae99SBarry Smith         s_t = bases[n5] + (i)*x_t;
170647c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
170747c6ae99SBarry Smith       }
170847c6ae99SBarry Smith     }
170947c6ae99SBarry Smith 
171047c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
171147c6ae99SBarry Smith       if (n6 >= 0) { /* left above */
171247c6ae99SBarry Smith         x_t = lx[n6 % m]*dof;
171347c6ae99SBarry Smith         /* y_t = ly[(n6/m)]; */
171447c6ae99SBarry Smith         s_t = bases[n6] + (i)*x_t - s_x;
171547c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
171647c6ae99SBarry Smith       }
171747c6ae99SBarry Smith       if (n7 >= 0) { /* directly above */
171847c6ae99SBarry Smith         x_t = x;
171947c6ae99SBarry Smith         /* y_t = ly[(n7/m)]; */
172047c6ae99SBarry Smith         s_t = bases[n7] + (i-1)*x_t;
172147c6ae99SBarry Smith         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
172247c6ae99SBarry Smith       }
172347c6ae99SBarry Smith       if (n8 >= 0) { /* right above */
172447c6ae99SBarry Smith         x_t = lx[n8 % m]*dof;
172547c6ae99SBarry Smith         /* y_t = ly[(n8/m)]; */
172647c6ae99SBarry Smith         s_t = bases[n8] + (i-1)*x_t;
172747c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
172847c6ae99SBarry Smith       }
172947c6ae99SBarry Smith     }
173047c6ae99SBarry Smith   }
173147c6ae99SBarry Smith   ierr = PetscFree2(bases,ldims);CHKERRQ(ierr);
173247c6ae99SBarry Smith 
173347c6ae99SBarry Smith   dd->m  = m;  dd->n  = n;
173447c6ae99SBarry Smith   dd->xs = xs; dd->xe = xe; dd->ys = ys; dd->ye = ye; dd->zs = 0; dd->ze = 1;
173547c6ae99SBarry Smith   dd->Xs = Xs; dd->Xe = Xe; dd->Ys = Ys; dd->Ye = Ye; dd->Zs = 0; dd->Ze = 1;
173647c6ae99SBarry Smith 
173747c6ae99SBarry Smith   ierr = VecDestroy(local);CHKERRQ(ierr);
173847c6ae99SBarry Smith   ierr = VecDestroy(global);CHKERRQ(ierr);
173947c6ae99SBarry Smith 
174047c6ae99SBarry Smith   dd->gtol      = gtol;
174147c6ae99SBarry Smith   dd->ltog      = ltog;
174247c6ae99SBarry Smith   dd->idx       = idx;
174347c6ae99SBarry Smith   dd->Nl        = nn;
174447c6ae99SBarry Smith   dd->base      = base;
17459a42bb27SBarry Smith   da->ops->view = DMView_DA_2d;
174647c6ae99SBarry Smith 
174747c6ae99SBarry Smith   /*
174847c6ae99SBarry Smith      Set the local to global ordering in the global vector, this allows use
174947c6ae99SBarry Smith      of VecSetValuesLocal().
175047c6ae99SBarry Smith   */
1751*1411c6eeSJed Brown   ierr = ISLocalToGlobalMappingCreate(comm,nn,idx,PETSC_OWN_POINTER,&da->ltogmap);CHKERRQ(ierr);
1752*1411c6eeSJed Brown   ierr = ISLocalToGlobalMappingBlock(da->ltogmap,dd->w,&da->ltogmapb);CHKERRQ(ierr);
1753*1411c6eeSJed Brown   ierr = PetscLogObjectParent(da,da->ltogmap);CHKERRQ(ierr);
175447c6ae99SBarry Smith 
175547c6ae99SBarry Smith   dd->ltol = PETSC_NULL;
175647c6ae99SBarry Smith   dd->ao   = PETSC_NULL;
175747c6ae99SBarry Smith 
175847c6ae99SBarry Smith   PetscFunctionReturn(0);
175947c6ae99SBarry Smith }
176047c6ae99SBarry Smith 
176147c6ae99SBarry Smith #undef __FUNCT__
1762aa219208SBarry Smith #define __FUNCT__ "DMDACreate2d"
176347c6ae99SBarry Smith /*@C
1764aa219208SBarry Smith    DMDACreate2d -  Creates an object that will manage the communication of  two-dimensional
176547c6ae99SBarry Smith    regular array data that is distributed across some processors.
176647c6ae99SBarry Smith 
176747c6ae99SBarry Smith    Collective on MPI_Comm
176847c6ae99SBarry Smith 
176947c6ae99SBarry Smith    Input Parameters:
177047c6ae99SBarry Smith +  comm - MPI communicator
177147c6ae99SBarry Smith .  wrap - type of periodicity should the array have.
1772aa219208SBarry Smith          Use one of DMDA_NONPERIODIC, DMDA_XPERIODIC, DMDA_YPERIODIC, or DMDA_XYPERIODIC.
1773aa219208SBarry Smith .  stencil_type - stencil type.  Use either DMDA_STENCIL_BOX or DMDA_STENCIL_STAR.
177447c6ae99SBarry 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
177547c6ae99SBarry Smith             from the command line with -da_grid_x <M> -da_grid_y <N>)
177647c6ae99SBarry Smith .  m,n - corresponding number of processors in each dimension
177747c6ae99SBarry Smith          (or PETSC_DECIDE to have calculated)
177847c6ae99SBarry Smith .  dof - number of degrees of freedom per node
177947c6ae99SBarry Smith .  s - stencil width
178047c6ae99SBarry Smith -  lx, ly - arrays containing the number of nodes in each cell along
178147c6ae99SBarry Smith            the x and y coordinates, or PETSC_NULL. If non-null, these
178247c6ae99SBarry Smith            must be of length as m and n, and the corresponding
178347c6ae99SBarry Smith            m and n cannot be PETSC_DECIDE. The sum of the lx[] entries
178447c6ae99SBarry Smith            must be M, and the sum of the ly[] entries must be N.
178547c6ae99SBarry Smith 
178647c6ae99SBarry Smith    Output Parameter:
178747c6ae99SBarry Smith .  da - the resulting distributed array object
178847c6ae99SBarry Smith 
178947c6ae99SBarry Smith    Options Database Key:
1790aa219208SBarry Smith +  -da_view - Calls DMView() at the conclusion of DMDACreate2d()
179147c6ae99SBarry Smith .  -da_grid_x <nx> - number of grid points in x direction, if M < 0
179247c6ae99SBarry Smith .  -da_grid_y <ny> - number of grid points in y direction, if N < 0
179347c6ae99SBarry Smith .  -da_processors_x <nx> - number of processors in x direction
179447c6ae99SBarry Smith .  -da_processors_y <ny> - number of processors in y direction
179547c6ae99SBarry Smith .  -da_refine_x - refinement ratio in x direction
179647c6ae99SBarry Smith -  -da_refine_y - refinement ratio in y direction
179747c6ae99SBarry Smith 
179847c6ae99SBarry Smith    Level: beginner
179947c6ae99SBarry Smith 
180047c6ae99SBarry Smith    Notes:
1801aa219208SBarry Smith    The stencil type DMDA_STENCIL_STAR with width 1 corresponds to the
1802aa219208SBarry Smith    standard 5-pt stencil, while DMDA_STENCIL_BOX with width 1 denotes
180347c6ae99SBarry Smith    the standard 9-pt stencil.
180447c6ae99SBarry Smith 
1805aa219208SBarry Smith    The array data itself is NOT stored in the DMDA, it is stored in Vec objects;
1806564755cdSBarry Smith    The appropriate vector objects can be obtained with calls to DMCreateGlobalVector()
1807564755cdSBarry Smith    and DMCreateLocalVector() and calls to VecDuplicate() if more are needed.
180847c6ae99SBarry Smith 
180947c6ae99SBarry Smith .keywords: distributed array, create, two-dimensional
181047c6ae99SBarry Smith 
1811aa219208SBarry Smith .seealso: DMDestroy(), DMView(), DMDACreate1d(), DMDACreate3d(), DMGlobalToLocalBegin(), DMDAGetRefinementFactor(),
1812aa219208SBarry Smith           DMGlobalToLocalEnd(), DMLocalToGlobalBegin(), DMDALocalToLocalBegin(), DMDALocalToLocalEnd(), DMDASetRefinementFactor(),
1813aa219208SBarry Smith           DMDAGetInfo(), DMCreateGlobalVector(), DMCreateLocalVector(), DMDACreateNaturalVector(), DMDALoad(), DMDAGetOwnershipRanges()
181447c6ae99SBarry Smith 
181547c6ae99SBarry Smith @*/
1816aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDACreate2d(MPI_Comm comm,DMDAPeriodicType wrap,DMDAStencilType stencil_type,
18179a42bb27SBarry Smith                           PetscInt M,PetscInt N,PetscInt m,PetscInt n,PetscInt dof,PetscInt s,const PetscInt lx[],const PetscInt ly[],DM *da)
181847c6ae99SBarry Smith {
181947c6ae99SBarry Smith   PetscErrorCode ierr;
182047c6ae99SBarry Smith 
182147c6ae99SBarry Smith   PetscFunctionBegin;
1822aa219208SBarry Smith   ierr = DMDACreate(comm, da);CHKERRQ(ierr);
1823aa219208SBarry Smith   ierr = DMDASetDim(*da, 2);CHKERRQ(ierr);
1824aa219208SBarry Smith   ierr = DMDASetSizes(*da, M, N, 1);CHKERRQ(ierr);
1825aa219208SBarry Smith   ierr = DMDASetNumProcs(*da, m, n, PETSC_DECIDE);CHKERRQ(ierr);
1826aa219208SBarry Smith   ierr = DMDASetPeriodicity(*da, wrap);CHKERRQ(ierr);
1827aa219208SBarry Smith   ierr = DMDASetDof(*da, dof);CHKERRQ(ierr);
1828aa219208SBarry Smith   ierr = DMDASetStencilType(*da, stencil_type);CHKERRQ(ierr);
1829aa219208SBarry Smith   ierr = DMDASetStencilWidth(*da, s);CHKERRQ(ierr);
1830aa219208SBarry Smith   ierr = DMDASetOwnershipRanges(*da, lx, ly, PETSC_NULL);CHKERRQ(ierr);
183147c6ae99SBarry Smith   /* This violates the behavior for other classes, but right now users expect negative dimensions to be handled this way */
18329a42bb27SBarry Smith   ierr = DMSetFromOptions(*da);CHKERRQ(ierr);
18339a42bb27SBarry Smith   ierr = DMSetUp(*da);CHKERRQ(ierr);
18347242296bSJed Brown   ierr = DMView_DA_Private(*da);CHKERRQ(ierr);
183547c6ae99SBarry Smith   PetscFunctionReturn(0);
183647c6ae99SBarry Smith }
1837