xref: /petsc/src/dm/impls/da/da2.c (revision 3da9ae135d9cba4ac09f57c869ee1dab666150e3)
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);
37*3da9ae13SJed Brown     } else {
38*3da9ae13SJed 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__
1715dff015aSBarry Smith #define __FUNCT__ "DMGetElements_DA_2d_P1"
1729a42bb27SBarry Smith PetscErrorCode DMGetElements_DA_2d_P1(DM da,PetscInt *n,const PetscInt *e[])
17347c6ae99SBarry Smith {
17447c6ae99SBarry Smith   PetscErrorCode ierr;
17547c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
17647c6ae99SBarry Smith   PetscInt       i,j,cnt,xs,xe = dd->xe,ys,ye = dd->ye,Xs = dd->Xs, Xe = dd->Xe, Ys = dd->Ys;
17747c6ae99SBarry Smith 
17847c6ae99SBarry Smith   PetscFunctionBegin;
17947c6ae99SBarry Smith   if (!dd->e) {
18047c6ae99SBarry Smith     if (dd->xs == Xs) xs = dd->xs; else xs = dd->xs - 1;
18147c6ae99SBarry Smith     if (dd->ys == Ys) ys = dd->ys; else ys = dd->ys - 1;
18247c6ae99SBarry Smith     dd->ne = 2*(xe - xs - 1)*(ye - ys - 1);
18347c6ae99SBarry Smith     ierr   = PetscMalloc((1 + 3*dd->ne)*sizeof(PetscInt),&dd->e);CHKERRQ(ierr);
18447c6ae99SBarry Smith     cnt    = 0;
18547c6ae99SBarry Smith     for (j=ys; j<ye-1; j++) {
18647c6ae99SBarry Smith       for (i=xs; i<xe-1; i++) {
18747c6ae99SBarry Smith         dd->e[cnt]   = i - Xs + (j - Ys)*(Xe - Xs);
18847c6ae99SBarry Smith         dd->e[cnt+1] = i - Xs + 1 + (j - Ys)*(Xe - Xs);
18947c6ae99SBarry Smith         dd->e[cnt+2] = i - Xs + (j - Ys + 1)*(Xe - Xs);
19047c6ae99SBarry Smith 
19147c6ae99SBarry Smith         dd->e[cnt+3] = i - Xs + 1 + (j - Ys + 1)*(Xe - Xs);
19247c6ae99SBarry Smith         dd->e[cnt+4] = i - Xs + (j - Ys + 1)*(Xe - Xs);
19347c6ae99SBarry Smith         dd->e[cnt+5] = i - Xs + 1 + (j - Ys)*(Xe - Xs);
19447c6ae99SBarry Smith         cnt += 6;
19547c6ae99SBarry Smith       }
19647c6ae99SBarry Smith     }
19747c6ae99SBarry Smith   }
19847c6ae99SBarry Smith   *n = dd->ne;
19947c6ae99SBarry Smith   *e = dd->e;
20047c6ae99SBarry Smith   PetscFunctionReturn(0);
20147c6ae99SBarry Smith }
20247c6ae99SBarry Smith 
20347c6ae99SBarry Smith #undef __FUNCT__
204aa219208SBarry Smith #define __FUNCT__ "DMDASetLocalFunction"
20547c6ae99SBarry Smith /*@C
206aa219208SBarry Smith        DMDASetLocalFunction - Caches in a DM a local function.
20747c6ae99SBarry Smith 
208aa219208SBarry Smith    Logically Collective on DMDA
20947c6ae99SBarry Smith 
21047c6ae99SBarry Smith    Input Parameter:
21147c6ae99SBarry Smith +  da - initial distributed array
21247c6ae99SBarry Smith -  lf - the local function
21347c6ae99SBarry Smith 
21447c6ae99SBarry Smith    Level: intermediate
21547c6ae99SBarry Smith 
21647c6ae99SBarry Smith    Notes: The routine SNESDAFormFunction() uses this the cached function to evaluate the user provided function.
21747c6ae99SBarry Smith 
21847c6ae99SBarry Smith .keywords:  distributed array, refine
21947c6ae99SBarry Smith 
220aa219208SBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalFunctioni()
22147c6ae99SBarry Smith @*/
222aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDASetLocalFunction(DM da,DMDALocalFunction1 lf)
22347c6ae99SBarry Smith {
22447c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
22547c6ae99SBarry Smith   PetscFunctionBegin;
22647c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
22747c6ae99SBarry Smith   dd->lf    = lf;
22847c6ae99SBarry Smith   PetscFunctionReturn(0);
22947c6ae99SBarry Smith }
23047c6ae99SBarry Smith 
23147c6ae99SBarry Smith #undef __FUNCT__
232aa219208SBarry Smith #define __FUNCT__ "DMDASetLocalFunctioni"
23347c6ae99SBarry Smith /*@C
234aa219208SBarry Smith        DMDASetLocalFunctioni - Caches in a DM a local function that evaluates a single component
23547c6ae99SBarry Smith 
236aa219208SBarry Smith    Logically Collective on DMDA
23747c6ae99SBarry Smith 
23847c6ae99SBarry Smith    Input Parameter:
23947c6ae99SBarry Smith +  da - initial distributed array
24047c6ae99SBarry Smith -  lfi - the local function
24147c6ae99SBarry Smith 
24247c6ae99SBarry Smith    Level: intermediate
24347c6ae99SBarry Smith 
24447c6ae99SBarry Smith .keywords:  distributed array, refine
24547c6ae99SBarry Smith 
246aa219208SBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalFunction()
24747c6ae99SBarry Smith @*/
248aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDASetLocalFunctioni(DM da,PetscErrorCode (*lfi)(DMDALocalInfo*,MatStencil*,void*,PetscScalar*,void*))
24947c6ae99SBarry Smith {
25047c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
25147c6ae99SBarry Smith   PetscFunctionBegin;
25247c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
25347c6ae99SBarry Smith   dd->lfi = lfi;
25447c6ae99SBarry Smith   PetscFunctionReturn(0);
25547c6ae99SBarry Smith }
25647c6ae99SBarry Smith 
25747c6ae99SBarry Smith #undef __FUNCT__
258aa219208SBarry Smith #define __FUNCT__ "DMDASetLocalFunctionib"
25947c6ae99SBarry Smith /*@C
260aa219208SBarry Smith        DMDASetLocalFunctionib - Caches in a DM a block local function that evaluates a single component
26147c6ae99SBarry Smith 
262aa219208SBarry Smith    Logically Collective on DMDA
26347c6ae99SBarry Smith 
26447c6ae99SBarry Smith    Input Parameter:
26547c6ae99SBarry Smith +  da - initial distributed array
26647c6ae99SBarry Smith -  lfi - the local function
26747c6ae99SBarry Smith 
26847c6ae99SBarry Smith    Level: intermediate
26947c6ae99SBarry Smith 
27047c6ae99SBarry Smith .keywords:  distributed array, refine
27147c6ae99SBarry Smith 
272aa219208SBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalFunction()
27347c6ae99SBarry Smith @*/
274aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDASetLocalFunctionib(DM da,PetscErrorCode (*lfi)(DMDALocalInfo*,MatStencil*,void*,PetscScalar*,void*))
27547c6ae99SBarry Smith {
27647c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
27747c6ae99SBarry Smith   PetscFunctionBegin;
27847c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
27947c6ae99SBarry Smith   dd->lfib = lfi;
28047c6ae99SBarry Smith   PetscFunctionReturn(0);
28147c6ae99SBarry Smith }
28247c6ae99SBarry Smith 
28347c6ae99SBarry Smith #undef __FUNCT__
284aa219208SBarry Smith #define __FUNCT__ "DMDASetLocalAdicFunction_Private"
285aa219208SBarry Smith PetscErrorCode DMDASetLocalAdicFunction_Private(DM da,DMDALocalFunction1 ad_lf)
28647c6ae99SBarry Smith {
28747c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
28847c6ae99SBarry Smith   PetscFunctionBegin;
28947c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
29047c6ae99SBarry Smith   dd->adic_lf = ad_lf;
29147c6ae99SBarry Smith   PetscFunctionReturn(0);
29247c6ae99SBarry Smith }
29347c6ae99SBarry Smith 
29447c6ae99SBarry Smith /*MC
295aa219208SBarry Smith        DMDASetLocalAdicFunctioni - Caches in a DM a local functioni computed by ADIC/ADIFOR
29647c6ae99SBarry Smith 
29747c6ae99SBarry Smith    Synopsis:
298aa219208SBarry Smith    PetscErrorCode DMDASetLocalAdicFunctioni(DM da,PetscInt (ad_lf*)(DMDALocalInfo*,MatStencil*,void*,void*,void*)
29947c6ae99SBarry Smith 
300aa219208SBarry Smith    Logically Collective on DMDA
30147c6ae99SBarry Smith 
30247c6ae99SBarry Smith    Input Parameter:
30347c6ae99SBarry Smith +  da - initial distributed array
30447c6ae99SBarry Smith -  ad_lfi - the local function as computed by ADIC/ADIFOR
30547c6ae99SBarry Smith 
30647c6ae99SBarry Smith    Level: intermediate
30747c6ae99SBarry Smith 
30847c6ae99SBarry Smith .keywords:  distributed array, refine
30947c6ae99SBarry Smith 
310aa219208SBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalFunction(),
311aa219208SBarry Smith           DMDASetLocalJacobian(), DMDASetLocalFunctioni()
31247c6ae99SBarry Smith M*/
31347c6ae99SBarry Smith 
31447c6ae99SBarry Smith #undef __FUNCT__
315aa219208SBarry Smith #define __FUNCT__ "DMDASetLocalAdicFunctioni_Private"
316aa219208SBarry Smith PetscErrorCode DMDASetLocalAdicFunctioni_Private(DM da,PetscErrorCode (*ad_lfi)(DMDALocalInfo*,MatStencil*,void*,void*,void*))
31747c6ae99SBarry Smith {
31847c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
31947c6ae99SBarry Smith   PetscFunctionBegin;
32047c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
32147c6ae99SBarry Smith   dd->adic_lfi = ad_lfi;
32247c6ae99SBarry Smith   PetscFunctionReturn(0);
32347c6ae99SBarry Smith }
32447c6ae99SBarry Smith 
32547c6ae99SBarry Smith /*MC
326aa219208SBarry Smith        DMDASetLocalAdicMFFunctioni - Caches in a DM a local functioni computed by ADIC/ADIFOR
32747c6ae99SBarry Smith 
32847c6ae99SBarry Smith    Synopsis:
329aa219208SBarry Smith    PetscErrorCode  DMDASetLocalAdicFunctioni(DM da,int (ad_lf*)(DMDALocalInfo*,MatStencil*,void*,void*,void*)
33047c6ae99SBarry Smith 
331aa219208SBarry Smith    Logically Collective on DMDA
33247c6ae99SBarry Smith 
33347c6ae99SBarry Smith    Input Parameter:
33447c6ae99SBarry Smith +  da - initial distributed array
33547c6ae99SBarry Smith -  admf_lfi - the local matrix-free function as computed by ADIC/ADIFOR
33647c6ae99SBarry Smith 
33747c6ae99SBarry Smith    Level: intermediate
33847c6ae99SBarry Smith 
33947c6ae99SBarry Smith .keywords:  distributed array, refine
34047c6ae99SBarry Smith 
341aa219208SBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalFunction(),
342aa219208SBarry Smith           DMDASetLocalJacobian(), DMDASetLocalFunctioni()
34347c6ae99SBarry Smith M*/
34447c6ae99SBarry Smith 
34547c6ae99SBarry Smith #undef __FUNCT__
346aa219208SBarry Smith #define __FUNCT__ "DMDASetLocalAdicMFFunctioni_Private"
347aa219208SBarry Smith PetscErrorCode DMDASetLocalAdicMFFunctioni_Private(DM da,PetscErrorCode (*admf_lfi)(DMDALocalInfo*,MatStencil*,void*,void*,void*))
34847c6ae99SBarry Smith {
34947c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
35047c6ae99SBarry Smith   PetscFunctionBegin;
35147c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
35247c6ae99SBarry Smith   dd->adicmf_lfi = admf_lfi;
35347c6ae99SBarry Smith   PetscFunctionReturn(0);
35447c6ae99SBarry Smith }
35547c6ae99SBarry Smith 
35647c6ae99SBarry Smith /*MC
357aa219208SBarry Smith        DMDASetLocalAdicFunctionib - Caches in a DM a block local functioni computed by ADIC/ADIFOR
35847c6ae99SBarry Smith 
35947c6ae99SBarry Smith    Synopsis:
360aa219208SBarry Smith    PetscErrorCode DMDASetLocalAdicFunctionib(DM da,PetscInt (ad_lf*)(DMDALocalInfo*,MatStencil*,void*,void*,void*)
36147c6ae99SBarry Smith 
362aa219208SBarry Smith    Logically Collective on DMDA
36347c6ae99SBarry Smith 
36447c6ae99SBarry Smith    Input Parameter:
36547c6ae99SBarry Smith +  da - initial distributed array
36647c6ae99SBarry Smith -  ad_lfi - the local function as computed by ADIC/ADIFOR
36747c6ae99SBarry Smith 
36847c6ae99SBarry Smith    Level: intermediate
36947c6ae99SBarry Smith 
37047c6ae99SBarry Smith .keywords:  distributed array, refine
37147c6ae99SBarry Smith 
372aa219208SBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalFunction(),
373aa219208SBarry Smith           DMDASetLocalJacobian(), DMDASetLocalFunctionib()
37447c6ae99SBarry Smith M*/
37547c6ae99SBarry Smith 
37647c6ae99SBarry Smith #undef __FUNCT__
377aa219208SBarry Smith #define __FUNCT__ "DMDASetLocalAdicFunctionib_Private"
378aa219208SBarry Smith PetscErrorCode DMDASetLocalAdicFunctionib_Private(DM da,PetscErrorCode (*ad_lfi)(DMDALocalInfo*,MatStencil*,void*,void*,void*))
37947c6ae99SBarry Smith {
38047c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
38147c6ae99SBarry Smith   PetscFunctionBegin;
38247c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
38347c6ae99SBarry Smith   dd->adic_lfib = ad_lfi;
38447c6ae99SBarry Smith   PetscFunctionReturn(0);
38547c6ae99SBarry Smith }
38647c6ae99SBarry Smith 
38747c6ae99SBarry Smith /*MC
388aa219208SBarry Smith        DMDASetLocalAdicMFFunctionib - Caches in a DM a block local functioni computed by ADIC/ADIFOR
38947c6ae99SBarry Smith 
39047c6ae99SBarry Smith    Synopsis:
391aa219208SBarry Smith    PetscErrorCode  DMDASetLocalAdicFunctionib(DM da,int (ad_lf*)(DMDALocalInfo*,MatStencil*,void*,void*,void*)
39247c6ae99SBarry Smith 
393aa219208SBarry Smith    Logically Collective on DMDA
39447c6ae99SBarry Smith 
39547c6ae99SBarry Smith    Input Parameter:
39647c6ae99SBarry Smith +  da - initial distributed array
39747c6ae99SBarry Smith -  admf_lfi - the local matrix-free function as computed by ADIC/ADIFOR
39847c6ae99SBarry Smith 
39947c6ae99SBarry Smith    Level: intermediate
40047c6ae99SBarry Smith 
40147c6ae99SBarry Smith .keywords:  distributed array, refine
40247c6ae99SBarry Smith 
403aa219208SBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalFunction(),
404aa219208SBarry Smith           DMDASetLocalJacobian(), DMDASetLocalFunctionib()
40547c6ae99SBarry Smith M*/
40647c6ae99SBarry Smith 
40747c6ae99SBarry Smith #undef __FUNCT__
408aa219208SBarry Smith #define __FUNCT__ "DMDASetLocalAdicMFFunctionib_Private"
409aa219208SBarry Smith PetscErrorCode DMDASetLocalAdicMFFunctionib_Private(DM da,PetscErrorCode (*admf_lfi)(DMDALocalInfo*,MatStencil*,void*,void*,void*))
41047c6ae99SBarry Smith {
41147c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
41247c6ae99SBarry Smith   PetscFunctionBegin;
41347c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
41447c6ae99SBarry Smith   dd->adicmf_lfib = admf_lfi;
41547c6ae99SBarry Smith   PetscFunctionReturn(0);
41647c6ae99SBarry Smith }
41747c6ae99SBarry Smith 
41847c6ae99SBarry Smith /*MC
419aa219208SBarry Smith        DMDASetLocalAdicMFFunction - Caches in a DM a local function computed by ADIC/ADIFOR
42047c6ae99SBarry Smith 
42147c6ae99SBarry Smith    Synopsis:
422aa219208SBarry Smith    PetscErrorCode DMDASetLocalAdicMFFunction(DM da,DMDALocalFunction1 ad_lf)
42347c6ae99SBarry Smith 
424aa219208SBarry Smith    Logically Collective on DMDA
42547c6ae99SBarry Smith 
42647c6ae99SBarry Smith    Input Parameter:
42747c6ae99SBarry Smith +  da - initial distributed array
42847c6ae99SBarry Smith -  ad_lf - the local function as computed by ADIC/ADIFOR
42947c6ae99SBarry Smith 
43047c6ae99SBarry Smith    Level: intermediate
43147c6ae99SBarry Smith 
43247c6ae99SBarry Smith .keywords:  distributed array, refine
43347c6ae99SBarry Smith 
434aa219208SBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalFunction(),
435aa219208SBarry Smith           DMDASetLocalJacobian()
43647c6ae99SBarry Smith M*/
43747c6ae99SBarry Smith 
43847c6ae99SBarry Smith #undef __FUNCT__
439aa219208SBarry Smith #define __FUNCT__ "DMDASetLocalAdicMFFunction_Private"
440aa219208SBarry Smith PetscErrorCode DMDASetLocalAdicMFFunction_Private(DM da,DMDALocalFunction1 ad_lf)
44147c6ae99SBarry Smith {
44247c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
44347c6ae99SBarry Smith   PetscFunctionBegin;
44447c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
44547c6ae99SBarry Smith   dd->adicmf_lf = ad_lf;
44647c6ae99SBarry Smith   PetscFunctionReturn(0);
44747c6ae99SBarry Smith }
44847c6ae99SBarry Smith 
44947c6ae99SBarry Smith /*@C
450aa219208SBarry Smith        DMDASetLocalJacobian - Caches in a DM a local Jacobian computation function
45147c6ae99SBarry Smith 
452aa219208SBarry Smith    Logically Collective on DMDA
45347c6ae99SBarry Smith 
45447c6ae99SBarry Smith 
45547c6ae99SBarry Smith    Input Parameter:
45647c6ae99SBarry Smith +  da - initial distributed array
45747c6ae99SBarry Smith -  lj - the local Jacobian
45847c6ae99SBarry Smith 
45947c6ae99SBarry Smith    Level: intermediate
46047c6ae99SBarry Smith 
46147c6ae99SBarry Smith    Notes: The routine SNESDAFormFunction() uses this the cached function to evaluate the user provided function.
46247c6ae99SBarry Smith 
46347c6ae99SBarry Smith .keywords:  distributed array, refine
46447c6ae99SBarry Smith 
465aa219208SBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalFunction()
46647c6ae99SBarry Smith @*/
46747c6ae99SBarry Smith #undef __FUNCT__
468aa219208SBarry Smith #define __FUNCT__ "DMDASetLocalJacobian"
469aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDASetLocalJacobian(DM da,DMDALocalFunction1 lj)
47047c6ae99SBarry Smith {
47147c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
47247c6ae99SBarry Smith   PetscFunctionBegin;
47347c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
47447c6ae99SBarry Smith   dd->lj    = lj;
47547c6ae99SBarry Smith   PetscFunctionReturn(0);
47647c6ae99SBarry Smith }
47747c6ae99SBarry Smith 
47847c6ae99SBarry Smith #undef __FUNCT__
479aa219208SBarry Smith #define __FUNCT__ "DMDAGetLocalFunction"
48047c6ae99SBarry Smith /*@C
481aa219208SBarry Smith        DMDAGetLocalFunction - Gets from a DM a local function and its ADIC/ADIFOR Jacobian
48247c6ae99SBarry Smith 
48347c6ae99SBarry Smith    Note Collective
48447c6ae99SBarry Smith 
48547c6ae99SBarry Smith    Input Parameter:
48647c6ae99SBarry Smith .  da - initial distributed array
48747c6ae99SBarry Smith 
48847c6ae99SBarry Smith    Output Parameter:
48947c6ae99SBarry Smith .  lf - the local function
49047c6ae99SBarry Smith 
49147c6ae99SBarry Smith    Level: intermediate
49247c6ae99SBarry Smith 
49347c6ae99SBarry Smith .keywords:  distributed array, refine
49447c6ae99SBarry Smith 
495aa219208SBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalJacobian(), DMDASetLocalFunction()
49647c6ae99SBarry Smith @*/
497aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDAGetLocalFunction(DM da,DMDALocalFunction1 *lf)
49847c6ae99SBarry Smith {
49947c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
50047c6ae99SBarry Smith   PetscFunctionBegin;
50147c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
50247c6ae99SBarry Smith   if (lf)       *lf = dd->lf;
50347c6ae99SBarry Smith   PetscFunctionReturn(0);
50447c6ae99SBarry Smith }
50547c6ae99SBarry Smith 
50647c6ae99SBarry Smith #undef __FUNCT__
507aa219208SBarry Smith #define __FUNCT__ "DMDAGetLocalJacobian"
50847c6ae99SBarry Smith /*@C
509aa219208SBarry Smith        DMDAGetLocalJacobian - Gets from a DM a local jacobian
51047c6ae99SBarry Smith 
51147c6ae99SBarry Smith    Not Collective
51247c6ae99SBarry Smith 
51347c6ae99SBarry Smith    Input Parameter:
51447c6ae99SBarry Smith .  da - initial distributed array
51547c6ae99SBarry Smith 
51647c6ae99SBarry Smith    Output Parameter:
51747c6ae99SBarry Smith .  lj - the local jacobian
51847c6ae99SBarry Smith 
51947c6ae99SBarry Smith    Level: intermediate
52047c6ae99SBarry Smith 
52147c6ae99SBarry Smith .keywords:  distributed array, refine
52247c6ae99SBarry Smith 
523aa219208SBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalJacobian()
52447c6ae99SBarry Smith @*/
525aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDAGetLocalJacobian(DM da,DMDALocalFunction1 *lj)
52647c6ae99SBarry Smith {
52747c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
52847c6ae99SBarry Smith   PetscFunctionBegin;
52947c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
53047c6ae99SBarry Smith   if (lj) *lj = dd->lj;
53147c6ae99SBarry Smith   PetscFunctionReturn(0);
53247c6ae99SBarry Smith }
53347c6ae99SBarry Smith 
53447c6ae99SBarry Smith #undef __FUNCT__
535aa219208SBarry Smith #define __FUNCT__ "DMDAFormFunction"
53647c6ae99SBarry Smith /*@
537aa219208SBarry Smith     DMDAFormFunction - Evaluates a user provided function on each processor that
538aa219208SBarry Smith         share a DMDA
53947c6ae99SBarry Smith 
54047c6ae99SBarry Smith    Input Parameters:
5419a42bb27SBarry Smith +    da - the DM that defines the grid
54247c6ae99SBarry Smith .    vu - input vector
54347c6ae99SBarry Smith .    vfu - output vector
54447c6ae99SBarry Smith -    w - any user data
54547c6ae99SBarry Smith 
54647c6ae99SBarry Smith     Notes: Does NOT do ghost updates on vu upon entry
54747c6ae99SBarry Smith 
548aa219208SBarry Smith            This should eventually replace DMDAFormFunction1
54947c6ae99SBarry Smith 
55047c6ae99SBarry Smith     Level: advanced
55147c6ae99SBarry Smith 
552aa219208SBarry Smith .seealso: DMDAComputeJacobian1WithAdic()
55347c6ae99SBarry Smith 
55447c6ae99SBarry Smith @*/
555aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDAFormFunction(DM da,PetscErrorCode (*lf)(void),Vec vu,Vec vfu,void *w)
55647c6ae99SBarry Smith {
55747c6ae99SBarry Smith   PetscErrorCode ierr;
55847c6ae99SBarry Smith   void           *u,*fu;
559aa219208SBarry Smith   DMDALocalInfo    info;
560aa219208SBarry Smith   PetscErrorCode (*f)(DMDALocalInfo*,void*,void*,void*) = (PetscErrorCode (*)(DMDALocalInfo*,void*,void*,void*))lf;
56147c6ae99SBarry Smith 
56247c6ae99SBarry Smith   PetscFunctionBegin;
563aa219208SBarry Smith   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
564aa219208SBarry Smith   ierr = DMDAVecGetArray(da,vu,&u);CHKERRQ(ierr);
565aa219208SBarry Smith   ierr = DMDAVecGetArray(da,vfu,&fu);CHKERRQ(ierr);
56647c6ae99SBarry Smith 
56747c6ae99SBarry Smith   ierr = (*f)(&info,u,fu,w);
56847c6ae99SBarry Smith   if (PetscExceptionValue(ierr)) {
569aa219208SBarry Smith     PetscErrorCode pierr = DMDAVecRestoreArray(da,vu,&u);CHKERRQ(pierr);
570aa219208SBarry Smith     pierr = DMDAVecRestoreArray(da,vfu,&fu);CHKERRQ(pierr);
57147c6ae99SBarry Smith   }
57247c6ae99SBarry Smith  CHKERRQ(ierr);
57347c6ae99SBarry Smith 
574aa219208SBarry Smith   ierr = DMDAVecRestoreArray(da,vu,&u);CHKERRQ(ierr);
575aa219208SBarry Smith   ierr = DMDAVecRestoreArray(da,vfu,&fu);CHKERRQ(ierr);
57647c6ae99SBarry Smith   PetscFunctionReturn(0);
57747c6ae99SBarry Smith }
57847c6ae99SBarry Smith 
57947c6ae99SBarry Smith #undef __FUNCT__
580aa219208SBarry Smith #define __FUNCT__ "DMDAFormFunctionLocal"
58147c6ae99SBarry Smith /*@C
582aa219208SBarry Smith    DMDAFormFunctionLocal - This is a universal function evaluation routine for
5839a42bb27SBarry Smith    a local DM function.
58447c6ae99SBarry Smith 
585aa219208SBarry Smith    Collective on DMDA
58647c6ae99SBarry Smith 
58747c6ae99SBarry Smith    Input Parameters:
5889a42bb27SBarry Smith +  da - the DM context
58947c6ae99SBarry Smith .  func - The local function
59047c6ae99SBarry Smith .  X - input vector
59147c6ae99SBarry Smith .  F - function vector
59247c6ae99SBarry Smith -  ctx - A user context
59347c6ae99SBarry Smith 
59447c6ae99SBarry Smith    Level: intermediate
59547c6ae99SBarry Smith 
596aa219208SBarry Smith .seealso: DMDASetLocalFunction(), DMDASetLocalJacobian(), DMDASetLocalAdicFunction(), DMDASetLocalAdicMFFunction(),
59747c6ae99SBarry Smith           SNESSetFunction(), SNESSetJacobian()
59847c6ae99SBarry Smith 
59947c6ae99SBarry Smith @*/
600aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDAFormFunctionLocal(DM da, DMDALocalFunction1 func, Vec X, Vec F, void *ctx)
60147c6ae99SBarry Smith {
60247c6ae99SBarry Smith   Vec            localX;
603aa219208SBarry Smith   DMDALocalInfo    info;
60447c6ae99SBarry Smith   void          *u;
60547c6ae99SBarry Smith   void          *fu;
60647c6ae99SBarry Smith   PetscErrorCode ierr;
60747c6ae99SBarry Smith 
60847c6ae99SBarry Smith   PetscFunctionBegin;
6099a42bb27SBarry Smith   ierr = DMGetLocalVector(da,&localX);CHKERRQ(ierr);
61047c6ae99SBarry Smith   /*
61147c6ae99SBarry Smith      Scatter ghost points to local vector, using the 2-step process
6129a42bb27SBarry Smith         DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
61347c6ae99SBarry Smith   */
6149a42bb27SBarry Smith   ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
6159a42bb27SBarry Smith   ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
616aa219208SBarry Smith   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
617aa219208SBarry Smith   ierr = DMDAVecGetArray(da,localX,&u);CHKERRQ(ierr);
618aa219208SBarry Smith   ierr = DMDAVecGetArray(da,F,&fu);CHKERRQ(ierr);
61947c6ae99SBarry Smith   ierr = (*func)(&info,u,fu,ctx);
62047c6ae99SBarry Smith   if (PetscExceptionValue(ierr)) {
621aa219208SBarry Smith     PetscErrorCode pierr = DMDAVecRestoreArray(da,localX,&u);CHKERRQ(pierr);
622aa219208SBarry Smith     pierr = DMDAVecRestoreArray(da,F,&fu);CHKERRQ(pierr);
62347c6ae99SBarry Smith   }
62447c6ae99SBarry Smith  CHKERRQ(ierr);
625aa219208SBarry Smith   ierr = DMDAVecRestoreArray(da,localX,&u);CHKERRQ(ierr);
626aa219208SBarry Smith   ierr = DMDAVecRestoreArray(da,F,&fu);CHKERRQ(ierr);
62747c6ae99SBarry Smith   if (PetscExceptionValue(ierr)) {
6289a42bb27SBarry Smith     PetscErrorCode pierr = DMRestoreLocalVector(da,&localX);CHKERRQ(pierr);
62947c6ae99SBarry Smith   }
63047c6ae99SBarry Smith  CHKERRQ(ierr);
6319a42bb27SBarry Smith   ierr = DMRestoreLocalVector(da,&localX);CHKERRQ(ierr);
63247c6ae99SBarry Smith   PetscFunctionReturn(0);
63347c6ae99SBarry Smith }
63447c6ae99SBarry Smith 
63547c6ae99SBarry Smith #undef __FUNCT__
636aa219208SBarry Smith #define __FUNCT__ "DMDAFormFunctionLocalGhost"
63747c6ae99SBarry Smith /*@C
638aa219208SBarry Smith    DMDAFormFunctionLocalGhost - This is a universal function evaluation routine for
6399a42bb27SBarry Smith    a local DM function, but the ghost values of the output are communicated and added.
64047c6ae99SBarry Smith 
641aa219208SBarry Smith    Collective on DMDA
64247c6ae99SBarry Smith 
64347c6ae99SBarry Smith    Input Parameters:
6449a42bb27SBarry Smith +  da - the DM context
64547c6ae99SBarry Smith .  func - The local function
64647c6ae99SBarry Smith .  X - input vector
64747c6ae99SBarry Smith .  F - function vector
64847c6ae99SBarry Smith -  ctx - A user context
64947c6ae99SBarry Smith 
65047c6ae99SBarry Smith    Level: intermediate
65147c6ae99SBarry Smith 
652aa219208SBarry Smith .seealso: DMDASetLocalFunction(), DMDASetLocalJacobian(), DMDASetLocalAdicFunction(), DMDASetLocalAdicMFFunction(),
65347c6ae99SBarry Smith           SNESSetFunction(), SNESSetJacobian()
65447c6ae99SBarry Smith 
65547c6ae99SBarry Smith @*/
656aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDAFormFunctionLocalGhost(DM da, DMDALocalFunction1 func, Vec X, Vec F, void *ctx)
65747c6ae99SBarry Smith {
65847c6ae99SBarry Smith   Vec            localX, localF;
659aa219208SBarry Smith   DMDALocalInfo    info;
66047c6ae99SBarry Smith   void          *u;
66147c6ae99SBarry Smith   void          *fu;
66247c6ae99SBarry Smith   PetscErrorCode ierr;
66347c6ae99SBarry Smith 
66447c6ae99SBarry Smith   PetscFunctionBegin;
6659a42bb27SBarry Smith   ierr = DMGetLocalVector(da,&localX);CHKERRQ(ierr);
6669a42bb27SBarry Smith   ierr = DMGetLocalVector(da,&localF);CHKERRQ(ierr);
66747c6ae99SBarry Smith   /*
66847c6ae99SBarry Smith      Scatter ghost points to local vector, using the 2-step process
6699a42bb27SBarry Smith         DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
67047c6ae99SBarry Smith   */
6719a42bb27SBarry Smith   ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
6729a42bb27SBarry Smith   ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
67347c6ae99SBarry Smith   ierr = VecSet(F, 0.0);CHKERRQ(ierr);
67447c6ae99SBarry Smith   ierr = VecSet(localF, 0.0);CHKERRQ(ierr);
675aa219208SBarry Smith   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
676aa219208SBarry Smith   ierr = DMDAVecGetArray(da,localX,&u);CHKERRQ(ierr);
677aa219208SBarry Smith   ierr = DMDAVecGetArray(da,localF,&fu);CHKERRQ(ierr);
67847c6ae99SBarry Smith   ierr = (*func)(&info,u,fu,ctx);
67947c6ae99SBarry Smith   if (PetscExceptionValue(ierr)) {
680aa219208SBarry Smith     PetscErrorCode pierr = DMDAVecRestoreArray(da,localX,&u);CHKERRQ(pierr);
681aa219208SBarry Smith     pierr = DMDAVecRestoreArray(da,localF,&fu);CHKERRQ(pierr);
68247c6ae99SBarry Smith   }
68347c6ae99SBarry Smith   CHKERRQ(ierr);
6849a42bb27SBarry Smith   ierr = DMLocalToGlobalBegin(da,localF,ADD_VALUES,F);CHKERRQ(ierr);
6859a42bb27SBarry Smith   ierr = DMLocalToGlobalEnd(da,localF,ADD_VALUES,F);CHKERRQ(ierr);
686aa219208SBarry Smith   ierr = DMDAVecRestoreArray(da,localX,&u);CHKERRQ(ierr);
687aa219208SBarry Smith   ierr = DMDAVecRestoreArray(da,localF,&fu);CHKERRQ(ierr);
68847c6ae99SBarry Smith   if (PetscExceptionValue(ierr)) {
6899a42bb27SBarry Smith     PetscErrorCode pierr = DMRestoreLocalVector(da,&localX);CHKERRQ(pierr);
6909a42bb27SBarry Smith   ierr = DMRestoreLocalVector(da,&localF);CHKERRQ(ierr);
69147c6ae99SBarry Smith   }
69247c6ae99SBarry Smith   CHKERRQ(ierr);
6939a42bb27SBarry Smith   ierr = DMRestoreLocalVector(da,&localX);CHKERRQ(ierr);
6949a42bb27SBarry Smith   ierr = DMRestoreLocalVector(da,&localF);CHKERRQ(ierr);
69547c6ae99SBarry Smith   PetscFunctionReturn(0);
69647c6ae99SBarry Smith }
69747c6ae99SBarry Smith 
69847c6ae99SBarry Smith #undef __FUNCT__
699aa219208SBarry Smith #define __FUNCT__ "DMDAFormFunction1"
70047c6ae99SBarry Smith /*@
701aa219208SBarry Smith     DMDAFormFunction1 - Evaluates a user provided function on each processor that
702aa219208SBarry Smith         share a DMDA
70347c6ae99SBarry Smith 
70447c6ae99SBarry Smith    Input Parameters:
7059a42bb27SBarry Smith +    da - the DM that defines the grid
70647c6ae99SBarry Smith .    vu - input vector
70747c6ae99SBarry Smith .    vfu - output vector
70847c6ae99SBarry Smith -    w - any user data
70947c6ae99SBarry Smith 
71047c6ae99SBarry Smith     Notes: Does NOT do ghost updates on vu upon entry
71147c6ae99SBarry Smith 
71247c6ae99SBarry Smith     Level: advanced
71347c6ae99SBarry Smith 
714aa219208SBarry Smith .seealso: DMDAComputeJacobian1WithAdic()
71547c6ae99SBarry Smith 
71647c6ae99SBarry Smith @*/
717aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDAFormFunction1(DM da,Vec vu,Vec vfu,void *w)
71847c6ae99SBarry Smith {
71947c6ae99SBarry Smith   PetscErrorCode ierr;
72047c6ae99SBarry Smith   void           *u,*fu;
721aa219208SBarry Smith   DMDALocalInfo    info;
72247c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
72347c6ae99SBarry Smith 
72447c6ae99SBarry Smith   PetscFunctionBegin;
72547c6ae99SBarry Smith 
726aa219208SBarry Smith   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
727aa219208SBarry Smith   ierr = DMDAVecGetArray(da,vu,&u);CHKERRQ(ierr);
728aa219208SBarry Smith   ierr = DMDAVecGetArray(da,vfu,&fu);CHKERRQ(ierr);
72947c6ae99SBarry Smith 
73047c6ae99SBarry Smith   CHKMEMQ;
73147c6ae99SBarry Smith   ierr = (*dd->lf)(&info,u,fu,w);
73247c6ae99SBarry Smith   if (PetscExceptionValue(ierr)) {
733aa219208SBarry Smith     PetscErrorCode pierr = DMDAVecRestoreArray(da,vu,&u);CHKERRQ(pierr);
734aa219208SBarry Smith     pierr = DMDAVecRestoreArray(da,vfu,&fu);CHKERRQ(pierr);
73547c6ae99SBarry Smith   }
73647c6ae99SBarry Smith   CHKERRQ(ierr);
73747c6ae99SBarry Smith   CHKMEMQ;
73847c6ae99SBarry Smith 
739aa219208SBarry Smith   ierr = DMDAVecRestoreArray(da,vu,&u);CHKERRQ(ierr);
740aa219208SBarry Smith   ierr = DMDAVecRestoreArray(da,vfu,&fu);CHKERRQ(ierr);
74147c6ae99SBarry Smith   PetscFunctionReturn(0);
74247c6ae99SBarry Smith }
74347c6ae99SBarry Smith 
74447c6ae99SBarry Smith #undef __FUNCT__
745aa219208SBarry Smith #define __FUNCT__ "DMDAFormFunctioniTest1"
746aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDAFormFunctioniTest1(DM da,void *w)
74747c6ae99SBarry Smith {
74847c6ae99SBarry Smith   Vec            vu,fu,fui;
74947c6ae99SBarry Smith   PetscErrorCode ierr;
75047c6ae99SBarry Smith   PetscInt       i,n;
75147c6ae99SBarry Smith   PetscScalar    *ui;
75247c6ae99SBarry Smith   PetscRandom    rnd;
75347c6ae99SBarry Smith   PetscReal      norm;
75447c6ae99SBarry Smith 
75547c6ae99SBarry Smith   PetscFunctionBegin;
7569a42bb27SBarry Smith   ierr = DMGetLocalVector(da,&vu);CHKERRQ(ierr);
75747c6ae99SBarry Smith   ierr = PetscRandomCreate(PETSC_COMM_SELF,&rnd);CHKERRQ(ierr);
75847c6ae99SBarry Smith   ierr = PetscRandomSetFromOptions(rnd);CHKERRQ(ierr);
75947c6ae99SBarry Smith   ierr = VecSetRandom(vu,rnd);CHKERRQ(ierr);
76047c6ae99SBarry Smith   ierr = PetscRandomDestroy(rnd);CHKERRQ(ierr);
76147c6ae99SBarry Smith 
7629a42bb27SBarry Smith   ierr = DMGetGlobalVector(da,&fu);CHKERRQ(ierr);
7639a42bb27SBarry Smith   ierr = DMGetGlobalVector(da,&fui);CHKERRQ(ierr);
76447c6ae99SBarry Smith 
765aa219208SBarry Smith   ierr = DMDAFormFunction1(da,vu,fu,w);CHKERRQ(ierr);
76647c6ae99SBarry Smith 
76747c6ae99SBarry Smith   ierr = VecGetArray(fui,&ui);CHKERRQ(ierr);
76847c6ae99SBarry Smith   ierr = VecGetLocalSize(fui,&n);CHKERRQ(ierr);
76947c6ae99SBarry Smith   for (i=0; i<n; i++) {
770aa219208SBarry Smith     ierr = DMDAFormFunctioni1(da,i,vu,ui+i,w);CHKERRQ(ierr);
77147c6ae99SBarry Smith   }
77247c6ae99SBarry Smith   ierr = VecRestoreArray(fui,&ui);CHKERRQ(ierr);
77347c6ae99SBarry Smith 
77447c6ae99SBarry Smith   ierr = VecAXPY(fui,-1.0,fu);CHKERRQ(ierr);
77547c6ae99SBarry Smith   ierr = VecNorm(fui,NORM_2,&norm);CHKERRQ(ierr);
77647c6ae99SBarry Smith   ierr = PetscPrintf(((PetscObject)da)->comm,"Norm of difference in vectors %G\n",norm);CHKERRQ(ierr);
77747c6ae99SBarry Smith   ierr = VecView(fu,0);CHKERRQ(ierr);
77847c6ae99SBarry Smith   ierr = VecView(fui,0);CHKERRQ(ierr);
77947c6ae99SBarry Smith 
7809a42bb27SBarry Smith   ierr = DMRestoreLocalVector(da,&vu);CHKERRQ(ierr);
7819a42bb27SBarry Smith   ierr = DMRestoreGlobalVector(da,&fu);CHKERRQ(ierr);
7829a42bb27SBarry Smith   ierr = DMRestoreGlobalVector(da,&fui);CHKERRQ(ierr);
78347c6ae99SBarry Smith   PetscFunctionReturn(0);
78447c6ae99SBarry Smith }
78547c6ae99SBarry Smith 
78647c6ae99SBarry Smith #undef __FUNCT__
787aa219208SBarry Smith #define __FUNCT__ "DMDAFormFunctioni1"
78847c6ae99SBarry Smith /*@
789aa219208SBarry Smith     DMDAFormFunctioni1 - Evaluates a user provided point-wise function
79047c6ae99SBarry Smith 
79147c6ae99SBarry Smith    Input Parameters:
7929a42bb27SBarry Smith +    da - the DM that defines the grid
79347c6ae99SBarry Smith .    i - the component of the function we wish to compute (must be local)
79447c6ae99SBarry Smith .    vu - input vector
79547c6ae99SBarry Smith .    vfu - output value
79647c6ae99SBarry Smith -    w - any user data
79747c6ae99SBarry Smith 
79847c6ae99SBarry Smith     Notes: Does NOT do ghost updates on vu upon entry
79947c6ae99SBarry Smith 
80047c6ae99SBarry Smith     Level: advanced
80147c6ae99SBarry Smith 
802aa219208SBarry Smith .seealso: DMDAComputeJacobian1WithAdic()
80347c6ae99SBarry Smith 
80447c6ae99SBarry Smith @*/
805aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDAFormFunctioni1(DM da,PetscInt i,Vec vu,PetscScalar *vfu,void *w)
80647c6ae99SBarry Smith {
80747c6ae99SBarry Smith   PetscErrorCode ierr;
80847c6ae99SBarry Smith   void           *u;
809aa219208SBarry Smith   DMDALocalInfo    info;
81047c6ae99SBarry Smith   MatStencil     stencil;
81147c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
81247c6ae99SBarry Smith 
81347c6ae99SBarry Smith   PetscFunctionBegin;
81447c6ae99SBarry Smith 
815aa219208SBarry Smith   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
816aa219208SBarry Smith   ierr = DMDAVecGetArray(da,vu,&u);CHKERRQ(ierr);
81747c6ae99SBarry Smith 
81847c6ae99SBarry Smith   /* figure out stencil value from i */
81947c6ae99SBarry Smith   stencil.c = i % info.dof;
82047c6ae99SBarry Smith   stencil.i = (i % (info.xm*info.dof))/info.dof;
82147c6ae99SBarry Smith   stencil.j = (i % (info.xm*info.ym*info.dof))/(info.xm*info.dof);
82247c6ae99SBarry Smith   stencil.k = i/(info.xm*info.ym*info.dof);
82347c6ae99SBarry Smith 
82447c6ae99SBarry Smith   ierr = (*dd->lfi)(&info,&stencil,u,vfu,w);CHKERRQ(ierr);
82547c6ae99SBarry Smith 
826aa219208SBarry Smith   ierr = DMDAVecRestoreArray(da,vu,&u);CHKERRQ(ierr);
82747c6ae99SBarry Smith   PetscFunctionReturn(0);
82847c6ae99SBarry Smith }
82947c6ae99SBarry Smith 
83047c6ae99SBarry Smith #undef __FUNCT__
831aa219208SBarry Smith #define __FUNCT__ "DMDAFormFunctionib1"
83247c6ae99SBarry Smith /*@
833aa219208SBarry Smith     DMDAFormFunctionib1 - Evaluates a user provided point-block function
83447c6ae99SBarry Smith 
83547c6ae99SBarry Smith    Input Parameters:
8369a42bb27SBarry Smith +    da - the DM that defines the grid
83747c6ae99SBarry Smith .    i - the component of the function we wish to compute (must be local)
83847c6ae99SBarry Smith .    vu - input vector
83947c6ae99SBarry Smith .    vfu - output value
84047c6ae99SBarry Smith -    w - any user data
84147c6ae99SBarry Smith 
84247c6ae99SBarry Smith     Notes: Does NOT do ghost updates on vu upon entry
84347c6ae99SBarry Smith 
84447c6ae99SBarry Smith     Level: advanced
84547c6ae99SBarry Smith 
846aa219208SBarry Smith .seealso: DMDAComputeJacobian1WithAdic()
84747c6ae99SBarry Smith 
84847c6ae99SBarry Smith @*/
849aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDAFormFunctionib1(DM da,PetscInt i,Vec vu,PetscScalar *vfu,void *w)
85047c6ae99SBarry Smith {
85147c6ae99SBarry Smith   PetscErrorCode ierr;
85247c6ae99SBarry Smith   void           *u;
853aa219208SBarry Smith   DMDALocalInfo    info;
85447c6ae99SBarry Smith   MatStencil     stencil;
85547c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
85647c6ae99SBarry Smith 
85747c6ae99SBarry Smith   PetscFunctionBegin;
858aa219208SBarry Smith   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
859aa219208SBarry Smith   ierr = DMDAVecGetArray(da,vu,&u);CHKERRQ(ierr);
86047c6ae99SBarry Smith 
86147c6ae99SBarry Smith   /* figure out stencil value from i */
86247c6ae99SBarry Smith   stencil.c = i % info.dof;
86347c6ae99SBarry Smith   if (stencil.c) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONG,"Point-block functions can only be called for the entire block");
86447c6ae99SBarry Smith   stencil.i = (i % (info.xm*info.dof))/info.dof;
86547c6ae99SBarry Smith   stencil.j = (i % (info.xm*info.ym*info.dof))/(info.xm*info.dof);
86647c6ae99SBarry Smith   stencil.k = i/(info.xm*info.ym*info.dof);
86747c6ae99SBarry Smith 
86847c6ae99SBarry Smith   ierr = (*dd->lfib)(&info,&stencil,u,vfu,w);CHKERRQ(ierr);
86947c6ae99SBarry Smith 
870aa219208SBarry Smith   ierr = DMDAVecRestoreArray(da,vu,&u);CHKERRQ(ierr);
87147c6ae99SBarry Smith   PetscFunctionReturn(0);
87247c6ae99SBarry Smith }
87347c6ae99SBarry Smith 
87447c6ae99SBarry Smith #if defined(new)
87547c6ae99SBarry Smith #undef __FUNCT__
876aa219208SBarry Smith #define __FUNCT__ "DMDAGetDiagonal_MFFD"
87747c6ae99SBarry Smith /*
878aa219208SBarry Smith   DMDAGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix where local
879aa219208SBarry Smith     function lives on a DMDA
88047c6ae99SBarry Smith 
88147c6ae99SBarry Smith         y ~= (F(u + ha) - F(u))/h,
88247c6ae99SBarry Smith   where F = nonlinear function, as set by SNESSetFunction()
88347c6ae99SBarry Smith         u = current iterate
88447c6ae99SBarry Smith         h = difference interval
88547c6ae99SBarry Smith */
886aa219208SBarry Smith PetscErrorCode DMDAGetDiagonal_MFFD(DM da,Vec U,Vec a)
88747c6ae99SBarry Smith {
88847c6ae99SBarry Smith   PetscScalar    h,*aa,*ww,v;
88947c6ae99SBarry Smith   PetscReal      epsilon = PETSC_SQRT_MACHINE_EPSILON,umin = 100.0*PETSC_SQRT_MACHINE_EPSILON;
89047c6ae99SBarry Smith   PetscErrorCode ierr;
89147c6ae99SBarry Smith   PetscInt       gI,nI;
89247c6ae99SBarry Smith   MatStencil     stencil;
893aa219208SBarry Smith   DMDALocalInfo    info;
89447c6ae99SBarry Smith 
89547c6ae99SBarry Smith   PetscFunctionBegin;
89647c6ae99SBarry Smith   ierr = (*ctx->func)(0,U,a,ctx->funcctx);CHKERRQ(ierr);
89747c6ae99SBarry Smith   ierr = (*ctx->funcisetbase)(U,ctx->funcctx);CHKERRQ(ierr);
89847c6ae99SBarry Smith 
89947c6ae99SBarry Smith   ierr = VecGetArray(U,&ww);CHKERRQ(ierr);
90047c6ae99SBarry Smith   ierr = VecGetArray(a,&aa);CHKERRQ(ierr);
90147c6ae99SBarry Smith 
90247c6ae99SBarry Smith   nI = 0;
90347c6ae99SBarry Smith     h  = ww[gI];
90447c6ae99SBarry Smith     if (h == 0.0) h = 1.0;
90547c6ae99SBarry Smith #if !defined(PETSC_USE_COMPLEX)
90647c6ae99SBarry Smith     if (h < umin && h >= 0.0)      h = umin;
90747c6ae99SBarry Smith     else if (h < 0.0 && h > -umin) h = -umin;
90847c6ae99SBarry Smith #else
90947c6ae99SBarry Smith     if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0)     h = umin;
91047c6ae99SBarry Smith     else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin;
91147c6ae99SBarry Smith #endif
91247c6ae99SBarry Smith     h     *= epsilon;
91347c6ae99SBarry Smith 
91447c6ae99SBarry Smith     ww[gI] += h;
91547c6ae99SBarry Smith     ierr          = (*ctx->funci)(i,w,&v,ctx->funcctx);CHKERRQ(ierr);
91647c6ae99SBarry Smith     aa[nI]  = (v - aa[nI])/h;
91747c6ae99SBarry Smith     ww[gI] -= h;
91847c6ae99SBarry Smith     nI++;
91947c6ae99SBarry Smith   }
92047c6ae99SBarry Smith   ierr = VecRestoreArray(U,&ww);CHKERRQ(ierr);
92147c6ae99SBarry Smith   ierr = VecRestoreArray(a,&aa);CHKERRQ(ierr);
92247c6ae99SBarry Smith   PetscFunctionReturn(0);
92347c6ae99SBarry Smith }
92447c6ae99SBarry Smith #endif
92547c6ae99SBarry Smith 
92647c6ae99SBarry Smith #if defined(PETSC_HAVE_ADIC)
92747c6ae99SBarry Smith EXTERN_C_BEGIN
92847c6ae99SBarry Smith #include "adic/ad_utils.h"
92947c6ae99SBarry Smith EXTERN_C_END
93047c6ae99SBarry Smith 
93147c6ae99SBarry Smith #undef __FUNCT__
932aa219208SBarry Smith #define __FUNCT__ "DMDAComputeJacobian1WithAdic"
93347c6ae99SBarry Smith /*@C
934aa219208SBarry Smith     DMDAComputeJacobian1WithAdic - Evaluates a adiC provided Jacobian function on each processor that
935aa219208SBarry Smith         share a DMDA
93647c6ae99SBarry Smith 
93747c6ae99SBarry Smith    Input Parameters:
9389a42bb27SBarry Smith +    da - the DM that defines the grid
93947c6ae99SBarry Smith .    vu - input vector (ghosted)
94047c6ae99SBarry Smith .    J - output matrix
94147c6ae99SBarry Smith -    w - any user data
94247c6ae99SBarry Smith 
94347c6ae99SBarry Smith    Level: advanced
94447c6ae99SBarry Smith 
94547c6ae99SBarry Smith     Notes: Does NOT do ghost updates on vu upon entry
94647c6ae99SBarry Smith 
947aa219208SBarry Smith .seealso: DMDAFormFunction1()
94847c6ae99SBarry Smith 
94947c6ae99SBarry Smith @*/
950aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDAComputeJacobian1WithAdic(DM da,Vec vu,Mat J,void *w)
95147c6ae99SBarry Smith {
95247c6ae99SBarry Smith   PetscErrorCode ierr;
95347c6ae99SBarry Smith   PetscInt       gtdof,tdof;
95447c6ae99SBarry Smith   PetscScalar    *ustart;
955aa219208SBarry Smith   DMDALocalInfo    info;
95647c6ae99SBarry Smith   void           *ad_u,*ad_f,*ad_ustart,*ad_fstart;
95747c6ae99SBarry Smith   ISColoring     iscoloring;
95847c6ae99SBarry Smith 
95947c6ae99SBarry Smith   PetscFunctionBegin;
960aa219208SBarry Smith   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
96147c6ae99SBarry Smith 
96247c6ae99SBarry Smith   PetscADResetIndep();
96347c6ae99SBarry Smith 
96447c6ae99SBarry Smith   /* get space for derivative objects.  */
965aa219208SBarry Smith   ierr = DMDAGetAdicArray(da,PETSC_TRUE,&ad_u,&ad_ustart,&gtdof);CHKERRQ(ierr);
966aa219208SBarry Smith   ierr = DMDAGetAdicArray(da,PETSC_FALSE,&ad_f,&ad_fstart,&tdof);CHKERRQ(ierr);
96747c6ae99SBarry Smith   ierr = VecGetArray(vu,&ustart);CHKERRQ(ierr);
96894013140SBarry Smith   ierr = DMGetColoring(da,IS_COLORING_GHOSTED,MATAIJ,&iscoloring);CHKERRQ(ierr);
96947c6ae99SBarry Smith 
97047c6ae99SBarry Smith   PetscADSetValueAndColor(ad_ustart,gtdof,iscoloring->colors,ustart);
97147c6ae99SBarry Smith 
97247c6ae99SBarry Smith   ierr = VecRestoreArray(vu,&ustart);CHKERRQ(ierr);
97347c6ae99SBarry Smith   ierr = ISColoringDestroy(iscoloring);CHKERRQ(ierr);
97447c6ae99SBarry Smith   ierr = PetscADIncrementTotalGradSize(iscoloring->n);CHKERRQ(ierr);
97547c6ae99SBarry Smith   PetscADSetIndepDone();
97647c6ae99SBarry Smith 
977aa219208SBarry Smith   ierr = PetscLogEventBegin(DMDA_LocalADFunction,0,0,0,0);CHKERRQ(ierr);
97847c6ae99SBarry Smith   ierr = (*dd->adic_lf)(&info,ad_u,ad_f,w);CHKERRQ(ierr);
979aa219208SBarry Smith   ierr = PetscLogEventEnd(DMDA_LocalADFunction,0,0,0,0);CHKERRQ(ierr);
98047c6ae99SBarry Smith 
98147c6ae99SBarry Smith   /* stick the values into the matrix */
98247c6ae99SBarry Smith   ierr = MatSetValuesAdic(J,(PetscScalar**)ad_fstart);CHKERRQ(ierr);
98347c6ae99SBarry Smith 
98447c6ae99SBarry Smith   /* return space for derivative objects.  */
985aa219208SBarry Smith   ierr = DMDARestoreAdicArray(da,PETSC_TRUE,&ad_u,&ad_ustart,&gtdof);CHKERRQ(ierr);
986aa219208SBarry Smith   ierr = DMDARestoreAdicArray(da,PETSC_FALSE,&ad_f,&ad_fstart,&tdof);CHKERRQ(ierr);
98747c6ae99SBarry Smith   PetscFunctionReturn(0);
98847c6ae99SBarry Smith }
98947c6ae99SBarry Smith 
99047c6ae99SBarry Smith #undef __FUNCT__
991aa219208SBarry Smith #define __FUNCT__ "DMDAMultiplyByJacobian1WithAdic"
99247c6ae99SBarry Smith /*@C
993aa219208SBarry Smith     DMDAMultiplyByJacobian1WithAdic - Applies an ADIC-provided Jacobian function to a vector on
994aa219208SBarry Smith     each processor that shares a DMDA.
99547c6ae99SBarry Smith 
99647c6ae99SBarry Smith     Input Parameters:
9979a42bb27SBarry Smith +   da - the DM that defines the grid
99847c6ae99SBarry Smith .   vu - Jacobian is computed at this point (ghosted)
99947c6ae99SBarry Smith .   v - product is done on this vector (ghosted)
100047c6ae99SBarry Smith .   fu - output vector = J(vu)*v (not ghosted)
100147c6ae99SBarry Smith -   w - any user data
100247c6ae99SBarry Smith 
100347c6ae99SBarry Smith     Notes:
100447c6ae99SBarry Smith     This routine does NOT do ghost updates on vu upon entry.
100547c6ae99SBarry Smith 
100647c6ae99SBarry Smith    Level: advanced
100747c6ae99SBarry Smith 
1008aa219208SBarry Smith .seealso: DMDAFormFunction1()
100947c6ae99SBarry Smith 
101047c6ae99SBarry Smith @*/
1011aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDAMultiplyByJacobian1WithAdic(DM da,Vec vu,Vec v,Vec f,void *w)
101247c6ae99SBarry Smith {
101347c6ae99SBarry Smith   PetscErrorCode ierr;
101447c6ae99SBarry Smith   PetscInt       i,gtdof,tdof;
101547c6ae99SBarry Smith   PetscScalar    *avu,*av,*af,*ad_vustart,*ad_fstart;
1016aa219208SBarry Smith   DMDALocalInfo    info;
101747c6ae99SBarry Smith   void           *ad_vu,*ad_f;
101847c6ae99SBarry Smith 
101947c6ae99SBarry Smith   PetscFunctionBegin;
1020aa219208SBarry Smith   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
102147c6ae99SBarry Smith 
102247c6ae99SBarry Smith   /* get space for derivative objects.  */
1023aa219208SBarry Smith   ierr = DMDAGetAdicMFArray(da,PETSC_TRUE,&ad_vu,&ad_vustart,&gtdof);CHKERRQ(ierr);
1024aa219208SBarry Smith   ierr = DMDAGetAdicMFArray(da,PETSC_FALSE,&ad_f,&ad_fstart,&tdof);CHKERRQ(ierr);
102547c6ae99SBarry Smith 
102647c6ae99SBarry Smith   /* copy input vector into derivative object */
102747c6ae99SBarry Smith   ierr = VecGetArray(vu,&avu);CHKERRQ(ierr);
102847c6ae99SBarry Smith   ierr = VecGetArray(v,&av);CHKERRQ(ierr);
102947c6ae99SBarry Smith   for (i=0; i<gtdof; i++) {
103047c6ae99SBarry Smith     ad_vustart[2*i]   = avu[i];
103147c6ae99SBarry Smith     ad_vustart[2*i+1] = av[i];
103247c6ae99SBarry Smith   }
103347c6ae99SBarry Smith   ierr = VecRestoreArray(vu,&avu);CHKERRQ(ierr);
103447c6ae99SBarry Smith   ierr = VecRestoreArray(v,&av);CHKERRQ(ierr);
103547c6ae99SBarry Smith 
103647c6ae99SBarry Smith   PetscADResetIndep();
103747c6ae99SBarry Smith   ierr = PetscADIncrementTotalGradSize(1);CHKERRQ(ierr);
103847c6ae99SBarry Smith   PetscADSetIndepDone();
103947c6ae99SBarry Smith 
104047c6ae99SBarry Smith   ierr = (*dd->adicmf_lf)(&info,ad_vu,ad_f,w);CHKERRQ(ierr);
104147c6ae99SBarry Smith 
104247c6ae99SBarry Smith   /* stick the values into the vector */
104347c6ae99SBarry Smith   ierr = VecGetArray(f,&af);CHKERRQ(ierr);
104447c6ae99SBarry Smith   for (i=0; i<tdof; i++) {
104547c6ae99SBarry Smith     af[i] = ad_fstart[2*i+1];
104647c6ae99SBarry Smith   }
104747c6ae99SBarry Smith   ierr = VecRestoreArray(f,&af);CHKERRQ(ierr);
104847c6ae99SBarry Smith 
104947c6ae99SBarry Smith   /* return space for derivative objects.  */
1050aa219208SBarry Smith   ierr = DMDARestoreAdicMFArray(da,PETSC_TRUE,&ad_vu,&ad_vustart,&gtdof);CHKERRQ(ierr);
1051aa219208SBarry Smith   ierr = DMDARestoreAdicMFArray(da,PETSC_FALSE,&ad_f,&ad_fstart,&tdof);CHKERRQ(ierr);
105247c6ae99SBarry Smith   PetscFunctionReturn(0);
105347c6ae99SBarry Smith }
105447c6ae99SBarry Smith #endif
105547c6ae99SBarry Smith 
105647c6ae99SBarry Smith #undef __FUNCT__
1057aa219208SBarry Smith #define __FUNCT__ "DMDAComputeJacobian1"
105847c6ae99SBarry Smith /*@
1059aa219208SBarry Smith     DMDAComputeJacobian1 - Evaluates a local Jacobian function on each processor that
1060aa219208SBarry Smith         share a DMDA
106147c6ae99SBarry Smith 
106247c6ae99SBarry Smith    Input Parameters:
10639a42bb27SBarry Smith +    da - the DM that defines the grid
106447c6ae99SBarry Smith .    vu - input vector (ghosted)
106547c6ae99SBarry Smith .    J - output matrix
106647c6ae99SBarry Smith -    w - any user data
106747c6ae99SBarry Smith 
106847c6ae99SBarry Smith     Notes: Does NOT do ghost updates on vu upon entry
106947c6ae99SBarry Smith 
107047c6ae99SBarry Smith     Level: advanced
107147c6ae99SBarry Smith 
1072aa219208SBarry Smith .seealso: DMDAFormFunction1()
107347c6ae99SBarry Smith 
107447c6ae99SBarry Smith @*/
1075aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDAComputeJacobian1(DM da,Vec vu,Mat J,void *w)
107647c6ae99SBarry Smith {
107747c6ae99SBarry Smith   PetscErrorCode ierr;
107847c6ae99SBarry Smith   void           *u;
1079aa219208SBarry Smith   DMDALocalInfo    info;
108047c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
108147c6ae99SBarry Smith 
108247c6ae99SBarry Smith   PetscFunctionBegin;
1083aa219208SBarry Smith   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
1084aa219208SBarry Smith   ierr = DMDAVecGetArray(da,vu,&u);CHKERRQ(ierr);
108547c6ae99SBarry Smith   ierr = (*dd->lj)(&info,u,J,w);CHKERRQ(ierr);
1086aa219208SBarry Smith   ierr = DMDAVecRestoreArray(da,vu,&u);CHKERRQ(ierr);
108747c6ae99SBarry Smith   PetscFunctionReturn(0);
108847c6ae99SBarry Smith }
108947c6ae99SBarry Smith 
109047c6ae99SBarry Smith 
109147c6ae99SBarry Smith #undef __FUNCT__
1092aa219208SBarry Smith #define __FUNCT__ "DMDAComputeJacobian1WithAdifor"
109347c6ae99SBarry Smith /*
1094aa219208SBarry Smith     DMDAComputeJacobian1WithAdifor - Evaluates a ADIFOR provided Jacobian local function on each processor that
1095aa219208SBarry Smith         share a DMDA
109647c6ae99SBarry Smith 
109747c6ae99SBarry Smith    Input Parameters:
10989a42bb27SBarry Smith +    da - the DM that defines the grid
109947c6ae99SBarry Smith .    vu - input vector (ghosted)
110047c6ae99SBarry Smith .    J - output matrix
110147c6ae99SBarry Smith -    w - any user data
110247c6ae99SBarry Smith 
110347c6ae99SBarry Smith     Notes: Does NOT do ghost updates on vu upon entry
110447c6ae99SBarry Smith 
1105aa219208SBarry Smith .seealso: DMDAFormFunction1()
110647c6ae99SBarry Smith 
110747c6ae99SBarry Smith */
1108aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDAComputeJacobian1WithAdifor(DM da,Vec vu,Mat J,void *w)
110947c6ae99SBarry Smith {
111047c6ae99SBarry Smith   PetscErrorCode  ierr;
111147c6ae99SBarry Smith   PetscInt        i,Nc,N;
111247c6ae99SBarry Smith   ISColoringValue *color;
1113aa219208SBarry Smith   DMDALocalInfo     info;
111447c6ae99SBarry Smith   PetscScalar     *u,*g_u,*g_f,*f = 0,*p_u;
111547c6ae99SBarry Smith   ISColoring      iscoloring;
111647c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
1117aa219208SBarry Smith   void            (*lf)(PetscInt*,DMDALocalInfo*,PetscScalar*,PetscScalar*,PetscInt*,PetscScalar*,PetscScalar*,PetscInt*,void*,PetscErrorCode*) =
1118aa219208SBarry Smith                   (void (*)(PetscInt*,DMDALocalInfo*,PetscScalar*,PetscScalar*,PetscInt*,PetscScalar*,PetscScalar*,PetscInt*,void*,PetscErrorCode*))*dd->adifor_lf;
111947c6ae99SBarry Smith 
112047c6ae99SBarry Smith   PetscFunctionBegin;
112194013140SBarry Smith   ierr = DMGetColoring(da,IS_COLORING_GHOSTED,MATAIJ,&iscoloring);CHKERRQ(ierr);
112247c6ae99SBarry Smith   Nc   = iscoloring->n;
1123aa219208SBarry Smith   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
112447c6ae99SBarry Smith   N    = info.gxm*info.gym*info.gzm*info.dof;
112547c6ae99SBarry Smith 
112647c6ae99SBarry Smith   /* get space for derivative objects.  */
112747c6ae99SBarry Smith   ierr  = PetscMalloc(Nc*info.gxm*info.gym*info.gzm*info.dof*sizeof(PetscScalar),&g_u);CHKERRQ(ierr);
112847c6ae99SBarry Smith   ierr  = PetscMemzero(g_u,Nc*info.gxm*info.gym*info.gzm*info.dof*sizeof(PetscScalar));CHKERRQ(ierr);
112947c6ae99SBarry Smith   p_u   = g_u;
113047c6ae99SBarry Smith   color = iscoloring->colors;
113147c6ae99SBarry Smith   for (i=0; i<N; i++) {
113247c6ae99SBarry Smith     p_u[*color++] = 1.0;
113347c6ae99SBarry Smith     p_u          += Nc;
113447c6ae99SBarry Smith   }
113547c6ae99SBarry Smith   ierr = ISColoringDestroy(iscoloring);CHKERRQ(ierr);
113647c6ae99SBarry 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);
113747c6ae99SBarry Smith 
113847c6ae99SBarry Smith   /* Seed the input array g_u with coloring information */
113947c6ae99SBarry Smith 
114047c6ae99SBarry Smith   ierr = VecGetArray(vu,&u);CHKERRQ(ierr);
114147c6ae99SBarry Smith   (lf)(&Nc,&info,u,g_u,&Nc,f,g_f,&Nc,w,&ierr);CHKERRQ(ierr);
114247c6ae99SBarry Smith   ierr = VecRestoreArray(vu,&u);CHKERRQ(ierr);
114347c6ae99SBarry Smith 
114447c6ae99SBarry Smith   /* stick the values into the matrix */
114547c6ae99SBarry Smith   /* PetscScalarView(Nc*info.xm*info.ym,g_f,0); */
114647c6ae99SBarry Smith   ierr = MatSetValuesAdifor(J,Nc,g_f);CHKERRQ(ierr);
114747c6ae99SBarry Smith 
114847c6ae99SBarry Smith   /* return space for derivative objects.  */
114947c6ae99SBarry Smith   ierr = PetscFree(g_u);CHKERRQ(ierr);
115047c6ae99SBarry Smith   ierr = PetscFree2(g_f,f);CHKERRQ(ierr);
115147c6ae99SBarry Smith   PetscFunctionReturn(0);
115247c6ae99SBarry Smith }
115347c6ae99SBarry Smith 
115447c6ae99SBarry Smith #undef __FUNCT__
1155aa219208SBarry Smith #define __FUNCT__ "DMDAFormJacobianLocal"
115647c6ae99SBarry Smith /*@C
1157aa219208SBarry Smith    DMDAFormjacobianLocal - This is a universal Jacobian evaluation routine for
11589a42bb27SBarry Smith    a local DM function.
115947c6ae99SBarry Smith 
1160aa219208SBarry Smith    Collective on DMDA
116147c6ae99SBarry Smith 
116247c6ae99SBarry Smith    Input Parameters:
11639a42bb27SBarry Smith +  da - the DM context
116447c6ae99SBarry Smith .  func - The local function
116547c6ae99SBarry Smith .  X - input vector
116647c6ae99SBarry Smith .  J - Jacobian matrix
116747c6ae99SBarry Smith -  ctx - A user context
116847c6ae99SBarry Smith 
116947c6ae99SBarry Smith    Level: intermediate
117047c6ae99SBarry Smith 
1171aa219208SBarry Smith .seealso: DMDASetLocalFunction(), DMDASetLocalJacobian(), DMDASetLocalAdicFunction(), DMDASetLocalAdicMFFunction(),
117247c6ae99SBarry Smith           SNESSetFunction(), SNESSetJacobian()
117347c6ae99SBarry Smith 
117447c6ae99SBarry Smith @*/
1175aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDAFormJacobianLocal(DM da, DMDALocalFunction1 func, Vec X, Mat J, void *ctx)
117647c6ae99SBarry Smith {
117747c6ae99SBarry Smith   Vec            localX;
1178aa219208SBarry Smith   DMDALocalInfo    info;
117947c6ae99SBarry Smith   void          *u;
118047c6ae99SBarry Smith   PetscErrorCode ierr;
118147c6ae99SBarry Smith 
118247c6ae99SBarry Smith   PetscFunctionBegin;
11839a42bb27SBarry Smith   ierr = DMGetLocalVector(da,&localX);CHKERRQ(ierr);
118447c6ae99SBarry Smith   /*
118547c6ae99SBarry Smith      Scatter ghost points to local vector, using the 2-step process
11869a42bb27SBarry Smith         DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
118747c6ae99SBarry Smith   */
11889a42bb27SBarry Smith   ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
11899a42bb27SBarry Smith   ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
1190aa219208SBarry Smith   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
1191aa219208SBarry Smith   ierr = DMDAVecGetArray(da,localX,&u);CHKERRQ(ierr);
119247c6ae99SBarry Smith   ierr = (*func)(&info,u,J,ctx);
119347c6ae99SBarry Smith   if (PetscExceptionValue(ierr)) {
1194aa219208SBarry Smith     PetscErrorCode pierr = DMDAVecRestoreArray(da,localX,&u);CHKERRQ(pierr);
119547c6ae99SBarry Smith   }
119647c6ae99SBarry Smith   CHKERRQ(ierr);
1197aa219208SBarry Smith   ierr = DMDAVecRestoreArray(da,localX,&u);CHKERRQ(ierr);
119847c6ae99SBarry Smith   if (PetscExceptionValue(ierr)) {
11999a42bb27SBarry Smith     PetscErrorCode pierr = DMRestoreLocalVector(da,&localX);CHKERRQ(pierr);
120047c6ae99SBarry Smith   }
120147c6ae99SBarry Smith   CHKERRQ(ierr);
12029a42bb27SBarry Smith   ierr = DMRestoreLocalVector(da,&localX);CHKERRQ(ierr);
120347c6ae99SBarry Smith   PetscFunctionReturn(0);
120447c6ae99SBarry Smith }
120547c6ae99SBarry Smith 
120647c6ae99SBarry Smith #undef __FUNCT__
1207aa219208SBarry Smith #define __FUNCT__ "DMDAMultiplyByJacobian1WithAD"
120847c6ae99SBarry Smith /*@C
1209aa219208SBarry Smith     DMDAMultiplyByJacobian1WithAD - Applies a Jacobian function supplied by ADIFOR or ADIC
1210aa219208SBarry Smith     to a vector on each processor that shares a DMDA.
121147c6ae99SBarry Smith 
121247c6ae99SBarry Smith    Input Parameters:
12139a42bb27SBarry Smith +    da - the DM that defines the grid
121447c6ae99SBarry Smith .    vu - Jacobian is computed at this point (ghosted)
121547c6ae99SBarry Smith .    v - product is done on this vector (ghosted)
121647c6ae99SBarry Smith .    fu - output vector = J(vu)*v (not ghosted)
121747c6ae99SBarry Smith -    w - any user data
121847c6ae99SBarry Smith 
121947c6ae99SBarry Smith     Notes:
122047c6ae99SBarry Smith     This routine does NOT do ghost updates on vu and v upon entry.
122147c6ae99SBarry Smith 
1222aa219208SBarry Smith     Automatically calls DMDAMultiplyByJacobian1WithAdifor() or DMDAMultiplyByJacobian1WithAdic()
1223aa219208SBarry Smith     depending on whether DMDASetLocalAdicMFFunction() or DMDASetLocalAdiforMFFunction() was called.
122447c6ae99SBarry Smith 
122547c6ae99SBarry Smith    Level: advanced
122647c6ae99SBarry Smith 
1227aa219208SBarry Smith .seealso: DMDAFormFunction1(), DMDAMultiplyByJacobian1WithAdifor(), DMDAMultiplyByJacobian1WithAdic()
122847c6ae99SBarry Smith 
122947c6ae99SBarry Smith @*/
1230aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDAMultiplyByJacobian1WithAD(DM da,Vec u,Vec v,Vec f,void *w)
123147c6ae99SBarry Smith {
123247c6ae99SBarry Smith   PetscErrorCode ierr;
123347c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
123447c6ae99SBarry Smith 
123547c6ae99SBarry Smith   PetscFunctionBegin;
123647c6ae99SBarry Smith   if (dd->adicmf_lf) {
123747c6ae99SBarry Smith #if defined(PETSC_HAVE_ADIC)
1238aa219208SBarry Smith     ierr = DMDAMultiplyByJacobian1WithAdic(da,u,v,f,w);CHKERRQ(ierr);
123947c6ae99SBarry Smith #else
124047c6ae99SBarry Smith     SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP_SYS,"Requires ADIC to be installed and cannot use complex numbers");
124147c6ae99SBarry Smith #endif
124247c6ae99SBarry Smith   } else if (dd->adiformf_lf) {
1243aa219208SBarry Smith     ierr = DMDAMultiplyByJacobian1WithAdifor(da,u,v,f,w);CHKERRQ(ierr);
124447c6ae99SBarry Smith   } else {
1245aa219208SBarry Smith     SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ORDER,"Must call DMDASetLocalAdiforMFFunction() or DMDASetLocalAdicMFFunction() before using");
124647c6ae99SBarry Smith   }
124747c6ae99SBarry Smith   PetscFunctionReturn(0);
124847c6ae99SBarry Smith }
124947c6ae99SBarry Smith 
125047c6ae99SBarry Smith 
125147c6ae99SBarry Smith #undef __FUNCT__
1252aa219208SBarry Smith #define __FUNCT__ "DMDAMultiplyByJacobian1WithAdifor"
125347c6ae99SBarry Smith /*@C
1254aa219208SBarry Smith     DMDAMultiplyByJacobian1WithAdifor - Applies a ADIFOR provided Jacobian function on each processor that
12559a42bb27SBarry Smith         share a DM to a vector
125647c6ae99SBarry Smith 
125747c6ae99SBarry Smith    Input Parameters:
12589a42bb27SBarry Smith +    da - the DM that defines the grid
125947c6ae99SBarry Smith .    vu - Jacobian is computed at this point (ghosted)
126047c6ae99SBarry Smith .    v - product is done on this vector (ghosted)
126147c6ae99SBarry Smith .    fu - output vector = J(vu)*v (not ghosted)
126247c6ae99SBarry Smith -    w - any user data
126347c6ae99SBarry Smith 
126447c6ae99SBarry Smith     Notes: Does NOT do ghost updates on vu and v upon entry
126547c6ae99SBarry Smith 
126647c6ae99SBarry Smith    Level: advanced
126747c6ae99SBarry Smith 
1268aa219208SBarry Smith .seealso: DMDAFormFunction1()
126947c6ae99SBarry Smith 
127047c6ae99SBarry Smith @*/
1271aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDAMultiplyByJacobian1WithAdifor(DM da,Vec u,Vec v,Vec f,void *w)
127247c6ae99SBarry Smith {
127347c6ae99SBarry Smith   PetscErrorCode ierr;
127447c6ae99SBarry Smith   PetscScalar    *au,*av,*af,*awork;
127547c6ae99SBarry Smith   Vec            work;
1276aa219208SBarry Smith   DMDALocalInfo    info;
127747c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
1278aa219208SBarry Smith   void           (*lf)(DMDALocalInfo*,PetscScalar*,PetscScalar*,PetscScalar*,PetscScalar*,void*,PetscErrorCode*) =
1279aa219208SBarry Smith                  (void (*)(DMDALocalInfo*,PetscScalar*,PetscScalar*,PetscScalar*,PetscScalar*,void*,PetscErrorCode*))*dd->adiformf_lf;
128047c6ae99SBarry Smith 
128147c6ae99SBarry Smith   PetscFunctionBegin;
1282aa219208SBarry Smith   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
128347c6ae99SBarry Smith 
12849a42bb27SBarry Smith   ierr = DMGetGlobalVector(da,&work);CHKERRQ(ierr);
128547c6ae99SBarry Smith   ierr = VecGetArray(u,&au);CHKERRQ(ierr);
128647c6ae99SBarry Smith   ierr = VecGetArray(v,&av);CHKERRQ(ierr);
128747c6ae99SBarry Smith   ierr = VecGetArray(f,&af);CHKERRQ(ierr);
128847c6ae99SBarry Smith   ierr = VecGetArray(work,&awork);CHKERRQ(ierr);
128947c6ae99SBarry Smith   (lf)(&info,au,av,awork,af,w,&ierr);CHKERRQ(ierr);
129047c6ae99SBarry Smith   ierr = VecRestoreArray(u,&au);CHKERRQ(ierr);
129147c6ae99SBarry Smith   ierr = VecRestoreArray(v,&av);CHKERRQ(ierr);
129247c6ae99SBarry Smith   ierr = VecRestoreArray(f,&af);CHKERRQ(ierr);
129347c6ae99SBarry Smith   ierr = VecRestoreArray(work,&awork);CHKERRQ(ierr);
12949a42bb27SBarry Smith   ierr = DMRestoreGlobalVector(da,&work);CHKERRQ(ierr);
129547c6ae99SBarry Smith 
129647c6ae99SBarry Smith   PetscFunctionReturn(0);
129747c6ae99SBarry Smith }
129847c6ae99SBarry Smith 
129947c6ae99SBarry Smith #undef __FUNCT__
13009a42bb27SBarry Smith #define __FUNCT__ "DMSetUp_DA_2D"
13019a42bb27SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMSetUp_DA_2D(DM da)
130247c6ae99SBarry Smith {
130347c6ae99SBarry Smith   DM_DA               *dd = (DM_DA*)da->data;
130447c6ae99SBarry Smith   const PetscInt       M            = dd->M;
130547c6ae99SBarry Smith   const PetscInt       N            = dd->N;
130647c6ae99SBarry Smith   PetscInt             m            = dd->m;
130747c6ae99SBarry Smith   PetscInt             n            = dd->n;
130847c6ae99SBarry Smith   const PetscInt       dof          = dd->w;
130947c6ae99SBarry Smith   const PetscInt       s            = dd->s;
1310aa219208SBarry Smith   const DMDAPeriodicType wrap         = dd->wrap;
1311aa219208SBarry Smith   const DMDAStencilType  stencil_type = dd->stencil_type;
131247c6ae99SBarry Smith   PetscInt            *lx           = dd->lx;
131347c6ae99SBarry Smith   PetscInt            *ly           = dd->ly;
131447c6ae99SBarry Smith   MPI_Comm             comm;
131547c6ae99SBarry Smith   PetscMPIInt          rank,size;
131647c6ae99SBarry Smith   PetscInt             xs,xe,ys,ye,x,y,Xs,Xe,Ys,Ye,start,end;
131747c6ae99SBarry Smith   PetscInt             up,down,left,i,n0,n1,n2,n3,n5,n6,n7,n8,*idx,nn;
131847c6ae99SBarry Smith   PetscInt             xbase,*bases,*ldims,j,x_t,y_t,s_t,base,count;
131947c6ae99SBarry Smith   PetscInt             s_x,s_y; /* s proportionalized to w */
132047c6ae99SBarry Smith   PetscInt             sn0 = 0,sn2 = 0,sn6 = 0,sn8 = 0;
132147c6ae99SBarry Smith   Vec                  local,global;
132247c6ae99SBarry Smith   VecScatter           ltog,gtol;
132347c6ae99SBarry Smith   IS                   to,from;
132447c6ae99SBarry Smith   PetscErrorCode       ierr;
132547c6ae99SBarry Smith 
132647c6ae99SBarry Smith   PetscFunctionBegin;
132747c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
132847c6ae99SBarry Smith 
132947c6ae99SBarry Smith   if (dof < 1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Must have 1 or more degrees of freedom per node: %D",dof);
133047c6ae99SBarry Smith   if (s < 0) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Stencil width cannot be negative: %D",s);
133147c6ae99SBarry Smith 
133247c6ae99SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
133347c6ae99SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
133447c6ae99SBarry Smith 
1335264ace61SBarry Smith   da->ops->getelements = DMGetElements_DA_2d_P1;
133647c6ae99SBarry Smith 
133747c6ae99SBarry Smith   dd->dim         = 2;
1338aa219208SBarry Smith   dd->elementtype = DMDA_ELEMENT_P1;
133947c6ae99SBarry Smith   ierr = PetscMalloc(dof*sizeof(char*),&dd->fieldname);CHKERRQ(ierr);
134047c6ae99SBarry Smith   ierr = PetscMemzero(dd->fieldname,dof*sizeof(char*));CHKERRQ(ierr);
134147c6ae99SBarry Smith 
134247c6ae99SBarry Smith   if (m != PETSC_DECIDE) {
134347c6ae99SBarry Smith     if (m < 1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in X direction: %D",m);
134447c6ae99SBarry Smith     else if (m > size) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in X direction: %D %d",m,size);
134547c6ae99SBarry Smith   }
134647c6ae99SBarry Smith   if (n != PETSC_DECIDE) {
134747c6ae99SBarry Smith     if (n < 1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Y direction: %D",n);
134847c6ae99SBarry Smith     else if (n > size) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Y direction: %D %d",n,size);
134947c6ae99SBarry Smith   }
135047c6ae99SBarry Smith 
135147c6ae99SBarry Smith   if (m == PETSC_DECIDE || n == PETSC_DECIDE) {
135247c6ae99SBarry Smith     if (n != PETSC_DECIDE) {
135347c6ae99SBarry Smith       m = size/n;
135447c6ae99SBarry Smith     } else if (m != PETSC_DECIDE) {
135547c6ae99SBarry Smith       n = size/m;
135647c6ae99SBarry Smith     } else {
135747c6ae99SBarry Smith       /* try for squarish distribution */
135847c6ae99SBarry Smith       m = (PetscInt)(0.5 + sqrt(((double)M)*((double)size)/((double)N)));
135947c6ae99SBarry Smith       if (!m) m = 1;
136047c6ae99SBarry Smith       while (m > 0) {
136147c6ae99SBarry Smith 	n = size/m;
136247c6ae99SBarry Smith 	if (m*n == size) break;
136347c6ae99SBarry Smith 	m--;
136447c6ae99SBarry Smith       }
136547c6ae99SBarry Smith       if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;}
136647c6ae99SBarry Smith     }
136747c6ae99SBarry 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 ");
136847c6ae99SBarry Smith   } else if (m*n != size) SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"Given Bad partition");
136947c6ae99SBarry Smith 
137047c6ae99SBarry Smith   if (M < m) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in x direction is too fine! %D %D",M,m);
137147c6ae99SBarry Smith   if (N < n) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in y direction is too fine! %D %D",N,n);
137247c6ae99SBarry Smith 
137347c6ae99SBarry Smith   /*
137447c6ae99SBarry Smith      Determine locally owned region
137547c6ae99SBarry Smith      xs is the first local node number, x is the number of local nodes
137647c6ae99SBarry Smith   */
137747c6ae99SBarry Smith   if (!lx) {
137847c6ae99SBarry Smith     ierr = PetscMalloc(m*sizeof(PetscInt), &dd->lx);CHKERRQ(ierr);
137947c6ae99SBarry Smith     lx = dd->lx;
138047c6ae99SBarry Smith     for (i=0; i<m; i++) {
138147c6ae99SBarry Smith       lx[i] = M/m + ((M % m) > i);
138247c6ae99SBarry Smith     }
138347c6ae99SBarry Smith   }
138447c6ae99SBarry Smith   x  = lx[rank % m];
138547c6ae99SBarry Smith   xs = 0;
138647c6ae99SBarry Smith   for (i=0; i<(rank % m); i++) {
138747c6ae99SBarry Smith     xs += lx[i];
138847c6ae99SBarry Smith   }
138947c6ae99SBarry Smith #if defined(PETSC_USE_DEBUG)
139047c6ae99SBarry Smith   left = xs;
139147c6ae99SBarry Smith   for (i=(rank % m); i<m; i++) {
139247c6ae99SBarry Smith     left += lx[i];
139347c6ae99SBarry Smith   }
139447c6ae99SBarry 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);
139547c6ae99SBarry Smith #endif
139647c6ae99SBarry Smith 
139747c6ae99SBarry Smith   /*
139847c6ae99SBarry Smith      Determine locally owned region
139947c6ae99SBarry Smith      ys is the first local node number, y is the number of local nodes
140047c6ae99SBarry Smith   */
140147c6ae99SBarry Smith   if (!ly) {
140247c6ae99SBarry Smith     ierr = PetscMalloc(n*sizeof(PetscInt), &dd->ly);CHKERRQ(ierr);
140347c6ae99SBarry Smith     ly = dd->ly;
140447c6ae99SBarry Smith     for (i=0; i<n; i++) {
140547c6ae99SBarry Smith       ly[i] = N/n + ((N % n) > i);
140647c6ae99SBarry Smith     }
140747c6ae99SBarry Smith   }
140847c6ae99SBarry Smith   y  = ly[rank/m];
140947c6ae99SBarry Smith   ys = 0;
141047c6ae99SBarry Smith   for (i=0; i<(rank/m); i++) {
141147c6ae99SBarry Smith     ys += ly[i];
141247c6ae99SBarry Smith   }
141347c6ae99SBarry Smith #if defined(PETSC_USE_DEBUG)
141447c6ae99SBarry Smith   left = ys;
141547c6ae99SBarry Smith   for (i=(rank/m); i<n; i++) {
141647c6ae99SBarry Smith     left += ly[i];
141747c6ae99SBarry Smith   }
141847c6ae99SBarry 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);
141947c6ae99SBarry Smith #endif
142047c6ae99SBarry Smith 
142147c6ae99SBarry 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);
142247c6ae99SBarry 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);
142347c6ae99SBarry Smith   xe = xs + x;
142447c6ae99SBarry Smith   ye = ys + y;
142547c6ae99SBarry Smith 
142647c6ae99SBarry Smith   /* determine ghost region */
142747c6ae99SBarry Smith   /* Assume No Periodicity */
142847c6ae99SBarry Smith   if (xs-s > 0) Xs = xs - s; else Xs = 0;
142947c6ae99SBarry Smith   if (ys-s > 0) Ys = ys - s; else Ys = 0;
143047c6ae99SBarry Smith   if (xe+s <= M) Xe = xe + s; else Xe = M;
143147c6ae99SBarry Smith   if (ye+s <= N) Ye = ye + s; else Ye = N;
143247c6ae99SBarry Smith 
143347c6ae99SBarry Smith   /* X Periodic */
1434aa219208SBarry Smith   if (DMDAXPeriodic(wrap)){
143547c6ae99SBarry Smith     Xs = xs - s;
143647c6ae99SBarry Smith     Xe = xe + s;
143747c6ae99SBarry Smith   }
143847c6ae99SBarry Smith 
143947c6ae99SBarry Smith   /* Y Periodic */
1440aa219208SBarry Smith   if (DMDAYPeriodic(wrap)){
144147c6ae99SBarry Smith     Ys = ys - s;
144247c6ae99SBarry Smith     Ye = ye + s;
144347c6ae99SBarry Smith   }
144447c6ae99SBarry Smith 
144547c6ae99SBarry Smith   /* Resize all X parameters to reflect w */
144647c6ae99SBarry Smith   x   *= dof;
144747c6ae99SBarry Smith   xs  *= dof;
144847c6ae99SBarry Smith   xe  *= dof;
144947c6ae99SBarry Smith   Xs  *= dof;
145047c6ae99SBarry Smith   Xe  *= dof;
145147c6ae99SBarry Smith   s_x = s*dof;
145247c6ae99SBarry Smith   s_y = s;
145347c6ae99SBarry Smith 
145447c6ae99SBarry Smith   /* determine starting point of each processor */
145547c6ae99SBarry Smith   nn    = x*y;
145647c6ae99SBarry Smith   ierr  = PetscMalloc2(size+1,PetscInt,&bases,size,PetscInt,&ldims);CHKERRQ(ierr);
145747c6ae99SBarry Smith   ierr  = MPI_Allgather(&nn,1,MPIU_INT,ldims,1,MPIU_INT,comm);CHKERRQ(ierr);
145847c6ae99SBarry Smith   bases[0] = 0;
145947c6ae99SBarry Smith   for (i=1; i<=size; i++) {
146047c6ae99SBarry Smith     bases[i] = ldims[i-1];
146147c6ae99SBarry Smith   }
146247c6ae99SBarry Smith   for (i=1; i<=size; i++) {
146347c6ae99SBarry Smith     bases[i] += bases[i-1];
146447c6ae99SBarry Smith   }
146547c6ae99SBarry Smith 
146647c6ae99SBarry Smith   /* allocate the base parallel and sequential vectors */
146747c6ae99SBarry Smith   dd->Nlocal = x*y;
146847c6ae99SBarry Smith   ierr = VecCreateMPIWithArray(comm,dd->Nlocal,PETSC_DECIDE,0,&global);CHKERRQ(ierr);
146947c6ae99SBarry Smith   ierr = VecSetBlockSize(global,dof);CHKERRQ(ierr);
147047c6ae99SBarry Smith   dd->nlocal = (Xe-Xs)*(Ye-Ys);
147147c6ae99SBarry Smith   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,dd->nlocal,0,&local);CHKERRQ(ierr);
147247c6ae99SBarry Smith   ierr = VecSetBlockSize(local,dof);CHKERRQ(ierr);
147347c6ae99SBarry Smith 
147447c6ae99SBarry Smith 
147547c6ae99SBarry Smith   /* generate appropriate vector scatters */
147647c6ae99SBarry Smith   /* local to global inserts non-ghost point region into global */
147747c6ae99SBarry Smith   ierr = VecGetOwnershipRange(global,&start,&end);CHKERRQ(ierr);
147847c6ae99SBarry Smith   ierr = ISCreateStride(comm,x*y,start,1,&to);CHKERRQ(ierr);
147947c6ae99SBarry Smith 
148047c6ae99SBarry Smith   left  = xs - Xs; down  = ys - Ys; up    = down + y;
148147c6ae99SBarry Smith   ierr = PetscMalloc(x*(up - down)*sizeof(PetscInt),&idx);CHKERRQ(ierr);
148247c6ae99SBarry Smith   count = 0;
148347c6ae99SBarry Smith   for (i=down; i<up; i++) {
148447c6ae99SBarry Smith     for (j=0; j<x/dof; j++) {
148547c6ae99SBarry Smith       idx[count++] = (left + i*(Xe-Xs) + j*dof)/dof;
148647c6ae99SBarry Smith     }
148747c6ae99SBarry Smith   }
148847c6ae99SBarry Smith   ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&from);CHKERRQ(ierr);
148947c6ae99SBarry Smith 
149047c6ae99SBarry Smith   ierr = VecScatterCreate(local,from,global,to,&ltog);CHKERRQ(ierr);
149147c6ae99SBarry Smith   ierr = PetscLogObjectParent(dd,ltog);CHKERRQ(ierr);
149247c6ae99SBarry Smith   ierr = ISDestroy(from);CHKERRQ(ierr);
149347c6ae99SBarry Smith   ierr = ISDestroy(to);CHKERRQ(ierr);
149447c6ae99SBarry Smith 
149547c6ae99SBarry Smith   /* global to local must include ghost points */
1496aa219208SBarry Smith   if (stencil_type == DMDA_STENCIL_BOX) {
149747c6ae99SBarry Smith     ierr = ISCreateStride(comm,(Xe-Xs)*(Ye-Ys),0,1,&to);CHKERRQ(ierr);
149847c6ae99SBarry Smith   } else {
149947c6ae99SBarry Smith     /* must drop into cross shape region */
150047c6ae99SBarry Smith     /*       ---------|
150147c6ae99SBarry Smith             |  top    |
150247c6ae99SBarry Smith          |---         ---|
150347c6ae99SBarry Smith          |   middle      |
150447c6ae99SBarry Smith          |               |
150547c6ae99SBarry Smith          ----         ----
150647c6ae99SBarry Smith             | bottom  |
150747c6ae99SBarry Smith             -----------
150847c6ae99SBarry Smith         Xs xs        xe  Xe */
150947c6ae99SBarry Smith     /* bottom */
151047c6ae99SBarry Smith     left  = xs - Xs; down = ys - Ys; up    = down + y;
151147c6ae99SBarry Smith     count = down*(xe-xs) + (up-down)*(Xe-Xs) + (Ye-Ys-up)*(xe-xs);
151247c6ae99SBarry Smith     ierr  = PetscMalloc(count*sizeof(PetscInt)/dof,&idx);CHKERRQ(ierr);
151347c6ae99SBarry Smith     count = 0;
151447c6ae99SBarry Smith     for (i=0; i<down; i++) {
151547c6ae99SBarry Smith       for (j=0; j<xe-xs; j += dof) {
151647c6ae99SBarry Smith         idx[count++] = left + i*(Xe-Xs) + j;
151747c6ae99SBarry Smith       }
151847c6ae99SBarry Smith     }
151947c6ae99SBarry Smith     /* middle */
152047c6ae99SBarry Smith     for (i=down; i<up; i++) {
152147c6ae99SBarry Smith       for (j=0; j<Xe-Xs; j += dof) {
152247c6ae99SBarry Smith         idx[count++] = i*(Xe-Xs) + j;
152347c6ae99SBarry Smith       }
152447c6ae99SBarry Smith     }
152547c6ae99SBarry Smith     /* top */
152647c6ae99SBarry Smith     for (i=up; i<Ye-Ys; i++) {
152747c6ae99SBarry Smith       for (j=0; j<xe-xs; j += dof) {
152847c6ae99SBarry Smith         idx[count++] = (left + i*(Xe-Xs) + j);
152947c6ae99SBarry Smith       }
153047c6ae99SBarry Smith     }
153147c6ae99SBarry Smith     for (i=0; i<count; i++) idx[i] = idx[i]/dof;
153247c6ae99SBarry Smith     ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr);
153347c6ae99SBarry Smith   }
153447c6ae99SBarry Smith 
153547c6ae99SBarry Smith 
153647c6ae99SBarry Smith   /* determine who lies on each side of us stored in    n6 n7 n8
153747c6ae99SBarry Smith                                                         n3    n5
153847c6ae99SBarry Smith                                                         n0 n1 n2
153947c6ae99SBarry Smith   */
154047c6ae99SBarry Smith 
154147c6ae99SBarry Smith   /* Assume the Non-Periodic Case */
154247c6ae99SBarry Smith   n1 = rank - m;
154347c6ae99SBarry Smith   if (rank % m) {
154447c6ae99SBarry Smith     n0 = n1 - 1;
154547c6ae99SBarry Smith   } else {
154647c6ae99SBarry Smith     n0 = -1;
154747c6ae99SBarry Smith   }
154847c6ae99SBarry Smith   if ((rank+1) % m) {
154947c6ae99SBarry Smith     n2 = n1 + 1;
155047c6ae99SBarry Smith     n5 = rank + 1;
155147c6ae99SBarry Smith     n8 = rank + m + 1; if (n8 >= m*n) n8 = -1;
155247c6ae99SBarry Smith   } else {
155347c6ae99SBarry Smith     n2 = -1; n5 = -1; n8 = -1;
155447c6ae99SBarry Smith   }
155547c6ae99SBarry Smith   if (rank % m) {
155647c6ae99SBarry Smith     n3 = rank - 1;
155747c6ae99SBarry Smith     n6 = n3 + m; if (n6 >= m*n) n6 = -1;
155847c6ae99SBarry Smith   } else {
155947c6ae99SBarry Smith     n3 = -1; n6 = -1;
156047c6ae99SBarry Smith   }
156147c6ae99SBarry Smith   n7 = rank + m; if (n7 >= m*n) n7 = -1;
156247c6ae99SBarry Smith 
156347c6ae99SBarry Smith 
156447c6ae99SBarry Smith   /* Modify for Periodic Cases */
1565aa219208SBarry Smith   if (wrap == DMDA_YPERIODIC) {  /* Handle Top and Bottom Sides */
156647c6ae99SBarry Smith     if (n1 < 0) n1 = rank + m * (n-1);
156747c6ae99SBarry Smith     if (n7 < 0) n7 = rank - m * (n-1);
156847c6ae99SBarry Smith     if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1;
156947c6ae99SBarry Smith     if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1;
157047c6ae99SBarry Smith     if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1;
157147c6ae99SBarry Smith     if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1;
1572aa219208SBarry Smith   } else if (wrap == DMDA_XPERIODIC) { /* Handle Left and Right Sides */
157347c6ae99SBarry Smith     if (n3 < 0) n3 = rank + (m-1);
157447c6ae99SBarry Smith     if (n5 < 0) n5 = rank - (m-1);
157547c6ae99SBarry Smith     if ((n1 >= 0) && (n0 < 0)) n0 = rank-1;
157647c6ae99SBarry Smith     if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1;
157747c6ae99SBarry Smith     if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1;
157847c6ae99SBarry Smith     if ((n7 >= 0) && (n8 < 0)) n8 = rank+1;
1579aa219208SBarry Smith   } else if (wrap == DMDA_XYPERIODIC) {
158047c6ae99SBarry Smith 
158147c6ae99SBarry Smith     /* Handle all four corners */
158247c6ae99SBarry Smith     if ((n6 < 0) && (n7 < 0) && (n3 < 0)) n6 = m-1;
158347c6ae99SBarry Smith     if ((n8 < 0) && (n7 < 0) && (n5 < 0)) n8 = 0;
158447c6ae99SBarry Smith     if ((n2 < 0) && (n5 < 0) && (n1 < 0)) n2 = size-m;
158547c6ae99SBarry Smith     if ((n0 < 0) && (n3 < 0) && (n1 < 0)) n0 = size-1;
158647c6ae99SBarry Smith 
158747c6ae99SBarry Smith     /* Handle Top and Bottom Sides */
158847c6ae99SBarry Smith     if (n1 < 0) n1 = rank + m * (n-1);
158947c6ae99SBarry Smith     if (n7 < 0) n7 = rank - m * (n-1);
159047c6ae99SBarry Smith     if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1;
159147c6ae99SBarry Smith     if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1;
159247c6ae99SBarry Smith     if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1;
159347c6ae99SBarry Smith     if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1;
159447c6ae99SBarry Smith 
159547c6ae99SBarry Smith     /* Handle Left and Right Sides */
159647c6ae99SBarry Smith     if (n3 < 0) n3 = rank + (m-1);
159747c6ae99SBarry Smith     if (n5 < 0) n5 = rank - (m-1);
159847c6ae99SBarry Smith     if ((n1 >= 0) && (n0 < 0)) n0 = rank-1;
159947c6ae99SBarry Smith     if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1;
160047c6ae99SBarry Smith     if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1;
160147c6ae99SBarry Smith     if ((n7 >= 0) && (n8 < 0)) n8 = rank+1;
160247c6ae99SBarry Smith   }
160347c6ae99SBarry Smith   ierr = PetscMalloc(9*sizeof(PetscInt),&dd->neighbors);CHKERRQ(ierr);
160447c6ae99SBarry Smith   dd->neighbors[0] = n0;
160547c6ae99SBarry Smith   dd->neighbors[1] = n1;
160647c6ae99SBarry Smith   dd->neighbors[2] = n2;
160747c6ae99SBarry Smith   dd->neighbors[3] = n3;
160847c6ae99SBarry Smith   dd->neighbors[4] = rank;
160947c6ae99SBarry Smith   dd->neighbors[5] = n5;
161047c6ae99SBarry Smith   dd->neighbors[6] = n6;
161147c6ae99SBarry Smith   dd->neighbors[7] = n7;
161247c6ae99SBarry Smith   dd->neighbors[8] = n8;
161347c6ae99SBarry Smith 
1614aa219208SBarry Smith   if (stencil_type == DMDA_STENCIL_STAR) {
161547c6ae99SBarry Smith     /* save corner processor numbers */
161647c6ae99SBarry Smith     sn0 = n0; sn2 = n2; sn6 = n6; sn8 = n8;
161747c6ae99SBarry Smith     n0 = n2 = n6 = n8 = -1;
161847c6ae99SBarry Smith   }
161947c6ae99SBarry Smith 
162047c6ae99SBarry Smith   ierr = PetscMalloc((x+2*s_x)*(y+2*s_y)*sizeof(PetscInt),&idx);CHKERRQ(ierr);
162147c6ae99SBarry Smith   ierr = PetscLogObjectMemory(da,(x+2*s_x)*(y+2*s_y)*sizeof(PetscInt));CHKERRQ(ierr);
162247c6ae99SBarry Smith   nn = 0;
162347c6ae99SBarry Smith 
162447c6ae99SBarry Smith   xbase = bases[rank];
162547c6ae99SBarry Smith   for (i=1; i<=s_y; i++) {
162647c6ae99SBarry Smith     if (n0 >= 0) { /* left below */
162747c6ae99SBarry Smith       x_t = lx[n0 % m]*dof;
162847c6ae99SBarry Smith       y_t = ly[(n0/m)];
162947c6ae99SBarry Smith       s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x;
163047c6ae99SBarry Smith       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
163147c6ae99SBarry Smith     }
163247c6ae99SBarry Smith     if (n1 >= 0) { /* directly below */
163347c6ae99SBarry Smith       x_t = x;
163447c6ae99SBarry Smith       y_t = ly[(n1/m)];
163547c6ae99SBarry Smith       s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t;
163647c6ae99SBarry Smith       for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
163747c6ae99SBarry Smith     }
163847c6ae99SBarry Smith     if (n2 >= 0) { /* right below */
163947c6ae99SBarry Smith       x_t = lx[n2 % m]*dof;
164047c6ae99SBarry Smith       y_t = ly[(n2/m)];
164147c6ae99SBarry Smith       s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t;
164247c6ae99SBarry Smith       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
164347c6ae99SBarry Smith     }
164447c6ae99SBarry Smith   }
164547c6ae99SBarry Smith 
164647c6ae99SBarry Smith   for (i=0; i<y; i++) {
164747c6ae99SBarry Smith     if (n3 >= 0) { /* directly left */
164847c6ae99SBarry Smith       x_t = lx[n3 % m]*dof;
164947c6ae99SBarry Smith       /* y_t = y; */
165047c6ae99SBarry Smith       s_t = bases[n3] + (i+1)*x_t - s_x;
165147c6ae99SBarry Smith       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
165247c6ae99SBarry Smith     }
165347c6ae99SBarry Smith 
165447c6ae99SBarry Smith     for (j=0; j<x; j++) { idx[nn++] = xbase++; } /* interior */
165547c6ae99SBarry Smith 
165647c6ae99SBarry Smith     if (n5 >= 0) { /* directly right */
165747c6ae99SBarry Smith       x_t = lx[n5 % m]*dof;
165847c6ae99SBarry Smith       /* y_t = y; */
165947c6ae99SBarry Smith       s_t = bases[n5] + (i)*x_t;
166047c6ae99SBarry Smith       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
166147c6ae99SBarry Smith     }
166247c6ae99SBarry Smith   }
166347c6ae99SBarry Smith 
166447c6ae99SBarry Smith   for (i=1; i<=s_y; i++) {
166547c6ae99SBarry Smith     if (n6 >= 0) { /* left above */
166647c6ae99SBarry Smith       x_t = lx[n6 % m]*dof;
166747c6ae99SBarry Smith       /* y_t = ly[(n6/m)]; */
166847c6ae99SBarry Smith       s_t = bases[n6] + (i)*x_t - s_x;
166947c6ae99SBarry Smith       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
167047c6ae99SBarry Smith     }
167147c6ae99SBarry Smith     if (n7 >= 0) { /* directly above */
167247c6ae99SBarry Smith       x_t = x;
167347c6ae99SBarry Smith       /* y_t = ly[(n7/m)]; */
167447c6ae99SBarry Smith       s_t = bases[n7] + (i-1)*x_t;
167547c6ae99SBarry Smith       for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
167647c6ae99SBarry Smith     }
167747c6ae99SBarry Smith     if (n8 >= 0) { /* right above */
167847c6ae99SBarry Smith       x_t = lx[n8 % m]*dof;
167947c6ae99SBarry Smith       /* y_t = ly[(n8/m)]; */
168047c6ae99SBarry Smith       s_t = bases[n8] + (i-1)*x_t;
168147c6ae99SBarry Smith       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
168247c6ae99SBarry Smith     }
168347c6ae99SBarry Smith   }
168447c6ae99SBarry Smith 
168547c6ae99SBarry Smith   base = bases[rank];
168647c6ae99SBarry Smith   {
168747c6ae99SBarry Smith     PetscInt nnn = nn/dof,*iidx;
168847c6ae99SBarry Smith     ierr = PetscMalloc(nnn*sizeof(PetscInt),&iidx);CHKERRQ(ierr);
168947c6ae99SBarry Smith     for (i=0; i<nnn; i++) {
169047c6ae99SBarry Smith       iidx[i] = idx[dof*i]/dof;
169147c6ae99SBarry Smith     }
169247c6ae99SBarry Smith     ierr = ISCreateBlock(comm,dof,nnn,iidx,PETSC_OWN_POINTER,&from);CHKERRQ(ierr);
169347c6ae99SBarry Smith   }
169447c6ae99SBarry Smith   ierr = VecScatterCreate(global,from,local,to,&gtol);CHKERRQ(ierr);
169547c6ae99SBarry Smith   ierr = PetscLogObjectParent(da,gtol);CHKERRQ(ierr);
169647c6ae99SBarry Smith   ierr = ISDestroy(to);CHKERRQ(ierr);
169747c6ae99SBarry Smith   ierr = ISDestroy(from);CHKERRQ(ierr);
169847c6ae99SBarry Smith 
1699aa219208SBarry Smith   if (stencil_type == DMDA_STENCIL_STAR) {
170047c6ae99SBarry Smith     /*
170147c6ae99SBarry Smith         Recompute the local to global mappings, this time keeping the
170247c6ae99SBarry Smith       information about the cross corner processor numbers.
170347c6ae99SBarry Smith     */
170447c6ae99SBarry Smith     n0 = sn0; n2 = sn2; n6 = sn6; n8 = sn8;
170547c6ae99SBarry Smith     nn = 0;
170647c6ae99SBarry Smith     xbase = bases[rank];
170747c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
170847c6ae99SBarry Smith       if (n0 >= 0) { /* left below */
170947c6ae99SBarry Smith         x_t = lx[n0 % m]*dof;
171047c6ae99SBarry Smith         y_t = ly[(n0/m)];
171147c6ae99SBarry Smith         s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x;
171247c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
171347c6ae99SBarry Smith       }
171447c6ae99SBarry Smith       if (n1 >= 0) { /* directly below */
171547c6ae99SBarry Smith         x_t = x;
171647c6ae99SBarry Smith         y_t = ly[(n1/m)];
171747c6ae99SBarry Smith         s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t;
171847c6ae99SBarry Smith         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
171947c6ae99SBarry Smith       }
172047c6ae99SBarry Smith       if (n2 >= 0) { /* right below */
172147c6ae99SBarry Smith         x_t = lx[n2 % m]*dof;
172247c6ae99SBarry Smith         y_t = ly[(n2/m)];
172347c6ae99SBarry Smith         s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t;
172447c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
172547c6ae99SBarry Smith       }
172647c6ae99SBarry Smith     }
172747c6ae99SBarry Smith 
172847c6ae99SBarry Smith     for (i=0; i<y; i++) {
172947c6ae99SBarry Smith       if (n3 >= 0) { /* directly left */
173047c6ae99SBarry Smith         x_t = lx[n3 % m]*dof;
173147c6ae99SBarry Smith         /* y_t = y; */
173247c6ae99SBarry Smith         s_t = bases[n3] + (i+1)*x_t - s_x;
173347c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
173447c6ae99SBarry Smith       }
173547c6ae99SBarry Smith 
173647c6ae99SBarry Smith       for (j=0; j<x; j++) { idx[nn++] = xbase++; } /* interior */
173747c6ae99SBarry Smith 
173847c6ae99SBarry Smith       if (n5 >= 0) { /* directly right */
173947c6ae99SBarry Smith         x_t = lx[n5 % m]*dof;
174047c6ae99SBarry Smith         /* y_t = y; */
174147c6ae99SBarry Smith         s_t = bases[n5] + (i)*x_t;
174247c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
174347c6ae99SBarry Smith       }
174447c6ae99SBarry Smith     }
174547c6ae99SBarry Smith 
174647c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
174747c6ae99SBarry Smith       if (n6 >= 0) { /* left above */
174847c6ae99SBarry Smith         x_t = lx[n6 % m]*dof;
174947c6ae99SBarry Smith         /* y_t = ly[(n6/m)]; */
175047c6ae99SBarry Smith         s_t = bases[n6] + (i)*x_t - s_x;
175147c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
175247c6ae99SBarry Smith       }
175347c6ae99SBarry Smith       if (n7 >= 0) { /* directly above */
175447c6ae99SBarry Smith         x_t = x;
175547c6ae99SBarry Smith         /* y_t = ly[(n7/m)]; */
175647c6ae99SBarry Smith         s_t = bases[n7] + (i-1)*x_t;
175747c6ae99SBarry Smith         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
175847c6ae99SBarry Smith       }
175947c6ae99SBarry Smith       if (n8 >= 0) { /* right above */
176047c6ae99SBarry Smith         x_t = lx[n8 % m]*dof;
176147c6ae99SBarry Smith         /* y_t = ly[(n8/m)]; */
176247c6ae99SBarry Smith         s_t = bases[n8] + (i-1)*x_t;
176347c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
176447c6ae99SBarry Smith       }
176547c6ae99SBarry Smith     }
176647c6ae99SBarry Smith   }
176747c6ae99SBarry Smith   ierr = PetscFree2(bases,ldims);CHKERRQ(ierr);
176847c6ae99SBarry Smith 
176947c6ae99SBarry Smith   dd->m  = m;  dd->n  = n;
177047c6ae99SBarry Smith   dd->xs = xs; dd->xe = xe; dd->ys = ys; dd->ye = ye; dd->zs = 0; dd->ze = 1;
177147c6ae99SBarry Smith   dd->Xs = Xs; dd->Xe = Xe; dd->Ys = Ys; dd->Ye = Ye; dd->Zs = 0; dd->Ze = 1;
177247c6ae99SBarry Smith 
177347c6ae99SBarry Smith   ierr = VecDestroy(local);CHKERRQ(ierr);
177447c6ae99SBarry Smith   ierr = VecDestroy(global);CHKERRQ(ierr);
177547c6ae99SBarry Smith 
177647c6ae99SBarry Smith   dd->gtol      = gtol;
177747c6ae99SBarry Smith   dd->ltog      = ltog;
177847c6ae99SBarry Smith   dd->idx       = idx;
177947c6ae99SBarry Smith   dd->Nl        = nn;
178047c6ae99SBarry Smith   dd->base      = base;
17819a42bb27SBarry Smith   da->ops->view = DMView_DA_2d;
178247c6ae99SBarry Smith 
178347c6ae99SBarry Smith   /*
178447c6ae99SBarry Smith      Set the local to global ordering in the global vector, this allows use
178547c6ae99SBarry Smith      of VecSetValuesLocal().
178647c6ae99SBarry Smith   */
178747c6ae99SBarry Smith   ierr = ISLocalToGlobalMappingCreate(comm,nn,idx,PETSC_OWN_POINTER,&dd->ltogmap);CHKERRQ(ierr);
178847c6ae99SBarry Smith   ierr = ISLocalToGlobalMappingBlock(dd->ltogmap,dd->w,&dd->ltogmapb);CHKERRQ(ierr);
178947c6ae99SBarry Smith   ierr = PetscLogObjectParent(da,dd->ltogmap);CHKERRQ(ierr);
179047c6ae99SBarry Smith 
179147c6ae99SBarry Smith   dd->ltol = PETSC_NULL;
179247c6ae99SBarry Smith   dd->ao   = PETSC_NULL;
179347c6ae99SBarry Smith 
179447c6ae99SBarry Smith   PetscFunctionReturn(0);
179547c6ae99SBarry Smith }
179647c6ae99SBarry Smith 
179747c6ae99SBarry Smith #undef __FUNCT__
1798aa219208SBarry Smith #define __FUNCT__ "DMDACreate2d"
179947c6ae99SBarry Smith /*@C
1800aa219208SBarry Smith    DMDACreate2d -  Creates an object that will manage the communication of  two-dimensional
180147c6ae99SBarry Smith    regular array data that is distributed across some processors.
180247c6ae99SBarry Smith 
180347c6ae99SBarry Smith    Collective on MPI_Comm
180447c6ae99SBarry Smith 
180547c6ae99SBarry Smith    Input Parameters:
180647c6ae99SBarry Smith +  comm - MPI communicator
180747c6ae99SBarry Smith .  wrap - type of periodicity should the array have.
1808aa219208SBarry Smith          Use one of DMDA_NONPERIODIC, DMDA_XPERIODIC, DMDA_YPERIODIC, or DMDA_XYPERIODIC.
1809aa219208SBarry Smith .  stencil_type - stencil type.  Use either DMDA_STENCIL_BOX or DMDA_STENCIL_STAR.
181047c6ae99SBarry 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
181147c6ae99SBarry Smith             from the command line with -da_grid_x <M> -da_grid_y <N>)
181247c6ae99SBarry Smith .  m,n - corresponding number of processors in each dimension
181347c6ae99SBarry Smith          (or PETSC_DECIDE to have calculated)
181447c6ae99SBarry Smith .  dof - number of degrees of freedom per node
181547c6ae99SBarry Smith .  s - stencil width
181647c6ae99SBarry Smith -  lx, ly - arrays containing the number of nodes in each cell along
181747c6ae99SBarry Smith            the x and y coordinates, or PETSC_NULL. If non-null, these
181847c6ae99SBarry Smith            must be of length as m and n, and the corresponding
181947c6ae99SBarry Smith            m and n cannot be PETSC_DECIDE. The sum of the lx[] entries
182047c6ae99SBarry Smith            must be M, and the sum of the ly[] entries must be N.
182147c6ae99SBarry Smith 
182247c6ae99SBarry Smith    Output Parameter:
182347c6ae99SBarry Smith .  da - the resulting distributed array object
182447c6ae99SBarry Smith 
182547c6ae99SBarry Smith    Options Database Key:
1826aa219208SBarry Smith +  -da_view - Calls DMView() at the conclusion of DMDACreate2d()
182747c6ae99SBarry Smith .  -da_grid_x <nx> - number of grid points in x direction, if M < 0
182847c6ae99SBarry Smith .  -da_grid_y <ny> - number of grid points in y direction, if N < 0
182947c6ae99SBarry Smith .  -da_processors_x <nx> - number of processors in x direction
183047c6ae99SBarry Smith .  -da_processors_y <ny> - number of processors in y direction
183147c6ae99SBarry Smith .  -da_refine_x - refinement ratio in x direction
183247c6ae99SBarry Smith -  -da_refine_y - refinement ratio in y direction
183347c6ae99SBarry Smith 
183447c6ae99SBarry Smith    Level: beginner
183547c6ae99SBarry Smith 
183647c6ae99SBarry Smith    Notes:
1837aa219208SBarry Smith    The stencil type DMDA_STENCIL_STAR with width 1 corresponds to the
1838aa219208SBarry Smith    standard 5-pt stencil, while DMDA_STENCIL_BOX with width 1 denotes
183947c6ae99SBarry Smith    the standard 9-pt stencil.
184047c6ae99SBarry Smith 
1841aa219208SBarry Smith    The array data itself is NOT stored in the DMDA, it is stored in Vec objects;
1842564755cdSBarry Smith    The appropriate vector objects can be obtained with calls to DMCreateGlobalVector()
1843564755cdSBarry Smith    and DMCreateLocalVector() and calls to VecDuplicate() if more are needed.
184447c6ae99SBarry Smith 
184547c6ae99SBarry Smith .keywords: distributed array, create, two-dimensional
184647c6ae99SBarry Smith 
1847aa219208SBarry Smith .seealso: DMDestroy(), DMView(), DMDACreate1d(), DMDACreate3d(), DMGlobalToLocalBegin(), DMDAGetRefinementFactor(),
1848aa219208SBarry Smith           DMGlobalToLocalEnd(), DMLocalToGlobalBegin(), DMDALocalToLocalBegin(), DMDALocalToLocalEnd(), DMDASetRefinementFactor(),
1849aa219208SBarry Smith           DMDAGetInfo(), DMCreateGlobalVector(), DMCreateLocalVector(), DMDACreateNaturalVector(), DMDALoad(), DMDAGetOwnershipRanges()
185047c6ae99SBarry Smith 
185147c6ae99SBarry Smith @*/
1852aa219208SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDACreate2d(MPI_Comm comm,DMDAPeriodicType wrap,DMDAStencilType stencil_type,
18539a42bb27SBarry Smith                           PetscInt M,PetscInt N,PetscInt m,PetscInt n,PetscInt dof,PetscInt s,const PetscInt lx[],const PetscInt ly[],DM *da)
185447c6ae99SBarry Smith {
185547c6ae99SBarry Smith   PetscErrorCode ierr;
185647c6ae99SBarry Smith 
185747c6ae99SBarry Smith   PetscFunctionBegin;
1858aa219208SBarry Smith   ierr = DMDACreate(comm, da);CHKERRQ(ierr);
1859aa219208SBarry Smith   ierr = DMDASetDim(*da, 2);CHKERRQ(ierr);
1860aa219208SBarry Smith   ierr = DMDASetSizes(*da, M, N, 1);CHKERRQ(ierr);
1861aa219208SBarry Smith   ierr = DMDASetNumProcs(*da, m, n, PETSC_DECIDE);CHKERRQ(ierr);
1862aa219208SBarry Smith   ierr = DMDASetPeriodicity(*da, wrap);CHKERRQ(ierr);
1863aa219208SBarry Smith   ierr = DMDASetDof(*da, dof);CHKERRQ(ierr);
1864aa219208SBarry Smith   ierr = DMDASetStencilType(*da, stencil_type);CHKERRQ(ierr);
1865aa219208SBarry Smith   ierr = DMDASetStencilWidth(*da, s);CHKERRQ(ierr);
1866aa219208SBarry Smith   ierr = DMDASetOwnershipRanges(*da, lx, ly, PETSC_NULL);CHKERRQ(ierr);
186747c6ae99SBarry Smith   /* This violates the behavior for other classes, but right now users expect negative dimensions to be handled this way */
18689a42bb27SBarry Smith   ierr = DMSetFromOptions(*da);CHKERRQ(ierr);
18699a42bb27SBarry Smith   ierr = DMSetUp(*da);CHKERRQ(ierr);
187047c6ae99SBarry Smith   PetscFunctionReturn(0);
187147c6ae99SBarry Smith }
1872