Actual source code: plate2.c

petsc-3.7.3 2016-08-01
Report Typos and Errors
  1: #include <petscdmda.h>
  2: #include <petsctao.h>

  4: static  char help[] =
  5: "This example demonstrates use of the TAO package to \n\
  6: solve an unconstrained minimization problem.  This example is based on a \n\
  7: problem from the MINPACK-2 test suite.  Given a rectangular 2-D domain, \n\
  8: boundary values along the edges of the domain, and a plate represented by \n\
  9: lower boundary conditions, the objective is to find the\n\
 10: surface with the minimal area that satisfies the boundary conditions.\n\
 11: The command line options are:\n\
 12:   -mx <xg>, where <xg> = number of grid points in the 1st coordinate direction\n\
 13:   -my <yg>, where <yg> = number of grid points in the 2nd coordinate direction\n\
 14:   -bmx <bxg>, where <bxg> = number of grid points under plate in 1st direction\n\
 15:   -bmy <byg>, where <byg> = number of grid points under plate in 2nd direction\n\
 16:   -bheight <ht>, where <ht> = height of the plate\n\
 17:   -start <st>, where <st> =0 for zero vector, <st> >0 for random start, and <st> <0 \n\
 18:                for an average of the boundary conditions\n\n";

 20: /*T
 21:    Concepts: TAO^Solving a bound constrained minimization problem
 22:    Routines: TaoCreate();
 23:    Routines: TaoSetType(); TaoSetObjectiveAndGradientRoutine();
 24:    Routines: TaoSetHessianRoutine();
 25:    Routines: TaoSetInitialVector();
 26:    Routines: TaoSetVariableBounds();
 27:    Routines: TaoSetFromOptions();
 28:    Routines: TaoSolve(); TaoView();
 29:    Routines: TaoDestroy();
 30:    Processors: n
 31: T*/


 34: /*
 35:    User-defined application context - contains data needed by the
 36:    application-provided call-back routines, FormFunctionGradient(),
 37:    FormHessian().
 38: */
 39: typedef struct {
 40:   /* problem parameters */
 41:   PetscReal      bheight;                  /* Height of plate under the surface */
 42:   PetscInt       mx, my;                   /* discretization in x, y directions */
 43:   PetscInt       bmx,bmy;                  /* Size of plate under the surface */
 44:   Vec            Bottom, Top, Left, Right; /* boundary values */

 46:   /* Working space */
 47:   Vec         localX, localV;           /* ghosted local vector */
 48:   DM          dm;                       /* distributed array data structure */
 49:   Mat         H;
 50: } AppCtx;

 52: /* -------- User-defined Routines --------- */

 54: static PetscErrorCode MSA_BoundaryConditions(AppCtx*);
 55: static PetscErrorCode MSA_InitialPoint(AppCtx*,Vec);
 56: static PetscErrorCode MSA_Plate(Vec,Vec,void*);
 57: PetscErrorCode FormFunctionGradient(Tao,Vec,PetscReal*,Vec,void*);
 58: PetscErrorCode FormHessian(Tao,Vec,Mat,Mat,void*);

 60: /* For testing matrix free submatrices */
 61: PetscErrorCode MatrixFreeHessian(Tao,Vec,Mat, Mat,void*);
 62: PetscErrorCode MyMatMult(Mat,Vec,Vec);

 66: int main( int argc, char **argv )
 67: {
 68:   PetscErrorCode         ierr;                 /* used to check for functions returning nonzeros */
 69:   PetscInt               Nx, Ny;               /* number of processors in x- and y- directions */
 70:   PetscInt               m, N;                 /* number of local and global elements in vectors */
 71:   Vec                    x,xl,xu;               /* solution vector  and bounds*/
 72:   PetscBool              flg;                /* A return variable when checking for user options */
 73:   Tao                    tao;                  /* Tao solver context */
 74:   ISLocalToGlobalMapping isltog;   /* local-to-global mapping object */
 75:   Mat                    H_shell;                  /* to test matrix-free submatrices */
 76:   AppCtx                 user;                 /* user-defined work context */

 78:   /* Initialize PETSc, TAO */
 79:   PetscInitialize( &argc, &argv,(char *)0,help );

 81:   /* Specify default dimension of the problem */
 82:   user.mx = 10; user.my = 10; user.bheight=0.1;

 84:   /* Check for any command line arguments that override defaults */
 85:   PetscOptionsGetInt(NULL,NULL,"-mx",&user.mx,&flg);
 86:   PetscOptionsGetInt(NULL,NULL,"-my",&user.my,&flg);
 87:   PetscOptionsGetReal(NULL,NULL,"-bheight",&user.bheight,&flg);

 89:   user.bmx = user.mx/2; user.bmy = user.my/2;
 90:   PetscOptionsGetInt(NULL,NULL,"-bmx",&user.bmx,&flg);
 91:   PetscOptionsGetInt(NULL,NULL,"-bmy",&user.bmy,&flg);

 93:   PetscPrintf(PETSC_COMM_WORLD,"\n---- Minimum Surface Area With Plate Problem -----\n");
 94:   PetscPrintf(PETSC_COMM_WORLD,"mx:%D, my:%D, bmx:%D, bmy:%D, height:%g\n",user.mx,user.my,user.bmx,user.bmy,(double)user.bheight);

 96:   /* Calculate any derived values from parameters */
 97:   N    = user.mx*user.my;

 99:   /* Let Petsc determine the dimensions of the local vectors */
100:   Nx = PETSC_DECIDE; Ny = PETSC_DECIDE;

102:   /*
103:      A two dimensional distributed array will help define this problem,
104:      which derives from an elliptic PDE on two dimensional domain.  From
105:      the distributed array, Create the vectors.
106:   */
107:   DMDACreate2d(MPI_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,
108:                       DMDA_STENCIL_BOX,user.mx,user.my,Nx,Ny,1,1,
109:                       NULL,NULL,&user.dm);

111:   /*
112:      Extract global and local vectors from DM; The local vectors are
113:      used solely as work space for the evaluation of the function,
114:      gradient, and Hessian.  Duplicate for remaining vectors that are
115:      the same types.
116:   */
117:   DMCreateGlobalVector(user.dm,&x); /* Solution */
118:   DMCreateLocalVector(user.dm,&user.localX);
119:   VecDuplicate(user.localX,&user.localV);

121:   VecDuplicate(x,&xl);
122:   VecDuplicate(x,&xu);



126:   /* The TAO code begins here */

128:   /*
129:      Create TAO solver and set desired solution method
130:      The method must either be TAOTRON or TAOBLMVM
131:      If TAOBLMVM is used, then hessian function is not called.
132:   */
133:   TaoCreate(PETSC_COMM_WORLD,&tao);
134:   TaoSetType(tao,TAOBLMVM);

136:   /* Set initial solution guess; */
137:   MSA_BoundaryConditions(&user);
138:   MSA_InitialPoint(&user,x);
139:   TaoSetInitialVector(tao,x);

141:   /* Set routines for function, gradient and hessian evaluation */
142:   TaoSetObjectiveAndGradientRoutine(tao,FormFunctionGradient,(void*) &user);

144:   VecGetLocalSize(x,&m);
145:   MatCreateAIJ(MPI_COMM_WORLD,m,m,N,N,7,NULL,3,NULL,&(user.H));
146:   MatSetOption(user.H,MAT_SYMMETRIC,PETSC_TRUE);

148:   DMGetLocalToGlobalMapping(user.dm,&isltog);
149:   MatSetLocalToGlobalMapping(user.H,isltog,isltog);
150:   PetscOptionsHasName(NULL,NULL,"-matrixfree",&flg);
151:   if (flg) {
152:       MatCreateShell(PETSC_COMM_WORLD,m,m,N,N,(void*)&user,&H_shell);
153:       MatShellSetOperation(H_shell,MATOP_MULT,(void(*)(void))MyMatMult);
154:       MatSetOption(H_shell,MAT_SYMMETRIC,PETSC_TRUE);
155:       TaoSetHessianRoutine(tao,H_shell,H_shell,MatrixFreeHessian,(void*)&user);
156:   } else {
157:       TaoSetHessianRoutine(tao,user.H,user.H,FormHessian,(void*)&user);
158:   }

160:   /* Set Variable bounds */
161:   MSA_Plate(xl,xu,(void*)&user);
162:   TaoSetVariableBounds(tao,xl,xu);

164:   /* Check for any tao command line options */
165:   TaoSetFromOptions(tao);

167:   /* SOLVE THE APPLICATION */
168:   TaoSolve(tao);

170:   TaoView(tao,PETSC_VIEWER_STDOUT_WORLD);

172:   /* Free TAO data structures */
173:   TaoDestroy(&tao);

175:   /* Free PETSc data structures */
176:   VecDestroy(&x);
177:   VecDestroy(&xl);
178:   VecDestroy(&xu);
179:   MatDestroy(&user.H);
180:   VecDestroy(&user.localX);
181:   VecDestroy(&user.localV);
182:   VecDestroy(&user.Bottom);
183:   VecDestroy(&user.Top);
184:   VecDestroy(&user.Left);
185:   VecDestroy(&user.Right);
186:   DMDestroy(&user.dm);
187:   if (flg) {
188:     MatDestroy(&H_shell);
189:   }
190:   PetscFinalize();
191:   return 0;
192: }

