Actual source code: jbearing3.c

  1: /*$Id: jbearing.c, v 1.1 2002/08/08 11:35 lopezca@mauddib.mcs.anl.gov $*/

  3: /* Program usage: mpirun -np <proc> jbearing [-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[] ="Pressure distribution in a Journal Bearing. \n\
 16: This example is based on the problem DPJB from the MINPACK-2 test suite.\n\
 17: This pressure journal bearing problem is an example of elliptic variational\n\
 18: problem defined over a two dimensional rectangle. By discretizing the domain \n\
 19: into triangular elements, the pressure surrounding the journal bearing is\n\
 20: defined as the minimum of a quadratic function whose variables are bounded\n\
 21: below by zero. The command line options are:\n\
 22:   -ecc <ecc>, where <ecc> = epsilon parameter\n\
 23:   -b <b>, where <b> = half the upper limit in the 2nd coordinate direction\n\
 24:   -mx <xg>, where <xg> = number of grid points in the 1st coordinate direction\n\
 25:   -my <yg>, where <yg> = number of grid points in the 2nd coordinate direction\n\
 26:   -nlevels <nlevels>, where <nlevels> = number of levels in multigrid\n\
 27:   -byelement, if computation is made by functions on rectangular elements\n\
 28:   -adic, if AD is used (AD is not used by default)\n\n";

 30: /*T
 31:    Concepts: TAO - Solving a bounded minimization problem
 32:    Routines: TaoInitialize(); TaoFinalize();
 33:    Routines: TaoCreate(); TaoDestroy();
 34:    Routines: DAApplicationCreate(); DAApplicationDestroy();
 35:    Routines: DAAppSetVariableBoundsRoutine;
 36:    Routines: DAAppSetElementObjectiveAndGradientRoutine();
 37:    Routines: DAAppSetElementHessianRoutine();
 38:    Routines: DAAppSetObjectiveAndGradientRoutine();
 39:    Routines: DAAppSetADElementFunctionGradient();
 40:    Routines: DAAppSetHessianRoutine();
 41:    Routines: TaoSetOptions();
 42:    Routines: TaoGetSolutionStatus(); TaoDAAppSolve();
 43:    Routines: DAAppSetBeforeMonitor(); DAAppSetAfterMonitor
 44:    Routines: DAAppGetSolution(); TaoView();
 45:    Routines: DAAppGetInterpolationMatrix();
 46:    Processors: n
 47: T*/
 48:  
 49: /*
 50:    User-defined application context - contains data needed by the
 51:    application-provided call-back routines.
 52: */

 54: int  ad_JBearLocalFunction(int[2] ,DERIV_TYPE[4], DERIV_TYPE *, void*);
 55: typedef struct {

 57:   InactiveDouble      *wq, *wl;      /* vectors with the parameters w_q(x) and w_l(x) */
 58:   InactiveDouble      hx, hy;        /* increment size in both directions */
 59:   InactiveDouble      area;          /* area of the triangles */

 61: } ADFGCtx;


 64: typedef struct {

 66:   PetscReal      ecc;           /* epsilon value */
 67:   PetscReal      b;             /* 0.5 * upper limit for 2nd variable */
 68:   double      *wq, *wl;      /* vectors with the parameters w_q(x) and w_l(x) */
 69:   double      hx, hy;        /* increment size in both directions */
 70:   double      area;          /* area of the triangles */

 72:   int         mx, my;        /* discretization including boundaries */

 74:   ADFGCtx     fgctx;         /* Used only when an ADIC generated gradient is used */

 76: } AppCtx;

 78: /* User-defined routines found in this file */
 79: static int AppCtxInitialize(void *ptr);
 80: static int FormInitialGuess(DA, Vec);

 82: static int JBearLocalFunctionGradient(int[2], double x[4], double *f, double g[4], void *ptr);
 83: static int JBearLocalHessian(int[2], double x[4], double H[4][4], void *ptr);

 85: static int WholeJBearFunctionGradient(TAO_APPLICATION,DA,Vec,double *,Vec,void*);
 86: static int WholeJBearHessian(TAO_APPLICATION,DA,Vec,Mat,void*);

 88: static int DASetBounds(TAO_APPLICATION, DA, Vec, Vec, void*);

 90: static int MyGridMonitorBefore(TAO_APPLICATION, DA, int, void *);
 91: static int MyGridMonitorAfter1(TAO_APPLICATION, DA, int, void *);

 95: int main( int argc, char **argv ) {

 97:   int             info,iter;             /* used to check for functions returning nonzeros */
 98:   int             nlevels;                                             /* multigrid levels */
 99:   int             Nx,Ny;
100:   double          ff,gnorm;
101:   DA              DAarray[20];
102:   Vec             X;
103:   PetscTruth      flg, PreLoad = PETSC_TRUE;                                      /* flags */
104:   AppCtx          user;                                       /* user-defined work context */
105:   TaoMethod       method = "tao_gpcg";                              /* minimization method */
106:   TAO_SOLVER      tao;                                        /* TAO_SOLVER solver context */
107:   TAO_APPLICATION JBearApp;                                    /* The PETSc application */
108:   TaoTerminateReason reason;
109:   KSP ksp;
110:   PC  pc;

112:   /* Initialize TAO */
113:   PetscInitialize(&argc, &argv, (char *)0, help);
114:   TaoInitialize(&argc, &argv, (char *)0, help);

116:   PreLoadBegin(PreLoad,"Solve");
117:   
118:   info = AppCtxInitialize((void*)&user); CHKERRQ(info);
119:   
120:   nlevels=5;
121:   info = PetscOptionsGetInt(PETSC_NULL,"-nlevels",&nlevels,&flg); CHKERRQ(info);
122:   if (PreLoadIt == 0) {
123:     nlevels = 1; user.mx = 11; user.my = 11; }

125:   PetscPrintf(MPI_COMM_WORLD,"\n---- Journal Bearing Problem -----\n\n");

127:   /* Let PETSc determine the vector distribution */
128:   Nx = PETSC_DECIDE; Ny = PETSC_DECIDE;

130:   /* Create distributed array (DA) to manage parallel grid and vectors  */
131:   info = DACreate2d(PETSC_COMM_WORLD,DA_NONPERIODIC,DA_STENCIL_BOX,user.mx,
132:                     user.my,Nx,Ny,1,1,PETSC_NULL,PETSC_NULL,&DAarray[0]); CHKERRQ(info);
133:   for (iter=1;iter<nlevels;iter++){
134:     info = DARefine(DAarray[iter-1],PETSC_COMM_WORLD,&DAarray[iter]); CHKERRQ(info);
135:   }

137:   /* Create TAO solver and set desired solution method */
138:   info = TaoCreate(MPI_COMM_WORLD,method,&tao); CHKERRQ(info);
139:   info = TaoApplicationCreate(PETSC_COMM_WORLD, &JBearApp); CHKERRQ(info);
140:   info = TaoAppSetDAApp(JBearApp,DAarray,nlevels); CHKERRQ(info);

142:   /* Sets routine bounds evaluation */
143:   info = DAAppSetVariableBoundsRoutine(JBearApp,DASetBounds,(void *)&user); CHKERRQ(info);

145:   info = PetscOptionsHasName(TAO_NULL, "-byelement", &flg); CHKERRQ(info);
146:   if (flg) {

148:     /* Sets routines for function and gradient evaluation, element by element */
149:     info = PetscOptionsHasName(TAO_NULL, "-adic", &flg); CHKERRQ(info);
150:     if (flg) {
151:       info = DAAppSetADElementFunctionGradient(JBearApp,ad_JBearLocalFunction,248,(void *)&user.fgctx); CHKERRQ(info);
152:     } else {
153:       info = DAAppSetElementObjectiveAndGradientRoutine(JBearApp,JBearLocalFunctionGradient,63,(void *)&user); CHKERRQ(info);
154:     }
155:     /* Sets routines for Hessian evaluation, element by element */
156:     info = DAAppSetElementHessianRoutine(JBearApp,JBearLocalHessian,16,(void*)&user); CHKERRQ(info);

158:   } else {

160:     /* Sets routines for function and gradient evaluation, all in one routine */
161:     info = DAAppSetObjectiveAndGradientRoutine(JBearApp,WholeJBearFunctionGradient,(void *)&user); CHKERRQ(info);

163:     /* Sets routines for Hessian evaluation, all in one routine */
164:     info = DAAppSetHessianRoutine(JBearApp,WholeJBearHessian,(void*)&user); CHKERRQ(info);    
165:     
166:   }

168:   info = DAAppSetBeforeMonitor(JBearApp,MyGridMonitorBefore,(void*)&user); CHKERRQ(info);
169:   info = DAAppSetAfterMonitor(JBearApp,MyGridMonitorAfter1,(void*)&user); CHKERRQ(info);
170:   info = PetscOptionsHasName(TAO_NULL,"-tao_monitor", &flg); CHKERRQ(info);
171:   if (flg){
172:     info = DAAppPrintInterpolationError(JBearApp); CHKERRQ(info);
173:     info = DAAppPrintStageTimes(JBearApp); CHKERRQ(info);
174:   }
175:   /* Check for any tao command line options */
176:   info = TaoAppSetRelativeTolerance(JBearApp,1.0e-6); CHKERRQ(info);
177:   info = TaoSetTolerances(tao,0,0,0,0); CHKERRQ(info);
178:   info = TaoSetGradientTolerances(tao,0,0,0); CHKERRQ(info);

180:   info = TaoAppGetKSP(JBearApp,&ksp);CHKERRQ(info);
181:   info = KSPSetType(ksp,KSPCG); CHKERRQ(info);
182:   info = KSPSetTolerances(ksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,10);CHKERRQ(info);
183:   info = KSPGetPC(ksp,&pc);CHKERRQ(info);
184:   info = PCSetType(pc,PCBJACOBI);CHKERRQ(info);

186:   info = TaoSetOptions(JBearApp,tao); CHKERRQ(info);

188:   info = DAAppGetSolution(JBearApp,0,&X); CHKERRQ(info);
189:   info = FormInitialGuess(DAarray[0],X); CHKERRQ(info);
190:   info = DAAppSetInitialSolution(JBearApp,X); CHKERRQ(info);
191:   /* SOLVE THE APPLICATION */
192:   info = TaoDAAppSolve(JBearApp,tao);  CHKERRQ(info);

194:   /* Get information on termination */
195:   info = TaoGetSolutionStatus(tao,&iter,&ff,&gnorm,0,0,&reason); CHKERRQ(info);
196:   if (reason <= 0 ){
197:     PetscPrintf(MPI_COMM_WORLD,"Try a different TAO method, adjust some parameters, or check the function evaluation routines\n");
198:     PetscPrintf(MPI_COMM_WORLD," Iterations: %d,  Function Value: %4.2e, Residual: %4.2e \n",iter,ff,gnorm);
199:   }

201:   info = PetscOptionsHasName(PETSC_NULL,"-view_sol",&flg); CHKERRQ(info);
202:   if (flg){
203:     info = DAAppGetSolution(JBearApp,nlevels-1,&X); CHKERRQ(info);
204:     info=VecView(X,PETSC_VIEWER_STDOUT_WORLD); CHKERRQ(info);
205:   }

207:   /*  To View TAO solver information */
208:   // info = TaoView(tao); CHKERRQ(info);

210:   /* Free TAO data structures */
211:   info = TaoDestroy(tao); CHKERRQ(info);
212:   info = TaoAppDestroy(JBearApp); CHKERRQ(info);

214:   /* Free PETSc data structures */
215:   for (iter=0;iter<nlevels;iter++){
216:     info = DADestroy(DAarray[iter]); CHKERRQ(info);
217:   }

219:   PreLoadEnd();

221:   /* Finalize TAO */
222:   TaoFinalize();
223:   PetscFinalize();

225:   return 0;
226: } /* main */



