Actual source code: plate2.c

petsc-master 2016-12-02
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 );if (ierr) return ierr;

 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,DMDA_STENCIL_BOX,user.mx,user.my,Nx,Ny,1,1,NULL,NULL,&user.dm);
108:   DMSetFromOptions(user.dm);
109:   DMSetUp(user.dm);
110:   /*
111:      Extract global and local vectors from DM; The local vectors are
112:      used solely as work space for the evaluation of the function,
113:      gradient, and Hessian.  Duplicate for remaining vectors that are
114:      the same types.
115:   */
116:   DMCreateGlobalVector(user.dm,&x); /* Solution */
117:   DMCreateLocalVector(user.dm,&user.localX);
118:   VecDuplicate(user.localX,&user.localV);

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

123:   /* The TAO code begins here */

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

133:   /* Set initial solution guess; */
134:   MSA_BoundaryConditions(&user);
135:   MSA_InitialPoint(&user,x);
136:   TaoSetInitialVector(tao,x);

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

141:   VecGetLocalSize(x,&m);
142:   MatCreateAIJ(MPI_COMM_WORLD,m,m,N,N,7,NULL,3,NULL,&(user.H));
143:   MatSetOption(user.H,MAT_SYMMETRIC,PETSC_TRUE);

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

157:   /* Set Variable bounds */
158:   MSA_Plate(xl,xu,(void*)&user);
159:   TaoSetVariableBounds(tao,xl,xu);

161:   /* Check for any tao command line options */
162:   TaoSetFromOptions(tao);

164:   /* SOLVE THE APPLICATION */
165:   TaoSolve(tao);

167:   TaoView(tao,PETSC_VIEWER_STDOUT_WORLD);

169:   /* Free TAO data structures */
170:   TaoDestroy(&tao);

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

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

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

200:     Output Parameters:
201: .   fcn     - the function value
202: .   G      - vector containing the newly evaluated gradient

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

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

234:   /* Scatter ghost points to local vector */
235:   DMGlobalToLocalBegin(user->dm,X,INSERT_VALUES,localX);
236:   DMGlobalToLocalEnd(user->dm,X,INSERT_VALUES,localX);

238:   /* Initialize vector to zero */
239:   VecSet(localG, zero);

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

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

254:       xc = x[row];
255:       xlt=xrb=xl=xr=xb=xt=xc;

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

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

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

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

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

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

301:       df1dxc = d1*hydhx;
302:       df2dxc = ( d1*hydhx + d4*hxdhy );
303:       df3dxc = d3*hxdhy;
304:       df4dxc = ( d2*hydhx + d3*hxdhy );
305:       df5dxc = d2*hydhx;
306:       df6dxc = d4*hxdhy;

308:       d1 *= rhx;
309:       d2 *= rhx;
310:       d3 *= rhy;
311:       d4 *= rhy;
312:       d5 *= rhy;
313:       d6 *= rhx;
314:       d7 *= rhy;
315:       d8 *= rhx;

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

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

326:       df1dxc /= f1;
327:       df2dxc /= f2;
328:       df3dxc /= f3;
329:       df4dxc /= f4;
330:       df5dxc /= f5;
331:       df6dxc /= f6;

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

335:     }
336:   }


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

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

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

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


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

393:   /* Scatter values to global vector */
394:   DMLocalToGlobalBegin(user->dm,localG,INSERT_VALUES,G);
395:   DMLocalToGlobalEnd(user->dm,localG,INSERT_VALUES,G);

397:   PetscLogFlops(70*xm*ym);

399:   return 0;
400: }

402: /* ------------------------------------------------------------------- */
405: /*
406:    FormHessian - Evaluates Hessian matrix.

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

413:    Output Parameters:
414: .  A    - Hessian matrix
415: .  B    - optionally different preconditioning matrix

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

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


459:   /* Set various matrix options */
460:   MatSetOption(Hessian,MAT_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE);