196: /*  FormFunctionGradient - Evaluates f(x) and gradient g(x).

198:     Input Parameters:
199: .   tao     - the Tao context
200: .   X      - input vector
201: .   userCtx - optional user-defined context, as set by TaoSetObjectiveAndGradientRoutine()

203:     Output Parameters:
204: .   fcn     - the function value
205: .   G      - vector containing the newly evaluated gradient

207:    Notes:
208:    In this case, we discretize the domain and Create triangles. The
209:    surface of each triangle is planar, whose surface area can be easily
210:    computed.  The total surface area is found by sweeping through the grid
211:    and computing the surface area of the two triangles that have their
212:    right angle at the grid point.  The diagonal line segments on the
213:    grid that define the triangles run from top left to lower right.
214:    The numbering of points starts at the lower left and runs left to
215:    right, then bottom to top.
216: */
217: PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *fcn, Vec G,void *userCtx)
218: {
219:   AppCtx         *user = (AppCtx *) userCtx;
221:   PetscInt       i,j,row;
222:   PetscInt       mx=user->mx, my=user->my;
223:   PetscInt       xs,xm,gxs,gxm,ys,ym,gys,gym;
224:   PetscReal      ft=0;
225:   PetscReal      zero=0.0;
226:   PetscReal      hx=1.0/(mx+1),hy=1.0/(my+1), hydhx=hy/hx, hxdhy=hx/hy, area=0.5*hx*hy;
227:   PetscReal      rhx=mx+1, rhy=my+1;
228:   PetscReal      f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb;
229:   PetscReal      df1dxc,df2dxc,df3dxc,df4dxc,df5dxc,df6dxc;
230:   PetscReal      *g, *x,*left,*right,*bottom,*top;
231:   Vec            localX = user->localX, localG = user->localV;

233:   /* Get local mesh boundaries */
234:   DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL);
235:   DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL);