230: /*----- The following two routines
231:    
232:   MyGridMonitorBefore    MyGridMonitorAfter

234:   help diplay info of iterations at every grid level -------*/

238: static int MyGridMonitorBefore(TAO_APPLICATION myapp, DA da, int level, void *ctx) {

240:   AppCtx *user = (AppCtx*)ctx;
241:   int info,mx,my;
242:   double t;

244:   info = DAGetInfo(da,PETSC_NULL,&mx,&my,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL,
245:                    PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL);CHKERRQ(info);
246:   user->mx = mx;
247:   user->my = my;

249:   user->hx = (2.0 * 3.14159265358979) / (user->mx - 1);
250:   user->hy = (2.0 * user->b) / (user->my - 1);
251:   user->area = 0.5 * user->hx * user->hy;

253:   user->wq = new double[user->mx];
254:   user->wl = new double[user->mx];
255:   for (int i=0; i<user->mx; i++) {
256:     t = 1.0 + user->ecc * cos(i*user->hx);
257:     user->wq[i] = t*t*t;
258:     user->wl[i] = user->ecc * sin(i*user->hx);
259:   }
260:   info = PetscLogFlops(8 + 7 * user->mx); CHKERRQ(info);

262:   user->fgctx.hx   = user->hx;
263:   user->fgctx.hy   = user->hy;
264:   user->fgctx.area = user->area;
265:   user->fgctx.wq = user->wq;
266:   user->fgctx.wl = user->wl;

268:   PetscPrintf(MPI_COMM_WORLD,"Grid: %d,    mx: %d     my: %d   \n",level,mx,my);

270:   return 0;
271: }