462:   /* Initialize matrix entries to zero */
463:   MatAssembled(Hessian,&assembled);
464:   if (assembled){MatZeroEntries(Hessian);}

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

470:   /* Scatter ghost points to local vector */
471:   DMGlobalToLocalBegin(user->dm,X,INSERT_VALUES,localX);
472:   DMGlobalToLocalEnd(user->dm,X,INSERT_VALUES,localX);

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

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

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

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

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

489:       xc = x[row];
490:       xlt=xrb=xl=xr=xb=xt=xc;

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

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

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

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

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


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

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


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

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

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

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

565:       k=0;
566:       if (j>0){
567:         v[k]=hb; col[k]=row - gxm; k++;
568:       }

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

574:       if (i>0){
575:         v[k]= hl; col[k]=row - 1; k++;
576:       }

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

580:       if (i < mx-1 ){
581:         v[k]= hr; col[k]=row+1; k++;
582:       }

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

588:       if (j < my-1 ){
589:         v[k]= ht; col[k] = row+gxm; k++;
590:       }

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

598:     }
599:   }

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

608:   /* Assemble the matrix */
609:   MatAssemblyBegin(Hessian,MAT_FINAL_ASSEMBLY);
610:   MatAssemblyEnd(Hessian,MAT_FINAL_ASSEMBLY);

612:   PetscLogFlops(199*xm*ym);
613:   return 0;
614: }

616: /* ------------------------------------------------------------------- */
619: /*
620:    MSA_BoundaryConditions -  Calculates the boundary conditions for
621:    the region.

623:    Input Parameter:
624: .  user - user-defined application context

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

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


649:   bsize=xm+2;
650:   lsize=ym+2;
651:   rsize=ym+2;
652:   tsize=xm+2;

654:   VecCreateMPI(MPI_COMM_WORLD,bsize,PETSC_DECIDE,&Bottom);
655:   VecCreateMPI(MPI_COMM_WORLD,tsize,PETSC_DECIDE,&Top);
656:   VecCreateMPI(MPI_COMM_WORLD,lsize,PETSC_DECIDE,&Left);
657:   VecCreateMPI(MPI_COMM_WORLD,rsize,PETSC_DECIDE,&Right);

659:   user->Top=Top;
660:   user->Left=Left;
661:   user->Bottom=Bottom;
662:   user->Right=Right;

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

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

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

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

724:   /* Scale the boundary if desired */

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

739:   PetscOptionsGetReal(NULL,NULL,"-left",&scl,&flg);
740:   if (flg){
741:     VecScale(Left, scl);
742:   }
743:   return 0;
744: }


747: /* ------------------------------------------------------------------- */
750: /*
751:    MSA_Plate -  Calculates an obstacle for surface to stretch over.

753:    Input Parameter:
754: .  user - user-defined application context

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

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

770:   user->bmy = PetscMax(0,user->bmy);user->bmy = PetscMin(my,user->bmy);
771:   user->bmx = PetscMax(0,user->bmx);user->bmx = PetscMin(mx,user->bmx);
772:   bmy=user->bmy; bmx=user->bmx;

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

776:   VecSet(XL, lb);
777:   VecSet(XU, ub);

779:   VecGetArray(XL,&xl);

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

809:   return 0;
810: }


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

819:    Input Parameters:
820: .  user - user-defined application context
821: .  X - vector for initial guess

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

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

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

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

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

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

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

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

873:     /* Restore vectors */
874:     VecRestoreArray(localX,&x);

876:     VecRestoreArray(user->Left,&left);
877:     VecRestoreArray(user->Top,&top);
878:     VecRestoreArray(user->Bottom,&bottom);
879:     VecRestoreArray(user->Right,&right);

881:     /* Scatter values into global vector */
882:     DMLocalToGlobalBegin(user->dm,localX,INSERT_VALUES,X);
883:     DMLocalToGlobalEnd(user->dm,localX,INSERT_VALUES,X);

885:   }
886:   return 0;
887: }

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