237:   /* Scatter ghost points to local vector */
238:   DMGlobalToLocalBegin(user->dm,X,INSERT_VALUES,localX);
239:   DMGlobalToLocalEnd(user->dm,X,INSERT_VALUES,localX);

241:   /* Initialize vector to zero */
242:   VecSet(localG, zero);

244:   /* Get pointers to vector data */
245:   VecGetArray(localX,&x);
246:   VecGetArray(localG,&g);
247:   VecGetArray(user->Top,&top);
248:   VecGetArray(user->Bottom,&bottom);
249:   VecGetArray(user->Left,&left);
250:   VecGetArray(user->Right,&right);

252:   /* Compute function over the locally owned part of the mesh */
253:   for (j=ys; j<ys+ym; j++){
254:     for (i=xs; i< xs+xm; i++){
255:       row=(j-gys)*gxm + (i-gxs);

257:       xc = x[row];
258:       xlt=xrb=xl=xr=xb=xt=xc;

260:       if (i==0){ /* left side */
261:         xl= left[j-ys+1];
262:         xlt = left[j-ys+2];
263:       } else {
264:         xl = x[row-1];
265:       }

267:       if (j==0){ /* bottom side */
268:         xb=bottom[i-xs+1];
269:         xrb = bottom[i-xs+2];
270:       } else {
271:         xb = x[row-gxm];
272:       }

274:       if (i+1 == gxs+gxm){ /* right side */
275:         xr=right[j-ys+1];
276:         xrb = right[j-ys];
277:       } else {
278:         xr = x[row+1];
279:       }

281:       if (j+1==gys+gym){ /* top side */
282:         xt=top[i-xs+1];
283:         xlt = top[i-xs];
284:       }else {
285:         xt = x[row+gxm];
286:       }

288:       if (i>gxs && j+1<gys+gym){
289:         xlt = x[row-1+gxm];
290:       }
291:       if (j>gys && i+1<gxs+gxm){
292:         xrb = x[row+1-gxm];
293:       }

295:       d1 = (xc-xl);
296:       d2 = (xc-xr);
297:       d3 = (xc-xt);
298:       d4 = (xc-xb);
299:       d5 = (xr-xrb);
300:       d6 = (xrb-xb);
301:       d7 = (xlt-xl);
302:       d8 = (xt-xlt);

304:       df1dxc = d1*hydhx;
305:       df2dxc = ( d1*hydhx + d4*hxdhy );
306:       df3dxc = d3*hxdhy;
307:       df4dxc = ( d2*hydhx + d3*hxdhy );
308:       df5dxc = d2*hydhx;
309:       df6dxc = d4*hxdhy;

311:       d1 *= rhx;
312:       d2 *= rhx;
313:       d3 *= rhy;
314:       d4 *= rhy;
315:       d5 *= rhy;
316:       d6 *= rhx;
317:       d7 *= rhy;
318:       d8 *= rhx;

320:       f1 = PetscSqrtScalar( 1.0 + d1*d1 + d7*d7);
321:       f2 = PetscSqrtScalar( 1.0 + d1*d1 + d4*d4);
322:       f3 = PetscSqrtScalar( 1.0 + d3*d3 + d8*d8);
323:       f4 = PetscSqrtScalar( 1.0 + d3*d3 + d2*d2);
324:       f5 = PetscSqrtScalar( 1.0 + d2*d2 + d5*d5);
325:       f6 = PetscSqrtScalar( 1.0 + d4*d4 + d6*d6);

327:       ft = ft + (f2 + f4);

329:       df1dxc /= f1;
330:       df2dxc /= f2;
331:       df3dxc /= f3;
332:       df4dxc /= f4;
333:       df5dxc /= f5;
334:       df6dxc /= f6;

336:       g[row] = (df1dxc+df2dxc+df3dxc+df4dxc+df5dxc+df6dxc ) * 0.5;

338:     }
339:   }