275: static int MyGridMonitorAfter1(TAO_APPLICATION myapp, DA da, int level, void *ctx){
276:   
277:   AppCtx *user = (AppCtx*)ctx;
278:   delete [] user->wq; delete [] user->wl;
279:   return 0;
280: }



284: /*------- USER-DEFINED: initialize the application context information -------*/

288: /*
289:   AppCtxInitialize - Sets initial values for the application context parameters

291:   Input:
292:     ptr - void user-defined application context

294:   Output:
295:     ptr - user-defined application context with the default or user-provided
296:              parameters
297: */
298: static int AppCtxInitialize(void *ptr) {

300:   AppCtx *user = (AppCtx*)ptr;
301:   PetscTruth    flg;            /* flag for PETSc calls */
302:   int info;

304:   /* Specify default parameters */
305:   user->ecc = 0.1;
306:   user->b = 10.0;
307:   user->mx = user->my = 11;

309:   /* Check for command line arguments that override defaults */
310:   info = PetscOptionsGetReal(TAO_NULL, "-ecc", &user->ecc, &flg); CHKERRQ(info);
311:   info = PetscOptionsGetReal(TAO_NULL, "-b", &user->b, &flg); CHKERRQ(info);
312:   info = PetscOptionsGetInt(TAO_NULL, "-mx", &user->mx, &flg); CHKERRQ(info);
313:   info = PetscOptionsGetInt(TAO_NULL, "-my", &user->my, &flg); CHKERRQ(info);

315:   return 0;
316: } /* AppCtxInitialize */


