Actual source code: eptorsion1.c
1: /*$Id$*/
3: /* Program usage: mpirun -np 1 eptorsion1 [-help] [all TAO options] */
5: /* ----------------------------------------------------------------------
7: Elastic-plastic torsion problem.
9: The elastic plastic torsion problem arises from the determination
10: of the stress field on an infinitely long cylindrical bar, which is
11: equivalent to the solution of the following problem:
13: min{ .5 * integral(||gradient(v(x))||^2 dx) - C * integral(v(x) dx)}
14:
15: where C is the torsion angle per unit length.
17: The multiprocessor version of this code is eptorsion2.c.
19: ---------------------------------------------------------------------- */
21: /*
22: Include "tao.h" so that we can use TAO solvers. Note that this
23: file automatically includes files for lower-level support, such as those
24: provided by the PETSc library:
25: petsc.h - base PETSc routines petscvec.h - vectors
26: petscsys.h - sysem routines petscmat.h - matrices
27: petscis.h - index sets petscksp.h - Krylov subspace methods
28: petscviewer.h - viewers petscpc.h - preconditioners
29: */
31: #include "tao.h"
34: static char help[]=
35: "Demonstrates use of the TAO package to solve \n\
36: unconstrained minimization problems on a single processor. This example \n\
37: is based on the Elastic-Plastic Torsion (dept) problem from the MINPACK-2 \n\
38: test suite.\n\
39: The command line options are:\n\
40: -mx <xg>, where <xg> = number of grid points in the 1st coordinate direction\n\
41: -my <yg>, where <yg> = number of grid points in the 2nd coordinate direction\n\
42: -par <param>, where <param> = angle of twist per unit length\n\n";
44: /*T
45: Concepts: TAO - Solving an unconstrained minimization problem
46: Routines: TaoInitialize(); TaoFinalize();
47: Routines: TaoApplicationCreate(); TaoAppDestroy();
48: Routines: TaoCreate(); TaoDestroy();
49: Routines: TaoAppSetObjectiveAndGradientRoutine();
50: Routines: TaoAppSetHessianMat(); TaoAppSetHessianRoutine();
51: Routines: TaoSetOptions();
52: Routines: TaoAppSetInitialSolutionVec();
53: Routines: TaoSolveApplication();
54: Routines: TaoGetSolutionStatus(); TaoAppGetKSP();
55: Processors: 1
56: T*/
58: /*
59: User-defined application context - contains data needed by the
60: application-provided call-back routines, FormFunction(),
61: FormGradient(), and FormHessian().
62: */
64: typedef struct {
65: PetscReal param; /* nonlinearity parameter */
66: int mx, my; /* discretization in x- and y-directions */
67: int ndim; /* problem dimension */
68: Vec s, y, xvec; /* work space for computing Hessian */
69: PetscReal hx, hy; /* mesh spacing in x- and y-directions */
70: } AppCtx;
72: /* -------- User-defined Routines --------- */
74: int FormInitialGuess(AppCtx*,Vec);
75: int FormFunction(TAO_APPLICATION,Vec,double*,void*);
76: int FormGradient(TAO_APPLICATION,Vec,Vec,void*);
77: int FormHessian(TAO_APPLICATION,Vec,Mat*,Mat*, MatStructure *,void*);
78: int HessianProductMat(Mat,Vec,Vec);
79: int HessianProduct(void*,Vec,Vec);
80: int MatrixFreeHessian(TAO_APPLICATION,Vec,Mat*,Mat*,MatStructure*,void*);
81: int FormFunctionGradient(TAO_APPLICATION,Vec,double *,Vec,void *);
85: int main(int argc,char **argv)
86: {
87: int info; /* used to check for functions returning nonzeros */
88: int mx=10; /* discretization in x-direction */
89: int my=10; /* discretization in y-direction */
90: Vec x; /* solution, gradient vectors */
91: PetscTruth flg; /* A return value when checking for use options */
92: TAO_SOLVER tao; /* TAO_SOLVER solver context */
93: TAO_APPLICATION eptorsionapp; /* The PETSc application */
94: Mat H; /* Hessian matrix */
95: int iter; /* iteration information */
96: double ff,gnorm;
97: TaoTerminateReason reason;
98: KSP ksp; /* PETSc Krylov subspace solver */
99: PC pc; /* PETSc preconditioner */
100: AppCtx user; /* application context */
101: int size; /* number of processes */
102: double one=1.0;
105: /* Initialize TAO,PETSc */
106: PetscInitialize(&argc,&argv,(char *)0,help);
107: TaoInitialize(&argc,&argv,(char *)0,help);
109: /* Optional: Check that only one processor is being used. */
110: info = MPI_Comm_size(MPI_COMM_WORLD,&size); CHKERRQ(info);
111: if (size >1) {
112: PetscPrintf(PETSC_COMM_SELF,"This example is intended for single processor use!\n");
113: PetscPrintf(PETSC_COMM_SELF,"Try the example eptorsion2!\n");
114: SETERRQ(1,"Incorrect number of processors");
115: }
117: /* Specify default parameters for the problem, check for command-line overrides */
118: user.param = 5.0;
119: info = PetscOptionsGetInt(TAO_NULL,"-my",&my,&flg); CHKERRQ(info);
120: info = PetscOptionsGetInt(TAO_NULL,"-mx",&mx,&flg); CHKERRQ(info);
121: info = PetscOptionsGetReal(TAO_NULL,"-par",&user.param,&flg); CHKERRQ(info);
124: PetscPrintf(PETSC_COMM_SELF,"\n---- Elastic-Plastic Torsion Problem -----\n");
125: PetscPrintf(PETSC_COMM_SELF,"mx: %d my: %d \n\n",mx,my);
126: user.ndim = mx * my; user.mx = mx; user.my = my;
128: user.hx = one/(mx+1); user.hy = one/(my+1);
131: /* Allocate vectors */
132: info = VecCreateSeq(MPI_COMM_SELF,user.ndim,&user.y); CHKERRQ(info);
133: info = VecDuplicate(user.y,&user.s); CHKERRQ(info);
134: info = VecDuplicate(user.y,&x); CHKERRQ(info);
136: /* The TAO code begins here */
138: /* Create TAO solver and set desired solution method */
139: info = TaoCreate(MPI_COMM_SELF,"tao_lmvm",&tao); CHKERRQ(info);
140: info = TaoApplicationCreate(PETSC_COMM_SELF,&eptorsionapp); CHKERRQ(info);
142: /* Set solution vector with an initial guess */
143: info = FormInitialGuess(&user,x); CHKERRQ(info);
144: info = TaoAppSetInitialSolutionVec(eptorsionapp,x); CHKERRQ(info);
146: /* Set routine for function and gradient evaluation */
147: info = TaoAppSetObjectiveAndGradientRoutine(eptorsionapp,FormFunctionGradient,(void *)&user); CHKERRQ(info);
149: /* From command line options, determine if using matrix-free hessian */
150: info = PetscOptionsHasName(TAO_NULL,"-my_tao_mf",&flg); CHKERRQ(info);
151: if (flg) {
152: info = MatCreateShell(MPI_COMM_SELF,user.ndim,user.ndim,user.ndim,
153: user.ndim,(void*)&user,&H); CHKERRQ(info);
154: info = MatShellSetOperation(H,MATOP_MULT,(void(*)())HessianProductMat); CHKERRQ
155: (info);
156: info = MatSetOption(H,MAT_SYMMETRIC); CHKERRQ(info);
158: info = TaoAppSetHessianRoutine(eptorsionapp,MatrixFreeHessian,(void *)&user); CHKERRQ(info);
159: info = TaoAppSetHessianMat(eptorsionapp,H,H); CHKERRQ(info);
161: /* Set null preconditioner. Alternatively, set user-provided
162: preconditioner or explicitly form preconditioning matrix */
163: info = TaoAppGetKSP(eptorsionapp,&ksp); CHKERRQ(info);
164: if (ksp){
165: info = KSPGetPC(ksp,&pc); CHKERRQ(info);
166: info = PCSetType(pc,PCNONE); CHKERRQ(info);
167: }
168: } else {
170: info = MatCreateSeqAIJ(MPI_COMM_SELF,user.ndim,user.ndim,5,TAO_NULL,&H); CHKERRQ(info);
171: info = MatSetOption(H,MAT_SYMMETRIC); CHKERRQ(info);
172: info = TaoAppSetHessianRoutine(eptorsionapp,FormHessian,(void *)&user); CHKERRQ(info);
173: info = TaoAppSetHessianMat(eptorsionapp,H,H); CHKERRQ(info);
174: }
176: /* Check for any TAO command line options */
177: info = TaoSetOptions(eptorsionapp,tao); CHKERRQ(info);
180: /* Modify the PETSc KSP structure */
181: info = TaoAppGetKSP(eptorsionapp,&ksp); CHKERRQ(info);
182: if (ksp) {
183: info = KSPSetType(ksp,KSPCG); CHKERRQ(info);
184: }
186: /* SOLVE THE APPLICATION */
187: info = TaoSolveApplication(eptorsionapp,tao); CHKERRQ(info);
188: if (ksp) {
189: KSPView(ksp,PETSC_VIEWER_STDOUT_SELF); CHKERRQ(info);
190: }
192: /*
193: To View TAO solver information use
194: info = TaoView(tao); CHKERRQ(info);
195: */
197: /* Get information on termination */
198: info = TaoGetSolutionStatus(tao,&iter,&ff,&gnorm,0,0,&reason); CHKERRQ(info);
199: if (reason <= 0){
200: PetscPrintf(MPI_COMM_WORLD,"Try a different TAO method, adjust some parameters, or check the function evaluation routines\n");
201: PetscPrintf(MPI_COMM_WORLD,"Iter: %d, f: %4.2e, residual: %4.2e\n",iter,ff,gnorm); CHKERRQ(info);
202: }
204: /* Free TAO data structures */
205: info = TaoDestroy(tao); CHKERRQ(info);
206: info = TaoAppDestroy(eptorsionapp); CHKERRQ(info);
208: /* Free PETSc data structures */
209: info = VecDestroy(user.s); CHKERRQ(info);
210: info = VecDestroy(user.y); CHKERRQ(info);
211: info = VecDestroy(x); CHKERRQ(info);
212: info = MatDestroy(H); CHKERRQ(info);
214: /* Finalize TAO, PETSc */
215: TaoFinalize();
216: PetscFinalize();
218: return 0;
219: }
221: /* ------------------------------------------------------------------- */
224: /*
225: FormInitialGuess - Computes an initial approximation to the solution.
227: Input Parameters:
228: . user - user-defined application context
229: . X - vector
231: Output Parameters:
232: . X - vector
233: */
234: int FormInitialGuess(AppCtx *user,Vec X)
235: {
236: double hx = user->hx, hy = user->hy, temp;
237: PetscScalar val;
238: int info, i, j, k, nx = user->mx, ny = user->my;
240: /* Compute initial guess */
241: for (j=0; j<ny; j++) {
242: temp = TaoMin(j+1,ny-j)*hy;
243: for (i=0; i<nx; i++) {
244: k = nx*j + i;
245: val = TaoMin((TaoMin(i+1,nx-i))*hx,temp);
246: info = VecSetValues(X,1,&k,&val,ADD_VALUES); CHKERRQ(info);
247: }
248: }
250: return 0;
251: }
252: /* ------------------------------------------------------------------- */
255: /*
256: FormFunctionGradient - Evaluates the function and corresponding gradient.
257:
258: Input Parameters:
259: tao - the TAO_APPLICATION context
260: X - the input vector
261: ptr - optional user-defined context, as set by TaoSetFunction()
263: Output Parameters:
264: f - the newly evaluated function
265: G - the newly evaluated gradient
266: */
267: int FormFunctionGradient(TAO_APPLICATION tao,Vec X,double *f,Vec G,void *ptr)
268: {
269: int info;
270: info = FormFunction(tao,X,f,ptr);CHKERRQ(info);
271: info = FormGradient(tao,X,G,ptr);CHKERRQ(info);
272: return 0;
273: }
274: /* ------------------------------------------------------------------- */
277: /*
278: FormFunction - Evaluates the function, f(X).
280: Input Parameters:
281: . taoapp - the TAO_APPLICATION context
282: . X - the input vector
283: . ptr - optional user-defined context, as set by TaoSetFunction()
285: Output Parameters:
286: . f - the newly evaluated function
287: */
288: int FormFunction(TAO_APPLICATION taoapp,Vec X,double *f,void *ptr)
289: {
290: AppCtx *user = (AppCtx *) ptr;
291: double hx = user->hx, hy = user->hy, area, three = 3.0, p5 = 0.5;
292: double zero = 0.0, vb, vl, vr, vt, dvdx, dvdy, flin = 0.0, fquad = 0.0;
293: double v, cdiv3 = user->param/three;
294: PetscScalar *x;
295: int info, nx = user->mx, ny = user->my, i, j, k;
297: /* Get pointer to vector data */
298: info = VecGetArray(X,&x); CHKERRQ(info);
300: /* Compute function contributions over the lower triangular elements */
301: for (j=-1; j<ny; j++) {
302: for (i=-1; i<nx; i++) {
303: k = nx*j + i;
304: v = zero;
305: vr = zero;
306: vt = zero;
307: if (i >= 0 && j >= 0) v = x[k];
308: if (i < nx-1 && j > -1) vr = x[k+1];
309: if (i > -1 && j < ny-1) vt = x[k+nx];
310: dvdx = (vr-v)/hx;
311: dvdy = (vt-v)/hy;
312: fquad += dvdx*dvdx + dvdy*dvdy;
313: flin -= cdiv3*(v+vr+vt);
314: }
315: }
317: /* Compute function contributions over the upper triangular elements */
318: for (j=0; j<=ny; j++) {
319: for (i=0; i<=nx; i++) {
320: k = nx*j + i;
321: vb = zero;
322: vl = zero;
323: v = zero;
324: if (i < nx && j > 0) vb = x[k-nx];
325: if (i > 0 && j < ny) vl = x[k-1];
326: if (i < nx && j < ny) v = x[k];
327: dvdx = (v-vl)/hx;
328: dvdy = (v-vb)/hy;
329: fquad = fquad + dvdx*dvdx + dvdy*dvdy;
330: flin = flin - cdiv3*(vb+vl+v);
331: }
332: }
334: /* Restore vector */
335: info = VecRestoreArray(X,&x); CHKERRQ(info);
337: /* Scale the function */
338: area = p5*hx*hy;
339: *f = area*(p5*fquad+flin);
341: info = PetscLogFlops(nx*ny*24); CHKERRQ(info);
342: return 0;
343: }
344: /* ------------------------------------------------------------------- */
347: /*
348: FormGradient - Evaluates the gradient, G(X).
350: Input Parameters:
351: . taoapp - the TAO_APPLICATION context
352: . X - input vector
353: . ptr - optional user-defined context
354:
355: Output Parameters:
356: . G - vector containing the newly evaluated gradient
357: */
358: int FormGradient(TAO_APPLICATION taoapp,Vec X,Vec G,void *ptr)
359: {
360: AppCtx *user = (AppCtx *) ptr;
361: PetscScalar zero=0.0, p5=0.5, three = 3.0, area, val, *x;
362: int info, nx = user->mx, ny = user->my, ind, i, j, k;
363: double hx = user->hx, hy = user->hy;
364: double vb, vl, vr, vt, dvdx, dvdy;
365: double v, cdiv3 = user->param/three;
367: /* Initialize gradient to zero */
368: info = VecSet(G, zero); CHKERRQ(info);
370: /* Get pointer to vector data */
371: info = VecGetArray(X,&x); CHKERRQ(info);
373: /* Compute gradient contributions over the lower triangular elements */
374: for (j=-1; j<ny; j++) {
375: for (i=-1; i<nx; i++) {
376: k = nx*j + i;
377: v = zero;
378: vr = zero;
379: vt = zero;
380: if (i >= 0 && j >= 0) v = x[k];
381: if (i < nx-1 && j > -1) vr = x[k+1];
382: if (i > -1 && j < ny-1) vt = x[k+nx];
383: dvdx = (vr-v)/hx;
384: dvdy = (vt-v)/hy;
385: if (i != -1 && j != -1) {
386: ind = k; val = - dvdx/hx - dvdy/hy - cdiv3;
387: info = VecSetValues(G,1,&ind,&val,ADD_VALUES); CHKERRQ(info);
388: }
389: if (i != nx-1 && j != -1) {
390: ind = k+1; val = dvdx/hx - cdiv3;
391: info = VecSetValues(G,1,&ind,&val,ADD_VALUES); CHKERRQ(info);
392: }
393: if (i != -1 && j != ny-1) {
394: ind = k+nx; val = dvdy/hy - cdiv3;
395: info = VecSetValues(G,1,&ind,&val,ADD_VALUES); CHKERRQ(info);
396: }
397: }
398: }
400: /* Compute gradient contributions over the upper triangular elements */
401: for (j=0; j<=ny; j++) {
402: for (i=0; i<=nx; i++) {
403: k = nx*j + i;
404: vb = zero;
405: vl = zero;
406: v = zero;
407: if (i < nx && j > 0) vb = x[k-nx];
408: if (i > 0 && j < ny) vl = x[k-1];
409: if (i < nx && j < ny) v = x[k];
410: dvdx = (v-vl)/hx;
411: dvdy = (v-vb)/hy;
412: if (i != nx && j != 0) {
413: ind = k-nx; val = - dvdy/hy - cdiv3;
414: info = VecSetValues(G,1,&ind,&val,ADD_VALUES); CHKERRQ(info);
415: }
416: if (i != 0 && j != ny) {
417: ind = k-1; val = - dvdx/hx - cdiv3;
418: info = VecSetValues(G,1,&ind,&val,ADD_VALUES); CHKERRQ(info);
419: }
420: if (i != nx && j != ny) {
421: ind = k; val = dvdx/hx + dvdy/hy - cdiv3;
422: info = VecSetValues(G,1,&ind,&val,ADD_VALUES); CHKERRQ(info);
423: }
424: }
425: }
427: /* Restore vector */
428: info = VecRestoreArray(X,&x); CHKERRQ(info);
430: /* Assemble gradient vector */
431: info = VecAssemblyBegin(G); CHKERRQ(info);
432: info = VecAssemblyEnd(G); CHKERRQ(info);
434: /* Scale the gradient */
435: area = p5*hx*hy;
436: info = VecScale(G, area); CHKERRQ(info);
437:
438: info = PetscLogFlops(nx*ny*24); CHKERRQ(info);
439: return 0;
440: }
442: /* ------------------------------------------------------------------- */
445: /*
446: FormHessian - Forms the Hessian matrix.
448: Input Parameters:
449: . taoapp - the TAO_APPLICATION context
450: . X - the input vector
451: . ptr - optional user-defined context, as set by TaoSetHessian()
452:
453: Output Parameters:
454: . H - Hessian matrix
455: . PrecH - optionally different preconditioning Hessian
456: . flag - flag indicating matrix structure
458: Notes:
459: This routine is intended simply as an example of the interface
460: to a Hessian evaluation routine. Since this example compute the
461: Hessian a column at a time, it is not particularly efficient and
462: is not recommended.
463: */
464: int FormHessian(TAO_APPLICATION taoapp,Vec X,Mat *HH,Mat *Hpre, MatStructure *flg, void *ptr)
465: {
466: AppCtx *user = (AppCtx *) ptr;
467: int i, j, info, ndim = user->ndim;
468: PetscScalar *y, zero = 0.0, one = 1.0;
469: Mat H=*HH;
470: *Hpre = H;
471: PetscTruth assembled;
473: /* Set location of vector */
474: user->xvec = X;
476: /* Initialize Hessian entries and work vector to zero */
477: info = MatAssembled(H,&assembled); CHKERRQ(info);
478: if (assembled){info = MatZeroEntries(H); CHKERRQ(info);}
480: info = VecSet(user->s, zero); CHKERRQ(info);
482: /* Loop over matrix columns to compute entries of the Hessian */
483: for (j=0; j<ndim; j++) {
485: info = VecSetValues(user->s,1,&j,&one,INSERT_VALUES); CHKERRQ(info);
486: info = VecAssemblyBegin(user->s); CHKERRQ(info);
487: info = VecAssemblyEnd(user->s); CHKERRQ(info);
489: info = HessianProduct(ptr,user->s,user->y); CHKERRQ(info);
491: info = VecSetValues(user->s,1,&j,&zero,INSERT_VALUES); CHKERRQ(info);
492: info = VecAssemblyBegin(user->s); CHKERRQ(info);
493: info = VecAssemblyEnd(user->s); CHKERRQ(info);
495: info = VecGetArray(user->y,&y); CHKERRQ(info);
496: for (i=0; i<ndim; i++) {
497: if (y[i] != zero) {
498: info = MatSetValues(H,1,&i,1,&j,&y[i],ADD_VALUES); CHKERRQ(info);
499: }
500: }
501: info = VecRestoreArray(user->y,&y); CHKERRQ(info);
503: }
505: *flg=SAME_NONZERO_PATTERN;
507: /* Assemble matrix */
508: info = MatAssemblyBegin(H,MAT_FINAL_ASSEMBLY); CHKERRQ(info);
509: info = MatAssemblyEnd(H,MAT_FINAL_ASSEMBLY); CHKERRQ(info);
511: return 0;
512: }
515: /* ------------------------------------------------------------------- */
518: /*
519: MatrixFreeHessian - Sets a pointer for use in computing Hessian-vector
520: products.
521:
522: Input Parameters:
523: . taoapp - the TAO_APPLICATION context
524: . X - the input vector
525: . ptr - optional user-defined context, as set by TaoSetHessian()
526:
527: Output Parameters:
528: . H - Hessian matrix
529: . PrecH - optionally different preconditioning Hessian
530: . flag - flag indicating matrix structure
531: */
532: int MatrixFreeHessian(TAO_APPLICATION taoapp,Vec X,Mat *H,Mat *PrecH,
533: MatStructure *flag,void *ptr)
534: {
535: AppCtx *user = (AppCtx *) ptr;
537: /* Sets location of vector for use in computing matrix-vector products
538: of the form H(X)*y */
540: user->xvec = X;
541: return 0;
542: }
544: /* ------------------------------------------------------------------- */
547: /*
548: HessianProductMat - Computes the matrix-vector product
549: y = mat*svec.
551: Input Parameters:
552: . mat - input matrix
553: . svec - input vector
555: Output Parameters:
556: . y - solution vector
557: */
558: int HessianProductMat(Mat mat,Vec svec,Vec y)
559: {
560: void *ptr;
561: MatShellGetContext(mat,&ptr);
562: HessianProduct(ptr,svec,y);
564: return 0;
565: }
567: /* ------------------------------------------------------------------- */
570: /*
571: Hessian Product - Computes the matrix-vector product:
572: y = f''(x)*svec.
574: Input Parameters
575: . ptr - pointer to the user-defined context
576: . svec - input vector
578: Output Parameters:
579: . y - product vector
580: */
581: int HessianProduct(void *ptr,Vec svec,Vec y)
582: {
583: AppCtx *user = (AppCtx *)ptr;
584: PetscScalar p5 = 0.5, zero = 0.0, one = 1.0, hx, hy, val, area, *x, *s;
585: double v, vb, vl, vr, vt, hxhx, hyhy;
586: int nx, ny, i, j, k, info, ind;
588: nx = user->mx;
589: ny = user->my;
590: hx = user->hx;
591: hy = user->hy;
592: hxhx = one/(hx*hx);
593: hyhy = one/(hy*hy);
595: /* Get pointers to vector data */
596: info = VecGetArray(user->xvec,&x); CHKERRQ(info);
597: info = VecGetArray(svec,&s); CHKERRQ(info);
599: /* Initialize product vector to zero */
600: info = VecSet(y, zero); CHKERRQ(info);
602: /* Compute f''(x)*s over the lower triangular elements */
603: for (j=-1; j<ny; j++) {
604: for (i=-1; i<nx; i++) {
605: k = nx*j + i;
606: v = zero;
607: vr = zero;
608: vt = zero;
609: if (i != -1 && j != -1) v = s[k];
610: if (i != nx-1 && j != -1) {
611: vr = s[k+1];
612: ind = k+1; val = hxhx*(vr-v);
613: info = VecSetValues(y,1,&ind,&val,ADD_VALUES); CHKERRQ(info);
614: }
615: if (i != -1 && j != ny-1) {
616: vt = s[k+nx];
617: ind = k+nx; val = hyhy*(vt-v);
618: info = VecSetValues(y,1,&ind,&val,ADD_VALUES); CHKERRQ(info);
619: }
620: if (i != -1 && j != -1) {
621: ind = k; val = hxhx*(v-vr) + hyhy*(v-vt);
622: info = VecSetValues(y,1,&ind,&val,ADD_VALUES); CHKERRQ(info);
623: }
624: }
625: }
626:
627: /* Compute f''(x)*s over the upper triangular elements */
628: for (j=0; j<=ny; j++) {
629: for (i=0; i<=nx; i++) {
630: k = nx*j + i;
631: v = zero;
632: vl = zero;
633: vb = zero;
634: if (i != nx && j != ny) v = s[k];
635: if (i != nx && j != 0) {
636: vb = s[k-nx];
637: ind = k-nx; val = hyhy*(vb-v);
638: info = VecSetValues(y,1,&ind,&val,ADD_VALUES); CHKERRQ(info);
639: }
640: if (i != 0 && j != ny) {
641: vl = s[k-1];
642: ind = k-1; val = hxhx*(vl-v);
643: info = VecSetValues(y,1,&ind,&val,ADD_VALUES); CHKERRQ(info);
644: }
645: if (i != nx && j != ny) {
646: ind = k; val = hxhx*(v-vl) + hyhy*(v-vb);
647: info = VecSetValues(y,1,&ind,&val,ADD_VALUES); CHKERRQ(info);
648: }
649: }
650: }
651: /* Restore vector data */
652: info = VecRestoreArray(svec,&s); CHKERRQ(info);
653: info = VecRestoreArray(user->xvec,&x); CHKERRQ(info);
655: /* Assemble vector */
656: info = VecAssemblyBegin(y); CHKERRQ(info);
657: info = VecAssemblyEnd(y); CHKERRQ(info);
658:
659: /* Scale resulting vector by area */
660: area = p5*hx*hy;
661: info = VecScale(y, area); CHKERRQ(info);
663: PetscLogFlops(nx*ny*18);
664:
665: return 0;
666: }