342:   /* Compute triangular areas along the border of the domain. */
343:   if (xs==0){ /* left side */
344:     for (j=ys; j<ys+ym; j++){
345:       d3=(left[j-ys+1] - left[j-ys+2])*rhy;
346:       d2=(left[j-ys+1] - x[(j-gys)*gxm])*rhx;
347:       ft = ft+PetscSqrtScalar( 1.0 + d3*d3 + d2*d2);
348:     }
349:   }
350:   if (ys==0){ /* bottom side */
351:     for (i=xs; i<xs+xm; i++){
352:       d2=(bottom[i+1-xs]-bottom[i-xs+2])*rhx;
353:       d3=(bottom[i-xs+1]-x[i-gxs])*rhy;
354:       ft = ft+PetscSqrtScalar( 1.0 + d3*d3 + d2*d2);
355:     }
356:   }

358:   if (xs+xm==mx){ /* right side */
359:     for (j=ys; j< ys+ym; j++){
360:       d1=(x[(j+1-gys)*gxm-1]-right[j-ys+1])*rhx;
361:       d4=(right[j-ys]-right[j-ys+1])*rhy;
362:       ft = ft+PetscSqrtScalar( 1.0 + d1*d1 + d4*d4);
363:     }
364:   }
365:   if (ys+ym==my){ /* top side */
366:     for (i=xs; i<xs+xm; i++){
367:       d1=(x[(gym-1)*gxm + i-gxs] - top[i-xs+1])*rhy;
368:       d4=(top[i-xs+1] - top[i-xs])*rhx;
369:       ft = ft+PetscSqrtScalar( 1.0 + d1*d1 + d4*d4);
370:     }
371:   }

373:   if (ys==0 && xs==0){
374:     d1=(left[0]-left[1])*rhy;
375:     d2=(bottom[0]-bottom[1])*rhx;
376:     ft +=PetscSqrtScalar( 1.0 + d1*d1 + d2*d2);
377:   }
378:   if (ys+ym == my && xs+xm == mx){
379:     d1=(right[ym+1] - right[ym])*rhy;
380:     d2=(top[xm+1] - top[xm])*rhx;
381:     ft +=PetscSqrtScalar( 1.0 + d1*d1 + d2*d2);
382:   }

384:   ft=ft*area;
385:   MPI_Allreduce(&ft,fcn,1,MPIU_REAL,MPIU_SUM,MPI_COMM_WORLD);


388:   /* Restore vectors */
389:   VecRestoreArray(localX,&x);
390:   VecRestoreArray(localG,&g);
391:   VecRestoreArray(user->Left,&left);
392:   VecRestoreArray(user->Top,&top);
393:   VecRestoreArray(user->Bottom,&bottom);
394:   VecRestoreArray(user->Right,&right);

396:   /* Scatter values to global vector */
397:   DMLocalToGlobalBegin(user->dm,localG,INSERT_VALUES,G);
398:   DMLocalToGlobalEnd(user->dm,localG,INSERT_VALUES,G);

400:   PetscLogFlops(70*xm*ym);

402:   return 0;
403: }