321: static int FormInitialGuess(DA da, Vec X)
322: {
323:   int    info, i, j, mx;
324:   int    xs, ys, xm, ym, xe, ye;
325:   PetscReal hx, val;
326:   double **x;

328:   /* Get local mesh boundaries */
329:   info = DAGetInfo(da,PETSC_NULL,&mx,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL,
330:                    PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL);CHKERRQ(info);
331:   hx = 2.0*4.0*atan(1.0)/(mx-1);

333:   info = DAGetCorners(da,&xs,&ys,PETSC_NULL,&xm,&ym,PETSC_NULL); CHKERRQ(info);
334:   xe = xs+xm; ye = ys+ym;

336:   info = DAVecGetArray(da, X, (void**)&x); CHKERRQ(info);
337:   /* Compute initial guess over locally owned part of mesh */
338:   for (j=ys; j<ye; j++) {  /*  for (j=0; j<my; j++) */
339:     for (i=xs; i<xe; i++) {  /*  for (i=0; i<mx; i++) */
340:       val = PetscMax(sin((i+1.0)*hx),0.0);
341:       x[j][i] = val;
342:       x[j][i] = 0;
343:     }
344:   }
345:   info = DAVecRestoreArray(da, X, (void**)&x); CHKERRQ(info);

347:   return 0;
348: }


