Actual source code: minsurf3.c
1: /*$Id: minsurf3.c,v 1.1 2002/05/13 19:32:57 benson Exp $*/
3: /* Program usage: mpirun -np <proc> minsurf3 [-help] [all TAO options] */
5: /* minsurf1: boundary values like in COPS example:
6: - parabola for top and bottom boundaries
7: - zero o.w. */
9: /*
10: Include "tao.h" so we can use TAO solvers.
11: petscda.h for distributed array
12: ad_deriv.h for AD gradient
13: */
15: #include "petscda.h"
16: #include "tao.h"
17: #include "taodaapplication.h"
19: static char help[] =
20: "Given a rectangular 2-D domain and boundary values along \n\
21: the edges of the domain, the objective is to find the surface with \n\
22: the minimal area that satisfies the boundary conditions.\n\
23: The command line options are:\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\
29: -bottom <b>, where <b> = bottom bound on the rectangular domain\n\
30: -top <t>, where <t> = bottom bound on the rectangular domain\n\
31: -left <l>, where <l> = left bound on the rectangular domain\n\
32: -right <r>, where <r> = right bound on the rectangular domain\n\
33: -cops <r>, if COPS boundaries should be use (default MINPACK bounds)\n\
34: -obst <obst>, where <obst> = 1 is obstacle is present, zero otherwise \n\n";
36: /*T
37: Concepts: TAO - Solving a bounded minimization problem
38: Routines: TaoInitialize(); TaoFinalize();
39: Routines: TaoCreate(); TaoDestroy();
40: Routines: DAAppSetVariableBoundsRoutine();
41: Routines: DAAppSetElementObjectiveAndGradientRoutine();
42: Routines: DAAppSetElementHessianRoutine();
43: Routines: DAAppSetObjectiveAndGradientRoutine();
44: Routines: DAAppSetADElementFunctionGradient();
45: Routines: DAAppSetHessianRoutine();
46: Routines: TaoSetOptions();
47: Routines: TaoAppGetSolutionStatus(); TaoDAAppSolve();
48: Routines: DAAppSetBeforeMonitor(); TaoView();
49: Routines: DAAppGetSolution();
50: Routines: DAAppGetInterpolationMatrix();
51: Processors: n
52: T*/
54: /*
55: User-defined application context - contains data needed by the
56: application-provided call-back routines.
57: */
60: /*
61: This structure is used only when an ADIC generated gradient is used.
62: An InactiveDouble type is a double
63: */
64: typedef struct {
65:
66: InactiveDouble hx, hy; /* increment size in both directions */
67: InactiveDouble area; /* area of the triangles */
69: } ADFGCtx;
70: int ad_MSurfLocalFunction(int[2], DERIV_TYPE[4], DERIV_TYPE*, void*);
72: typedef struct {
73: PetscReal b, t, l, r; /* domain boundaries */
74: PetscReal bheight; /* height of obstacle */
75: PetscReal fx, fy; /* relative size of obstacle */
76: double hx, hy; /* increment size in both directions */
77: double area; /* area of the triangles */
78: int mx, my; /* discretization including boundaries */
79: int bmx, bmy; /* discretization of obstacle */
80: ADFGCtx fgctx; /* Used only when an ADIC generated gradient is used */
81: } AppCtx;
83: /* User-defined routines */
84: static int AppCtxInitialize(void *);
86: static int MSurfLocalFunctionGradient(int[2], double[4], double *, double[4], void *);
87: static int WholeMSurfFunctionGradient(TAO_APPLICATION,DA,Vec,double *,Vec,void*);
89: static int MSurfLocalHessian(int[2], double[4], double[4][4], void *);
90: static int WholeMSurfHessian(TAO_APPLICATION,DA,Vec,Mat,void*);
92: static int COPS_Bounds(TAO_APPLICATION, DA, Vec, Vec, void*);
93: static int MINPACK_Bounds(TAO_APPLICATION, DA, Vec, Vec, void *);
95: static int MyGridMonitorBefore(TAO_APPLICATION, DA, int, void *);
100: int main( int argc, char **argv ) {
102: int info; /* used to check for functions returning nonzeros */
103: int iter,nlevels;
104: int Nx,Ny;
105: double ff,gnorm;
106: DA DAarray[20];
107: Vec X;
108: PetscTruth flg, PreLoad = PETSC_TRUE; /* flags */
109: TaoMethod method = "tao_bnls"; /* minimization method */
110: TaoTerminateReason reason;
111: TAO_SOLVER tao; /* TAO_SOLVER solver context */
112: TAO_APPLICATION MSurfApp; /* The PETSc application */
113: AppCtx user; /* user-defined work context */
114: KSP ksp;
115: PC pc;
117: /* Initialize TAO */
118: PetscInitialize(&argc, &argv, (char *)0, help);
119: TaoInitialize(&argc, &argv, (char *)0, help);
121: PreLoadBegin(PreLoad,"Solve");
123: info = AppCtxInitialize((void*)&user); CHKERRQ(info);
125: nlevels=4;
126: info = PetscOptionsGetInt(PETSC_NULL,"-nlevels",&nlevels,&flg); CHKERRQ(info);
127: if (PreLoadIt == 0) {
128: nlevels = 1; user.mx = 11; user.my = 11;
129: }
131: PetscPrintf(MPI_COMM_WORLD,"\n---- Minimal Surface Problem (simple boundary) -----\n\n");
133: /* Let PETSc determine the vector distribution */
134: Nx = PETSC_DECIDE; Ny = PETSC_DECIDE;
136: /* Create distributed array (DA) to manage parallel grid and vectors */
137: info = DACreate2d(PETSC_COMM_WORLD,DA_NONPERIODIC,DA_STENCIL_BOX,user.mx,
138: user.my,Nx,Ny,1,1,PETSC_NULL,PETSC_NULL,&DAarray[0]); CHKERRQ(info);
139: for (iter=1;iter<nlevels;iter++){
140: info = DARefine(DAarray[iter-1],PETSC_COMM_WORLD,&DAarray[iter]); CHKERRQ(info);
141: }
143: /* Create TAO solver and set desired solution method */
144: info = TaoCreate(MPI_COMM_WORLD,method,&tao); CHKERRQ(info);
145: info = TaoApplicationCreate(PETSC_COMM_WORLD, &MSurfApp); CHKERRQ(info);
146: info = TaoAppSetDAApp(MSurfApp,DAarray,nlevels); CHKERRQ(info);
148: /* Sets routines for function, gradient and bounds evaluation */
149: info = PetscOptionsHasName(TAO_NULL, "-cops", &flg); CHKERRQ(info);
150: if (flg){
151: info = DAAppSetVariableBoundsRoutine(MSurfApp,COPS_Bounds,(void *)&user); CHKERRQ(info);
152: } else {
153: info = DAAppSetVariableBoundsRoutine(MSurfApp,MINPACK_Bounds,(void *)&user); CHKERRQ(info);
154: }
156: info = PetscOptionsHasName(TAO_NULL, "-byelement", &flg); CHKERRQ(info);
157: if (flg) {
159: /* Sets routines for function and gradient evaluation, element by element */
160: info = PetscOptionsHasName(TAO_NULL, "-adic", &flg); CHKERRQ(info);
161: if (flg) {
162: info = DAAppSetADElementFunctionGradient(MSurfApp,ad_MSurfLocalFunction,150,(void *)&user.fgctx); CHKERRQ(info);
163: } else {
164: info = DAAppSetElementObjectiveAndGradientRoutine(MSurfApp,MSurfLocalFunctionGradient,36,(void *)&user); CHKERRQ(info);
165: }
167: /* Sets routines for Hessian evaluation, element by element */
168: info = DAAppSetElementHessianRoutine(MSurfApp,MSurfLocalHessian,87,(void*)&user); CHKERRQ(info);
170: } else {
172: /* Sets routines for function and gradient evaluation, all in one routine */
173: info = DAAppSetObjectiveAndGradientRoutine(MSurfApp,WholeMSurfFunctionGradient,(void *)&user); CHKERRQ(info);
175: /* Sets routines for Hessian evaluation, all in one routine */
176: info = DAAppSetHessianRoutine(MSurfApp,WholeMSurfHessian,(void*)&user); CHKERRQ(info);
177:
178: }
180: info = DAAppSetBeforeMonitor(MSurfApp,MyGridMonitorBefore,(void*)&user); CHKERRQ(info);
181: info = PetscOptionsHasName(TAO_NULL,"-tao_monitor", &flg); CHKERRQ(info);
182: if (flg){
183: info = DAAppPrintInterpolationError(MSurfApp); CHKERRQ(info);
184: info = DAAppPrintStageTimes(MSurfApp); CHKERRQ(info);
185: }
187: info = TaoAppSetRelativeTolerance(MSurfApp,1.0e-6); CHKERRQ(info);
188: info = TaoSetTolerances(tao,0,0,0,0); CHKERRQ(info);
189: info = TaoSetGradientTolerances(tao,0,0,0); CHKERRQ(info);
191: info = TaoAppGetKSP(MSurfApp,&ksp);CHKERRQ(info);
192: info = KSPSetType(ksp,KSPCG); CHKERRQ(info);
193: info = KSPSetTolerances(ksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,100);CHKERRQ(info);
194: info = KSPGetPC(ksp,&pc);CHKERRQ(info);
195: info = PCSetType(pc,PCBJACOBI);CHKERRQ(info);
197: /* Check for any tao command line options */
198: info = TaoSetOptions(MSurfApp,tao); CHKERRQ(info);
200: /* SOLVE THE APPLICATION */
201: info = TaoDAAppSolve(MSurfApp,tao); CHKERRQ(info);
203: /* Get information on termination */
204: info = TaoGetSolutionStatus(tao,&iter,&ff,&gnorm,0,0,&reason); CHKERRQ(info);
205: if (reason <= 0 ){
206: PetscPrintf(MPI_COMM_WORLD,"Try a different TAO method, adjust some parameters, or check the function evaluation routines\n");
207: PetscPrintf(MPI_COMM_WORLD," Iterations: %d, Function Value: %4.2e, Residual: %4.2e \n",iter,ff,gnorm);
208: }
210: info = PetscOptionsHasName(PETSC_NULL,"-view_sol",&flg); CHKERRQ(info);
211: if (flg){
212: info = DAAppGetSolution(MSurfApp,nlevels-1,&X); CHKERRQ(info);
213: info=VecView(X,PETSC_VIEWER_STDOUT_WORLD); CHKERRQ(info);
214: }
216: /* To View TAO solver information */
217: // info = TaoView(tao); CHKERRQ(info);
219: /* Free TAO data structures */
220: info = TaoDestroy(tao); CHKERRQ(info);
221: info = TaoAppDestroy(MSurfApp); CHKERRQ(info);
223: /* Free PETSc data structures */
224: for (iter=0;iter<nlevels;iter++){
225: info = DADestroy(DAarray[iter]); CHKERRQ(info);
226: }
228: PreLoadEnd();
230: /* Finalize TAO */
231: TaoFinalize();
232: PetscFinalize();
234: return 0;
235: } /* main */
239: /*----- The following two routines
240: MyGridMonitorBefore MyGridMonitorAfter
241: help diplay info of iterations at every grid level
242: */
246: static int MyGridMonitorBefore(TAO_APPLICATION myapp, DA da, int level, void *ctx) {
248: AppCtx *user = (AppCtx*)ctx;
249: int info,mx,my;
251: info = DAGetInfo(da,PETSC_NULL,&mx,&my,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL,
252: PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL);CHKERRQ(info);
253: user->mx = mx;
254: user->my = my;
255: user->hx = (user->r - user->l) / (user->mx - 1);
256: user->hy = (user->t - user->b) / (user->my - 1);
257: user->area = 0.5 * user->hx * user->hy;
258: user->fgctx.hx = user->hx;
259: user->fgctx.hy = user->hy;
260: user->fgctx.area = user->area;
262: user->bmx=(int)((mx+1)*user->fx); user->bmy=(int)((my+1)*user->fy);
264: PetscPrintf(MPI_COMM_WORLD,"Grid: %d, mx: %d my: %d \n",level,mx,my);
266: return 0;
267: }
271: /*
272: AppCtxInitialize - Sets initial values for the application context parameters
274: Input:
275: ptr - void user-defined application context
277: Output:
278: ptr - user-defined application context with the default or user-provided
279: parameters
280: */
281: static int AppCtxInitialize(void *ptr) {
283: AppCtx *user = (AppCtx*)ptr;
284: PetscTruth flg;
285: int info;
287: /* Specify default parameters */
288: user->mx = user->my = 11;
289: user->b = -0.5;
290: user->t = 0.5;
291: user->l = -0.5;
292: user->r = 0.5;
293: user->fx=0.5;
294: user->fy=0.5;
295: user->bheight=0.0;
297: /* Check for command line arguments that override defaults */
298: info = PetscOptionsGetInt(TAO_NULL, "-mx", &user->mx, &flg); CHKERRQ(info);
299: info = PetscOptionsGetInt(TAO_NULL, "-my", &user->my, &flg); CHKERRQ(info);
300: info = PetscOptionsGetReal(TAO_NULL, "-bottom", &user->b, &flg); CHKERRQ(info);
301: info = PetscOptionsGetReal(TAO_NULL, "-top", &user->t, &flg); CHKERRQ(info);
302: info = PetscOptionsGetReal(TAO_NULL, "-left", &user->l, &flg); CHKERRQ(info);
303: info = PetscOptionsGetReal(TAO_NULL, "-right", &user->r, &flg); CHKERRQ(info);
304: info = PetscOptionsGetReal(PETSC_NULL,"-bmx",&user->fx,&flg); CHKERRQ(info);
305: info = PetscOptionsGetReal(PETSC_NULL,"-bmy",&user->fy,&flg); CHKERRQ(info);
306: info = PetscOptionsGetReal(PETSC_NULL,"-bheight",&user->bheight,&flg); CHKERRQ(info);
308: user->hx = (user->r - user->l) / (user->mx - 1);
309: user->hy = (user->t - user->b) / (user->my - 1);
310: user->area = 0.5 * user->hx * user->hy;
311: info = PetscLogFlops(8); CHKERRQ(info);
313: return 0;
315: } /* AppCtxInitialize */
318: /*------- USER-DEFINED: set the upper and lower bounds for the variables -------*/
322: /*
323: FormBounds - Forms bounds on the variables
325: Input:
326: user - user-defined application context
328: Output:
329: XL - vector of lower bounds
330: XU - vector of upper bounds
332: */
333: static int COPS_Bounds(TAO_APPLICATION tao, DA da, Vec XL, Vec XU, void *ptr) {
335: AppCtx *user = (AppCtx*)ptr;
336: PetscTruth flg;
337: int i, j, mx, my, info, xs, xm, ys, ym;
338: double lb = -TAO_INFINITY;
339: double ub = TAO_INFINITY;
340: double **xl, **xu;
341: double xi, xi1, xi2;
342: double cx, cy, radx, rady, height = 1.0;
344: info = DAGetInfo(da,PETSC_NULL,&mx,&my,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL,
345: PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL);CHKERRQ(info);
346: user->hx = (user->r - user->l) / (mx - 1);
347: user->hy = (user->t - user->b) / (my - 1);
348: user->area = 0.5 * user->hx * user->hy;
350: info = DAVecGetArray(da, XL, (void**)&xl); CHKERRQ(info);
351: info = DAVecGetArray(da, XU, (void**)&xu); CHKERRQ(info);
352: info = DAGetCorners(da, &xs, &ys, TAO_NULL, &xm, &ym, TAO_NULL); CHKERRQ(info);
354: /* Compute default bounds */
355: for (j = ys; j < ys+ym; j++){
356: for (i = xs; i < xs+xm; i++){
358: if (j == 0 || j == my - 1) {
359: xi = user->l + i * user->hx;
360: xl[j][i] = xu[j][i] = -4 * (xi - user->l) * (xi - user->r);
361: } else if (i == 0 || i == mx - 1) {
362: xl[j][i] = xu[j][i] = 0.0;
363: } else {
364: xl[j][i] = lb;
365: xu[j][i] = ub;
366: }
368: }
369: }
371: /* Adjust lower bound if obstacle is present */
372: info = PetscOptionsHasName(PETSC_NULL, "-obst", &flg); CHKERRQ(info);
373: if (flg) {
374: radx = (user->r - user->l) * 0.25;
375: cx = user->l + 2.0 * radx;
376: rady = (user->t - user->b) * 0.25;
377: cy = user->b + 2.0 * rady;
378: for (j = ys; j < ys+ym; j++){
379: for (i = xs; i < xs+xm; i++){
380: xi1 = user->l + i * user->hx;
381: xi2 = user->b + j * user->hy;
382: if ( fabs(xi1 - cx) <= radx && fabs(xi2 - cy) <= rady ) {
383: xl[j][i] = height;
384: }
385: }
386: }
387: info = PetscLogFlops(8 + xm * ym * 6); CHKERRQ(info);
388: }
390: info = DAVecRestoreArray(da, XL, (void**)&xl); CHKERRQ(info);
391: info = DAVecRestoreArray(da, XU, (void**)&xu); CHKERRQ(info);
393: info = PetscLogFlops(12 * ym + 6); CHKERRQ(info);
395: return 0;
397: } /* DAGetBounds2d */
402: static int MINPACK_Bounds(TAO_APPLICATION tao, DA da, Vec XL,Vec XU, void *user1){
404: AppCtx *user=(AppCtx*)user1;
405: int i,j,k,limit=0,info,maxits=5;
406: int xs,ys,xm,ym;
407: int mx, my, bmy, bmx;
408: int row=0, bsize=0, lsize=0, tsize=0, rsize=0;
409: double bheight=user->bheight;
410: double one=1.0, two=2.0, three=3.0, tol=1e-10;
411: double fnorm,det,hx,hy,xt=0,yt=0;
412: double u1,u2,nf1,nf2,njac11,njac12,njac21,njac22;
413: double b=-0.5, t=0.5, l=-0.5, r=0.5;
414: PetscScalar scl,*xl,*xu, lb=TAO_NINFINITY, ub=TAO_INFINITY;
415: PetscTruth flg;
416: double **xll;
418: info = DAGetInfo(da,PETSC_NULL,&mx,&my,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL,
419: PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL);CHKERRQ(info);
421: bheight=user->bheight,lb=-1000; ub=1000;
422: bmx=user->bmx; bmy=user->bmy;
424: info = DAGetCorners(da,&xs,&ys,PETSC_NULL,&xm,&ym,PETSC_NULL); CHKERRQ(info);
426: info = VecSet(XL, lb); CHKERRQ(info);
427: info = VecSet(XU, ub); CHKERRQ(info);
429: /*
430: Get pointers to vector data
431: */
432: info = DAVecGetArray(da,XL,(void**)&xll); CHKERRQ(info);
434: if (ys==0) bsize=xm;
435: if (xs==0) lsize=ym;
436: if (xs+xm==mx) rsize=ym;
437: if (ys+ym==my) tsize=xm;
438:
439: hx= (r-l)/(mx-1); hy=(t-b)/(my-1);
441: /* Compute the optional lower box */
442: for (i=xs; i< xs+xm; i++){
443: for (j=ys; j<ys+ym; j++){
444:
445: if (i>=(mx-bmx)/2 && i<mx-(mx-bmx)/2 && j>=(my-bmy)/2 && j<my-(my-bmy)/2 ){
446: xll[j][i] = bheight;
447: }
449: }
450: }
452: info = DAVecRestoreArray(da,XL,(void**)&xll); CHKERRQ(info);
454: /* Boundary Values */
455: info = VecGetArray(XL,&xl); CHKERRQ(info);
456: info = VecGetArray(XU,&xu); CHKERRQ(info);
458: for (j=0; j<4; j++){
459: if (j==0){
460: yt=b;
461: xt=l+hx*xs;
462: limit=bsize;
463: scl=1.0;
464: info = PetscOptionsGetReal(PETSC_NULL,"-bottom",&scl,&flg);
465: } else if (j==1){
466: yt=t;
467: xt=l+hx*xs;
468: limit=tsize;
469: scl=1.0;
470: info = PetscOptionsGetReal(PETSC_NULL,"-top",&scl,&flg);
471: } else if (j==2){
472: yt=b+hy*ys;
473: xt=l;
474: limit=lsize;
475: scl=1.0;
476: info = PetscOptionsGetReal(PETSC_NULL,"-left",&scl,&flg);
477: } else if (j==3){
478: yt=b+hy*ys;
479: xt=r;
480: limit=rsize;
481: scl=1.0;
482: info = PetscOptionsGetReal(PETSC_NULL,"-right",&scl,&flg);
483: }
485: for (i=0; i<limit; i++){
486: u1=xt;
487: u2=-yt;
488: for (k=0; k<maxits; k++){
489: nf1=u1 + u1*u2*u2 - u1*u1*u1/three-xt;
490: nf2=-u2 - u1*u1*u2 + u2*u2*u2/three-yt;
491: fnorm=sqrt(nf1*nf1+nf2*nf2);
492: if (fnorm <= tol) break;
493: njac11=one+u2*u2-u1*u1;
494: njac12=two*u1*u2;
495: njac21=-two*u1*u2;
496: njac22=-one - u1*u1 + u2*u2;
497: det = njac11*njac22-njac21*njac12;
498: u1 = u1-(njac22*nf1-njac12*nf2)/det;
499: u2 = u2-(njac11*nf2-njac21*nf1)/det;
500: }
502: if (j==0){
503: row = i;
504: } else if (j==1){
505: row = (ym-1)*xm+i;
506: } else if (j==2){
507: row = (xm*i);
508: } else if (j==3){
509: row = xm*(i+1)-1;
510: }
511:
512: xl[row]=(u1*u1-u2*u2)*scl;
513: xu[row]=(u1*u1-u2*u2)*scl;
515: if (j==0 || j==1) {
516: xt=xt+hx;
517: } else if (j==2 || j==3){
518: yt=yt+hy;
519: }
520:
521: }
522:
523: }
525: info = VecRestoreArray(XL,&xl); CHKERRQ(info);
526: info = VecRestoreArray(XU,&xu); CHKERRQ(info);
528: return 0;
529: }
531: /*------- USER-DEFINED: routine to evaluate the function and gradient
532: at a local (rectangular element) level -------*/
536: /*
537: MSurfLocalFunctionGradient - Evaluates function and gradient over the
538: local rectangular element
540: Input:
541: coor - vector with the indices of the position of current element
542: in the first, second and third directions
543: x - current point (values over the current rectangular element)
544: ptr - user-defined application context
546: Output:
547: f - value of the objective funtion at the local rectangular element
548: g - gradient of the local function
550: */
551: static int MSurfLocalFunctionGradient(int coor[2], double x[4], double *f, double g[4], void *ptr) {
553: AppCtx *user = (AppCtx*)ptr;
555: double hx, hy, area;
556: double dvdx, dvdy, flow, fup;
557: double d1,d2;
559: hx = user->hx;
560: hy = user->hy;
561: area = user->area;
563: /* lower triangle contribution */
564: dvdx = (x[0] - x[1]);
565: dvdy = (x[0] - x[2]);
566: g[1] = (dvdx * (hy/hx))/(-2);
567: g[2] = (dvdy * (hx/hy))/(-2);
568: dvdx /= hx;
569: dvdy /= hy;
570: flow = sqrt( 1 + dvdx * dvdx + dvdy * dvdy );
571: *f=flow;
572: g[1] /= flow;
573: g[2] /= flow;
574: g[0] = -(g[1]+g[2]);
576: /* upper triangle contribution */
577: dvdx = (x[3] - x[2]);
578: dvdy = (x[3] - x[1]);
579: d1 = (dvdy*(hx/hy))/(-2);
580: d2 = (dvdx*(hy/hx))/(-2);
581: dvdx /= hx;
582: dvdy /= hy;
583: fup = sqrt( 1 + dvdx * dvdx + dvdy * dvdy );
584: *f += fup;
585: g[1] += d1/fup;
586: g[2] += d2/fup;
587: g[3] = -(d1+d2)/fup;
589: *f *= area;
591: return 0;
592: } /* MSurfLocalFunctionGradient */
598: /*
599: MSurfLocalHessian - Computes the Hessian of the local (partial) function
600: defined over the current rectangle
602: Input:
603: coor - vector with the indices of the position of current element
604: in the first, second and third directions
605: x - current local solution (over the rectangle only)
606: ptr - user-defined application context
608: Output:
609: H - Hessian matrix of the local function (wrt the four
610: points of the rectangle only)
612: */
613: static int MSurfLocalHessian(int coor[2], double x[4], double H[4][4],void *ptr) {
615: AppCtx *user = (AppCtx*)ptr;
616: double hx, hy, area, byhxhx, byhyhy;
617: double dvdx, dvdy, flow, fup;
618: double areadivf, areadivf3;
620: hx = user->hx;
621: hy = user->hy;
622: area = user->area;
623:
624: byhxhx = 1.0 / (hx * hx);
625: byhyhy = 1.0 / (hy * hy);
627: /* 0 is 0,0; 1 is 1,0; 2 is 0,1; 3 is 1,1 */
628: dvdx = (x[0] - x[1]) / hx; /* lower triangle contribution */
629: dvdy = (x[0] - x[2]) / hy;
630: flow = sqrt( 1 + dvdx * dvdx + dvdy * dvdy );
631: dvdx = dvdx / hx;
632: dvdy = dvdy / hy;
633: areadivf = area / flow;
634: areadivf3 = areadivf / (flow * flow);
635: H[0][0] = areadivf * (byhxhx + byhyhy) - areadivf3 * (dvdx + dvdy) * (dvdx + dvdy);
636: H[0][1] = areadivf * (-byhxhx) + areadivf3 * (dvdx + dvdy) * (dvdx);
637: H[0][2] = areadivf * (-byhyhy) + areadivf3 * (dvdx + dvdy) * (dvdy);
638: H[0][3] = 0.0;
639: H[1][1] = areadivf * byhxhx - areadivf3 * dvdx * dvdx;
640: H[1][2] = areadivf3 * (-dvdx) * dvdy;
641: H[2][2] = areadivf * byhyhy - areadivf3 * dvdy * dvdy;
643: /* upper triangle contribution */
644: dvdx = (x[3] - x[2]) / hx;
645: dvdy = (x[3] - x[1]) / hy;
646: fup = sqrt( 1 + dvdx * dvdx + dvdy * dvdy );
647: dvdx = dvdx / hx;
648: dvdy = dvdy / hy;
649: areadivf = area / fup;
650: areadivf3 = areadivf / (fup * fup);
651: H[1][1] += areadivf * byhyhy - areadivf3 * dvdy * dvdy;
652: H[1][2] += areadivf3 * (-dvdy) * dvdx;
653: H[2][2] += areadivf * byhxhx - areadivf3 * (dvdx * dvdx);
654: H[1][3] = areadivf * (-byhyhy) + areadivf3 * (dvdx + dvdy) * dvdy;
655: H[2][3] = areadivf * (-byhxhx) + areadivf3 * (dvdx + dvdy) * dvdx;
656: H[3][3] = areadivf * (byhxhx + byhyhy) - areadivf3 * (dvdx + dvdy) * (dvdx + dvdy);
658: H[1][0] = H[0][1];
659: H[2][0] = H[0][2];
660: H[3][0] = H[0][3];
661: H[2][1] = H[1][2];
662: H[3][1] = H[1][3];
663: H[3][2] = H[2][3];
665: return 0;
667: } /* MSurfLocalHessian */
670: /*------- USER-DEFINED: routine to evaluate the function
671: and gradient at the whole grid -------*/
675: /*
676: WholeMSurfFunctionGradient - Evaluates function and gradient over the
677: whole grid
679: Input:
680: daapplication - TAO application object
681: da - distributed array
682: X - the current point, at which the function and gradient are evaluated
683: ptr - user-defined application context
685: Output:
686: f - value of the objective funtion at X
687: G - gradient at X
688: */
689: static int WholeMSurfFunctionGradient(TAO_APPLICATION daapplication, DA da, Vec X, double *f, Vec G, void *ptr) {
691: AppCtx *user = (AppCtx*)ptr;
692: Vec localX, localG;
693: int info, i, j;
694: int xs, xm, gxs, gxm, xe, ys, ym, gys, gym, ye;
695: double **x, **g;
696: double floc = 0.0;
697: PetscScalar zero = 0.0;
699: double hx, hy, area;
700: double dvdx, dvdy, flow, fup;
701: double areadivf;
703: hx = user->hx;
704: hy = user->hy;
705: area = user->area;
707: info = DAGetLocalVector(da, &localX); CHKERRQ(info);
708: info = DAGetLocalVector(da, &localG); CHKERRQ(info);
709: info = VecSet(G, zero); CHKERRQ(info);
710: info = VecSet(localG, zero); CHKERRQ(info);
712: info = DAGlobalToLocalBegin(da, X, INSERT_VALUES, localX); CHKERRQ(info);
713: info = DAGlobalToLocalEnd(da, X, INSERT_VALUES, localX); CHKERRQ(info);
715: info = DAVecGetArray(da, localX, (void**)&x); CHKERRQ(info);
716: info = DAVecGetArray(da, localG, (void**)&g); CHKERRQ(info);
718: info = DAGetCorners(da, &xs, &ys, TAO_NULL, &xm, &ym, TAO_NULL); CHKERRQ(info);
719: info = DAGetGhostCorners(da, &gxs, &gys, TAO_NULL, &gxm, &gym, TAO_NULL); CHKERRQ(info);
721: xe = gxs + gxm - 1;
722: ye = gys + gym - 1;
723: for (j = ys; j < ye; j++) {
724: for (i = xs; i < xe; i++) {
726: /* lower triangle contribution */
727: dvdx = (x[j][i] - x[j][i+1]) / hx;
728: dvdy = (x[j][i] - x[j+1][i]) / hy;
729: flow = sqrt( 1 + dvdx * dvdx + dvdy * dvdy );
730: areadivf = area / flow;
731: g[j][i] += (dvdx / hx + dvdy / hy) * areadivf;
732: g[j][i+1] += (-dvdx / hx) * areadivf;
733: g[j+1][i] += (-dvdy / hy) * areadivf;
735: /* upper triangle contribution */
736: dvdx = (x[j+1][i+1] - x[j+1][i]) / hx;
737: dvdy = (x[j+1][i+1] - x[j][i+1]) / hy;
738: fup = sqrt( 1 + dvdx * dvdx + dvdy * dvdy );
739: areadivf = area / fup;
740: g[j][i+1] += (-dvdy / hy) * areadivf;
741: g[j+1][i] += (-dvdx / hx) * areadivf;
742: g[j+1][i+1] += (dvdx / hx + dvdy / hy) * areadivf;
744: floc += area * (flow + fup);
746: }
747: }
749: info = MPI_Allreduce(&floc, f, 1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD); CHKERRQ(info);
751: info = DAVecRestoreArray(da, localX, (void**)&x); CHKERRQ(info);
752: info = DAVecRestoreArray(da, localG, (void**)&g); CHKERRQ(info);
754: info = DALocalToGlobalBegin(da, localG, G); CHKERRQ(info);
755: info = DALocalToGlobalEnd(da, localG, G); CHKERRQ(info);
757: info = DARestoreLocalVector(da, &localX); CHKERRQ(info);
758: info = DARestoreLocalVector(da, &localG); CHKERRQ(info);
760: info = PetscLogFlops((xe-xs) * (ye-ys) * 42); CHKERRQ(info);
762: return 0;
763: } /* WholeMSurfFunctionGradient */
766: /*------- USER-DEFINED: routine to evaluate the Hessian
767: at the whole grid -------*/
771: /*
772: WholeMSurfHessian - Evaluates Hessian over the whole grid
774: Input:
775: daapplication - TAO application object
776: da - distributed array
777: X - the current point, at which the function and gradient are evaluated
778: ptr - user-defined application context
780: Output:
781: H - Hessian at X
782: */
783: static int WholeMSurfHessian(TAO_APPLICATION daapplication, DA da, Vec X, Mat H, void *ptr) {
785: AppCtx *user = (AppCtx*)ptr;
786: Vec localX;
787: int info, i, j, ind[4];
788: int xs, xm, gxs, gxm, xe, ys, ym, gys, gym, ye;
789: double smallH[4][4];
790: double **x;
792: double hx, hy, area, byhxhx, byhyhy;
793: double dvdx, dvdy, flow, fup;
794: double areadivf, areadivf3;
795: PetscTruth assembled;
797: hx = user->hx;
798: hy = user->hy;
799: area = user->area;
800:
801: byhxhx = 1.0 / (hx * hx);
802: byhyhy = 1.0 / (hy * hy);
804: info = DAGetLocalVector(da, &localX); CHKERRQ(info);
805: info = MatAssembled(H,&assembled); CHKERRQ(info);
806: if (assembled){info = MatZeroEntries(H); CHKERRQ(info);}
808: info = DAGlobalToLocalBegin(da, X, INSERT_VALUES, localX); CHKERRQ(info);
809: info = DAGlobalToLocalEnd(da, X, INSERT_VALUES, localX); CHKERRQ(info);
811: info = DAVecGetArray(da, localX, (void**)&x); CHKERRQ(info);
813: info = DAGetCorners(da, &xs, &ys, TAO_NULL, &xm, &ym, TAO_NULL); CHKERRQ(info);
814: info = DAGetGhostCorners(da, &gxs, &gys, TAO_NULL, &gxm, &gym, TAO_NULL); CHKERRQ(info);
816: xe = gxs + gxm - 1;
817: ye = gys + gym - 1;
818: for (j = ys; j < ye; j++) {
819: for (i = xs; i < xe; i++) {
821: /* 0 is 0,0; 1 is 1,0; 2 is 0,1; 3 is 1,1 */
822: dvdx = (x[j][i] - x[j][i+1]) / hx; /* lower triangle contribution */
823: dvdy = (x[j][i] - x[j+1][i]) / hy;
824: flow = sqrt( 1 + dvdx * dvdx + dvdy * dvdy );
825: dvdx = dvdx / hx;
826: dvdy = dvdy / hy;
827: areadivf = area / flow;
828: areadivf3 = areadivf / (flow * flow);
829: smallH[0][0] = areadivf * (byhxhx + byhyhy) - areadivf3 * (dvdx + dvdy) * (dvdx + dvdy);
830: smallH[0][1] = areadivf * (-byhxhx) + areadivf3 * (dvdx + dvdy) * (dvdx);
831: smallH[0][2] = areadivf * (-byhyhy) + areadivf3 * (dvdx + dvdy) * (dvdy);
832: smallH[0][3] = 0.0;
833: smallH[1][1] = areadivf * byhxhx - areadivf3 * dvdx * dvdx;
834: smallH[1][2] = areadivf3 * (-dvdx) * dvdy;
835: smallH[2][2] = areadivf * byhyhy - areadivf3 * dvdy * dvdy;
837: /* upper triangle contribution */
838: dvdx = (x[j+1][i+1] - x[j+1][i]) / hx;
839: dvdy = (x[j+1][i+1] - x[j][i+1]) / hy;
840: fup = sqrt( 1 + dvdx * dvdx + dvdy * dvdy );
841: dvdx = dvdx / hx;
842: dvdy = dvdy / hy;
843: areadivf = area / fup;
844: areadivf3 = areadivf / (fup * fup);
845: smallH[1][1] += areadivf * byhyhy - areadivf3 * dvdy * dvdy;
846: smallH[1][2] += areadivf3 * (-dvdy) * dvdx;
847: smallH[2][2] += areadivf * byhxhx - areadivf3 * (dvdx * dvdx);
848: smallH[1][3] = areadivf * (-byhyhy) + areadivf3 * (dvdx + dvdy) * dvdy;
849: smallH[2][3] = areadivf * (-byhxhx) + areadivf3 * (dvdx + dvdy) * dvdx;
850: smallH[3][3] = areadivf * (byhxhx + byhyhy) - areadivf3 * (dvdx + dvdy) * (dvdx + dvdy);
852: smallH[1][0] = smallH[0][1];
853: smallH[2][0] = smallH[0][2];
854: smallH[3][0] = smallH[0][3];
855: smallH[2][1] = smallH[1][2];
856: smallH[3][1] = smallH[1][3];
857: smallH[3][2] = smallH[2][3];
859: ind[0] = (j-gys) * gxm + (i-gxs);
860: ind[1] = ind[0] + 1;
861: ind[2] = ind[0] + gxm;
862: ind[3] = ind[2] + 1;
863: info = MatSetValuesLocal(H,4,ind,4,ind,(PetscScalar*)smallH,ADD_VALUES); CHKERRQ(info);
865: }
866: }
868: info = DAVecRestoreArray(da, localX, (void**)&x); CHKERRQ(info);
870: info = MatAssemblyBegin(H, MAT_FINAL_ASSEMBLY); CHKERRQ(info);
871: info = MatAssemblyEnd(H, MAT_FINAL_ASSEMBLY); CHKERRQ(info);
872: info = MatSetOption(H, MAT_SYMMETRIC); CHKERRQ(info);
874: info = DARestoreLocalVector(da, &localX); CHKERRQ(info);
876: info = PetscLogFlops((xe-xs) * (ye-ys) * 83 + 4); CHKERRQ(info);
877: return 0;
879: } /* WholeMSurfHessian */