405: /* ------------------------------------------------------------------- */
408: /*
409:    FormHessian - Evaluates Hessian matrix.

411:    Input Parameters:
412: .  tao  - the Tao context
413: .  x    - input vector
414: .  ptr  - optional user-defined context, as set by TaoSetHessianRoutine()

416:    Output Parameters:
417: .  A    - Hessian matrix
418: .  B    - optionally different preconditioning matrix

420:    Notes:
421:    Due to mesh point reordering with DMs, we must always work
422:    with the local mesh points, and then transform them to the new
423:    global numbering with the local-to-global mapping.  We cannot work
424:    directly with the global numbers for the original uniprocessor mesh!

426:    Two methods are available for imposing this transformation
427:    when setting matrix entries:
428:      (A) MatSetValuesLocal(), using the local ordering (including
429:          ghost points!)
430:          - Do the following two steps once, before calling TaoSolve()
431:            - Use DMGetISLocalToGlobalMapping() to extract the
432:              local-to-global map from the DM
433:            - Associate this map with the matrix by calling
434:              MatSetLocalToGlobalMapping()
435:          - Then set matrix entries using the local ordering
436:            by calling MatSetValuesLocal()
437:      (B) MatSetValues(), using the global ordering
438:          - Use DMGetGlobalIndices() to extract the local-to-global map
439:          - Then apply this map explicitly yourself
440:          - Set matrix entries using the global ordering by calling
441:            MatSetValues()
442:    Option (A) seems cleaner/easier in many cases, and is the procedure
443:    used in this example.
444: */
445: PetscErrorCode FormHessian(Tao tao,Vec X,Mat Hptr, Mat Hessian, void *ptr)
446: {
448:   AppCtx         *user = (AppCtx *) ptr;
449:   PetscInt       i,j,k,row;
450:   PetscInt       mx=user->mx, my=user->my;
451:   PetscInt       xs,xm,gxs,gxm,ys,ym,gys,gym,col[7];
452:   PetscReal      hx=1.0/(mx+1), hy=1.0/(my+1), hydhx=hy/hx, hxdhy=hx/hy;
453:   PetscReal      rhx=mx+1, rhy=my+1;
454:   PetscReal      f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb;
455:   PetscReal      hl,hr,ht,hb,hc,htl,hbr;
456:   PetscReal      *x,*left,*right,*bottom,*top;
457:   PetscReal      v[7];
458:   Vec            localX = user->localX;
459:   PetscBool      assembled;


462:   /* Set various matrix options */
463:   MatSetOption(Hessian,MAT_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE);

465:   /* Initialize matrix entries to zero */
466:   MatAssembled(Hessian,&assembled);
467:   if (assembled){MatZeroEntries(Hessian);}

469:   /* Get local mesh boundaries */
470:   DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL);
471:   DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL);

473:   /* Scatter ghost points to local vector */
474:   DMGlobalToLocalBegin(user->dm,X,INSERT_VALUES,localX);
475:   DMGlobalToLocalEnd(user->dm,X,INSERT_VALUES,localX);

477:   /* Get pointers to vector data */
478:   VecGetArray(localX,&x);
479:   VecGetArray(user->Top,&top);
480:   VecGetArray(user->Bottom,&bottom);
481:   VecGetArray(user->Left,&left);
482:   VecGetArray(user->Right,&right);

484:   /* Compute Hessian over the locally owned part of the mesh */