351: /*------- USER-DEFINED: set the upper and lower bounds for the variables  -------*/
354: /*
355:   FormBounds - Forms bounds on the variables

357:   Input:
358:     user - user-defined application context

360:   Output:
361:     XL - vector of lower bounds
362:     XU - vector of upper bounds

364: */
365: static int DASetBounds(TAO_APPLICATION daapplication, DA da, Vec XL, Vec XU, void *ptr) {

367:   AppCtx *user = (AppCtx*)ptr;
368:   int i, j, info, mx, my;
369:   int xs, xm, ys, ym;
370:   double **xl, **xu;

372:   mx = user->mx;
373:   my = user->my;

375:   info = DAVecGetArray(da, XL, (void**)&xl); CHKERRQ(info);
376:   info = DAVecGetArray(da, XU, (void**)&xu); CHKERRQ(info);
377:   info = DAGetCorners(da, &xs, &ys, TAO_NULL, &xm, &ym, TAO_NULL); CHKERRQ(info);

379:   for (j = ys; j < ys+ym; j++){
380:     for (i = xs; i < xs+xm; i++){
381:       xl[j][i] = 0.0;
382:       if (i == 0 || j == 0 || i == mx - 1 || j == my - 1) {
383:         xu[j][i] = 0.0;
384:       } else {
385:         xu[j][i] = TAO_INFINITY;
386:       }
387:     }
388:   }

390:   info = DAVecRestoreArray(da, XL, (void**)&xl); CHKERRQ(info);
391:   info = DAVecRestoreArray(da, XU, (void**)&xu); CHKERRQ(info);

393:   return 0;

395: } /* DASetBounds */



401: /*
402:   JBearLocalFunctionGradient - Evaluates function and gradient over the 
403:       local rectangular element

405:   Input:
406:     coor - vector with the indices of the position of current element
407:              in the first, second and third directions
408:     x - current point (values over the current rectangular element)
409:     ptr - user-defined application context

411:   Output:
412:     f - value of the objective funtion at the local rectangular element
413:     g - gradient of the local function

415: */
416: static int JBearLocalFunctionGradient(int coor[2], double x[4], double *f, double g[4], void *ptr) {

418:   AppCtx *user = (AppCtx*)ptr;

420:   double avgWq, sqGrad, avgWV, fl, fu;
421:   double hx, hy, area, aread3, *wq, *wl;
422:   double dvdx, dvdy;
423:   int i;

425:   hx = user->hx;
426:   hy = user->hy;
427:   area = user->area;
428:   aread3 = area / 3.0;
429:   wq = user->wq;
430:   wl = user->wl;
431:   i = coor[0];

433:   /* lower triangle contribution */
434:   dvdx = (x[0] - x[1]) / hx;
435:   dvdy = (x[0] - x[2]) / hy;
436:   sqGrad = dvdx * dvdx + dvdy * dvdy;
437:   avgWq = (2.0 * wq[i] + wq[i+1]) / 6.0;
438:   avgWV = (wl[i]*x[0] + wl[i+1]*x[1] + wl[i]*x[2]) / 3.0;
439:   fl = avgWq * sqGrad - avgWV;

441:   dvdx = dvdx * hy * avgWq;
442:   dvdy = dvdy * hx * avgWq;
443:   g[0] = ( dvdx + dvdy ) - wl[i] * aread3;
444:   g[1] = ( -dvdx ) - wl[i+1] * aread3;
445:   g[2] = ( -dvdy ) - wl[i] * aread3;

447:   /* upper triangle contribution */
448:   dvdx = (x[3] - x[2]) / hx; 
449:   dvdy = (x[3] - x[1]) / hy;
450:   sqGrad = dvdx * dvdx + dvdy * dvdy;
451:   avgWq = (2.0 * wq[i+1] + wq[i]) / 6.0;
452:   avgWV = (wl[i+1]*x[1] + wl[i]*x[2] + wl[i+1]*x[3]) / 3.0;
453:   fu = avgWq * sqGrad - avgWV;

455:   dvdx = dvdx * hy * avgWq;
456:   dvdy = dvdy * hx * avgWq;
457:   g[1] += (-dvdy) - wl[i+1] * aread3;
458:   g[2] +=  (-dvdx) - wl[i] * aread3;
459:   g[3] = ( dvdx + dvdy ) - wl[i+1] * aread3;

461:   *f = area * (fl + fu);

463:   return 0;
464: } /* JBearLocalFunctionGradient */


