xref: /petsc/src/dm/impls/da/da2.c (revision 5dff015ad9e8ad435245d6dadb2d7c7a046f27b8)
147c6ae99SBarry Smith #define PETSCDM_DLL
247c6ae99SBarry Smith 
347c6ae99SBarry Smith #include "private/daimpl.h"    /*I   "petscda.h"   I*/
447c6ae99SBarry Smith 
547c6ae99SBarry Smith #undef __FUNCT__
647c6ae99SBarry Smith #define __FUNCT__ "DAView_2d"
747c6ae99SBarry Smith PetscErrorCode DAView_2d(DA da,PetscViewer viewer)
847c6ae99SBarry Smith {
947c6ae99SBarry Smith   PetscErrorCode ierr;
1047c6ae99SBarry Smith   PetscMPIInt    rank;
1147c6ae99SBarry Smith   PetscBool      iascii,isdraw;
1247c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
1347c6ae99SBarry Smith 
1447c6ae99SBarry Smith   PetscFunctionBegin;
1547c6ae99SBarry Smith   ierr = MPI_Comm_rank(((PetscObject)da)->comm,&rank);CHKERRQ(ierr);
1647c6ae99SBarry Smith 
1747c6ae99SBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1847c6ae99SBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
1947c6ae99SBarry Smith   if (iascii) {
2047c6ae99SBarry Smith     PetscViewerFormat format;
2147c6ae99SBarry Smith 
2247c6ae99SBarry Smith     ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr);
2347c6ae99SBarry Smith     if (format != PETSC_VIEWER_ASCII_VTK && format != PETSC_VIEWER_ASCII_VTK_CELL) {
2447c6ae99SBarry Smith       DALocalInfo info;
2547c6ae99SBarry Smith       ierr = DAGetLocalInfo(da,&info);CHKERRQ(ierr);
2647c6ae99SBarry 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);
2747c6ae99SBarry 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);
2847c6ae99SBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
2947c6ae99SBarry Smith     }
3047c6ae99SBarry Smith   } else if (isdraw) {
3147c6ae99SBarry Smith     PetscDraw       draw;
3247c6ae99SBarry Smith     double     ymin = -1*dd->s-1,ymax = dd->N+dd->s;
3347c6ae99SBarry Smith     double     xmin = -1*dd->s-1,xmax = dd->M+dd->s;
3447c6ae99SBarry Smith     double     x,y;
3547c6ae99SBarry Smith     PetscInt   base,*idx;
3647c6ae99SBarry Smith     char       node[10];
3747c6ae99SBarry Smith     PetscBool  isnull;
3847c6ae99SBarry Smith 
3947c6ae99SBarry Smith     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
4047c6ae99SBarry Smith     ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
4147c6ae99SBarry Smith     if (!dd->coordinates) {
4247c6ae99SBarry Smith       ierr = PetscDrawSetCoordinates(draw,xmin,ymin,xmax,ymax);CHKERRQ(ierr);
4347c6ae99SBarry Smith     }
4447c6ae99SBarry Smith     ierr = PetscDrawSynchronizedClear(draw);CHKERRQ(ierr);
4547c6ae99SBarry Smith 
4647c6ae99SBarry Smith     /* first processor draw all node lines */
4747c6ae99SBarry Smith     if (!rank) {
4847c6ae99SBarry Smith       ymin = 0.0; ymax = dd->N - 1;
4947c6ae99SBarry Smith       for (xmin=0; xmin<dd->M; xmin++) {
5047c6ae99SBarry Smith         ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_BLACK);CHKERRQ(ierr);
5147c6ae99SBarry Smith       }
5247c6ae99SBarry Smith       xmin = 0.0; xmax = dd->M - 1;
5347c6ae99SBarry Smith       for (ymin=0; ymin<dd->N; ymin++) {
5447c6ae99SBarry Smith         ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_BLACK);CHKERRQ(ierr);
5547c6ae99SBarry Smith       }
5647c6ae99SBarry Smith     }
5747c6ae99SBarry Smith     ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr);
5847c6ae99SBarry Smith     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
5947c6ae99SBarry Smith 
6047c6ae99SBarry Smith     /* draw my box */
6147c6ae99SBarry Smith     ymin = dd->ys; ymax = dd->ye - 1; xmin = dd->xs/dd->w;
6247c6ae99SBarry Smith     xmax =(dd->xe-1)/dd->w;
6347c6ae99SBarry Smith     ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_RED);CHKERRQ(ierr);
6447c6ae99SBarry Smith     ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_RED);CHKERRQ(ierr);
6547c6ae99SBarry Smith     ierr = PetscDrawLine(draw,xmin,ymax,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr);
6647c6ae99SBarry Smith     ierr = PetscDrawLine(draw,xmax,ymin,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr);
6747c6ae99SBarry Smith 
6847c6ae99SBarry Smith     /* put in numbers */
6947c6ae99SBarry Smith     base = (dd->base)/dd->w;
7047c6ae99SBarry Smith     for (y=ymin; y<=ymax; y++) {
7147c6ae99SBarry Smith       for (x=xmin; x<=xmax; x++) {
7247c6ae99SBarry Smith         sprintf(node,"%d",(int)base++);
7347c6ae99SBarry Smith         ierr = PetscDrawString(draw,x,y,PETSC_DRAW_BLACK,node);CHKERRQ(ierr);
7447c6ae99SBarry Smith       }
7547c6ae99SBarry Smith     }
7647c6ae99SBarry Smith 
7747c6ae99SBarry Smith     ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr);
7847c6ae99SBarry Smith     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
7947c6ae99SBarry Smith     /* overlay ghost numbers, useful for error checking */
8047c6ae99SBarry Smith     /* put in numbers */
8147c6ae99SBarry Smith 
8247c6ae99SBarry Smith     base = 0; idx = dd->idx;
8347c6ae99SBarry Smith     ymin = dd->Ys; ymax = dd->Ye; xmin = dd->Xs; xmax = dd->Xe;
8447c6ae99SBarry Smith     for (y=ymin; y<ymax; y++) {
8547c6ae99SBarry Smith       for (x=xmin; x<xmax; x++) {
8647c6ae99SBarry Smith         if ((base % dd->w) == 0) {
8747c6ae99SBarry Smith           sprintf(node,"%d",(int)(idx[base]/dd->w));
8847c6ae99SBarry Smith           ierr = PetscDrawString(draw,x/dd->w,y,PETSC_DRAW_BLUE,node);CHKERRQ(ierr);
8947c6ae99SBarry Smith         }
9047c6ae99SBarry Smith         base++;
9147c6ae99SBarry Smith       }
9247c6ae99SBarry Smith     }
9347c6ae99SBarry Smith     ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr);
9447c6ae99SBarry Smith     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
9547c6ae99SBarry Smith   } else {
9647c6ae99SBarry Smith     SETERRQ1(((PetscObject)viewer)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for DA2d",((PetscObject)viewer)->type_name);
9747c6ae99SBarry Smith   }
9847c6ae99SBarry Smith   PetscFunctionReturn(0);
9947c6ae99SBarry Smith }
10047c6ae99SBarry Smith 
10147c6ae99SBarry Smith /*
10247c6ae99SBarry Smith       M is number of grid points
10347c6ae99SBarry Smith       m is number of processors
10447c6ae99SBarry Smith 
10547c6ae99SBarry Smith */
10647c6ae99SBarry Smith #undef __FUNCT__
10747c6ae99SBarry Smith #define __FUNCT__ "DASplitComm2d"
10847c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DASplitComm2d(MPI_Comm comm,PetscInt M,PetscInt N,PetscInt sw,MPI_Comm *outcomm)
10947c6ae99SBarry Smith {
11047c6ae99SBarry Smith   PetscErrorCode ierr;
11147c6ae99SBarry Smith   PetscInt       m,n = 0,x = 0,y = 0;
11247c6ae99SBarry Smith   PetscMPIInt    size,csize,rank;
11347c6ae99SBarry Smith 
11447c6ae99SBarry Smith   PetscFunctionBegin;
11547c6ae99SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
11647c6ae99SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
11747c6ae99SBarry Smith 
11847c6ae99SBarry Smith   csize = 4*size;
11947c6ae99SBarry Smith   do {
12047c6ae99SBarry 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);
12147c6ae99SBarry Smith     csize   = csize/4;
12247c6ae99SBarry Smith 
12347c6ae99SBarry Smith     m = (PetscInt)(0.5 + sqrt(((double)M)*((double)csize)/((double)N)));
12447c6ae99SBarry Smith     if (!m) m = 1;
12547c6ae99SBarry Smith     while (m > 0) {
12647c6ae99SBarry Smith       n = csize/m;
12747c6ae99SBarry Smith       if (m*n == csize) break;
12847c6ae99SBarry Smith       m--;
12947c6ae99SBarry Smith     }
13047c6ae99SBarry Smith     if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;}
13147c6ae99SBarry Smith 
13247c6ae99SBarry Smith     x = M/m + ((M % m) > ((csize-1) % m));
13347c6ae99SBarry Smith     y = (N + (csize-1)/m)/n;
13447c6ae99SBarry Smith   } while ((x < 4 || y < 4) && csize > 1);
13547c6ae99SBarry Smith   if (size != csize) {
13647c6ae99SBarry Smith     MPI_Group    entire_group,sub_group;
13747c6ae99SBarry Smith     PetscMPIInt  i,*groupies;
13847c6ae99SBarry Smith 
13947c6ae99SBarry Smith     ierr = MPI_Comm_group(comm,&entire_group);CHKERRQ(ierr);
14047c6ae99SBarry Smith     ierr = PetscMalloc(csize*sizeof(PetscInt),&groupies);CHKERRQ(ierr);
14147c6ae99SBarry Smith     for (i=0; i<csize; i++) {
14247c6ae99SBarry Smith       groupies[i] = (rank/csize)*csize + i;
14347c6ae99SBarry Smith     }
14447c6ae99SBarry Smith     ierr = MPI_Group_incl(entire_group,csize,groupies,&sub_group);CHKERRQ(ierr);
14547c6ae99SBarry Smith     ierr = PetscFree(groupies);CHKERRQ(ierr);
14647c6ae99SBarry Smith     ierr = MPI_Comm_create(comm,sub_group,outcomm);CHKERRQ(ierr);
14747c6ae99SBarry Smith     ierr = MPI_Group_free(&entire_group);CHKERRQ(ierr);
14847c6ae99SBarry Smith     ierr = MPI_Group_free(&sub_group);CHKERRQ(ierr);
14947c6ae99SBarry Smith     ierr = PetscInfo1(0,"DASplitComm2d:Creating redundant coarse problems of size %d\n",csize);CHKERRQ(ierr);
15047c6ae99SBarry Smith   } else {
15147c6ae99SBarry Smith     *outcomm = comm;
15247c6ae99SBarry Smith   }
15347c6ae99SBarry Smith   PetscFunctionReturn(0);
15447c6ae99SBarry Smith }
15547c6ae99SBarry Smith 
15647c6ae99SBarry Smith #undef __FUNCT__
157*5dff015aSBarry Smith #define __FUNCT__ "DMGetElements_DA_2d_P1"
158*5dff015aSBarry Smith PetscErrorCode DMGetElements_DA_2d_P1(DA da,PetscInt *n,const PetscInt *e[])
15947c6ae99SBarry Smith {
16047c6ae99SBarry Smith   PetscErrorCode ierr;
16147c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
16247c6ae99SBarry Smith   PetscInt       i,j,cnt,xs,xe = dd->xe,ys,ye = dd->ye,Xs = dd->Xs, Xe = dd->Xe, Ys = dd->Ys;
16347c6ae99SBarry Smith 
16447c6ae99SBarry Smith   PetscFunctionBegin;
16547c6ae99SBarry Smith   if (!dd->e) {
16647c6ae99SBarry Smith     if (dd->xs == Xs) xs = dd->xs; else xs = dd->xs - 1;
16747c6ae99SBarry Smith     if (dd->ys == Ys) ys = dd->ys; else ys = dd->ys - 1;
16847c6ae99SBarry Smith     dd->ne = 2*(xe - xs - 1)*(ye - ys - 1);
16947c6ae99SBarry Smith     ierr   = PetscMalloc((1 + 3*dd->ne)*sizeof(PetscInt),&dd->e);CHKERRQ(ierr);
17047c6ae99SBarry Smith     cnt    = 0;
17147c6ae99SBarry Smith     for (j=ys; j<ye-1; j++) {
17247c6ae99SBarry Smith       for (i=xs; i<xe-1; i++) {
17347c6ae99SBarry Smith         dd->e[cnt]   = i - Xs + (j - Ys)*(Xe - Xs);
17447c6ae99SBarry Smith         dd->e[cnt+1] = i - Xs + 1 + (j - Ys)*(Xe - Xs);
17547c6ae99SBarry Smith         dd->e[cnt+2] = i - Xs + (j - Ys + 1)*(Xe - Xs);
17647c6ae99SBarry Smith 
17747c6ae99SBarry Smith         dd->e[cnt+3] = i - Xs + 1 + (j - Ys + 1)*(Xe - Xs);
17847c6ae99SBarry Smith         dd->e[cnt+4] = i - Xs + (j - Ys + 1)*(Xe - Xs);
17947c6ae99SBarry Smith         dd->e[cnt+5] = i - Xs + 1 + (j - Ys)*(Xe - Xs);
18047c6ae99SBarry Smith         cnt += 6;
18147c6ae99SBarry Smith       }
18247c6ae99SBarry Smith     }
18347c6ae99SBarry Smith   }
18447c6ae99SBarry Smith   *n = dd->ne;
18547c6ae99SBarry Smith   *e = dd->e;
18647c6ae99SBarry Smith   PetscFunctionReturn(0);
18747c6ae99SBarry Smith }
18847c6ae99SBarry Smith 
18947c6ae99SBarry Smith #undef __FUNCT__
19047c6ae99SBarry Smith #define __FUNCT__ "DASetLocalFunction"
19147c6ae99SBarry Smith /*@C
19247c6ae99SBarry Smith        DASetLocalFunction - Caches in a DA a local function.
19347c6ae99SBarry Smith 
19447c6ae99SBarry Smith    Logically Collective on DA
19547c6ae99SBarry Smith 
19647c6ae99SBarry Smith    Input Parameter:
19747c6ae99SBarry Smith +  da - initial distributed array
19847c6ae99SBarry Smith -  lf - the local function
19947c6ae99SBarry Smith 
20047c6ae99SBarry Smith    Level: intermediate
20147c6ae99SBarry Smith 
20247c6ae99SBarry Smith    Notes: The routine SNESDAFormFunction() uses this the cached function to evaluate the user provided function.
20347c6ae99SBarry Smith 
20447c6ae99SBarry Smith .keywords:  distributed array, refine
20547c6ae99SBarry Smith 
20647c6ae99SBarry Smith .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DAGetLocalFunction(), DASetLocalFunctioni()
20747c6ae99SBarry Smith @*/
20847c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DASetLocalFunction(DA da,DALocalFunction1 lf)
20947c6ae99SBarry Smith {
21047c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
21147c6ae99SBarry Smith   PetscFunctionBegin;
21247c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
21347c6ae99SBarry Smith   dd->lf    = lf;
21447c6ae99SBarry Smith   PetscFunctionReturn(0);
21547c6ae99SBarry Smith }
21647c6ae99SBarry Smith 
21747c6ae99SBarry Smith #undef __FUNCT__
21847c6ae99SBarry Smith #define __FUNCT__ "DASetLocalFunctioni"
21947c6ae99SBarry Smith /*@C
22047c6ae99SBarry Smith        DASetLocalFunctioni - Caches in a DA a local function that evaluates a single component
22147c6ae99SBarry Smith 
22247c6ae99SBarry Smith    Logically Collective on DA
22347c6ae99SBarry Smith 
22447c6ae99SBarry Smith    Input Parameter:
22547c6ae99SBarry Smith +  da - initial distributed array
22647c6ae99SBarry Smith -  lfi - the local function
22747c6ae99SBarry Smith 
22847c6ae99SBarry Smith    Level: intermediate
22947c6ae99SBarry Smith 
23047c6ae99SBarry Smith .keywords:  distributed array, refine
23147c6ae99SBarry Smith 
23247c6ae99SBarry Smith .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DAGetLocalFunction(), DASetLocalFunction()
23347c6ae99SBarry Smith @*/
23447c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DASetLocalFunctioni(DA da,PetscErrorCode (*lfi)(DALocalInfo*,MatStencil*,void*,PetscScalar*,void*))
23547c6ae99SBarry Smith {
23647c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
23747c6ae99SBarry Smith   PetscFunctionBegin;
23847c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
23947c6ae99SBarry Smith   dd->lfi = lfi;
24047c6ae99SBarry Smith   PetscFunctionReturn(0);
24147c6ae99SBarry Smith }
24247c6ae99SBarry Smith 
24347c6ae99SBarry Smith #undef __FUNCT__
24447c6ae99SBarry Smith #define __FUNCT__ "DASetLocalFunctionib"
24547c6ae99SBarry Smith /*@C
24647c6ae99SBarry Smith        DASetLocalFunctionib - Caches in a DA a block local function that evaluates a single component
24747c6ae99SBarry Smith 
24847c6ae99SBarry Smith    Logically Collective on DA
24947c6ae99SBarry Smith 
25047c6ae99SBarry Smith    Input Parameter:
25147c6ae99SBarry Smith +  da - initial distributed array
25247c6ae99SBarry Smith -  lfi - the local function
25347c6ae99SBarry Smith 
25447c6ae99SBarry Smith    Level: intermediate
25547c6ae99SBarry Smith 
25647c6ae99SBarry Smith .keywords:  distributed array, refine
25747c6ae99SBarry Smith 
25847c6ae99SBarry Smith .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DAGetLocalFunction(), DASetLocalFunction()
25947c6ae99SBarry Smith @*/
26047c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DASetLocalFunctionib(DA da,PetscErrorCode (*lfi)(DALocalInfo*,MatStencil*,void*,PetscScalar*,void*))
26147c6ae99SBarry Smith {
26247c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
26347c6ae99SBarry Smith   PetscFunctionBegin;
26447c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
26547c6ae99SBarry Smith   dd->lfib = lfi;
26647c6ae99SBarry Smith   PetscFunctionReturn(0);
26747c6ae99SBarry Smith }
26847c6ae99SBarry Smith 
26947c6ae99SBarry Smith #undef __FUNCT__
27047c6ae99SBarry Smith #define __FUNCT__ "DASetLocalAdicFunction_Private"
27147c6ae99SBarry Smith PetscErrorCode DASetLocalAdicFunction_Private(DA da,DALocalFunction1 ad_lf)
27247c6ae99SBarry Smith {
27347c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
27447c6ae99SBarry Smith   PetscFunctionBegin;
27547c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
27647c6ae99SBarry Smith   dd->adic_lf = ad_lf;
27747c6ae99SBarry Smith   PetscFunctionReturn(0);
27847c6ae99SBarry Smith }
27947c6ae99SBarry Smith 
28047c6ae99SBarry Smith /*MC
28147c6ae99SBarry Smith        DASetLocalAdicFunctioni - Caches in a DA a local functioni computed by ADIC/ADIFOR
28247c6ae99SBarry Smith 
28347c6ae99SBarry Smith    Synopsis:
28447c6ae99SBarry Smith    PetscErrorCode DASetLocalAdicFunctioni(DA da,PetscInt (ad_lf*)(DALocalInfo*,MatStencil*,void*,void*,void*)
28547c6ae99SBarry Smith 
28647c6ae99SBarry Smith    Logically Collective on DA
28747c6ae99SBarry Smith 
28847c6ae99SBarry Smith    Input Parameter:
28947c6ae99SBarry Smith +  da - initial distributed array
29047c6ae99SBarry Smith -  ad_lfi - the local function as computed by ADIC/ADIFOR
29147c6ae99SBarry Smith 
29247c6ae99SBarry Smith    Level: intermediate
29347c6ae99SBarry Smith 
29447c6ae99SBarry Smith .keywords:  distributed array, refine
29547c6ae99SBarry Smith 
29647c6ae99SBarry Smith .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DAGetLocalFunction(), DASetLocalFunction(),
29747c6ae99SBarry Smith           DASetLocalJacobian(), DASetLocalFunctioni()
29847c6ae99SBarry Smith M*/
29947c6ae99SBarry Smith 
30047c6ae99SBarry Smith #undef __FUNCT__
30147c6ae99SBarry Smith #define __FUNCT__ "DASetLocalAdicFunctioni_Private"
30247c6ae99SBarry Smith PetscErrorCode DASetLocalAdicFunctioni_Private(DA da,PetscErrorCode (*ad_lfi)(DALocalInfo*,MatStencil*,void*,void*,void*))
30347c6ae99SBarry Smith {
30447c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
30547c6ae99SBarry Smith   PetscFunctionBegin;
30647c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
30747c6ae99SBarry Smith   dd->adic_lfi = ad_lfi;
30847c6ae99SBarry Smith   PetscFunctionReturn(0);
30947c6ae99SBarry Smith }
31047c6ae99SBarry Smith 
31147c6ae99SBarry Smith /*MC
31247c6ae99SBarry Smith        DASetLocalAdicMFFunctioni - Caches in a DA a local functioni computed by ADIC/ADIFOR
31347c6ae99SBarry Smith 
31447c6ae99SBarry Smith    Synopsis:
31547c6ae99SBarry Smith    PetscErrorCode  DASetLocalAdicFunctioni(DA da,int (ad_lf*)(DALocalInfo*,MatStencil*,void*,void*,void*)
31647c6ae99SBarry Smith 
31747c6ae99SBarry Smith    Logically Collective on DA
31847c6ae99SBarry Smith 
31947c6ae99SBarry Smith    Input Parameter:
32047c6ae99SBarry Smith +  da - initial distributed array
32147c6ae99SBarry Smith -  admf_lfi - the local matrix-free function as computed by ADIC/ADIFOR
32247c6ae99SBarry Smith 
32347c6ae99SBarry Smith    Level: intermediate
32447c6ae99SBarry Smith 
32547c6ae99SBarry Smith .keywords:  distributed array, refine
32647c6ae99SBarry Smith 
32747c6ae99SBarry Smith .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DAGetLocalFunction(), DASetLocalFunction(),
32847c6ae99SBarry Smith           DASetLocalJacobian(), DASetLocalFunctioni()
32947c6ae99SBarry Smith M*/
33047c6ae99SBarry Smith 
33147c6ae99SBarry Smith #undef __FUNCT__
33247c6ae99SBarry Smith #define __FUNCT__ "DASetLocalAdicMFFunctioni_Private"
33347c6ae99SBarry Smith PetscErrorCode DASetLocalAdicMFFunctioni_Private(DA da,PetscErrorCode (*admf_lfi)(DALocalInfo*,MatStencil*,void*,void*,void*))
33447c6ae99SBarry Smith {
33547c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
33647c6ae99SBarry Smith   PetscFunctionBegin;
33747c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
33847c6ae99SBarry Smith   dd->adicmf_lfi = admf_lfi;
33947c6ae99SBarry Smith   PetscFunctionReturn(0);
34047c6ae99SBarry Smith }
34147c6ae99SBarry Smith 
34247c6ae99SBarry Smith /*MC
34347c6ae99SBarry Smith        DASetLocalAdicFunctionib - Caches in a DA a block local functioni computed by ADIC/ADIFOR
34447c6ae99SBarry Smith 
34547c6ae99SBarry Smith    Synopsis:
34647c6ae99SBarry Smith    PetscErrorCode DASetLocalAdicFunctionib(DA da,PetscInt (ad_lf*)(DALocalInfo*,MatStencil*,void*,void*,void*)
34747c6ae99SBarry Smith 
34847c6ae99SBarry Smith    Logically Collective on DA
34947c6ae99SBarry Smith 
35047c6ae99SBarry Smith    Input Parameter:
35147c6ae99SBarry Smith +  da - initial distributed array
35247c6ae99SBarry Smith -  ad_lfi - the local function as computed by ADIC/ADIFOR
35347c6ae99SBarry Smith 
35447c6ae99SBarry Smith    Level: intermediate
35547c6ae99SBarry Smith 
35647c6ae99SBarry Smith .keywords:  distributed array, refine
35747c6ae99SBarry Smith 
35847c6ae99SBarry Smith .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DAGetLocalFunction(), DASetLocalFunction(),
35947c6ae99SBarry Smith           DASetLocalJacobian(), DASetLocalFunctionib()
36047c6ae99SBarry Smith M*/
36147c6ae99SBarry Smith 
36247c6ae99SBarry Smith #undef __FUNCT__
36347c6ae99SBarry Smith #define __FUNCT__ "DASetLocalAdicFunctionib_Private"
36447c6ae99SBarry Smith PetscErrorCode DASetLocalAdicFunctionib_Private(DA da,PetscErrorCode (*ad_lfi)(DALocalInfo*,MatStencil*,void*,void*,void*))
36547c6ae99SBarry Smith {
36647c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
36747c6ae99SBarry Smith   PetscFunctionBegin;
36847c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
36947c6ae99SBarry Smith   dd->adic_lfib = ad_lfi;
37047c6ae99SBarry Smith   PetscFunctionReturn(0);
37147c6ae99SBarry Smith }
37247c6ae99SBarry Smith 
37347c6ae99SBarry Smith /*MC
37447c6ae99SBarry Smith        DASetLocalAdicMFFunctionib - Caches in a DA a block local functioni computed by ADIC/ADIFOR
37547c6ae99SBarry Smith 
37647c6ae99SBarry Smith    Synopsis:
37747c6ae99SBarry Smith    PetscErrorCode  DASetLocalAdicFunctionib(DA da,int (ad_lf*)(DALocalInfo*,MatStencil*,void*,void*,void*)
37847c6ae99SBarry Smith 
37947c6ae99SBarry Smith    Logically Collective on DA
38047c6ae99SBarry Smith 
38147c6ae99SBarry Smith    Input Parameter:
38247c6ae99SBarry Smith +  da - initial distributed array
38347c6ae99SBarry Smith -  admf_lfi - the local matrix-free function as computed by ADIC/ADIFOR
38447c6ae99SBarry Smith 
38547c6ae99SBarry Smith    Level: intermediate
38647c6ae99SBarry Smith 
38747c6ae99SBarry Smith .keywords:  distributed array, refine
38847c6ae99SBarry Smith 
38947c6ae99SBarry Smith .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DAGetLocalFunction(), DASetLocalFunction(),
39047c6ae99SBarry Smith           DASetLocalJacobian(), DASetLocalFunctionib()
39147c6ae99SBarry Smith M*/
39247c6ae99SBarry Smith 
39347c6ae99SBarry Smith #undef __FUNCT__
39447c6ae99SBarry Smith #define __FUNCT__ "DASetLocalAdicMFFunctionib_Private"
39547c6ae99SBarry Smith PetscErrorCode DASetLocalAdicMFFunctionib_Private(DA da,PetscErrorCode (*admf_lfi)(DALocalInfo*,MatStencil*,void*,void*,void*))
39647c6ae99SBarry Smith {
39747c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
39847c6ae99SBarry Smith   PetscFunctionBegin;
39947c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
40047c6ae99SBarry Smith   dd->adicmf_lfib = admf_lfi;
40147c6ae99SBarry Smith   PetscFunctionReturn(0);
40247c6ae99SBarry Smith }
40347c6ae99SBarry Smith 
40447c6ae99SBarry Smith /*MC
40547c6ae99SBarry Smith        DASetLocalAdicMFFunction - Caches in a DA a local function computed by ADIC/ADIFOR
40647c6ae99SBarry Smith 
40747c6ae99SBarry Smith    Synopsis:
40847c6ae99SBarry Smith    PetscErrorCode DASetLocalAdicMFFunction(DA da,DALocalFunction1 ad_lf)
40947c6ae99SBarry Smith 
41047c6ae99SBarry Smith    Logically Collective on DA
41147c6ae99SBarry Smith 
41247c6ae99SBarry Smith    Input Parameter:
41347c6ae99SBarry Smith +  da - initial distributed array
41447c6ae99SBarry Smith -  ad_lf - the local function as computed by ADIC/ADIFOR
41547c6ae99SBarry Smith 
41647c6ae99SBarry Smith    Level: intermediate
41747c6ae99SBarry Smith 
41847c6ae99SBarry Smith .keywords:  distributed array, refine
41947c6ae99SBarry Smith 
42047c6ae99SBarry Smith .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DAGetLocalFunction(), DASetLocalFunction(),
42147c6ae99SBarry Smith           DASetLocalJacobian()
42247c6ae99SBarry Smith M*/
42347c6ae99SBarry Smith 
42447c6ae99SBarry Smith #undef __FUNCT__
42547c6ae99SBarry Smith #define __FUNCT__ "DASetLocalAdicMFFunction_Private"
42647c6ae99SBarry Smith PetscErrorCode DASetLocalAdicMFFunction_Private(DA da,DALocalFunction1 ad_lf)
42747c6ae99SBarry Smith {
42847c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
42947c6ae99SBarry Smith   PetscFunctionBegin;
43047c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
43147c6ae99SBarry Smith   dd->adicmf_lf = ad_lf;
43247c6ae99SBarry Smith   PetscFunctionReturn(0);
43347c6ae99SBarry Smith }
43447c6ae99SBarry Smith 
43547c6ae99SBarry Smith /*@C
43647c6ae99SBarry Smith        DASetLocalJacobian - Caches in a DA a local Jacobian computation function
43747c6ae99SBarry Smith 
43847c6ae99SBarry Smith    Logically Collective on DA
43947c6ae99SBarry Smith 
44047c6ae99SBarry Smith 
44147c6ae99SBarry Smith    Input Parameter:
44247c6ae99SBarry Smith +  da - initial distributed array
44347c6ae99SBarry Smith -  lj - the local Jacobian
44447c6ae99SBarry Smith 
44547c6ae99SBarry Smith    Level: intermediate
44647c6ae99SBarry Smith 
44747c6ae99SBarry Smith    Notes: The routine SNESDAFormFunction() uses this the cached function to evaluate the user provided function.
44847c6ae99SBarry Smith 
44947c6ae99SBarry Smith .keywords:  distributed array, refine
45047c6ae99SBarry Smith 
45147c6ae99SBarry Smith .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DAGetLocalFunction(), DASetLocalFunction()
45247c6ae99SBarry Smith @*/
45347c6ae99SBarry Smith #undef __FUNCT__
45447c6ae99SBarry Smith #define __FUNCT__ "DASetLocalJacobian"
45547c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DASetLocalJacobian(DA da,DALocalFunction1 lj)
45647c6ae99SBarry Smith {
45747c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
45847c6ae99SBarry Smith   PetscFunctionBegin;
45947c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
46047c6ae99SBarry Smith   dd->lj    = lj;
46147c6ae99SBarry Smith   PetscFunctionReturn(0);
46247c6ae99SBarry Smith }
46347c6ae99SBarry Smith 
46447c6ae99SBarry Smith #undef __FUNCT__
46547c6ae99SBarry Smith #define __FUNCT__ "DAGetLocalFunction"
46647c6ae99SBarry Smith /*@C
46747c6ae99SBarry Smith        DAGetLocalFunction - Gets from a DA a local function and its ADIC/ADIFOR Jacobian
46847c6ae99SBarry Smith 
46947c6ae99SBarry Smith    Note Collective
47047c6ae99SBarry Smith 
47147c6ae99SBarry Smith    Input Parameter:
47247c6ae99SBarry Smith .  da - initial distributed array
47347c6ae99SBarry Smith 
47447c6ae99SBarry Smith    Output Parameter:
47547c6ae99SBarry Smith .  lf - the local function
47647c6ae99SBarry Smith 
47747c6ae99SBarry Smith    Level: intermediate
47847c6ae99SBarry Smith 
47947c6ae99SBarry Smith .keywords:  distributed array, refine
48047c6ae99SBarry Smith 
48147c6ae99SBarry Smith .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DAGetLocalJacobian(), DASetLocalFunction()
48247c6ae99SBarry Smith @*/
48347c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DAGetLocalFunction(DA da,DALocalFunction1 *lf)
48447c6ae99SBarry Smith {
48547c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
48647c6ae99SBarry Smith   PetscFunctionBegin;
48747c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
48847c6ae99SBarry Smith   if (lf)       *lf = dd->lf;
48947c6ae99SBarry Smith   PetscFunctionReturn(0);
49047c6ae99SBarry Smith }
49147c6ae99SBarry Smith 
49247c6ae99SBarry Smith #undef __FUNCT__
49347c6ae99SBarry Smith #define __FUNCT__ "DAGetLocalJacobian"
49447c6ae99SBarry Smith /*@C
49547c6ae99SBarry Smith        DAGetLocalJacobian - Gets from a DA a local jacobian
49647c6ae99SBarry Smith 
49747c6ae99SBarry Smith    Not Collective
49847c6ae99SBarry Smith 
49947c6ae99SBarry Smith    Input Parameter:
50047c6ae99SBarry Smith .  da - initial distributed array
50147c6ae99SBarry Smith 
50247c6ae99SBarry Smith    Output Parameter:
50347c6ae99SBarry Smith .  lj - the local jacobian
50447c6ae99SBarry Smith 
50547c6ae99SBarry Smith    Level: intermediate
50647c6ae99SBarry Smith 
50747c6ae99SBarry Smith .keywords:  distributed array, refine
50847c6ae99SBarry Smith 
50947c6ae99SBarry Smith .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DAGetLocalFunction(), DASetLocalJacobian()
51047c6ae99SBarry Smith @*/
51147c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DAGetLocalJacobian(DA da,DALocalFunction1 *lj)
51247c6ae99SBarry Smith {
51347c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
51447c6ae99SBarry Smith   PetscFunctionBegin;
51547c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
51647c6ae99SBarry Smith   if (lj) *lj = dd->lj;
51747c6ae99SBarry Smith   PetscFunctionReturn(0);
51847c6ae99SBarry Smith }
51947c6ae99SBarry Smith 
52047c6ae99SBarry Smith #undef __FUNCT__
52147c6ae99SBarry Smith #define __FUNCT__ "DAFormFunction"
52247c6ae99SBarry Smith /*@
52347c6ae99SBarry Smith     DAFormFunction - Evaluates a user provided function on each processor that
52447c6ae99SBarry Smith         share a DA
52547c6ae99SBarry Smith 
52647c6ae99SBarry Smith    Input Parameters:
52747c6ae99SBarry Smith +    da - the DA that defines the grid
52847c6ae99SBarry Smith .    vu - input vector
52947c6ae99SBarry Smith .    vfu - output vector
53047c6ae99SBarry Smith -    w - any user data
53147c6ae99SBarry Smith 
53247c6ae99SBarry Smith     Notes: Does NOT do ghost updates on vu upon entry
53347c6ae99SBarry Smith 
53447c6ae99SBarry Smith            This should eventually replace DAFormFunction1
53547c6ae99SBarry Smith 
53647c6ae99SBarry Smith     Level: advanced
53747c6ae99SBarry Smith 
53847c6ae99SBarry Smith .seealso: DAComputeJacobian1WithAdic()
53947c6ae99SBarry Smith 
54047c6ae99SBarry Smith @*/
54147c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DAFormFunction(DA da,PetscErrorCode (*lf)(void),Vec vu,Vec vfu,void *w)
54247c6ae99SBarry Smith {
54347c6ae99SBarry Smith   PetscErrorCode ierr;
54447c6ae99SBarry Smith   void           *u,*fu;
54547c6ae99SBarry Smith   DALocalInfo    info;
54647c6ae99SBarry Smith   PetscErrorCode (*f)(DALocalInfo*,void*,void*,void*) = (PetscErrorCode (*)(DALocalInfo*,void*,void*,void*))lf;
54747c6ae99SBarry Smith 
54847c6ae99SBarry Smith   PetscFunctionBegin;
54947c6ae99SBarry Smith   ierr = DAGetLocalInfo(da,&info);CHKERRQ(ierr);
55047c6ae99SBarry Smith   ierr = DAVecGetArray(da,vu,&u);CHKERRQ(ierr);
55147c6ae99SBarry Smith   ierr = DAVecGetArray(da,vfu,&fu);CHKERRQ(ierr);
55247c6ae99SBarry Smith 
55347c6ae99SBarry Smith   ierr = (*f)(&info,u,fu,w);
55447c6ae99SBarry Smith   if (PetscExceptionValue(ierr)) {
55547c6ae99SBarry Smith     PetscErrorCode pierr = DAVecRestoreArray(da,vu,&u);CHKERRQ(pierr);
55647c6ae99SBarry Smith     pierr = DAVecRestoreArray(da,vfu,&fu);CHKERRQ(pierr);
55747c6ae99SBarry Smith   }
55847c6ae99SBarry Smith  CHKERRQ(ierr);
55947c6ae99SBarry Smith 
56047c6ae99SBarry Smith   ierr = DAVecRestoreArray(da,vu,&u);CHKERRQ(ierr);
56147c6ae99SBarry Smith   ierr = DAVecRestoreArray(da,vfu,&fu);CHKERRQ(ierr);
56247c6ae99SBarry Smith   PetscFunctionReturn(0);
56347c6ae99SBarry Smith }
56447c6ae99SBarry Smith 
56547c6ae99SBarry Smith #undef __FUNCT__
56647c6ae99SBarry Smith #define __FUNCT__ "DAFormFunctionLocal"
56747c6ae99SBarry Smith /*@C
56847c6ae99SBarry Smith    DAFormFunctionLocal - This is a universal function evaluation routine for
56947c6ae99SBarry Smith    a local DA function.
57047c6ae99SBarry Smith 
57147c6ae99SBarry Smith    Collective on DA
57247c6ae99SBarry Smith 
57347c6ae99SBarry Smith    Input Parameters:
57447c6ae99SBarry Smith +  da - the DA context
57547c6ae99SBarry Smith .  func - The local function
57647c6ae99SBarry Smith .  X - input vector
57747c6ae99SBarry Smith .  F - function vector
57847c6ae99SBarry Smith -  ctx - A user context
57947c6ae99SBarry Smith 
58047c6ae99SBarry Smith    Level: intermediate
58147c6ae99SBarry Smith 
58247c6ae99SBarry Smith .seealso: DASetLocalFunction(), DASetLocalJacobian(), DASetLocalAdicFunction(), DASetLocalAdicMFFunction(),
58347c6ae99SBarry Smith           SNESSetFunction(), SNESSetJacobian()
58447c6ae99SBarry Smith 
58547c6ae99SBarry Smith @*/
58647c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DAFormFunctionLocal(DA da, DALocalFunction1 func, Vec X, Vec F, void *ctx)
58747c6ae99SBarry Smith {
58847c6ae99SBarry Smith   Vec            localX;
58947c6ae99SBarry Smith   DALocalInfo    info;
59047c6ae99SBarry Smith   void          *u;
59147c6ae99SBarry Smith   void          *fu;
59247c6ae99SBarry Smith   PetscErrorCode ierr;
59347c6ae99SBarry Smith 
59447c6ae99SBarry Smith   PetscFunctionBegin;
59547c6ae99SBarry Smith   ierr = DAGetLocalVector(da,&localX);CHKERRQ(ierr);
59647c6ae99SBarry Smith   /*
59747c6ae99SBarry Smith      Scatter ghost points to local vector, using the 2-step process
59847c6ae99SBarry Smith         DAGlobalToLocalBegin(), DAGlobalToLocalEnd().
59947c6ae99SBarry Smith   */
60047c6ae99SBarry Smith   ierr = DAGlobalToLocalBegin(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
60147c6ae99SBarry Smith   ierr = DAGlobalToLocalEnd(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
60247c6ae99SBarry Smith   ierr = DAGetLocalInfo(da,&info);CHKERRQ(ierr);
60347c6ae99SBarry Smith   ierr = DAVecGetArray(da,localX,&u);CHKERRQ(ierr);
60447c6ae99SBarry Smith   ierr = DAVecGetArray(da,F,&fu);CHKERRQ(ierr);
60547c6ae99SBarry Smith   ierr = (*func)(&info,u,fu,ctx);
60647c6ae99SBarry Smith   if (PetscExceptionValue(ierr)) {
60747c6ae99SBarry Smith     PetscErrorCode pierr = DAVecRestoreArray(da,localX,&u);CHKERRQ(pierr);
60847c6ae99SBarry Smith     pierr = DAVecRestoreArray(da,F,&fu);CHKERRQ(pierr);
60947c6ae99SBarry Smith   }
61047c6ae99SBarry Smith  CHKERRQ(ierr);
61147c6ae99SBarry Smith   ierr = DAVecRestoreArray(da,localX,&u);CHKERRQ(ierr);
61247c6ae99SBarry Smith   ierr = DAVecRestoreArray(da,F,&fu);CHKERRQ(ierr);
61347c6ae99SBarry Smith   if (PetscExceptionValue(ierr)) {
61447c6ae99SBarry Smith     PetscErrorCode pierr = DARestoreLocalVector(da,&localX);CHKERRQ(pierr);
61547c6ae99SBarry Smith   }
61647c6ae99SBarry Smith  CHKERRQ(ierr);
61747c6ae99SBarry Smith   ierr = DARestoreLocalVector(da,&localX);CHKERRQ(ierr);
61847c6ae99SBarry Smith   PetscFunctionReturn(0);
61947c6ae99SBarry Smith }
62047c6ae99SBarry Smith 
62147c6ae99SBarry Smith #undef __FUNCT__
62247c6ae99SBarry Smith #define __FUNCT__ "DAFormFunctionLocalGhost"
62347c6ae99SBarry Smith /*@C
62447c6ae99SBarry Smith    DAFormFunctionLocalGhost - This is a universal function evaluation routine for
62547c6ae99SBarry Smith    a local DA function, but the ghost values of the output are communicated and added.
62647c6ae99SBarry Smith 
62747c6ae99SBarry Smith    Collective on DA
62847c6ae99SBarry Smith 
62947c6ae99SBarry Smith    Input Parameters:
63047c6ae99SBarry Smith +  da - the DA context
63147c6ae99SBarry Smith .  func - The local function
63247c6ae99SBarry Smith .  X - input vector
63347c6ae99SBarry Smith .  F - function vector
63447c6ae99SBarry Smith -  ctx - A user context
63547c6ae99SBarry Smith 
63647c6ae99SBarry Smith    Level: intermediate
63747c6ae99SBarry Smith 
63847c6ae99SBarry Smith .seealso: DASetLocalFunction(), DASetLocalJacobian(), DASetLocalAdicFunction(), DASetLocalAdicMFFunction(),
63947c6ae99SBarry Smith           SNESSetFunction(), SNESSetJacobian()
64047c6ae99SBarry Smith 
64147c6ae99SBarry Smith @*/
64247c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DAFormFunctionLocalGhost(DA da, DALocalFunction1 func, Vec X, Vec F, void *ctx)
64347c6ae99SBarry Smith {
64447c6ae99SBarry Smith   Vec            localX, localF;
64547c6ae99SBarry Smith   DALocalInfo    info;
64647c6ae99SBarry Smith   void          *u;
64747c6ae99SBarry Smith   void          *fu;
64847c6ae99SBarry Smith   PetscErrorCode ierr;
64947c6ae99SBarry Smith 
65047c6ae99SBarry Smith   PetscFunctionBegin;
65147c6ae99SBarry Smith   ierr = DAGetLocalVector(da,&localX);CHKERRQ(ierr);
65247c6ae99SBarry Smith   ierr = DAGetLocalVector(da,&localF);CHKERRQ(ierr);
65347c6ae99SBarry Smith   /*
65447c6ae99SBarry Smith      Scatter ghost points to local vector, using the 2-step process
65547c6ae99SBarry Smith         DAGlobalToLocalBegin(), DAGlobalToLocalEnd().
65647c6ae99SBarry Smith   */
65747c6ae99SBarry Smith   ierr = DAGlobalToLocalBegin(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
65847c6ae99SBarry Smith   ierr = DAGlobalToLocalEnd(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
65947c6ae99SBarry Smith   ierr = VecSet(F, 0.0);CHKERRQ(ierr);
66047c6ae99SBarry Smith   ierr = VecSet(localF, 0.0);CHKERRQ(ierr);
66147c6ae99SBarry Smith   ierr = DAGetLocalInfo(da,&info);CHKERRQ(ierr);
66247c6ae99SBarry Smith   ierr = DAVecGetArray(da,localX,&u);CHKERRQ(ierr);
66347c6ae99SBarry Smith   ierr = DAVecGetArray(da,localF,&fu);CHKERRQ(ierr);
66447c6ae99SBarry Smith   ierr = (*func)(&info,u,fu,ctx);
66547c6ae99SBarry Smith   if (PetscExceptionValue(ierr)) {
66647c6ae99SBarry Smith     PetscErrorCode pierr = DAVecRestoreArray(da,localX,&u);CHKERRQ(pierr);
66747c6ae99SBarry Smith     pierr = DAVecRestoreArray(da,localF,&fu);CHKERRQ(pierr);
66847c6ae99SBarry Smith   }
66947c6ae99SBarry Smith   CHKERRQ(ierr);
67047c6ae99SBarry Smith   ierr = DALocalToGlobalBegin(da,localF,F);CHKERRQ(ierr);
67147c6ae99SBarry Smith   ierr = DALocalToGlobalEnd(da,localF,F);CHKERRQ(ierr);
67247c6ae99SBarry Smith   ierr = DAVecRestoreArray(da,localX,&u);CHKERRQ(ierr);
67347c6ae99SBarry Smith   ierr = DAVecRestoreArray(da,localF,&fu);CHKERRQ(ierr);
67447c6ae99SBarry Smith   if (PetscExceptionValue(ierr)) {
67547c6ae99SBarry Smith     PetscErrorCode pierr = DARestoreLocalVector(da,&localX);CHKERRQ(pierr);
67647c6ae99SBarry Smith   ierr = DARestoreLocalVector(da,&localF);CHKERRQ(ierr);
67747c6ae99SBarry Smith   }
67847c6ae99SBarry Smith   CHKERRQ(ierr);
67947c6ae99SBarry Smith   ierr = DARestoreLocalVector(da,&localX);CHKERRQ(ierr);
68047c6ae99SBarry Smith   ierr = DARestoreLocalVector(da,&localF);CHKERRQ(ierr);
68147c6ae99SBarry Smith   PetscFunctionReturn(0);
68247c6ae99SBarry Smith }
68347c6ae99SBarry Smith 
68447c6ae99SBarry Smith #undef __FUNCT__
68547c6ae99SBarry Smith #define __FUNCT__ "DAFormFunction1"
68647c6ae99SBarry Smith /*@
68747c6ae99SBarry Smith     DAFormFunction1 - Evaluates a user provided function on each processor that
68847c6ae99SBarry Smith         share a DA
68947c6ae99SBarry Smith 
69047c6ae99SBarry Smith    Input Parameters:
69147c6ae99SBarry Smith +    da - the DA that defines the grid
69247c6ae99SBarry Smith .    vu - input vector
69347c6ae99SBarry Smith .    vfu - output vector
69447c6ae99SBarry Smith -    w - any user data
69547c6ae99SBarry Smith 
69647c6ae99SBarry Smith     Notes: Does NOT do ghost updates on vu upon entry
69747c6ae99SBarry Smith 
69847c6ae99SBarry Smith     Level: advanced
69947c6ae99SBarry Smith 
70047c6ae99SBarry Smith .seealso: DAComputeJacobian1WithAdic()
70147c6ae99SBarry Smith 
70247c6ae99SBarry Smith @*/
70347c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DAFormFunction1(DA da,Vec vu,Vec vfu,void *w)
70447c6ae99SBarry Smith {
70547c6ae99SBarry Smith   PetscErrorCode ierr;
70647c6ae99SBarry Smith   void           *u,*fu;
70747c6ae99SBarry Smith   DALocalInfo    info;
70847c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
70947c6ae99SBarry Smith 
71047c6ae99SBarry Smith   PetscFunctionBegin;
71147c6ae99SBarry Smith 
71247c6ae99SBarry Smith   ierr = DAGetLocalInfo(da,&info);CHKERRQ(ierr);
71347c6ae99SBarry Smith   ierr = DAVecGetArray(da,vu,&u);CHKERRQ(ierr);
71447c6ae99SBarry Smith   ierr = DAVecGetArray(da,vfu,&fu);CHKERRQ(ierr);
71547c6ae99SBarry Smith 
71647c6ae99SBarry Smith   CHKMEMQ;
71747c6ae99SBarry Smith   ierr = (*dd->lf)(&info,u,fu,w);
71847c6ae99SBarry Smith   if (PetscExceptionValue(ierr)) {
71947c6ae99SBarry Smith     PetscErrorCode pierr = DAVecRestoreArray(da,vu,&u);CHKERRQ(pierr);
72047c6ae99SBarry Smith     pierr = DAVecRestoreArray(da,vfu,&fu);CHKERRQ(pierr);
72147c6ae99SBarry Smith   }
72247c6ae99SBarry Smith   CHKERRQ(ierr);
72347c6ae99SBarry Smith   CHKMEMQ;
72447c6ae99SBarry Smith 
72547c6ae99SBarry Smith   ierr = DAVecRestoreArray(da,vu,&u);CHKERRQ(ierr);
72647c6ae99SBarry Smith   ierr = DAVecRestoreArray(da,vfu,&fu);CHKERRQ(ierr);
72747c6ae99SBarry Smith   PetscFunctionReturn(0);
72847c6ae99SBarry Smith }
72947c6ae99SBarry Smith 
73047c6ae99SBarry Smith #undef __FUNCT__
73147c6ae99SBarry Smith #define __FUNCT__ "DAFormFunctioniTest1"
73247c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DAFormFunctioniTest1(DA da,void *w)
73347c6ae99SBarry Smith {
73447c6ae99SBarry Smith   Vec            vu,fu,fui;
73547c6ae99SBarry Smith   PetscErrorCode ierr;
73647c6ae99SBarry Smith   PetscInt       i,n;
73747c6ae99SBarry Smith   PetscScalar    *ui;
73847c6ae99SBarry Smith   PetscRandom    rnd;
73947c6ae99SBarry Smith   PetscReal      norm;
74047c6ae99SBarry Smith 
74147c6ae99SBarry Smith   PetscFunctionBegin;
74247c6ae99SBarry Smith   ierr = DAGetLocalVector(da,&vu);CHKERRQ(ierr);
74347c6ae99SBarry Smith   ierr = PetscRandomCreate(PETSC_COMM_SELF,&rnd);CHKERRQ(ierr);
74447c6ae99SBarry Smith   ierr = PetscRandomSetFromOptions(rnd);CHKERRQ(ierr);
74547c6ae99SBarry Smith   ierr = VecSetRandom(vu,rnd);CHKERRQ(ierr);
74647c6ae99SBarry Smith   ierr = PetscRandomDestroy(rnd);CHKERRQ(ierr);
74747c6ae99SBarry Smith 
74847c6ae99SBarry Smith   ierr = DAGetGlobalVector(da,&fu);CHKERRQ(ierr);
74947c6ae99SBarry Smith   ierr = DAGetGlobalVector(da,&fui);CHKERRQ(ierr);
75047c6ae99SBarry Smith 
75147c6ae99SBarry Smith   ierr = DAFormFunction1(da,vu,fu,w);CHKERRQ(ierr);
75247c6ae99SBarry Smith 
75347c6ae99SBarry Smith   ierr = VecGetArray(fui,&ui);CHKERRQ(ierr);
75447c6ae99SBarry Smith   ierr = VecGetLocalSize(fui,&n);CHKERRQ(ierr);
75547c6ae99SBarry Smith   for (i=0; i<n; i++) {
75647c6ae99SBarry Smith     ierr = DAFormFunctioni1(da,i,vu,ui+i,w);CHKERRQ(ierr);
75747c6ae99SBarry Smith   }
75847c6ae99SBarry Smith   ierr = VecRestoreArray(fui,&ui);CHKERRQ(ierr);
75947c6ae99SBarry Smith 
76047c6ae99SBarry Smith   ierr = VecAXPY(fui,-1.0,fu);CHKERRQ(ierr);
76147c6ae99SBarry Smith   ierr = VecNorm(fui,NORM_2,&norm);CHKERRQ(ierr);
76247c6ae99SBarry Smith   ierr = PetscPrintf(((PetscObject)da)->comm,"Norm of difference in vectors %G\n",norm);CHKERRQ(ierr);
76347c6ae99SBarry Smith   ierr = VecView(fu,0);CHKERRQ(ierr);
76447c6ae99SBarry Smith   ierr = VecView(fui,0);CHKERRQ(ierr);
76547c6ae99SBarry Smith 
76647c6ae99SBarry Smith   ierr = DARestoreLocalVector(da,&vu);CHKERRQ(ierr);
76747c6ae99SBarry Smith   ierr = DARestoreGlobalVector(da,&fu);CHKERRQ(ierr);
76847c6ae99SBarry Smith   ierr = DARestoreGlobalVector(da,&fui);CHKERRQ(ierr);
76947c6ae99SBarry Smith   PetscFunctionReturn(0);
77047c6ae99SBarry Smith }
77147c6ae99SBarry Smith 
77247c6ae99SBarry Smith #undef __FUNCT__
77347c6ae99SBarry Smith #define __FUNCT__ "DAFormFunctioni1"
77447c6ae99SBarry Smith /*@
77547c6ae99SBarry Smith     DAFormFunctioni1 - Evaluates a user provided point-wise function
77647c6ae99SBarry Smith 
77747c6ae99SBarry Smith    Input Parameters:
77847c6ae99SBarry Smith +    da - the DA that defines the grid
77947c6ae99SBarry Smith .    i - the component of the function we wish to compute (must be local)
78047c6ae99SBarry Smith .    vu - input vector
78147c6ae99SBarry Smith .    vfu - output value
78247c6ae99SBarry Smith -    w - any user data
78347c6ae99SBarry Smith 
78447c6ae99SBarry Smith     Notes: Does NOT do ghost updates on vu upon entry
78547c6ae99SBarry Smith 
78647c6ae99SBarry Smith     Level: advanced
78747c6ae99SBarry Smith 
78847c6ae99SBarry Smith .seealso: DAComputeJacobian1WithAdic()
78947c6ae99SBarry Smith 
79047c6ae99SBarry Smith @*/
79147c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DAFormFunctioni1(DA da,PetscInt i,Vec vu,PetscScalar *vfu,void *w)
79247c6ae99SBarry Smith {
79347c6ae99SBarry Smith   PetscErrorCode ierr;
79447c6ae99SBarry Smith   void           *u;
79547c6ae99SBarry Smith   DALocalInfo    info;
79647c6ae99SBarry Smith   MatStencil     stencil;
79747c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
79847c6ae99SBarry Smith 
79947c6ae99SBarry Smith   PetscFunctionBegin;
80047c6ae99SBarry Smith 
80147c6ae99SBarry Smith   ierr = DAGetLocalInfo(da,&info);CHKERRQ(ierr);
80247c6ae99SBarry Smith   ierr = DAVecGetArray(da,vu,&u);CHKERRQ(ierr);
80347c6ae99SBarry Smith 
80447c6ae99SBarry Smith   /* figure out stencil value from i */
80547c6ae99SBarry Smith   stencil.c = i % info.dof;
80647c6ae99SBarry Smith   stencil.i = (i % (info.xm*info.dof))/info.dof;
80747c6ae99SBarry Smith   stencil.j = (i % (info.xm*info.ym*info.dof))/(info.xm*info.dof);
80847c6ae99SBarry Smith   stencil.k = i/(info.xm*info.ym*info.dof);
80947c6ae99SBarry Smith 
81047c6ae99SBarry Smith   ierr = (*dd->lfi)(&info,&stencil,u,vfu,w);CHKERRQ(ierr);
81147c6ae99SBarry Smith 
81247c6ae99SBarry Smith   ierr = DAVecRestoreArray(da,vu,&u);CHKERRQ(ierr);
81347c6ae99SBarry Smith   PetscFunctionReturn(0);
81447c6ae99SBarry Smith }
81547c6ae99SBarry Smith 
81647c6ae99SBarry Smith #undef __FUNCT__
81747c6ae99SBarry Smith #define __FUNCT__ "DAFormFunctionib1"
81847c6ae99SBarry Smith /*@
81947c6ae99SBarry Smith     DAFormFunctionib1 - Evaluates a user provided point-block function
82047c6ae99SBarry Smith 
82147c6ae99SBarry Smith    Input Parameters:
82247c6ae99SBarry Smith +    da - the DA that defines the grid
82347c6ae99SBarry Smith .    i - the component of the function we wish to compute (must be local)
82447c6ae99SBarry Smith .    vu - input vector
82547c6ae99SBarry Smith .    vfu - output value
82647c6ae99SBarry Smith -    w - any user data
82747c6ae99SBarry Smith 
82847c6ae99SBarry Smith     Notes: Does NOT do ghost updates on vu upon entry
82947c6ae99SBarry Smith 
83047c6ae99SBarry Smith     Level: advanced
83147c6ae99SBarry Smith 
83247c6ae99SBarry Smith .seealso: DAComputeJacobian1WithAdic()
83347c6ae99SBarry Smith 
83447c6ae99SBarry Smith @*/
83547c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DAFormFunctionib1(DA da,PetscInt i,Vec vu,PetscScalar *vfu,void *w)
83647c6ae99SBarry Smith {
83747c6ae99SBarry Smith   PetscErrorCode ierr;
83847c6ae99SBarry Smith   void           *u;
83947c6ae99SBarry Smith   DALocalInfo    info;
84047c6ae99SBarry Smith   MatStencil     stencil;
84147c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
84247c6ae99SBarry Smith 
84347c6ae99SBarry Smith   PetscFunctionBegin;
84447c6ae99SBarry Smith   ierr = DAGetLocalInfo(da,&info);CHKERRQ(ierr);
84547c6ae99SBarry Smith   ierr = DAVecGetArray(da,vu,&u);CHKERRQ(ierr);
84647c6ae99SBarry Smith 
84747c6ae99SBarry Smith   /* figure out stencil value from i */
84847c6ae99SBarry Smith   stencil.c = i % info.dof;
84947c6ae99SBarry Smith   if (stencil.c) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONG,"Point-block functions can only be called for the entire block");
85047c6ae99SBarry Smith   stencil.i = (i % (info.xm*info.dof))/info.dof;
85147c6ae99SBarry Smith   stencil.j = (i % (info.xm*info.ym*info.dof))/(info.xm*info.dof);
85247c6ae99SBarry Smith   stencil.k = i/(info.xm*info.ym*info.dof);
85347c6ae99SBarry Smith 
85447c6ae99SBarry Smith   ierr = (*dd->lfib)(&info,&stencil,u,vfu,w);CHKERRQ(ierr);
85547c6ae99SBarry Smith 
85647c6ae99SBarry Smith   ierr = DAVecRestoreArray(da,vu,&u);CHKERRQ(ierr);
85747c6ae99SBarry Smith   PetscFunctionReturn(0);
85847c6ae99SBarry Smith }
85947c6ae99SBarry Smith 
86047c6ae99SBarry Smith #if defined(new)
86147c6ae99SBarry Smith #undef __FUNCT__
86247c6ae99SBarry Smith #define __FUNCT__ "DAGetDiagonal_MFFD"
86347c6ae99SBarry Smith /*
86447c6ae99SBarry Smith   DAGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix where local
86547c6ae99SBarry Smith     function lives on a DA
86647c6ae99SBarry Smith 
86747c6ae99SBarry Smith         y ~= (F(u + ha) - F(u))/h,
86847c6ae99SBarry Smith   where F = nonlinear function, as set by SNESSetFunction()
86947c6ae99SBarry Smith         u = current iterate
87047c6ae99SBarry Smith         h = difference interval
87147c6ae99SBarry Smith */
87247c6ae99SBarry Smith PetscErrorCode DAGetDiagonal_MFFD(DA da,Vec U,Vec a)
87347c6ae99SBarry Smith {
87447c6ae99SBarry Smith   PetscScalar    h,*aa,*ww,v;
87547c6ae99SBarry Smith   PetscReal      epsilon = PETSC_SQRT_MACHINE_EPSILON,umin = 100.0*PETSC_SQRT_MACHINE_EPSILON;
87647c6ae99SBarry Smith   PetscErrorCode ierr;
87747c6ae99SBarry Smith   PetscInt       gI,nI;
87847c6ae99SBarry Smith   MatStencil     stencil;
87947c6ae99SBarry Smith   DALocalInfo    info;
88047c6ae99SBarry Smith 
88147c6ae99SBarry Smith   PetscFunctionBegin;
88247c6ae99SBarry Smith   ierr = (*ctx->func)(0,U,a,ctx->funcctx);CHKERRQ(ierr);
88347c6ae99SBarry Smith   ierr = (*ctx->funcisetbase)(U,ctx->funcctx);CHKERRQ(ierr);
88447c6ae99SBarry Smith 
88547c6ae99SBarry Smith   ierr = VecGetArray(U,&ww);CHKERRQ(ierr);
88647c6ae99SBarry Smith   ierr = VecGetArray(a,&aa);CHKERRQ(ierr);
88747c6ae99SBarry Smith 
88847c6ae99SBarry Smith   nI = 0;
88947c6ae99SBarry Smith     h  = ww[gI];
89047c6ae99SBarry Smith     if (h == 0.0) h = 1.0;
89147c6ae99SBarry Smith #if !defined(PETSC_USE_COMPLEX)
89247c6ae99SBarry Smith     if (h < umin && h >= 0.0)      h = umin;
89347c6ae99SBarry Smith     else if (h < 0.0 && h > -umin) h = -umin;
89447c6ae99SBarry Smith #else
89547c6ae99SBarry Smith     if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0)     h = umin;
89647c6ae99SBarry Smith     else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin;
89747c6ae99SBarry Smith #endif
89847c6ae99SBarry Smith     h     *= epsilon;
89947c6ae99SBarry Smith 
90047c6ae99SBarry Smith     ww[gI] += h;
90147c6ae99SBarry Smith     ierr          = (*ctx->funci)(i,w,&v,ctx->funcctx);CHKERRQ(ierr);
90247c6ae99SBarry Smith     aa[nI]  = (v - aa[nI])/h;
90347c6ae99SBarry Smith     ww[gI] -= h;
90447c6ae99SBarry Smith     nI++;
90547c6ae99SBarry Smith   }
90647c6ae99SBarry Smith   ierr = VecRestoreArray(U,&ww);CHKERRQ(ierr);
90747c6ae99SBarry Smith   ierr = VecRestoreArray(a,&aa);CHKERRQ(ierr);
90847c6ae99SBarry Smith   PetscFunctionReturn(0);
90947c6ae99SBarry Smith }
91047c6ae99SBarry Smith #endif
91147c6ae99SBarry Smith 
91247c6ae99SBarry Smith #if defined(PETSC_HAVE_ADIC)
91347c6ae99SBarry Smith EXTERN_C_BEGIN
91447c6ae99SBarry Smith #include "adic/ad_utils.h"
91547c6ae99SBarry Smith EXTERN_C_END
91647c6ae99SBarry Smith 
91747c6ae99SBarry Smith #undef __FUNCT__
91847c6ae99SBarry Smith #define __FUNCT__ "DAComputeJacobian1WithAdic"
91947c6ae99SBarry Smith /*@C
92047c6ae99SBarry Smith     DAComputeJacobian1WithAdic - Evaluates a adiC provided Jacobian function on each processor that
92147c6ae99SBarry Smith         share a DA
92247c6ae99SBarry Smith 
92347c6ae99SBarry Smith    Input Parameters:
92447c6ae99SBarry Smith +    da - the DA that defines the grid
92547c6ae99SBarry Smith .    vu - input vector (ghosted)
92647c6ae99SBarry Smith .    J - output matrix
92747c6ae99SBarry Smith -    w - any user data
92847c6ae99SBarry Smith 
92947c6ae99SBarry Smith    Level: advanced
93047c6ae99SBarry Smith 
93147c6ae99SBarry Smith     Notes: Does NOT do ghost updates on vu upon entry
93247c6ae99SBarry Smith 
93347c6ae99SBarry Smith .seealso: DAFormFunction1()
93447c6ae99SBarry Smith 
93547c6ae99SBarry Smith @*/
93647c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DAComputeJacobian1WithAdic(DA da,Vec vu,Mat J,void *w)
93747c6ae99SBarry Smith {
93847c6ae99SBarry Smith   PetscErrorCode ierr;
93947c6ae99SBarry Smith   PetscInt       gtdof,tdof;
94047c6ae99SBarry Smith   PetscScalar    *ustart;
94147c6ae99SBarry Smith   DALocalInfo    info;
94247c6ae99SBarry Smith   void           *ad_u,*ad_f,*ad_ustart,*ad_fstart;
94347c6ae99SBarry Smith   ISColoring     iscoloring;
94447c6ae99SBarry Smith 
94547c6ae99SBarry Smith   PetscFunctionBegin;
94647c6ae99SBarry Smith   ierr = DAGetLocalInfo(da,&info);CHKERRQ(ierr);
94747c6ae99SBarry Smith 
94847c6ae99SBarry Smith   PetscADResetIndep();
94947c6ae99SBarry Smith 
95047c6ae99SBarry Smith   /* get space for derivative objects.  */
95147c6ae99SBarry Smith   ierr = DAGetAdicArray(da,PETSC_TRUE,&ad_u,&ad_ustart,&gtdof);CHKERRQ(ierr);
95247c6ae99SBarry Smith   ierr = DAGetAdicArray(da,PETSC_FALSE,&ad_f,&ad_fstart,&tdof);CHKERRQ(ierr);
95347c6ae99SBarry Smith   ierr = VecGetArray(vu,&ustart);CHKERRQ(ierr);
95447c6ae99SBarry Smith   ierr = DAGetColoring(da,IS_COLORING_GHOSTED,MATAIJ,&iscoloring);CHKERRQ(ierr);
95547c6ae99SBarry Smith 
95647c6ae99SBarry Smith   PetscADSetValueAndColor(ad_ustart,gtdof,iscoloring->colors,ustart);
95747c6ae99SBarry Smith 
95847c6ae99SBarry Smith   ierr = VecRestoreArray(vu,&ustart);CHKERRQ(ierr);
95947c6ae99SBarry Smith   ierr = ISColoringDestroy(iscoloring);CHKERRQ(ierr);
96047c6ae99SBarry Smith   ierr = PetscADIncrementTotalGradSize(iscoloring->n);CHKERRQ(ierr);
96147c6ae99SBarry Smith   PetscADSetIndepDone();
96247c6ae99SBarry Smith 
96347c6ae99SBarry Smith   ierr = PetscLogEventBegin(DA_LocalADFunction,0,0,0,0);CHKERRQ(ierr);
96447c6ae99SBarry Smith   ierr = (*dd->adic_lf)(&info,ad_u,ad_f,w);CHKERRQ(ierr);
96547c6ae99SBarry Smith   ierr = PetscLogEventEnd(DA_LocalADFunction,0,0,0,0);CHKERRQ(ierr);
96647c6ae99SBarry Smith 
96747c6ae99SBarry Smith   /* stick the values into the matrix */
96847c6ae99SBarry Smith   ierr = MatSetValuesAdic(J,(PetscScalar**)ad_fstart);CHKERRQ(ierr);
96947c6ae99SBarry Smith 
97047c6ae99SBarry Smith   /* return space for derivative objects.  */
97147c6ae99SBarry Smith   ierr = DARestoreAdicArray(da,PETSC_TRUE,&ad_u,&ad_ustart,&gtdof);CHKERRQ(ierr);
97247c6ae99SBarry Smith   ierr = DARestoreAdicArray(da,PETSC_FALSE,&ad_f,&ad_fstart,&tdof);CHKERRQ(ierr);
97347c6ae99SBarry Smith   PetscFunctionReturn(0);
97447c6ae99SBarry Smith }
97547c6ae99SBarry Smith 
97647c6ae99SBarry Smith #undef __FUNCT__
97747c6ae99SBarry Smith #define __FUNCT__ "DAMultiplyByJacobian1WithAdic"
97847c6ae99SBarry Smith /*@C
97947c6ae99SBarry Smith     DAMultiplyByJacobian1WithAdic - Applies an ADIC-provided Jacobian function to a vector on
98047c6ae99SBarry Smith     each processor that shares a DA.
98147c6ae99SBarry Smith 
98247c6ae99SBarry Smith     Input Parameters:
98347c6ae99SBarry Smith +   da - the DA that defines the grid
98447c6ae99SBarry Smith .   vu - Jacobian is computed at this point (ghosted)
98547c6ae99SBarry Smith .   v - product is done on this vector (ghosted)
98647c6ae99SBarry Smith .   fu - output vector = J(vu)*v (not ghosted)
98747c6ae99SBarry Smith -   w - any user data
98847c6ae99SBarry Smith 
98947c6ae99SBarry Smith     Notes:
99047c6ae99SBarry Smith     This routine does NOT do ghost updates on vu upon entry.
99147c6ae99SBarry Smith 
99247c6ae99SBarry Smith    Level: advanced
99347c6ae99SBarry Smith 
99447c6ae99SBarry Smith .seealso: DAFormFunction1()
99547c6ae99SBarry Smith 
99647c6ae99SBarry Smith @*/
99747c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DAMultiplyByJacobian1WithAdic(DA da,Vec vu,Vec v,Vec f,void *w)
99847c6ae99SBarry Smith {
99947c6ae99SBarry Smith   PetscErrorCode ierr;
100047c6ae99SBarry Smith   PetscInt       i,gtdof,tdof;
100147c6ae99SBarry Smith   PetscScalar    *avu,*av,*af,*ad_vustart,*ad_fstart;
100247c6ae99SBarry Smith   DALocalInfo    info;
100347c6ae99SBarry Smith   void           *ad_vu,*ad_f;
100447c6ae99SBarry Smith 
100547c6ae99SBarry Smith   PetscFunctionBegin;
100647c6ae99SBarry Smith   ierr = DAGetLocalInfo(da,&info);CHKERRQ(ierr);
100747c6ae99SBarry Smith 
100847c6ae99SBarry Smith   /* get space for derivative objects.  */
100947c6ae99SBarry Smith   ierr = DAGetAdicMFArray(da,PETSC_TRUE,&ad_vu,&ad_vustart,&gtdof);CHKERRQ(ierr);
101047c6ae99SBarry Smith   ierr = DAGetAdicMFArray(da,PETSC_FALSE,&ad_f,&ad_fstart,&tdof);CHKERRQ(ierr);
101147c6ae99SBarry Smith 
101247c6ae99SBarry Smith   /* copy input vector into derivative object */
101347c6ae99SBarry Smith   ierr = VecGetArray(vu,&avu);CHKERRQ(ierr);
101447c6ae99SBarry Smith   ierr = VecGetArray(v,&av);CHKERRQ(ierr);
101547c6ae99SBarry Smith   for (i=0; i<gtdof; i++) {
101647c6ae99SBarry Smith     ad_vustart[2*i]   = avu[i];
101747c6ae99SBarry Smith     ad_vustart[2*i+1] = av[i];
101847c6ae99SBarry Smith   }
101947c6ae99SBarry Smith   ierr = VecRestoreArray(vu,&avu);CHKERRQ(ierr);
102047c6ae99SBarry Smith   ierr = VecRestoreArray(v,&av);CHKERRQ(ierr);
102147c6ae99SBarry Smith 
102247c6ae99SBarry Smith   PetscADResetIndep();
102347c6ae99SBarry Smith   ierr = PetscADIncrementTotalGradSize(1);CHKERRQ(ierr);
102447c6ae99SBarry Smith   PetscADSetIndepDone();
102547c6ae99SBarry Smith 
102647c6ae99SBarry Smith   ierr = (*dd->adicmf_lf)(&info,ad_vu,ad_f,w);CHKERRQ(ierr);
102747c6ae99SBarry Smith 
102847c6ae99SBarry Smith   /* stick the values into the vector */
102947c6ae99SBarry Smith   ierr = VecGetArray(f,&af);CHKERRQ(ierr);
103047c6ae99SBarry Smith   for (i=0; i<tdof; i++) {
103147c6ae99SBarry Smith     af[i] = ad_fstart[2*i+1];
103247c6ae99SBarry Smith   }
103347c6ae99SBarry Smith   ierr = VecRestoreArray(f,&af);CHKERRQ(ierr);
103447c6ae99SBarry Smith 
103547c6ae99SBarry Smith   /* return space for derivative objects.  */
103647c6ae99SBarry Smith   ierr = DARestoreAdicMFArray(da,PETSC_TRUE,&ad_vu,&ad_vustart,&gtdof);CHKERRQ(ierr);
103747c6ae99SBarry Smith   ierr = DARestoreAdicMFArray(da,PETSC_FALSE,&ad_f,&ad_fstart,&tdof);CHKERRQ(ierr);
103847c6ae99SBarry Smith   PetscFunctionReturn(0);
103947c6ae99SBarry Smith }
104047c6ae99SBarry Smith #endif
104147c6ae99SBarry Smith 
104247c6ae99SBarry Smith #undef __FUNCT__
104347c6ae99SBarry Smith #define __FUNCT__ "DAComputeJacobian1"
104447c6ae99SBarry Smith /*@
104547c6ae99SBarry Smith     DAComputeJacobian1 - Evaluates a local Jacobian function on each processor that
104647c6ae99SBarry Smith         share a DA
104747c6ae99SBarry Smith 
104847c6ae99SBarry Smith    Input Parameters:
104947c6ae99SBarry Smith +    da - the DA that defines the grid
105047c6ae99SBarry Smith .    vu - input vector (ghosted)
105147c6ae99SBarry Smith .    J - output matrix
105247c6ae99SBarry Smith -    w - any user data
105347c6ae99SBarry Smith 
105447c6ae99SBarry Smith     Notes: Does NOT do ghost updates on vu upon entry
105547c6ae99SBarry Smith 
105647c6ae99SBarry Smith     Level: advanced
105747c6ae99SBarry Smith 
105847c6ae99SBarry Smith .seealso: DAFormFunction1()
105947c6ae99SBarry Smith 
106047c6ae99SBarry Smith @*/
106147c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DAComputeJacobian1(DA da,Vec vu,Mat J,void *w)
106247c6ae99SBarry Smith {
106347c6ae99SBarry Smith   PetscErrorCode ierr;
106447c6ae99SBarry Smith   void           *u;
106547c6ae99SBarry Smith   DALocalInfo    info;
106647c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
106747c6ae99SBarry Smith 
106847c6ae99SBarry Smith   PetscFunctionBegin;
106947c6ae99SBarry Smith   ierr = DAGetLocalInfo(da,&info);CHKERRQ(ierr);
107047c6ae99SBarry Smith   ierr = DAVecGetArray(da,vu,&u);CHKERRQ(ierr);
107147c6ae99SBarry Smith   ierr = (*dd->lj)(&info,u,J,w);CHKERRQ(ierr);
107247c6ae99SBarry Smith   ierr = DAVecRestoreArray(da,vu,&u);CHKERRQ(ierr);
107347c6ae99SBarry Smith   PetscFunctionReturn(0);
107447c6ae99SBarry Smith }
107547c6ae99SBarry Smith 
107647c6ae99SBarry Smith 
107747c6ae99SBarry Smith #undef __FUNCT__
107847c6ae99SBarry Smith #define __FUNCT__ "DAComputeJacobian1WithAdifor"
107947c6ae99SBarry Smith /*
108047c6ae99SBarry Smith     DAComputeJacobian1WithAdifor - Evaluates a ADIFOR provided Jacobian local function on each processor that
108147c6ae99SBarry Smith         share a DA
108247c6ae99SBarry Smith 
108347c6ae99SBarry Smith    Input Parameters:
108447c6ae99SBarry Smith +    da - the DA that defines the grid
108547c6ae99SBarry Smith .    vu - input vector (ghosted)
108647c6ae99SBarry Smith .    J - output matrix
108747c6ae99SBarry Smith -    w - any user data
108847c6ae99SBarry Smith 
108947c6ae99SBarry Smith     Notes: Does NOT do ghost updates on vu upon entry
109047c6ae99SBarry Smith 
109147c6ae99SBarry Smith .seealso: DAFormFunction1()
109247c6ae99SBarry Smith 
109347c6ae99SBarry Smith */
109447c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DAComputeJacobian1WithAdifor(DA da,Vec vu,Mat J,void *w)
109547c6ae99SBarry Smith {
109647c6ae99SBarry Smith   PetscErrorCode  ierr;
109747c6ae99SBarry Smith   PetscInt        i,Nc,N;
109847c6ae99SBarry Smith   ISColoringValue *color;
109947c6ae99SBarry Smith   DALocalInfo     info;
110047c6ae99SBarry Smith   PetscScalar     *u,*g_u,*g_f,*f = 0,*p_u;
110147c6ae99SBarry Smith   ISColoring      iscoloring;
110247c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
110347c6ae99SBarry Smith   void            (*lf)(PetscInt*,DALocalInfo*,PetscScalar*,PetscScalar*,PetscInt*,PetscScalar*,PetscScalar*,PetscInt*,void*,PetscErrorCode*) =
110447c6ae99SBarry Smith                   (void (*)(PetscInt*,DALocalInfo*,PetscScalar*,PetscScalar*,PetscInt*,PetscScalar*,PetscScalar*,PetscInt*,void*,PetscErrorCode*))*dd->adifor_lf;
110547c6ae99SBarry Smith 
110647c6ae99SBarry Smith   PetscFunctionBegin;
110747c6ae99SBarry Smith   ierr = DAGetColoring(da,IS_COLORING_GHOSTED,MATAIJ,&iscoloring);CHKERRQ(ierr);
110847c6ae99SBarry Smith   Nc   = iscoloring->n;
110947c6ae99SBarry Smith   ierr = DAGetLocalInfo(da,&info);CHKERRQ(ierr);
111047c6ae99SBarry Smith   N    = info.gxm*info.gym*info.gzm*info.dof;
111147c6ae99SBarry Smith 
111247c6ae99SBarry Smith   /* get space for derivative objects.  */
111347c6ae99SBarry Smith   ierr  = PetscMalloc(Nc*info.gxm*info.gym*info.gzm*info.dof*sizeof(PetscScalar),&g_u);CHKERRQ(ierr);
111447c6ae99SBarry Smith   ierr  = PetscMemzero(g_u,Nc*info.gxm*info.gym*info.gzm*info.dof*sizeof(PetscScalar));CHKERRQ(ierr);
111547c6ae99SBarry Smith   p_u   = g_u;
111647c6ae99SBarry Smith   color = iscoloring->colors;
111747c6ae99SBarry Smith   for (i=0; i<N; i++) {
111847c6ae99SBarry Smith     p_u[*color++] = 1.0;
111947c6ae99SBarry Smith     p_u          += Nc;
112047c6ae99SBarry Smith   }
112147c6ae99SBarry Smith   ierr = ISColoringDestroy(iscoloring);CHKERRQ(ierr);
112247c6ae99SBarry 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);
112347c6ae99SBarry Smith 
112447c6ae99SBarry Smith   /* Seed the input array g_u with coloring information */
112547c6ae99SBarry Smith 
112647c6ae99SBarry Smith   ierr = VecGetArray(vu,&u);CHKERRQ(ierr);
112747c6ae99SBarry Smith   (lf)(&Nc,&info,u,g_u,&Nc,f,g_f,&Nc,w,&ierr);CHKERRQ(ierr);
112847c6ae99SBarry Smith   ierr = VecRestoreArray(vu,&u);CHKERRQ(ierr);
112947c6ae99SBarry Smith 
113047c6ae99SBarry Smith   /* stick the values into the matrix */
113147c6ae99SBarry Smith   /* PetscScalarView(Nc*info.xm*info.ym,g_f,0); */
113247c6ae99SBarry Smith   ierr = MatSetValuesAdifor(J,Nc,g_f);CHKERRQ(ierr);
113347c6ae99SBarry Smith 
113447c6ae99SBarry Smith   /* return space for derivative objects.  */
113547c6ae99SBarry Smith   ierr = PetscFree(g_u);CHKERRQ(ierr);
113647c6ae99SBarry Smith   ierr = PetscFree2(g_f,f);CHKERRQ(ierr);
113747c6ae99SBarry Smith   PetscFunctionReturn(0);
113847c6ae99SBarry Smith }
113947c6ae99SBarry Smith 
114047c6ae99SBarry Smith #undef __FUNCT__
114147c6ae99SBarry Smith #define __FUNCT__ "DAFormJacobianLocal"
114247c6ae99SBarry Smith /*@C
114347c6ae99SBarry Smith    DAFormjacobianLocal - This is a universal Jacobian evaluation routine for
114447c6ae99SBarry Smith    a local DA function.
114547c6ae99SBarry Smith 
114647c6ae99SBarry Smith    Collective on DA
114747c6ae99SBarry Smith 
114847c6ae99SBarry Smith    Input Parameters:
114947c6ae99SBarry Smith +  da - the DA context
115047c6ae99SBarry Smith .  func - The local function
115147c6ae99SBarry Smith .  X - input vector
115247c6ae99SBarry Smith .  J - Jacobian matrix
115347c6ae99SBarry Smith -  ctx - A user context
115447c6ae99SBarry Smith 
115547c6ae99SBarry Smith    Level: intermediate
115647c6ae99SBarry Smith 
115747c6ae99SBarry Smith .seealso: DASetLocalFunction(), DASetLocalJacobian(), DASetLocalAdicFunction(), DASetLocalAdicMFFunction(),
115847c6ae99SBarry Smith           SNESSetFunction(), SNESSetJacobian()
115947c6ae99SBarry Smith 
116047c6ae99SBarry Smith @*/
116147c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DAFormJacobianLocal(DA da, DALocalFunction1 func, Vec X, Mat J, void *ctx)
116247c6ae99SBarry Smith {
116347c6ae99SBarry Smith   Vec            localX;
116447c6ae99SBarry Smith   DALocalInfo    info;
116547c6ae99SBarry Smith   void          *u;
116647c6ae99SBarry Smith   PetscErrorCode ierr;
116747c6ae99SBarry Smith 
116847c6ae99SBarry Smith   PetscFunctionBegin;
116947c6ae99SBarry Smith   ierr = DAGetLocalVector(da,&localX);CHKERRQ(ierr);
117047c6ae99SBarry Smith   /*
117147c6ae99SBarry Smith      Scatter ghost points to local vector, using the 2-step process
117247c6ae99SBarry Smith         DAGlobalToLocalBegin(), DAGlobalToLocalEnd().
117347c6ae99SBarry Smith   */
117447c6ae99SBarry Smith   ierr = DAGlobalToLocalBegin(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
117547c6ae99SBarry Smith   ierr = DAGlobalToLocalEnd(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
117647c6ae99SBarry Smith   ierr = DAGetLocalInfo(da,&info);CHKERRQ(ierr);
117747c6ae99SBarry Smith   ierr = DAVecGetArray(da,localX,&u);CHKERRQ(ierr);
117847c6ae99SBarry Smith   ierr = (*func)(&info,u,J,ctx);
117947c6ae99SBarry Smith   if (PetscExceptionValue(ierr)) {
118047c6ae99SBarry Smith     PetscErrorCode pierr = DAVecRestoreArray(da,localX,&u);CHKERRQ(pierr);
118147c6ae99SBarry Smith   }
118247c6ae99SBarry Smith   CHKERRQ(ierr);
118347c6ae99SBarry Smith   ierr = DAVecRestoreArray(da,localX,&u);CHKERRQ(ierr);
118447c6ae99SBarry Smith   if (PetscExceptionValue(ierr)) {
118547c6ae99SBarry Smith     PetscErrorCode pierr = DARestoreLocalVector(da,&localX);CHKERRQ(pierr);
118647c6ae99SBarry Smith   }
118747c6ae99SBarry Smith   CHKERRQ(ierr);
118847c6ae99SBarry Smith   ierr = DARestoreLocalVector(da,&localX);CHKERRQ(ierr);
118947c6ae99SBarry Smith   PetscFunctionReturn(0);
119047c6ae99SBarry Smith }
119147c6ae99SBarry Smith 
119247c6ae99SBarry Smith #undef __FUNCT__
119347c6ae99SBarry Smith #define __FUNCT__ "DAMultiplyByJacobian1WithAD"
119447c6ae99SBarry Smith /*@C
119547c6ae99SBarry Smith     DAMultiplyByJacobian1WithAD - Applies a Jacobian function supplied by ADIFOR or ADIC
119647c6ae99SBarry Smith     to a vector on each processor that shares a DA.
119747c6ae99SBarry Smith 
119847c6ae99SBarry Smith    Input Parameters:
119947c6ae99SBarry Smith +    da - the DA that defines the grid
120047c6ae99SBarry Smith .    vu - Jacobian is computed at this point (ghosted)
120147c6ae99SBarry Smith .    v - product is done on this vector (ghosted)
120247c6ae99SBarry Smith .    fu - output vector = J(vu)*v (not ghosted)
120347c6ae99SBarry Smith -    w - any user data
120447c6ae99SBarry Smith 
120547c6ae99SBarry Smith     Notes:
120647c6ae99SBarry Smith     This routine does NOT do ghost updates on vu and v upon entry.
120747c6ae99SBarry Smith 
120847c6ae99SBarry Smith     Automatically calls DAMultiplyByJacobian1WithAdifor() or DAMultiplyByJacobian1WithAdic()
120947c6ae99SBarry Smith     depending on whether DASetLocalAdicMFFunction() or DASetLocalAdiforMFFunction() was called.
121047c6ae99SBarry Smith 
121147c6ae99SBarry Smith    Level: advanced
121247c6ae99SBarry Smith 
121347c6ae99SBarry Smith .seealso: DAFormFunction1(), DAMultiplyByJacobian1WithAdifor(), DAMultiplyByJacobian1WithAdic()
121447c6ae99SBarry Smith 
121547c6ae99SBarry Smith @*/
121647c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DAMultiplyByJacobian1WithAD(DA da,Vec u,Vec v,Vec f,void *w)
121747c6ae99SBarry Smith {
121847c6ae99SBarry Smith   PetscErrorCode ierr;
121947c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
122047c6ae99SBarry Smith 
122147c6ae99SBarry Smith   PetscFunctionBegin;
122247c6ae99SBarry Smith   if (dd->adicmf_lf) {
122347c6ae99SBarry Smith #if defined(PETSC_HAVE_ADIC)
122447c6ae99SBarry Smith     ierr = DAMultiplyByJacobian1WithAdic(da,u,v,f,w);CHKERRQ(ierr);
122547c6ae99SBarry Smith #else
122647c6ae99SBarry Smith     SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP_SYS,"Requires ADIC to be installed and cannot use complex numbers");
122747c6ae99SBarry Smith #endif
122847c6ae99SBarry Smith   } else if (dd->adiformf_lf) {
122947c6ae99SBarry Smith     ierr = DAMultiplyByJacobian1WithAdifor(da,u,v,f,w);CHKERRQ(ierr);
123047c6ae99SBarry Smith   } else {
123147c6ae99SBarry Smith     SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ORDER,"Must call DASetLocalAdiforMFFunction() or DASetLocalAdicMFFunction() before using");
123247c6ae99SBarry Smith   }
123347c6ae99SBarry Smith   PetscFunctionReturn(0);
123447c6ae99SBarry Smith }
123547c6ae99SBarry Smith 
123647c6ae99SBarry Smith 
123747c6ae99SBarry Smith #undef __FUNCT__
123847c6ae99SBarry Smith #define __FUNCT__ "DAMultiplyByJacobian1WithAdifor"
123947c6ae99SBarry Smith /*@C
124047c6ae99SBarry Smith     DAMultiplyByJacobian1WithAdifor - Applies a ADIFOR provided Jacobian function on each processor that
124147c6ae99SBarry Smith         share a DA to a vector
124247c6ae99SBarry Smith 
124347c6ae99SBarry Smith    Input Parameters:
124447c6ae99SBarry Smith +    da - the DA that defines the grid
124547c6ae99SBarry Smith .    vu - Jacobian is computed at this point (ghosted)
124647c6ae99SBarry Smith .    v - product is done on this vector (ghosted)
124747c6ae99SBarry Smith .    fu - output vector = J(vu)*v (not ghosted)
124847c6ae99SBarry Smith -    w - any user data
124947c6ae99SBarry Smith 
125047c6ae99SBarry Smith     Notes: Does NOT do ghost updates on vu and v upon entry
125147c6ae99SBarry Smith 
125247c6ae99SBarry Smith    Level: advanced
125347c6ae99SBarry Smith 
125447c6ae99SBarry Smith .seealso: DAFormFunction1()
125547c6ae99SBarry Smith 
125647c6ae99SBarry Smith @*/
125747c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DAMultiplyByJacobian1WithAdifor(DA da,Vec u,Vec v,Vec f,void *w)
125847c6ae99SBarry Smith {
125947c6ae99SBarry Smith   PetscErrorCode ierr;
126047c6ae99SBarry Smith   PetscScalar    *au,*av,*af,*awork;
126147c6ae99SBarry Smith   Vec            work;
126247c6ae99SBarry Smith   DALocalInfo    info;
126347c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
126447c6ae99SBarry Smith   void           (*lf)(DALocalInfo*,PetscScalar*,PetscScalar*,PetscScalar*,PetscScalar*,void*,PetscErrorCode*) =
126547c6ae99SBarry Smith                  (void (*)(DALocalInfo*,PetscScalar*,PetscScalar*,PetscScalar*,PetscScalar*,void*,PetscErrorCode*))*dd->adiformf_lf;
126647c6ae99SBarry Smith 
126747c6ae99SBarry Smith   PetscFunctionBegin;
126847c6ae99SBarry Smith   ierr = DAGetLocalInfo(da,&info);CHKERRQ(ierr);
126947c6ae99SBarry Smith 
127047c6ae99SBarry Smith   ierr = DAGetGlobalVector(da,&work);CHKERRQ(ierr);
127147c6ae99SBarry Smith   ierr = VecGetArray(u,&au);CHKERRQ(ierr);
127247c6ae99SBarry Smith   ierr = VecGetArray(v,&av);CHKERRQ(ierr);
127347c6ae99SBarry Smith   ierr = VecGetArray(f,&af);CHKERRQ(ierr);
127447c6ae99SBarry Smith   ierr = VecGetArray(work,&awork);CHKERRQ(ierr);
127547c6ae99SBarry Smith   (lf)(&info,au,av,awork,af,w,&ierr);CHKERRQ(ierr);
127647c6ae99SBarry Smith   ierr = VecRestoreArray(u,&au);CHKERRQ(ierr);
127747c6ae99SBarry Smith   ierr = VecRestoreArray(v,&av);CHKERRQ(ierr);
127847c6ae99SBarry Smith   ierr = VecRestoreArray(f,&af);CHKERRQ(ierr);
127947c6ae99SBarry Smith   ierr = VecRestoreArray(work,&awork);CHKERRQ(ierr);
128047c6ae99SBarry Smith   ierr = DARestoreGlobalVector(da,&work);CHKERRQ(ierr);
128147c6ae99SBarry Smith 
128247c6ae99SBarry Smith   PetscFunctionReturn(0);
128347c6ae99SBarry Smith }
128447c6ae99SBarry Smith 
128547c6ae99SBarry Smith EXTERN_C_BEGIN
128647c6ae99SBarry Smith #undef __FUNCT__
128747c6ae99SBarry Smith #define __FUNCT__ "DACreate_2D"
128847c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DACreate_2D(DA da)
128947c6ae99SBarry Smith {
129047c6ae99SBarry Smith   DM_DA               *dd = (DM_DA*)da->data;
129147c6ae99SBarry Smith   const PetscInt       dim          = dd->dim;
129247c6ae99SBarry Smith   const PetscInt       M            = dd->M;
129347c6ae99SBarry Smith   const PetscInt       N            = dd->N;
129447c6ae99SBarry Smith   PetscInt             m            = dd->m;
129547c6ae99SBarry Smith   PetscInt             n            = dd->n;
129647c6ae99SBarry Smith   const PetscInt       dof          = dd->w;
129747c6ae99SBarry Smith   const PetscInt       s            = dd->s;
129847c6ae99SBarry Smith   const DAPeriodicType wrap         = dd->wrap;
129947c6ae99SBarry Smith   const DAStencilType  stencil_type = dd->stencil_type;
130047c6ae99SBarry Smith   PetscInt            *lx           = dd->lx;
130147c6ae99SBarry Smith   PetscInt            *ly           = dd->ly;
130247c6ae99SBarry Smith   MPI_Comm             comm;
130347c6ae99SBarry Smith   PetscMPIInt          rank,size;
130447c6ae99SBarry Smith   PetscInt             xs,xe,ys,ye,x,y,Xs,Xe,Ys,Ye,start,end;
130547c6ae99SBarry Smith   PetscInt             up,down,left,i,n0,n1,n2,n3,n5,n6,n7,n8,*idx,nn;
130647c6ae99SBarry Smith   PetscInt             xbase,*bases,*ldims,j,x_t,y_t,s_t,base,count;
130747c6ae99SBarry Smith   PetscInt             s_x,s_y; /* s proportionalized to w */
130847c6ae99SBarry Smith   PetscInt             sn0 = 0,sn2 = 0,sn6 = 0,sn8 = 0;
130947c6ae99SBarry Smith   Vec                  local,global;
131047c6ae99SBarry Smith   VecScatter           ltog,gtol;
131147c6ae99SBarry Smith   IS                   to,from;
131247c6ae99SBarry Smith   PetscErrorCode       ierr;
131347c6ae99SBarry Smith 
131447c6ae99SBarry Smith   PetscFunctionBegin;
131547c6ae99SBarry Smith #ifndef PETSC_USE_DYNAMIC_LIBRARIES
131647c6ae99SBarry Smith   ierr = DMInitializePackage(PETSC_NULL);CHKERRQ(ierr);
131747c6ae99SBarry Smith #endif
131847c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
131947c6ae99SBarry Smith 
132047c6ae99SBarry Smith   if (dim != PETSC_DECIDE && dim != 2) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Dimension should be 2: %D",dim);
132147c6ae99SBarry Smith   if (dof < 1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Must have 1 or more degrees of freedom per node: %D",dof);
132247c6ae99SBarry Smith   if (s < 0) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Stencil width cannot be negative: %D",s);
132347c6ae99SBarry Smith 
132447c6ae99SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
132547c6ae99SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
132647c6ae99SBarry Smith 
1327*5dff015aSBarry Smith   da->ops->getelements = DMGetElements_DM_2d_P1;
132847c6ae99SBarry Smith 
132947c6ae99SBarry Smith   dd->dim         = 2;
133047c6ae99SBarry Smith   dd->elementtype = DA_ELEMENT_P1;
133147c6ae99SBarry Smith   ierr = PetscMalloc(dof*sizeof(char*),&dd->fieldname);CHKERRQ(ierr);
133247c6ae99SBarry Smith   ierr = PetscMemzero(dd->fieldname,dof*sizeof(char*));CHKERRQ(ierr);
133347c6ae99SBarry Smith 
133447c6ae99SBarry Smith   if (m != PETSC_DECIDE) {
133547c6ae99SBarry Smith     if (m < 1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in X direction: %D",m);
133647c6ae99SBarry Smith     else if (m > size) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in X direction: %D %d",m,size);
133747c6ae99SBarry Smith   }
133847c6ae99SBarry Smith   if (n != PETSC_DECIDE) {
133947c6ae99SBarry Smith     if (n < 1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Y direction: %D",n);
134047c6ae99SBarry Smith     else if (n > size) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Y direction: %D %d",n,size);
134147c6ae99SBarry Smith   }
134247c6ae99SBarry Smith 
134347c6ae99SBarry Smith   if (m == PETSC_DECIDE || n == PETSC_DECIDE) {
134447c6ae99SBarry Smith     if (n != PETSC_DECIDE) {
134547c6ae99SBarry Smith       m = size/n;
134647c6ae99SBarry Smith     } else if (m != PETSC_DECIDE) {
134747c6ae99SBarry Smith       n = size/m;
134847c6ae99SBarry Smith     } else {
134947c6ae99SBarry Smith       /* try for squarish distribution */
135047c6ae99SBarry Smith       m = (PetscInt)(0.5 + sqrt(((double)M)*((double)size)/((double)N)));
135147c6ae99SBarry Smith       if (!m) m = 1;
135247c6ae99SBarry Smith       while (m > 0) {
135347c6ae99SBarry Smith 	n = size/m;
135447c6ae99SBarry Smith 	if (m*n == size) break;
135547c6ae99SBarry Smith 	m--;
135647c6ae99SBarry Smith       }
135747c6ae99SBarry Smith       if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;}
135847c6ae99SBarry Smith     }
135947c6ae99SBarry 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 ");
136047c6ae99SBarry Smith   } else if (m*n != size) SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"Given Bad partition");
136147c6ae99SBarry Smith 
136247c6ae99SBarry Smith   if (M < m) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in x direction is too fine! %D %D",M,m);
136347c6ae99SBarry Smith   if (N < n) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in y direction is too fine! %D %D",N,n);
136447c6ae99SBarry Smith 
136547c6ae99SBarry Smith   /*
136647c6ae99SBarry Smith      Determine locally owned region
136747c6ae99SBarry Smith      xs is the first local node number, x is the number of local nodes
136847c6ae99SBarry Smith   */
136947c6ae99SBarry Smith   if (!lx) {
137047c6ae99SBarry Smith     ierr = PetscMalloc(m*sizeof(PetscInt), &dd->lx);CHKERRQ(ierr);
137147c6ae99SBarry Smith     lx = dd->lx;
137247c6ae99SBarry Smith     for (i=0; i<m; i++) {
137347c6ae99SBarry Smith       lx[i] = M/m + ((M % m) > i);
137447c6ae99SBarry Smith     }
137547c6ae99SBarry Smith   }
137647c6ae99SBarry Smith   x  = lx[rank % m];
137747c6ae99SBarry Smith   xs = 0;
137847c6ae99SBarry Smith   for (i=0; i<(rank % m); i++) {
137947c6ae99SBarry Smith     xs += lx[i];
138047c6ae99SBarry Smith   }
138147c6ae99SBarry Smith #if defined(PETSC_USE_DEBUG)
138247c6ae99SBarry Smith   left = xs;
138347c6ae99SBarry Smith   for (i=(rank % m); i<m; i++) {
138447c6ae99SBarry Smith     left += lx[i];
138547c6ae99SBarry Smith   }
138647c6ae99SBarry 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);
138747c6ae99SBarry Smith #endif
138847c6ae99SBarry Smith 
138947c6ae99SBarry Smith   /*
139047c6ae99SBarry Smith      Determine locally owned region
139147c6ae99SBarry Smith      ys is the first local node number, y is the number of local nodes
139247c6ae99SBarry Smith   */
139347c6ae99SBarry Smith   if (!ly) {
139447c6ae99SBarry Smith     ierr = PetscMalloc(n*sizeof(PetscInt), &dd->ly);CHKERRQ(ierr);
139547c6ae99SBarry Smith     ly = dd->ly;
139647c6ae99SBarry Smith     for (i=0; i<n; i++) {
139747c6ae99SBarry Smith       ly[i] = N/n + ((N % n) > i);
139847c6ae99SBarry Smith     }
139947c6ae99SBarry Smith   }
140047c6ae99SBarry Smith   y  = ly[rank/m];
140147c6ae99SBarry Smith   ys = 0;
140247c6ae99SBarry Smith   for (i=0; i<(rank/m); i++) {
140347c6ae99SBarry Smith     ys += ly[i];
140447c6ae99SBarry Smith   }
140547c6ae99SBarry Smith #if defined(PETSC_USE_DEBUG)
140647c6ae99SBarry Smith   left = ys;
140747c6ae99SBarry Smith   for (i=(rank/m); i<n; i++) {
140847c6ae99SBarry Smith     left += ly[i];
140947c6ae99SBarry Smith   }
141047c6ae99SBarry 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);
141147c6ae99SBarry Smith #endif
141247c6ae99SBarry Smith 
141347c6ae99SBarry 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);
141447c6ae99SBarry 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);
141547c6ae99SBarry Smith   xe = xs + x;
141647c6ae99SBarry Smith   ye = ys + y;
141747c6ae99SBarry Smith 
141847c6ae99SBarry Smith   /* determine ghost region */
141947c6ae99SBarry Smith   /* Assume No Periodicity */
142047c6ae99SBarry Smith   if (xs-s > 0) Xs = xs - s; else Xs = 0;
142147c6ae99SBarry Smith   if (ys-s > 0) Ys = ys - s; else Ys = 0;
142247c6ae99SBarry Smith   if (xe+s <= M) Xe = xe + s; else Xe = M;
142347c6ae99SBarry Smith   if (ye+s <= N) Ye = ye + s; else Ye = N;
142447c6ae99SBarry Smith 
142547c6ae99SBarry Smith   /* X Periodic */
142647c6ae99SBarry Smith   if (DAXPeriodic(wrap)){
142747c6ae99SBarry Smith     Xs = xs - s;
142847c6ae99SBarry Smith     Xe = xe + s;
142947c6ae99SBarry Smith   }
143047c6ae99SBarry Smith 
143147c6ae99SBarry Smith   /* Y Periodic */
143247c6ae99SBarry Smith   if (DAYPeriodic(wrap)){
143347c6ae99SBarry Smith     Ys = ys - s;
143447c6ae99SBarry Smith     Ye = ye + s;
143547c6ae99SBarry Smith   }
143647c6ae99SBarry Smith 
143747c6ae99SBarry Smith   /* Resize all X parameters to reflect w */
143847c6ae99SBarry Smith   x   *= dof;
143947c6ae99SBarry Smith   xs  *= dof;
144047c6ae99SBarry Smith   xe  *= dof;
144147c6ae99SBarry Smith   Xs  *= dof;
144247c6ae99SBarry Smith   Xe  *= dof;
144347c6ae99SBarry Smith   s_x = s*dof;
144447c6ae99SBarry Smith   s_y = s;
144547c6ae99SBarry Smith 
144647c6ae99SBarry Smith   /* determine starting point of each processor */
144747c6ae99SBarry Smith   nn    = x*y;
144847c6ae99SBarry Smith   ierr  = PetscMalloc2(size+1,PetscInt,&bases,size,PetscInt,&ldims);CHKERRQ(ierr);
144947c6ae99SBarry Smith   ierr  = MPI_Allgather(&nn,1,MPIU_INT,ldims,1,MPIU_INT,comm);CHKERRQ(ierr);
145047c6ae99SBarry Smith   bases[0] = 0;
145147c6ae99SBarry Smith   for (i=1; i<=size; i++) {
145247c6ae99SBarry Smith     bases[i] = ldims[i-1];
145347c6ae99SBarry Smith   }
145447c6ae99SBarry Smith   for (i=1; i<=size; i++) {
145547c6ae99SBarry Smith     bases[i] += bases[i-1];
145647c6ae99SBarry Smith   }
145747c6ae99SBarry Smith 
145847c6ae99SBarry Smith   /* allocate the base parallel and sequential vectors */
145947c6ae99SBarry Smith   dd->Nlocal = x*y;
146047c6ae99SBarry Smith   ierr = VecCreateMPIWithArray(comm,dd->Nlocal,PETSC_DECIDE,0,&global);CHKERRQ(ierr);
146147c6ae99SBarry Smith   ierr = VecSetBlockSize(global,dof);CHKERRQ(ierr);
146247c6ae99SBarry Smith   dd->nlocal = (Xe-Xs)*(Ye-Ys);
146347c6ae99SBarry Smith   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,dd->nlocal,0,&local);CHKERRQ(ierr);
146447c6ae99SBarry Smith   ierr = VecSetBlockSize(local,dof);CHKERRQ(ierr);
146547c6ae99SBarry Smith 
146647c6ae99SBarry Smith 
146747c6ae99SBarry Smith   /* generate appropriate vector scatters */
146847c6ae99SBarry Smith   /* local to global inserts non-ghost point region into global */
146947c6ae99SBarry Smith   ierr = VecGetOwnershipRange(global,&start,&end);CHKERRQ(ierr);
147047c6ae99SBarry Smith   ierr = ISCreateStride(comm,x*y,start,1,&to);CHKERRQ(ierr);
147147c6ae99SBarry Smith 
147247c6ae99SBarry Smith   left  = xs - Xs; down  = ys - Ys; up    = down + y;
147347c6ae99SBarry Smith   ierr = PetscMalloc(x*(up - down)*sizeof(PetscInt),&idx);CHKERRQ(ierr);
147447c6ae99SBarry Smith   count = 0;
147547c6ae99SBarry Smith   for (i=down; i<up; i++) {
147647c6ae99SBarry Smith     for (j=0; j<x/dof; j++) {
147747c6ae99SBarry Smith       idx[count++] = (left + i*(Xe-Xs) + j*dof)/dof;
147847c6ae99SBarry Smith     }
147947c6ae99SBarry Smith   }
148047c6ae99SBarry Smith   ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&from);CHKERRQ(ierr);
148147c6ae99SBarry Smith 
148247c6ae99SBarry Smith   ierr = VecScatterCreate(local,from,global,to,&ltog);CHKERRQ(ierr);
148347c6ae99SBarry Smith   ierr = PetscLogObjectParent(dd,ltog);CHKERRQ(ierr);
148447c6ae99SBarry Smith   ierr = ISDestroy(from);CHKERRQ(ierr);
148547c6ae99SBarry Smith   ierr = ISDestroy(to);CHKERRQ(ierr);
148647c6ae99SBarry Smith 
148747c6ae99SBarry Smith   /* global to local must include ghost points */
148847c6ae99SBarry Smith   if (stencil_type == DA_STENCIL_BOX) {
148947c6ae99SBarry Smith     ierr = ISCreateStride(comm,(Xe-Xs)*(Ye-Ys),0,1,&to);CHKERRQ(ierr);
149047c6ae99SBarry Smith   } else {
149147c6ae99SBarry Smith     /* must drop into cross shape region */
149247c6ae99SBarry Smith     /*       ---------|
149347c6ae99SBarry Smith             |  top    |
149447c6ae99SBarry Smith          |---         ---|
149547c6ae99SBarry Smith          |   middle      |
149647c6ae99SBarry Smith          |               |
149747c6ae99SBarry Smith          ----         ----
149847c6ae99SBarry Smith             | bottom  |
149947c6ae99SBarry Smith             -----------
150047c6ae99SBarry Smith         Xs xs        xe  Xe */
150147c6ae99SBarry Smith     /* bottom */
150247c6ae99SBarry Smith     left  = xs - Xs; down = ys - Ys; up    = down + y;
150347c6ae99SBarry Smith     count = down*(xe-xs) + (up-down)*(Xe-Xs) + (Ye-Ys-up)*(xe-xs);
150447c6ae99SBarry Smith     ierr  = PetscMalloc(count*sizeof(PetscInt)/dof,&idx);CHKERRQ(ierr);
150547c6ae99SBarry Smith     count = 0;
150647c6ae99SBarry Smith     for (i=0; i<down; i++) {
150747c6ae99SBarry Smith       for (j=0; j<xe-xs; j += dof) {
150847c6ae99SBarry Smith         idx[count++] = left + i*(Xe-Xs) + j;
150947c6ae99SBarry Smith       }
151047c6ae99SBarry Smith     }
151147c6ae99SBarry Smith     /* middle */
151247c6ae99SBarry Smith     for (i=down; i<up; i++) {
151347c6ae99SBarry Smith       for (j=0; j<Xe-Xs; j += dof) {
151447c6ae99SBarry Smith         idx[count++] = i*(Xe-Xs) + j;
151547c6ae99SBarry Smith       }
151647c6ae99SBarry Smith     }
151747c6ae99SBarry Smith     /* top */
151847c6ae99SBarry Smith     for (i=up; i<Ye-Ys; i++) {
151947c6ae99SBarry Smith       for (j=0; j<xe-xs; j += dof) {
152047c6ae99SBarry Smith         idx[count++] = (left + i*(Xe-Xs) + j);
152147c6ae99SBarry Smith       }
152247c6ae99SBarry Smith     }
152347c6ae99SBarry Smith     for (i=0; i<count; i++) idx[i] = idx[i]/dof;
152447c6ae99SBarry Smith     ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr);
152547c6ae99SBarry Smith   }
152647c6ae99SBarry Smith 
152747c6ae99SBarry Smith 
152847c6ae99SBarry Smith   /* determine who lies on each side of us stored in    n6 n7 n8
152947c6ae99SBarry Smith                                                         n3    n5
153047c6ae99SBarry Smith                                                         n0 n1 n2
153147c6ae99SBarry Smith   */
153247c6ae99SBarry Smith 
153347c6ae99SBarry Smith   /* Assume the Non-Periodic Case */
153447c6ae99SBarry Smith   n1 = rank - m;
153547c6ae99SBarry Smith   if (rank % m) {
153647c6ae99SBarry Smith     n0 = n1 - 1;
153747c6ae99SBarry Smith   } else {
153847c6ae99SBarry Smith     n0 = -1;
153947c6ae99SBarry Smith   }
154047c6ae99SBarry Smith   if ((rank+1) % m) {
154147c6ae99SBarry Smith     n2 = n1 + 1;
154247c6ae99SBarry Smith     n5 = rank + 1;
154347c6ae99SBarry Smith     n8 = rank + m + 1; if (n8 >= m*n) n8 = -1;
154447c6ae99SBarry Smith   } else {
154547c6ae99SBarry Smith     n2 = -1; n5 = -1; n8 = -1;
154647c6ae99SBarry Smith   }
154747c6ae99SBarry Smith   if (rank % m) {
154847c6ae99SBarry Smith     n3 = rank - 1;
154947c6ae99SBarry Smith     n6 = n3 + m; if (n6 >= m*n) n6 = -1;
155047c6ae99SBarry Smith   } else {
155147c6ae99SBarry Smith     n3 = -1; n6 = -1;
155247c6ae99SBarry Smith   }
155347c6ae99SBarry Smith   n7 = rank + m; if (n7 >= m*n) n7 = -1;
155447c6ae99SBarry Smith 
155547c6ae99SBarry Smith 
155647c6ae99SBarry Smith   /* Modify for Periodic Cases */
155747c6ae99SBarry Smith   if (wrap == DA_YPERIODIC) {  /* Handle Top and Bottom Sides */
155847c6ae99SBarry Smith     if (n1 < 0) n1 = rank + m * (n-1);
155947c6ae99SBarry Smith     if (n7 < 0) n7 = rank - m * (n-1);
156047c6ae99SBarry Smith     if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1;
156147c6ae99SBarry Smith     if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1;
156247c6ae99SBarry Smith     if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1;
156347c6ae99SBarry Smith     if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1;
156447c6ae99SBarry Smith   } else if (wrap == DA_XPERIODIC) { /* Handle Left and Right Sides */
156547c6ae99SBarry Smith     if (n3 < 0) n3 = rank + (m-1);
156647c6ae99SBarry Smith     if (n5 < 0) n5 = rank - (m-1);
156747c6ae99SBarry Smith     if ((n1 >= 0) && (n0 < 0)) n0 = rank-1;
156847c6ae99SBarry Smith     if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1;
156947c6ae99SBarry Smith     if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1;
157047c6ae99SBarry Smith     if ((n7 >= 0) && (n8 < 0)) n8 = rank+1;
157147c6ae99SBarry Smith   } else if (wrap == DA_XYPERIODIC) {
157247c6ae99SBarry Smith 
157347c6ae99SBarry Smith     /* Handle all four corners */
157447c6ae99SBarry Smith     if ((n6 < 0) && (n7 < 0) && (n3 < 0)) n6 = m-1;
157547c6ae99SBarry Smith     if ((n8 < 0) && (n7 < 0) && (n5 < 0)) n8 = 0;
157647c6ae99SBarry Smith     if ((n2 < 0) && (n5 < 0) && (n1 < 0)) n2 = size-m;
157747c6ae99SBarry Smith     if ((n0 < 0) && (n3 < 0) && (n1 < 0)) n0 = size-1;
157847c6ae99SBarry Smith 
157947c6ae99SBarry Smith     /* Handle Top and Bottom Sides */
158047c6ae99SBarry Smith     if (n1 < 0) n1 = rank + m * (n-1);
158147c6ae99SBarry Smith     if (n7 < 0) n7 = rank - m * (n-1);
158247c6ae99SBarry Smith     if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1;
158347c6ae99SBarry Smith     if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1;
158447c6ae99SBarry Smith     if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1;
158547c6ae99SBarry Smith     if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1;
158647c6ae99SBarry Smith 
158747c6ae99SBarry Smith     /* Handle Left and Right Sides */
158847c6ae99SBarry Smith     if (n3 < 0) n3 = rank + (m-1);
158947c6ae99SBarry Smith     if (n5 < 0) n5 = rank - (m-1);
159047c6ae99SBarry Smith     if ((n1 >= 0) && (n0 < 0)) n0 = rank-1;
159147c6ae99SBarry Smith     if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1;
159247c6ae99SBarry Smith     if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1;
159347c6ae99SBarry Smith     if ((n7 >= 0) && (n8 < 0)) n8 = rank+1;
159447c6ae99SBarry Smith   }
159547c6ae99SBarry Smith   ierr = PetscMalloc(9*sizeof(PetscInt),&dd->neighbors);CHKERRQ(ierr);
159647c6ae99SBarry Smith   dd->neighbors[0] = n0;
159747c6ae99SBarry Smith   dd->neighbors[1] = n1;
159847c6ae99SBarry Smith   dd->neighbors[2] = n2;
159947c6ae99SBarry Smith   dd->neighbors[3] = n3;
160047c6ae99SBarry Smith   dd->neighbors[4] = rank;
160147c6ae99SBarry Smith   dd->neighbors[5] = n5;
160247c6ae99SBarry Smith   dd->neighbors[6] = n6;
160347c6ae99SBarry Smith   dd->neighbors[7] = n7;
160447c6ae99SBarry Smith   dd->neighbors[8] = n8;
160547c6ae99SBarry Smith 
160647c6ae99SBarry Smith   if (stencil_type == DA_STENCIL_STAR) {
160747c6ae99SBarry Smith     /* save corner processor numbers */
160847c6ae99SBarry Smith     sn0 = n0; sn2 = n2; sn6 = n6; sn8 = n8;
160947c6ae99SBarry Smith     n0 = n2 = n6 = n8 = -1;
161047c6ae99SBarry Smith   }
161147c6ae99SBarry Smith 
161247c6ae99SBarry Smith   ierr = PetscMalloc((x+2*s_x)*(y+2*s_y)*sizeof(PetscInt),&idx);CHKERRQ(ierr);
161347c6ae99SBarry Smith   ierr = PetscLogObjectMemory(da,(x+2*s_x)*(y+2*s_y)*sizeof(PetscInt));CHKERRQ(ierr);
161447c6ae99SBarry Smith   nn = 0;
161547c6ae99SBarry Smith 
161647c6ae99SBarry Smith   xbase = bases[rank];
161747c6ae99SBarry Smith   for (i=1; i<=s_y; i++) {
161847c6ae99SBarry Smith     if (n0 >= 0) { /* left below */
161947c6ae99SBarry Smith       x_t = lx[n0 % m]*dof;
162047c6ae99SBarry Smith       y_t = ly[(n0/m)];
162147c6ae99SBarry Smith       s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x;
162247c6ae99SBarry Smith       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
162347c6ae99SBarry Smith     }
162447c6ae99SBarry Smith     if (n1 >= 0) { /* directly below */
162547c6ae99SBarry Smith       x_t = x;
162647c6ae99SBarry Smith       y_t = ly[(n1/m)];
162747c6ae99SBarry Smith       s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t;
162847c6ae99SBarry Smith       for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
162947c6ae99SBarry Smith     }
163047c6ae99SBarry Smith     if (n2 >= 0) { /* right below */
163147c6ae99SBarry Smith       x_t = lx[n2 % m]*dof;
163247c6ae99SBarry Smith       y_t = ly[(n2/m)];
163347c6ae99SBarry Smith       s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t;
163447c6ae99SBarry Smith       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
163547c6ae99SBarry Smith     }
163647c6ae99SBarry Smith   }
163747c6ae99SBarry Smith 
163847c6ae99SBarry Smith   for (i=0; i<y; i++) {
163947c6ae99SBarry Smith     if (n3 >= 0) { /* directly left */
164047c6ae99SBarry Smith       x_t = lx[n3 % m]*dof;
164147c6ae99SBarry Smith       /* y_t = y; */
164247c6ae99SBarry Smith       s_t = bases[n3] + (i+1)*x_t - s_x;
164347c6ae99SBarry Smith       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
164447c6ae99SBarry Smith     }
164547c6ae99SBarry Smith 
164647c6ae99SBarry Smith     for (j=0; j<x; j++) { idx[nn++] = xbase++; } /* interior */
164747c6ae99SBarry Smith 
164847c6ae99SBarry Smith     if (n5 >= 0) { /* directly right */
164947c6ae99SBarry Smith       x_t = lx[n5 % m]*dof;
165047c6ae99SBarry Smith       /* y_t = y; */
165147c6ae99SBarry Smith       s_t = bases[n5] + (i)*x_t;
165247c6ae99SBarry Smith       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
165347c6ae99SBarry Smith     }
165447c6ae99SBarry Smith   }
165547c6ae99SBarry Smith 
165647c6ae99SBarry Smith   for (i=1; i<=s_y; i++) {
165747c6ae99SBarry Smith     if (n6 >= 0) { /* left above */
165847c6ae99SBarry Smith       x_t = lx[n6 % m]*dof;
165947c6ae99SBarry Smith       /* y_t = ly[(n6/m)]; */
166047c6ae99SBarry Smith       s_t = bases[n6] + (i)*x_t - s_x;
166147c6ae99SBarry Smith       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
166247c6ae99SBarry Smith     }
166347c6ae99SBarry Smith     if (n7 >= 0) { /* directly above */
166447c6ae99SBarry Smith       x_t = x;
166547c6ae99SBarry Smith       /* y_t = ly[(n7/m)]; */
166647c6ae99SBarry Smith       s_t = bases[n7] + (i-1)*x_t;
166747c6ae99SBarry Smith       for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
166847c6ae99SBarry Smith     }
166947c6ae99SBarry Smith     if (n8 >= 0) { /* right above */
167047c6ae99SBarry Smith       x_t = lx[n8 % m]*dof;
167147c6ae99SBarry Smith       /* y_t = ly[(n8/m)]; */
167247c6ae99SBarry Smith       s_t = bases[n8] + (i-1)*x_t;
167347c6ae99SBarry Smith       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
167447c6ae99SBarry Smith     }
167547c6ae99SBarry Smith   }
167647c6ae99SBarry Smith 
167747c6ae99SBarry Smith   base = bases[rank];
167847c6ae99SBarry Smith   {
167947c6ae99SBarry Smith     PetscInt nnn = nn/dof,*iidx;
168047c6ae99SBarry Smith     ierr = PetscMalloc(nnn*sizeof(PetscInt),&iidx);CHKERRQ(ierr);
168147c6ae99SBarry Smith     for (i=0; i<nnn; i++) {
168247c6ae99SBarry Smith       iidx[i] = idx[dof*i]/dof;
168347c6ae99SBarry Smith     }
168447c6ae99SBarry Smith     ierr = ISCreateBlock(comm,dof,nnn,iidx,PETSC_OWN_POINTER,&from);CHKERRQ(ierr);
168547c6ae99SBarry Smith   }
168647c6ae99SBarry Smith   ierr = VecScatterCreate(global,from,local,to,&gtol);CHKERRQ(ierr);
168747c6ae99SBarry Smith   ierr = PetscLogObjectParent(da,gtol);CHKERRQ(ierr);
168847c6ae99SBarry Smith   ierr = ISDestroy(to);CHKERRQ(ierr);
168947c6ae99SBarry Smith   ierr = ISDestroy(from);CHKERRQ(ierr);
169047c6ae99SBarry Smith 
169147c6ae99SBarry Smith   if (stencil_type == DA_STENCIL_STAR) {
169247c6ae99SBarry Smith     /*
169347c6ae99SBarry Smith         Recompute the local to global mappings, this time keeping the
169447c6ae99SBarry Smith       information about the cross corner processor numbers.
169547c6ae99SBarry Smith     */
169647c6ae99SBarry Smith     n0 = sn0; n2 = sn2; n6 = sn6; n8 = sn8;
169747c6ae99SBarry Smith     nn = 0;
169847c6ae99SBarry Smith     xbase = bases[rank];
169947c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
170047c6ae99SBarry Smith       if (n0 >= 0) { /* left below */
170147c6ae99SBarry Smith         x_t = lx[n0 % m]*dof;
170247c6ae99SBarry Smith         y_t = ly[(n0/m)];
170347c6ae99SBarry Smith         s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x;
170447c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
170547c6ae99SBarry Smith       }
170647c6ae99SBarry Smith       if (n1 >= 0) { /* directly below */
170747c6ae99SBarry Smith         x_t = x;
170847c6ae99SBarry Smith         y_t = ly[(n1/m)];
170947c6ae99SBarry Smith         s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t;
171047c6ae99SBarry Smith         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
171147c6ae99SBarry Smith       }
171247c6ae99SBarry Smith       if (n2 >= 0) { /* right below */
171347c6ae99SBarry Smith         x_t = lx[n2 % m]*dof;
171447c6ae99SBarry Smith         y_t = ly[(n2/m)];
171547c6ae99SBarry Smith         s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t;
171647c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
171747c6ae99SBarry Smith       }
171847c6ae99SBarry Smith     }
171947c6ae99SBarry Smith 
172047c6ae99SBarry Smith     for (i=0; i<y; i++) {
172147c6ae99SBarry Smith       if (n3 >= 0) { /* directly left */
172247c6ae99SBarry Smith         x_t = lx[n3 % m]*dof;
172347c6ae99SBarry Smith         /* y_t = y; */
172447c6ae99SBarry Smith         s_t = bases[n3] + (i+1)*x_t - s_x;
172547c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
172647c6ae99SBarry Smith       }
172747c6ae99SBarry Smith 
172847c6ae99SBarry Smith       for (j=0; j<x; j++) { idx[nn++] = xbase++; } /* interior */
172947c6ae99SBarry Smith 
173047c6ae99SBarry Smith       if (n5 >= 0) { /* directly right */
173147c6ae99SBarry Smith         x_t = lx[n5 % m]*dof;
173247c6ae99SBarry Smith         /* y_t = y; */
173347c6ae99SBarry Smith         s_t = bases[n5] + (i)*x_t;
173447c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
173547c6ae99SBarry Smith       }
173647c6ae99SBarry Smith     }
173747c6ae99SBarry Smith 
173847c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
173947c6ae99SBarry Smith       if (n6 >= 0) { /* left above */
174047c6ae99SBarry Smith         x_t = lx[n6 % m]*dof;
174147c6ae99SBarry Smith         /* y_t = ly[(n6/m)]; */
174247c6ae99SBarry Smith         s_t = bases[n6] + (i)*x_t - s_x;
174347c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
174447c6ae99SBarry Smith       }
174547c6ae99SBarry Smith       if (n7 >= 0) { /* directly above */
174647c6ae99SBarry Smith         x_t = x;
174747c6ae99SBarry Smith         /* y_t = ly[(n7/m)]; */
174847c6ae99SBarry Smith         s_t = bases[n7] + (i-1)*x_t;
174947c6ae99SBarry Smith         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
175047c6ae99SBarry Smith       }
175147c6ae99SBarry Smith       if (n8 >= 0) { /* right above */
175247c6ae99SBarry Smith         x_t = lx[n8 % m]*dof;
175347c6ae99SBarry Smith         /* y_t = ly[(n8/m)]; */
175447c6ae99SBarry Smith         s_t = bases[n8] + (i-1)*x_t;
175547c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
175647c6ae99SBarry Smith       }
175747c6ae99SBarry Smith     }
175847c6ae99SBarry Smith   }
175947c6ae99SBarry Smith   ierr = PetscFree2(bases,ldims);CHKERRQ(ierr);
176047c6ae99SBarry Smith 
176147c6ae99SBarry Smith   dd->m  = m;  dd->n  = n;
176247c6ae99SBarry Smith   dd->xs = xs; dd->xe = xe; dd->ys = ys; dd->ye = ye; dd->zs = 0; dd->ze = 1;
176347c6ae99SBarry Smith   dd->Xs = Xs; dd->Xe = Xe; dd->Ys = Ys; dd->Ye = Ye; dd->Zs = 0; dd->Ze = 1;
176447c6ae99SBarry Smith 
176547c6ae99SBarry Smith   ierr = VecDestroy(local);CHKERRQ(ierr);
176647c6ae99SBarry Smith   ierr = VecDestroy(global);CHKERRQ(ierr);
176747c6ae99SBarry Smith 
176847c6ae99SBarry Smith   dd->gtol      = gtol;
176947c6ae99SBarry Smith   dd->ltog      = ltog;
177047c6ae99SBarry Smith   dd->idx       = idx;
177147c6ae99SBarry Smith   dd->Nl        = nn;
177247c6ae99SBarry Smith   dd->base      = base;
177347c6ae99SBarry Smith   da->ops->view = DAView_2d;
177447c6ae99SBarry Smith 
177547c6ae99SBarry Smith   /*
177647c6ae99SBarry Smith      Set the local to global ordering in the global vector, this allows use
177747c6ae99SBarry Smith      of VecSetValuesLocal().
177847c6ae99SBarry Smith   */
177947c6ae99SBarry Smith   ierr = ISLocalToGlobalMappingCreate(comm,nn,idx,PETSC_OWN_POINTER,&dd->ltogmap);CHKERRQ(ierr);
178047c6ae99SBarry Smith   ierr = ISLocalToGlobalMappingBlock(dd->ltogmap,dd->w,&dd->ltogmapb);CHKERRQ(ierr);
178147c6ae99SBarry Smith   ierr = PetscLogObjectParent(da,dd->ltogmap);CHKERRQ(ierr);
178247c6ae99SBarry Smith 
178347c6ae99SBarry Smith   dd->ltol = PETSC_NULL;
178447c6ae99SBarry Smith   dd->ao   = PETSC_NULL;
178547c6ae99SBarry Smith 
178647c6ae99SBarry Smith   PetscFunctionReturn(0);
178747c6ae99SBarry Smith }
178847c6ae99SBarry Smith EXTERN_C_END
178947c6ae99SBarry Smith 
179047c6ae99SBarry Smith #undef __FUNCT__
179147c6ae99SBarry Smith #define __FUNCT__ "DACreate2d"
179247c6ae99SBarry Smith /*@C
179347c6ae99SBarry Smith    DACreate2d -  Creates an object that will manage the communication of  two-dimensional
179447c6ae99SBarry Smith    regular array data that is distributed across some processors.
179547c6ae99SBarry Smith 
179647c6ae99SBarry Smith    Collective on MPI_Comm
179747c6ae99SBarry Smith 
179847c6ae99SBarry Smith    Input Parameters:
179947c6ae99SBarry Smith +  comm - MPI communicator
180047c6ae99SBarry Smith .  wrap - type of periodicity should the array have.
180147c6ae99SBarry Smith          Use one of DA_NONPERIODIC, DA_XPERIODIC, DA_YPERIODIC, or DA_XYPERIODIC.
180247c6ae99SBarry Smith .  stencil_type - stencil type.  Use either DA_STENCIL_BOX or DA_STENCIL_STAR.
180347c6ae99SBarry 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
180447c6ae99SBarry Smith             from the command line with -da_grid_x <M> -da_grid_y <N>)
180547c6ae99SBarry Smith .  m,n - corresponding number of processors in each dimension
180647c6ae99SBarry Smith          (or PETSC_DECIDE to have calculated)
180747c6ae99SBarry Smith .  dof - number of degrees of freedom per node
180847c6ae99SBarry Smith .  s - stencil width
180947c6ae99SBarry Smith -  lx, ly - arrays containing the number of nodes in each cell along
181047c6ae99SBarry Smith            the x and y coordinates, or PETSC_NULL. If non-null, these
181147c6ae99SBarry Smith            must be of length as m and n, and the corresponding
181247c6ae99SBarry Smith            m and n cannot be PETSC_DECIDE. The sum of the lx[] entries
181347c6ae99SBarry Smith            must be M, and the sum of the ly[] entries must be N.
181447c6ae99SBarry Smith 
181547c6ae99SBarry Smith    Output Parameter:
181647c6ae99SBarry Smith .  da - the resulting distributed array object
181747c6ae99SBarry Smith 
181847c6ae99SBarry Smith    Options Database Key:
181947c6ae99SBarry Smith +  -da_view - Calls DAView() at the conclusion of DACreate2d()
182047c6ae99SBarry Smith .  -da_grid_x <nx> - number of grid points in x direction, if M < 0
182147c6ae99SBarry Smith .  -da_grid_y <ny> - number of grid points in y direction, if N < 0
182247c6ae99SBarry Smith .  -da_processors_x <nx> - number of processors in x direction
182347c6ae99SBarry Smith .  -da_processors_y <ny> - number of processors in y direction
182447c6ae99SBarry Smith .  -da_refine_x - refinement ratio in x direction
182547c6ae99SBarry Smith -  -da_refine_y - refinement ratio in y direction
182647c6ae99SBarry Smith 
182747c6ae99SBarry Smith    Level: beginner
182847c6ae99SBarry Smith 
182947c6ae99SBarry Smith    Notes:
183047c6ae99SBarry Smith    The stencil type DA_STENCIL_STAR with width 1 corresponds to the
183147c6ae99SBarry Smith    standard 5-pt stencil, while DA_STENCIL_BOX with width 1 denotes
183247c6ae99SBarry Smith    the standard 9-pt stencil.
183347c6ae99SBarry Smith 
183447c6ae99SBarry Smith    The array data itself is NOT stored in the DA, it is stored in Vec objects;
183547c6ae99SBarry Smith    The appropriate vector objects can be obtained with calls to DACreateGlobalVector()
183647c6ae99SBarry Smith    and DACreateLocalVector() and calls to VecDuplicate() if more are needed.
183747c6ae99SBarry Smith 
183847c6ae99SBarry Smith .keywords: distributed array, create, two-dimensional
183947c6ae99SBarry Smith 
184047c6ae99SBarry Smith .seealso: DADestroy(), DAView(), DACreate1d(), DACreate3d(), DAGlobalToLocalBegin(), DAGetRefinementFactor(),
184147c6ae99SBarry Smith           DAGlobalToLocalEnd(), DALocalToGlobal(), DALocalToLocalBegin(), DALocalToLocalEnd(), DASetRefinementFactor(),
184247c6ae99SBarry Smith           DAGetInfo(), DACreateGlobalVector(), DACreateLocalVector(), DACreateNaturalVector(), DALoad(), DAView(), DAGetOwnershipRanges()
184347c6ae99SBarry Smith 
184447c6ae99SBarry Smith @*/
184547c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DACreate2d(MPI_Comm comm,DAPeriodicType wrap,DAStencilType stencil_type,
184647c6ae99SBarry Smith                           PetscInt M,PetscInt N,PetscInt m,PetscInt n,PetscInt dof,PetscInt s,const PetscInt lx[],const PetscInt ly[],DA *da)
184747c6ae99SBarry Smith {
184847c6ae99SBarry Smith   PetscErrorCode ierr;
184947c6ae99SBarry Smith 
185047c6ae99SBarry Smith   PetscFunctionBegin;
185147c6ae99SBarry Smith   ierr = DACreate(comm, da);CHKERRQ(ierr);
185247c6ae99SBarry Smith   ierr = DASetDim(*da, 2);CHKERRQ(ierr);
185347c6ae99SBarry Smith   ierr = DASetSizes(*da, M, N, 1);CHKERRQ(ierr);
185447c6ae99SBarry Smith   ierr = DASetNumProcs(*da, m, n, PETSC_DECIDE);CHKERRQ(ierr);
185547c6ae99SBarry Smith   ierr = DASetPeriodicity(*da, wrap);CHKERRQ(ierr);
185647c6ae99SBarry Smith   ierr = DASetDof(*da, dof);CHKERRQ(ierr);
185747c6ae99SBarry Smith   ierr = DASetStencilType(*da, stencil_type);CHKERRQ(ierr);
185847c6ae99SBarry Smith   ierr = DASetStencilWidth(*da, s);CHKERRQ(ierr);
185947c6ae99SBarry Smith   ierr = DASetOwnershipRanges(*da, lx, ly, PETSC_NULL);CHKERRQ(ierr);
186047c6ae99SBarry Smith   /* This violates the behavior for other classes, but right now users expect negative dimensions to be handled this way */
186147c6ae99SBarry Smith   ierr = DASetFromOptions(*da);CHKERRQ(ierr);
186247c6ae99SBarry Smith   ierr = DASetType(*da, DA2D);CHKERRQ(ierr);
186347c6ae99SBarry Smith   PetscFunctionReturn(0);
186447c6ae99SBarry Smith }
1865