486:   for (i=xs; i< xs+xm; i++){

488:     for (j=ys; j<ys+ym; j++){

490:       row=(j-gys)*gxm + (i-gxs);

492:       xc = x[row];
493:       xlt=xrb=xl=xr=xb=xt=xc;

495:       /* Left side */
496:       if (i==gxs){
497:         xl= left[j-ys+1];
498:         xlt = left[j-ys+2];
499:       } else {
500:         xl = x[row-1];
501:       }

503:       if (j==gys){
504:         xb=bottom[i-xs+1];
505:         xrb = bottom[i-xs+2];
506:       } else {
507:         xb = x[row-gxm];
508:       }

510:       if (i+1 == gxs+gxm){
511:         xr=right[j-ys+1];
512:         xrb = right[j-ys];
513:       } else {
514:         xr = x[row+1];
515:       }

517:       if (j+1==gys+gym){
518:         xt=top[i-xs+1];
519:         xlt = top[i-xs];
520:       }else {
521:         xt = x[row+gxm];
522:       }

524:       if (i>gxs && j+1<gys+gym){
525:         xlt = x[row-1+gxm];
526:       }
527:       if (j>gys && i+1<gxs+gxm){
528:         xrb = x[row+1-gxm];
529:       }


532:       d1 = (xc-xl)*rhx;
533:       d2 = (xc-xr)*rhx;
534:       d3 = (xc-xt)*rhy;
535:       d4 = (xc-xb)*rhy;
536:       d5 = (xrb-xr)*rhy;
537:       d6 = (xrb-xb)*rhx;
538:       d7 = (xlt-xl)*rhy;
539:       d8 = (xlt-xt)*rhx;

541:       f1 = PetscSqrtScalar( 1.0 + d1*d1 + d7*d7);
542:       f2 = PetscSqrtScalar( 1.0 + d1*d1 + d4*d4);
543:       f3 = PetscSqrtScalar( 1.0 + d3*d3 + d8*d8);
544:       f4 = PetscSqrtScalar( 1.0 + d3*d3 + d2*d2);
545:       f5 = PetscSqrtScalar( 1.0 + d2*d2 + d5*d5);
546:       f6 = PetscSqrtScalar( 1.0 + d4*d4 + d6*d6);


549:       hl = (-hydhx*(1.0+d7*d7)+d1*d7)/(f1*f1*f1)+
550:         (-hydhx*(1.0+d4*d4)+d1*d4)/(f2*f2*f2);
551:       hr = (-hydhx*(1.0+d5*d5)+d2*d5)/(f5*f5*f5)+
552:         (-hydhx*(1.0+d3*d3)+d2*d3)/(f4*f4*f4);
553:       ht = (-hxdhy*(1.0+d8*d8)+d3*d8)/(f3*f3*f3)+
554:         (-hxdhy*(1.0+d2*d2)+d2*d3)/(f4*f4*f4);
555:       hb = (-hxdhy*(1.0+d6*d6)+d4*d6)/(f6*f6*f6)+
556:         (-hxdhy*(1.0+d1*d1)+d1*d4)/(f2*f2*f2);

558:       hbr = -d2*d5/(f5*f5*f5) - d4*d6/(f6*f6*f6);
559:       htl = -d1*d7/(f1*f1*f1) - d3*d8/(f3*f3*f3);

561:       hc = hydhx*(1.0+d7*d7)/(f1*f1*f1) + hxdhy*(1.0+d8*d8)/(f3*f3*f3) +
562:         hydhx*(1.0+d5*d5)/(f5*f5*f5) + hxdhy*(1.0+d6*d6)/(f6*f6*f6) +
563:         (hxdhy*(1.0+d1*d1)+hydhx*(1.0+d4*d4)-2*d1*d4)/(f2*f2*f2) +
564:         (hxdhy*(1.0+d2*d2)+hydhx*(1.0+d3*d3)-2*d2*d3)/(f4*f4*f4);

566:       hl*=0.5; hr*=0.5; ht*=0.5; hb*=0.5; hbr*=0.5; htl*=0.5;  hc*=0.5;

568:       k=0;
569:       if (j>0){
570:         v[k]=hb; col[k]=row - gxm; k++;
571:       }

573:       if (j>0 && i < mx -1){
574:         v[k]=hbr; col[k]=row - gxm+1; k++;
575:       }

577:       if (i>0){
578:         v[k]= hl; col[k]=row - 1; k++;
579:       }

581:       v[k]= hc; col[k]=row; k++;

583:       if (i < mx-1 ){
584:         v[k]= hr; col[k]=row+1; k++;
585:       }

587:       if (i>0 && j < my-1 ){
588:         v[k]= htl; col[k] = row+gxm-1; k++;
589:       }

591:       if (j < my-1 ){
592:         v[k]= ht; col[k] = row+gxm; k++;
593:       }

595:       /*
596:          Set matrix values using local numbering, which was defined
597:          earlier, in the main routine.
598:       */
599:       MatSetValuesLocal(Hessian,1,&row,k,col,v,INSERT_VALUES);

601:     }
602:   }

604:   /* Restore vectors */
605:   VecRestoreArray(localX,&x);
606:   VecRestoreArray(user->Left,&left);
607:   VecRestoreArray(user->Top,&top);
608:   VecRestoreArray(user->Bottom,&bottom);
609:   VecRestoreArray(user->Right,&right);

611:   /* Assemble the matrix */
612:   MatAssemblyBegin(Hessian,MAT_FINAL_ASSEMBLY);
613:   MatAssemblyEnd(Hessian,MAT_FINAL_ASSEMBLY);

615:   PetscLogFlops(199*xm*ym);
616:   return 0;
617: }

619: /* ------------------------------------------------------------------- */
622: /*
623:    MSA_BoundaryConditions -  Calculates the boundary conditions for
624:    the region.

626:    Input Parameter:
627: .  user - user-defined application context

629:    Output Parameter:
630: .  user - user-defined application context
631: */
632: static PetscErrorCode MSA_BoundaryConditions(AppCtx * user)
633: {
634:   int        ierr;
635:   PetscInt   i,j,k,maxits=5,limit=0;
636:   PetscInt   xs,ys,xm,ym,gxs,gys,gxm,gym;
637:   PetscInt   mx=user->mx,my=user->my;
638:   PetscInt   bsize=0, lsize=0, tsize=0, rsize=0;
639:   PetscReal  one=1.0, two=2.0, three=3.0, scl=1.0, tol=1e-10;
640:   PetscReal  fnorm,det,hx,hy,xt=0,yt=0;
641:   PetscReal  u1,u2,nf1,nf2,njac11,njac12,njac21,njac22;
642:   PetscReal  b=-0.5, t=0.5, l=-0.5, r=0.5;
643:   PetscReal  *boundary;
644:   PetscBool  flg;
645:   Vec        Bottom,Top,Right,Left;

647:   /* Get local mesh boundaries */
648:   DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL);
649:   DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL);


652:   bsize=xm+2;
653:   lsize=ym+2;
654:   rsize=ym+2;
655:   tsize=xm+2;

