xref: /petsc/src/dm/impls/da/da2.c (revision 47c6ae997ffd1b2afd66b6474dff5950ae8613d1)
1*47c6ae99SBarry Smith #define PETSCDM_DLL
2*47c6ae99SBarry Smith 
3*47c6ae99SBarry Smith #include "private/daimpl.h"    /*I   "petscda.h"   I*/
4*47c6ae99SBarry Smith 
5*47c6ae99SBarry Smith #undef __FUNCT__
6*47c6ae99SBarry Smith #define __FUNCT__ "DAView_2d"
7*47c6ae99SBarry Smith PetscErrorCode DAView_2d(DA da,PetscViewer viewer)
8*47c6ae99SBarry Smith {
9*47c6ae99SBarry Smith   PetscErrorCode ierr;
10*47c6ae99SBarry Smith   PetscMPIInt    rank;
11*47c6ae99SBarry Smith   PetscBool      iascii,isdraw;
12*47c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
13*47c6ae99SBarry Smith 
14*47c6ae99SBarry Smith   PetscFunctionBegin;
15*47c6ae99SBarry Smith   ierr = MPI_Comm_rank(((PetscObject)da)->comm,&rank);CHKERRQ(ierr);
16*47c6ae99SBarry Smith 
17*47c6ae99SBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
18*47c6ae99SBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
19*47c6ae99SBarry Smith   if (iascii) {
20*47c6ae99SBarry Smith     PetscViewerFormat format;
21*47c6ae99SBarry Smith 
22*47c6ae99SBarry Smith     ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr);
23*47c6ae99SBarry Smith     if (format != PETSC_VIEWER_ASCII_VTK && format != PETSC_VIEWER_ASCII_VTK_CELL) {
24*47c6ae99SBarry Smith       DALocalInfo info;
25*47c6ae99SBarry Smith       ierr = DAGetLocalInfo(da,&info);CHKERRQ(ierr);
26*47c6ae99SBarry 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);
27*47c6ae99SBarry 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);
28*47c6ae99SBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
29*47c6ae99SBarry Smith     }
30*47c6ae99SBarry Smith   } else if (isdraw) {
31*47c6ae99SBarry Smith     PetscDraw       draw;
32*47c6ae99SBarry Smith     double     ymin = -1*dd->s-1,ymax = dd->N+dd->s;
33*47c6ae99SBarry Smith     double     xmin = -1*dd->s-1,xmax = dd->M+dd->s;
34*47c6ae99SBarry Smith     double     x,y;
35*47c6ae99SBarry Smith     PetscInt   base,*idx;
36*47c6ae99SBarry Smith     char       node[10];
37*47c6ae99SBarry Smith     PetscBool  isnull;
38*47c6ae99SBarry Smith 
39*47c6ae99SBarry Smith     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
40*47c6ae99SBarry Smith     ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
41*47c6ae99SBarry Smith     if (!dd->coordinates) {
42*47c6ae99SBarry Smith       ierr = PetscDrawSetCoordinates(draw,xmin,ymin,xmax,ymax);CHKERRQ(ierr);
43*47c6ae99SBarry Smith     }
44*47c6ae99SBarry Smith     ierr = PetscDrawSynchronizedClear(draw);CHKERRQ(ierr);
45*47c6ae99SBarry Smith 
46*47c6ae99SBarry Smith     /* first processor draw all node lines */
47*47c6ae99SBarry Smith     if (!rank) {
48*47c6ae99SBarry Smith       ymin = 0.0; ymax = dd->N - 1;
49*47c6ae99SBarry Smith       for (xmin=0; xmin<dd->M; xmin++) {
50*47c6ae99SBarry Smith         ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_BLACK);CHKERRQ(ierr);
51*47c6ae99SBarry Smith       }
52*47c6ae99SBarry Smith       xmin = 0.0; xmax = dd->M - 1;
53*47c6ae99SBarry Smith       for (ymin=0; ymin<dd->N; ymin++) {
54*47c6ae99SBarry Smith         ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_BLACK);CHKERRQ(ierr);
55*47c6ae99SBarry Smith       }
56*47c6ae99SBarry Smith     }
57*47c6ae99SBarry Smith     ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr);
58*47c6ae99SBarry Smith     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
59*47c6ae99SBarry Smith 
60*47c6ae99SBarry Smith     /* draw my box */
61*47c6ae99SBarry Smith     ymin = dd->ys; ymax = dd->ye - 1; xmin = dd->xs/dd->w;
62*47c6ae99SBarry Smith     xmax =(dd->xe-1)/dd->w;
63*47c6ae99SBarry Smith     ierr = PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_RED);CHKERRQ(ierr);
64*47c6ae99SBarry Smith     ierr = PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_RED);CHKERRQ(ierr);
65*47c6ae99SBarry Smith     ierr = PetscDrawLine(draw,xmin,ymax,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr);
66*47c6ae99SBarry Smith     ierr = PetscDrawLine(draw,xmax,ymin,xmax,ymax,PETSC_DRAW_RED);CHKERRQ(ierr);
67*47c6ae99SBarry Smith 
68*47c6ae99SBarry Smith     /* put in numbers */
69*47c6ae99SBarry Smith     base = (dd->base)/dd->w;
70*47c6ae99SBarry Smith     for (y=ymin; y<=ymax; y++) {
71*47c6ae99SBarry Smith       for (x=xmin; x<=xmax; x++) {
72*47c6ae99SBarry Smith         sprintf(node,"%d",(int)base++);
73*47c6ae99SBarry Smith         ierr = PetscDrawString(draw,x,y,PETSC_DRAW_BLACK,node);CHKERRQ(ierr);
74*47c6ae99SBarry Smith       }
75*47c6ae99SBarry Smith     }
76*47c6ae99SBarry Smith 
77*47c6ae99SBarry Smith     ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr);
78*47c6ae99SBarry Smith     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
79*47c6ae99SBarry Smith     /* overlay ghost numbers, useful for error checking */
80*47c6ae99SBarry Smith     /* put in numbers */
81*47c6ae99SBarry Smith 
82*47c6ae99SBarry Smith     base = 0; idx = dd->idx;
83*47c6ae99SBarry Smith     ymin = dd->Ys; ymax = dd->Ye; xmin = dd->Xs; xmax = dd->Xe;
84*47c6ae99SBarry Smith     for (y=ymin; y<ymax; y++) {
85*47c6ae99SBarry Smith       for (x=xmin; x<xmax; x++) {
86*47c6ae99SBarry Smith         if ((base % dd->w) == 0) {
87*47c6ae99SBarry Smith           sprintf(node,"%d",(int)(idx[base]/dd->w));
88*47c6ae99SBarry Smith           ierr = PetscDrawString(draw,x/dd->w,y,PETSC_DRAW_BLUE,node);CHKERRQ(ierr);
89*47c6ae99SBarry Smith         }
90*47c6ae99SBarry Smith         base++;
91*47c6ae99SBarry Smith       }
92*47c6ae99SBarry Smith     }
93*47c6ae99SBarry Smith     ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr);
94*47c6ae99SBarry Smith     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
95*47c6ae99SBarry Smith   } else {
96*47c6ae99SBarry Smith     SETERRQ1(((PetscObject)viewer)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for DA2d",((PetscObject)viewer)->type_name);
97*47c6ae99SBarry Smith   }
98*47c6ae99SBarry Smith   PetscFunctionReturn(0);
99*47c6ae99SBarry Smith }
100*47c6ae99SBarry Smith 
101*47c6ae99SBarry Smith /*
102*47c6ae99SBarry Smith       M is number of grid points
103*47c6ae99SBarry Smith       m is number of processors
104*47c6ae99SBarry Smith 
105*47c6ae99SBarry Smith */
106*47c6ae99SBarry Smith #undef __FUNCT__
107*47c6ae99SBarry Smith #define __FUNCT__ "DASplitComm2d"
108*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DASplitComm2d(MPI_Comm comm,PetscInt M,PetscInt N,PetscInt sw,MPI_Comm *outcomm)
109*47c6ae99SBarry Smith {
110*47c6ae99SBarry Smith   PetscErrorCode ierr;
111*47c6ae99SBarry Smith   PetscInt       m,n = 0,x = 0,y = 0;
112*47c6ae99SBarry Smith   PetscMPIInt    size,csize,rank;
113*47c6ae99SBarry Smith 
114*47c6ae99SBarry Smith   PetscFunctionBegin;
115*47c6ae99SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
116*47c6ae99SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
117*47c6ae99SBarry Smith 
118*47c6ae99SBarry Smith   csize = 4*size;
119*47c6ae99SBarry Smith   do {
120*47c6ae99SBarry 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);
121*47c6ae99SBarry Smith     csize   = csize/4;
122*47c6ae99SBarry Smith 
123*47c6ae99SBarry Smith     m = (PetscInt)(0.5 + sqrt(((double)M)*((double)csize)/((double)N)));
124*47c6ae99SBarry Smith     if (!m) m = 1;
125*47c6ae99SBarry Smith     while (m > 0) {
126*47c6ae99SBarry Smith       n = csize/m;
127*47c6ae99SBarry Smith       if (m*n == csize) break;
128*47c6ae99SBarry Smith       m--;
129*47c6ae99SBarry Smith     }
130*47c6ae99SBarry Smith     if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;}
131*47c6ae99SBarry Smith 
132*47c6ae99SBarry Smith     x = M/m + ((M % m) > ((csize-1) % m));
133*47c6ae99SBarry Smith     y = (N + (csize-1)/m)/n;
134*47c6ae99SBarry Smith   } while ((x < 4 || y < 4) && csize > 1);
135*47c6ae99SBarry Smith   if (size != csize) {
136*47c6ae99SBarry Smith     MPI_Group    entire_group,sub_group;
137*47c6ae99SBarry Smith     PetscMPIInt  i,*groupies;
138*47c6ae99SBarry Smith 
139*47c6ae99SBarry Smith     ierr = MPI_Comm_group(comm,&entire_group);CHKERRQ(ierr);
140*47c6ae99SBarry Smith     ierr = PetscMalloc(csize*sizeof(PetscInt),&groupies);CHKERRQ(ierr);
141*47c6ae99SBarry Smith     for (i=0; i<csize; i++) {
142*47c6ae99SBarry Smith       groupies[i] = (rank/csize)*csize + i;
143*47c6ae99SBarry Smith     }
144*47c6ae99SBarry Smith     ierr = MPI_Group_incl(entire_group,csize,groupies,&sub_group);CHKERRQ(ierr);
145*47c6ae99SBarry Smith     ierr = PetscFree(groupies);CHKERRQ(ierr);
146*47c6ae99SBarry Smith     ierr = MPI_Comm_create(comm,sub_group,outcomm);CHKERRQ(ierr);
147*47c6ae99SBarry Smith     ierr = MPI_Group_free(&entire_group);CHKERRQ(ierr);
148*47c6ae99SBarry Smith     ierr = MPI_Group_free(&sub_group);CHKERRQ(ierr);
149*47c6ae99SBarry Smith     ierr = PetscInfo1(0,"DASplitComm2d:Creating redundant coarse problems of size %d\n",csize);CHKERRQ(ierr);
150*47c6ae99SBarry Smith   } else {
151*47c6ae99SBarry Smith     *outcomm = comm;
152*47c6ae99SBarry Smith   }
153*47c6ae99SBarry Smith   PetscFunctionReturn(0);
154*47c6ae99SBarry Smith }
155*47c6ae99SBarry Smith 
156*47c6ae99SBarry Smith #undef __FUNCT__
157*47c6ae99SBarry Smith #define __FUNCT__ "DAGetElements_2d_P1"
158*47c6ae99SBarry Smith PetscErrorCode DAGetElements_2d_P1(DA da,PetscInt *n,const PetscInt *e[])
159*47c6ae99SBarry Smith {
160*47c6ae99SBarry Smith   PetscErrorCode ierr;
161*47c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
162*47c6ae99SBarry Smith   PetscInt       i,j,cnt,xs,xe = dd->xe,ys,ye = dd->ye,Xs = dd->Xs, Xe = dd->Xe, Ys = dd->Ys;
163*47c6ae99SBarry Smith 
164*47c6ae99SBarry Smith   PetscFunctionBegin;
165*47c6ae99SBarry Smith   if (!dd->e) {
166*47c6ae99SBarry Smith     if (dd->xs == Xs) xs = dd->xs; else xs = dd->xs - 1;
167*47c6ae99SBarry Smith     if (dd->ys == Ys) ys = dd->ys; else ys = dd->ys - 1;
168*47c6ae99SBarry Smith     dd->ne = 2*(xe - xs - 1)*(ye - ys - 1);
169*47c6ae99SBarry Smith     ierr   = PetscMalloc((1 + 3*dd->ne)*sizeof(PetscInt),&dd->e);CHKERRQ(ierr);
170*47c6ae99SBarry Smith     cnt    = 0;
171*47c6ae99SBarry Smith     for (j=ys; j<ye-1; j++) {
172*47c6ae99SBarry Smith       for (i=xs; i<xe-1; i++) {
173*47c6ae99SBarry Smith         dd->e[cnt]   = i - Xs + (j - Ys)*(Xe - Xs);
174*47c6ae99SBarry Smith         dd->e[cnt+1] = i - Xs + 1 + (j - Ys)*(Xe - Xs);
175*47c6ae99SBarry Smith         dd->e[cnt+2] = i - Xs + (j - Ys + 1)*(Xe - Xs);
176*47c6ae99SBarry Smith 
177*47c6ae99SBarry Smith         dd->e[cnt+3] = i - Xs + 1 + (j - Ys + 1)*(Xe - Xs);
178*47c6ae99SBarry Smith         dd->e[cnt+4] = i - Xs + (j - Ys + 1)*(Xe - Xs);
179*47c6ae99SBarry Smith         dd->e[cnt+5] = i - Xs + 1 + (j - Ys)*(Xe - Xs);
180*47c6ae99SBarry Smith         cnt += 6;
181*47c6ae99SBarry Smith       }
182*47c6ae99SBarry Smith     }
183*47c6ae99SBarry Smith   }
184*47c6ae99SBarry Smith   *n = dd->ne;
185*47c6ae99SBarry Smith   *e = dd->e;
186*47c6ae99SBarry Smith   PetscFunctionReturn(0);
187*47c6ae99SBarry Smith }
188*47c6ae99SBarry Smith 
189*47c6ae99SBarry Smith #undef __FUNCT__
190*47c6ae99SBarry Smith #define __FUNCT__ "DASetLocalFunction"
191*47c6ae99SBarry Smith /*@C
192*47c6ae99SBarry Smith        DASetLocalFunction - Caches in a DA a local function.
193*47c6ae99SBarry Smith 
194*47c6ae99SBarry Smith    Logically Collective on DA
195*47c6ae99SBarry Smith 
196*47c6ae99SBarry Smith    Input Parameter:
197*47c6ae99SBarry Smith +  da - initial distributed array
198*47c6ae99SBarry Smith -  lf - the local function
199*47c6ae99SBarry Smith 
200*47c6ae99SBarry Smith    Level: intermediate
201*47c6ae99SBarry Smith 
202*47c6ae99SBarry Smith    Notes: The routine SNESDAFormFunction() uses this the cached function to evaluate the user provided function.
203*47c6ae99SBarry Smith 
204*47c6ae99SBarry Smith .keywords:  distributed array, refine
205*47c6ae99SBarry Smith 
206*47c6ae99SBarry Smith .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DAGetLocalFunction(), DASetLocalFunctioni()
207*47c6ae99SBarry Smith @*/
208*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DASetLocalFunction(DA da,DALocalFunction1 lf)
209*47c6ae99SBarry Smith {
210*47c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
211*47c6ae99SBarry Smith   PetscFunctionBegin;
212*47c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
213*47c6ae99SBarry Smith   dd->lf    = lf;
214*47c6ae99SBarry Smith   PetscFunctionReturn(0);
215*47c6ae99SBarry Smith }
216*47c6ae99SBarry Smith 
217*47c6ae99SBarry Smith #undef __FUNCT__
218*47c6ae99SBarry Smith #define __FUNCT__ "DASetLocalFunctioni"
219*47c6ae99SBarry Smith /*@C
220*47c6ae99SBarry Smith        DASetLocalFunctioni - Caches in a DA a local function that evaluates a single component
221*47c6ae99SBarry Smith 
222*47c6ae99SBarry Smith    Logically Collective on DA
223*47c6ae99SBarry Smith 
224*47c6ae99SBarry Smith    Input Parameter:
225*47c6ae99SBarry Smith +  da - initial distributed array
226*47c6ae99SBarry Smith -  lfi - the local function
227*47c6ae99SBarry Smith 
228*47c6ae99SBarry Smith    Level: intermediate
229*47c6ae99SBarry Smith 
230*47c6ae99SBarry Smith .keywords:  distributed array, refine
231*47c6ae99SBarry Smith 
232*47c6ae99SBarry Smith .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DAGetLocalFunction(), DASetLocalFunction()
233*47c6ae99SBarry Smith @*/
234*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DASetLocalFunctioni(DA da,PetscErrorCode (*lfi)(DALocalInfo*,MatStencil*,void*,PetscScalar*,void*))
235*47c6ae99SBarry Smith {
236*47c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
237*47c6ae99SBarry Smith   PetscFunctionBegin;
238*47c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
239*47c6ae99SBarry Smith   dd->lfi = lfi;
240*47c6ae99SBarry Smith   PetscFunctionReturn(0);
241*47c6ae99SBarry Smith }
242*47c6ae99SBarry Smith 
243*47c6ae99SBarry Smith #undef __FUNCT__
244*47c6ae99SBarry Smith #define __FUNCT__ "DASetLocalFunctionib"
245*47c6ae99SBarry Smith /*@C
246*47c6ae99SBarry Smith        DASetLocalFunctionib - Caches in a DA a block local function that evaluates a single component
247*47c6ae99SBarry Smith 
248*47c6ae99SBarry Smith    Logically Collective on DA
249*47c6ae99SBarry Smith 
250*47c6ae99SBarry Smith    Input Parameter:
251*47c6ae99SBarry Smith +  da - initial distributed array
252*47c6ae99SBarry Smith -  lfi - the local function
253*47c6ae99SBarry Smith 
254*47c6ae99SBarry Smith    Level: intermediate
255*47c6ae99SBarry Smith 
256*47c6ae99SBarry Smith .keywords:  distributed array, refine
257*47c6ae99SBarry Smith 
258*47c6ae99SBarry Smith .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DAGetLocalFunction(), DASetLocalFunction()
259*47c6ae99SBarry Smith @*/
260*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DASetLocalFunctionib(DA da,PetscErrorCode (*lfi)(DALocalInfo*,MatStencil*,void*,PetscScalar*,void*))
261*47c6ae99SBarry Smith {
262*47c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
263*47c6ae99SBarry Smith   PetscFunctionBegin;
264*47c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
265*47c6ae99SBarry Smith   dd->lfib = lfi;
266*47c6ae99SBarry Smith   PetscFunctionReturn(0);
267*47c6ae99SBarry Smith }
268*47c6ae99SBarry Smith 
269*47c6ae99SBarry Smith #undef __FUNCT__
270*47c6ae99SBarry Smith #define __FUNCT__ "DASetLocalAdicFunction_Private"
271*47c6ae99SBarry Smith PetscErrorCode DASetLocalAdicFunction_Private(DA da,DALocalFunction1 ad_lf)
272*47c6ae99SBarry Smith {
273*47c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
274*47c6ae99SBarry Smith   PetscFunctionBegin;
275*47c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
276*47c6ae99SBarry Smith   dd->adic_lf = ad_lf;
277*47c6ae99SBarry Smith   PetscFunctionReturn(0);
278*47c6ae99SBarry Smith }
279*47c6ae99SBarry Smith 
280*47c6ae99SBarry Smith /*MC
281*47c6ae99SBarry Smith        DASetLocalAdicFunctioni - Caches in a DA a local functioni computed by ADIC/ADIFOR
282*47c6ae99SBarry Smith 
283*47c6ae99SBarry Smith    Synopsis:
284*47c6ae99SBarry Smith    PetscErrorCode DASetLocalAdicFunctioni(DA da,PetscInt (ad_lf*)(DALocalInfo*,MatStencil*,void*,void*,void*)
285*47c6ae99SBarry Smith 
286*47c6ae99SBarry Smith    Logically Collective on DA
287*47c6ae99SBarry Smith 
288*47c6ae99SBarry Smith    Input Parameter:
289*47c6ae99SBarry Smith +  da - initial distributed array
290*47c6ae99SBarry Smith -  ad_lfi - the local function as computed by ADIC/ADIFOR
291*47c6ae99SBarry Smith 
292*47c6ae99SBarry Smith    Level: intermediate
293*47c6ae99SBarry Smith 
294*47c6ae99SBarry Smith .keywords:  distributed array, refine
295*47c6ae99SBarry Smith 
296*47c6ae99SBarry Smith .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DAGetLocalFunction(), DASetLocalFunction(),
297*47c6ae99SBarry Smith           DASetLocalJacobian(), DASetLocalFunctioni()
298*47c6ae99SBarry Smith M*/
299*47c6ae99SBarry Smith 
300*47c6ae99SBarry Smith #undef __FUNCT__
301*47c6ae99SBarry Smith #define __FUNCT__ "DASetLocalAdicFunctioni_Private"
302*47c6ae99SBarry Smith PetscErrorCode DASetLocalAdicFunctioni_Private(DA da,PetscErrorCode (*ad_lfi)(DALocalInfo*,MatStencil*,void*,void*,void*))
303*47c6ae99SBarry Smith {
304*47c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
305*47c6ae99SBarry Smith   PetscFunctionBegin;
306*47c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
307*47c6ae99SBarry Smith   dd->adic_lfi = ad_lfi;
308*47c6ae99SBarry Smith   PetscFunctionReturn(0);
309*47c6ae99SBarry Smith }
310*47c6ae99SBarry Smith 
311*47c6ae99SBarry Smith /*MC
312*47c6ae99SBarry Smith        DASetLocalAdicMFFunctioni - Caches in a DA a local functioni computed by ADIC/ADIFOR
313*47c6ae99SBarry Smith 
314*47c6ae99SBarry Smith    Synopsis:
315*47c6ae99SBarry Smith    PetscErrorCode  DASetLocalAdicFunctioni(DA da,int (ad_lf*)(DALocalInfo*,MatStencil*,void*,void*,void*)
316*47c6ae99SBarry Smith 
317*47c6ae99SBarry Smith    Logically Collective on DA
318*47c6ae99SBarry Smith 
319*47c6ae99SBarry Smith    Input Parameter:
320*47c6ae99SBarry Smith +  da - initial distributed array
321*47c6ae99SBarry Smith -  admf_lfi - the local matrix-free function as computed by ADIC/ADIFOR
322*47c6ae99SBarry Smith 
323*47c6ae99SBarry Smith    Level: intermediate
324*47c6ae99SBarry Smith 
325*47c6ae99SBarry Smith .keywords:  distributed array, refine
326*47c6ae99SBarry Smith 
327*47c6ae99SBarry Smith .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DAGetLocalFunction(), DASetLocalFunction(),
328*47c6ae99SBarry Smith           DASetLocalJacobian(), DASetLocalFunctioni()
329*47c6ae99SBarry Smith M*/
330*47c6ae99SBarry Smith 
331*47c6ae99SBarry Smith #undef __FUNCT__
332*47c6ae99SBarry Smith #define __FUNCT__ "DASetLocalAdicMFFunctioni_Private"
333*47c6ae99SBarry Smith PetscErrorCode DASetLocalAdicMFFunctioni_Private(DA da,PetscErrorCode (*admf_lfi)(DALocalInfo*,MatStencil*,void*,void*,void*))
334*47c6ae99SBarry Smith {
335*47c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
336*47c6ae99SBarry Smith   PetscFunctionBegin;
337*47c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
338*47c6ae99SBarry Smith   dd->adicmf_lfi = admf_lfi;
339*47c6ae99SBarry Smith   PetscFunctionReturn(0);
340*47c6ae99SBarry Smith }
341*47c6ae99SBarry Smith 
342*47c6ae99SBarry Smith /*MC
343*47c6ae99SBarry Smith        DASetLocalAdicFunctionib - Caches in a DA a block local functioni computed by ADIC/ADIFOR
344*47c6ae99SBarry Smith 
345*47c6ae99SBarry Smith    Synopsis:
346*47c6ae99SBarry Smith    PetscErrorCode DASetLocalAdicFunctionib(DA da,PetscInt (ad_lf*)(DALocalInfo*,MatStencil*,void*,void*,void*)
347*47c6ae99SBarry Smith 
348*47c6ae99SBarry Smith    Logically Collective on DA
349*47c6ae99SBarry Smith 
350*47c6ae99SBarry Smith    Input Parameter:
351*47c6ae99SBarry Smith +  da - initial distributed array
352*47c6ae99SBarry Smith -  ad_lfi - the local function as computed by ADIC/ADIFOR
353*47c6ae99SBarry Smith 
354*47c6ae99SBarry Smith    Level: intermediate
355*47c6ae99SBarry Smith 
356*47c6ae99SBarry Smith .keywords:  distributed array, refine
357*47c6ae99SBarry Smith 
358*47c6ae99SBarry Smith .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DAGetLocalFunction(), DASetLocalFunction(),
359*47c6ae99SBarry Smith           DASetLocalJacobian(), DASetLocalFunctionib()
360*47c6ae99SBarry Smith M*/
361*47c6ae99SBarry Smith 
362*47c6ae99SBarry Smith #undef __FUNCT__
363*47c6ae99SBarry Smith #define __FUNCT__ "DASetLocalAdicFunctionib_Private"
364*47c6ae99SBarry Smith PetscErrorCode DASetLocalAdicFunctionib_Private(DA da,PetscErrorCode (*ad_lfi)(DALocalInfo*,MatStencil*,void*,void*,void*))
365*47c6ae99SBarry Smith {
366*47c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
367*47c6ae99SBarry Smith   PetscFunctionBegin;
368*47c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
369*47c6ae99SBarry Smith   dd->adic_lfib = ad_lfi;
370*47c6ae99SBarry Smith   PetscFunctionReturn(0);
371*47c6ae99SBarry Smith }
372*47c6ae99SBarry Smith 
373*47c6ae99SBarry Smith /*MC
374*47c6ae99SBarry Smith        DASetLocalAdicMFFunctionib - Caches in a DA a block local functioni computed by ADIC/ADIFOR
375*47c6ae99SBarry Smith 
376*47c6ae99SBarry Smith    Synopsis:
377*47c6ae99SBarry Smith    PetscErrorCode  DASetLocalAdicFunctionib(DA da,int (ad_lf*)(DALocalInfo*,MatStencil*,void*,void*,void*)
378*47c6ae99SBarry Smith 
379*47c6ae99SBarry Smith    Logically Collective on DA
380*47c6ae99SBarry Smith 
381*47c6ae99SBarry Smith    Input Parameter:
382*47c6ae99SBarry Smith +  da - initial distributed array
383*47c6ae99SBarry Smith -  admf_lfi - the local matrix-free function as computed by ADIC/ADIFOR
384*47c6ae99SBarry Smith 
385*47c6ae99SBarry Smith    Level: intermediate
386*47c6ae99SBarry Smith 
387*47c6ae99SBarry Smith .keywords:  distributed array, refine
388*47c6ae99SBarry Smith 
389*47c6ae99SBarry Smith .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DAGetLocalFunction(), DASetLocalFunction(),
390*47c6ae99SBarry Smith           DASetLocalJacobian(), DASetLocalFunctionib()
391*47c6ae99SBarry Smith M*/
392*47c6ae99SBarry Smith 
393*47c6ae99SBarry Smith #undef __FUNCT__
394*47c6ae99SBarry Smith #define __FUNCT__ "DASetLocalAdicMFFunctionib_Private"
395*47c6ae99SBarry Smith PetscErrorCode DASetLocalAdicMFFunctionib_Private(DA da,PetscErrorCode (*admf_lfi)(DALocalInfo*,MatStencil*,void*,void*,void*))
396*47c6ae99SBarry Smith {
397*47c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
398*47c6ae99SBarry Smith   PetscFunctionBegin;
399*47c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
400*47c6ae99SBarry Smith   dd->adicmf_lfib = admf_lfi;
401*47c6ae99SBarry Smith   PetscFunctionReturn(0);
402*47c6ae99SBarry Smith }
403*47c6ae99SBarry Smith 
404*47c6ae99SBarry Smith /*MC
405*47c6ae99SBarry Smith        DASetLocalAdicMFFunction - Caches in a DA a local function computed by ADIC/ADIFOR
406*47c6ae99SBarry Smith 
407*47c6ae99SBarry Smith    Synopsis:
408*47c6ae99SBarry Smith    PetscErrorCode DASetLocalAdicMFFunction(DA da,DALocalFunction1 ad_lf)
409*47c6ae99SBarry Smith 
410*47c6ae99SBarry Smith    Logically Collective on DA
411*47c6ae99SBarry Smith 
412*47c6ae99SBarry Smith    Input Parameter:
413*47c6ae99SBarry Smith +  da - initial distributed array
414*47c6ae99SBarry Smith -  ad_lf - the local function as computed by ADIC/ADIFOR
415*47c6ae99SBarry Smith 
416*47c6ae99SBarry Smith    Level: intermediate
417*47c6ae99SBarry Smith 
418*47c6ae99SBarry Smith .keywords:  distributed array, refine
419*47c6ae99SBarry Smith 
420*47c6ae99SBarry Smith .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DAGetLocalFunction(), DASetLocalFunction(),
421*47c6ae99SBarry Smith           DASetLocalJacobian()
422*47c6ae99SBarry Smith M*/
423*47c6ae99SBarry Smith 
424*47c6ae99SBarry Smith #undef __FUNCT__
425*47c6ae99SBarry Smith #define __FUNCT__ "DASetLocalAdicMFFunction_Private"
426*47c6ae99SBarry Smith PetscErrorCode DASetLocalAdicMFFunction_Private(DA da,DALocalFunction1 ad_lf)
427*47c6ae99SBarry Smith {
428*47c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
429*47c6ae99SBarry Smith   PetscFunctionBegin;
430*47c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
431*47c6ae99SBarry Smith   dd->adicmf_lf = ad_lf;
432*47c6ae99SBarry Smith   PetscFunctionReturn(0);
433*47c6ae99SBarry Smith }
434*47c6ae99SBarry Smith 
435*47c6ae99SBarry Smith /*@C
436*47c6ae99SBarry Smith        DASetLocalJacobian - Caches in a DA a local Jacobian computation function
437*47c6ae99SBarry Smith 
438*47c6ae99SBarry Smith    Logically Collective on DA
439*47c6ae99SBarry Smith 
440*47c6ae99SBarry Smith 
441*47c6ae99SBarry Smith    Input Parameter:
442*47c6ae99SBarry Smith +  da - initial distributed array
443*47c6ae99SBarry Smith -  lj - the local Jacobian
444*47c6ae99SBarry Smith 
445*47c6ae99SBarry Smith    Level: intermediate
446*47c6ae99SBarry Smith 
447*47c6ae99SBarry Smith    Notes: The routine SNESDAFormFunction() uses this the cached function to evaluate the user provided function.
448*47c6ae99SBarry Smith 
449*47c6ae99SBarry Smith .keywords:  distributed array, refine
450*47c6ae99SBarry Smith 
451*47c6ae99SBarry Smith .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DAGetLocalFunction(), DASetLocalFunction()
452*47c6ae99SBarry Smith @*/
453*47c6ae99SBarry Smith #undef __FUNCT__
454*47c6ae99SBarry Smith #define __FUNCT__ "DASetLocalJacobian"
455*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DASetLocalJacobian(DA da,DALocalFunction1 lj)
456*47c6ae99SBarry Smith {
457*47c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
458*47c6ae99SBarry Smith   PetscFunctionBegin;
459*47c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
460*47c6ae99SBarry Smith   dd->lj    = lj;
461*47c6ae99SBarry Smith   PetscFunctionReturn(0);
462*47c6ae99SBarry Smith }
463*47c6ae99SBarry Smith 
464*47c6ae99SBarry Smith #undef __FUNCT__
465*47c6ae99SBarry Smith #define __FUNCT__ "DAGetLocalFunction"
466*47c6ae99SBarry Smith /*@C
467*47c6ae99SBarry Smith        DAGetLocalFunction - Gets from a DA a local function and its ADIC/ADIFOR Jacobian
468*47c6ae99SBarry Smith 
469*47c6ae99SBarry Smith    Note Collective
470*47c6ae99SBarry Smith 
471*47c6ae99SBarry Smith    Input Parameter:
472*47c6ae99SBarry Smith .  da - initial distributed array
473*47c6ae99SBarry Smith 
474*47c6ae99SBarry Smith    Output Parameter:
475*47c6ae99SBarry Smith .  lf - the local function
476*47c6ae99SBarry Smith 
477*47c6ae99SBarry Smith    Level: intermediate
478*47c6ae99SBarry Smith 
479*47c6ae99SBarry Smith .keywords:  distributed array, refine
480*47c6ae99SBarry Smith 
481*47c6ae99SBarry Smith .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DAGetLocalJacobian(), DASetLocalFunction()
482*47c6ae99SBarry Smith @*/
483*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DAGetLocalFunction(DA da,DALocalFunction1 *lf)
484*47c6ae99SBarry Smith {
485*47c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
486*47c6ae99SBarry Smith   PetscFunctionBegin;
487*47c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
488*47c6ae99SBarry Smith   if (lf)       *lf = dd->lf;
489*47c6ae99SBarry Smith   PetscFunctionReturn(0);
490*47c6ae99SBarry Smith }
491*47c6ae99SBarry Smith 
492*47c6ae99SBarry Smith #undef __FUNCT__
493*47c6ae99SBarry Smith #define __FUNCT__ "DAGetLocalJacobian"
494*47c6ae99SBarry Smith /*@C
495*47c6ae99SBarry Smith        DAGetLocalJacobian - Gets from a DA a local jacobian
496*47c6ae99SBarry Smith 
497*47c6ae99SBarry Smith    Not Collective
498*47c6ae99SBarry Smith 
499*47c6ae99SBarry Smith    Input Parameter:
500*47c6ae99SBarry Smith .  da - initial distributed array
501*47c6ae99SBarry Smith 
502*47c6ae99SBarry Smith    Output Parameter:
503*47c6ae99SBarry Smith .  lj - the local jacobian
504*47c6ae99SBarry Smith 
505*47c6ae99SBarry Smith    Level: intermediate
506*47c6ae99SBarry Smith 
507*47c6ae99SBarry Smith .keywords:  distributed array, refine
508*47c6ae99SBarry Smith 
509*47c6ae99SBarry Smith .seealso: DACreate1d(), DACreate2d(), DACreate3d(), DADestroy(), DAGetLocalFunction(), DASetLocalJacobian()
510*47c6ae99SBarry Smith @*/
511*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DAGetLocalJacobian(DA da,DALocalFunction1 *lj)
512*47c6ae99SBarry Smith {
513*47c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
514*47c6ae99SBarry Smith   PetscFunctionBegin;
515*47c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
516*47c6ae99SBarry Smith   if (lj) *lj = dd->lj;
517*47c6ae99SBarry Smith   PetscFunctionReturn(0);
518*47c6ae99SBarry Smith }
519*47c6ae99SBarry Smith 
520*47c6ae99SBarry Smith #undef __FUNCT__
521*47c6ae99SBarry Smith #define __FUNCT__ "DAFormFunction"
522*47c6ae99SBarry Smith /*@
523*47c6ae99SBarry Smith     DAFormFunction - Evaluates a user provided function on each processor that
524*47c6ae99SBarry Smith         share a DA
525*47c6ae99SBarry Smith 
526*47c6ae99SBarry Smith    Input Parameters:
527*47c6ae99SBarry Smith +    da - the DA that defines the grid
528*47c6ae99SBarry Smith .    vu - input vector
529*47c6ae99SBarry Smith .    vfu - output vector
530*47c6ae99SBarry Smith -    w - any user data
531*47c6ae99SBarry Smith 
532*47c6ae99SBarry Smith     Notes: Does NOT do ghost updates on vu upon entry
533*47c6ae99SBarry Smith 
534*47c6ae99SBarry Smith            This should eventually replace DAFormFunction1
535*47c6ae99SBarry Smith 
536*47c6ae99SBarry Smith     Level: advanced
537*47c6ae99SBarry Smith 
538*47c6ae99SBarry Smith .seealso: DAComputeJacobian1WithAdic()
539*47c6ae99SBarry Smith 
540*47c6ae99SBarry Smith @*/
541*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DAFormFunction(DA da,PetscErrorCode (*lf)(void),Vec vu,Vec vfu,void *w)
542*47c6ae99SBarry Smith {
543*47c6ae99SBarry Smith   PetscErrorCode ierr;
544*47c6ae99SBarry Smith   void           *u,*fu;
545*47c6ae99SBarry Smith   DALocalInfo    info;
546*47c6ae99SBarry Smith   PetscErrorCode (*f)(DALocalInfo*,void*,void*,void*) = (PetscErrorCode (*)(DALocalInfo*,void*,void*,void*))lf;
547*47c6ae99SBarry Smith 
548*47c6ae99SBarry Smith   PetscFunctionBegin;
549*47c6ae99SBarry Smith   ierr = DAGetLocalInfo(da,&info);CHKERRQ(ierr);
550*47c6ae99SBarry Smith   ierr = DAVecGetArray(da,vu,&u);CHKERRQ(ierr);
551*47c6ae99SBarry Smith   ierr = DAVecGetArray(da,vfu,&fu);CHKERRQ(ierr);
552*47c6ae99SBarry Smith 
553*47c6ae99SBarry Smith   ierr = (*f)(&info,u,fu,w);
554*47c6ae99SBarry Smith   if (PetscExceptionValue(ierr)) {
555*47c6ae99SBarry Smith     PetscErrorCode pierr = DAVecRestoreArray(da,vu,&u);CHKERRQ(pierr);
556*47c6ae99SBarry Smith     pierr = DAVecRestoreArray(da,vfu,&fu);CHKERRQ(pierr);
557*47c6ae99SBarry Smith   }
558*47c6ae99SBarry Smith  CHKERRQ(ierr);
559*47c6ae99SBarry Smith 
560*47c6ae99SBarry Smith   ierr = DAVecRestoreArray(da,vu,&u);CHKERRQ(ierr);
561*47c6ae99SBarry Smith   ierr = DAVecRestoreArray(da,vfu,&fu);CHKERRQ(ierr);
562*47c6ae99SBarry Smith   PetscFunctionReturn(0);
563*47c6ae99SBarry Smith }
564*47c6ae99SBarry Smith 
565*47c6ae99SBarry Smith #undef __FUNCT__
566*47c6ae99SBarry Smith #define __FUNCT__ "DAFormFunctionLocal"
567*47c6ae99SBarry Smith /*@C
568*47c6ae99SBarry Smith    DAFormFunctionLocal - This is a universal function evaluation routine for
569*47c6ae99SBarry Smith    a local DA function.
570*47c6ae99SBarry Smith 
571*47c6ae99SBarry Smith    Collective on DA
572*47c6ae99SBarry Smith 
573*47c6ae99SBarry Smith    Input Parameters:
574*47c6ae99SBarry Smith +  da - the DA context
575*47c6ae99SBarry Smith .  func - The local function
576*47c6ae99SBarry Smith .  X - input vector
577*47c6ae99SBarry Smith .  F - function vector
578*47c6ae99SBarry Smith -  ctx - A user context
579*47c6ae99SBarry Smith 
580*47c6ae99SBarry Smith    Level: intermediate
581*47c6ae99SBarry Smith 
582*47c6ae99SBarry Smith .seealso: DASetLocalFunction(), DASetLocalJacobian(), DASetLocalAdicFunction(), DASetLocalAdicMFFunction(),
583*47c6ae99SBarry Smith           SNESSetFunction(), SNESSetJacobian()
584*47c6ae99SBarry Smith 
585*47c6ae99SBarry Smith @*/
586*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DAFormFunctionLocal(DA da, DALocalFunction1 func, Vec X, Vec F, void *ctx)
587*47c6ae99SBarry Smith {
588*47c6ae99SBarry Smith   Vec            localX;
589*47c6ae99SBarry Smith   DALocalInfo    info;
590*47c6ae99SBarry Smith   void          *u;
591*47c6ae99SBarry Smith   void          *fu;
592*47c6ae99SBarry Smith   PetscErrorCode ierr;
593*47c6ae99SBarry Smith 
594*47c6ae99SBarry Smith   PetscFunctionBegin;
595*47c6ae99SBarry Smith   ierr = DAGetLocalVector(da,&localX);CHKERRQ(ierr);
596*47c6ae99SBarry Smith   /*
597*47c6ae99SBarry Smith      Scatter ghost points to local vector, using the 2-step process
598*47c6ae99SBarry Smith         DAGlobalToLocalBegin(), DAGlobalToLocalEnd().
599*47c6ae99SBarry Smith   */
600*47c6ae99SBarry Smith   ierr = DAGlobalToLocalBegin(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
601*47c6ae99SBarry Smith   ierr = DAGlobalToLocalEnd(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
602*47c6ae99SBarry Smith   ierr = DAGetLocalInfo(da,&info);CHKERRQ(ierr);
603*47c6ae99SBarry Smith   ierr = DAVecGetArray(da,localX,&u);CHKERRQ(ierr);
604*47c6ae99SBarry Smith   ierr = DAVecGetArray(da,F,&fu);CHKERRQ(ierr);
605*47c6ae99SBarry Smith   ierr = (*func)(&info,u,fu,ctx);
606*47c6ae99SBarry Smith   if (PetscExceptionValue(ierr)) {
607*47c6ae99SBarry Smith     PetscErrorCode pierr = DAVecRestoreArray(da,localX,&u);CHKERRQ(pierr);
608*47c6ae99SBarry Smith     pierr = DAVecRestoreArray(da,F,&fu);CHKERRQ(pierr);
609*47c6ae99SBarry Smith   }
610*47c6ae99SBarry Smith  CHKERRQ(ierr);
611*47c6ae99SBarry Smith   ierr = DAVecRestoreArray(da,localX,&u);CHKERRQ(ierr);
612*47c6ae99SBarry Smith   ierr = DAVecRestoreArray(da,F,&fu);CHKERRQ(ierr);
613*47c6ae99SBarry Smith   if (PetscExceptionValue(ierr)) {
614*47c6ae99SBarry Smith     PetscErrorCode pierr = DARestoreLocalVector(da,&localX);CHKERRQ(pierr);
615*47c6ae99SBarry Smith   }
616*47c6ae99SBarry Smith  CHKERRQ(ierr);
617*47c6ae99SBarry Smith   ierr = DARestoreLocalVector(da,&localX);CHKERRQ(ierr);
618*47c6ae99SBarry Smith   PetscFunctionReturn(0);
619*47c6ae99SBarry Smith }
620*47c6ae99SBarry Smith 
621*47c6ae99SBarry Smith #undef __FUNCT__
622*47c6ae99SBarry Smith #define __FUNCT__ "DAFormFunctionLocalGhost"
623*47c6ae99SBarry Smith /*@C
624*47c6ae99SBarry Smith    DAFormFunctionLocalGhost - This is a universal function evaluation routine for
625*47c6ae99SBarry Smith    a local DA function, but the ghost values of the output are communicated and added.
626*47c6ae99SBarry Smith 
627*47c6ae99SBarry Smith    Collective on DA
628*47c6ae99SBarry Smith 
629*47c6ae99SBarry Smith    Input Parameters:
630*47c6ae99SBarry Smith +  da - the DA context
631*47c6ae99SBarry Smith .  func - The local function
632*47c6ae99SBarry Smith .  X - input vector
633*47c6ae99SBarry Smith .  F - function vector
634*47c6ae99SBarry Smith -  ctx - A user context
635*47c6ae99SBarry Smith 
636*47c6ae99SBarry Smith    Level: intermediate
637*47c6ae99SBarry Smith 
638*47c6ae99SBarry Smith .seealso: DASetLocalFunction(), DASetLocalJacobian(), DASetLocalAdicFunction(), DASetLocalAdicMFFunction(),
639*47c6ae99SBarry Smith           SNESSetFunction(), SNESSetJacobian()
640*47c6ae99SBarry Smith 
641*47c6ae99SBarry Smith @*/
642*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DAFormFunctionLocalGhost(DA da, DALocalFunction1 func, Vec X, Vec F, void *ctx)
643*47c6ae99SBarry Smith {
644*47c6ae99SBarry Smith   Vec            localX, localF;
645*47c6ae99SBarry Smith   DALocalInfo    info;
646*47c6ae99SBarry Smith   void          *u;
647*47c6ae99SBarry Smith   void          *fu;
648*47c6ae99SBarry Smith   PetscErrorCode ierr;
649*47c6ae99SBarry Smith 
650*47c6ae99SBarry Smith   PetscFunctionBegin;
651*47c6ae99SBarry Smith   ierr = DAGetLocalVector(da,&localX);CHKERRQ(ierr);
652*47c6ae99SBarry Smith   ierr = DAGetLocalVector(da,&localF);CHKERRQ(ierr);
653*47c6ae99SBarry Smith   /*
654*47c6ae99SBarry Smith      Scatter ghost points to local vector, using the 2-step process
655*47c6ae99SBarry Smith         DAGlobalToLocalBegin(), DAGlobalToLocalEnd().
656*47c6ae99SBarry Smith   */
657*47c6ae99SBarry Smith   ierr = DAGlobalToLocalBegin(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
658*47c6ae99SBarry Smith   ierr = DAGlobalToLocalEnd(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
659*47c6ae99SBarry Smith   ierr = VecSet(F, 0.0);CHKERRQ(ierr);
660*47c6ae99SBarry Smith   ierr = VecSet(localF, 0.0);CHKERRQ(ierr);
661*47c6ae99SBarry Smith   ierr = DAGetLocalInfo(da,&info);CHKERRQ(ierr);
662*47c6ae99SBarry Smith   ierr = DAVecGetArray(da,localX,&u);CHKERRQ(ierr);
663*47c6ae99SBarry Smith   ierr = DAVecGetArray(da,localF,&fu);CHKERRQ(ierr);
664*47c6ae99SBarry Smith   ierr = (*func)(&info,u,fu,ctx);
665*47c6ae99SBarry Smith   if (PetscExceptionValue(ierr)) {
666*47c6ae99SBarry Smith     PetscErrorCode pierr = DAVecRestoreArray(da,localX,&u);CHKERRQ(pierr);
667*47c6ae99SBarry Smith     pierr = DAVecRestoreArray(da,localF,&fu);CHKERRQ(pierr);
668*47c6ae99SBarry Smith   }
669*47c6ae99SBarry Smith   CHKERRQ(ierr);
670*47c6ae99SBarry Smith   ierr = DALocalToGlobalBegin(da,localF,F);CHKERRQ(ierr);
671*47c6ae99SBarry Smith   ierr = DALocalToGlobalEnd(da,localF,F);CHKERRQ(ierr);
672*47c6ae99SBarry Smith   ierr = DAVecRestoreArray(da,localX,&u);CHKERRQ(ierr);
673*47c6ae99SBarry Smith   ierr = DAVecRestoreArray(da,localF,&fu);CHKERRQ(ierr);
674*47c6ae99SBarry Smith   if (PetscExceptionValue(ierr)) {
675*47c6ae99SBarry Smith     PetscErrorCode pierr = DARestoreLocalVector(da,&localX);CHKERRQ(pierr);
676*47c6ae99SBarry Smith   ierr = DARestoreLocalVector(da,&localF);CHKERRQ(ierr);
677*47c6ae99SBarry Smith   }
678*47c6ae99SBarry Smith   CHKERRQ(ierr);
679*47c6ae99SBarry Smith   ierr = DARestoreLocalVector(da,&localX);CHKERRQ(ierr);
680*47c6ae99SBarry Smith   ierr = DARestoreLocalVector(da,&localF);CHKERRQ(ierr);
681*47c6ae99SBarry Smith   PetscFunctionReturn(0);
682*47c6ae99SBarry Smith }
683*47c6ae99SBarry Smith 
684*47c6ae99SBarry Smith #undef __FUNCT__
685*47c6ae99SBarry Smith #define __FUNCT__ "DAFormFunction1"
686*47c6ae99SBarry Smith /*@
687*47c6ae99SBarry Smith     DAFormFunction1 - Evaluates a user provided function on each processor that
688*47c6ae99SBarry Smith         share a DA
689*47c6ae99SBarry Smith 
690*47c6ae99SBarry Smith    Input Parameters:
691*47c6ae99SBarry Smith +    da - the DA that defines the grid
692*47c6ae99SBarry Smith .    vu - input vector
693*47c6ae99SBarry Smith .    vfu - output vector
694*47c6ae99SBarry Smith -    w - any user data
695*47c6ae99SBarry Smith 
696*47c6ae99SBarry Smith     Notes: Does NOT do ghost updates on vu upon entry
697*47c6ae99SBarry Smith 
698*47c6ae99SBarry Smith     Level: advanced
699*47c6ae99SBarry Smith 
700*47c6ae99SBarry Smith .seealso: DAComputeJacobian1WithAdic()
701*47c6ae99SBarry Smith 
702*47c6ae99SBarry Smith @*/
703*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DAFormFunction1(DA da,Vec vu,Vec vfu,void *w)
704*47c6ae99SBarry Smith {
705*47c6ae99SBarry Smith   PetscErrorCode ierr;
706*47c6ae99SBarry Smith   void           *u,*fu;
707*47c6ae99SBarry Smith   DALocalInfo    info;
708*47c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
709*47c6ae99SBarry Smith 
710*47c6ae99SBarry Smith   PetscFunctionBegin;
711*47c6ae99SBarry Smith 
712*47c6ae99SBarry Smith   ierr = DAGetLocalInfo(da,&info);CHKERRQ(ierr);
713*47c6ae99SBarry Smith   ierr = DAVecGetArray(da,vu,&u);CHKERRQ(ierr);
714*47c6ae99SBarry Smith   ierr = DAVecGetArray(da,vfu,&fu);CHKERRQ(ierr);
715*47c6ae99SBarry Smith 
716*47c6ae99SBarry Smith   CHKMEMQ;
717*47c6ae99SBarry Smith   ierr = (*dd->lf)(&info,u,fu,w);
718*47c6ae99SBarry Smith   if (PetscExceptionValue(ierr)) {
719*47c6ae99SBarry Smith     PetscErrorCode pierr = DAVecRestoreArray(da,vu,&u);CHKERRQ(pierr);
720*47c6ae99SBarry Smith     pierr = DAVecRestoreArray(da,vfu,&fu);CHKERRQ(pierr);
721*47c6ae99SBarry Smith   }
722*47c6ae99SBarry Smith   CHKERRQ(ierr);
723*47c6ae99SBarry Smith   CHKMEMQ;
724*47c6ae99SBarry Smith 
725*47c6ae99SBarry Smith   ierr = DAVecRestoreArray(da,vu,&u);CHKERRQ(ierr);
726*47c6ae99SBarry Smith   ierr = DAVecRestoreArray(da,vfu,&fu);CHKERRQ(ierr);
727*47c6ae99SBarry Smith   PetscFunctionReturn(0);
728*47c6ae99SBarry Smith }
729*47c6ae99SBarry Smith 
730*47c6ae99SBarry Smith #undef __FUNCT__
731*47c6ae99SBarry Smith #define __FUNCT__ "DAFormFunctioniTest1"
732*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DAFormFunctioniTest1(DA da,void *w)
733*47c6ae99SBarry Smith {
734*47c6ae99SBarry Smith   Vec            vu,fu,fui;
735*47c6ae99SBarry Smith   PetscErrorCode ierr;
736*47c6ae99SBarry Smith   PetscInt       i,n;
737*47c6ae99SBarry Smith   PetscScalar    *ui;
738*47c6ae99SBarry Smith   PetscRandom    rnd;
739*47c6ae99SBarry Smith   PetscReal      norm;
740*47c6ae99SBarry Smith 
741*47c6ae99SBarry Smith   PetscFunctionBegin;
742*47c6ae99SBarry Smith   ierr = DAGetLocalVector(da,&vu);CHKERRQ(ierr);
743*47c6ae99SBarry Smith   ierr = PetscRandomCreate(PETSC_COMM_SELF,&rnd);CHKERRQ(ierr);
744*47c6ae99SBarry Smith   ierr = PetscRandomSetFromOptions(rnd);CHKERRQ(ierr);
745*47c6ae99SBarry Smith   ierr = VecSetRandom(vu,rnd);CHKERRQ(ierr);
746*47c6ae99SBarry Smith   ierr = PetscRandomDestroy(rnd);CHKERRQ(ierr);
747*47c6ae99SBarry Smith 
748*47c6ae99SBarry Smith   ierr = DAGetGlobalVector(da,&fu);CHKERRQ(ierr);
749*47c6ae99SBarry Smith   ierr = DAGetGlobalVector(da,&fui);CHKERRQ(ierr);
750*47c6ae99SBarry Smith 
751*47c6ae99SBarry Smith   ierr = DAFormFunction1(da,vu,fu,w);CHKERRQ(ierr);
752*47c6ae99SBarry Smith 
753*47c6ae99SBarry Smith   ierr = VecGetArray(fui,&ui);CHKERRQ(ierr);
754*47c6ae99SBarry Smith   ierr = VecGetLocalSize(fui,&n);CHKERRQ(ierr);
755*47c6ae99SBarry Smith   for (i=0; i<n; i++) {
756*47c6ae99SBarry Smith     ierr = DAFormFunctioni1(da,i,vu,ui+i,w);CHKERRQ(ierr);
757*47c6ae99SBarry Smith   }
758*47c6ae99SBarry Smith   ierr = VecRestoreArray(fui,&ui);CHKERRQ(ierr);
759*47c6ae99SBarry Smith 
760*47c6ae99SBarry Smith   ierr = VecAXPY(fui,-1.0,fu);CHKERRQ(ierr);
761*47c6ae99SBarry Smith   ierr = VecNorm(fui,NORM_2,&norm);CHKERRQ(ierr);
762*47c6ae99SBarry Smith   ierr = PetscPrintf(((PetscObject)da)->comm,"Norm of difference in vectors %G\n",norm);CHKERRQ(ierr);
763*47c6ae99SBarry Smith   ierr = VecView(fu,0);CHKERRQ(ierr);
764*47c6ae99SBarry Smith   ierr = VecView(fui,0);CHKERRQ(ierr);
765*47c6ae99SBarry Smith 
766*47c6ae99SBarry Smith   ierr = DARestoreLocalVector(da,&vu);CHKERRQ(ierr);
767*47c6ae99SBarry Smith   ierr = DARestoreGlobalVector(da,&fu);CHKERRQ(ierr);
768*47c6ae99SBarry Smith   ierr = DARestoreGlobalVector(da,&fui);CHKERRQ(ierr);
769*47c6ae99SBarry Smith   PetscFunctionReturn(0);
770*47c6ae99SBarry Smith }
771*47c6ae99SBarry Smith 
772*47c6ae99SBarry Smith #undef __FUNCT__
773*47c6ae99SBarry Smith #define __FUNCT__ "DAFormFunctioni1"
774*47c6ae99SBarry Smith /*@
775*47c6ae99SBarry Smith     DAFormFunctioni1 - Evaluates a user provided point-wise function
776*47c6ae99SBarry Smith 
777*47c6ae99SBarry Smith    Input Parameters:
778*47c6ae99SBarry Smith +    da - the DA that defines the grid
779*47c6ae99SBarry Smith .    i - the component of the function we wish to compute (must be local)
780*47c6ae99SBarry Smith .    vu - input vector
781*47c6ae99SBarry Smith .    vfu - output value
782*47c6ae99SBarry Smith -    w - any user data
783*47c6ae99SBarry Smith 
784*47c6ae99SBarry Smith     Notes: Does NOT do ghost updates on vu upon entry
785*47c6ae99SBarry Smith 
786*47c6ae99SBarry Smith     Level: advanced
787*47c6ae99SBarry Smith 
788*47c6ae99SBarry Smith .seealso: DAComputeJacobian1WithAdic()
789*47c6ae99SBarry Smith 
790*47c6ae99SBarry Smith @*/
791*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DAFormFunctioni1(DA da,PetscInt i,Vec vu,PetscScalar *vfu,void *w)
792*47c6ae99SBarry Smith {
793*47c6ae99SBarry Smith   PetscErrorCode ierr;
794*47c6ae99SBarry Smith   void           *u;
795*47c6ae99SBarry Smith   DALocalInfo    info;
796*47c6ae99SBarry Smith   MatStencil     stencil;
797*47c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
798*47c6ae99SBarry Smith 
799*47c6ae99SBarry Smith   PetscFunctionBegin;
800*47c6ae99SBarry Smith 
801*47c6ae99SBarry Smith   ierr = DAGetLocalInfo(da,&info);CHKERRQ(ierr);
802*47c6ae99SBarry Smith   ierr = DAVecGetArray(da,vu,&u);CHKERRQ(ierr);
803*47c6ae99SBarry Smith 
804*47c6ae99SBarry Smith   /* figure out stencil value from i */
805*47c6ae99SBarry Smith   stencil.c = i % info.dof;
806*47c6ae99SBarry Smith   stencil.i = (i % (info.xm*info.dof))/info.dof;
807*47c6ae99SBarry Smith   stencil.j = (i % (info.xm*info.ym*info.dof))/(info.xm*info.dof);
808*47c6ae99SBarry Smith   stencil.k = i/(info.xm*info.ym*info.dof);
809*47c6ae99SBarry Smith 
810*47c6ae99SBarry Smith   ierr = (*dd->lfi)(&info,&stencil,u,vfu,w);CHKERRQ(ierr);
811*47c6ae99SBarry Smith 
812*47c6ae99SBarry Smith   ierr = DAVecRestoreArray(da,vu,&u);CHKERRQ(ierr);
813*47c6ae99SBarry Smith   PetscFunctionReturn(0);
814*47c6ae99SBarry Smith }
815*47c6ae99SBarry Smith 
816*47c6ae99SBarry Smith #undef __FUNCT__
817*47c6ae99SBarry Smith #define __FUNCT__ "DAFormFunctionib1"
818*47c6ae99SBarry Smith /*@
819*47c6ae99SBarry Smith     DAFormFunctionib1 - Evaluates a user provided point-block function
820*47c6ae99SBarry Smith 
821*47c6ae99SBarry Smith    Input Parameters:
822*47c6ae99SBarry Smith +    da - the DA that defines the grid
823*47c6ae99SBarry Smith .    i - the component of the function we wish to compute (must be local)
824*47c6ae99SBarry Smith .    vu - input vector
825*47c6ae99SBarry Smith .    vfu - output value
826*47c6ae99SBarry Smith -    w - any user data
827*47c6ae99SBarry Smith 
828*47c6ae99SBarry Smith     Notes: Does NOT do ghost updates on vu upon entry
829*47c6ae99SBarry Smith 
830*47c6ae99SBarry Smith     Level: advanced
831*47c6ae99SBarry Smith 
832*47c6ae99SBarry Smith .seealso: DAComputeJacobian1WithAdic()
833*47c6ae99SBarry Smith 
834*47c6ae99SBarry Smith @*/
835*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DAFormFunctionib1(DA da,PetscInt i,Vec vu,PetscScalar *vfu,void *w)
836*47c6ae99SBarry Smith {
837*47c6ae99SBarry Smith   PetscErrorCode ierr;
838*47c6ae99SBarry Smith   void           *u;
839*47c6ae99SBarry Smith   DALocalInfo    info;
840*47c6ae99SBarry Smith   MatStencil     stencil;
841*47c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
842*47c6ae99SBarry Smith 
843*47c6ae99SBarry Smith   PetscFunctionBegin;
844*47c6ae99SBarry Smith   ierr = DAGetLocalInfo(da,&info);CHKERRQ(ierr);
845*47c6ae99SBarry Smith   ierr = DAVecGetArray(da,vu,&u);CHKERRQ(ierr);
846*47c6ae99SBarry Smith 
847*47c6ae99SBarry Smith   /* figure out stencil value from i */
848*47c6ae99SBarry Smith   stencil.c = i % info.dof;
849*47c6ae99SBarry Smith   if (stencil.c) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONG,"Point-block functions can only be called for the entire block");
850*47c6ae99SBarry Smith   stencil.i = (i % (info.xm*info.dof))/info.dof;
851*47c6ae99SBarry Smith   stencil.j = (i % (info.xm*info.ym*info.dof))/(info.xm*info.dof);
852*47c6ae99SBarry Smith   stencil.k = i/(info.xm*info.ym*info.dof);
853*47c6ae99SBarry Smith 
854*47c6ae99SBarry Smith   ierr = (*dd->lfib)(&info,&stencil,u,vfu,w);CHKERRQ(ierr);
855*47c6ae99SBarry Smith 
856*47c6ae99SBarry Smith   ierr = DAVecRestoreArray(da,vu,&u);CHKERRQ(ierr);
857*47c6ae99SBarry Smith   PetscFunctionReturn(0);
858*47c6ae99SBarry Smith }
859*47c6ae99SBarry Smith 
860*47c6ae99SBarry Smith #if defined(new)
861*47c6ae99SBarry Smith #undef __FUNCT__
862*47c6ae99SBarry Smith #define __FUNCT__ "DAGetDiagonal_MFFD"
863*47c6ae99SBarry Smith /*
864*47c6ae99SBarry Smith   DAGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix where local
865*47c6ae99SBarry Smith     function lives on a DA
866*47c6ae99SBarry Smith 
867*47c6ae99SBarry Smith         y ~= (F(u + ha) - F(u))/h,
868*47c6ae99SBarry Smith   where F = nonlinear function, as set by SNESSetFunction()
869*47c6ae99SBarry Smith         u = current iterate
870*47c6ae99SBarry Smith         h = difference interval
871*47c6ae99SBarry Smith */
872*47c6ae99SBarry Smith PetscErrorCode DAGetDiagonal_MFFD(DA da,Vec U,Vec a)
873*47c6ae99SBarry Smith {
874*47c6ae99SBarry Smith   PetscScalar    h,*aa,*ww,v;
875*47c6ae99SBarry Smith   PetscReal      epsilon = PETSC_SQRT_MACHINE_EPSILON,umin = 100.0*PETSC_SQRT_MACHINE_EPSILON;
876*47c6ae99SBarry Smith   PetscErrorCode ierr;
877*47c6ae99SBarry Smith   PetscInt       gI,nI;
878*47c6ae99SBarry Smith   MatStencil     stencil;
879*47c6ae99SBarry Smith   DALocalInfo    info;
880*47c6ae99SBarry Smith 
881*47c6ae99SBarry Smith   PetscFunctionBegin;
882*47c6ae99SBarry Smith   ierr = (*ctx->func)(0,U,a,ctx->funcctx);CHKERRQ(ierr);
883*47c6ae99SBarry Smith   ierr = (*ctx->funcisetbase)(U,ctx->funcctx);CHKERRQ(ierr);
884*47c6ae99SBarry Smith 
885*47c6ae99SBarry Smith   ierr = VecGetArray(U,&ww);CHKERRQ(ierr);
886*47c6ae99SBarry Smith   ierr = VecGetArray(a,&aa);CHKERRQ(ierr);
887*47c6ae99SBarry Smith 
888*47c6ae99SBarry Smith   nI = 0;
889*47c6ae99SBarry Smith     h  = ww[gI];
890*47c6ae99SBarry Smith     if (h == 0.0) h = 1.0;
891*47c6ae99SBarry Smith #if !defined(PETSC_USE_COMPLEX)
892*47c6ae99SBarry Smith     if (h < umin && h >= 0.0)      h = umin;
893*47c6ae99SBarry Smith     else if (h < 0.0 && h > -umin) h = -umin;
894*47c6ae99SBarry Smith #else
895*47c6ae99SBarry Smith     if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0)     h = umin;
896*47c6ae99SBarry Smith     else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin;
897*47c6ae99SBarry Smith #endif
898*47c6ae99SBarry Smith     h     *= epsilon;
899*47c6ae99SBarry Smith 
900*47c6ae99SBarry Smith     ww[gI] += h;
901*47c6ae99SBarry Smith     ierr          = (*ctx->funci)(i,w,&v,ctx->funcctx);CHKERRQ(ierr);
902*47c6ae99SBarry Smith     aa[nI]  = (v - aa[nI])/h;
903*47c6ae99SBarry Smith     ww[gI] -= h;
904*47c6ae99SBarry Smith     nI++;
905*47c6ae99SBarry Smith   }
906*47c6ae99SBarry Smith   ierr = VecRestoreArray(U,&ww);CHKERRQ(ierr);
907*47c6ae99SBarry Smith   ierr = VecRestoreArray(a,&aa);CHKERRQ(ierr);
908*47c6ae99SBarry Smith   PetscFunctionReturn(0);
909*47c6ae99SBarry Smith }
910*47c6ae99SBarry Smith #endif
911*47c6ae99SBarry Smith 
912*47c6ae99SBarry Smith #if defined(PETSC_HAVE_ADIC)
913*47c6ae99SBarry Smith EXTERN_C_BEGIN
914*47c6ae99SBarry Smith #include "adic/ad_utils.h"
915*47c6ae99SBarry Smith EXTERN_C_END
916*47c6ae99SBarry Smith 
917*47c6ae99SBarry Smith #undef __FUNCT__
918*47c6ae99SBarry Smith #define __FUNCT__ "DAComputeJacobian1WithAdic"
919*47c6ae99SBarry Smith /*@C
920*47c6ae99SBarry Smith     DAComputeJacobian1WithAdic - Evaluates a adiC provided Jacobian function on each processor that
921*47c6ae99SBarry Smith         share a DA
922*47c6ae99SBarry Smith 
923*47c6ae99SBarry Smith    Input Parameters:
924*47c6ae99SBarry Smith +    da - the DA that defines the grid
925*47c6ae99SBarry Smith .    vu - input vector (ghosted)
926*47c6ae99SBarry Smith .    J - output matrix
927*47c6ae99SBarry Smith -    w - any user data
928*47c6ae99SBarry Smith 
929*47c6ae99SBarry Smith    Level: advanced
930*47c6ae99SBarry Smith 
931*47c6ae99SBarry Smith     Notes: Does NOT do ghost updates on vu upon entry
932*47c6ae99SBarry Smith 
933*47c6ae99SBarry Smith .seealso: DAFormFunction1()
934*47c6ae99SBarry Smith 
935*47c6ae99SBarry Smith @*/
936*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DAComputeJacobian1WithAdic(DA da,Vec vu,Mat J,void *w)
937*47c6ae99SBarry Smith {
938*47c6ae99SBarry Smith   PetscErrorCode ierr;
939*47c6ae99SBarry Smith   PetscInt       gtdof,tdof;
940*47c6ae99SBarry Smith   PetscScalar    *ustart;
941*47c6ae99SBarry Smith   DALocalInfo    info;
942*47c6ae99SBarry Smith   void           *ad_u,*ad_f,*ad_ustart,*ad_fstart;
943*47c6ae99SBarry Smith   ISColoring     iscoloring;
944*47c6ae99SBarry Smith 
945*47c6ae99SBarry Smith   PetscFunctionBegin;
946*47c6ae99SBarry Smith   ierr = DAGetLocalInfo(da,&info);CHKERRQ(ierr);
947*47c6ae99SBarry Smith 
948*47c6ae99SBarry Smith   PetscADResetIndep();
949*47c6ae99SBarry Smith 
950*47c6ae99SBarry Smith   /* get space for derivative objects.  */
951*47c6ae99SBarry Smith   ierr = DAGetAdicArray(da,PETSC_TRUE,&ad_u,&ad_ustart,&gtdof);CHKERRQ(ierr);
952*47c6ae99SBarry Smith   ierr = DAGetAdicArray(da,PETSC_FALSE,&ad_f,&ad_fstart,&tdof);CHKERRQ(ierr);
953*47c6ae99SBarry Smith   ierr = VecGetArray(vu,&ustart);CHKERRQ(ierr);
954*47c6ae99SBarry Smith   ierr = DAGetColoring(da,IS_COLORING_GHOSTED,MATAIJ,&iscoloring);CHKERRQ(ierr);
955*47c6ae99SBarry Smith 
956*47c6ae99SBarry Smith   PetscADSetValueAndColor(ad_ustart,gtdof,iscoloring->colors,ustart);
957*47c6ae99SBarry Smith 
958*47c6ae99SBarry Smith   ierr = VecRestoreArray(vu,&ustart);CHKERRQ(ierr);
959*47c6ae99SBarry Smith   ierr = ISColoringDestroy(iscoloring);CHKERRQ(ierr);
960*47c6ae99SBarry Smith   ierr = PetscADIncrementTotalGradSize(iscoloring->n);CHKERRQ(ierr);
961*47c6ae99SBarry Smith   PetscADSetIndepDone();
962*47c6ae99SBarry Smith 
963*47c6ae99SBarry Smith   ierr = PetscLogEventBegin(DA_LocalADFunction,0,0,0,0);CHKERRQ(ierr);
964*47c6ae99SBarry Smith   ierr = (*dd->adic_lf)(&info,ad_u,ad_f,w);CHKERRQ(ierr);
965*47c6ae99SBarry Smith   ierr = PetscLogEventEnd(DA_LocalADFunction,0,0,0,0);CHKERRQ(ierr);
966*47c6ae99SBarry Smith 
967*47c6ae99SBarry Smith   /* stick the values into the matrix */
968*47c6ae99SBarry Smith   ierr = MatSetValuesAdic(J,(PetscScalar**)ad_fstart);CHKERRQ(ierr);
969*47c6ae99SBarry Smith 
970*47c6ae99SBarry Smith   /* return space for derivative objects.  */
971*47c6ae99SBarry Smith   ierr = DARestoreAdicArray(da,PETSC_TRUE,&ad_u,&ad_ustart,&gtdof);CHKERRQ(ierr);
972*47c6ae99SBarry Smith   ierr = DARestoreAdicArray(da,PETSC_FALSE,&ad_f,&ad_fstart,&tdof);CHKERRQ(ierr);
973*47c6ae99SBarry Smith   PetscFunctionReturn(0);
974*47c6ae99SBarry Smith }
975*47c6ae99SBarry Smith 
976*47c6ae99SBarry Smith #undef __FUNCT__
977*47c6ae99SBarry Smith #define __FUNCT__ "DAMultiplyByJacobian1WithAdic"
978*47c6ae99SBarry Smith /*@C
979*47c6ae99SBarry Smith     DAMultiplyByJacobian1WithAdic - Applies an ADIC-provided Jacobian function to a vector on
980*47c6ae99SBarry Smith     each processor that shares a DA.
981*47c6ae99SBarry Smith 
982*47c6ae99SBarry Smith     Input Parameters:
983*47c6ae99SBarry Smith +   da - the DA that defines the grid
984*47c6ae99SBarry Smith .   vu - Jacobian is computed at this point (ghosted)
985*47c6ae99SBarry Smith .   v - product is done on this vector (ghosted)
986*47c6ae99SBarry Smith .   fu - output vector = J(vu)*v (not ghosted)
987*47c6ae99SBarry Smith -   w - any user data
988*47c6ae99SBarry Smith 
989*47c6ae99SBarry Smith     Notes:
990*47c6ae99SBarry Smith     This routine does NOT do ghost updates on vu upon entry.
991*47c6ae99SBarry Smith 
992*47c6ae99SBarry Smith    Level: advanced
993*47c6ae99SBarry Smith 
994*47c6ae99SBarry Smith .seealso: DAFormFunction1()
995*47c6ae99SBarry Smith 
996*47c6ae99SBarry Smith @*/
997*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DAMultiplyByJacobian1WithAdic(DA da,Vec vu,Vec v,Vec f,void *w)
998*47c6ae99SBarry Smith {
999*47c6ae99SBarry Smith   PetscErrorCode ierr;
1000*47c6ae99SBarry Smith   PetscInt       i,gtdof,tdof;
1001*47c6ae99SBarry Smith   PetscScalar    *avu,*av,*af,*ad_vustart,*ad_fstart;
1002*47c6ae99SBarry Smith   DALocalInfo    info;
1003*47c6ae99SBarry Smith   void           *ad_vu,*ad_f;
1004*47c6ae99SBarry Smith 
1005*47c6ae99SBarry Smith   PetscFunctionBegin;
1006*47c6ae99SBarry Smith   ierr = DAGetLocalInfo(da,&info);CHKERRQ(ierr);
1007*47c6ae99SBarry Smith 
1008*47c6ae99SBarry Smith   /* get space for derivative objects.  */
1009*47c6ae99SBarry Smith   ierr = DAGetAdicMFArray(da,PETSC_TRUE,&ad_vu,&ad_vustart,&gtdof);CHKERRQ(ierr);
1010*47c6ae99SBarry Smith   ierr = DAGetAdicMFArray(da,PETSC_FALSE,&ad_f,&ad_fstart,&tdof);CHKERRQ(ierr);
1011*47c6ae99SBarry Smith 
1012*47c6ae99SBarry Smith   /* copy input vector into derivative object */
1013*47c6ae99SBarry Smith   ierr = VecGetArray(vu,&avu);CHKERRQ(ierr);
1014*47c6ae99SBarry Smith   ierr = VecGetArray(v,&av);CHKERRQ(ierr);
1015*47c6ae99SBarry Smith   for (i=0; i<gtdof; i++) {
1016*47c6ae99SBarry Smith     ad_vustart[2*i]   = avu[i];
1017*47c6ae99SBarry Smith     ad_vustart[2*i+1] = av[i];
1018*47c6ae99SBarry Smith   }
1019*47c6ae99SBarry Smith   ierr = VecRestoreArray(vu,&avu);CHKERRQ(ierr);
1020*47c6ae99SBarry Smith   ierr = VecRestoreArray(v,&av);CHKERRQ(ierr);
1021*47c6ae99SBarry Smith 
1022*47c6ae99SBarry Smith   PetscADResetIndep();
1023*47c6ae99SBarry Smith   ierr = PetscADIncrementTotalGradSize(1);CHKERRQ(ierr);
1024*47c6ae99SBarry Smith   PetscADSetIndepDone();
1025*47c6ae99SBarry Smith 
1026*47c6ae99SBarry Smith   ierr = (*dd->adicmf_lf)(&info,ad_vu,ad_f,w);CHKERRQ(ierr);
1027*47c6ae99SBarry Smith 
1028*47c6ae99SBarry Smith   /* stick the values into the vector */
1029*47c6ae99SBarry Smith   ierr = VecGetArray(f,&af);CHKERRQ(ierr);
1030*47c6ae99SBarry Smith   for (i=0; i<tdof; i++) {
1031*47c6ae99SBarry Smith     af[i] = ad_fstart[2*i+1];
1032*47c6ae99SBarry Smith   }
1033*47c6ae99SBarry Smith   ierr = VecRestoreArray(f,&af);CHKERRQ(ierr);
1034*47c6ae99SBarry Smith 
1035*47c6ae99SBarry Smith   /* return space for derivative objects.  */
1036*47c6ae99SBarry Smith   ierr = DARestoreAdicMFArray(da,PETSC_TRUE,&ad_vu,&ad_vustart,&gtdof);CHKERRQ(ierr);
1037*47c6ae99SBarry Smith   ierr = DARestoreAdicMFArray(da,PETSC_FALSE,&ad_f,&ad_fstart,&tdof);CHKERRQ(ierr);
1038*47c6ae99SBarry Smith   PetscFunctionReturn(0);
1039*47c6ae99SBarry Smith }
1040*47c6ae99SBarry Smith #endif
1041*47c6ae99SBarry Smith 
1042*47c6ae99SBarry Smith #undef __FUNCT__
1043*47c6ae99SBarry Smith #define __FUNCT__ "DAComputeJacobian1"
1044*47c6ae99SBarry Smith /*@
1045*47c6ae99SBarry Smith     DAComputeJacobian1 - Evaluates a local Jacobian function on each processor that
1046*47c6ae99SBarry Smith         share a DA
1047*47c6ae99SBarry Smith 
1048*47c6ae99SBarry Smith    Input Parameters:
1049*47c6ae99SBarry Smith +    da - the DA that defines the grid
1050*47c6ae99SBarry Smith .    vu - input vector (ghosted)
1051*47c6ae99SBarry Smith .    J - output matrix
1052*47c6ae99SBarry Smith -    w - any user data
1053*47c6ae99SBarry Smith 
1054*47c6ae99SBarry Smith     Notes: Does NOT do ghost updates on vu upon entry
1055*47c6ae99SBarry Smith 
1056*47c6ae99SBarry Smith     Level: advanced
1057*47c6ae99SBarry Smith 
1058*47c6ae99SBarry Smith .seealso: DAFormFunction1()
1059*47c6ae99SBarry Smith 
1060*47c6ae99SBarry Smith @*/
1061*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DAComputeJacobian1(DA da,Vec vu,Mat J,void *w)
1062*47c6ae99SBarry Smith {
1063*47c6ae99SBarry Smith   PetscErrorCode ierr;
1064*47c6ae99SBarry Smith   void           *u;
1065*47c6ae99SBarry Smith   DALocalInfo    info;
1066*47c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
1067*47c6ae99SBarry Smith 
1068*47c6ae99SBarry Smith   PetscFunctionBegin;
1069*47c6ae99SBarry Smith   ierr = DAGetLocalInfo(da,&info);CHKERRQ(ierr);
1070*47c6ae99SBarry Smith   ierr = DAVecGetArray(da,vu,&u);CHKERRQ(ierr);
1071*47c6ae99SBarry Smith   ierr = (*dd->lj)(&info,u,J,w);CHKERRQ(ierr);
1072*47c6ae99SBarry Smith   ierr = DAVecRestoreArray(da,vu,&u);CHKERRQ(ierr);
1073*47c6ae99SBarry Smith   PetscFunctionReturn(0);
1074*47c6ae99SBarry Smith }
1075*47c6ae99SBarry Smith 
1076*47c6ae99SBarry Smith 
1077*47c6ae99SBarry Smith #undef __FUNCT__
1078*47c6ae99SBarry Smith #define __FUNCT__ "DAComputeJacobian1WithAdifor"
1079*47c6ae99SBarry Smith /*
1080*47c6ae99SBarry Smith     DAComputeJacobian1WithAdifor - Evaluates a ADIFOR provided Jacobian local function on each processor that
1081*47c6ae99SBarry Smith         share a DA
1082*47c6ae99SBarry Smith 
1083*47c6ae99SBarry Smith    Input Parameters:
1084*47c6ae99SBarry Smith +    da - the DA that defines the grid
1085*47c6ae99SBarry Smith .    vu - input vector (ghosted)
1086*47c6ae99SBarry Smith .    J - output matrix
1087*47c6ae99SBarry Smith -    w - any user data
1088*47c6ae99SBarry Smith 
1089*47c6ae99SBarry Smith     Notes: Does NOT do ghost updates on vu upon entry
1090*47c6ae99SBarry Smith 
1091*47c6ae99SBarry Smith .seealso: DAFormFunction1()
1092*47c6ae99SBarry Smith 
1093*47c6ae99SBarry Smith */
1094*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DAComputeJacobian1WithAdifor(DA da,Vec vu,Mat J,void *w)
1095*47c6ae99SBarry Smith {
1096*47c6ae99SBarry Smith   PetscErrorCode  ierr;
1097*47c6ae99SBarry Smith   PetscInt        i,Nc,N;
1098*47c6ae99SBarry Smith   ISColoringValue *color;
1099*47c6ae99SBarry Smith   DALocalInfo     info;
1100*47c6ae99SBarry Smith   PetscScalar     *u,*g_u,*g_f,*f = 0,*p_u;
1101*47c6ae99SBarry Smith   ISColoring      iscoloring;
1102*47c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
1103*47c6ae99SBarry Smith   void            (*lf)(PetscInt*,DALocalInfo*,PetscScalar*,PetscScalar*,PetscInt*,PetscScalar*,PetscScalar*,PetscInt*,void*,PetscErrorCode*) =
1104*47c6ae99SBarry Smith                   (void (*)(PetscInt*,DALocalInfo*,PetscScalar*,PetscScalar*,PetscInt*,PetscScalar*,PetscScalar*,PetscInt*,void*,PetscErrorCode*))*dd->adifor_lf;
1105*47c6ae99SBarry Smith 
1106*47c6ae99SBarry Smith   PetscFunctionBegin;
1107*47c6ae99SBarry Smith   ierr = DAGetColoring(da,IS_COLORING_GHOSTED,MATAIJ,&iscoloring);CHKERRQ(ierr);
1108*47c6ae99SBarry Smith   Nc   = iscoloring->n;
1109*47c6ae99SBarry Smith   ierr = DAGetLocalInfo(da,&info);CHKERRQ(ierr);
1110*47c6ae99SBarry Smith   N    = info.gxm*info.gym*info.gzm*info.dof;
1111*47c6ae99SBarry Smith 
1112*47c6ae99SBarry Smith   /* get space for derivative objects.  */
1113*47c6ae99SBarry Smith   ierr  = PetscMalloc(Nc*info.gxm*info.gym*info.gzm*info.dof*sizeof(PetscScalar),&g_u);CHKERRQ(ierr);
1114*47c6ae99SBarry Smith   ierr  = PetscMemzero(g_u,Nc*info.gxm*info.gym*info.gzm*info.dof*sizeof(PetscScalar));CHKERRQ(ierr);
1115*47c6ae99SBarry Smith   p_u   = g_u;
1116*47c6ae99SBarry Smith   color = iscoloring->colors;
1117*47c6ae99SBarry Smith   for (i=0; i<N; i++) {
1118*47c6ae99SBarry Smith     p_u[*color++] = 1.0;
1119*47c6ae99SBarry Smith     p_u          += Nc;
1120*47c6ae99SBarry Smith   }
1121*47c6ae99SBarry Smith   ierr = ISColoringDestroy(iscoloring);CHKERRQ(ierr);
1122*47c6ae99SBarry 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);
1123*47c6ae99SBarry Smith 
1124*47c6ae99SBarry Smith   /* Seed the input array g_u with coloring information */
1125*47c6ae99SBarry Smith 
1126*47c6ae99SBarry Smith   ierr = VecGetArray(vu,&u);CHKERRQ(ierr);
1127*47c6ae99SBarry Smith   (lf)(&Nc,&info,u,g_u,&Nc,f,g_f,&Nc,w,&ierr);CHKERRQ(ierr);
1128*47c6ae99SBarry Smith   ierr = VecRestoreArray(vu,&u);CHKERRQ(ierr);
1129*47c6ae99SBarry Smith 
1130*47c6ae99SBarry Smith   /* stick the values into the matrix */
1131*47c6ae99SBarry Smith   /* PetscScalarView(Nc*info.xm*info.ym,g_f,0); */
1132*47c6ae99SBarry Smith   ierr = MatSetValuesAdifor(J,Nc,g_f);CHKERRQ(ierr);
1133*47c6ae99SBarry Smith 
1134*47c6ae99SBarry Smith   /* return space for derivative objects.  */
1135*47c6ae99SBarry Smith   ierr = PetscFree(g_u);CHKERRQ(ierr);
1136*47c6ae99SBarry Smith   ierr = PetscFree2(g_f,f);CHKERRQ(ierr);
1137*47c6ae99SBarry Smith   PetscFunctionReturn(0);
1138*47c6ae99SBarry Smith }
1139*47c6ae99SBarry Smith 
1140*47c6ae99SBarry Smith #undef __FUNCT__
1141*47c6ae99SBarry Smith #define __FUNCT__ "DAFormJacobianLocal"
1142*47c6ae99SBarry Smith /*@C
1143*47c6ae99SBarry Smith    DAFormjacobianLocal - This is a universal Jacobian evaluation routine for
1144*47c6ae99SBarry Smith    a local DA function.
1145*47c6ae99SBarry Smith 
1146*47c6ae99SBarry Smith    Collective on DA
1147*47c6ae99SBarry Smith 
1148*47c6ae99SBarry Smith    Input Parameters:
1149*47c6ae99SBarry Smith +  da - the DA context
1150*47c6ae99SBarry Smith .  func - The local function
1151*47c6ae99SBarry Smith .  X - input vector
1152*47c6ae99SBarry Smith .  J - Jacobian matrix
1153*47c6ae99SBarry Smith -  ctx - A user context
1154*47c6ae99SBarry Smith 
1155*47c6ae99SBarry Smith    Level: intermediate
1156*47c6ae99SBarry Smith 
1157*47c6ae99SBarry Smith .seealso: DASetLocalFunction(), DASetLocalJacobian(), DASetLocalAdicFunction(), DASetLocalAdicMFFunction(),
1158*47c6ae99SBarry Smith           SNESSetFunction(), SNESSetJacobian()
1159*47c6ae99SBarry Smith 
1160*47c6ae99SBarry Smith @*/
1161*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DAFormJacobianLocal(DA da, DALocalFunction1 func, Vec X, Mat J, void *ctx)
1162*47c6ae99SBarry Smith {
1163*47c6ae99SBarry Smith   Vec            localX;
1164*47c6ae99SBarry Smith   DALocalInfo    info;
1165*47c6ae99SBarry Smith   void          *u;
1166*47c6ae99SBarry Smith   PetscErrorCode ierr;
1167*47c6ae99SBarry Smith 
1168*47c6ae99SBarry Smith   PetscFunctionBegin;
1169*47c6ae99SBarry Smith   ierr = DAGetLocalVector(da,&localX);CHKERRQ(ierr);
1170*47c6ae99SBarry Smith   /*
1171*47c6ae99SBarry Smith      Scatter ghost points to local vector, using the 2-step process
1172*47c6ae99SBarry Smith         DAGlobalToLocalBegin(), DAGlobalToLocalEnd().
1173*47c6ae99SBarry Smith   */
1174*47c6ae99SBarry Smith   ierr = DAGlobalToLocalBegin(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
1175*47c6ae99SBarry Smith   ierr = DAGlobalToLocalEnd(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
1176*47c6ae99SBarry Smith   ierr = DAGetLocalInfo(da,&info);CHKERRQ(ierr);
1177*47c6ae99SBarry Smith   ierr = DAVecGetArray(da,localX,&u);CHKERRQ(ierr);
1178*47c6ae99SBarry Smith   ierr = (*func)(&info,u,J,ctx);
1179*47c6ae99SBarry Smith   if (PetscExceptionValue(ierr)) {
1180*47c6ae99SBarry Smith     PetscErrorCode pierr = DAVecRestoreArray(da,localX,&u);CHKERRQ(pierr);
1181*47c6ae99SBarry Smith   }
1182*47c6ae99SBarry Smith   CHKERRQ(ierr);
1183*47c6ae99SBarry Smith   ierr = DAVecRestoreArray(da,localX,&u);CHKERRQ(ierr);
1184*47c6ae99SBarry Smith   if (PetscExceptionValue(ierr)) {
1185*47c6ae99SBarry Smith     PetscErrorCode pierr = DARestoreLocalVector(da,&localX);CHKERRQ(pierr);
1186*47c6ae99SBarry Smith   }
1187*47c6ae99SBarry Smith   CHKERRQ(ierr);
1188*47c6ae99SBarry Smith   ierr = DARestoreLocalVector(da,&localX);CHKERRQ(ierr);
1189*47c6ae99SBarry Smith   PetscFunctionReturn(0);
1190*47c6ae99SBarry Smith }
1191*47c6ae99SBarry Smith 
1192*47c6ae99SBarry Smith #undef __FUNCT__
1193*47c6ae99SBarry Smith #define __FUNCT__ "DAMultiplyByJacobian1WithAD"
1194*47c6ae99SBarry Smith /*@C
1195*47c6ae99SBarry Smith     DAMultiplyByJacobian1WithAD - Applies a Jacobian function supplied by ADIFOR or ADIC
1196*47c6ae99SBarry Smith     to a vector on each processor that shares a DA.
1197*47c6ae99SBarry Smith 
1198*47c6ae99SBarry Smith    Input Parameters:
1199*47c6ae99SBarry Smith +    da - the DA that defines the grid
1200*47c6ae99SBarry Smith .    vu - Jacobian is computed at this point (ghosted)
1201*47c6ae99SBarry Smith .    v - product is done on this vector (ghosted)
1202*47c6ae99SBarry Smith .    fu - output vector = J(vu)*v (not ghosted)
1203*47c6ae99SBarry Smith -    w - any user data
1204*47c6ae99SBarry Smith 
1205*47c6ae99SBarry Smith     Notes:
1206*47c6ae99SBarry Smith     This routine does NOT do ghost updates on vu and v upon entry.
1207*47c6ae99SBarry Smith 
1208*47c6ae99SBarry Smith     Automatically calls DAMultiplyByJacobian1WithAdifor() or DAMultiplyByJacobian1WithAdic()
1209*47c6ae99SBarry Smith     depending on whether DASetLocalAdicMFFunction() or DASetLocalAdiforMFFunction() was called.
1210*47c6ae99SBarry Smith 
1211*47c6ae99SBarry Smith    Level: advanced
1212*47c6ae99SBarry Smith 
1213*47c6ae99SBarry Smith .seealso: DAFormFunction1(), DAMultiplyByJacobian1WithAdifor(), DAMultiplyByJacobian1WithAdic()
1214*47c6ae99SBarry Smith 
1215*47c6ae99SBarry Smith @*/
1216*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DAMultiplyByJacobian1WithAD(DA da,Vec u,Vec v,Vec f,void *w)
1217*47c6ae99SBarry Smith {
1218*47c6ae99SBarry Smith   PetscErrorCode ierr;
1219*47c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
1220*47c6ae99SBarry Smith 
1221*47c6ae99SBarry Smith   PetscFunctionBegin;
1222*47c6ae99SBarry Smith   if (dd->adicmf_lf) {
1223*47c6ae99SBarry Smith #if defined(PETSC_HAVE_ADIC)
1224*47c6ae99SBarry Smith     ierr = DAMultiplyByJacobian1WithAdic(da,u,v,f,w);CHKERRQ(ierr);
1225*47c6ae99SBarry Smith #else
1226*47c6ae99SBarry Smith     SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP_SYS,"Requires ADIC to be installed and cannot use complex numbers");
1227*47c6ae99SBarry Smith #endif
1228*47c6ae99SBarry Smith   } else if (dd->adiformf_lf) {
1229*47c6ae99SBarry Smith     ierr = DAMultiplyByJacobian1WithAdifor(da,u,v,f,w);CHKERRQ(ierr);
1230*47c6ae99SBarry Smith   } else {
1231*47c6ae99SBarry Smith     SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ORDER,"Must call DASetLocalAdiforMFFunction() or DASetLocalAdicMFFunction() before using");
1232*47c6ae99SBarry Smith   }
1233*47c6ae99SBarry Smith   PetscFunctionReturn(0);
1234*47c6ae99SBarry Smith }
1235*47c6ae99SBarry Smith 
1236*47c6ae99SBarry Smith 
1237*47c6ae99SBarry Smith #undef __FUNCT__
1238*47c6ae99SBarry Smith #define __FUNCT__ "DAMultiplyByJacobian1WithAdifor"
1239*47c6ae99SBarry Smith /*@C
1240*47c6ae99SBarry Smith     DAMultiplyByJacobian1WithAdifor - Applies a ADIFOR provided Jacobian function on each processor that
1241*47c6ae99SBarry Smith         share a DA to a vector
1242*47c6ae99SBarry Smith 
1243*47c6ae99SBarry Smith    Input Parameters:
1244*47c6ae99SBarry Smith +    da - the DA that defines the grid
1245*47c6ae99SBarry Smith .    vu - Jacobian is computed at this point (ghosted)
1246*47c6ae99SBarry Smith .    v - product is done on this vector (ghosted)
1247*47c6ae99SBarry Smith .    fu - output vector = J(vu)*v (not ghosted)
1248*47c6ae99SBarry Smith -    w - any user data
1249*47c6ae99SBarry Smith 
1250*47c6ae99SBarry Smith     Notes: Does NOT do ghost updates on vu and v upon entry
1251*47c6ae99SBarry Smith 
1252*47c6ae99SBarry Smith    Level: advanced
1253*47c6ae99SBarry Smith 
1254*47c6ae99SBarry Smith .seealso: DAFormFunction1()
1255*47c6ae99SBarry Smith 
1256*47c6ae99SBarry Smith @*/
1257*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DAMultiplyByJacobian1WithAdifor(DA da,Vec u,Vec v,Vec f,void *w)
1258*47c6ae99SBarry Smith {
1259*47c6ae99SBarry Smith   PetscErrorCode ierr;
1260*47c6ae99SBarry Smith   PetscScalar    *au,*av,*af,*awork;
1261*47c6ae99SBarry Smith   Vec            work;
1262*47c6ae99SBarry Smith   DALocalInfo    info;
1263*47c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
1264*47c6ae99SBarry Smith   void           (*lf)(DALocalInfo*,PetscScalar*,PetscScalar*,PetscScalar*,PetscScalar*,void*,PetscErrorCode*) =
1265*47c6ae99SBarry Smith                  (void (*)(DALocalInfo*,PetscScalar*,PetscScalar*,PetscScalar*,PetscScalar*,void*,PetscErrorCode*))*dd->adiformf_lf;
1266*47c6ae99SBarry Smith 
1267*47c6ae99SBarry Smith   PetscFunctionBegin;
1268*47c6ae99SBarry Smith   ierr = DAGetLocalInfo(da,&info);CHKERRQ(ierr);
1269*47c6ae99SBarry Smith 
1270*47c6ae99SBarry Smith   ierr = DAGetGlobalVector(da,&work);CHKERRQ(ierr);
1271*47c6ae99SBarry Smith   ierr = VecGetArray(u,&au);CHKERRQ(ierr);
1272*47c6ae99SBarry Smith   ierr = VecGetArray(v,&av);CHKERRQ(ierr);
1273*47c6ae99SBarry Smith   ierr = VecGetArray(f,&af);CHKERRQ(ierr);
1274*47c6ae99SBarry Smith   ierr = VecGetArray(work,&awork);CHKERRQ(ierr);
1275*47c6ae99SBarry Smith   (lf)(&info,au,av,awork,af,w,&ierr);CHKERRQ(ierr);
1276*47c6ae99SBarry Smith   ierr = VecRestoreArray(u,&au);CHKERRQ(ierr);
1277*47c6ae99SBarry Smith   ierr = VecRestoreArray(v,&av);CHKERRQ(ierr);
1278*47c6ae99SBarry Smith   ierr = VecRestoreArray(f,&af);CHKERRQ(ierr);
1279*47c6ae99SBarry Smith   ierr = VecRestoreArray(work,&awork);CHKERRQ(ierr);
1280*47c6ae99SBarry Smith   ierr = DARestoreGlobalVector(da,&work);CHKERRQ(ierr);
1281*47c6ae99SBarry Smith 
1282*47c6ae99SBarry Smith   PetscFunctionReturn(0);
1283*47c6ae99SBarry Smith }
1284*47c6ae99SBarry Smith 
1285*47c6ae99SBarry Smith EXTERN_C_BEGIN
1286*47c6ae99SBarry Smith #undef __FUNCT__
1287*47c6ae99SBarry Smith #define __FUNCT__ "DACreate_2D"
1288*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DACreate_2D(DA da)
1289*47c6ae99SBarry Smith {
1290*47c6ae99SBarry Smith   DM_DA               *dd = (DM_DA*)da->data;
1291*47c6ae99SBarry Smith   const PetscInt       dim          = dd->dim;
1292*47c6ae99SBarry Smith   const PetscInt       M            = dd->M;
1293*47c6ae99SBarry Smith   const PetscInt       N            = dd->N;
1294*47c6ae99SBarry Smith   PetscInt             m            = dd->m;
1295*47c6ae99SBarry Smith   PetscInt             n            = dd->n;
1296*47c6ae99SBarry Smith   const PetscInt       dof          = dd->w;
1297*47c6ae99SBarry Smith   const PetscInt       s            = dd->s;
1298*47c6ae99SBarry Smith   const DAPeriodicType wrap         = dd->wrap;
1299*47c6ae99SBarry Smith   const DAStencilType  stencil_type = dd->stencil_type;
1300*47c6ae99SBarry Smith   PetscInt            *lx           = dd->lx;
1301*47c6ae99SBarry Smith   PetscInt            *ly           = dd->ly;
1302*47c6ae99SBarry Smith   MPI_Comm             comm;
1303*47c6ae99SBarry Smith   PetscMPIInt          rank,size;
1304*47c6ae99SBarry Smith   PetscInt             xs,xe,ys,ye,x,y,Xs,Xe,Ys,Ye,start,end;
1305*47c6ae99SBarry Smith   PetscInt             up,down,left,i,n0,n1,n2,n3,n5,n6,n7,n8,*idx,nn;
1306*47c6ae99SBarry Smith   PetscInt             xbase,*bases,*ldims,j,x_t,y_t,s_t,base,count;
1307*47c6ae99SBarry Smith   PetscInt             s_x,s_y; /* s proportionalized to w */
1308*47c6ae99SBarry Smith   PetscInt             sn0 = 0,sn2 = 0,sn6 = 0,sn8 = 0;
1309*47c6ae99SBarry Smith   Vec                  local,global;
1310*47c6ae99SBarry Smith   VecScatter           ltog,gtol;
1311*47c6ae99SBarry Smith   IS                   to,from;
1312*47c6ae99SBarry Smith   PetscErrorCode       ierr;
1313*47c6ae99SBarry Smith 
1314*47c6ae99SBarry Smith   PetscFunctionBegin;
1315*47c6ae99SBarry Smith #ifndef PETSC_USE_DYNAMIC_LIBRARIES
1316*47c6ae99SBarry Smith   ierr = DMInitializePackage(PETSC_NULL);CHKERRQ(ierr);
1317*47c6ae99SBarry Smith #endif
1318*47c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
1319*47c6ae99SBarry Smith 
1320*47c6ae99SBarry Smith   if (dim != PETSC_DECIDE && dim != 2) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Dimension should be 2: %D",dim);
1321*47c6ae99SBarry Smith   if (dof < 1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Must have 1 or more degrees of freedom per node: %D",dof);
1322*47c6ae99SBarry Smith   if (s < 0) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Stencil width cannot be negative: %D",s);
1323*47c6ae99SBarry Smith 
1324*47c6ae99SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1325*47c6ae99SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1326*47c6ae99SBarry Smith 
1327*47c6ae99SBarry Smith   da->ops->getelements = DAGetElements_2d_P1;
1328*47c6ae99SBarry Smith 
1329*47c6ae99SBarry Smith   dd->dim         = 2;
1330*47c6ae99SBarry Smith   dd->elementtype = DA_ELEMENT_P1;
1331*47c6ae99SBarry Smith   ierr = PetscMalloc(dof*sizeof(char*),&dd->fieldname);CHKERRQ(ierr);
1332*47c6ae99SBarry Smith   ierr = PetscMemzero(dd->fieldname,dof*sizeof(char*));CHKERRQ(ierr);
1333*47c6ae99SBarry Smith 
1334*47c6ae99SBarry Smith   if (m != PETSC_DECIDE) {
1335*47c6ae99SBarry Smith     if (m < 1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in X direction: %D",m);
1336*47c6ae99SBarry Smith     else if (m > size) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in X direction: %D %d",m,size);
1337*47c6ae99SBarry Smith   }
1338*47c6ae99SBarry Smith   if (n != PETSC_DECIDE) {
1339*47c6ae99SBarry Smith     if (n < 1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Y direction: %D",n);
1340*47c6ae99SBarry Smith     else if (n > size) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Y direction: %D %d",n,size);
1341*47c6ae99SBarry Smith   }
1342*47c6ae99SBarry Smith 
1343*47c6ae99SBarry Smith   if (m == PETSC_DECIDE || n == PETSC_DECIDE) {
1344*47c6ae99SBarry Smith     if (n != PETSC_DECIDE) {
1345*47c6ae99SBarry Smith       m = size/n;
1346*47c6ae99SBarry Smith     } else if (m != PETSC_DECIDE) {
1347*47c6ae99SBarry Smith       n = size/m;
1348*47c6ae99SBarry Smith     } else {
1349*47c6ae99SBarry Smith       /* try for squarish distribution */
1350*47c6ae99SBarry Smith       m = (PetscInt)(0.5 + sqrt(((double)M)*((double)size)/((double)N)));
1351*47c6ae99SBarry Smith       if (!m) m = 1;
1352*47c6ae99SBarry Smith       while (m > 0) {
1353*47c6ae99SBarry Smith 	n = size/m;
1354*47c6ae99SBarry Smith 	if (m*n == size) break;
1355*47c6ae99SBarry Smith 	m--;
1356*47c6ae99SBarry Smith       }
1357*47c6ae99SBarry Smith       if (M > N && m < n) {PetscInt _m = m; m = n; n = _m;}
1358*47c6ae99SBarry Smith     }
1359*47c6ae99SBarry 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 ");
1360*47c6ae99SBarry Smith   } else if (m*n != size) SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"Given Bad partition");
1361*47c6ae99SBarry Smith 
1362*47c6ae99SBarry Smith   if (M < m) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in x direction is too fine! %D %D",M,m);
1363*47c6ae99SBarry Smith   if (N < n) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"Partition in y direction is too fine! %D %D",N,n);
1364*47c6ae99SBarry Smith 
1365*47c6ae99SBarry Smith   /*
1366*47c6ae99SBarry Smith      Determine locally owned region
1367*47c6ae99SBarry Smith      xs is the first local node number, x is the number of local nodes
1368*47c6ae99SBarry Smith   */
1369*47c6ae99SBarry Smith   if (!lx) {
1370*47c6ae99SBarry Smith     ierr = PetscMalloc(m*sizeof(PetscInt), &dd->lx);CHKERRQ(ierr);
1371*47c6ae99SBarry Smith     lx = dd->lx;
1372*47c6ae99SBarry Smith     for (i=0; i<m; i++) {
1373*47c6ae99SBarry Smith       lx[i] = M/m + ((M % m) > i);
1374*47c6ae99SBarry Smith     }
1375*47c6ae99SBarry Smith   }
1376*47c6ae99SBarry Smith   x  = lx[rank % m];
1377*47c6ae99SBarry Smith   xs = 0;
1378*47c6ae99SBarry Smith   for (i=0; i<(rank % m); i++) {
1379*47c6ae99SBarry Smith     xs += lx[i];
1380*47c6ae99SBarry Smith   }
1381*47c6ae99SBarry Smith #if defined(PETSC_USE_DEBUG)
1382*47c6ae99SBarry Smith   left = xs;
1383*47c6ae99SBarry Smith   for (i=(rank % m); i<m; i++) {
1384*47c6ae99SBarry Smith     left += lx[i];
1385*47c6ae99SBarry Smith   }
1386*47c6ae99SBarry 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);
1387*47c6ae99SBarry Smith #endif
1388*47c6ae99SBarry Smith 
1389*47c6ae99SBarry Smith   /*
1390*47c6ae99SBarry Smith      Determine locally owned region
1391*47c6ae99SBarry Smith      ys is the first local node number, y is the number of local nodes
1392*47c6ae99SBarry Smith   */
1393*47c6ae99SBarry Smith   if (!ly) {
1394*47c6ae99SBarry Smith     ierr = PetscMalloc(n*sizeof(PetscInt), &dd->ly);CHKERRQ(ierr);
1395*47c6ae99SBarry Smith     ly = dd->ly;
1396*47c6ae99SBarry Smith     for (i=0; i<n; i++) {
1397*47c6ae99SBarry Smith       ly[i] = N/n + ((N % n) > i);
1398*47c6ae99SBarry Smith     }
1399*47c6ae99SBarry Smith   }
1400*47c6ae99SBarry Smith   y  = ly[rank/m];
1401*47c6ae99SBarry Smith   ys = 0;
1402*47c6ae99SBarry Smith   for (i=0; i<(rank/m); i++) {
1403*47c6ae99SBarry Smith     ys += ly[i];
1404*47c6ae99SBarry Smith   }
1405*47c6ae99SBarry Smith #if defined(PETSC_USE_DEBUG)
1406*47c6ae99SBarry Smith   left = ys;
1407*47c6ae99SBarry Smith   for (i=(rank/m); i<n; i++) {
1408*47c6ae99SBarry Smith     left += ly[i];
1409*47c6ae99SBarry Smith   }
1410*47c6ae99SBarry 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);
1411*47c6ae99SBarry Smith #endif
1412*47c6ae99SBarry Smith 
1413*47c6ae99SBarry Smith   if (x < s) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local x-width of domain x %D is smaller than stencil width s %D",x,s);
1414*47c6ae99SBarry Smith   if (y < s) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local y-width of domain y %D is smaller than stencil width s %D",y,s);
1415*47c6ae99SBarry Smith   xe = xs + x;
1416*47c6ae99SBarry Smith   ye = ys + y;
1417*47c6ae99SBarry Smith 
1418*47c6ae99SBarry Smith   /* determine ghost region */
1419*47c6ae99SBarry Smith   /* Assume No Periodicity */
1420*47c6ae99SBarry Smith   if (xs-s > 0) Xs = xs - s; else Xs = 0;
1421*47c6ae99SBarry Smith   if (ys-s > 0) Ys = ys - s; else Ys = 0;
1422*47c6ae99SBarry Smith   if (xe+s <= M) Xe = xe + s; else Xe = M;
1423*47c6ae99SBarry Smith   if (ye+s <= N) Ye = ye + s; else Ye = N;
1424*47c6ae99SBarry Smith 
1425*47c6ae99SBarry Smith   /* X Periodic */
1426*47c6ae99SBarry Smith   if (DAXPeriodic(wrap)){
1427*47c6ae99SBarry Smith     Xs = xs - s;
1428*47c6ae99SBarry Smith     Xe = xe + s;
1429*47c6ae99SBarry Smith   }
1430*47c6ae99SBarry Smith 
1431*47c6ae99SBarry Smith   /* Y Periodic */
1432*47c6ae99SBarry Smith   if (DAYPeriodic(wrap)){
1433*47c6ae99SBarry Smith     Ys = ys - s;
1434*47c6ae99SBarry Smith     Ye = ye + s;
1435*47c6ae99SBarry Smith   }
1436*47c6ae99SBarry Smith 
1437*47c6ae99SBarry Smith   /* Resize all X parameters to reflect w */
1438*47c6ae99SBarry Smith   x   *= dof;
1439*47c6ae99SBarry Smith   xs  *= dof;
1440*47c6ae99SBarry Smith   xe  *= dof;
1441*47c6ae99SBarry Smith   Xs  *= dof;
1442*47c6ae99SBarry Smith   Xe  *= dof;
1443*47c6ae99SBarry Smith   s_x = s*dof;
1444*47c6ae99SBarry Smith   s_y = s;
1445*47c6ae99SBarry Smith 
1446*47c6ae99SBarry Smith   /* determine starting point of each processor */
1447*47c6ae99SBarry Smith   nn    = x*y;
1448*47c6ae99SBarry Smith   ierr  = PetscMalloc2(size+1,PetscInt,&bases,size,PetscInt,&ldims);CHKERRQ(ierr);
1449*47c6ae99SBarry Smith   ierr  = MPI_Allgather(&nn,1,MPIU_INT,ldims,1,MPIU_INT,comm);CHKERRQ(ierr);
1450*47c6ae99SBarry Smith   bases[0] = 0;
1451*47c6ae99SBarry Smith   for (i=1; i<=size; i++) {
1452*47c6ae99SBarry Smith     bases[i] = ldims[i-1];
1453*47c6ae99SBarry Smith   }
1454*47c6ae99SBarry Smith   for (i=1; i<=size; i++) {
1455*47c6ae99SBarry Smith     bases[i] += bases[i-1];
1456*47c6ae99SBarry Smith   }
1457*47c6ae99SBarry Smith 
1458*47c6ae99SBarry Smith   /* allocate the base parallel and sequential vectors */
1459*47c6ae99SBarry Smith   dd->Nlocal = x*y;
1460*47c6ae99SBarry Smith   ierr = VecCreateMPIWithArray(comm,dd->Nlocal,PETSC_DECIDE,0,&global);CHKERRQ(ierr);
1461*47c6ae99SBarry Smith   ierr = VecSetBlockSize(global,dof);CHKERRQ(ierr);
1462*47c6ae99SBarry Smith   dd->nlocal = (Xe-Xs)*(Ye-Ys);
1463*47c6ae99SBarry Smith   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,dd->nlocal,0,&local);CHKERRQ(ierr);
1464*47c6ae99SBarry Smith   ierr = VecSetBlockSize(local,dof);CHKERRQ(ierr);
1465*47c6ae99SBarry Smith 
1466*47c6ae99SBarry Smith 
1467*47c6ae99SBarry Smith   /* generate appropriate vector scatters */
1468*47c6ae99SBarry Smith   /* local to global inserts non-ghost point region into global */
1469*47c6ae99SBarry Smith   ierr = VecGetOwnershipRange(global,&start,&end);CHKERRQ(ierr);
1470*47c6ae99SBarry Smith   ierr = ISCreateStride(comm,x*y,start,1,&to);CHKERRQ(ierr);
1471*47c6ae99SBarry Smith 
1472*47c6ae99SBarry Smith   left  = xs - Xs; down  = ys - Ys; up    = down + y;
1473*47c6ae99SBarry Smith   ierr = PetscMalloc(x*(up - down)*sizeof(PetscInt),&idx);CHKERRQ(ierr);
1474*47c6ae99SBarry Smith   count = 0;
1475*47c6ae99SBarry Smith   for (i=down; i<up; i++) {
1476*47c6ae99SBarry Smith     for (j=0; j<x/dof; j++) {
1477*47c6ae99SBarry Smith       idx[count++] = (left + i*(Xe-Xs) + j*dof)/dof;
1478*47c6ae99SBarry Smith     }
1479*47c6ae99SBarry Smith   }
1480*47c6ae99SBarry Smith   ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&from);CHKERRQ(ierr);
1481*47c6ae99SBarry Smith 
1482*47c6ae99SBarry Smith   ierr = VecScatterCreate(local,from,global,to,&ltog);CHKERRQ(ierr);
1483*47c6ae99SBarry Smith   ierr = PetscLogObjectParent(dd,ltog);CHKERRQ(ierr);
1484*47c6ae99SBarry Smith   ierr = ISDestroy(from);CHKERRQ(ierr);
1485*47c6ae99SBarry Smith   ierr = ISDestroy(to);CHKERRQ(ierr);
1486*47c6ae99SBarry Smith 
1487*47c6ae99SBarry Smith   /* global to local must include ghost points */
1488*47c6ae99SBarry Smith   if (stencil_type == DA_STENCIL_BOX) {
1489*47c6ae99SBarry Smith     ierr = ISCreateStride(comm,(Xe-Xs)*(Ye-Ys),0,1,&to);CHKERRQ(ierr);
1490*47c6ae99SBarry Smith   } else {
1491*47c6ae99SBarry Smith     /* must drop into cross shape region */
1492*47c6ae99SBarry Smith     /*       ---------|
1493*47c6ae99SBarry Smith             |  top    |
1494*47c6ae99SBarry Smith          |---         ---|
1495*47c6ae99SBarry Smith          |   middle      |
1496*47c6ae99SBarry Smith          |               |
1497*47c6ae99SBarry Smith          ----         ----
1498*47c6ae99SBarry Smith             | bottom  |
1499*47c6ae99SBarry Smith             -----------
1500*47c6ae99SBarry Smith         Xs xs        xe  Xe */
1501*47c6ae99SBarry Smith     /* bottom */
1502*47c6ae99SBarry Smith     left  = xs - Xs; down = ys - Ys; up    = down + y;
1503*47c6ae99SBarry Smith     count = down*(xe-xs) + (up-down)*(Xe-Xs) + (Ye-Ys-up)*(xe-xs);
1504*47c6ae99SBarry Smith     ierr  = PetscMalloc(count*sizeof(PetscInt)/dof,&idx);CHKERRQ(ierr);
1505*47c6ae99SBarry Smith     count = 0;
1506*47c6ae99SBarry Smith     for (i=0; i<down; i++) {
1507*47c6ae99SBarry Smith       for (j=0; j<xe-xs; j += dof) {
1508*47c6ae99SBarry Smith         idx[count++] = left + i*(Xe-Xs) + j;
1509*47c6ae99SBarry Smith       }
1510*47c6ae99SBarry Smith     }
1511*47c6ae99SBarry Smith     /* middle */
1512*47c6ae99SBarry Smith     for (i=down; i<up; i++) {
1513*47c6ae99SBarry Smith       for (j=0; j<Xe-Xs; j += dof) {
1514*47c6ae99SBarry Smith         idx[count++] = i*(Xe-Xs) + j;
1515*47c6ae99SBarry Smith       }
1516*47c6ae99SBarry Smith     }
1517*47c6ae99SBarry Smith     /* top */
1518*47c6ae99SBarry Smith     for (i=up; i<Ye-Ys; i++) {
1519*47c6ae99SBarry Smith       for (j=0; j<xe-xs; j += dof) {
1520*47c6ae99SBarry Smith         idx[count++] = (left + i*(Xe-Xs) + j);
1521*47c6ae99SBarry Smith       }
1522*47c6ae99SBarry Smith     }
1523*47c6ae99SBarry Smith     for (i=0; i<count; i++) idx[i] = idx[i]/dof;
1524*47c6ae99SBarry Smith     ierr = ISCreateBlock(comm,dof,count,idx,PETSC_OWN_POINTER,&to);CHKERRQ(ierr);
1525*47c6ae99SBarry Smith   }
1526*47c6ae99SBarry Smith 
1527*47c6ae99SBarry Smith 
1528*47c6ae99SBarry Smith   /* determine who lies on each side of us stored in    n6 n7 n8
1529*47c6ae99SBarry Smith                                                         n3    n5
1530*47c6ae99SBarry Smith                                                         n0 n1 n2
1531*47c6ae99SBarry Smith   */
1532*47c6ae99SBarry Smith 
1533*47c6ae99SBarry Smith   /* Assume the Non-Periodic Case */
1534*47c6ae99SBarry Smith   n1 = rank - m;
1535*47c6ae99SBarry Smith   if (rank % m) {
1536*47c6ae99SBarry Smith     n0 = n1 - 1;
1537*47c6ae99SBarry Smith   } else {
1538*47c6ae99SBarry Smith     n0 = -1;
1539*47c6ae99SBarry Smith   }
1540*47c6ae99SBarry Smith   if ((rank+1) % m) {
1541*47c6ae99SBarry Smith     n2 = n1 + 1;
1542*47c6ae99SBarry Smith     n5 = rank + 1;
1543*47c6ae99SBarry Smith     n8 = rank + m + 1; if (n8 >= m*n) n8 = -1;
1544*47c6ae99SBarry Smith   } else {
1545*47c6ae99SBarry Smith     n2 = -1; n5 = -1; n8 = -1;
1546*47c6ae99SBarry Smith   }
1547*47c6ae99SBarry Smith   if (rank % m) {
1548*47c6ae99SBarry Smith     n3 = rank - 1;
1549*47c6ae99SBarry Smith     n6 = n3 + m; if (n6 >= m*n) n6 = -1;
1550*47c6ae99SBarry Smith   } else {
1551*47c6ae99SBarry Smith     n3 = -1; n6 = -1;
1552*47c6ae99SBarry Smith   }
1553*47c6ae99SBarry Smith   n7 = rank + m; if (n7 >= m*n) n7 = -1;
1554*47c6ae99SBarry Smith 
1555*47c6ae99SBarry Smith 
1556*47c6ae99SBarry Smith   /* Modify for Periodic Cases */
1557*47c6ae99SBarry Smith   if (wrap == DA_YPERIODIC) {  /* Handle Top and Bottom Sides */
1558*47c6ae99SBarry Smith     if (n1 < 0) n1 = rank + m * (n-1);
1559*47c6ae99SBarry Smith     if (n7 < 0) n7 = rank - m * (n-1);
1560*47c6ae99SBarry Smith     if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1;
1561*47c6ae99SBarry Smith     if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1;
1562*47c6ae99SBarry Smith     if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1;
1563*47c6ae99SBarry Smith     if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1;
1564*47c6ae99SBarry Smith   } else if (wrap == DA_XPERIODIC) { /* Handle Left and Right Sides */
1565*47c6ae99SBarry Smith     if (n3 < 0) n3 = rank + (m-1);
1566*47c6ae99SBarry Smith     if (n5 < 0) n5 = rank - (m-1);
1567*47c6ae99SBarry Smith     if ((n1 >= 0) && (n0 < 0)) n0 = rank-1;
1568*47c6ae99SBarry Smith     if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1;
1569*47c6ae99SBarry Smith     if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1;
1570*47c6ae99SBarry Smith     if ((n7 >= 0) && (n8 < 0)) n8 = rank+1;
1571*47c6ae99SBarry Smith   } else if (wrap == DA_XYPERIODIC) {
1572*47c6ae99SBarry Smith 
1573*47c6ae99SBarry Smith     /* Handle all four corners */
1574*47c6ae99SBarry Smith     if ((n6 < 0) && (n7 < 0) && (n3 < 0)) n6 = m-1;
1575*47c6ae99SBarry Smith     if ((n8 < 0) && (n7 < 0) && (n5 < 0)) n8 = 0;
1576*47c6ae99SBarry Smith     if ((n2 < 0) && (n5 < 0) && (n1 < 0)) n2 = size-m;
1577*47c6ae99SBarry Smith     if ((n0 < 0) && (n3 < 0) && (n1 < 0)) n0 = size-1;
1578*47c6ae99SBarry Smith 
1579*47c6ae99SBarry Smith     /* Handle Top and Bottom Sides */
1580*47c6ae99SBarry Smith     if (n1 < 0) n1 = rank + m * (n-1);
1581*47c6ae99SBarry Smith     if (n7 < 0) n7 = rank - m * (n-1);
1582*47c6ae99SBarry Smith     if ((n3 >= 0) && (n0 < 0)) n0 = size - m + rank - 1;
1583*47c6ae99SBarry Smith     if ((n3 >= 0) && (n6 < 0)) n6 = (rank%m)-1;
1584*47c6ae99SBarry Smith     if ((n5 >= 0) && (n2 < 0)) n2 = size - m + rank + 1;
1585*47c6ae99SBarry Smith     if ((n5 >= 0) && (n8 < 0)) n8 = (rank%m)+1;
1586*47c6ae99SBarry Smith 
1587*47c6ae99SBarry Smith     /* Handle Left and Right Sides */
1588*47c6ae99SBarry Smith     if (n3 < 0) n3 = rank + (m-1);
1589*47c6ae99SBarry Smith     if (n5 < 0) n5 = rank - (m-1);
1590*47c6ae99SBarry Smith     if ((n1 >= 0) && (n0 < 0)) n0 = rank-1;
1591*47c6ae99SBarry Smith     if ((n1 >= 0) && (n2 < 0)) n2 = rank-2*m+1;
1592*47c6ae99SBarry Smith     if ((n7 >= 0) && (n6 < 0)) n6 = rank+2*m-1;
1593*47c6ae99SBarry Smith     if ((n7 >= 0) && (n8 < 0)) n8 = rank+1;
1594*47c6ae99SBarry Smith   }
1595*47c6ae99SBarry Smith   ierr = PetscMalloc(9*sizeof(PetscInt),&dd->neighbors);CHKERRQ(ierr);
1596*47c6ae99SBarry Smith   dd->neighbors[0] = n0;
1597*47c6ae99SBarry Smith   dd->neighbors[1] = n1;
1598*47c6ae99SBarry Smith   dd->neighbors[2] = n2;
1599*47c6ae99SBarry Smith   dd->neighbors[3] = n3;
1600*47c6ae99SBarry Smith   dd->neighbors[4] = rank;
1601*47c6ae99SBarry Smith   dd->neighbors[5] = n5;
1602*47c6ae99SBarry Smith   dd->neighbors[6] = n6;
1603*47c6ae99SBarry Smith   dd->neighbors[7] = n7;
1604*47c6ae99SBarry Smith   dd->neighbors[8] = n8;
1605*47c6ae99SBarry Smith 
1606*47c6ae99SBarry Smith   if (stencil_type == DA_STENCIL_STAR) {
1607*47c6ae99SBarry Smith     /* save corner processor numbers */
1608*47c6ae99SBarry Smith     sn0 = n0; sn2 = n2; sn6 = n6; sn8 = n8;
1609*47c6ae99SBarry Smith     n0 = n2 = n6 = n8 = -1;
1610*47c6ae99SBarry Smith   }
1611*47c6ae99SBarry Smith 
1612*47c6ae99SBarry Smith   ierr = PetscMalloc((x+2*s_x)*(y+2*s_y)*sizeof(PetscInt),&idx);CHKERRQ(ierr);
1613*47c6ae99SBarry Smith   ierr = PetscLogObjectMemory(da,(x+2*s_x)*(y+2*s_y)*sizeof(PetscInt));CHKERRQ(ierr);
1614*47c6ae99SBarry Smith   nn = 0;
1615*47c6ae99SBarry Smith 
1616*47c6ae99SBarry Smith   xbase = bases[rank];
1617*47c6ae99SBarry Smith   for (i=1; i<=s_y; i++) {
1618*47c6ae99SBarry Smith     if (n0 >= 0) { /* left below */
1619*47c6ae99SBarry Smith       x_t = lx[n0 % m]*dof;
1620*47c6ae99SBarry Smith       y_t = ly[(n0/m)];
1621*47c6ae99SBarry Smith       s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x;
1622*47c6ae99SBarry Smith       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1623*47c6ae99SBarry Smith     }
1624*47c6ae99SBarry Smith     if (n1 >= 0) { /* directly below */
1625*47c6ae99SBarry Smith       x_t = x;
1626*47c6ae99SBarry Smith       y_t = ly[(n1/m)];
1627*47c6ae99SBarry Smith       s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t;
1628*47c6ae99SBarry Smith       for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1629*47c6ae99SBarry Smith     }
1630*47c6ae99SBarry Smith     if (n2 >= 0) { /* right below */
1631*47c6ae99SBarry Smith       x_t = lx[n2 % m]*dof;
1632*47c6ae99SBarry Smith       y_t = ly[(n2/m)];
1633*47c6ae99SBarry Smith       s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t;
1634*47c6ae99SBarry Smith       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1635*47c6ae99SBarry Smith     }
1636*47c6ae99SBarry Smith   }
1637*47c6ae99SBarry Smith 
1638*47c6ae99SBarry Smith   for (i=0; i<y; i++) {
1639*47c6ae99SBarry Smith     if (n3 >= 0) { /* directly left */
1640*47c6ae99SBarry Smith       x_t = lx[n3 % m]*dof;
1641*47c6ae99SBarry Smith       /* y_t = y; */
1642*47c6ae99SBarry Smith       s_t = bases[n3] + (i+1)*x_t - s_x;
1643*47c6ae99SBarry Smith       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1644*47c6ae99SBarry Smith     }
1645*47c6ae99SBarry Smith 
1646*47c6ae99SBarry Smith     for (j=0; j<x; j++) { idx[nn++] = xbase++; } /* interior */
1647*47c6ae99SBarry Smith 
1648*47c6ae99SBarry Smith     if (n5 >= 0) { /* directly right */
1649*47c6ae99SBarry Smith       x_t = lx[n5 % m]*dof;
1650*47c6ae99SBarry Smith       /* y_t = y; */
1651*47c6ae99SBarry Smith       s_t = bases[n5] + (i)*x_t;
1652*47c6ae99SBarry Smith       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1653*47c6ae99SBarry Smith     }
1654*47c6ae99SBarry Smith   }
1655*47c6ae99SBarry Smith 
1656*47c6ae99SBarry Smith   for (i=1; i<=s_y; i++) {
1657*47c6ae99SBarry Smith     if (n6 >= 0) { /* left above */
1658*47c6ae99SBarry Smith       x_t = lx[n6 % m]*dof;
1659*47c6ae99SBarry Smith       /* y_t = ly[(n6/m)]; */
1660*47c6ae99SBarry Smith       s_t = bases[n6] + (i)*x_t - s_x;
1661*47c6ae99SBarry Smith       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1662*47c6ae99SBarry Smith     }
1663*47c6ae99SBarry Smith     if (n7 >= 0) { /* directly above */
1664*47c6ae99SBarry Smith       x_t = x;
1665*47c6ae99SBarry Smith       /* y_t = ly[(n7/m)]; */
1666*47c6ae99SBarry Smith       s_t = bases[n7] + (i-1)*x_t;
1667*47c6ae99SBarry Smith       for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1668*47c6ae99SBarry Smith     }
1669*47c6ae99SBarry Smith     if (n8 >= 0) { /* right above */
1670*47c6ae99SBarry Smith       x_t = lx[n8 % m]*dof;
1671*47c6ae99SBarry Smith       /* y_t = ly[(n8/m)]; */
1672*47c6ae99SBarry Smith       s_t = bases[n8] + (i-1)*x_t;
1673*47c6ae99SBarry Smith       for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1674*47c6ae99SBarry Smith     }
1675*47c6ae99SBarry Smith   }
1676*47c6ae99SBarry Smith 
1677*47c6ae99SBarry Smith   base = bases[rank];
1678*47c6ae99SBarry Smith   {
1679*47c6ae99SBarry Smith     PetscInt nnn = nn/dof,*iidx;
1680*47c6ae99SBarry Smith     ierr = PetscMalloc(nnn*sizeof(PetscInt),&iidx);CHKERRQ(ierr);
1681*47c6ae99SBarry Smith     for (i=0; i<nnn; i++) {
1682*47c6ae99SBarry Smith       iidx[i] = idx[dof*i]/dof;
1683*47c6ae99SBarry Smith     }
1684*47c6ae99SBarry Smith     ierr = ISCreateBlock(comm,dof,nnn,iidx,PETSC_OWN_POINTER,&from);CHKERRQ(ierr);
1685*47c6ae99SBarry Smith   }
1686*47c6ae99SBarry Smith   ierr = VecScatterCreate(global,from,local,to,&gtol);CHKERRQ(ierr);
1687*47c6ae99SBarry Smith   ierr = PetscLogObjectParent(da,gtol);CHKERRQ(ierr);
1688*47c6ae99SBarry Smith   ierr = ISDestroy(to);CHKERRQ(ierr);
1689*47c6ae99SBarry Smith   ierr = ISDestroy(from);CHKERRQ(ierr);
1690*47c6ae99SBarry Smith 
1691*47c6ae99SBarry Smith   if (stencil_type == DA_STENCIL_STAR) {
1692*47c6ae99SBarry Smith     /*
1693*47c6ae99SBarry Smith         Recompute the local to global mappings, this time keeping the
1694*47c6ae99SBarry Smith       information about the cross corner processor numbers.
1695*47c6ae99SBarry Smith     */
1696*47c6ae99SBarry Smith     n0 = sn0; n2 = sn2; n6 = sn6; n8 = sn8;
1697*47c6ae99SBarry Smith     nn = 0;
1698*47c6ae99SBarry Smith     xbase = bases[rank];
1699*47c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
1700*47c6ae99SBarry Smith       if (n0 >= 0) { /* left below */
1701*47c6ae99SBarry Smith         x_t = lx[n0 % m]*dof;
1702*47c6ae99SBarry Smith         y_t = ly[(n0/m)];
1703*47c6ae99SBarry Smith         s_t = bases[n0] + x_t*y_t - (s_y-i)*x_t - s_x;
1704*47c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1705*47c6ae99SBarry Smith       }
1706*47c6ae99SBarry Smith       if (n1 >= 0) { /* directly below */
1707*47c6ae99SBarry Smith         x_t = x;
1708*47c6ae99SBarry Smith         y_t = ly[(n1/m)];
1709*47c6ae99SBarry Smith         s_t = bases[n1] + x_t*y_t - (s_y+1-i)*x_t;
1710*47c6ae99SBarry Smith         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1711*47c6ae99SBarry Smith       }
1712*47c6ae99SBarry Smith       if (n2 >= 0) { /* right below */
1713*47c6ae99SBarry Smith         x_t = lx[n2 % m]*dof;
1714*47c6ae99SBarry Smith         y_t = ly[(n2/m)];
1715*47c6ae99SBarry Smith         s_t = bases[n2] + x_t*y_t - (s_y+1-i)*x_t;
1716*47c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1717*47c6ae99SBarry Smith       }
1718*47c6ae99SBarry Smith     }
1719*47c6ae99SBarry Smith 
1720*47c6ae99SBarry Smith     for (i=0; i<y; i++) {
1721*47c6ae99SBarry Smith       if (n3 >= 0) { /* directly left */
1722*47c6ae99SBarry Smith         x_t = lx[n3 % m]*dof;
1723*47c6ae99SBarry Smith         /* y_t = y; */
1724*47c6ae99SBarry Smith         s_t = bases[n3] + (i+1)*x_t - s_x;
1725*47c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1726*47c6ae99SBarry Smith       }
1727*47c6ae99SBarry Smith 
1728*47c6ae99SBarry Smith       for (j=0; j<x; j++) { idx[nn++] = xbase++; } /* interior */
1729*47c6ae99SBarry Smith 
1730*47c6ae99SBarry Smith       if (n5 >= 0) { /* directly right */
1731*47c6ae99SBarry Smith         x_t = lx[n5 % m]*dof;
1732*47c6ae99SBarry Smith         /* y_t = y; */
1733*47c6ae99SBarry Smith         s_t = bases[n5] + (i)*x_t;
1734*47c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1735*47c6ae99SBarry Smith       }
1736*47c6ae99SBarry Smith     }
1737*47c6ae99SBarry Smith 
1738*47c6ae99SBarry Smith     for (i=1; i<=s_y; i++) {
1739*47c6ae99SBarry Smith       if (n6 >= 0) { /* left above */
1740*47c6ae99SBarry Smith         x_t = lx[n6 % m]*dof;
1741*47c6ae99SBarry Smith         /* y_t = ly[(n6/m)]; */
1742*47c6ae99SBarry Smith         s_t = bases[n6] + (i)*x_t - s_x;
1743*47c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1744*47c6ae99SBarry Smith       }
1745*47c6ae99SBarry Smith       if (n7 >= 0) { /* directly above */
1746*47c6ae99SBarry Smith         x_t = x;
1747*47c6ae99SBarry Smith         /* y_t = ly[(n7/m)]; */
1748*47c6ae99SBarry Smith         s_t = bases[n7] + (i-1)*x_t;
1749*47c6ae99SBarry Smith         for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1750*47c6ae99SBarry Smith       }
1751*47c6ae99SBarry Smith       if (n8 >= 0) { /* right above */
1752*47c6ae99SBarry Smith         x_t = lx[n8 % m]*dof;
1753*47c6ae99SBarry Smith         /* y_t = ly[(n8/m)]; */
1754*47c6ae99SBarry Smith         s_t = bases[n8] + (i-1)*x_t;
1755*47c6ae99SBarry Smith         for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1756*47c6ae99SBarry Smith       }
1757*47c6ae99SBarry Smith     }
1758*47c6ae99SBarry Smith   }
1759*47c6ae99SBarry Smith   ierr = PetscFree2(bases,ldims);CHKERRQ(ierr);
1760*47c6ae99SBarry Smith 
1761*47c6ae99SBarry Smith   dd->m  = m;  dd->n  = n;
1762*47c6ae99SBarry Smith   dd->xs = xs; dd->xe = xe; dd->ys = ys; dd->ye = ye; dd->zs = 0; dd->ze = 1;
1763*47c6ae99SBarry Smith   dd->Xs = Xs; dd->Xe = Xe; dd->Ys = Ys; dd->Ye = Ye; dd->Zs = 0; dd->Ze = 1;
1764*47c6ae99SBarry Smith 
1765*47c6ae99SBarry Smith   ierr = VecDestroy(local);CHKERRQ(ierr);
1766*47c6ae99SBarry Smith   ierr = VecDestroy(global);CHKERRQ(ierr);
1767*47c6ae99SBarry Smith 
1768*47c6ae99SBarry Smith   dd->gtol      = gtol;
1769*47c6ae99SBarry Smith   dd->ltog      = ltog;
1770*47c6ae99SBarry Smith   dd->idx       = idx;
1771*47c6ae99SBarry Smith   dd->Nl        = nn;
1772*47c6ae99SBarry Smith   dd->base      = base;
1773*47c6ae99SBarry Smith   da->ops->view = DAView_2d;
1774*47c6ae99SBarry Smith 
1775*47c6ae99SBarry Smith   /*
1776*47c6ae99SBarry Smith      Set the local to global ordering in the global vector, this allows use
1777*47c6ae99SBarry Smith      of VecSetValuesLocal().
1778*47c6ae99SBarry Smith   */
1779*47c6ae99SBarry Smith   ierr = ISLocalToGlobalMappingCreate(comm,nn,idx,PETSC_OWN_POINTER,&dd->ltogmap);CHKERRQ(ierr);
1780*47c6ae99SBarry Smith   ierr = ISLocalToGlobalMappingBlock(dd->ltogmap,dd->w,&dd->ltogmapb);CHKERRQ(ierr);
1781*47c6ae99SBarry Smith   ierr = PetscLogObjectParent(da,dd->ltogmap);CHKERRQ(ierr);
1782*47c6ae99SBarry Smith 
1783*47c6ae99SBarry Smith   dd->ltol = PETSC_NULL;
1784*47c6ae99SBarry Smith   dd->ao   = PETSC_NULL;
1785*47c6ae99SBarry Smith 
1786*47c6ae99SBarry Smith   PetscFunctionReturn(0);
1787*47c6ae99SBarry Smith }
1788*47c6ae99SBarry Smith EXTERN_C_END
1789*47c6ae99SBarry Smith 
1790*47c6ae99SBarry Smith #undef __FUNCT__
1791*47c6ae99SBarry Smith #define __FUNCT__ "DACreate2d"
1792*47c6ae99SBarry Smith /*@C
1793*47c6ae99SBarry Smith    DACreate2d -  Creates an object that will manage the communication of  two-dimensional
1794*47c6ae99SBarry Smith    regular array data that is distributed across some processors.
1795*47c6ae99SBarry Smith 
1796*47c6ae99SBarry Smith    Collective on MPI_Comm
1797*47c6ae99SBarry Smith 
1798*47c6ae99SBarry Smith    Input Parameters:
1799*47c6ae99SBarry Smith +  comm - MPI communicator
1800*47c6ae99SBarry Smith .  wrap - type of periodicity should the array have.
1801*47c6ae99SBarry Smith          Use one of DA_NONPERIODIC, DA_XPERIODIC, DA_YPERIODIC, or DA_XYPERIODIC.
1802*47c6ae99SBarry Smith .  stencil_type - stencil type.  Use either DA_STENCIL_BOX or DA_STENCIL_STAR.
1803*47c6ae99SBarry 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
1804*47c6ae99SBarry Smith             from the command line with -da_grid_x <M> -da_grid_y <N>)
1805*47c6ae99SBarry Smith .  m,n - corresponding number of processors in each dimension
1806*47c6ae99SBarry Smith          (or PETSC_DECIDE to have calculated)
1807*47c6ae99SBarry Smith .  dof - number of degrees of freedom per node
1808*47c6ae99SBarry Smith .  s - stencil width
1809*47c6ae99SBarry Smith -  lx, ly - arrays containing the number of nodes in each cell along
1810*47c6ae99SBarry Smith            the x and y coordinates, or PETSC_NULL. If non-null, these
1811*47c6ae99SBarry Smith            must be of length as m and n, and the corresponding
1812*47c6ae99SBarry Smith            m and n cannot be PETSC_DECIDE. The sum of the lx[] entries
1813*47c6ae99SBarry Smith            must be M, and the sum of the ly[] entries must be N.
1814*47c6ae99SBarry Smith 
1815*47c6ae99SBarry Smith    Output Parameter:
1816*47c6ae99SBarry Smith .  da - the resulting distributed array object
1817*47c6ae99SBarry Smith 
1818*47c6ae99SBarry Smith    Options Database Key:
1819*47c6ae99SBarry Smith +  -da_view - Calls DAView() at the conclusion of DACreate2d()
1820*47c6ae99SBarry Smith .  -da_grid_x <nx> - number of grid points in x direction, if M < 0
1821*47c6ae99SBarry Smith .  -da_grid_y <ny> - number of grid points in y direction, if N < 0
1822*47c6ae99SBarry Smith .  -da_processors_x <nx> - number of processors in x direction
1823*47c6ae99SBarry Smith .  -da_processors_y <ny> - number of processors in y direction
1824*47c6ae99SBarry Smith .  -da_refine_x - refinement ratio in x direction
1825*47c6ae99SBarry Smith -  -da_refine_y - refinement ratio in y direction
1826*47c6ae99SBarry Smith 
1827*47c6ae99SBarry Smith    Level: beginner
1828*47c6ae99SBarry Smith 
1829*47c6ae99SBarry Smith    Notes:
1830*47c6ae99SBarry Smith    The stencil type DA_STENCIL_STAR with width 1 corresponds to the
1831*47c6ae99SBarry Smith    standard 5-pt stencil, while DA_STENCIL_BOX with width 1 denotes
1832*47c6ae99SBarry Smith    the standard 9-pt stencil.
1833*47c6ae99SBarry Smith 
1834*47c6ae99SBarry Smith    The array data itself is NOT stored in the DA, it is stored in Vec objects;
1835*47c6ae99SBarry Smith    The appropriate vector objects can be obtained with calls to DACreateGlobalVector()
1836*47c6ae99SBarry Smith    and DACreateLocalVector() and calls to VecDuplicate() if more are needed.
1837*47c6ae99SBarry Smith 
1838*47c6ae99SBarry Smith .keywords: distributed array, create, two-dimensional
1839*47c6ae99SBarry Smith 
1840*47c6ae99SBarry Smith .seealso: DADestroy(), DAView(), DACreate1d(), DACreate3d(), DAGlobalToLocalBegin(), DAGetRefinementFactor(),
1841*47c6ae99SBarry Smith           DAGlobalToLocalEnd(), DALocalToGlobal(), DALocalToLocalBegin(), DALocalToLocalEnd(), DASetRefinementFactor(),
1842*47c6ae99SBarry Smith           DAGetInfo(), DACreateGlobalVector(), DACreateLocalVector(), DACreateNaturalVector(), DALoad(), DAView(), DAGetOwnershipRanges()
1843*47c6ae99SBarry Smith 
1844*47c6ae99SBarry Smith @*/
1845*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DACreate2d(MPI_Comm comm,DAPeriodicType wrap,DAStencilType stencil_type,
1846*47c6ae99SBarry Smith                           PetscInt M,PetscInt N,PetscInt m,PetscInt n,PetscInt dof,PetscInt s,const PetscInt lx[],const PetscInt ly[],DA *da)
1847*47c6ae99SBarry Smith {
1848*47c6ae99SBarry Smith   PetscErrorCode ierr;
1849*47c6ae99SBarry Smith 
1850*47c6ae99SBarry Smith   PetscFunctionBegin;
1851*47c6ae99SBarry Smith   ierr = DACreate(comm, da);CHKERRQ(ierr);
1852*47c6ae99SBarry Smith   ierr = DASetDim(*da, 2);CHKERRQ(ierr);
1853*47c6ae99SBarry Smith   ierr = DASetSizes(*da, M, N, 1);CHKERRQ(ierr);
1854*47c6ae99SBarry Smith   ierr = DASetNumProcs(*da, m, n, PETSC_DECIDE);CHKERRQ(ierr);
1855*47c6ae99SBarry Smith   ierr = DASetPeriodicity(*da, wrap);CHKERRQ(ierr);
1856*47c6ae99SBarry Smith   ierr = DASetDof(*da, dof);CHKERRQ(ierr);
1857*47c6ae99SBarry Smith   ierr = DASetStencilType(*da, stencil_type);CHKERRQ(ierr);
1858*47c6ae99SBarry Smith   ierr = DASetStencilWidth(*da, s);CHKERRQ(ierr);
1859*47c6ae99SBarry Smith   ierr = DASetOwnershipRanges(*da, lx, ly, PETSC_NULL);CHKERRQ(ierr);
1860*47c6ae99SBarry Smith   /* This violates the behavior for other classes, but right now users expect negative dimensions to be handled this way */
1861*47c6ae99SBarry Smith   ierr = DASetFromOptions(*da);CHKERRQ(ierr);
1862*47c6ae99SBarry Smith   ierr = DASetType(*da, DA2D);CHKERRQ(ierr);
1863*47c6ae99SBarry Smith   PetscFunctionReturn(0);
1864*47c6ae99SBarry Smith }
1865