Actual source code: combustion3.c
1: /*$Id: combustion.c,v 1.1 2002/08/08 9:39 lopezca@mauddib.mcs.anl.gov $*/
3: /* Program usage: mpirun -np <proc> combustion [-help] [all TAO options] */
5: /*
6: Include "tao.h" so we can use TAO solvers.
7: petscda.h for distributed array
8: ad_deriv.h for AD gradient
9: */
11: #include "petscda.h"
12: #include "tao.h"
13: #include "taodaapplication.h"
15: static char help[] = "Steady-State Combustion.\n\
16: We solve the Steady-State Combustion problem (MINPACK-2 test suite) in a 2D\n\
17: rectangular domain, using distributed arrays (DAs) to partition the parallel grid.\n\
18: The command line options include:\n\
19: -mx <xg>, where <xg> = number of grid points in the 1st coordinate direction\n\
20: -my <yg>, where <yg> = number of grid points in the 2nd coordinate direction\n\
21: -nlevels <nlevels>, where <nlevels> = number of levels in multigrid\n\
22: -byelement, if computation is made by functions on rectangular elements\n\
23: -adic, if AD is used (AD is not used by default)\n\
24: -par <parameter>, where <parameter> indicates the problem's nonlinearity\n\
25: parameter lambda (0 <= par <= 6.81)\n\n";
27: /*T
28: Concepts: TAO - Solving a bounded minimization problem
29: Routines: TaoInitialize(); TaoFinalize();
30: Routines: TaoCreate(); TaoDestroy();
31: Routines: DAApplicationCreate(); DAApplicationDestroy();
32: Routines: DAAppSetVariableBoundsRoutine();
33: Routines: DAAppSetElementObjectiveAndGradientRoutine();
34: Routines: DAAppSetElementHessianRoutine();
35: Routines: DAAppSetObjectiveAndGradientRoutine();
36: Routines: DAAppSetADElementFunctionGradient();
37: Routines: DAAppSetHessianRoutine();
38: Routines: TaoAppSetOptions();
39: Routines: TaoGetSolutionStatus(); TaoDAAppSolve();
40: Routines: DAAppSetMonitor(); TaoView();
41: Routines: DAAppGetSolution();
42: Routines: DAAppGetInterpolationMatrix();
43: Processors: n
44: T*/
46: /*
47: User-defined application context - contains data needed by the
48: application-provided call-back routines.
49: */
51: typedef struct {
52: InactiveDouble param;
53: InactiveDouble hx, hy; /* increment size in both directions */
54: InactiveDouble area; /* area of the triangles */
55: } ADFGCtx;
57: typedef struct {
58: PetscReal param; /* nonlinearity parameter */
59: double hx, hy, area; /* increments and area of the triangle */
60: int mx, my; /* discretization including boundaries */
61: ADFGCtx fgctx; /* Used only when an ADIC generated gradient is used */
62: } AppCtx;
63: int ad_CombLocalFunction(int[2], DERIV_TYPE[4], DERIV_TYPE*, void*);
65: /* User-defined routines foun in this file */
66: static int AppCtxInitialize(void *ptr);
67: static int FormInitialGuess(DA, Vec, AppCtx*);
69: static int CombLocalFunctionGradient(int[3], double x[4], double *f, double g[4], void *ptr);
70: static int WholeCombFunctionGradient(TAO_APPLICATION,DA,Vec,double *,Vec,void*);
72: static int CombLocalHessian(int[3], double x[4], double H[4][4], void *ptr);
73: static int WholeCombHessian(TAO_APPLICATION,DA,Vec,Mat,void*);
75: static int DAFixBoundary(TAO_APPLICATION, DA, Vec, Vec, void*);
77: static int MyGridMonitorBefore(TAO_APPLICATION, DA, int, void *);
82: int main( int argc, char **argv ) {
84: int info,iter; /* used to check for functions returning nonzeros */
85: int Nx,Ny;
86: int nlevels; /* multigrid levels */
87: double ff,gnorm;
88: DA DAarray[20];
89: Vec X;
90: KSP ksp;
91: PetscTruth flg, PreLoad = PETSC_TRUE; /* flags */
92: AppCtx user; /* user-defined work context */
93: TaoMethod method = "tao_tron"; /* minimization method */
94: TAO_SOLVER tao; /* TAO_SOLVER solver context */
95: TAO_APPLICATION CombApp; /* The PETSc application */
96: TaoTerminateReason reason;
98: /* Initialize TAO */
99: PetscInitialize(&argc, &argv, (char *)0, help);
100: TaoInitialize(&argc, &argv, (char *)0, help);
102: PreLoadBegin(PreLoad,"Solve");
103:
104: info = AppCtxInitialize((void*)&user); CHKERRQ(info);
106: nlevels=5;
107: info = PetscOptionsGetInt(PETSC_NULL,"-nlevels",&nlevels,&flg); CHKERRQ(info);
108: if (PreLoadIt == 0) {
109: nlevels = 1; user.mx = 11; user.my = 11;}
111: PetscPrintf(MPI_COMM_WORLD,"\n---- Steady-State Combustion Problem -----\n\n");
113: /* Let PETSc determine the vector distribution */
114: Nx = PETSC_DECIDE; Ny = PETSC_DECIDE;
116: /* Create distributed array (DA) to manage parallel grid and vectors */
117: info = DACreate2d(PETSC_COMM_WORLD,DA_NONPERIODIC,DA_STENCIL_BOX,user.mx,
118: user.my,Nx,Ny,1,1,PETSC_NULL,PETSC_NULL,&DAarray[0]); CHKERRQ(info);
119: for (iter=1;iter<nlevels;iter++){
120: info = DARefine(DAarray[iter-1],PETSC_COMM_WORLD,&DAarray[iter]); CHKERRQ(info);
121: }
123: /* Create TAO solver and set desired solution method */
124: info = TaoCreate(MPI_COMM_WORLD,method,&tao); CHKERRQ(info);
126: info = TaoApplicationCreate(PETSC_COMM_WORLD, &CombApp); CHKERRQ(info);
127: info = TaoAppSetDAApp(CombApp,DAarray,nlevels); CHKERRQ(info);
129: /* Sets routines for function, gradient and bounds evaluation */
130: info = DAAppSetVariableBoundsRoutine(CombApp,DAFixBoundary,(void *)&user); CHKERRQ(info);
132: info = PetscOptionsHasName(TAO_NULL, "-byelement", &flg); CHKERRQ(info);
133: if (flg) {
135: /* Sets routines for function and gradient evaluation, element by element */
136: info = PetscOptionsHasName(TAO_NULL, "-adic", &flg); CHKERRQ(info);
137: if (flg) {
138: info = DAAppSetADElementFunctionGradient(CombApp,ad_CombLocalFunction,228,(void *)&user.fgctx); CHKERRQ(info);
139: } else {
140: info = DAAppSetElementObjectiveAndGradientRoutine(CombApp,CombLocalFunctionGradient,51,(void *)&user); CHKERRQ(info);
141: }
142: /* Sets routines for Hessian evaluation, element by element */
143: info = DAAppSetElementHessianRoutine(CombApp,CombLocalHessian,21,(void*)&user); CHKERRQ(info);
145: } else {
147: /* Sets routines for function and gradient evaluation, all in one routine */
148: info = DAAppSetObjectiveAndGradientRoutine(CombApp,WholeCombFunctionGradient,(void *)&user); CHKERRQ(info);
150: /* Sets routines for Hessian evaluation, all in one routine */
151: info = DAAppSetHessianRoutine(CombApp,WholeCombHessian,(void*)&user); CHKERRQ(info);
152:
153: }
155: info = DAAppSetBeforeMonitor(CombApp,MyGridMonitorBefore,(void*)&user); CHKERRQ(info);
156: info = DAAppPrintInterpolationError(CombApp); CHKERRQ(info);
157: // info = DAAppPrintStageTimes(CombApp); CHKERRQ(info);
159: info = TaoAppSetRelativeTolerance(CombApp,1.0e-6); CHKERRQ(info);
160: info = TaoSetTolerances(tao,0,0,0,0); CHKERRQ(info);
161: info = TaoSetGradientTolerances(tao,0,0,0); CHKERRQ(info);
163: info = TaoAppGetKSP(CombApp,&ksp);CHKERRQ(info);
164: info = KSPSetType(ksp,KSPCG); CHKERRQ(info);
165: info = KSPSetTolerances(ksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,100);CHKERRQ(info);
167: /* Check for any tao command line options */
168: info = TaoSetOptions(CombApp, tao); CHKERRQ(info);
170: info = DAAppGetSolution(CombApp,0,&X); CHKERRQ(info);
171: info = FormInitialGuess(DAarray[0],X,&user); CHKERRQ(info);
172: info = DAAppSetInitialSolution(CombApp,X); CHKERRQ(info);
174: /* SOLVE THE APPLICATION */
175: info = TaoDAAppSolve(CombApp, tao); CHKERRQ(info);
177: /* Get information on termination */
178: info = TaoGetSolutionStatus(tao,&iter,&ff,&gnorm,0,0,&reason); CHKERRQ(info);
179: if (reason <= 0 ){
180: PetscPrintf(MPI_COMM_WORLD,"Try a different TAO method, adjust some parameters, or check the function evaluation routines\n");
181: PetscPrintf(MPI_COMM_WORLD," Iterations: %d, Function Value: %4.2e, Residual: %4.2e \n",iter,ff,gnorm);
182: }
184: info = PetscOptionsHasName(PETSC_NULL,"-view_sol",&flg); CHKERRQ(info);
185: if (flg){
186: info = DAAppGetSolution(CombApp,nlevels-1,&X); CHKERRQ(info);
187: info=VecView(X,PETSC_VIEWER_STDOUT_WORLD); CHKERRQ(info);
188: }
190: /* To View TAO solver information */
191: // info = TaoView(tao); CHKERRQ(info);
193: /* Free TAO data structures */
194: info = TaoDestroy(tao); CHKERRQ(info);
195: info = TaoAppDestroy(CombApp); CHKERRQ(info);
197: /* Free PETSc data structures */
198: for (iter=0;iter<nlevels;iter++){
199: info = DADestroy(DAarray[iter]); CHKERRQ(info);
200: }
202: PreLoadEnd();
204: /* Finalize TAO */
205: TaoFinalize();
206: PetscFinalize();
208: return 0;
209: } /* main */
213: /*----- The following two routines
214:
215: MyGridMonitorBefore MyGridMonitorAfter
217: help diplay info of iterations at every grid level -------*/
221: static int MyGridMonitorBefore(TAO_APPLICATION myapp, DA da, int level, void *ctx) {
223: AppCtx *user = (AppCtx*)ctx;
224: int info,mx,my;
226: info = DAGetInfo(da,PETSC_NULL,&mx,&my,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL,
227: PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL);CHKERRQ(info);
228: user->mx = mx;
229: user->my = my;
230: user->hx = 1.0 / (user->mx - 1);
231: user->hy = 1.0 / (user->my - 1);
232: user->area = 0.5 * user->hx * user->hy;
233: user->fgctx.hx = user->hx;
234: user->fgctx.hy = user->hy;
235: user->fgctx.area = user->area;
236: user->fgctx.param = user->param;
238: PetscPrintf(MPI_COMM_WORLD,"Grid: %d, mx: %d my: %d \n",level,mx,my);
239: return 0;
240: }
242: /*------- USER-DEFINED: initialize the application context information -------*/
246: /*
247: AppCtxInitialize - Sets initial values for the application context parameters
249: Input:
250: ptr - void user-defined application context
252: Output:
253: ptr - user-defined application context with the default or user-provided
254: parameters
255: */
256: static int AppCtxInitialize(void *ptr) {
258: AppCtx *user = (AppCtx*)ptr;
259: PetscReal LambdaMax = 6.81, LambdaMin = 0.0; /* bounds on parameter lambda */
260: PetscTruth flg; /* flag for PETSc calls */
261: int info;
263: /* Specify dimension of the problem */
264: user->param = 5.0;
265: user->mx = 11;
266: user->my = 11;
268: /* Check for any command line arguments that override defaults */
269: info = PetscOptionsGetReal(TAO_NULL, "-par", &user->param, &flg); CHKERRQ(info);
270: if (user->param >= LambdaMax || user->param <= LambdaMin) {
271: SETERRQ(1,"Lambda is out of range.");
272: }
273: info = PetscOptionsGetInt(PETSC_NULL,"-mx",&user->mx,&flg); CHKERRQ(info);
274: info = PetscOptionsGetInt(PETSC_NULL,"-my",&user->my,&flg); CHKERRQ(info);
276: user->hx = 1.0 / (user->mx - 1);
277: user->hy = 1.0 / (user->my - 1);
278: user->area = 0.5 * user->hx * user->hy;
279: info = PetscLogFlops(6); CHKERRQ(info);
281: return 0;
282: } /* AppCtxInitialize */
288: static int FormInitialGuess(DA da, Vec X, AppCtx *ctx)
289: {
290: int info, i, j, mx, my;
291: int xs, ys, xm, ym, xe, ye;
292: PetscReal hx, hy, temp, val, lambda;
293: double **x;
295: lambda = ctx->param;
296: lambda = lambda/(lambda+1.0);
298: /* Get local mesh boundaries */
299: info = DAGetInfo(da,PETSC_NULL,&mx,&my,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL,
300: PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL);CHKERRQ(info);
301: hx = 1.0/(mx-1); hy = 1.0/(my-1);
303: info = DAGetCorners(da,&xs,&ys,PETSC_NULL,&xm,&ym,PETSC_NULL); CHKERRQ(info);
304: xe = xs+xm; ye = ys+ym;
306: info = DAVecGetArray(da, X, (void**)&x); CHKERRQ(info);
307: /* Compute initial guess over locally owned part of mesh */
308: for (j=ys; j<ye; j++) { /* for (j=0; j<my; j++) */
309: temp = PetscMin(j+1,my-j)*hy;
310: for (i=xs; i<xe; i++) { /* for (i=0; i<mx; i++) */
311: val = lambda*sqrt(PetscMin((PetscMin(i+1,mx-i))*hx,temp));
312: x[j][i] = val;
313: }
314: }
315: info = DAVecRestoreArray(da, X, (void**)&x); CHKERRQ(info);
317: return 0;
318: }
321: /*------- USER-DEFINED: set the upper and lower bounds for the variables -------*/
325: /*
326: FormBounds - Forms bounds on the variables
328: Input:
329: user - user-defined application context
331: Output:
332: XL - vector of lower bounds
333: XU - vector of upper bounds
335: */
336: static int DAFixBoundary(TAO_APPLICATION daapplication, DA da, Vec XL, Vec XU, void *ptr)
337: {
338: AppCtx *user = (AppCtx*)ptr;
339: int i, j, mx, my, info, xs, xm, ys, ym;
340: double lb = -TAO_INFINITY;
341: double ub = TAO_INFINITY;
342: double **xl, **xu;
344: mx = user->mx; /* number of points including the boundary */
345: my = user->my;
347: info = DAVecGetArray(da, XL, (void**)&xl); CHKERRQ(info);
348: info = DAVecGetArray(da, XU, (void**)&xu); CHKERRQ(info);
349: info = DAGetCorners(da, &xs, &ys, TAO_NULL, &xm, &ym, TAO_NULL); CHKERRQ(info);
351: /* Compute initial guess over locally owned part of the grid */
352: for (j = ys; j < ys+ym; j++){
353: for (i = xs; i < xs+xm; i++){
354: if (i == 0 || j == 0 || i == mx - 1 || j == my - 1) {
355: xl[j][i] = xu[j][i] = 0.0;
356: } else {
357: xl[j][i] = lb;
358: xu[j][i] = ub;
359: }
360: }
361: }
363: info = DAVecRestoreArray(da, XL, (void**)&xl); CHKERRQ(info);
364: info = DAVecRestoreArray(da, XU, (void**)&xu); CHKERRQ(info);
366: return 0;
367: } /* DAFixBoundary */
370: /*------- USER-DEFINED: routine to evaluate the function and gradient
371: at a local (rectangular element) level -------*/
375: /*
376: CombLocalFunctionGradient - Evaluates function and gradient over the
377: local rectangular element
379: Input:
380: coor - vector with the indices of the position of current element
381: in the first, second and third directions
382: x - current point (values over the current rectangular element)
383: ptr - user-defined application context
385: Output:
386: f - value of the objective funtion at the local rectangular element
387: g - gradient of the local function
389: */
390: static int CombLocalFunctionGradient(int coor[2], double x[4], double *f, double g[4], void *ptr) {
392: AppCtx *user = (AppCtx*)ptr;
394: double lambdad3, hx, hy, area;
395: double fquad, fexp, dvdx, dvdy;
397: lambdad3 = user->param / 3.0;
398: hx = user->hx;
399: hy = user->hy;
400: area = user->area;
402: /* lower triangle contribution */
403: dvdx = (x[0] - x[1]) / hx;
404: dvdy = (x[0] - x[2]) / hy;
405: fquad = dvdx * dvdx + dvdy * dvdy;
406: fexp = exp(x[0]) + exp(x[1]) + exp(x[2]);
408: dvdx = 0.5 * dvdx * hy;
409: dvdy = 0.5 * dvdy * hx;
410: g[0] = dvdx + dvdy - exp(x[0]) * area * lambdad3;
411: g[1] = -dvdx - 2.0 * exp(x[1]) * area * lambdad3;
412: g[2] = -dvdy - 2.0 * exp(x[2]) * area * lambdad3;
414: /* upper triangle contribution */
415: dvdx = (x[3] - x[2]) / hx;
416: dvdy = (x[3] - x[1]) / hy;
417: fquad += dvdx * dvdx + dvdy * dvdy;
418: fexp += exp(x[1]) + exp(x[2]) + exp(x[3]);
420: dvdx = 0.5 * dvdx * hy;
421: dvdy = 0.5 * dvdy * hx;
422: g[1] += -dvdy;
423: g[2] += -dvdx;
424: g[3] = dvdx + dvdy - exp(x[3]) * area * lambdad3;
427: *f = area * (0.5 * fquad - lambdad3 * fexp);
429: return 0;
430: } /* CombLocalFunctionGradient */
433: /*------- USER-DEFINED: routine to evaluate the Hessian
434: at a local (rectangular element) level -------*/
438: /*
439: CombLocalHessian - Computes the Hessian of the local (partial) function
440: defined over the current rectangle
442: Input:
443: coor - vector with the indices of the position of current element
444: in the first, second and third directions
445: x - current local solution (over the rectangle only)
446: ptr - user-defined application context
448: Output:
449: H - Hessian matrix of the local function (wrt the four
450: points of the rectangle only)
452: */
453: static int CombLocalHessian(int coor[2], double x[4], double H[4][4],void *ptr) {
455: AppCtx *user = (AppCtx*)ptr;
456: double hx, hy, lambdad3, area, dxdy, dydx;
457: double diagxy, dexp, bandxy, bandyx;
459: hx = user->hx;
460: hy = user->hy;
461: lambdad3 = user->param / 3.0;
462: area = user->area;
463: dxdy = hx/hy;
464: dydx = hy/hx;
465: diagxy = 0.5 * (dxdy + dydx);
466: bandxy = -0.5 * dxdy;
467: bandyx = -0.5 * dydx;
469: /* Hessian contribution at 0,0 */
470: dexp = exp(x[0]) * area * lambdad3;
471: H[0][0] = diagxy - dexp;
472: H[0][1] = H[1][0] = bandyx;
473: H[0][2] = H[2][0] = bandxy;
474: H[0][3] = H[3][0] = 0.0;
476: /* Hessian contribution at 1,0 */
477: dexp = exp(x[1]) * area * 2.0 * lambdad3;
478: H[1][1] = diagxy - dexp;
479: H[1][2] = H[2][1] = 0.0;
480: H[1][3] = H[3][1] = bandxy;
482: /* Hessian contribution at 0,1 */
483: dexp = exp(x[2]) * area * 2.0 * lambdad3;
484: H[2][2] = diagxy - dexp;
485: H[2][3] = H[3][2] = bandyx;
487: /* Hessian contribution at 1,1 */
488: dexp = exp(x[3]) * area * lambdad3;
489: H[3][3] = diagxy - dexp;
491: return 0;
492: } /* CombLocalHessian */
495: /*------- USER-DEFINED: routine to evaluate the function
496: and gradient at the whole grid -------*/
500: /*
501: WholeCombFunctionGradient - Evaluates function and gradient over the
502: whole grid
504: Input:
505: daapplication - TAO application object
506: da - distributed array
507: X - the current point, at which the function and gradient are evaluated
508: ptr - user-defined application context
510: Output:
511: f - value of the objective funtion at X
512: G - gradient at X
513: */
514: static int WholeCombFunctionGradient(TAO_APPLICATION daapplication, DA da, Vec X, double *f, Vec G, void *ptr) {
516: AppCtx *user = (AppCtx*)ptr;
517: Vec localX, localG;
518: int info, i, j;
519: int xs, xm, gxs, gxm, xe, ys, ym, gys, gym, ye;
520: double **x, **g;
521: double floc = 0.0;
522: PetscScalar zero = 0.0;
524: double lambdad3, hx, hy, area;
525: double fquad, fexp, dvdx, dvdy;
527: lambdad3 = user->param / 3.0;
528: hx = user->hx;
529: hy = user->hy;
530: area = user->area;
532: info = DAGetLocalVector(da, &localX); CHKERRQ(info);
533: info = DAGetLocalVector(da, &localG); CHKERRQ(info);
534: info = VecSet(G, zero); CHKERRQ(info);
535: info = VecSet(localG, zero); CHKERRQ(info);
537: info = DAGlobalToLocalBegin(da, X, INSERT_VALUES, localX); CHKERRQ(info);
538: info = DAGlobalToLocalEnd(da, X, INSERT_VALUES, localX); CHKERRQ(info);
540: info = DAVecGetArray(da, localX, (void**)&x); CHKERRQ(info);
541: info = DAVecGetArray(da, localG, (void**)&g); CHKERRQ(info);
543: info = DAGetCorners(da, &xs, &ys, TAO_NULL, &xm, &ym, TAO_NULL); CHKERRQ(info);
544: info = DAGetGhostCorners(da, &gxs, &gys, TAO_NULL, &gxm, &gym, TAO_NULL); CHKERRQ(info);
546: xe = gxs + gxm - 1;
547: ye = gys + gym - 1;
548: for (j = ys; j < ye; j++) {
549: for (i = xs; i < xe; i++) {
551: /* lower triangle contribution */
552: dvdx = (x[j][i] - x[j][i+1]) / hx;
553: dvdy = (x[j][i] - x[j+1][i]) / hy;
554: fquad = dvdx * dvdx + dvdy * dvdy;
555: fexp = exp(x[j][i]) + exp(x[j][i+1]) + exp(x[j+1][i]);
557: dvdx = 0.5 * dvdx * hy;
558: dvdy = 0.5 * dvdy * hx;
559: g[j][i] += dvdx + dvdy - exp(x[j][i]) * area * lambdad3;
560: g[j][i+1] += -dvdx - 2.0 * exp(x[j][i+1]) * area * lambdad3;
561: g[j+1][i] += -dvdy - 2.0 * exp(x[j+1][i]) * area * lambdad3;
563: /* upper triangle contribution */
564: dvdx = (x[j+1][i+1] - x[j+1][i]) / hx;
565: dvdy = (x[j+1][i+1] - x[j][i+1]) / hy;
566: fquad += dvdx * dvdx + dvdy * dvdy;
567: fexp += exp(x[j][i+1]) + exp(x[j+1][i]) + exp(x[j+1][i+1]);
569: dvdx = 0.5 * dvdx * hy;
570: dvdy = 0.5 * dvdy * hx;
571: g[j][i+1] += -dvdy;
572: g[j+1][i] += -dvdx;
573: g[j+1][i+1] += dvdx + dvdy - exp(x[j+1][i+1]) * area * lambdad3;
575: floc += area * (0.5 * fquad - fexp * lambdad3);
577: }
578: }
580: info = MPI_Allreduce(&floc, f, 1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD); CHKERRQ(info);
582: info = DAVecRestoreArray(da, localX, (void**)&x); CHKERRQ(info);
583: info = DAVecRestoreArray(da, localG, (void**)&g); CHKERRQ(info);
585: info = DALocalToGlobalBegin(da, localG, G); CHKERRQ(info);
586: info = DALocalToGlobalEnd(da, localG, G); CHKERRQ(info);
588: info = DARestoreLocalVector(da, &localX); CHKERRQ(info);
589: info = DARestoreLocalVector(da, &localG); CHKERRQ(info);
591: info = PetscLogFlops((xe-xs) * (ye-ys) * 55 + 1); CHKERRQ(info);
592: return 0;
594: } /* WholeCombFunctionGradient */
597: /*------- USER-DEFINED: routine to evaluate the Hessian
598: at the whole grid -------*/
601: /*
602: WholeCombHessian - Evaluates Hessian over the whole grid
604: Input:
605: daapplication - TAO application object
606: da - distributed array
607: X - the current point, at which the function and gradient are evaluated
608: ptr - user-defined application context
610: Output:
611: H - Hessian at X
612: */
613: static int WholeCombHessian(TAO_APPLICATION daapplication, DA da, Vec X, Mat H, void *ptr) {
615: AppCtx *user = (AppCtx*)ptr;
616: Vec localX;
617: int info, i, j, ind[4];
618: int xs, xm, gxs, gxm, xe, ys, ym, gys, gym, ye;
619: double smallH[4][4];
620: double **x;
622: double hx, hy, lambdad3, area, dxdy, dydx;
623: double diagxy, dexp, bandxy, bandyx;
624: PetscTruth assembled;
627: hx = user->hx;
628: hy = user->hy;
629: lambdad3 = user->param / 3.0;
630: area = user->area;
631: dxdy = hx/hy;
632: dydx = hy/hx;
633: diagxy = 0.5 * (dxdy + dydx);
634: bandxy = -0.5 * dxdy;
635: bandyx = -0.5 * dydx;
637: info = DAGetLocalVector(da, &localX); CHKERRQ(info);
638: info = MatAssembled(H,&assembled); CHKERRQ(info);
639: if (assembled){info = MatZeroEntries(H); CHKERRQ(info);}
641: info = DAGlobalToLocalBegin(da, X, INSERT_VALUES, localX); CHKERRQ(info);
642: info = DAGlobalToLocalEnd(da, X, INSERT_VALUES, localX); CHKERRQ(info);
644: info = DAVecGetArray(da, localX, (void**)&x); CHKERRQ(info);
646: info = DAGetCorners(da, &xs, &ys, TAO_NULL, &xm, &ym, TAO_NULL); CHKERRQ(info);
647: info = DAGetGhostCorners(da, &gxs, &gys, TAO_NULL, &gxm, &gym, TAO_NULL); CHKERRQ(info);
649: xe = gxs + gxm - 1;
650: ye = gys + gym - 1;
651: for (j = ys; j < ye; j++) {
652: for (i = xs; i < xe; i++) {
654: /* Hessian contribution at 0,0 */
655: dexp = exp(x[j][i]) * area * lambdad3;
656: smallH[0][0] = diagxy - dexp;
657: smallH[0][1] = smallH[1][0] = bandyx;
658: smallH[0][2] = smallH[2][0] = bandxy;
659: smallH[0][3] = smallH[3][0] = 0.0;
661: /* Hessian contribution at 1,0 */
662: dexp = exp(x[j][i+1]) * area * 2.0 * lambdad3;
663: smallH[1][1] = diagxy - dexp;
664: smallH[1][2] = smallH[2][1] = 0.0;
665: smallH[1][3] = smallH[3][1] = bandxy;
667: /* Hessian contribution at 0,1 */
668: dexp = exp(x[j+1][i]) * area * 2.0 * lambdad3;
669: smallH[2][2] = diagxy - dexp;
670: smallH[2][3] = smallH[3][2] = bandyx;
672: /* Hessian contribution at 1,1 */
673: dexp = exp(x[j+1][i+1]) * area * lambdad3;
674: smallH[3][3] = diagxy - dexp;
676: ind[0] = (j-gys) * gxm + (i-gxs);
677: ind[1] = ind[0] + 1;
678: ind[2] = ind[0] + gxm;
679: ind[3] = ind[2] + 1;
680: info = MatSetValuesLocal(H,4,ind,4,ind,(PetscScalar*)smallH,ADD_VALUES); CHKERRQ(info);
682: }
683: }
685: info = DAVecRestoreArray(da, localX, (void**)&x); CHKERRQ(info);
687: info = MatAssemblyBegin(H, MAT_FINAL_ASSEMBLY); CHKERRQ(info);
688: info = MatAssemblyEnd(H, MAT_FINAL_ASSEMBLY); CHKERRQ(info);
689: info = MatSetOption(H, MAT_SYMMETRIC); CHKERRQ(info);
691: info = DARestoreLocalVector(da, &localX); CHKERRQ(info);
693: info = PetscLogFlops((xe-xs) * (ye-ys) * 14 + 7); CHKERRQ(info);
694: return 0;
696: } /* WholeCombHessian */