657:   VecCreateMPI(MPI_COMM_WORLD,bsize,PETSC_DECIDE,&Bottom);
658:   VecCreateMPI(MPI_COMM_WORLD,tsize,PETSC_DECIDE,&Top);
659:   VecCreateMPI(MPI_COMM_WORLD,lsize,PETSC_DECIDE,&Left);
660:   VecCreateMPI(MPI_COMM_WORLD,rsize,PETSC_DECIDE,&Right);

662:   user->Top=Top;
663:   user->Left=Left;
664:   user->Bottom=Bottom;
665:   user->Right=Right;

667:   hx= (r-l)/(mx+1); hy=(t-b)/(my+1);

669:   for (j=0; j<4; j++){
670:     if (j==0){
671:       yt=b;
672:       xt=l+hx*xs;
673:       limit=bsize;
674:       VecGetArray(Bottom,&boundary);
675:     } else if (j==1){
676:       yt=t;
677:       xt=l+hx*xs;
678:       limit=tsize;
679:       VecGetArray(Top,&boundary);
680:     } else if (j==2){
681:       yt=b+hy*ys;
682:       xt=l;
683:       limit=lsize;
684:       VecGetArray(Left,&boundary);
685:     } else if (j==3){
686:       yt=b+hy*ys;
687:       xt=r;
688:       limit=rsize;
689:       VecGetArray(Right,&boundary);
690:     }

692:     for (i=0; i<limit; i++){
693:       u1=xt;
694:       u2=-yt;
695:       for (k=0; k<maxits; k++){
696:         nf1=u1 + u1*u2*u2 - u1*u1*u1/three-xt;
697:         nf2=-u2 - u1*u1*u2 + u2*u2*u2/three-yt;
698:         fnorm=PetscSqrtScalar(nf1*nf1+nf2*nf2);
699:         if (fnorm <= tol) break;
700:         njac11=one+u2*u2-u1*u1;
701:         njac12=two*u1*u2;
702:         njac21=-two*u1*u2;
703:         njac22=-one - u1*u1 + u2*u2;
704:         det = njac11*njac22-njac21*njac12;
705:         u1 = u1-(njac22*nf1-njac12*nf2)/det;
706:         u2 = u2-(njac11*nf2-njac21*nf1)/det;
707:       }

709:       boundary[i]=u1*u1-u2*u2;
710:       if (j==0 || j==1) {
711:         xt=xt+hx;
712:       } else if (j==2 || j==3){
713:         yt=yt+hy;
714:       }
715:     }
716:     if (j==0){
717:       VecRestoreArray(Bottom,&boundary);
718:     } else if (j==1){
719:       VecRestoreArray(Top,&boundary);
720:     } else if (j==2){
721:       VecRestoreArray(Left,&boundary);
722:     } else if (j==3){
723:       VecRestoreArray(Right,&boundary);
724:     }
725:   }

727:   /* Scale the boundary if desired */

729:   PetscOptionsGetReal(NULL,NULL,"-bottom",&scl,&flg);
730:   if (flg){
731:     VecScale(Bottom, scl);
732:   }
733:   PetscOptionsGetReal(NULL,NULL,"-top",&scl,&flg);
734:   if (flg){
735:     VecScale(Top, scl);
736:   }
737:   PetscOptionsGetReal(NULL,NULL,"-right",&scl,&flg);
738:   if (flg){
739:     VecScale(Right, scl);
740:   }

742:   PetscOptionsGetReal(NULL,NULL,"-left",&scl,&flg);
743:   if (flg){
744:     VecScale(Left, scl);
745:   }
746:   return 0;
747: }


750: /* ------------------------------------------------------------------- */
753: /*
754:    MSA_Plate -  Calculates an obstacle for surface to stretch over.

756:    Input Parameter:
757: .  user - user-defined application context

759:    Output Parameter:
760: .  user - user-defined application context
761: */
762: static PetscErrorCode MSA_Plate(Vec XL,Vec XU,void *ctx){

764:   AppCtx         *user=(AppCtx *)ctx;
766:   PetscInt       i,j,row;
767:   PetscInt       xs,ys,xm,ym;
768:   PetscInt       mx=user->mx, my=user->my, bmy, bmx;
769:   PetscReal      t1,t2,t3;
770:   PetscReal      *xl, lb=PETSC_NINFINITY, ub=PETSC_INFINITY;
771:   PetscBool      cylinder;

773:   user->bmy = PetscMax(0,user->bmy);user->bmy = PetscMin(my,user->bmy);
774:   user->bmx = PetscMax(0,user->bmx);user->bmx = PetscMin(mx,user->bmx);
775:   bmy=user->bmy, bmx=user->bmx;

777:   DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL);

779:   VecSet(XL, lb);
780:   VecSet(XU, ub);

782:   VecGetArray(XL,&xl);