467: /*------- USER-DEFINED: routine to evaluate the Hessian
468:            at a local (rectangular element) level       -------*/
471: /*
472:   JBearLocalHessian - Computes the Hessian of the local (partial) function
473:          defined over the current rectangle

475:   Input:
476:     coor - vector with the indices of the position of current element
477:              in the first, second and third directions
478:     x - current local solution (over the rectangle only)
479:     ptr - user-defined application context

481:   Output:
482:     H - Hessian matrix of the local function (wrt the four
483:            points of the rectangle only)

485: */
486: static int JBearLocalHessian(int coor[2], double x[4], double H[4][4],void *ptr) {

488:   AppCtx *user = (AppCtx*)ptr;
489:   double wql, wqu;
490:   double hx, hy, dydx, dxdy;
491:   double *wq;
492:   double wqldydx,wqldxdy,wqudydx,wqudxdy;
493:   int i;

495:   hx = user->hx;
496:   hy = user->hy;
497:   wq = user->wq;
498:   i = coor[0];

500:   dxdy = hx / hy;
501:   dydx = hy / hx;
502:   wql = (2.0 * wq[i] + wq[i+1]) / 6.0;
503:   wqu = (wq[i] + 2.0 * wq[i+1]) / 6.0;
504:   wqldydx = wql * dydx;
505:   wqldxdy = wql * dxdy;
506:   wqudydx = wqu * dydx;
507:   wqudxdy = wqu * dxdy;

509:           /* Hessian contribution at 0,0 */
510:   H[0][0] = wqldxdy + wqldydx;
511:   H[0][1] = H[1][0] = -wqldydx;
512:   H[0][2] = H[2][0] = -wqldxdy;
513:   H[0][3] = H[3][0] = 0.0;

515:           /* Hessian contribution at 1,0 */
516:   H[1][1] = wqldydx + wqudxdy;
517:   H[1][2] = H[2][1] = 0.0;
518:   H[1][3] = H[3][1] = -wqudxdy; 

520:           /* Hessian contribution at 0,1 */
521:   H[2][2] = wqldxdy + wqudydx;
522:   H[2][3] = H[3][2] = -wqudydx; 

524:           /* Hessian contribution at 1,1 */
525:   H[3][3] = wqudydx + wqudxdy;

527:   return 0;

529: } /* JBearLocalHessian */


532: /*------- USER-DEFINED: routine to evaluate the function 
533:           and gradient at the whole grid             -------*/

