xref: /petsc/src/dm/impls/da/da2.c (revision b45d2f2cb7e031d9c0de5873eca80614ca7b863b)
19a42bb27SBarry Smith 
2*b45d2f2cSJed Brown #include <petsc-private/daimpl.h>    /*I   "petscdmda.h"   I*/
347c6ae99SBarry Smith 
447c6ae99SBarry Smith #undef __FUNCT__
59a42bb27SBarry Smith #define __FUNCT__ "DMView_DA_2d"
69a42bb27SBarry Smith PetscErrorCode DMView_DA_2d(DM da,PetscViewer viewer)
747c6ae99SBarry Smith {
847c6ae99SBarry Smith   PetscErrorCode ierr;
947c6ae99SBarry Smith   PetscMPIInt    rank;
109a42bb27SBarry Smith   PetscBool      iascii,isdraw,isbinary;
1147c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
129a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
139a42bb27SBarry Smith   PetscBool      ismatlab;
149a42bb27SBarry Smith #endif
1547c6ae99SBarry Smith 
1647c6ae99SBarry Smith   PetscFunctionBegin;
1747c6ae99SBarry Smith   ierr = MPI_Comm_rank(((PetscObject)da)->comm,&rank);CHKERRQ(ierr);
1847c6ae99SBarry Smith 
1947c6ae99SBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
2047c6ae99SBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
219a42bb27SBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
229a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
239a42bb27SBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERMATLAB,&ismatlab);CHKERRQ(ierr);
249a42bb27SBarry Smith #endif
2547c6ae99SBarry Smith   if (iascii) {
2647c6ae99SBarry Smith     PetscViewerFormat format;
2747c6ae99SBarry Smith 
2847c6ae99SBarry Smith     ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr);
2947c6ae99SBarry Smith     if (format != PETSC_VIEWER_ASCII_VTK && format != PETSC_VIEWER_ASCII_VTK_CELL) {
30aa219208SBarry Smith       DMDALocalInfo info;
31aa219208SBarry Smith       ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
327b23a99aSBarry Smith       ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
3347c6ae99SBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Processor [%d] M %D N %D m %D n %D w %D s %D\n",rank,dd->M,dd->N,dd->m,dd->n,dd->w,dd->s);CHKERRQ(ierr);
3447c6ae99SBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"X range of indices: %D %D, Y range of indices: %D %D\n",info.xs,info.xs+info.xm,info.ys,info.ys+info.ym);CHKERRQ(ierr);
3547c6ae99SBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
367b23a99aSBarry Smith       ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr);
373da9ae13SJed Brown     } else {
383da9ae13SJed Brown       ierr = DMView_DA_VTK(da,viewer);CHKERRQ(ierr);
3947c6ae99SBarry Smith     }
4047c6ae99SBarry Smith   } else if (isdraw) {
4147c6ae99SBarry Smith     PetscDraw  draw;
4247c6ae99SBarry Smith     double     ymin = -1*dd->s-1,ymax = dd->N+dd->s;
4347c6ae99SBarry Smith     double     xmin = -1*dd->s-1,xmax = dd->M+dd->s;
4447c6ae99SBarry Smith     double     x,y;
4547c6ae99SBarry Smith     PetscInt   base,*idx;
4647c6ae99SBarry Smith     char       node[10];
4747c6ae99SBarry Smith     PetscBool  isnull;
4847c6ae99SBarry Smith 
4947c6ae99SBarry Smith     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
5047c6ae99SBarry Smith     ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
5147c6ae99SBarry Smith     if (!dd->coordinates) {
5247c6ae99SBarry Smith       ierr = PetscDrawSetCoordinates(draw,xmin,ymin,xmax,ymax);CHKERRQ(ierr);
5347c6ae99SBarry Smith     }
5447c6ae99SBarry Smith     ierr = PetscDrawSynchronizedClear(draw);CHKERRQ(ierr);
5547c6ae99SBarry Smith 
5647c6ae99SBarry Smith     /* first processor draw all node lines */
5747c6ae99SBarry Smith     if (!rank) {
5847c6ae99SBarry Smith       ymin = 0.0; ymax = dd->N - 1;
5947c6ae99SBarry Smith       for (xmin=0; xmin<dd->M; xmin++) {
6047c6ae99SBarry Smith         ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_BLACK);CHKERRQ(ierr);
6147c6ae99SBarry Smith       }
6247c6ae99SBarry Smith       xmin = 0.0; xmax = dd->M - 1;
6347c6ae99SBarry Smith       for (ymin=0; ymin<dd->N; ymin++) {
6447c6ae99SBarry Smith         ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_BLACK);CHKERRQ(ierr);
6547c6ae99SBarry Smith       }
6647c6ae99SBarry Smith     }
6747c6ae99SBarry Smith     ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr);
6847c6ae99SBarry Smith     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
6947c6ae99SBarry Smith 
7047c6ae99SBarry Smith     /* draw my box */
7147c6ae99SBarry Smith     ymin = dd->ys; ymax = dd->ye - 1; xmin = dd->xs/dd->w;
7247c6ae99SBarry Smith     xmax =(dd->xe-1)/dd->w;
7347c6ae99SBarry Smith     ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_RED);CHKERRQ(ierr);
7447c6ae99SBarry Smith     ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_RED);CHKERRQ(ierr);
7547c6ae99SBarry Smith     ierr = PetscDrawLine(draw,xmin,ymax,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr);
7647c6ae99SBarry Smith     ierr = PetscDrawLine(draw,xmax,ymin,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr);
7747c6ae99SBarry Smith 
7847c6ae99SBarry Smith     /* put in numbers */
7947c6ae99SBarry Smith     base = (dd->base)/dd->w;
8047c6ae99SBarry Smith     for (y=ymin; y<=ymax; y++) {
8147c6ae99SBarry Smith       for (x=xmin; x<=xmax; x++) {
8247c6ae99SBarry Smith         sprintf(node,"%d",(int)base++);
8347c6ae99SBarry Smith         ierr = PetscDrawString(draw,x,y,PETSC_DRAW_BLACK,node);CHKERRQ(ierr);
8447c6ae99SBarry Smith       }
8547c6ae99SBarry Smith     }
8647c6ae99SBarry Smith 
8747c6ae99SBarry Smith     ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr);
8847c6ae99SBarry Smith     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
8947c6ae99SBarry Smith     /* overlay ghost numbers, useful for error checking */
9047c6ae99SBarry Smith     /* put in numbers */
9147c6ae99SBarry Smith 
9247c6ae99SBarry Smith     base = 0; idx = dd->idx;
9347c6ae99SBarry Smith     ymin = dd->Ys; ymax = dd->Ye; xmin = dd->Xs; xmax = dd->Xe;
9447c6ae99SBarry Smith     for (y=ymin; y<ymax; y++) {
9547c6ae99SBarry Smith       for (x=xmin; x<xmax; x++) {
9647c6ae99SBarry Smith         if ((base % dd->w) == 0) {
9747c6ae99SBarry Smith           sprintf(node,"%d",(int)(idx[base]/dd->w));
9847c6ae99SBarry Smith           ierr = PetscDrawString(draw,x/dd->w,y,PETSC_DRAW_BLUE,node);CHKERRQ(ierr);
9947c6ae99SBarry Smith         }
10047c6ae99SBarry Smith         base++;
10147c6ae99SBarry Smith       }
10247c6ae99SBarry Smith     }
10347c6ae99SBarry Smith     ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr);
10447c6ae99SBarry Smith     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
1059a42bb27SBarry Smith   } else if (isbinary){
1069a42bb27SBarry Smith     ierr = DMView_DA_Binary(da,viewer);CHKERRQ(ierr);
1079a42bb27SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
1089a42bb27SBarry Smith   } else if (ismatlab) {
1099a42bb27SBarry Smith     ierr = DMView_DA_Matlab(da,viewer);CHKERRQ(ierr);
1109a42bb27SBarry Smith #endif
111aa219208SBarry Smith   } else SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for DMDA 1d",((PetscObject)viewer)->type_name);
11247c6ae99SBarry Smith   PetscFunctionReturn(0);
11347c6ae99SBarry Smith }
11447c6ae99SBarry Smith 
11547c6ae99SBarry Smith /*
11647c6ae99SBarry Smith       M is number of grid points
11747c6ae99SBarry Smith       m is number of processors
11847c6ae99SBarry Smith 
11947c6ae99SBarry Smith */
12047c6ae99SBarry Smith #undef __FUNCT__
121aa219208SBarry Smith #define __FUNCT__ "DMDASplitComm2d"
1227087cfbeSBarry Smith PetscErrorCode  DMDASplitComm2d(MPI_Comm comm,PetscInt M,PetscInt N,PetscInt sw,MPI_Comm *outcomm)
12347c6ae99SBarry Smith {
12447c6ae99SBarry Smith   PetscErrorCode ierr;
12547c6ae99SBarry Smith   PetscInt       m,n = 0,x = 0,y = 0;
12647c6ae99SBarry Smith   PetscMPIInt    size,csize,rank;
12747c6ae99SBarry Smith 
12847c6ae99SBarry Smith   PetscFunctionBegin;
12947c6ae99SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
13047c6ae99SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
13147c6ae99SBarry Smith 
13247c6ae99SBarry Smith   csize = 4*size;
13347c6ae99SBarry Smith   do {
13447c6ae99SBarry Smith     if (csize % 4) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Cannot split communicator of size %d tried %d %D %D",size,csize,x,y);
13547c6ae99SBarry Smith     csize   = csize/4;
13647c6ae99SBarry Smith 
13747c6ae99SBarry Smith     m = (PetscInt)(0.5 + sqrt(((double)M)*((double)csize)/((double)N)));
13847c6ae99SBarry Smith     if (!m) m = 1;
13947c6ae99SBarry Smith     while (m > 0) {
14047c6ae99SBarry Smith       n = csize/m;
14147c6ae99SBarry Smith       if (m*n == csize) break;
14247c6ae99SBarry Smith       m--;
14347c6ae99SBarry Smith     }
14447c6ae99SBarry Smith     if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;}
14547c6ae99SBarry Smith 
14647c6ae99SBarry Smith     x = M/m + ((M % m) > ((csize-1) % m));
14747c6ae99SBarry Smith     y = (N + (csize-1)/m)/n;
14847c6ae99SBarry Smith   } while ((x < 4 || y < 4) && csize > 1);
14947c6ae99SBarry Smith   if (size != csize) {
15047c6ae99SBarry Smith     MPI_Group    entire_group,sub_group;
15147c6ae99SBarry Smith     PetscMPIInt  i,*groupies;
15247c6ae99SBarry Smith 
15347c6ae99SBarry Smith     ierr = MPI_Comm_group(comm,&entire_group);CHKERRQ(ierr);
15447c6ae99SBarry Smith     ierr = PetscMalloc(csize*sizeof(PetscInt),&groupies);CHKERRQ(ierr);
15547c6ae99SBarry Smith     for (i=0; i<csize; i++) {
15647c6ae99SBarry Smith       groupies[i] = (rank/csize)*csize + i;
15747c6ae99SBarry Smith     }
15847c6ae99SBarry Smith     ierr = MPI_Group_incl(entire_group,csize,groupies,&sub_group);CHKERRQ(ierr);
15947c6ae99SBarry Smith     ierr = PetscFree(groupies);CHKERRQ(ierr);
16047c6ae99SBarry Smith     ierr = MPI_Comm_create(comm,sub_group,outcomm);CHKERRQ(ierr);
16147c6ae99SBarry Smith     ierr = MPI_Group_free(&entire_group);CHKERRQ(ierr);
16247c6ae99SBarry Smith     ierr = MPI_Group_free(&sub_group);CHKERRQ(ierr);
163aa219208SBarry Smith     ierr = PetscInfo1(0,"DMDASplitComm2d:Creating redundant coarse problems of size %d\n",csize);CHKERRQ(ierr);
16447c6ae99SBarry Smith   } else {
16547c6ae99SBarry Smith     *outcomm = comm;
16647c6ae99SBarry Smith   }
16747c6ae99SBarry Smith   PetscFunctionReturn(0);
16847c6ae99SBarry Smith }
16947c6ae99SBarry Smith 
17047c6ae99SBarry Smith #undef __FUNCT__
171644e2e5bSBarry Smith #define __FUNCT__ "DMDAFunction"
172644e2e5bSBarry Smith static PetscErrorCode DMDAFunction(DM dm,Vec x,Vec F)
173644e2e5bSBarry Smith {
174644e2e5bSBarry Smith   PetscErrorCode ierr;
175644e2e5bSBarry Smith   Vec            localX;
176644e2e5bSBarry Smith 
177644e2e5bSBarry Smith   PetscFunctionBegin;
178644e2e5bSBarry Smith   ierr = DMGetLocalVector(dm,&localX);CHKERRQ(ierr);
179644e2e5bSBarry Smith   ierr = DMGlobalToLocalBegin(dm,x,INSERT_VALUES,localX);CHKERRQ(ierr);
180644e2e5bSBarry Smith   ierr = DMGlobalToLocalEnd(dm,x,INSERT_VALUES,localX);CHKERRQ(ierr);
181837ad101SBarry Smith   ierr = DMDAComputeFunction1(dm,localX,F,dm->ctx);CHKERRQ(ierr);
182644e2e5bSBarry Smith   ierr = DMRestoreLocalVector(dm,&localX);CHKERRQ(ierr);
183644e2e5bSBarry Smith   PetscFunctionReturn(0);
184644e2e5bSBarry Smith }
185644e2e5bSBarry Smith 
186644e2e5bSBarry Smith #undef __FUNCT__
187aa219208SBarry Smith #define __FUNCT__ "DMDASetLocalFunction"
18847c6ae99SBarry Smith /*@C
189aa219208SBarry Smith        DMDASetLocalFunction - Caches in a DM a local function.
19047c6ae99SBarry Smith 
191aa219208SBarry Smith    Logically Collective on DMDA
19247c6ae99SBarry Smith 
19347c6ae99SBarry Smith    Input Parameter:
19447c6ae99SBarry Smith +  da - initial distributed array
19547c6ae99SBarry Smith -  lf - the local function
19647c6ae99SBarry Smith 
19747c6ae99SBarry Smith    Level: intermediate
19847c6ae99SBarry Smith 
19931e16157SBarry Smith    Notes:
20031e16157SBarry Smith       If you used SNESSetFunction(snes,r,SNESDMDAComputeFunction,ctx); then the context passed to your function is the ctx set here.
20131e16157SBarry Smith 
20231e16157SBarry Smith       If you use SNESSetDM() and did not use SNESSetFunction() then the context passed to your function is the context set with DMSetApplicationContext()
20331e16157SBarry Smith 
20431e16157SBarry Smith    Developer Notes: It is possibly confusing which context is passed to the user function, it would be nice to unify them somehow.
20547c6ae99SBarry Smith 
20647c6ae99SBarry Smith .keywords:  distributed array, refine
20747c6ae99SBarry Smith 
208aa219208SBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalFunctioni()
20947c6ae99SBarry Smith @*/
2107087cfbeSBarry Smith PetscErrorCode  DMDASetLocalFunction(DM da,DMDALocalFunction1 lf)
21147c6ae99SBarry Smith {
212644e2e5bSBarry Smith   PetscErrorCode ierr;
21347c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
214644e2e5bSBarry Smith 
21547c6ae99SBarry Smith   PetscFunctionBegin;
21647c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
217644e2e5bSBarry Smith   ierr = DMSetFunction(da,DMDAFunction);CHKERRQ(ierr);
21847c6ae99SBarry Smith   dd->lf       = lf;
21947c6ae99SBarry Smith   PetscFunctionReturn(0);
22047c6ae99SBarry Smith }
22147c6ae99SBarry Smith 
22247c6ae99SBarry Smith #undef __FUNCT__
223aa219208SBarry Smith #define __FUNCT__ "DMDASetLocalFunctioni"
22447c6ae99SBarry Smith /*@C
225aa219208SBarry Smith        DMDASetLocalFunctioni - Caches in a DM a local function that evaluates a single component
22647c6ae99SBarry Smith 
227aa219208SBarry Smith    Logically Collective on DMDA
22847c6ae99SBarry Smith 
22947c6ae99SBarry Smith    Input Parameter:
23047c6ae99SBarry Smith +  da - initial distributed array
23147c6ae99SBarry Smith -  lfi - the local function
23247c6ae99SBarry Smith 
23347c6ae99SBarry Smith    Level: intermediate
23447c6ae99SBarry Smith 
23547c6ae99SBarry Smith .keywords:  distributed array, refine
23647c6ae99SBarry Smith 
237aa219208SBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalFunction()
23847c6ae99SBarry Smith @*/
2397087cfbeSBarry Smith PetscErrorCode  DMDASetLocalFunctioni(DM da,PetscErrorCode (*lfi)(DMDALocalInfo*,MatStencil*,void*,PetscScalar*,void*))
24047c6ae99SBarry Smith {
24147c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
24247c6ae99SBarry Smith   PetscFunctionBegin;
24347c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
24447c6ae99SBarry Smith   dd->lfi = lfi;
24547c6ae99SBarry Smith   PetscFunctionReturn(0);
24647c6ae99SBarry Smith }
24747c6ae99SBarry Smith 
24847c6ae99SBarry Smith #undef __FUNCT__
249aa219208SBarry Smith #define __FUNCT__ "DMDASetLocalFunctionib"
25047c6ae99SBarry Smith /*@C
251aa219208SBarry Smith        DMDASetLocalFunctionib - Caches in a DM a block local function that evaluates a single component
25247c6ae99SBarry Smith 
253aa219208SBarry Smith    Logically Collective on DMDA
25447c6ae99SBarry Smith 
25547c6ae99SBarry Smith    Input Parameter:
25647c6ae99SBarry Smith +  da - initial distributed array
25747c6ae99SBarry Smith -  lfi - the local function
25847c6ae99SBarry Smith 
25947c6ae99SBarry Smith    Level: intermediate
26047c6ae99SBarry Smith 
26147c6ae99SBarry Smith .keywords:  distributed array, refine
26247c6ae99SBarry Smith 
263aa219208SBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalFunction()
26447c6ae99SBarry Smith @*/
2657087cfbeSBarry Smith PetscErrorCode  DMDASetLocalFunctionib(DM da,PetscErrorCode (*lfi)(DMDALocalInfo*,MatStencil*,void*,PetscScalar*,void*))
26647c6ae99SBarry Smith {
26747c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
26847c6ae99SBarry Smith   PetscFunctionBegin;
26947c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
27047c6ae99SBarry Smith   dd->lfib = lfi;
27147c6ae99SBarry Smith   PetscFunctionReturn(0);
27247c6ae99SBarry Smith }
27347c6ae99SBarry Smith 
27447c6ae99SBarry Smith #undef __FUNCT__
275aa219208SBarry Smith #define __FUNCT__ "DMDASetLocalAdicFunction_Private"
276aa219208SBarry Smith PetscErrorCode DMDASetLocalAdicFunction_Private(DM da,DMDALocalFunction1 ad_lf)
27747c6ae99SBarry Smith {
27847c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
27947c6ae99SBarry Smith   PetscFunctionBegin;
28047c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
28147c6ae99SBarry Smith   dd->adic_lf = ad_lf;
28247c6ae99SBarry Smith   PetscFunctionReturn(0);
28347c6ae99SBarry Smith }
28447c6ae99SBarry Smith 
28547c6ae99SBarry Smith /*MC
286aa219208SBarry Smith        DMDASetLocalAdicFunctioni - Caches in a DM a local functioni computed by ADIC/ADIFOR
28747c6ae99SBarry Smith 
28847c6ae99SBarry Smith    Synopsis:
289aa219208SBarry Smith    PetscErrorCode DMDASetLocalAdicFunctioni(DM da,PetscInt (ad_lf*)(DMDALocalInfo*,MatStencil*,void*,void*,void*)
29047c6ae99SBarry Smith 
291aa219208SBarry Smith    Logically Collective on DMDA
29247c6ae99SBarry Smith 
29347c6ae99SBarry Smith    Input Parameter:
29447c6ae99SBarry Smith +  da - initial distributed array
29547c6ae99SBarry Smith -  ad_lfi - the local function as computed by ADIC/ADIFOR
29647c6ae99SBarry Smith 
29747c6ae99SBarry Smith    Level: intermediate
29847c6ae99SBarry Smith 
29947c6ae99SBarry Smith .keywords:  distributed array, refine
30047c6ae99SBarry Smith 
301aa219208SBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalFunction(),
302aa219208SBarry Smith           DMDASetLocalJacobian(), DMDASetLocalFunctioni()
30347c6ae99SBarry Smith M*/
30447c6ae99SBarry Smith 
30547c6ae99SBarry Smith #undef __FUNCT__
306aa219208SBarry Smith #define __FUNCT__ "DMDASetLocalAdicFunctioni_Private"
307aa219208SBarry Smith PetscErrorCode DMDASetLocalAdicFunctioni_Private(DM da,PetscErrorCode (*ad_lfi)(DMDALocalInfo*,MatStencil*,void*,void*,void*))
30847c6ae99SBarry Smith {
30947c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
31047c6ae99SBarry Smith   PetscFunctionBegin;
31147c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
31247c6ae99SBarry Smith   dd->adic_lfi = ad_lfi;
31347c6ae99SBarry Smith   PetscFunctionReturn(0);
31447c6ae99SBarry Smith }
31547c6ae99SBarry Smith 
31647c6ae99SBarry Smith /*MC
317aa219208SBarry Smith        DMDASetLocalAdicMFFunctioni - Caches in a DM a local functioni computed by ADIC/ADIFOR
31847c6ae99SBarry Smith 
31947c6ae99SBarry Smith    Synopsis:
320aa219208SBarry Smith    PetscErrorCode  DMDASetLocalAdicFunctioni(DM da,int (ad_lf*)(DMDALocalInfo*,MatStencil*,void*,void*,void*)
32147c6ae99SBarry Smith 
322aa219208SBarry Smith    Logically Collective on DMDA
32347c6ae99SBarry Smith 
32447c6ae99SBarry Smith    Input Parameter:
32547c6ae99SBarry Smith +  da - initial distributed array
32647c6ae99SBarry Smith -  admf_lfi - the local matrix-free function as computed by ADIC/ADIFOR
32747c6ae99SBarry Smith 
32847c6ae99SBarry Smith    Level: intermediate
32947c6ae99SBarry Smith 
33047c6ae99SBarry Smith .keywords:  distributed array, refine
33147c6ae99SBarry Smith 
332aa219208SBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalFunction(),
333aa219208SBarry Smith           DMDASetLocalJacobian(), DMDASetLocalFunctioni()
33447c6ae99SBarry Smith M*/
33547c6ae99SBarry Smith 
33647c6ae99SBarry Smith #undef __FUNCT__
337aa219208SBarry Smith #define __FUNCT__ "DMDASetLocalAdicMFFunctioni_Private"
338aa219208SBarry Smith PetscErrorCode DMDASetLocalAdicMFFunctioni_Private(DM da,PetscErrorCode (*admf_lfi)(DMDALocalInfo*,MatStencil*,void*,void*,void*))
33947c6ae99SBarry Smith {
34047c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
34147c6ae99SBarry Smith   PetscFunctionBegin;
34247c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
34347c6ae99SBarry Smith   dd->adicmf_lfi = admf_lfi;
34447c6ae99SBarry Smith   PetscFunctionReturn(0);
34547c6ae99SBarry Smith }
34647c6ae99SBarry Smith 
34747c6ae99SBarry Smith /*MC
348aa219208SBarry Smith        DMDASetLocalAdicFunctionib - Caches in a DM a block local functioni computed by ADIC/ADIFOR
34947c6ae99SBarry Smith 
35047c6ae99SBarry Smith    Synopsis:
351aa219208SBarry Smith    PetscErrorCode DMDASetLocalAdicFunctionib(DM da,PetscInt (ad_lf*)(DMDALocalInfo*,MatStencil*,void*,void*,void*)
35247c6ae99SBarry Smith 
353aa219208SBarry Smith    Logically Collective on DMDA
35447c6ae99SBarry Smith 
35547c6ae99SBarry Smith    Input Parameter:
35647c6ae99SBarry Smith +  da - initial distributed array
35747c6ae99SBarry Smith -  ad_lfi - the local function as computed by ADIC/ADIFOR
35847c6ae99SBarry Smith 
35947c6ae99SBarry Smith    Level: intermediate
36047c6ae99SBarry Smith 
36147c6ae99SBarry Smith .keywords:  distributed array, refine
36247c6ae99SBarry Smith 
363aa219208SBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalFunction(),
364aa219208SBarry Smith           DMDASetLocalJacobian(), DMDASetLocalFunctionib()
36547c6ae99SBarry Smith M*/
36647c6ae99SBarry Smith 
36747c6ae99SBarry Smith #undef __FUNCT__
368aa219208SBarry Smith #define __FUNCT__ "DMDASetLocalAdicFunctionib_Private"
369aa219208SBarry Smith PetscErrorCode DMDASetLocalAdicFunctionib_Private(DM da,PetscErrorCode (*ad_lfi)(DMDALocalInfo*,MatStencil*,void*,void*,void*))
37047c6ae99SBarry Smith {
37147c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
37247c6ae99SBarry Smith   PetscFunctionBegin;
37347c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
37447c6ae99SBarry Smith   dd->adic_lfib = ad_lfi;
37547c6ae99SBarry Smith   PetscFunctionReturn(0);
37647c6ae99SBarry Smith }
37747c6ae99SBarry Smith 
37847c6ae99SBarry Smith /*MC
379aa219208SBarry Smith        DMDASetLocalAdicMFFunctionib - Caches in a DM a block local functioni computed by ADIC/ADIFOR
38047c6ae99SBarry Smith 
38147c6ae99SBarry Smith    Synopsis:
382aa219208SBarry Smith    PetscErrorCode  DMDASetLocalAdicFunctionib(DM da,int (ad_lf*)(DMDALocalInfo*,MatStencil*,void*,void*,void*)
38347c6ae99SBarry Smith 
384aa219208SBarry Smith    Logically Collective on DMDA
38547c6ae99SBarry Smith 
38647c6ae99SBarry Smith    Input Parameter:
38747c6ae99SBarry Smith +  da - initial distributed array
38847c6ae99SBarry Smith -  admf_lfi - the local matrix-free function as computed by ADIC/ADIFOR
38947c6ae99SBarry Smith 
39047c6ae99SBarry Smith    Level: intermediate
39147c6ae99SBarry Smith 
39247c6ae99SBarry Smith .keywords:  distributed array, refine
39347c6ae99SBarry Smith 
394aa219208SBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalFunction(),
395aa219208SBarry Smith           DMDASetLocalJacobian(), DMDASetLocalFunctionib()
39647c6ae99SBarry Smith M*/
39747c6ae99SBarry Smith 
39847c6ae99SBarry Smith #undef __FUNCT__
399aa219208SBarry Smith #define __FUNCT__ "DMDASetLocalAdicMFFunctionib_Private"
400aa219208SBarry Smith PetscErrorCode DMDASetLocalAdicMFFunctionib_Private(DM da,PetscErrorCode (*admf_lfi)(DMDALocalInfo*,MatStencil*,void*,void*,void*))
40147c6ae99SBarry Smith {
40247c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
40347c6ae99SBarry Smith   PetscFunctionBegin;
40447c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
40547c6ae99SBarry Smith   dd->adicmf_lfib = admf_lfi;
40647c6ae99SBarry Smith   PetscFunctionReturn(0);
40747c6ae99SBarry Smith }
40847c6ae99SBarry Smith 
40947c6ae99SBarry Smith /*MC
410aa219208SBarry Smith        DMDASetLocalAdicMFFunction - Caches in a DM a local function computed by ADIC/ADIFOR
41147c6ae99SBarry Smith 
41247c6ae99SBarry Smith    Synopsis:
413aa219208SBarry Smith    PetscErrorCode DMDASetLocalAdicMFFunction(DM da,DMDALocalFunction1 ad_lf)
41447c6ae99SBarry Smith 
415aa219208SBarry Smith    Logically Collective on DMDA
41647c6ae99SBarry Smith 
41747c6ae99SBarry Smith    Input Parameter:
41847c6ae99SBarry Smith +  da - initial distributed array
41947c6ae99SBarry Smith -  ad_lf - the local function as computed by ADIC/ADIFOR
42047c6ae99SBarry Smith 
42147c6ae99SBarry Smith    Level: intermediate
42247c6ae99SBarry Smith 
42347c6ae99SBarry Smith .keywords:  distributed array, refine
42447c6ae99SBarry Smith 
425aa219208SBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalFunction(),
426aa219208SBarry Smith           DMDASetLocalJacobian()
42747c6ae99SBarry Smith M*/
42847c6ae99SBarry Smith 
42947c6ae99SBarry Smith #undef __FUNCT__
430aa219208SBarry Smith #define __FUNCT__ "DMDASetLocalAdicMFFunction_Private"
431aa219208SBarry Smith PetscErrorCode DMDASetLocalAdicMFFunction_Private(DM da,DMDALocalFunction1 ad_lf)
43247c6ae99SBarry Smith {
43347c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
43447c6ae99SBarry Smith   PetscFunctionBegin;
43547c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
43647c6ae99SBarry Smith   dd->adicmf_lf = ad_lf;
43747c6ae99SBarry Smith   PetscFunctionReturn(0);
43847c6ae99SBarry Smith }
43947c6ae99SBarry Smith 
440644e2e5bSBarry Smith #undef __FUNCT__
4412533e041SBarry Smith #define __FUNCT__ "DMDAJacobianDefaultLocal"
4422533e041SBarry Smith PetscErrorCode DMDAJacobianLocal(DM dm,Vec x,Mat A,Mat B, MatStructure *str)
4432533e041SBarry Smith {
4442533e041SBarry Smith   PetscErrorCode ierr;
4452533e041SBarry Smith   Vec            localX;
4462533e041SBarry Smith 
4472533e041SBarry Smith   PetscFunctionBegin;
4482533e041SBarry Smith   ierr = DMGetLocalVector(dm,&localX);CHKERRQ(ierr);
4492533e041SBarry Smith   ierr = DMGlobalToLocalBegin(dm,x,INSERT_VALUES,localX);CHKERRQ(ierr);
4502533e041SBarry Smith   ierr = DMGlobalToLocalEnd(dm,x,INSERT_VALUES,localX);CHKERRQ(ierr);
4512533e041SBarry Smith   ierr = MatFDColoringApply(B,dm->fd,localX,str,dm);CHKERRQ(ierr);
4522533e041SBarry Smith   ierr = DMRestoreLocalVector(dm,&localX);CHKERRQ(ierr);
4532533e041SBarry Smith   /* Assemble true Jacobian; if it is different */
4542533e041SBarry Smith   if (A != B) {
4552533e041SBarry Smith     ierr  = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
4562533e041SBarry Smith     ierr  = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
4572533e041SBarry Smith   }
4582533e041SBarry Smith   ierr  = MatSetOption(B,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
4592533e041SBarry Smith   *str = SAME_NONZERO_PATTERN;
4602533e041SBarry Smith   PetscFunctionReturn(0);
4612533e041SBarry Smith }
4622533e041SBarry Smith 
4632533e041SBarry Smith 
4642533e041SBarry Smith #undef __FUNCT__
465644e2e5bSBarry Smith #define __FUNCT__ "DMDAJacobian"
466644e2e5bSBarry Smith static PetscErrorCode DMDAJacobian(DM dm,Vec x,Mat A,Mat B, MatStructure *str)
467644e2e5bSBarry Smith {
468644e2e5bSBarry Smith   PetscErrorCode ierr;
469644e2e5bSBarry Smith   Vec            localX;
470644e2e5bSBarry Smith 
471644e2e5bSBarry Smith   PetscFunctionBegin;
472644e2e5bSBarry Smith   ierr = DMGetLocalVector(dm,&localX);CHKERRQ(ierr);
473644e2e5bSBarry Smith   ierr = DMGlobalToLocalBegin(dm,x,INSERT_VALUES,localX);CHKERRQ(ierr);
474644e2e5bSBarry Smith   ierr = DMGlobalToLocalEnd(dm,x,INSERT_VALUES,localX);CHKERRQ(ierr);
475644e2e5bSBarry Smith   ierr = DMDAComputeJacobian1(dm,localX,B,dm->ctx);CHKERRQ(ierr);
476644e2e5bSBarry Smith   ierr = DMRestoreLocalVector(dm,&localX);CHKERRQ(ierr);
477644e2e5bSBarry Smith   /* Assemble true Jacobian; if it is different */
478644e2e5bSBarry Smith   if (A != B) {
479644e2e5bSBarry Smith     ierr  = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
480644e2e5bSBarry Smith     ierr  = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
481644e2e5bSBarry Smith   }
482644e2e5bSBarry Smith   ierr  = MatSetOption(B,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
483644e2e5bSBarry Smith   *str = SAME_NONZERO_PATTERN;
484644e2e5bSBarry Smith   PetscFunctionReturn(0);
485644e2e5bSBarry Smith }
486644e2e5bSBarry Smith 
48747c6ae99SBarry Smith /*@C
488aa219208SBarry Smith        DMDASetLocalJacobian - Caches in a DM a local Jacobian computation function
48947c6ae99SBarry Smith 
490aa219208SBarry Smith    Logically Collective on DMDA
49147c6ae99SBarry Smith 
49247c6ae99SBarry Smith 
49347c6ae99SBarry Smith    Input Parameter:
49447c6ae99SBarry Smith +  da - initial distributed array
49547c6ae99SBarry Smith -  lj - the local Jacobian
49647c6ae99SBarry Smith 
49747c6ae99SBarry Smith    Level: intermediate
49847c6ae99SBarry Smith 
499837ad101SBarry Smith    Notes: The routine SNESDMDAComputeFunction() uses this the cached function to evaluate the user provided function.
50047c6ae99SBarry Smith 
50147c6ae99SBarry Smith .keywords:  distributed array, refine
50247c6ae99SBarry Smith 
503aa219208SBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalFunction()
50447c6ae99SBarry Smith @*/
50547c6ae99SBarry Smith #undef __FUNCT__
506aa219208SBarry Smith #define __FUNCT__ "DMDASetLocalJacobian"
5077087cfbeSBarry Smith PetscErrorCode  DMDASetLocalJacobian(DM da,DMDALocalFunction1 lj)
50847c6ae99SBarry Smith {
509644e2e5bSBarry Smith   PetscErrorCode ierr;
51047c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
511644e2e5bSBarry Smith 
51247c6ae99SBarry Smith   PetscFunctionBegin;
51347c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
514644e2e5bSBarry Smith   ierr = DMSetJacobian(da,DMDAJacobian);CHKERRQ(ierr);
51547c6ae99SBarry Smith   dd->lj    = lj;
51647c6ae99SBarry Smith   PetscFunctionReturn(0);
51747c6ae99SBarry Smith }
51847c6ae99SBarry Smith 
51947c6ae99SBarry Smith #undef __FUNCT__
520aa219208SBarry Smith #define __FUNCT__ "DMDAGetLocalFunction"
52147c6ae99SBarry Smith /*@C
522aa219208SBarry Smith        DMDAGetLocalFunction - Gets from a DM a local function and its ADIC/ADIFOR Jacobian
52347c6ae99SBarry Smith 
52447c6ae99SBarry Smith    Note Collective
52547c6ae99SBarry Smith 
52647c6ae99SBarry Smith    Input Parameter:
52747c6ae99SBarry Smith .  da - initial distributed array
52847c6ae99SBarry Smith 
52947c6ae99SBarry Smith    Output Parameter:
53047c6ae99SBarry Smith .  lf - the local function
53147c6ae99SBarry Smith 
53247c6ae99SBarry Smith    Level: intermediate
53347c6ae99SBarry Smith 
53447c6ae99SBarry Smith .keywords:  distributed array, refine
53547c6ae99SBarry Smith 
536aa219208SBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalJacobian(), DMDASetLocalFunction()
53747c6ae99SBarry Smith @*/
5387087cfbeSBarry Smith PetscErrorCode  DMDAGetLocalFunction(DM da,DMDALocalFunction1 *lf)
53947c6ae99SBarry Smith {
54047c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
54147c6ae99SBarry Smith   PetscFunctionBegin;
54247c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
54347c6ae99SBarry Smith   if (lf) *lf = dd->lf;
54447c6ae99SBarry Smith   PetscFunctionReturn(0);
54547c6ae99SBarry Smith }
54647c6ae99SBarry Smith 
54747c6ae99SBarry Smith #undef __FUNCT__
548aa219208SBarry Smith #define __FUNCT__ "DMDAGetLocalJacobian"
54947c6ae99SBarry Smith /*@C
550aa219208SBarry Smith        DMDAGetLocalJacobian - Gets from a DM a local jacobian
55147c6ae99SBarry Smith 
55247c6ae99SBarry Smith    Not Collective
55347c6ae99SBarry Smith 
55447c6ae99SBarry Smith    Input Parameter:
55547c6ae99SBarry Smith .  da - initial distributed array
55647c6ae99SBarry Smith 
55747c6ae99SBarry Smith    Output Parameter:
55847c6ae99SBarry Smith .  lj - the local jacobian
55947c6ae99SBarry Smith 
56047c6ae99SBarry Smith    Level: intermediate
56147c6ae99SBarry Smith 
56247c6ae99SBarry Smith .keywords:  distributed array, refine
56347c6ae99SBarry Smith 
564aa219208SBarry Smith .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalJacobian()
56547c6ae99SBarry Smith @*/
5667087cfbeSBarry Smith PetscErrorCode  DMDAGetLocalJacobian(DM da,DMDALocalFunction1 *lj)
56747c6ae99SBarry Smith {
56847c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
56947c6ae99SBarry Smith   PetscFunctionBegin;
57047c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
57147c6ae99SBarry Smith   if (lj) *lj = dd->lj;
57247c6ae99SBarry Smith   PetscFunctionReturn(0);
57347c6ae99SBarry Smith }
57447c6ae99SBarry Smith 
57547c6ae99SBarry Smith #undef __FUNCT__
576837ad101SBarry Smith #define __FUNCT__ "DMDAComputeFunction"
57747c6ae99SBarry Smith /*@
578837ad101SBarry Smith     DMDAComputeFunction - Evaluates a user provided function on each processor that
579aa219208SBarry Smith         share a DMDA
58047c6ae99SBarry Smith 
58147c6ae99SBarry Smith    Input Parameters:
5829a42bb27SBarry Smith +    da - the DM that defines the grid
58347c6ae99SBarry Smith .    vu - input vector
58447c6ae99SBarry Smith .    vfu - output vector
58547c6ae99SBarry Smith -    w - any user data
58647c6ae99SBarry Smith 
58747c6ae99SBarry Smith     Notes: Does NOT do ghost updates on vu upon entry
58847c6ae99SBarry Smith 
589837ad101SBarry Smith            This should eventually replace DMDAComputeFunction1
59047c6ae99SBarry Smith 
59147c6ae99SBarry Smith     Level: advanced
59247c6ae99SBarry Smith 
593aa219208SBarry Smith .seealso: DMDAComputeJacobian1WithAdic()
59447c6ae99SBarry Smith 
59547c6ae99SBarry Smith @*/
596837ad101SBarry Smith PetscErrorCode  DMDAComputeFunction(DM da,PetscErrorCode (*lf)(void),Vec vu,Vec vfu,void *w)
59747c6ae99SBarry Smith {
59847c6ae99SBarry Smith   PetscErrorCode ierr;
59947c6ae99SBarry Smith   void           *u,*fu;
600aa219208SBarry Smith   DMDALocalInfo  info;
601aa219208SBarry Smith   PetscErrorCode (*f)(DMDALocalInfo*,void*,void*,void*) = (PetscErrorCode (*)(DMDALocalInfo*,void*,void*,void*))lf;
60247c6ae99SBarry Smith 
60347c6ae99SBarry Smith   PetscFunctionBegin;
604aa219208SBarry Smith   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
605aa219208SBarry Smith   ierr = DMDAVecGetArray(da,vu,&u);CHKERRQ(ierr);
606aa219208SBarry Smith   ierr = DMDAVecGetArray(da,vfu,&fu);CHKERRQ(ierr);
60747c6ae99SBarry Smith 
60839d508bbSBarry Smith   ierr = (*f)(&info,u,fu,w);CHKERRQ(ierr);
60947c6ae99SBarry Smith 
610aa219208SBarry Smith   ierr = DMDAVecRestoreArray(da,vu,&u);CHKERRQ(ierr);
611aa219208SBarry Smith   ierr = DMDAVecRestoreArray(da,vfu,&fu);CHKERRQ(ierr);
61247c6ae99SBarry Smith   PetscFunctionReturn(0);
61347c6ae99SBarry Smith }
61447c6ae99SBarry Smith 
61547c6ae99SBarry Smith #undef __FUNCT__
616837ad101SBarry Smith #define __FUNCT__ "DMDAComputeFunctionLocal"
61747c6ae99SBarry Smith /*@C
618837ad101SBarry Smith    DMDAComputeFunctionLocal - This is a universal function evaluation routine for
6199a42bb27SBarry Smith    a local DM function.
62047c6ae99SBarry Smith 
621aa219208SBarry Smith    Collective on DMDA
62247c6ae99SBarry Smith 
62347c6ae99SBarry Smith    Input Parameters:
6249a42bb27SBarry Smith +  da - the DM context
62547c6ae99SBarry Smith .  func - The local function
62647c6ae99SBarry Smith .  X - input vector
62747c6ae99SBarry Smith .  F - function vector
62847c6ae99SBarry Smith -  ctx - A user context
62947c6ae99SBarry Smith 
63047c6ae99SBarry Smith    Level: intermediate
63147c6ae99SBarry Smith 
632aa219208SBarry Smith .seealso: DMDASetLocalFunction(), DMDASetLocalJacobian(), DMDASetLocalAdicFunction(), DMDASetLocalAdicMFFunction(),
63347c6ae99SBarry Smith           SNESSetFunction(), SNESSetJacobian()
63447c6ae99SBarry Smith 
63547c6ae99SBarry Smith @*/
636837ad101SBarry Smith PetscErrorCode  DMDAComputeFunctionLocal(DM da, DMDALocalFunction1 func, Vec X, Vec F, void *ctx)
63747c6ae99SBarry Smith {
63847c6ae99SBarry Smith   Vec            localX;
639aa219208SBarry Smith   DMDALocalInfo  info;
64047c6ae99SBarry Smith   void           *u;
64147c6ae99SBarry Smith   void           *fu;
64247c6ae99SBarry Smith   PetscErrorCode ierr;
64347c6ae99SBarry Smith 
64447c6ae99SBarry Smith   PetscFunctionBegin;
6459a42bb27SBarry Smith   ierr = DMGetLocalVector(da,&localX);CHKERRQ(ierr);
64647c6ae99SBarry Smith   /*
64747c6ae99SBarry Smith      Scatter ghost points to local vector, using the 2-step process
6489a42bb27SBarry Smith         DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
64947c6ae99SBarry Smith   */
6509a42bb27SBarry Smith   ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
6519a42bb27SBarry Smith   ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
652aa219208SBarry Smith   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
653aa219208SBarry Smith   ierr = DMDAVecGetArray(da,localX,&u);CHKERRQ(ierr);
654aa219208SBarry Smith   ierr = DMDAVecGetArray(da,F,&fu);CHKERRQ(ierr);
65539d508bbSBarry Smith   ierr = (*func)(&info,u,fu,ctx);CHKERRQ(ierr);
656aa219208SBarry Smith   ierr = DMDAVecRestoreArray(da,localX,&u);CHKERRQ(ierr);
657aa219208SBarry Smith   ierr = DMDAVecRestoreArray(da,F,&fu);CHKERRQ(ierr);
6589a42bb27SBarry Smith   ierr = DMRestoreLocalVector(da,&localX);CHKERRQ(ierr);
65947c6ae99SBarry Smith   PetscFunctionReturn(0);
66047c6ae99SBarry Smith }
66147c6ae99SBarry Smith 
66247c6ae99SBarry Smith #undef __FUNCT__
663837ad101SBarry Smith #define __FUNCT__ "DMDAComputeFunctionLocalGhost"
66447c6ae99SBarry Smith /*@C
665837ad101SBarry Smith    DMDAComputeFunctionLocalGhost - This is a universal function evaluation routine for
6669a42bb27SBarry Smith    a local DM function, but the ghost values of the output are communicated and added.
66747c6ae99SBarry Smith 
668aa219208SBarry Smith    Collective on DMDA
66947c6ae99SBarry Smith 
67047c6ae99SBarry Smith    Input Parameters:
6719a42bb27SBarry Smith +  da - the DM context
67247c6ae99SBarry Smith .  func - The local function
67347c6ae99SBarry Smith .  X - input vector
67447c6ae99SBarry Smith .  F - function vector
67547c6ae99SBarry Smith -  ctx - A user context
67647c6ae99SBarry Smith 
67747c6ae99SBarry Smith    Level: intermediate
67847c6ae99SBarry Smith 
679aa219208SBarry Smith .seealso: DMDASetLocalFunction(), DMDASetLocalJacobian(), DMDASetLocalAdicFunction(), DMDASetLocalAdicMFFunction(),
68047c6ae99SBarry Smith           SNESSetFunction(), SNESSetJacobian()
68147c6ae99SBarry Smith 
68247c6ae99SBarry Smith @*/
683837ad101SBarry Smith PetscErrorCode  DMDAComputeFunctionLocalGhost(DM da, DMDALocalFunction1 func, Vec X, Vec F, void *ctx)
68447c6ae99SBarry Smith {
68547c6ae99SBarry Smith   Vec            localX, localF;
686aa219208SBarry Smith   DMDALocalInfo  info;
68747c6ae99SBarry Smith   void           *u;
68847c6ae99SBarry Smith   void           *fu;
68947c6ae99SBarry Smith   PetscErrorCode ierr;
69047c6ae99SBarry Smith 
69147c6ae99SBarry Smith   PetscFunctionBegin;
6929a42bb27SBarry Smith   ierr = DMGetLocalVector(da,&localX);CHKERRQ(ierr);
6939a42bb27SBarry Smith   ierr = DMGetLocalVector(da,&localF);CHKERRQ(ierr);
69447c6ae99SBarry Smith   /*
69547c6ae99SBarry Smith      Scatter ghost points to local vector, using the 2-step process
6969a42bb27SBarry Smith         DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
69747c6ae99SBarry Smith   */
6989a42bb27SBarry Smith   ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
6999a42bb27SBarry Smith   ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
70047c6ae99SBarry Smith   ierr = VecSet(F, 0.0);CHKERRQ(ierr);
70147c6ae99SBarry Smith   ierr = VecSet(localF, 0.0);CHKERRQ(ierr);
702aa219208SBarry Smith   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
703aa219208SBarry Smith   ierr = DMDAVecGetArray(da,localX,&u);CHKERRQ(ierr);
704aa219208SBarry Smith   ierr = DMDAVecGetArray(da,localF,&fu);CHKERRQ(ierr);
70539d508bbSBarry Smith   ierr = (*func)(&info,u,fu,ctx);CHKERRQ(ierr);
7069a42bb27SBarry Smith   ierr = DMLocalToGlobalBegin(da,localF,ADD_VALUES,F);CHKERRQ(ierr);
7079a42bb27SBarry Smith   ierr = DMLocalToGlobalEnd(da,localF,ADD_VALUES,F);CHKERRQ(ierr);
708aa219208SBarry Smith   ierr = DMDAVecRestoreArray(da,localX,&u);CHKERRQ(ierr);
709aa219208SBarry Smith   ierr = DMDAVecRestoreArray(da,localF,&fu);CHKERRQ(ierr);
7109a42bb27SBarry Smith   ierr = DMRestoreLocalVector(da,&localX);CHKERRQ(ierr);
7119a42bb27SBarry Smith   ierr = DMRestoreLocalVector(da,&localF);CHKERRQ(ierr);
71247c6ae99SBarry Smith   PetscFunctionReturn(0);
71347c6ae99SBarry Smith }
71447c6ae99SBarry Smith 
71547c6ae99SBarry Smith #undef __FUNCT__
716837ad101SBarry Smith #define __FUNCT__ "DMDAComputeFunction1"
71747c6ae99SBarry Smith /*@
718837ad101SBarry Smith     DMDAComputeFunction1 - Evaluates a user provided function on each processor that
719aa219208SBarry Smith         share a DMDA
72047c6ae99SBarry Smith 
72147c6ae99SBarry Smith    Input Parameters:
7229a42bb27SBarry Smith +    da - the DM that defines the grid
72347c6ae99SBarry Smith .    vu - input vector
72447c6ae99SBarry Smith .    vfu - output vector
72547c6ae99SBarry Smith -    w - any user data
72647c6ae99SBarry Smith 
72747c6ae99SBarry Smith     Notes: Does NOT do ghost updates on vu upon entry
72847c6ae99SBarry Smith 
72947c6ae99SBarry Smith     Level: advanced
73047c6ae99SBarry Smith 
731aa219208SBarry Smith .seealso: DMDAComputeJacobian1WithAdic()
73247c6ae99SBarry Smith 
73347c6ae99SBarry Smith @*/
734837ad101SBarry Smith PetscErrorCode  DMDAComputeFunction1(DM da,Vec vu,Vec vfu,void *w)
73547c6ae99SBarry Smith {
73647c6ae99SBarry Smith   PetscErrorCode ierr;
73747c6ae99SBarry Smith   void           *u,*fu;
738aa219208SBarry Smith   DMDALocalInfo  info;
73947c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
74047c6ae99SBarry Smith 
74147c6ae99SBarry Smith   PetscFunctionBegin;
7428b0a5094SBarry Smith   if (!dd->lf) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_NULL,"DMDASetLocalFunction() never called");
743aa219208SBarry Smith   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
744aa219208SBarry Smith   ierr = DMDAVecGetArray(da,vu,&u);CHKERRQ(ierr);
745aa219208SBarry Smith   ierr = DMDAVecGetArray(da,vfu,&fu);CHKERRQ(ierr);
74647c6ae99SBarry Smith 
74747c6ae99SBarry Smith   CHKMEMQ;
74839d508bbSBarry Smith   ierr = (*dd->lf)(&info,u,fu,w);CHKERRQ(ierr);
74947c6ae99SBarry Smith   CHKMEMQ;
75047c6ae99SBarry Smith 
751aa219208SBarry Smith   ierr = DMDAVecRestoreArray(da,vu,&u);CHKERRQ(ierr);
752aa219208SBarry Smith   ierr = DMDAVecRestoreArray(da,vfu,&fu);CHKERRQ(ierr);
75347c6ae99SBarry Smith   PetscFunctionReturn(0);
75447c6ae99SBarry Smith }
75547c6ae99SBarry Smith 
75647c6ae99SBarry Smith #undef __FUNCT__
757837ad101SBarry Smith #define __FUNCT__ "DMDAComputeFunctioniTest1"
758837ad101SBarry Smith PetscErrorCode  DMDAComputeFunctioniTest1(DM da,void *w)
75947c6ae99SBarry Smith {
76047c6ae99SBarry Smith   Vec            vu,fu,fui;
76147c6ae99SBarry Smith   PetscErrorCode ierr;
76247c6ae99SBarry Smith   PetscInt       i,n;
76347c6ae99SBarry Smith   PetscScalar    *ui;
76447c6ae99SBarry Smith   PetscRandom    rnd;
76547c6ae99SBarry Smith   PetscReal      norm;
76647c6ae99SBarry Smith 
76747c6ae99SBarry Smith   PetscFunctionBegin;
7689a42bb27SBarry Smith   ierr = DMGetLocalVector(da,&vu);CHKERRQ(ierr);
76947c6ae99SBarry Smith   ierr = PetscRandomCreate(PETSC_COMM_SELF,&rnd);CHKERRQ(ierr);
77047c6ae99SBarry Smith   ierr = PetscRandomSetFromOptions(rnd);CHKERRQ(ierr);
77147c6ae99SBarry Smith   ierr = VecSetRandom(vu,rnd);CHKERRQ(ierr);
772fcfd50ebSBarry Smith   ierr = PetscRandomDestroy(&rnd);CHKERRQ(ierr);
77347c6ae99SBarry Smith 
7749a42bb27SBarry Smith   ierr = DMGetGlobalVector(da,&fu);CHKERRQ(ierr);
7759a42bb27SBarry Smith   ierr = DMGetGlobalVector(da,&fui);CHKERRQ(ierr);
77647c6ae99SBarry Smith 
777837ad101SBarry Smith   ierr = DMDAComputeFunction1(da,vu,fu,w);CHKERRQ(ierr);
77847c6ae99SBarry Smith 
77947c6ae99SBarry Smith   ierr = VecGetArray(fui,&ui);CHKERRQ(ierr);
78047c6ae99SBarry Smith   ierr = VecGetLocalSize(fui,&n);CHKERRQ(ierr);
78147c6ae99SBarry Smith   for (i=0; i<n; i++) {
782837ad101SBarry Smith     ierr = DMDAComputeFunctioni1(da,i,vu,ui+i,w);CHKERRQ(ierr);
78347c6ae99SBarry Smith   }
78447c6ae99SBarry Smith   ierr = VecRestoreArray(fui,&ui);CHKERRQ(ierr);
78547c6ae99SBarry Smith 
78647c6ae99SBarry Smith   ierr = VecAXPY(fui,-1.0,fu);CHKERRQ(ierr);
78747c6ae99SBarry Smith   ierr = VecNorm(fui,NORM_2,&norm);CHKERRQ(ierr);
78847c6ae99SBarry Smith   ierr = PetscPrintf(((PetscObject)da)->comm,"Norm of difference in vectors %G\n",norm);CHKERRQ(ierr);
78947c6ae99SBarry Smith   ierr = VecView(fu,0);CHKERRQ(ierr);
79047c6ae99SBarry Smith   ierr = VecView(fui,0);CHKERRQ(ierr);
79147c6ae99SBarry Smith 
7929a42bb27SBarry Smith   ierr = DMRestoreLocalVector(da,&vu);CHKERRQ(ierr);
7939a42bb27SBarry Smith   ierr = DMRestoreGlobalVector(da,&fu);CHKERRQ(ierr);
7949a42bb27SBarry Smith   ierr = DMRestoreGlobalVector(da,&fui);CHKERRQ(ierr);
79547c6ae99SBarry Smith   PetscFunctionReturn(0);
79647c6ae99SBarry Smith }
79747c6ae99SBarry Smith 
79847c6ae99SBarry Smith #undef __FUNCT__
799837ad101SBarry Smith #define __FUNCT__ "DMDAComputeFunctioni1"
80047c6ae99SBarry Smith /*@
801837ad101SBarry Smith     DMDAComputeFunctioni1 - Evaluates a user provided point-wise function
80247c6ae99SBarry Smith 
80347c6ae99SBarry Smith    Input Parameters:
8049a42bb27SBarry Smith +    da - the DM that defines the grid
80547c6ae99SBarry Smith .    i - the component of the function we wish to compute (must be local)
80647c6ae99SBarry Smith .    vu - input vector
80747c6ae99SBarry Smith .    vfu - output value
80847c6ae99SBarry Smith -    w - any user data
80947c6ae99SBarry Smith 
81047c6ae99SBarry Smith     Notes: Does NOT do ghost updates on vu upon entry
81147c6ae99SBarry Smith 
81247c6ae99SBarry Smith     Level: advanced
81347c6ae99SBarry Smith 
814aa219208SBarry Smith .seealso: DMDAComputeJacobian1WithAdic()
81547c6ae99SBarry Smith 
81647c6ae99SBarry Smith @*/
817837ad101SBarry Smith PetscErrorCode  DMDAComputeFunctioni1(DM da,PetscInt i,Vec vu,PetscScalar *vfu,void *w)
81847c6ae99SBarry Smith {
81947c6ae99SBarry Smith   PetscErrorCode ierr;
82047c6ae99SBarry Smith   void           *u;
821aa219208SBarry Smith   DMDALocalInfo  info;
82247c6ae99SBarry Smith   MatStencil     stencil;
82347c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
82447c6ae99SBarry Smith 
82547c6ae99SBarry Smith   PetscFunctionBegin;
82647c6ae99SBarry Smith 
827aa219208SBarry Smith   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
828aa219208SBarry Smith   ierr = DMDAVecGetArray(da,vu,&u);CHKERRQ(ierr);
82947c6ae99SBarry Smith 
83047c6ae99SBarry Smith   /* figure out stencil value from i */
83147c6ae99SBarry Smith   stencil.c = i % info.dof;
83247c6ae99SBarry Smith   stencil.i = (i % (info.xm*info.dof))/info.dof;
83347c6ae99SBarry Smith   stencil.j = (i % (info.xm*info.ym*info.dof))/(info.xm*info.dof);
83447c6ae99SBarry Smith   stencil.k = i/(info.xm*info.ym*info.dof);
83547c6ae99SBarry Smith 
83647c6ae99SBarry Smith   ierr = (*dd->lfi)(&info,&stencil,u,vfu,w);CHKERRQ(ierr);
83747c6ae99SBarry Smith 
838aa219208SBarry Smith   ierr = DMDAVecRestoreArray(da,vu,&u);CHKERRQ(ierr);
83947c6ae99SBarry Smith   PetscFunctionReturn(0);
84047c6ae99SBarry Smith }
84147c6ae99SBarry Smith 
84247c6ae99SBarry Smith #undef __FUNCT__
843837ad101SBarry Smith #define __FUNCT__ "DMDAComputeFunctionib1"
84447c6ae99SBarry Smith /*@
845837ad101SBarry Smith     DMDAComputeFunctionib1 - Evaluates a user provided point-block function
84647c6ae99SBarry Smith 
84747c6ae99SBarry Smith    Input Parameters:
8489a42bb27SBarry Smith +    da - the DM that defines the grid
84947c6ae99SBarry Smith .    i - the component of the function we wish to compute (must be local)
85047c6ae99SBarry Smith .    vu - input vector
85147c6ae99SBarry Smith .    vfu - output value
85247c6ae99SBarry Smith -    w - any user data
85347c6ae99SBarry Smith 
85447c6ae99SBarry Smith     Notes: Does NOT do ghost updates on vu upon entry
85547c6ae99SBarry Smith 
85647c6ae99SBarry Smith     Level: advanced
85747c6ae99SBarry Smith 
858aa219208SBarry Smith .seealso: DMDAComputeJacobian1WithAdic()
85947c6ae99SBarry Smith 
86047c6ae99SBarry Smith @*/
861837ad101SBarry Smith PetscErrorCode  DMDAComputeFunctionib1(DM da,PetscInt i,Vec vu,PetscScalar *vfu,void *w)
86247c6ae99SBarry Smith {
86347c6ae99SBarry Smith   PetscErrorCode ierr;
86447c6ae99SBarry Smith   void           *u;
865aa219208SBarry Smith   DMDALocalInfo  info;
86647c6ae99SBarry Smith   MatStencil     stencil;
86747c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
86847c6ae99SBarry Smith 
86947c6ae99SBarry Smith   PetscFunctionBegin;
870aa219208SBarry Smith   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
871aa219208SBarry Smith   ierr = DMDAVecGetArray(da,vu,&u);CHKERRQ(ierr);
87247c6ae99SBarry Smith 
87347c6ae99SBarry Smith   /* figure out stencil value from i */
87447c6ae99SBarry Smith   stencil.c = i % info.dof;
87547c6ae99SBarry Smith   if (stencil.c) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONG,"Point-block functions can only be called for the entire block");
87647c6ae99SBarry Smith   stencil.i = (i % (info.xm*info.dof))/info.dof;
87747c6ae99SBarry Smith   stencil.j = (i % (info.xm*info.ym*info.dof))/(info.xm*info.dof);
87847c6ae99SBarry Smith   stencil.k = i/(info.xm*info.ym*info.dof);
87947c6ae99SBarry Smith 
88047c6ae99SBarry Smith   ierr = (*dd->lfib)(&info,&stencil,u,vfu,w);CHKERRQ(ierr);
88147c6ae99SBarry Smith 
882aa219208SBarry Smith   ierr = DMDAVecRestoreArray(da,vu,&u);CHKERRQ(ierr);
88347c6ae99SBarry Smith   PetscFunctionReturn(0);
88447c6ae99SBarry Smith }
88547c6ae99SBarry Smith 
88647c6ae99SBarry Smith #if defined(new)
88747c6ae99SBarry Smith #undef __FUNCT__
888aa219208SBarry Smith #define __FUNCT__ "DMDAGetDiagonal_MFFD"
88947c6ae99SBarry Smith /*
890aa219208SBarry Smith   DMDAGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix where local
891aa219208SBarry Smith     function lives on a DMDA
89247c6ae99SBarry Smith 
89347c6ae99SBarry Smith         y ~= (F(u + ha) - F(u))/h,
89447c6ae99SBarry Smith   where F = nonlinear function, as set by SNESSetFunction()
89547c6ae99SBarry Smith         u = current iterate
89647c6ae99SBarry Smith         h = difference interval
89747c6ae99SBarry Smith */
898aa219208SBarry Smith PetscErrorCode DMDAGetDiagonal_MFFD(DM da,Vec U,Vec a)
89947c6ae99SBarry Smith {
90047c6ae99SBarry Smith   PetscScalar    h,*aa,*ww,v;
90147c6ae99SBarry Smith   PetscReal      epsilon = PETSC_SQRT_MACHINE_EPSILON,umin = 100.0*PETSC_SQRT_MACHINE_EPSILON;
90247c6ae99SBarry Smith   PetscErrorCode ierr;
90347c6ae99SBarry Smith   PetscInt       gI,nI;
90447c6ae99SBarry Smith   MatStencil     stencil;
905aa219208SBarry Smith   DMDALocalInfo  info;
90647c6ae99SBarry Smith 
90747c6ae99SBarry Smith   PetscFunctionBegin;
90847c6ae99SBarry Smith   ierr = (*ctx->func)(0,U,a,ctx->funcctx);CHKERRQ(ierr);
90947c6ae99SBarry Smith   ierr = (*ctx->funcisetbase)(U,ctx->funcctx);CHKERRQ(ierr);
91047c6ae99SBarry Smith 
91147c6ae99SBarry Smith   ierr = VecGetArray(U,&ww);CHKERRQ(ierr);
91247c6ae99SBarry Smith   ierr = VecGetArray(a,&aa);CHKERRQ(ierr);
91347c6ae99SBarry Smith 
91447c6ae99SBarry Smith   nI = 0;
91547c6ae99SBarry Smith     h  = ww[gI];
91647c6ae99SBarry Smith     if (h == 0.0) h = 1.0;
91747c6ae99SBarry Smith #if !defined(PETSC_USE_COMPLEX)
91847c6ae99SBarry Smith     if (h < umin && h >= 0.0)      h = umin;
91947c6ae99SBarry Smith     else if (h < 0.0 && h > -umin) h = -umin;
92047c6ae99SBarry Smith #else
92147c6ae99SBarry Smith     if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0)     h = umin;
92247c6ae99SBarry Smith     else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin;
92347c6ae99SBarry Smith #endif
92447c6ae99SBarry Smith     h     *= epsilon;
92547c6ae99SBarry Smith 
92647c6ae99SBarry Smith     ww[gI] += h;
92747c6ae99SBarry Smith     ierr          = (*ctx->funci)(i,w,&v,ctx->funcctx);CHKERRQ(ierr);
92847c6ae99SBarry Smith     aa[nI]  = (v - aa[nI])/h;
92947c6ae99SBarry Smith     ww[gI] -= h;
93047c6ae99SBarry Smith     nI++;
93147c6ae99SBarry Smith   }
93247c6ae99SBarry Smith   ierr = VecRestoreArray(U,&ww);CHKERRQ(ierr);
93347c6ae99SBarry Smith   ierr = VecRestoreArray(a,&aa);CHKERRQ(ierr);
93447c6ae99SBarry Smith   PetscFunctionReturn(0);
93547c6ae99SBarry Smith }
93647c6ae99SBarry Smith #endif
93747c6ae99SBarry Smith 
93847c6ae99SBarry Smith #if defined(PETSC_HAVE_ADIC)
93947c6ae99SBarry Smith EXTERN_C_BEGIN
940c6db04a5SJed Brown #include <adic/ad_utils.h>
94147c6ae99SBarry Smith EXTERN_C_END
94247c6ae99SBarry Smith 
94347c6ae99SBarry Smith #undef __FUNCT__
944aa219208SBarry Smith #define __FUNCT__ "DMDAComputeJacobian1WithAdic"
94547c6ae99SBarry Smith /*@C
946aa219208SBarry Smith     DMDAComputeJacobian1WithAdic - Evaluates a adiC provided Jacobian function on each processor that
947aa219208SBarry Smith         share a DMDA
94847c6ae99SBarry Smith 
94947c6ae99SBarry Smith    Input Parameters:
9509a42bb27SBarry Smith +    da - the DM that defines the grid
95147c6ae99SBarry Smith .    vu - input vector (ghosted)
95247c6ae99SBarry Smith .    J - output matrix
95347c6ae99SBarry Smith -    w - any user data
95447c6ae99SBarry Smith 
95547c6ae99SBarry Smith    Level: advanced
95647c6ae99SBarry Smith 
95747c6ae99SBarry Smith     Notes: Does NOT do ghost updates on vu upon entry
95847c6ae99SBarry Smith 
959837ad101SBarry Smith .seealso: DMDAComputeFunction1()
96047c6ae99SBarry Smith 
96147c6ae99SBarry Smith @*/
9627087cfbeSBarry Smith PetscErrorCode  DMDAComputeJacobian1WithAdic(DM da,Vec vu,Mat J,void *w)
96347c6ae99SBarry Smith {
96447c6ae99SBarry Smith   PetscErrorCode ierr;
96547c6ae99SBarry Smith   PetscInt       gtdof,tdof;
96647c6ae99SBarry Smith   PetscScalar    *ustart;
967aa219208SBarry Smith   DMDALocalInfo  info;
96847c6ae99SBarry Smith   void           *ad_u,*ad_f,*ad_ustart,*ad_fstart;
96947c6ae99SBarry Smith   ISColoring     iscoloring;
97047c6ae99SBarry Smith 
97147c6ae99SBarry Smith   PetscFunctionBegin;
972aa219208SBarry Smith   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
97347c6ae99SBarry Smith 
97447c6ae99SBarry Smith   PetscADResetIndep();
97547c6ae99SBarry Smith 
97647c6ae99SBarry Smith   /* get space for derivative objects.  */
977aa219208SBarry Smith   ierr = DMDAGetAdicArray(da,PETSC_TRUE,&ad_u,&ad_ustart,&gtdof);CHKERRQ(ierr);
978aa219208SBarry Smith   ierr = DMDAGetAdicArray(da,PETSC_FALSE,&ad_f,&ad_fstart,&tdof);CHKERRQ(ierr);
97947c6ae99SBarry Smith   ierr = VecGetArray(vu,&ustart);CHKERRQ(ierr);
980e727c939SJed Brown   ierr = DMCreateColoring(da,IS_COLORING_GHOSTED,MATAIJ,&iscoloring);CHKERRQ(ierr);
98147c6ae99SBarry Smith 
98247c6ae99SBarry Smith   PetscADSetValueAndColor(ad_ustart,gtdof,iscoloring->colors,ustart);
98347c6ae99SBarry Smith 
98447c6ae99SBarry Smith   ierr = VecRestoreArray(vu,&ustart);CHKERRQ(ierr);
985fcfd50ebSBarry Smith   ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr);
98647c6ae99SBarry Smith   ierr = PetscADIncrementTotalGradSize(iscoloring->n);CHKERRQ(ierr);
98747c6ae99SBarry Smith   PetscADSetIndepDone();
98847c6ae99SBarry Smith 
989aa219208SBarry Smith   ierr = PetscLogEventBegin(DMDA_LocalADFunction,0,0,0,0);CHKERRQ(ierr);
99047c6ae99SBarry Smith   ierr = (*dd->adic_lf)(&info,ad_u,ad_f,w);CHKERRQ(ierr);
991aa219208SBarry Smith   ierr = PetscLogEventEnd(DMDA_LocalADFunction,0,0,0,0);CHKERRQ(ierr);
99247c6ae99SBarry Smith 
99347c6ae99SBarry Smith   /* stick the values into the matrix */
99447c6ae99SBarry Smith   ierr = MatSetValuesAdic(J,(PetscScalar**)ad_fstart);CHKERRQ(ierr);
99547c6ae99SBarry Smith 
99647c6ae99SBarry Smith   /* return space for derivative objects.  */
997aa219208SBarry Smith   ierr = DMDARestoreAdicArray(da,PETSC_TRUE,&ad_u,&ad_ustart,&gtdof);CHKERRQ(ierr);
998aa219208SBarry Smith   ierr = DMDARestoreAdicArray(da,PETSC_FALSE,&ad_f,&ad_fstart,&tdof);CHKERRQ(ierr);
99947c6ae99SBarry Smith   PetscFunctionReturn(0);
100047c6ae99SBarry Smith }
100147c6ae99SBarry Smith 
100247c6ae99SBarry Smith #undef __FUNCT__
1003aa219208SBarry Smith #define __FUNCT__ "DMDAMultiplyByJacobian1WithAdic"
100447c6ae99SBarry Smith /*@C
1005aa219208SBarry Smith     DMDAMultiplyByJacobian1WithAdic - Applies an ADIC-provided Jacobian function to a vector on
1006aa219208SBarry Smith     each processor that shares a DMDA.
100747c6ae99SBarry Smith 
100847c6ae99SBarry Smith     Input Parameters:
10099a42bb27SBarry Smith +   da - the DM that defines the grid
101047c6ae99SBarry Smith .   vu - Jacobian is computed at this point (ghosted)
101147c6ae99SBarry Smith .   v - product is done on this vector (ghosted)
101247c6ae99SBarry Smith .   fu - output vector = J(vu)*v (not ghosted)
101347c6ae99SBarry Smith -   w - any user data
101447c6ae99SBarry Smith 
101547c6ae99SBarry Smith     Notes:
101647c6ae99SBarry Smith     This routine does NOT do ghost updates on vu upon entry.
101747c6ae99SBarry Smith 
101847c6ae99SBarry Smith    Level: advanced
101947c6ae99SBarry Smith 
1020837ad101SBarry Smith .seealso: DMDAComputeFunction1()
102147c6ae99SBarry Smith 
102247c6ae99SBarry Smith @*/
10237087cfbeSBarry Smith PetscErrorCode  DMDAMultiplyByJacobian1WithAdic(DM da,Vec vu,Vec v,Vec f,void *w)
102447c6ae99SBarry Smith {
102547c6ae99SBarry Smith   PetscErrorCode ierr;
102647c6ae99SBarry Smith   PetscInt       i,gtdof,tdof;
102747c6ae99SBarry Smith   PetscScalar    *avu,*av,*af,*ad_vustart,*ad_fstart;
1028aa219208SBarry Smith   DMDALocalInfo  info;
102947c6ae99SBarry Smith   void           *ad_vu,*ad_f;
103047c6ae99SBarry Smith 
103147c6ae99SBarry Smith   PetscFunctionBegin;
1032aa219208SBarry Smith   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
103347c6ae99SBarry Smith 
103447c6ae99SBarry Smith   /* get space for derivative objects.  */
1035aa219208SBarry Smith   ierr = DMDAGetAdicMFArray(da,PETSC_TRUE,&ad_vu,&ad_vustart,&gtdof);CHKERRQ(ierr);
1036aa219208SBarry Smith   ierr = DMDAGetAdicMFArray(da,PETSC_FALSE,&ad_f,&ad_fstart,&tdof);CHKERRQ(ierr);
103747c6ae99SBarry Smith 
103847c6ae99SBarry Smith   /* copy input vector into derivative object */
103947c6ae99SBarry Smith   ierr = VecGetArray(vu,&avu);CHKERRQ(ierr);
104047c6ae99SBarry Smith   ierr = VecGetArray(v,&av);CHKERRQ(ierr);
104147c6ae99SBarry Smith   for (i=0; i<gtdof; i++) {
104247c6ae99SBarry Smith     ad_vustart[2*i]   = avu[i];
104347c6ae99SBarry Smith     ad_vustart[2*i+1] = av[i];
104447c6ae99SBarry Smith   }
104547c6ae99SBarry Smith   ierr = VecRestoreArray(vu,&avu);CHKERRQ(ierr);
104647c6ae99SBarry Smith   ierr = VecRestoreArray(v,&av);CHKERRQ(ierr);
104747c6ae99SBarry Smith 
104847c6ae99SBarry Smith   PetscADResetIndep();
104947c6ae99SBarry Smith   ierr = PetscADIncrementTotalGradSize(1);CHKERRQ(ierr);
105047c6ae99SBarry Smith   PetscADSetIndepDone();
105147c6ae99SBarry Smith 
105247c6ae99SBarry Smith   ierr = (*dd->adicmf_lf)(&info,ad_vu,ad_f,w);CHKERRQ(ierr);
105347c6ae99SBarry Smith 
105447c6ae99SBarry Smith   /* stick the values into the vector */
105547c6ae99SBarry Smith   ierr = VecGetArray(f,&af);CHKERRQ(ierr);
105647c6ae99SBarry Smith   for (i=0; i<tdof; i++) {
105747c6ae99SBarry Smith     af[i] = ad_fstart[2*i+1];
105847c6ae99SBarry Smith   }
105947c6ae99SBarry Smith   ierr = VecRestoreArray(f,&af);CHKERRQ(ierr);
106047c6ae99SBarry Smith 
106147c6ae99SBarry Smith   /* return space for derivative objects.  */
1062aa219208SBarry Smith   ierr = DMDARestoreAdicMFArray(da,PETSC_TRUE,&ad_vu,&ad_vustart,&gtdof);CHKERRQ(ierr);
1063aa219208SBarry Smith   ierr = DMDARestoreAdicMFArray(da,PETSC_FALSE,&ad_f,&ad_fstart,&tdof);CHKERRQ(ierr);
106447c6ae99SBarry Smith   PetscFunctionReturn(0);
106547c6ae99SBarry Smith }
106647c6ae99SBarry Smith #endif
106747c6ae99SBarry Smith 
106847c6ae99SBarry Smith #undef __FUNCT__
1069aa219208SBarry Smith #define __FUNCT__ "DMDAComputeJacobian1"
107047c6ae99SBarry Smith /*@
1071aa219208SBarry Smith     DMDAComputeJacobian1 - Evaluates a local Jacobian function on each processor that
1072aa219208SBarry Smith         share a DMDA
107347c6ae99SBarry Smith 
107447c6ae99SBarry Smith    Input Parameters:
10759a42bb27SBarry Smith +    da - the DM that defines the grid
107647c6ae99SBarry Smith .    vu - input vector (ghosted)
107747c6ae99SBarry Smith .    J - output matrix
107847c6ae99SBarry Smith -    w - any user data
107947c6ae99SBarry Smith 
108047c6ae99SBarry Smith     Notes: Does NOT do ghost updates on vu upon entry
108147c6ae99SBarry Smith 
108247c6ae99SBarry Smith     Level: advanced
108347c6ae99SBarry Smith 
1084837ad101SBarry Smith .seealso: DMDAComputeFunction1()
108547c6ae99SBarry Smith 
108647c6ae99SBarry Smith @*/
10877087cfbeSBarry Smith PetscErrorCode  DMDAComputeJacobian1(DM da,Vec vu,Mat J,void *w)
108847c6ae99SBarry Smith {
108947c6ae99SBarry Smith   PetscErrorCode ierr;
109047c6ae99SBarry Smith   void           *u;
1091aa219208SBarry Smith   DMDALocalInfo  info;
109247c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
109347c6ae99SBarry Smith 
109447c6ae99SBarry Smith   PetscFunctionBegin;
1095aa219208SBarry Smith   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
1096aa219208SBarry Smith   ierr = DMDAVecGetArray(da,vu,&u);CHKERRQ(ierr);
109747c6ae99SBarry Smith   ierr = (*dd->lj)(&info,u,J,w);CHKERRQ(ierr);
1098aa219208SBarry Smith   ierr = DMDAVecRestoreArray(da,vu,&u);CHKERRQ(ierr);
109947c6ae99SBarry Smith   PetscFunctionReturn(0);
110047c6ae99SBarry Smith }
110147c6ae99SBarry Smith 
110247c6ae99SBarry Smith 
110347c6ae99SBarry Smith #undef __FUNCT__
1104aa219208SBarry Smith #define __FUNCT__ "DMDAComputeJacobian1WithAdifor"
110547c6ae99SBarry Smith /*
1106aa219208SBarry Smith     DMDAComputeJacobian1WithAdifor - Evaluates a ADIFOR provided Jacobian local function on each processor that
1107aa219208SBarry Smith         share a DMDA
110847c6ae99SBarry Smith 
110947c6ae99SBarry Smith    Input Parameters:
11109a42bb27SBarry Smith +    da - the DM that defines the grid
111147c6ae99SBarry Smith .    vu - input vector (ghosted)
111247c6ae99SBarry Smith .    J - output matrix
111347c6ae99SBarry Smith -    w - any user data
111447c6ae99SBarry Smith 
111547c6ae99SBarry Smith     Notes: Does NOT do ghost updates on vu upon entry
111647c6ae99SBarry Smith 
1117837ad101SBarry Smith .seealso: DMDAComputeFunction1()
111847c6ae99SBarry Smith 
111947c6ae99SBarry Smith */
11207087cfbeSBarry Smith PetscErrorCode  DMDAComputeJacobian1WithAdifor(DM da,Vec vu,Mat J,void *w)
112147c6ae99SBarry Smith {
112247c6ae99SBarry Smith   PetscErrorCode  ierr;
112347c6ae99SBarry Smith   PetscInt        i,Nc,N;
112447c6ae99SBarry Smith   ISColoringValue *color;
1125aa219208SBarry Smith   DMDALocalInfo   info;
112647c6ae99SBarry Smith   PetscScalar     *u,*g_u,*g_f,*f = 0,*p_u;
112747c6ae99SBarry Smith   ISColoring      iscoloring;
112847c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
1129aa219208SBarry Smith   void            (*lf)(PetscInt*,DMDALocalInfo*,PetscScalar*,PetscScalar*,PetscInt*,PetscScalar*,PetscScalar*,PetscInt*,void*,PetscErrorCode*) =
1130aa219208SBarry Smith                   (void (*)(PetscInt*,DMDALocalInfo*,PetscScalar*,PetscScalar*,PetscInt*,PetscScalar*,PetscScalar*,PetscInt*,void*,PetscErrorCode*))*dd->adifor_lf;
113147c6ae99SBarry Smith 
113247c6ae99SBarry Smith   PetscFunctionBegin;
1133e727c939SJed Brown   ierr = DMCreateColoring(da,IS_COLORING_GHOSTED,MATAIJ,&iscoloring);CHKERRQ(ierr);
113447c6ae99SBarry Smith   Nc   = iscoloring->n;
1135aa219208SBarry Smith   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
113647c6ae99SBarry Smith   N    = info.gxm*info.gym*info.gzm*info.dof;
113747c6ae99SBarry Smith 
113847c6ae99SBarry Smith   /* get space for derivative objects.  */
113947c6ae99SBarry Smith   ierr  = PetscMalloc(Nc*info.gxm*info.gym*info.gzm*info.dof*sizeof(PetscScalar),&g_u);CHKERRQ(ierr);
114047c6ae99SBarry Smith   ierr  = PetscMemzero(g_u,Nc*info.gxm*info.gym*info.gzm*info.dof*sizeof(PetscScalar));CHKERRQ(ierr);
114147c6ae99SBarry Smith   p_u   = g_u;
114247c6ae99SBarry Smith   color = iscoloring->colors;
114347c6ae99SBarry Smith   for (i=0; i<N; i++) {
114447c6ae99SBarry Smith     p_u[*color++] = 1.0;
114547c6ae99SBarry Smith     p_u          += Nc;
114647c6ae99SBarry Smith   }
1147fcfd50ebSBarry Smith   ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr);
114847c6ae99SBarry 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);
114947c6ae99SBarry Smith 
115047c6ae99SBarry Smith   /* Seed the input array g_u with coloring information */
115147c6ae99SBarry Smith 
115247c6ae99SBarry Smith   ierr = VecGetArray(vu,&u);CHKERRQ(ierr);
115347c6ae99SBarry Smith   (lf)(&Nc,&info,u,g_u,&Nc,f,g_f,&Nc,w,&ierr);CHKERRQ(ierr);
115447c6ae99SBarry Smith   ierr = VecRestoreArray(vu,&u);CHKERRQ(ierr);
115547c6ae99SBarry Smith 
115647c6ae99SBarry Smith   /* stick the values into the matrix */
115747c6ae99SBarry Smith   /* PetscScalarView(Nc*info.xm*info.ym,g_f,0); */
115847c6ae99SBarry Smith   ierr = MatSetValuesAdifor(J,Nc,g_f);CHKERRQ(ierr);
115947c6ae99SBarry Smith 
116047c6ae99SBarry Smith   /* return space for derivative objects.  */
116147c6ae99SBarry Smith   ierr = PetscFree(g_u);CHKERRQ(ierr);
116247c6ae99SBarry Smith   ierr = PetscFree2(g_f,f);CHKERRQ(ierr);
116347c6ae99SBarry Smith   PetscFunctionReturn(0);
116447c6ae99SBarry Smith }
116547c6ae99SBarry Smith 
116647c6ae99SBarry Smith #undef __FUNCT__
1167aa219208SBarry Smith #define __FUNCT__ "DMDAFormJacobianLocal"
116847c6ae99SBarry Smith /*@C
1169aa219208SBarry Smith    DMDAFormjacobianLocal - This is a universal Jacobian evaluation routine for
11709a42bb27SBarry Smith    a local DM function.
117147c6ae99SBarry Smith 
1172aa219208SBarry Smith    Collective on DMDA
117347c6ae99SBarry Smith 
117447c6ae99SBarry Smith    Input Parameters:
11759a42bb27SBarry Smith +  da - the DM context
117647c6ae99SBarry Smith .  func - The local function
117747c6ae99SBarry Smith .  X - input vector
117847c6ae99SBarry Smith .  J - Jacobian matrix
117947c6ae99SBarry Smith -  ctx - A user context
118047c6ae99SBarry Smith 
118147c6ae99SBarry Smith    Level: intermediate
118247c6ae99SBarry Smith 
1183aa219208SBarry Smith .seealso: DMDASetLocalFunction(), DMDASetLocalJacobian(), DMDASetLocalAdicFunction(), DMDASetLocalAdicMFFunction(),
118447c6ae99SBarry Smith           SNESSetFunction(), SNESSetJacobian()
118547c6ae99SBarry Smith 
118647c6ae99SBarry Smith @*/
11877087cfbeSBarry Smith PetscErrorCode  DMDAFormJacobianLocal(DM da, DMDALocalFunction1 func, Vec X, Mat J, void *ctx)
118847c6ae99SBarry Smith {
118947c6ae99SBarry Smith   Vec            localX;
1190aa219208SBarry Smith   DMDALocalInfo  info;
119147c6ae99SBarry Smith   void           *u;
119247c6ae99SBarry Smith   PetscErrorCode ierr;
119347c6ae99SBarry Smith 
119447c6ae99SBarry Smith   PetscFunctionBegin;
11959a42bb27SBarry Smith   ierr = DMGetLocalVector(da,&localX);CHKERRQ(ierr);
119647c6ae99SBarry Smith   /*
119747c6ae99SBarry Smith      Scatter ghost points to local vector, using the 2-step process
11989a42bb27SBarry Smith         DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
119947c6ae99SBarry Smith   */
12009a42bb27SBarry Smith   ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
12019a42bb27SBarry Smith   ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
1202aa219208SBarry Smith   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
1203aa219208SBarry Smith   ierr = DMDAVecGetArray(da,localX,&u);CHKERRQ(ierr);
120439d508bbSBarry Smith   ierr = (*func)(&info,u,J,ctx);CHKERRQ(ierr);
1205aa219208SBarry Smith   ierr = DMDAVecRestoreArray(da,localX,&u);CHKERRQ(ierr);
12069a42bb27SBarry Smith   ierr = DMRestoreLocalVector(da,&localX);CHKERRQ(ierr);
120747c6ae99SBarry Smith   PetscFunctionReturn(0);
120847c6ae99SBarry Smith }
120947c6ae99SBarry Smith 
121047c6ae99SBarry Smith #undef __FUNCT__
1211aa219208SBarry Smith #define __FUNCT__ "DMDAMultiplyByJacobian1WithAD"
121247c6ae99SBarry Smith /*@C
1213aa219208SBarry Smith     DMDAMultiplyByJacobian1WithAD - Applies a Jacobian function supplied by ADIFOR or ADIC
1214aa219208SBarry Smith     to a vector on each processor that shares a DMDA.
121547c6ae99SBarry Smith 
121647c6ae99SBarry Smith    Input Parameters:
12179a42bb27SBarry Smith +    da - the DM that defines the grid
121847c6ae99SBarry Smith .    vu - Jacobian is computed at this point (ghosted)
121947c6ae99SBarry Smith .    v - product is done on this vector (ghosted)
122047c6ae99SBarry Smith .    fu - output vector = J(vu)*v (not ghosted)
122147c6ae99SBarry Smith -    w - any user data
122247c6ae99SBarry Smith 
122347c6ae99SBarry Smith     Notes:
122447c6ae99SBarry Smith     This routine does NOT do ghost updates on vu and v upon entry.
122547c6ae99SBarry Smith 
1226aa219208SBarry Smith     Automatically calls DMDAMultiplyByJacobian1WithAdifor() or DMDAMultiplyByJacobian1WithAdic()
1227aa219208SBarry Smith     depending on whether DMDASetLocalAdicMFFunction() or DMDASetLocalAdiforMFFunction() was called.
122847c6ae99SBarry Smith 
122947c6ae99SBarry Smith    Level: advanced
123047c6ae99SBarry Smith 
1231837ad101SBarry Smith .seealso: DMDAComputeFunction1(), DMDAMultiplyByJacobian1WithAdifor(), DMDAMultiplyByJacobian1WithAdic()
123247c6ae99SBarry Smith 
123347c6ae99SBarry Smith @*/
12347087cfbeSBarry Smith PetscErrorCode  DMDAMultiplyByJacobian1WithAD(DM da,Vec u,Vec v,Vec f,void *w)
123547c6ae99SBarry Smith {
123647c6ae99SBarry Smith   PetscErrorCode ierr;
123747c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
123847c6ae99SBarry Smith 
123947c6ae99SBarry Smith   PetscFunctionBegin;
124047c6ae99SBarry Smith   if (dd->adicmf_lf) {
124147c6ae99SBarry Smith #if defined(PETSC_HAVE_ADIC)
1242aa219208SBarry Smith     ierr = DMDAMultiplyByJacobian1WithAdic(da,u,v,f,w);CHKERRQ(ierr);
124347c6ae99SBarry Smith #else
124447c6ae99SBarry Smith     SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP_SYS,"Requires ADIC to be installed and cannot use complex numbers");
124547c6ae99SBarry Smith #endif
124647c6ae99SBarry Smith   } else if (dd->adiformf_lf) {
1247aa219208SBarry Smith     ierr = DMDAMultiplyByJacobian1WithAdifor(da,u,v,f,w);CHKERRQ(ierr);
124847c6ae99SBarry Smith   } else {
1249aa219208SBarry Smith     SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ORDER,"Must call DMDASetLocalAdiforMFFunction() or DMDASetLocalAdicMFFunction() before using");
125047c6ae99SBarry Smith   }
125147c6ae99SBarry Smith   PetscFunctionReturn(0);
125247c6ae99SBarry Smith }
125347c6ae99SBarry Smith 
125447c6ae99SBarry Smith 
125547c6ae99SBarry Smith #undef __FUNCT__
1256aa219208SBarry Smith #define __FUNCT__ "DMDAMultiplyByJacobian1WithAdifor"
125747c6ae99SBarry Smith /*@C
1258aa219208SBarry Smith     DMDAMultiplyByJacobian1WithAdifor - Applies a ADIFOR provided Jacobian function on each processor that
12599a42bb27SBarry Smith         share a DM to a vector
126047c6ae99SBarry Smith 
126147c6ae99SBarry Smith    Input Parameters:
12629a42bb27SBarry Smith +    da - the DM that defines the grid
126347c6ae99SBarry Smith .    vu - Jacobian is computed at this point (ghosted)
126447c6ae99SBarry Smith .    v - product is done on this vector (ghosted)
126547c6ae99SBarry Smith .    fu - output vector = J(vu)*v (not ghosted)
126647c6ae99SBarry Smith -    w - any user data
126747c6ae99SBarry Smith 
126847c6ae99SBarry Smith     Notes: Does NOT do ghost updates on vu and v upon entry
126947c6ae99SBarry Smith 
127047c6ae99SBarry Smith    Level: advanced
127147c6ae99SBarry Smith 
1272837ad101SBarry Smith .seealso: DMDAComputeFunction1()
127347c6ae99SBarry Smith 
127447c6ae99SBarry Smith @*/
12757087cfbeSBarry Smith PetscErrorCode  DMDAMultiplyByJacobian1WithAdifor(DM da,Vec u,Vec v,Vec f,void *w)
127647c6ae99SBarry Smith {
127747c6ae99SBarry Smith   PetscErrorCode ierr;
127847c6ae99SBarry Smith   PetscScalar    *au,*av,*af,*awork;
127947c6ae99SBarry Smith   Vec            work;
1280aa219208SBarry Smith   DMDALocalInfo  info;
128147c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
1282aa219208SBarry Smith   void           (*lf)(DMDALocalInfo*,PetscScalar*,PetscScalar*,PetscScalar*,PetscScalar*,void*,PetscErrorCode*) =
1283aa219208SBarry Smith                  (void (*)(DMDALocalInfo*,PetscScalar*,PetscScalar*,PetscScalar*,PetscScalar*,void*,PetscErrorCode*))*dd->adiformf_lf;
128447c6ae99SBarry Smith 
128547c6ae99SBarry Smith   PetscFunctionBegin;
1286aa219208SBarry Smith   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
128747c6ae99SBarry Smith 
12889a42bb27SBarry Smith   ierr = DMGetGlobalVector(da,&work);CHKERRQ(ierr);
128947c6ae99SBarry Smith   ierr = VecGetArray(u,&au);CHKERRQ(ierr);
129047c6ae99SBarry Smith   ierr = VecGetArray(v,&av);CHKERRQ(ierr);
129147c6ae99SBarry Smith   ierr = VecGetArray(f,&af);CHKERRQ(ierr);
129247c6ae99SBarry Smith   ierr = VecGetArray(work,&awork);CHKERRQ(ierr);
129347c6ae99SBarry Smith   (lf)(&info,au,av,awork,af,w,&ierr);CHKERRQ(ierr);
129447c6ae99SBarry Smith   ierr = VecRestoreArray(u,&au);CHKERRQ(ierr);
129547c6ae99SBarry Smith   ierr = VecRestoreArray(v,&av);CHKERRQ(ierr);
129647c6ae99SBarry Smith   ierr = VecRestoreArray(f,&af);CHKERRQ(ierr);
129747c6ae99SBarry Smith   ierr = VecRestoreArray(work,&awork);CHKERRQ(ierr);
12989a42bb27SBarry Smith   ierr = DMRestoreGlobalVector(da,&work);CHKERRQ(ierr);
129947c6ae99SBarry Smith 
130047c6ae99SBarry Smith   PetscFunctionReturn(0);
130147c6ae99SBarry Smith }
130247c6ae99SBarry Smith 
130347c6ae99SBarry Smith #undef __FUNCT__
13049a42bb27SBarry Smith #define __FUNCT__ "DMSetUp_DA_2D"
13057087cfbeSBarry Smith PetscErrorCode  DMSetUp_DA_2D(DM da)
130647c6ae99SBarry Smith {
130747c6ae99SBarry Smith   DM_DA                  *dd = (DM_DA*)da->data;
130847c6ae99SBarry Smith   const PetscInt         M            = dd->M;
130947c6ae99SBarry Smith   const PetscInt         N            = dd->N;
131047c6ae99SBarry Smith   PetscInt               m            = dd->m;
131147c6ae99SBarry Smith   PetscInt               n            = dd->n;
131247c6ae99SBarry Smith   const PetscInt         dof          = dd->w;
131347c6ae99SBarry Smith   const PetscInt         s            = dd->s;
13141321219cSEthan Coon   const DMDABoundaryType bx         = dd->bx;
13151321219cSEthan Coon   const DMDABoundaryType by         = dd->by;
1316aa219208SBarry Smith   const DMDAStencilType  stencil_type = dd->stencil_type;
131747c6ae99SBarry Smith   PetscInt               *lx           = dd->lx;
131847c6ae99SBarry Smith   PetscInt               *ly           = dd->ly;
131947c6ae99SBarry Smith   MPI_Comm               comm;
132047c6ae99SBarry Smith   PetscMPIInt            rank,size;
1321ce00eea3SSatish Balay   PetscInt               xs,xe,ys,ye,x,y,Xs,Xe,Ys,Ye,start,end,IXs,IXe,IYs,IYe;
1322ce00eea3SSatish Balay   PetscInt               up,down,left,right,i,n0,n1,n2,n3,n5,n6,n7,n8,*idx,nn,*idx_cpy;
1323ce00eea3SSatish Balay   const PetscInt         *idx_full;
1324db87c5ecSEthan Coon   PetscInt               xbase,*bases,*ldims,j,x_t,y_t,s_t,base,count;
132547c6ae99SBarry Smith   PetscInt               s_x,s_y; /* s proportionalized to w */
132647c6ae99SBarry Smith   PetscInt               sn0 = 0,sn2 = 0,sn6 = 0,sn8 = 0;
132747c6ae99SBarry Smith   Vec                    local,global;
132847c6ae99SBarry Smith   VecScatter             ltog,gtol;
1329ce00eea3SSatish Balay   IS                     to,from,ltogis;
133047c6ae99SBarry Smith   PetscErrorCode         ierr;
133147c6ae99SBarry Smith 
133247c6ae99SBarry Smith   PetscFunctionBegin;
13333855c12bSBarry Smith   if (dof < 1) SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Must have 1 or more degrees of freedom per node: %D",dof);
13343855c12bSBarry Smith   if (s < 0) SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Stencil width cannot be negative: %D",s);
133547c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
13363855c12bSBarry Smith #if !defined(PETSC_USE_64BIT_INDICES)
13373855c12bSBarry Smith   if (((Petsc64bitInt) M)*((Petsc64bitInt) N)*((Petsc64bitInt) dof) > (Petsc64bitInt) PETSC_MPI_INT_MAX) SETERRQ3(comm,PETSC_ERR_INT_OVERFLOW,"Mesh of %D by %D by %D (dof) is too large for 32 bit indices",M,N,dof);
13383855c12bSBarry Smith #endif
133947c6ae99SBarry Smith 
134047c6ae99SBarry Smith   if (dof < 1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Must have 1 or more degrees of freedom per node: %D",dof);
134147c6ae99SBarry Smith   if (s < 0) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Stencil width cannot be negative: %D",s);
134247c6ae99SBarry Smith 
134347c6ae99SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
134447c6ae99SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
134547c6ae99SBarry Smith 
134647c6ae99SBarry Smith   dd->dim         = 2;
134747c6ae99SBarry Smith   ierr = PetscMalloc(dof*sizeof(char*),&dd->fieldname);CHKERRQ(ierr);
134847c6ae99SBarry Smith   ierr = PetscMemzero(dd->fieldname,dof*sizeof(char*));CHKERRQ(ierr);
134947c6ae99SBarry Smith 
135047c6ae99SBarry Smith   if (m != PETSC_DECIDE) {
135147c6ae99SBarry Smith     if (m < 1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in X direction: %D",m);
135247c6ae99SBarry Smith     else if (m > size) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in X direction: %D %d",m,size);
135347c6ae99SBarry Smith   }
135447c6ae99SBarry Smith   if (n != PETSC_DECIDE) {
135547c6ae99SBarry Smith     if (n < 1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Y direction: %D",n);
135647c6ae99SBarry Smith     else if (n > size) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Y direction: %D %d",n,size);
135747c6ae99SBarry Smith   }
135847c6ae99SBarry Smith 
135947c6ae99SBarry Smith   if (m == PETSC_DECIDE || n == PETSC_DECIDE) {
136047c6ae99SBarry Smith     if (n != PETSC_DECIDE) {
136147c6ae99SBarry Smith       m = size/n;
136247c6ae99SBarry Smith     } else if (m != PETSC_DECIDE) {
136347c6ae99SBarry Smith       n = size/m;
136447c6ae99SBarry Smith     } else {
136547c6ae99SBarry Smith       /* try for squarish distribution */
136647c6ae99SBarry Smith       m = (PetscInt)(0.5 + sqrt(((double)M)*((double)size)/((double)N)));
136747c6ae99SBarry Smith       if (!m) m = 1;
136847c6ae99SBarry Smith       while (m > 0) {
136947c6ae99SBarry Smith 	n = size/m;
137047c6ae99SBarry Smith 	if (m*n == size) break;
137147c6ae99SBarry Smith 	m--;
137247c6ae99SBarry Smith       }
137347c6ae99SBarry Smith       if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;}
137447c6ae99SBarry Smith     }
137547c6ae99SBarry 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 ");
137647c6ae99SBarry Smith   } else if (m*n != size) SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"Given Bad partition");
137747c6ae99SBarry Smith 
137847c6ae99SBarry Smith   if (M < m) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in x direction is too fine! %D %D",M,m);
137947c6ae99SBarry Smith   if (N < n) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in y direction is too fine! %D %D",N,n);
138047c6ae99SBarry Smith 
138147c6ae99SBarry Smith   /*
138247c6ae99SBarry Smith      Determine locally owned region
138347c6ae99SBarry Smith      xs is the first local node number, x is the number of local nodes
138447c6ae99SBarry Smith   */
138547c6ae99SBarry Smith   if (!lx) {
138647c6ae99SBarry Smith     ierr = PetscMalloc(m*sizeof(PetscInt), &dd->lx);CHKERRQ(ierr);
138747c6ae99SBarry Smith     lx = dd->lx;
138847c6ae99SBarry Smith     for (i=0; i<m; i++) {
138947c6ae99SBarry Smith       lx[i] = M/m + ((M % m) > i);
139047c6ae99SBarry Smith     }
139147c6ae99SBarry Smith   }
139247c6ae99SBarry Smith   x  = lx[rank % m];
139347c6ae99SBarry Smith   xs = 0;
139447c6ae99SBarry Smith   for (i=0; i<(rank % m); i++) {
139547c6ae99SBarry Smith     xs += lx[i];
139647c6ae99SBarry Smith   }
139747c6ae99SBarry Smith #if defined(PETSC_USE_DEBUG)
139847c6ae99SBarry Smith   left = xs;
139947c6ae99SBarry Smith   for (i=(rank % m); i<m; i++) {
140047c6ae99SBarry Smith     left += lx[i];
140147c6ae99SBarry Smith   }
140247c6ae99SBarry 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);
140347c6ae99SBarry Smith #endif
140447c6ae99SBarry Smith 
140547c6ae99SBarry Smith   /*
140647c6ae99SBarry Smith      Determine locally owned region
140747c6ae99SBarry Smith      ys is the first local node number, y is the number of local nodes
140847c6ae99SBarry Smith   */
140947c6ae99SBarry Smith   if (!ly) {
141047c6ae99SBarry Smith     ierr = PetscMalloc(n*sizeof(PetscInt), &dd->ly);CHKERRQ(ierr);
141147c6ae99SBarry Smith     ly = dd->ly;
141247c6ae99SBarry Smith     for (i=0; i<n; i++) {
141347c6ae99SBarry Smith       ly[i] = N/n + ((N % n) > i);
141447c6ae99SBarry Smith     }
141547c6ae99SBarry Smith   }
141647c6ae99SBarry Smith   y  = ly[rank/m];
141747c6ae99SBarry Smith   ys = 0;
141847c6ae99SBarry Smith   for (i=0; i<(rank/m); i++) {
141947c6ae99SBarry Smith     ys += ly[i];
142047c6ae99SBarry Smith   }
142147c6ae99SBarry Smith #if defined(PETSC_USE_DEBUG)
142247c6ae99SBarry Smith   left = ys;
142347c6ae99SBarry Smith   for (i=(rank/m); i<n; i++) {
142447c6ae99SBarry Smith     left += ly[i];
142547c6ae99SBarry Smith   }
142647c6ae99SBarry 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);
142747c6ae99SBarry Smith #endif
142847c6ae99SBarry Smith 
1429bcea557cSEthan Coon   /*
1430bcea557cSEthan Coon    check if the scatter requires more than one process neighbor or wraps around
1431bcea557cSEthan Coon    the domain more than once
1432bcea557cSEthan Coon   */
1433bcea557cSEthan Coon   if ((x < s) && ((m > 1) || (bx == DMDA_BOUNDARY_PERIODIC))) {
1434bcea557cSEthan Coon     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local x-width of domain x %D is smaller than stencil width s %D",x,s);
1435bcea557cSEthan Coon   }
1436bcea557cSEthan Coon   if ((y < s) && ((n > 1) || (by == DMDA_BOUNDARY_PERIODIC))) {
1437bcea557cSEthan Coon     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local y-width of domain y %D is smaller than stencil width s %D",y,s);
1438bcea557cSEthan Coon   }
143947c6ae99SBarry Smith   xe = xs + x;
144047c6ae99SBarry Smith   ye = ys + y;
144147c6ae99SBarry Smith 
1442ce00eea3SSatish Balay   /* determine ghost region (Xs) and region scattered into (IXs)  */
144347c6ae99SBarry Smith   /* Assume No Periodicity */
1444ce00eea3SSatish Balay   if (xs-s > 0) { Xs = xs - s; IXs = xs - s; } else { Xs = 0; IXs = 0; }
1445ce00eea3SSatish Balay   if (xe+s <= M) { Xe = xe + s; IXe = xe + s; } else { Xe = M; IXe = M; }
1446ce00eea3SSatish Balay   if (ys-s > 0) { Ys = ys - s; IYs = ys - s; } else { Ys = 0; IYs = 0; }
1447ce00eea3SSatish Balay   if (ye+s <= N) { Ye = ye + s; IYe = ye + s; } else { Ye = N; IYe = N; }
144847c6ae99SBarry Smith 
1449ce00eea3SSatish Balay   /* fix for periodicity/ghosted */
14501321219cSEthan Coon   if (bx) { Xs = xs - s; Xe = xe + s; }
14511321219cSEthan Coon   if (bx == DMDA_BOUNDARY_PERIODIC) { IXs = xs - s; IXe = xe + s; }
14521321219cSEthan Coon   if (by) { Ys = ys - s; Ye = ye + s; }
14531321219cSEthan Coon   if (by == DMDA_BOUNDARY_PERIODIC) { IYs = ys - s; IYe = ye + s; }
145447c6ae99SBarry Smith 
145547c6ae99SBarry Smith   /* Resize all X parameters to reflect w */
1456ce00eea3SSatish Balay   s_x = s;
145747c6ae99SBarry Smith   s_y = s;
145847c6ae99SBarry Smith 
145947c6ae99SBarry Smith   /* determine starting point of each processor */
146047c6ae99SBarry Smith   nn    = x*y;
146147c6ae99SBarry Smith   ierr  = PetscMalloc2(size+1,PetscInt,&bases,size,PetscInt,&ldims);CHKERRQ(ierr);
146247c6ae99SBarry Smith   ierr  = MPI_Allgather(&nn,1,MPIU_INT,ldims,1,MPIU_INT,comm);CHKERRQ(ierr);
146347c6ae99SBarry Smith   bases[0] = 0;
146447c6ae99SBarry Smith   for (i=1; i<=size; i++) {
146547c6ae99SBarry Smith     bases[i] = ldims[i-1];
146647c6ae99SBarry Smith   }
146747c6ae99SBarry Smith   for (i=1; i<=size; i++) {
146847c6ae99SBarry Smith     bases[i] += bases[i-1];
146947c6ae99SBarry Smith   }
1470ce00eea3SSatish Balay   base = bases[rank]*dof;
147147c6ae99SBarry Smith 
147247c6ae99SBarry Smith   /* allocate the base parallel and sequential vectors */
1473ce00eea3SSatish Balay   dd->Nlocal = x*y*dof;
147447c6ae99SBarry Smith   ierr = VecCreateMPIWithArray(comm,dd->Nlocal,PETSC_DECIDE,0,&global);CHKERRQ(ierr);
147547c6ae99SBarry Smith   ierr = VecSetBlockSize(global,dof);CHKERRQ(ierr);
1476ce00eea3SSatish Balay   dd->nlocal = (Xe-Xs)*(Ye-Ys)*dof;
147747c6ae99SBarry Smith   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,dd->nlocal,0,&local);CHKERRQ(ierr);
147847c6ae99SBarry Smith   ierr = VecSetBlockSize(local,dof);CHKERRQ(ierr);
147947c6ae99SBarry Smith 
148047c6ae99SBarry Smith   /* generate appropriate vector scatters */
148147c6ae99SBarry Smith   /* local to global inserts non-ghost point region into global */
148247c6ae99SBarry Smith   ierr = VecGetOwnershipRange(global,&start,&end);CHKERRQ(ierr);
1483ce00eea3SSatish Balay   ierr = ISCreateStride(comm,x*y*dof,start,1,&to);CHKERRQ(ierr);
148447c6ae99SBarry Smith 
1485db87c5ecSEthan Coon   count = x*y;
1486ce00eea3SSatish Balay   ierr = PetscMalloc(x*y*sizeof(PetscInt),&idx);CHKERRQ(ierr);
1487ce00eea3SSatish Balay   left = xs - Xs; right = left + x;
1488ce00eea3SSatish Balay   down = ys - Ys; up = down + y;
148947c6ae99SBarry Smith   count = 0;
149047c6ae99SBarry Smith   for (i=down; i<up; i++) {
1491ce00eea3SSatish Balay     for (j=left; j<right; j++) {
1492ce00eea3SSatish Balay       idx[count++] = i*(Xe-Xs) + j;
149347c6ae99SBarry Smith     }
149447c6ae99SBarry Smith   }
149547c6ae99SBarry Smith 
1496ce00eea3SSatish Balay   ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&from);CHKERRQ(ierr);
149747c6ae99SBarry Smith   ierr = VecScatterCreate(local,from,global,to,&ltog);CHKERRQ(ierr);
149847c6ae99SBarry Smith   ierr = PetscLogObjectParent(dd,ltog);CHKERRQ(ierr);
1499fcfd50ebSBarry Smith   ierr = ISDestroy(&from);CHKERRQ(ierr);
1500fcfd50ebSBarry Smith   ierr = ISDestroy(&to);CHKERRQ(ierr);
150147c6ae99SBarry Smith 
1502ce00eea3SSatish Balay   /* global to local must include ghost points within the domain,
1503ce00eea3SSatish Balay      but not ghost points outside the domain that aren't periodic */
1504aa219208SBarry Smith   if (stencil_type == DMDA_STENCIL_BOX) {
1505db87c5ecSEthan Coon     count = (IXe-IXs)*(IYe-IYs);
1506db87c5ecSEthan Coon     ierr  = PetscMalloc(count*sizeof(PetscInt),&idx);CHKERRQ(ierr);
1507ce00eea3SSatish Balay 
1508ce00eea3SSatish Balay     left = IXs - Xs; right = left + (IXe-IXs);
1509ce00eea3SSatish Balay     down = IYs - Ys; up = down + (IYe-IYs);
1510ce00eea3SSatish Balay     count = 0;
1511ce00eea3SSatish Balay     for (i=down; i<up; i++) {
1512ce00eea3SSatish Balay       for (j=left; j<right; j++) {
1513ce00eea3SSatish Balay         idx[count++] = j + i*(Xe-Xs);
1514ce00eea3SSatish Balay       }
1515ce00eea3SSatish Balay     }
1516ce00eea3SSatish Balay     ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr);
1517ce00eea3SSatish Balay 
151847c6ae99SBarry Smith   } else {
151947c6ae99SBarry Smith     /* must drop into cross shape region */
152047c6ae99SBarry Smith     /*       ---------|
152147c6ae99SBarry Smith             |  top    |
1522ce00eea3SSatish Balay          |---         ---| up
152347c6ae99SBarry Smith          |   middle      |
152447c6ae99SBarry Smith          |               |
1525ce00eea3SSatish Balay          ----         ---- down
152647c6ae99SBarry Smith             | bottom  |
152747c6ae99SBarry Smith             -----------
152847c6ae99SBarry Smith          Xs xs        xe Xe */
1529db87c5ecSEthan Coon     count = (ys-IYs)*x + y*(IXe-IXs) + (IYe-ye)*x;
1530db87c5ecSEthan Coon     ierr  = PetscMalloc(count*sizeof(PetscInt),&idx);CHKERRQ(ierr);
1531ce00eea3SSatish Balay 
1532ce00eea3SSatish Balay     left = xs - Xs; right = left + x;
1533ce00eea3SSatish Balay     down = ys - Ys; up = down + y;
153447c6ae99SBarry Smith     count = 0;
1535ce00eea3SSatish Balay     /* bottom */
1536ce00eea3SSatish Balay     for (i=(IYs-Ys); i<down; i++) {
1537ce00eea3SSatish Balay       for (j=left; j<right; j++) {
1538ce00eea3SSatish Balay         idx[count++] = j + i*(Xe-Xs);
153947c6ae99SBarry Smith       }
154047c6ae99SBarry Smith     }
154147c6ae99SBarry Smith     /* middle */
154247c6ae99SBarry Smith     for (i=down; i<up; i++) {
1543ce00eea3SSatish Balay       for (j=(IXs-Xs); j<(IXe-Xs); j++) {
1544ce00eea3SSatish Balay         idx[count++] = j + i*(Xe-Xs);
154547c6ae99SBarry Smith       }
154647c6ae99SBarry Smith     }
154747c6ae99SBarry Smith     /* top */
1548ce00eea3SSatish Balay     for (i=up; i<up+IYe-ye; i++) {
1549ce00eea3SSatish Balay       for (j=left; j<right; j++) {
1550ce00eea3SSatish Balay         idx[count++] = j + i*(Xe-Xs);
155147c6ae99SBarry Smith       }
155247c6ae99SBarry Smith     }
155347c6ae99SBarry Smith     ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr);
155447c6ae99SBarry Smith   }
155547c6ae99SBarry Smith 
155647c6ae99SBarry Smith 
155747c6ae99SBarry Smith   /* determine who lies on each side of us stored in    n6 n7 n8
155847c6ae99SBarry Smith                                                         n3    n5
155947c6ae99SBarry Smith                                                         n0 n1 n2
156047c6ae99SBarry Smith   */
156147c6ae99SBarry Smith 
156247c6ae99SBarry Smith   /* Assume the Non-Periodic Case */
156347c6ae99SBarry Smith   n1 = rank - m;
156447c6ae99SBarry Smith   if (rank % m) {
156547c6ae99SBarry Smith     n0 = n1 - 1;
156647c6ae99SBarry Smith   } else {
156747c6ae99SBarry Smith     n0 = -1;
156847c6ae99SBarry Smith   }
156947c6ae99SBarry Smith   if ((rank+1) % m) {
157047c6ae99SBarry Smith     n2 = n1 + 1;
157147c6ae99SBarry Smith     n5 = rank + 1;
157247c6ae99SBarry Smith     n8 = rank + m + 1; if (n8 >= m*n) n8 = -1;
157347c6ae99SBarry Smith   } else {
157447c6ae99SBarry Smith     n2 = -1; n5 = -1; n8 = -1;
157547c6ae99SBarry Smith   }
157647c6ae99SBarry Smith   if (rank % m) {
157747c6ae99SBarry Smith     n3 = rank - 1;
157847c6ae99SBarry Smith     n6 = n3 + m; if (n6 >= m*n) n6 = -1;
157947c6ae99SBarry Smith   } else {
158047c6ae99SBarry Smith     n3 = -1; n6 = -1;
158147c6ae99SBarry Smith   }
158247c6ae99SBarry Smith   n7 = rank + m; if (n7 >= m*n) n7 = -1;
158347c6ae99SBarry Smith 
15841321219cSEthan Coon   if (bx == DMDA_BOUNDARY_PERIODIC && by == DMDA_BOUNDARY_PERIODIC) {
158547c6ae99SBarry Smith   /* Modify for Periodic Cases */
158647c6ae99SBarry Smith     /* Handle all four corners */
158747c6ae99SBarry Smith     if ((n6 < 0) && (n7 < 0) && (n3 < 0)) n6 = m-1;
158847c6ae99SBarry Smith     if ((n8 < 0) && (n7 < 0) && (n5 < 0)) n8 = 0;
158947c6ae99SBarry Smith     if ((n2 < 0) && (n5 < 0) && (n1 < 0)) n2 = size-m;
159047c6ae99SBarry Smith     if ((n0 < 0) && (n3 < 0) && (n1 < 0)) n0 = size-1;
159147c6ae99SBarry Smith 
159247c6ae99SBarry Smith     /* Handle Top and Bottom Sides */
159347c6ae99SBarry Smith     if (n1 < 0) n1 = rank + m * (n-1);
159447c6ae99SBarry Smith     if (n7 < 0) n7 = rank - m * (n-1);
159547c6ae99SBarry Smith     if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1;
159647c6ae99SBarry Smith     if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1;
159747c6ae99SBarry Smith     if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1;
159847c6ae99SBarry Smith     if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1;
159947c6ae99SBarry Smith 
160047c6ae99SBarry Smith     /* Handle Left and Right Sides */
160147c6ae99SBarry Smith     if (n3 < 0) n3 = rank + (m-1);
160247c6ae99SBarry Smith     if (n5 < 0) n5 = rank - (m-1);
160347c6ae99SBarry Smith     if ((n1 >= 0) && (n0 < 0)) n0 = rank-1;
160447c6ae99SBarry Smith     if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1;
160547c6ae99SBarry Smith     if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1;
160647c6ae99SBarry Smith     if ((n7 >= 0) && (n8 < 0)) n8 = rank+1;
16071321219cSEthan Coon   } else if (by == DMDA_BOUNDARY_PERIODIC) {  /* Handle Top and Bottom Sides */
1608ce00eea3SSatish Balay     if (n1 < 0) n1 = rank + m * (n-1);
1609ce00eea3SSatish Balay     if (n7 < 0) n7 = rank - m * (n-1);
1610ce00eea3SSatish Balay     if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1;
1611ce00eea3SSatish Balay     if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1;
1612ce00eea3SSatish Balay     if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1;
1613ce00eea3SSatish Balay     if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1;
16141321219cSEthan Coon   } else if (bx == DMDA_BOUNDARY_PERIODIC) { /* Handle Left and Right Sides */
1615ce00eea3SSatish Balay     if (n3 < 0) n3 = rank + (m-1);
1616ce00eea3SSatish Balay     if (n5 < 0) n5 = rank - (m-1);
1617ce00eea3SSatish Balay     if ((n1 >= 0) && (n0 < 0)) n0 = rank-1;
1618ce00eea3SSatish Balay     if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1;
1619ce00eea3SSatish Balay     if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1;
1620ce00eea3SSatish Balay     if ((n7 >= 0) && (n8 < 0)) n8 = rank+1;
162147c6ae99SBarry Smith   }
1622ce00eea3SSatish Balay 
162347c6ae99SBarry Smith   ierr = PetscMalloc(9*sizeof(PetscInt),&dd->neighbors);CHKERRQ(ierr);
162447c6ae99SBarry Smith   dd->neighbors[0] = n0;
162547c6ae99SBarry Smith   dd->neighbors[1] = n1;
162647c6ae99SBarry Smith   dd->neighbors[2] = n2;
162747c6ae99SBarry Smith   dd->neighbors[3] = n3;
162847c6ae99SBarry Smith   dd->neighbors[4] = rank;
162947c6ae99SBarry Smith   dd->neighbors[5] = n5;
163047c6ae99SBarry Smith   dd->neighbors[6] = n6;
163147c6ae99SBarry Smith   dd->neighbors[7] = n7;
163247c6ae99SBarry Smith   dd->neighbors[8] = n8;
163347c6ae99SBarry Smith 
1634aa219208SBarry Smith   if (stencil_type == DMDA_STENCIL_STAR) {
163547c6ae99SBarry Smith     /* save corner processor numbers */
163647c6ae99SBarry Smith     sn0 = n0; sn2 = n2; sn6 = n6; sn8 = n8;
163747c6ae99SBarry Smith     n0 = n2 = n6 = n8 = -1;
163847c6ae99SBarry Smith   }
163947c6ae99SBarry Smith 
1640ce00eea3SSatish Balay   ierr = PetscMalloc((Xe-Xs)*(Ye-Ys)*sizeof(PetscInt),&idx);CHKERRQ(ierr);
1641ce00eea3SSatish Balay   ierr = PetscLogObjectMemory(da,(Xe-Xs)*(Ye-Ys)*sizeof(PetscInt));CHKERRQ(ierr);
164247c6ae99SBarry Smith 
1643ce00eea3SSatish Balay   nn = 0;
164447c6ae99SBarry Smith   xbase = bases[rank];
164547c6ae99SBarry Smith   for (i=1; i<=s_y; i++) {
164647c6ae99SBarry Smith     if (n0 >= 0) { /* left below */
1647ce00eea3SSatish Balay       x_t = lx[n0 % m];
164847c6ae99SBarry Smith       y_t = ly[(n0/m)];
164947c6ae99SBarry Smith       s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x;
165047c6ae99SBarry Smith       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
165147c6ae99SBarry Smith     }
165247c6ae99SBarry Smith     if (n1 >= 0) { /* directly below */
165347c6ae99SBarry Smith       x_t = x;
165447c6ae99SBarry Smith       y_t = ly[(n1/m)];
165547c6ae99SBarry Smith       s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t;
165647c6ae99SBarry Smith       for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
165747c6ae99SBarry Smith     }
165847c6ae99SBarry Smith     if (n2 >= 0) { /* right below */
1659ce00eea3SSatish Balay       x_t = lx[n2 % m];
166047c6ae99SBarry Smith       y_t = ly[(n2/m)];
166147c6ae99SBarry Smith       s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t;
166247c6ae99SBarry Smith       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
166347c6ae99SBarry Smith     }
166447c6ae99SBarry Smith   }
166547c6ae99SBarry Smith 
166647c6ae99SBarry Smith   for (i=0; i<y; i++) {
166747c6ae99SBarry Smith     if (n3 >= 0) { /* directly left */
1668ce00eea3SSatish Balay       x_t = lx[n3 % m];
166947c6ae99SBarry Smith       /* y_t = y; */
167047c6ae99SBarry Smith       s_t = bases[n3] + (i+1)*x_t - s_x;
167147c6ae99SBarry Smith       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
167247c6ae99SBarry Smith     }
167347c6ae99SBarry Smith 
167447c6ae99SBarry Smith     for (j=0; j<x; j++) { idx[nn++] = xbase++; } /* interior */
167547c6ae99SBarry Smith 
167647c6ae99SBarry Smith     if (n5 >= 0) { /* directly right */
1677ce00eea3SSatish Balay       x_t = lx[n5 % m];
167847c6ae99SBarry Smith       /* y_t = y; */
167947c6ae99SBarry Smith       s_t = bases[n5] + (i)*x_t;
168047c6ae99SBarry Smith       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
168147c6ae99SBarry Smith     }
168247c6ae99SBarry Smith   }
168347c6ae99SBarry Smith 
168447c6ae99SBarry Smith   for (i=1; i<=s_y; i++) {
168547c6ae99SBarry Smith     if (n6 >= 0) { /* left above */
1686ce00eea3SSatish Balay       x_t = lx[n6 % m];
168747c6ae99SBarry Smith       /* y_t = ly[(n6/m)]; */
168847c6ae99SBarry Smith       s_t = bases[n6] + (i)*x_t - s_x;
168947c6ae99SBarry Smith       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
169047c6ae99SBarry Smith     }
169147c6ae99SBarry Smith     if (n7 >= 0) { /* directly above */
169247c6ae99SBarry Smith       x_t = x;
169347c6ae99SBarry Smith       /* y_t = ly[(n7/m)]; */
169447c6ae99SBarry Smith       s_t = bases[n7] + (i-1)*x_t;
169547c6ae99SBarry Smith       for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
169647c6ae99SBarry Smith     }
169747c6ae99SBarry Smith     if (n8 >= 0) { /* right above */
1698ce00eea3SSatish Balay       x_t = lx[n8 % m];
169947c6ae99SBarry Smith       /* y_t = ly[(n8/m)]; */
170047c6ae99SBarry Smith       s_t = bases[n8] + (i-1)*x_t;
170147c6ae99SBarry Smith       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
170247c6ae99SBarry Smith     }
170347c6ae99SBarry Smith   }
170447c6ae99SBarry Smith 
1705ce00eea3SSatish Balay   ierr = ISCreateBlock(comm,dof,nn,idx,PETSC_COPY_VALUES,&from);CHKERRQ(ierr);
170647c6ae99SBarry Smith   ierr = VecScatterCreate(global,from,local,to,&gtol);CHKERRQ(ierr);
170747c6ae99SBarry Smith   ierr = PetscLogObjectParent(da,gtol);CHKERRQ(ierr);
1708fcfd50ebSBarry Smith   ierr = ISDestroy(&to);CHKERRQ(ierr);
1709fcfd50ebSBarry Smith   ierr = ISDestroy(&from);CHKERRQ(ierr);
171047c6ae99SBarry Smith 
1711aa219208SBarry Smith   if (stencil_type == DMDA_STENCIL_STAR) {
1712ce00eea3SSatish Balay     n0 = sn0; n2 = sn2; n6 = sn6; n8 = sn8;
1713ce00eea3SSatish Balay   }
1714ce00eea3SSatish Balay 
1715ce00eea3SSatish Balay   if ((stencil_type == DMDA_STENCIL_STAR) ||
17161321219cSEthan Coon       (bx && bx != DMDA_BOUNDARY_PERIODIC) ||
17171321219cSEthan Coon       (by && by != DMDA_BOUNDARY_PERIODIC)) {
171847c6ae99SBarry Smith     /*
171947c6ae99SBarry Smith         Recompute the local to global mappings, this time keeping the
1720ce00eea3SSatish Balay       information about the cross corner processor numbers and any ghosted
1721ce00eea3SSatish Balay       but not periodic indices.
172247c6ae99SBarry Smith     */
172347c6ae99SBarry Smith     nn = 0;
172447c6ae99SBarry Smith     xbase = bases[rank];
172547c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
172647c6ae99SBarry Smith       if (n0 >= 0) { /* left below */
1727ce00eea3SSatish Balay         x_t = lx[n0 % m];
172847c6ae99SBarry Smith         y_t = ly[(n0/m)];
172947c6ae99SBarry Smith         s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x;
173047c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1731ce00eea3SSatish Balay       } else if (xs-Xs > 0 && ys-Ys > 0) {
1732ce00eea3SSatish Balay         for (j=0; j<s_x; j++) { idx[nn++] = -1;}
173347c6ae99SBarry Smith       }
173447c6ae99SBarry Smith       if (n1 >= 0) { /* directly below */
173547c6ae99SBarry Smith         x_t = x;
173647c6ae99SBarry Smith         y_t = ly[(n1/m)];
173747c6ae99SBarry Smith         s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t;
173847c6ae99SBarry Smith         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1739ce00eea3SSatish Balay       } else if (ys-Ys > 0) {
1740ce00eea3SSatish Balay         for (j=0; j<x; j++) { idx[nn++] = -1;}
174147c6ae99SBarry Smith       }
174247c6ae99SBarry Smith       if (n2 >= 0) { /* right below */
1743ce00eea3SSatish Balay         x_t = lx[n2 % m];
174447c6ae99SBarry Smith         y_t = ly[(n2/m)];
174547c6ae99SBarry Smith         s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t;
174647c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1747ce00eea3SSatish Balay       } else if (Xe-xe> 0 && ys-Ys > 0) {
1748ce00eea3SSatish Balay         for (j=0; j<s_x; j++) { idx[nn++] = -1;}
174947c6ae99SBarry Smith       }
175047c6ae99SBarry Smith     }
175147c6ae99SBarry Smith 
175247c6ae99SBarry Smith     for (i=0; i<y; i++) {
175347c6ae99SBarry Smith       if (n3 >= 0) { /* directly left */
1754ce00eea3SSatish Balay         x_t = lx[n3 % m];
175547c6ae99SBarry Smith         /* y_t = y; */
175647c6ae99SBarry Smith         s_t = bases[n3] + (i+1)*x_t - s_x;
175747c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1758ce00eea3SSatish Balay       } else if (xs-Xs > 0) {
1759ce00eea3SSatish Balay         for (j=0; j<s_x; j++) { idx[nn++] = -1;}
176047c6ae99SBarry Smith       }
176147c6ae99SBarry Smith 
176247c6ae99SBarry Smith       for (j=0; j<x; j++) { idx[nn++] = xbase++; } /* interior */
176347c6ae99SBarry Smith 
176447c6ae99SBarry Smith       if (n5 >= 0) { /* directly right */
1765ce00eea3SSatish Balay         x_t = lx[n5 % m];
176647c6ae99SBarry Smith         /* y_t = y; */
176747c6ae99SBarry Smith         s_t = bases[n5] + (i)*x_t;
176847c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1769ce00eea3SSatish Balay       } else if (Xe-xe > 0) {
1770ce00eea3SSatish Balay         for (j=0; j<s_x; j++) { idx[nn++] = -1;}
177147c6ae99SBarry Smith       }
177247c6ae99SBarry Smith     }
177347c6ae99SBarry Smith 
177447c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
177547c6ae99SBarry Smith       if (n6 >= 0) { /* left above */
1776ce00eea3SSatish Balay         x_t = lx[n6 % m];
177747c6ae99SBarry Smith         /* y_t = ly[(n6/m)]; */
177847c6ae99SBarry Smith         s_t = bases[n6] + (i)*x_t - s_x;
177947c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1780ce00eea3SSatish Balay       } else if (xs-Xs > 0 && Ye-ye > 0) {
1781ce00eea3SSatish Balay         for (j=0; j<s_x; j++) { idx[nn++] = -1;}
178247c6ae99SBarry Smith       }
178347c6ae99SBarry Smith       if (n7 >= 0) { /* directly above */
178447c6ae99SBarry Smith         x_t = x;
178547c6ae99SBarry Smith         /* y_t = ly[(n7/m)]; */
178647c6ae99SBarry Smith         s_t = bases[n7] + (i-1)*x_t;
178747c6ae99SBarry Smith         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1788ce00eea3SSatish Balay       } else if (Ye-ye > 0) {
1789ce00eea3SSatish Balay         for (j=0; j<x; j++) { idx[nn++] = -1;}
179047c6ae99SBarry Smith       }
179147c6ae99SBarry Smith       if (n8 >= 0) { /* right above */
1792ce00eea3SSatish Balay         x_t = lx[n8 % m];
179347c6ae99SBarry Smith         /* y_t = ly[(n8/m)]; */
179447c6ae99SBarry Smith         s_t = bases[n8] + (i-1)*x_t;
179547c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1796ce00eea3SSatish Balay       } else if (Xe-xe > 0 && Ye-ye > 0) {
1797ce00eea3SSatish Balay         for (j=0; j<s_x; j++) { idx[nn++] = -1;}
179847c6ae99SBarry Smith       }
179947c6ae99SBarry Smith     }
180047c6ae99SBarry Smith   }
1801ce00eea3SSatish Balay   /*
1802ce00eea3SSatish Balay      Set the local to global ordering in the global vector, this allows use
1803ce00eea3SSatish Balay      of VecSetValuesLocal().
1804ce00eea3SSatish Balay   */
1805ce00eea3SSatish Balay   ierr = ISCreateBlock(comm,dof,nn,idx,PETSC_OWN_POINTER,&ltogis);CHKERRQ(ierr);
1806ce00eea3SSatish Balay   ierr = PetscMalloc(nn*dof*sizeof(PetscInt),&idx_cpy);CHKERRQ(ierr);
1807db87c5ecSEthan Coon   ierr = PetscLogObjectMemory(da,nn*dof*sizeof(PetscInt));CHKERRQ(ierr);
1808ce00eea3SSatish Balay   ierr = ISGetIndices(ltogis, &idx_full);
1809ce00eea3SSatish Balay   ierr = PetscMemcpy(idx_cpy,idx_full,nn*dof*sizeof(PetscInt));CHKERRQ(ierr);
1810ce00eea3SSatish Balay   ierr = ISRestoreIndices(ltogis, &idx_full);
1811ce00eea3SSatish Balay   ierr = ISLocalToGlobalMappingCreateIS(ltogis,&da->ltogmap);CHKERRQ(ierr);
1812ce00eea3SSatish Balay   ierr = PetscLogObjectParent(da,da->ltogmap);CHKERRQ(ierr);
1813fcfd50ebSBarry Smith   ierr = ISDestroy(&ltogis);CHKERRQ(ierr);
1814ce00eea3SSatish Balay   ierr = ISLocalToGlobalMappingBlock(da->ltogmap,dd->w,&da->ltogmapb);CHKERRQ(ierr);
1815ce00eea3SSatish Balay   ierr = PetscLogObjectParent(da,da->ltogmap);CHKERRQ(ierr);
181647c6ae99SBarry Smith 
1817ce00eea3SSatish Balay   ierr = PetscFree2(bases,ldims);CHKERRQ(ierr);
181847c6ae99SBarry Smith   dd->m  = m;  dd->n  = n;
1819ce00eea3SSatish Balay   /* note petsc expects xs/xe/Xs/Xe to be multiplied by #dofs in many places */
1820ce00eea3SSatish Balay   dd->xs = xs*dof; dd->xe = xe*dof; dd->ys = ys; dd->ye = ye; dd->zs = 0; dd->ze = 1;
1821ce00eea3SSatish Balay   dd->Xs = Xs*dof; dd->Xe = Xe*dof; dd->Ys = Ys; dd->Ye = Ye; dd->Zs = 0; dd->Ze = 1;
182247c6ae99SBarry Smith 
1823fcfd50ebSBarry Smith   ierr = VecDestroy(&local);CHKERRQ(ierr);
1824fcfd50ebSBarry Smith   ierr = VecDestroy(&global);CHKERRQ(ierr);
182547c6ae99SBarry Smith 
182647c6ae99SBarry Smith   dd->gtol      = gtol;
182747c6ae99SBarry Smith   dd->ltog      = ltog;
1828ce00eea3SSatish Balay   dd->idx       = idx_cpy;
1829ce00eea3SSatish Balay   dd->Nl        = nn*dof;
183047c6ae99SBarry Smith   dd->base      = base;
18319a42bb27SBarry Smith   da->ops->view = DMView_DA_2d;
183247c6ae99SBarry Smith   dd->ltol = PETSC_NULL;
183347c6ae99SBarry Smith   dd->ao   = PETSC_NULL;
183447c6ae99SBarry Smith 
183547c6ae99SBarry Smith   PetscFunctionReturn(0);
183647c6ae99SBarry Smith }
183747c6ae99SBarry Smith 
183847c6ae99SBarry Smith #undef __FUNCT__
1839aa219208SBarry Smith #define __FUNCT__ "DMDACreate2d"
184047c6ae99SBarry Smith /*@C
1841aa219208SBarry Smith    DMDACreate2d -  Creates an object that will manage the communication of  two-dimensional
184247c6ae99SBarry Smith    regular array data that is distributed across some processors.
184347c6ae99SBarry Smith 
184447c6ae99SBarry Smith    Collective on MPI_Comm
184547c6ae99SBarry Smith 
184647c6ae99SBarry Smith    Input Parameters:
184747c6ae99SBarry Smith +  comm - MPI communicator
18481321219cSEthan Coon .  bx,by - type of ghost nodes the array have.
18491321219cSEthan Coon          Use one of DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_GHOSTED, DMDA_BOUNDARY_PERIODIC.
1850aa219208SBarry Smith .  stencil_type - stencil type.  Use either DMDA_STENCIL_BOX or DMDA_STENCIL_STAR.
185147c6ae99SBarry 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
185247c6ae99SBarry Smith             from the command line with -da_grid_x <M> -da_grid_y <N>)
185347c6ae99SBarry Smith .  m,n - corresponding number of processors in each dimension
185447c6ae99SBarry Smith          (or PETSC_DECIDE to have calculated)
185547c6ae99SBarry Smith .  dof - number of degrees of freedom per node
185647c6ae99SBarry Smith .  s - stencil width
185747c6ae99SBarry Smith -  lx, ly - arrays containing the number of nodes in each cell along
185847c6ae99SBarry Smith            the x and y coordinates, or PETSC_NULL. If non-null, these
185947c6ae99SBarry Smith            must be of length as m and n, and the corresponding
186047c6ae99SBarry Smith            m and n cannot be PETSC_DECIDE. The sum of the lx[] entries
186147c6ae99SBarry Smith            must be M, and the sum of the ly[] entries must be N.
186247c6ae99SBarry Smith 
186347c6ae99SBarry Smith    Output Parameter:
186447c6ae99SBarry Smith .  da - the resulting distributed array object
186547c6ae99SBarry Smith 
186647c6ae99SBarry Smith    Options Database Key:
1867aa219208SBarry Smith +  -da_view - Calls DMView() at the conclusion of DMDACreate2d()
186847c6ae99SBarry Smith .  -da_grid_x <nx> - number of grid points in x direction, if M < 0
186947c6ae99SBarry Smith .  -da_grid_y <ny> - number of grid points in y direction, if N < 0
187047c6ae99SBarry Smith .  -da_processors_x <nx> - number of processors in x direction
187147c6ae99SBarry Smith .  -da_processors_y <ny> - number of processors in y direction
1872e0f5d30fSBarry Smith .  -da_refine_x <rx> - refinement ratio in x direction
1873e0f5d30fSBarry Smith .  -da_refine_y <ry> - refinement ratio in y direction
1874e0f5d30fSBarry Smith -  -da_refine <n> - refine the DMDA n times before creating, if M or N < 0
1875e0f5d30fSBarry Smith 
187647c6ae99SBarry Smith 
187747c6ae99SBarry Smith    Level: beginner
187847c6ae99SBarry Smith 
187947c6ae99SBarry Smith    Notes:
1880aa219208SBarry Smith    The stencil type DMDA_STENCIL_STAR with width 1 corresponds to the
1881aa219208SBarry Smith    standard 5-pt stencil, while DMDA_STENCIL_BOX with width 1 denotes
188247c6ae99SBarry Smith    the standard 9-pt stencil.
188347c6ae99SBarry Smith 
1884aa219208SBarry Smith    The array data itself is NOT stored in the DMDA, it is stored in Vec objects;
1885564755cdSBarry Smith    The appropriate vector objects can be obtained with calls to DMCreateGlobalVector()
1886564755cdSBarry Smith    and DMCreateLocalVector() and calls to VecDuplicate() if more are needed.
188747c6ae99SBarry Smith 
188847c6ae99SBarry Smith .keywords: distributed array, create, two-dimensional
188947c6ae99SBarry Smith 
1890aa219208SBarry Smith .seealso: DMDestroy(), DMView(), DMDACreate1d(), DMDACreate3d(), DMGlobalToLocalBegin(), DMDAGetRefinementFactor(),
1891aa219208SBarry Smith           DMGlobalToLocalEnd(), DMLocalToGlobalBegin(), DMDALocalToLocalBegin(), DMDALocalToLocalEnd(), DMDASetRefinementFactor(),
1892d461ba97SBarry Smith           DMDAGetInfo(), DMCreateGlobalVector(), DMCreateLocalVector(), DMDACreateNaturalVector(), DMLoad(), DMDAGetOwnershipRanges()
189347c6ae99SBarry Smith 
189447c6ae99SBarry Smith @*/
1895fe16a2e9SBarry Smith 
18961321219cSEthan Coon PetscErrorCode  DMDACreate2d(MPI_Comm comm,DMDABoundaryType bx,DMDABoundaryType by,DMDAStencilType stencil_type,
18979a42bb27SBarry Smith                           PetscInt M,PetscInt N,PetscInt m,PetscInt n,PetscInt dof,PetscInt s,const PetscInt lx[],const PetscInt ly[],DM *da)
189847c6ae99SBarry Smith {
189947c6ae99SBarry Smith   PetscErrorCode ierr;
190047c6ae99SBarry Smith 
190147c6ae99SBarry Smith   PetscFunctionBegin;
1902aa219208SBarry Smith   ierr = DMDACreate(comm, da);CHKERRQ(ierr);
1903aa219208SBarry Smith   ierr = DMDASetDim(*da, 2);CHKERRQ(ierr);
1904aa219208SBarry Smith   ierr = DMDASetSizes(*da, M, N, 1);CHKERRQ(ierr);
1905aa219208SBarry Smith   ierr = DMDASetNumProcs(*da, m, n, PETSC_DECIDE);CHKERRQ(ierr);
1906755f258dSLisandro Dalcin   ierr = DMDASetBoundaryType(*da, bx, by, DMDA_BOUNDARY_NONE);CHKERRQ(ierr);
1907aa219208SBarry Smith   ierr = DMDASetDof(*da, dof);CHKERRQ(ierr);
1908aa219208SBarry Smith   ierr = DMDASetStencilType(*da, stencil_type);CHKERRQ(ierr);
1909aa219208SBarry Smith   ierr = DMDASetStencilWidth(*da, s);CHKERRQ(ierr);
1910aa219208SBarry Smith   ierr = DMDASetOwnershipRanges(*da, lx, ly, PETSC_NULL);CHKERRQ(ierr);
191147c6ae99SBarry Smith   /* This violates the behavior for other classes, but right now users expect negative dimensions to be handled this way */
19129a42bb27SBarry Smith   ierr = DMSetFromOptions(*da);CHKERRQ(ierr);
19139a42bb27SBarry Smith   ierr = DMSetUp(*da);CHKERRQ(ierr);
19147242296bSJed Brown   ierr = DMView_DA_Private(*da);CHKERRQ(ierr);
191547c6ae99SBarry Smith   PetscFunctionReturn(0);
191647c6ae99SBarry Smith }
1917