784:   PetscOptionsHasName(NULL,NULL,"-cylinder",&cylinder);
785:   /* Compute the optional lower box */
786:   if (cylinder){
787:     for (i=xs; i< xs+xm; i++){
788:       for (j=ys; j<ys+ym; j++){
789:         row=(j-ys)*xm + (i-xs);
790:         t1=(2.0*i-mx)*bmy;
791:         t2=(2.0*j-my)*bmx;
792:         t3=bmx*bmx*bmy*bmy;
793:         if ( t1*t1 + t2*t2 <= t3 ){
794:           xl[row] = user->bheight;
795:         }
796:       }
797:     }
798:   } else {
799:     /* Compute the optional lower box */
800:     for (i=xs; i< xs+xm; i++){
801:       for (j=ys; j<ys+ym; j++){
802:         row=(j-ys)*xm + (i-xs);
803:         if (i>=(mx-bmx)/2 && i<mx-(mx-bmx)/2 &&
804:             j>=(my-bmy)/2 && j<my-(my-bmy)/2 ){
805:           xl[row] = user->bheight;
806:         }
807:       }
808:     }
809:   }
810:     VecRestoreArray(XL,&xl);

812:   return 0;
813: }


816: /* ------------------------------------------------------------------- */
819: /*
820:    MSA_InitialPoint - Calculates the initial guess in one of three ways.

822:    Input Parameters:
823: .  user - user-defined application context
824: .  X - vector for initial guess

826:    Output Parameters:
827: .  X - newly computed initial guess
828: */
829: static PetscErrorCode MSA_InitialPoint(AppCtx * user, Vec X)
830: {
832:   PetscInt       start=-1,i,j;
833:   PetscReal      zero=0.0;
834:   PetscBool      flg;

836:   PetscOptionsGetInt(NULL,NULL,"-start",&start,&flg);
837:   if (flg && start==0){ /* The zero vector is reasonable */
838:     VecSet(X, zero);
839:   } else if (flg && start>0){ /* Try a random start between -0.5 and 0.5 */
840:     PetscRandom rctx;  PetscReal np5=-0.5;

842:     PetscRandomCreate(MPI_COMM_WORLD,&rctx);
843:     for (i=0; i<start; i++){
844:       VecSetRandom(X, rctx);
845:     }
846:     PetscRandomDestroy(&rctx);
847:     VecShift(X, np5);

849:   } else { /* Take an average of the boundary conditions */

851:     PetscInt row,xs,xm,gxs,gxm,ys,ym,gys,gym;
852:     PetscInt mx=user->mx,my=user->my;
853:     PetscReal *x,*left,*right,*bottom,*top;
854:     Vec    localX = user->localX;

856:     /* Get local mesh boundaries */
857:     DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL);
858:     DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL);

860:     /* Get pointers to vector data */
861:     VecGetArray(user->Top,&top);
862:     VecGetArray(user->Bottom,&bottom);
863:     VecGetArray(user->Left,&left);
864:     VecGetArray(user->Right,&right);

866:     VecGetArray(localX,&x);
867:     /* Perform local computations */
868:     for (j=ys; j<ys+ym; j++){
869:       for (i=xs; i< xs+xm; i++){
870:         row=(j-gys)*gxm + (i-gxs);
871:         x[row] = ( (j+1)*bottom[i-xs+1]/my + (my-j+1)*top[i-xs+1]/(my+2)+
872:                    (i+1)*left[j-ys+1]/mx + (mx-i+1)*right[j-ys+1]/(mx+2))/2.0;
873:       }
874:     }

876:     /* Restore vectors */
877:     VecRestoreArray(localX,&x);

879:     VecRestoreArray(user->Left,&left);
880:     VecRestoreArray(user->Top,&top);
881:     VecRestoreArray(user->Bottom,&bottom);
882:     VecRestoreArray(user->Right,&right);

884:     /* Scatter values into global vector */
885:     DMLocalToGlobalBegin(user->dm,localX,INSERT_VALUES,X);
886:     DMLocalToGlobalEnd(user->dm,localX,INSERT_VALUES,X);

888:   }
889:   return 0;
890: }

892: /* For testing matrix free submatrices */
895: PetscErrorCode MatrixFreeHessian(Tao tao, Vec x, Mat H, Mat Hpre, void *ptr)
896: {
898:   AppCtx         *user = (AppCtx*)ptr;
900:   FormHessian(tao,x,user->H,user->H,ptr);
901:   return(0);
902: }
905: PetscErrorCode MyMatMult(Mat H_shell, Vec X, Vec Y)
906: {
908:   void           *ptr;
909:   AppCtx         *user;
911:   MatShellGetContext(H_shell,&ptr);
912:   user = (AppCtx*)ptr;
913:   MatMult(user->H,X,Y);
914:   return(0);
915: }