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: PetscInt mx, my; /* discretization in x- and y-directions */
67: PetscInt 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: PetscInt mx=10; /* discretization in x-direction */
89: PetscInt 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: TaoInt 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(PETSC_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(PETSC_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(PETSC_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,PETSC_TRUE); 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(PETSC_COMM_SELF,user.ndim,user.ndim,5,TAO_NULL,&H); CHKERRQ(info);
171: info = MatSetOption(H,MAT_SYMMETRIC,PETSC_TRUE); CHKERRQ(info);
173: info = TaoAppSetHessianRoutine(eptorsionapp,FormHessian,(void *)&user); CHKERRQ(info);
174: info = TaoAppSetHessianMat(eptorsionapp,H,H); CHKERRQ(info);
175: }
177: /* Check for any TAO command line options */
178: info = TaoSetOptions(eptorsionapp,tao); CHKERRQ(info);
181: /* Modify the PETSc KSP structure */
182: info = TaoAppGetKSP(eptorsionapp,&ksp); CHKERRQ(info);
183: if (ksp) {
184: info = KSPSetType(ksp,KSPCG); CHKERRQ(info);
185: }
187: /* SOLVE THE APPLICATION */
188: info = TaoSolveApplication(eptorsionapp,tao); CHKERRQ(info);
189: if (ksp) {
190: KSPView(ksp,PETSC_VIEWER_STDOUT_SELF); CHKERRQ(info);
191: }
193: /*
194: To View TAO solver information use
195: info = TaoView(tao); CHKERRQ(info);
196: */
198: /* Get information on termination */
199: info = TaoGetSolutionStatus(tao,&iter,&ff,&gnorm,0,0,&reason); CHKERRQ(info);
200: if (reason <= 0){
201: PetscPrintf(PETSC_COMM_WORLD,"Try a different TAO method, adjust some parameters, or check the function evaluation routines\n");
202: PetscPrintf(PETSC_COMM_WORLD,"Iter: %d, f: %4.2e, residual: %4.2e\n",iter,ff,gnorm); CHKERRQ(info);
203: }
205: /* Free TAO data structures */
206: info = TaoDestroy(tao); CHKERRQ(info);
207: info = TaoAppDestroy(eptorsionapp); CHKERRQ(info);
209: /* Free PETSc data structures */
210: info = VecDestroy(user.s); CHKERRQ(info);
211: info = VecDestroy(user.y); CHKERRQ(info);
212: info = VecDestroy(x); CHKERRQ(info);
213: info = MatDestroy(H); CHKERRQ(info);
215: /* Finalize TAO, PETSc */
216: TaoFinalize();
217: PetscFinalize();
219: return 0;
220: }
222: /* ------------------------------------------------------------------- */
225: /*
226: FormInitialGuess - Computes an initial approximation to the solution.
228: Input Parameters:
229: . user - user-defined application context
230: . X - vector
232: Output Parameters:
233: . X - vector
234: */
235: int FormInitialGuess(AppCtx *user,Vec X)
236: {
237: double hx = user->hx, hy = user->hy, temp;
238: PetscScalar val;
239: int info;
240: PetscInt i, j, k, nx = user->mx, ny = user->my;
242: /* Compute initial guess */
243: for (j=0; j<ny; j++) {
244: temp = TaoMin(j+1,ny-j)*hy;
245: for (i=0; i<nx; i++) {
246: k = nx*j + i;
247: val = TaoMin((TaoMin(i+1,nx-i))*hx,temp);
248: info = VecSetValues(X,1,&k,&val,ADD_VALUES); CHKERRQ(info);
249: }
250: }
252: return 0;
253: }
254: /* ------------------------------------------------------------------- */
257: /*
258: FormFunctionGradient - Evaluates the function and corresponding gradient.
259:
260: Input Parameters:
261: tao - the TAO_APPLICATION context
262: X - the input vector
263: ptr - optional user-defined context, as set by TaoSetFunction()
265: Output Parameters:
266: f - the newly evaluated function
267: G - the newly evaluated gradient
268: */
269: int FormFunctionGradient(TAO_APPLICATION tao,Vec X,double *f,Vec G,void *ptr)
270: {
271: int info;
272: info = FormFunction(tao,X,f,ptr);CHKERRQ(info);
273: info = FormGradient(tao,X,G,ptr);CHKERRQ(info);
274: return 0;
275: }
276: /* ------------------------------------------------------------------- */
279: /*
280: FormFunction - Evaluates the function, f(X).
282: Input Parameters:
283: . taoapp - the TAO_APPLICATION context
284: . X - the input vector
285: . ptr - optional user-defined context, as set by TaoSetFunction()
287: Output Parameters:
288: . f - the newly evaluated function
289: */
290: int FormFunction(TAO_APPLICATION taoapp,Vec X,double *f,void *ptr)
291: {
292: AppCtx *user = (AppCtx *) ptr;
293: double hx = user->hx, hy = user->hy, area, three = 3.0, p5 = 0.5;
294: double zero = 0.0, vb, vl, vr, vt, dvdx, dvdy, flin = 0.0, fquad = 0.0;
295: double v, cdiv3 = user->param/three;
296: PetscScalar *x;
297: int info;
298: PetscInt nx = user->mx, ny = user->my, i, j, k;
300: /* Get pointer to vector data */
301: info = VecGetArray(X,&x); CHKERRQ(info);
303: /* Compute function contributions over the lower triangular elements */
304: for (j=-1; j<ny; j++) {
305: for (i=-1; i<nx; i++) {
306: k = nx*j + i;
307: v = zero;
308: vr = zero;
309: vt = zero;
310: if (i >= 0 && j >= 0) v = x[k];
311: if (i < nx-1 && j > -1) vr = x[k+1];
312: if (i > -1 && j < ny-1) vt = x[k+nx];
313: dvdx = (vr-v)/hx;
314: dvdy = (vt-v)/hy;
315: fquad += dvdx*dvdx + dvdy*dvdy;
316: flin -= cdiv3*(v+vr+vt);
317: }
318: }
320: /* Compute function contributions over the upper triangular elements */
321: for (j=0; j<=ny; j++) {
322: for (i=0; i<=nx; i++) {
323: k = nx*j + i;
324: vb = zero;
325: vl = zero;
326: v = zero;
327: if (i < nx && j > 0) vb = x[k-nx];
328: if (i > 0 && j < ny) vl = x[k-1];
329: if (i < nx && j < ny) v = x[k];
330: dvdx = (v-vl)/hx;
331: dvdy = (v-vb)/hy;
332: fquad = fquad + dvdx*dvdx + dvdy*dvdy;
333: flin = flin - cdiv3*(vb+vl+v);
334: }
335: }
337: /* Restore vector */
338: info = VecRestoreArray(X,&x); CHKERRQ(info);
340: /* Scale the function */
341: area = p5*hx*hy;
342: *f = area*(p5*fquad+flin);
344: info = PetscLogFlops(nx*ny*24); CHKERRQ(info);
345: return 0;
346: }
347: /* ------------------------------------------------------------------- */
350: /*
351: FormGradient - Evaluates the gradient, G(X).
353: Input Parameters:
354: . taoapp - the TAO_APPLICATION context
355: . X - input vector
356: . ptr - optional user-defined context
357:
358: Output Parameters:
359: . G - vector containing the newly evaluated gradient
360: */
361: int FormGradient(TAO_APPLICATION taoapp,Vec X,Vec G,void *ptr)
362: {
363: AppCtx *user = (AppCtx *) ptr;
364: PetscScalar zero=0.0, p5=0.5, three = 3.0, area, val, *x;
365: int info;
366: PetscInt nx = user->mx, ny = user->my, ind, i, j, k;
367: double hx = user->hx, hy = user->hy;
368: double vb, vl, vr, vt, dvdx, dvdy;
369: double v, cdiv3 = user->param/three;
371: /* Initialize gradient to zero */
372: info = VecSet(G, zero); CHKERRQ(info);
374: /* Get pointer to vector data */
375: info = VecGetArray(X,&x); CHKERRQ(info);
377: /* Compute gradient contributions over the lower triangular elements */
378: for (j=-1; j<ny; j++) {
379: for (i=-1; i<nx; i++) {
380: k = nx*j + i;
381: v = zero;
382: vr = zero;
383: vt = zero;
384: if (i >= 0 && j >= 0) v = x[k];
385: if (i < nx-1 && j > -1) vr = x[k+1];
386: if (i > -1 && j < ny-1) vt = x[k+nx];
387: dvdx = (vr-v)/hx;
388: dvdy = (vt-v)/hy;
389: if (i != -1 && j != -1) {
390: ind = k; val = - dvdx/hx - dvdy/hy - cdiv3;
391: info = VecSetValues(G,1,&ind,&val,ADD_VALUES); CHKERRQ(info);
392: }
393: if (i != nx-1 && j != -1) {
394: ind = k+1; val = dvdx/hx - cdiv3;
395: info = VecSetValues(G,1,&ind,&val,ADD_VALUES); CHKERRQ(info);
396: }
397: if (i != -1 && j != ny-1) {
398: ind = k+nx; val = dvdy/hy - cdiv3;
399: info = VecSetValues(G,1,&ind,&val,ADD_VALUES); CHKERRQ(info);
400: }
401: }
402: }
404: /* Compute gradient contributions over the upper triangular elements */
405: for (j=0; j<=ny; j++) {
406: for (i=0; i<=nx; i++) {
407: k = nx*j + i;
408: vb = zero;
409: vl = zero;
410: v = zero;
411: if (i < nx && j > 0) vb = x[k-nx];
412: if (i > 0 && j < ny) vl = x[k-1];
413: if (i < nx && j < ny) v = x[k];
414: dvdx = (v-vl)/hx;
415: dvdy = (v-vb)/hy;
416: if (i != nx && j != 0) {
417: ind = k-nx; val = - dvdy/hy - cdiv3;
418: info = VecSetValues(G,1,&ind,&val,ADD_VALUES); CHKERRQ(info);
419: }
420: if (i != 0 && j != ny) {
421: ind = k-1; val = - dvdx/hx - cdiv3;
422: info = VecSetValues(G,1,&ind,&val,ADD_VALUES); CHKERRQ(info);
423: }
424: if (i != nx && j != ny) {
425: ind = k; val = dvdx/hx + dvdy/hy - cdiv3;
426: info = VecSetValues(G,1,&ind,&val,ADD_VALUES); CHKERRQ(info);
427: }
428: }
429: }
431: /* Restore vector */
432: info = VecRestoreArray(X,&x); CHKERRQ(info);
434: /* Assemble gradient vector */
435: info = VecAssemblyBegin(G); CHKERRQ(info);
436: info = VecAssemblyEnd(G); CHKERRQ(info);
438: /* Scale the gradient */
439: area = p5*hx*hy;
440: info = VecScale(G, area); CHKERRQ(info);
441:
442: info = PetscLogFlops(nx*ny*24); CHKERRQ(info);
443: return 0;
444: }
446: /* ------------------------------------------------------------------- */
449: /*
450: FormHessian - Forms the Hessian matrix.
452: Input Parameters:
453: . taoapp - the TAO_APPLICATION context
454: . X - the input vector
455: . ptr - optional user-defined context, as set by TaoSetHessian()
456:
457: Output Parameters:
458: . H - Hessian matrix
459: . PrecH - optionally different preconditioning Hessian
460: . flag - flag indicating matrix structure
462: Notes:
463: This routine is intended simply as an example of the interface
464: to a Hessian evaluation routine. Since this example compute the
465: Hessian a column at a time, it is not particularly efficient and
466: is not recommended.
467: */
468: int FormHessian(TAO_APPLICATION taoapp,Vec X,Mat *HH,Mat *Hpre, MatStructure *flg, void *ptr)
469: {
470: AppCtx *user = (AppCtx *) ptr;
471: int info;
472: PetscInt i,j, ndim = user->ndim;
473: PetscScalar *y, zero = 0.0, one = 1.0;
474: Mat H=*HH;
475: *Hpre = H;
476: PetscTruth assembled;
478: /* Set location of vector */
479: user->xvec = X;
481: /* Initialize Hessian entries and work vector to zero */
482: info = MatAssembled(H,&assembled); CHKERRQ(info);
483: if (assembled){info = MatZeroEntries(H); CHKERRQ(info);}
485: info = VecSet(user->s, zero); CHKERRQ(info);
487: /* Loop over matrix columns to compute entries of the Hessian */
488: for (j=0; j<ndim; j++) {
490: info = VecSetValues(user->s,1,&j,&one,INSERT_VALUES); CHKERRQ(info);
491: info = VecAssemblyBegin(user->s); CHKERRQ(info);
492: info = VecAssemblyEnd(user->s); CHKERRQ(info);
494: info = HessianProduct(ptr,user->s,user->y); CHKERRQ(info);
496: info = VecSetValues(user->s,1,&j,&zero,INSERT_VALUES); CHKERRQ(info);
497: info = VecAssemblyBegin(user->s); CHKERRQ(info);
498: info = VecAssemblyEnd(user->s); CHKERRQ(info);
500: info = VecGetArray(user->y,&y); CHKERRQ(info);
501: for (i=0; i<ndim; i++) {
502: if (y[i] != zero) {
503: info = MatSetValues(H,1,&i,1,&j,&y[i],ADD_VALUES); CHKERRQ(info);
504: }
505: }
506: info = VecRestoreArray(user->y,&y); CHKERRQ(info);
508: }
510: *flg=SAME_NONZERO_PATTERN;
512: /* Assemble matrix */
513: info = MatAssemblyBegin(H,MAT_FINAL_ASSEMBLY); CHKERRQ(info);
514: info = MatAssemblyEnd(H,MAT_FINAL_ASSEMBLY); CHKERRQ(info);
516: return 0;
517: }
520: /* ------------------------------------------------------------------- */
523: /*
524: MatrixFreeHessian - Sets a pointer for use in computing Hessian-vector
525: products.
526:
527: Input Parameters:
528: . taoapp - the TAO_APPLICATION context
529: . X - the input vector
530: . ptr - optional user-defined context, as set by TaoSetHessian()
531:
532: Output Parameters:
533: . H - Hessian matrix
534: . PrecH - optionally different preconditioning Hessian
535: . flag - flag indicating matrix structure
536: */
537: int MatrixFreeHessian(TAO_APPLICATION taoapp,Vec X,Mat *H,Mat *PrecH,
538: MatStructure *flag,void *ptr)
539: {
540: AppCtx *user = (AppCtx *) ptr;
542: /* Sets location of vector for use in computing matrix-vector products
543: of the form H(X)*y */
545: user->xvec = X;
546: return 0;
547: }
549: /* ------------------------------------------------------------------- */
552: /*
553: HessianProductMat - Computes the matrix-vector product
554: y = mat*svec.
556: Input Parameters:
557: . mat - input matrix
558: . svec - input vector
560: Output Parameters:
561: . y - solution vector
562: */
563: int HessianProductMat(Mat mat,Vec svec,Vec y)
564: {
565: void *ptr;
566: MatShellGetContext(mat,&ptr);
567: HessianProduct(ptr,svec,y);
569: return 0;
570: }
572: /* ------------------------------------------------------------------- */
575: /*
576: Hessian Product - Computes the matrix-vector product:
577: y = f''(x)*svec.
579: Input Parameters
580: . ptr - pointer to the user-defined context
581: . svec - input vector
583: Output Parameters:
584: . y - product vector
585: */
586: int HessianProduct(void *ptr,Vec svec,Vec y)
587: {
588: AppCtx *user = (AppCtx *)ptr;
589: PetscScalar p5 = 0.5, zero = 0.0, one = 1.0, hx, hy, val, area, *x, *s;
590: double v, vb, vl, vr, vt, hxhx, hyhy;
591: int info;
592: PetscInt nx, ny, i, j, k, ind;
594: nx = user->mx;
595: ny = user->my;
596: hx = user->hx;
597: hy = user->hy;
598: hxhx = one/(hx*hx);
599: hyhy = one/(hy*hy);
601: /* Get pointers to vector data */
602: info = VecGetArray(user->xvec,&x); CHKERRQ(info);
603: info = VecGetArray(svec,&s); CHKERRQ(info);
605: /* Initialize product vector to zero */
606: info = VecSet(y, zero); CHKERRQ(info);
608: /* Compute f''(x)*s over the lower triangular elements */
609: for (j=-1; j<ny; j++) {
610: for (i=-1; i<nx; i++) {
611: k = nx*j + i;
612: v = zero;
613: vr = zero;
614: vt = zero;
615: if (i != -1 && j != -1) v = s[k];
616: if (i != nx-1 && j != -1) {
617: vr = s[k+1];
618: ind = k+1; val = hxhx*(vr-v);
619: info = VecSetValues(y,1,&ind,&val,ADD_VALUES); CHKERRQ(info);
620: }
621: if (i != -1 && j != ny-1) {
622: vt = s[k+nx];
623: ind = k+nx; val = hyhy*(vt-v);
624: info = VecSetValues(y,1,&ind,&val,ADD_VALUES); CHKERRQ(info);
625: }
626: if (i != -1 && j != -1) {
627: ind = k; val = hxhx*(v-vr) + hyhy*(v-vt);
628: info = VecSetValues(y,1,&ind,&val,ADD_VALUES); CHKERRQ(info);
629: }
630: }
631: }
632:
633: /* Compute f''(x)*s over the upper triangular elements */
634: for (j=0; j<=ny; j++) {
635: for (i=0; i<=nx; i++) {
636: k = nx*j + i;
637: v = zero;
638: vl = zero;
639: vb = zero;
640: if (i != nx && j != ny) v = s[k];
641: if (i != nx && j != 0) {
642: vb = s[k-nx];
643: ind = k-nx; val = hyhy*(vb-v);
644: info = VecSetValues(y,1,&ind,&val,ADD_VALUES); CHKERRQ(info);
645: }
646: if (i != 0 && j != ny) {
647: vl = s[k-1];
648: ind = k-1; val = hxhx*(vl-v);
649: info = VecSetValues(y,1,&ind,&val,ADD_VALUES); CHKERRQ(info);
650: }
651: if (i != nx && j != ny) {
652: ind = k; val = hxhx*(v-vl) + hyhy*(v-vb);
653: info = VecSetValues(y,1,&ind,&val,ADD_VALUES); CHKERRQ(info);
654: }
655: }
656: }
657: /* Restore vector data */
658: info = VecRestoreArray(svec,&s); CHKERRQ(info);
659: info = VecRestoreArray(user->xvec,&x); CHKERRQ(info);
661: /* Assemble vector */
662: info = VecAssemblyBegin(y); CHKERRQ(info);
663: info = VecAssemblyEnd(y); CHKERRQ(info);
664:
665: /* Scale resulting vector by area */
666: area = p5*hx*hy;
667: info = VecScale(y, area); CHKERRQ(info);
669: info = PetscLogFlops(nx*ny*18); CHKERRQ(info);
670:
671: return 0;
672: }