537: /*
538:   WholeJBearFunctionGradient - Evaluates function and gradient over the 
539:       whole grid

541:   Input:
542:     daapplication - TAO application object
543:     da  - distributed array
544:     X   - the current point, at which the function and gradient are evaluated
545:     ptr - user-defined application context

547:   Output:
548:     f - value of the objective funtion at X
549:     G - gradient at X
550: */
551: static int WholeJBearFunctionGradient(TAO_APPLICATION daapplication, DA da, Vec X, double *f, Vec G, void *ptr) {

553:   AppCtx *user = (AppCtx*)ptr;
554:   Vec localX, localG;
555:   int info, i, j;
556:   int xs, xm, gxs, gxm, xe, ys, ym, gys, gym, ye;
557:   double **x, **g;
558:   double floc = 0.0;
559:   PetscScalar zero = 0.0;

561:   double avgWq, sqGrad, avgWV, fl, fu;
562:   double hx, hy, area, aread3, *wq, *wl;
563:   double dvdx, dvdy;

565:   hx = user->hx;
566:   hy = user->hy;
567:   area = user->area;
568:   aread3= area/3.0;
569:   wq = user->wq;
570:   wl = user->wl;

572:   info = DAGetLocalVector(da, &localX); CHKERRQ(info);
573:   info = DAGetLocalVector(da, &localG); CHKERRQ(info);
574:   info = VecSet(G, zero); CHKERRQ(info);
575:   info = VecSet(localG, zero); CHKERRQ(info);

577:   info = DAGlobalToLocalBegin(da, X, INSERT_VALUES, localX); CHKERRQ(info);
578:   info = DAGlobalToLocalEnd(da, X, INSERT_VALUES, localX); CHKERRQ(info);

580:   info = DAVecGetArray(da, localX, (void**)&x); CHKERRQ(info);
581:   info = DAVecGetArray(da, localG, (void**)&g); CHKERRQ(info);

583:   info = DAGetCorners(da, &xs, &ys, TAO_NULL, &xm, &ym, TAO_NULL); CHKERRQ(info);
584:   info = DAGetGhostCorners(da, &gxs, &gys, TAO_NULL, &gxm, &gym, TAO_NULL); CHKERRQ(info);

586:   xe = gxs + gxm - 1;
587:   ye = gys + gym - 1;
588:   for (j = ys; j < ye; j++) {
589:     for (i = xs; i < xe; i++) {

591:       /* lower triangle contribution */
592:       dvdx = (x[j][i] - x[j][i+1]) / hx;
593:       dvdy = (x[j][i] - x[j+1][i]) / hy;
594:       sqGrad = dvdx * dvdx + dvdy * dvdy;
595:       avgWq = (2.0 * wq[i] + wq[i+1]) / 6.0;
596:       avgWV = (wl[i]*x[j][i] + wl[i+1]*x[j][i+1] + wl[i]*x[j+1][i]) / 3.0;
597:       fl = avgWq * sqGrad - avgWV;

599:       dvdx = dvdx * hy * avgWq;
600:       dvdy = dvdy * hx * avgWq;
601:       g[j][i] +=  ( dvdx + dvdy ) - wl[i] * aread3;
602:       g[j][i+1] += ( -dvdx ) - wl[i+1] * aread3;
603:       g[j+1][i] += ( -dvdy ) - wl[i] * aread3;

605:       /* upper triangle contribution */
606:       dvdx = (x[j+1][i+1] - x[j+1][i]) / hx;
607:       dvdy = (x[j+1][i+1] - x[j][i+1]) / hy;
608:       sqGrad = dvdx * dvdx + dvdy * dvdy;
609:       avgWq = (2.0 * wq[i+1] + wq[i]) / 6.0;
610:       avgWV = (wl[i+1]*x[j][i+1] + wl[i]*x[j+1][i] + wl[i+1]*x[j+1][i+1]) / 3.0;
611:       fu = avgWq * sqGrad - avgWV;

613:       dvdx = dvdx * hy * avgWq;
614:       dvdy = dvdy * hx * avgWq;
615:       g[j][i+1] += (-dvdy) - wl[i+1] * aread3;
616:       g[j+1][i] +=  (-dvdx) - wl[i] * aread3;
617:       g[j+1][i+1] += ( dvdx + dvdy ) - wl[i+1] * aread3;

619:       floc += area * (fl + fu);

621:     }
622:   }

624:   info = MPI_Allreduce(&floc, f, 1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD); CHKERRQ(info);

626:   info = DAVecRestoreArray(da, localX, (void**)&x); CHKERRQ(info);
627:   info = DAVecRestoreArray(da, localG, (void**)&g); CHKERRQ(info);

629:   info = DALocalToGlobalBegin(da, localG, G); CHKERRQ(info);
630:   info = DALocalToGlobalEnd(da, localG, G); CHKERRQ(info);

632:   info = DARestoreLocalVector(da, &localX); CHKERRQ(info);
633:   info = DARestoreLocalVector(da, &localG); CHKERRQ(info);

635:   info = PetscLogFlops((xe-xs) * (ye-ys) * 67 + 1); CHKERRQ(info);
636:   return 0;
637: } /* WholeJBearFunctionGradient */


642: /*
643:   WholeJBearHessian - Evaluates Hessian over the whole grid

645:   Input:
646:     daapplication - TAO application object
647:     da  - distributed array
648:     X   - the current point, at which the function and gradient are evaluated
649:     ptr - user-defined application context

651:   Output:
652:     H - Hessian at X
653: */
654: static int WholeJBearHessian(TAO_APPLICATION daapplication, DA da, Vec X, Mat H, void *ptr) {

656:   AppCtx *user = (AppCtx*)ptr;
657:   int info, i, j, ind[4];
658:   int xs, xm, gxs, gxm, xe, ys, ym, gys, gym, ye;
659:   double smallH[4][4];
660:   double wql, wqu;
661:   double dydx, dxdy;
662:   double *wq;
663:   double wqldydx,wqldxdy,wqudydx,wqudxdy;
664:   PetscTruth assembled;

666:   wq = user->wq;
667:   dydx = user->hy / user->hx;
668:   dxdy = user->hx / user->hy;

670:   info = MatAssembled(H,&assembled); CHKERRQ(info);
671:   if (assembled){info = MatZeroEntries(H);  CHKERRQ(info);}


674:   info = DAGetCorners(da, &xs, &ys, TAO_NULL, &xm, &ym, TAO_NULL); CHKERRQ(info);
675:   info = DAGetGhostCorners(da, &gxs, &gys, TAO_NULL, &gxm, &gym, TAO_NULL); CHKERRQ(info);

677:   xe = gxs + gxm - 1;
678:   ye = gys + gym - 1;
679:   for (j = ys; j < ye; j++) {
680:     for (i = xs; i < xe; i++) {

682:       wql = (2.0 * wq[i] + wq[i+1]) / 6.0;
683:       wqu = (wq[i] + 2.0 * wq[i+1]) / 6.0;

685:       wqldydx = wql * dydx;
686:       wqldxdy = wql * dxdy;
687:       wqudydx = wqu * dydx;
688:       wqudxdy = wqu * dxdy;

690:           /* Hessian contribution at 0,0 */
691:       smallH[0][0] = wqldxdy + wqldydx;
692:       smallH[0][1] = smallH[1][0] = -wqldydx;
693:       smallH[0][2] = smallH[2][0] = -wqldxdy;
694:       smallH[0][3] = smallH[3][0] = 0.0;

696:           /* Hessian contribution at 1,0 */
697:       smallH[1][1] = (wqldydx + wqudxdy);
698:       smallH[1][2] = smallH[2][1] = 0.0;
699:       smallH[1][3] = smallH[3][1] = -wqudxdy;

701:           /* Hessian contribution at 0,1 */
702:       smallH[2][2] = (wqldxdy + wqudydx);
703:       smallH[2][3] = smallH[3][2] = -wqudydx;

705:           /* Hessian contribution at 1,1 */
706:       smallH[3][3] = wqudydx + wqudxdy;

708:       ind[0] = (j-gys) * gxm + (i-gxs);
709:       ind[1] = ind[0] + 1;
710:       ind[2] = ind[0] + gxm;
711:       ind[3] = ind[2] + 1;
712:       info = MatSetValuesLocal(H,4,ind,4,ind,(PetscScalar*)smallH,ADD_VALUES); CHKERRQ(info);

714:     }
715:   }


718:   info = MatAssemblyBegin(H, MAT_FINAL_ASSEMBLY); CHKERRQ(info);
719:   info = MatAssemblyEnd(H, MAT_FINAL_ASSEMBLY); CHKERRQ(info);
720:   info = MatSetOption(H, MAT_SYMMETRIC); CHKERRQ(info);


723:   info = PetscLogFlops((xe-xs) * (ye-ys) * 14 + 2); CHKERRQ(info);
724:   return 0;

726: } /* WholeJBearHessian */