xref: /petsc/src/tao/bound/tutorials/plate2.c (revision 5cab5458055e6544d97095d04e76587ba3d30732)
1 #include <petscdmda.h>
2 #include <petsctao.h>
3 
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";
19 
20 /*
21    User-defined application context - contains data needed by the
22    application-provided call-back routines, FormFunctionGradient(),
23    FormHessian().
24 */
25 typedef struct {
26   /* problem parameters */
27   PetscReal      bheight;                  /* Height of plate under the surface */
28   PetscInt       mx, my;                   /* discretization in x, y directions */
29   PetscInt       bmx,bmy;                  /* Size of plate under the surface */
30   Vec            Bottom, Top, Left, Right; /* boundary values */
31 
32   /* Working space */
33   Vec         localX, localV;           /* ghosted local vector */
34   DM          dm;                       /* distributed array data structure */
35   Mat         H;
36 } AppCtx;
37 
38 /* -------- User-defined Routines --------- */
39 
40 static PetscErrorCode MSA_BoundaryConditions(AppCtx*);
41 static PetscErrorCode MSA_InitialPoint(AppCtx*,Vec);
42 static PetscErrorCode MSA_Plate(Vec,Vec,void*);
43 PetscErrorCode FormFunctionGradient(Tao,Vec,PetscReal*,Vec,void*);
44 PetscErrorCode FormHessian(Tao,Vec,Mat,Mat,void*);
45 
46 /* For testing matrix free submatrices */
47 PetscErrorCode MatrixFreeHessian(Tao,Vec,Mat, Mat,void*);
48 PetscErrorCode MyMatMult(Mat,Vec,Vec);
49 
50 int main(int argc, char **argv)
51 {
52   PetscInt               Nx, Ny;               /* number of processors in x- and y- directions */
53   PetscInt               m, N;                 /* number of local and global elements in vectors */
54   Vec                    x,xl,xu;               /* solution vector  and bounds*/
55   PetscBool              flg;                /* A return variable when checking for user options */
56   Tao                    tao;                  /* Tao solver context */
57   ISLocalToGlobalMapping isltog;   /* local-to-global mapping object */
58   Mat                    H_shell;                  /* to test matrix-free submatrices */
59   AppCtx                 user;                 /* user-defined work context */
60 
61   /* Initialize PETSc, TAO */
62   PetscFunctionBeginUser;
63   PetscCall(PetscInitialize(&argc, &argv,(char *)0,help));
64 
65   /* Specify default dimension of the problem */
66   user.mx = 10; user.my = 10; user.bheight=0.1;
67 
68   /* Check for any command line arguments that override defaults */
69   PetscCall(PetscOptionsGetInt(NULL,NULL,"-mx",&user.mx,&flg));
70   PetscCall(PetscOptionsGetInt(NULL,NULL,"-my",&user.my,&flg));
71   PetscCall(PetscOptionsGetReal(NULL,NULL,"-bheight",&user.bheight,&flg));
72 
73   user.bmx = user.mx/2; user.bmy = user.my/2;
74   PetscCall(PetscOptionsGetInt(NULL,NULL,"-bmx",&user.bmx,&flg));
75   PetscCall(PetscOptionsGetInt(NULL,NULL,"-bmy",&user.bmy,&flg));
76 
77   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n---- Minimum Surface Area With Plate Problem -----\n"));
78   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"mx:%" PetscInt_FMT ", my:%" PetscInt_FMT ", bmx:%" PetscInt_FMT ", bmy:%" PetscInt_FMT ", height:%g\n",user.mx,user.my,user.bmx,user.bmy,(double)user.bheight));
79 
80   /* Calculate any derived values from parameters */
81   N    = user.mx*user.my;
82 
83   /* Let Petsc determine the dimensions of the local vectors */
84   Nx = PETSC_DECIDE; Ny = PETSC_DECIDE;
85 
86   /*
87      A two dimensional distributed array will help define this problem,
88      which derives from an elliptic PDE on two dimensional domain.  From
89      the distributed array, Create the vectors.
90   */
91   PetscCall(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));
92   PetscCall(DMSetFromOptions(user.dm));
93   PetscCall(DMSetUp(user.dm));
94   /*
95      Extract global and local vectors from DM; The local vectors are
96      used solely as work space for the evaluation of the function,
97      gradient, and Hessian.  Duplicate for remaining vectors that are
98      the same types.
99   */
100   PetscCall(DMCreateGlobalVector(user.dm,&x)); /* Solution */
101   PetscCall(DMCreateLocalVector(user.dm,&user.localX));
102   PetscCall(VecDuplicate(user.localX,&user.localV));
103 
104   PetscCall(VecDuplicate(x,&xl));
105   PetscCall(VecDuplicate(x,&xu));
106 
107   /* The TAO code begins here */
108 
109   /*
110      Create TAO solver and set desired solution method
111      The method must either be TAOTRON or TAOBLMVM
112      If TAOBLMVM is used, then hessian function is not called.
113   */
114   PetscCall(TaoCreate(PETSC_COMM_WORLD,&tao));
115   PetscCall(TaoSetType(tao,TAOBLMVM));
116 
117   /* Set initial solution guess; */
118   PetscCall(MSA_BoundaryConditions(&user));
119   PetscCall(MSA_InitialPoint(&user,x));
120   PetscCall(TaoSetSolution(tao,x));
121 
122   /* Set routines for function, gradient and hessian evaluation */
123   PetscCall(TaoSetObjectiveAndGradient(tao,NULL,FormFunctionGradient,(void*) &user));
124 
125   PetscCall(VecGetLocalSize(x,&m));
126   PetscCall(MatCreateAIJ(MPI_COMM_WORLD,m,m,N,N,7,NULL,3,NULL,&(user.H)));
127   PetscCall(MatSetOption(user.H,MAT_SYMMETRIC,PETSC_TRUE));
128 
129   PetscCall(DMGetLocalToGlobalMapping(user.dm,&isltog));
130   PetscCall(MatSetLocalToGlobalMapping(user.H,isltog,isltog));
131   PetscCall(PetscOptionsHasName(NULL,NULL,"-matrixfree",&flg));
132   if (flg) {
133       PetscCall(MatCreateShell(PETSC_COMM_WORLD,m,m,N,N,(void*)&user,&H_shell));
134       PetscCall(MatShellSetOperation(H_shell,MATOP_MULT,(void(*)(void))MyMatMult));
135       PetscCall(MatSetOption(H_shell,MAT_SYMMETRIC,PETSC_TRUE));
136       PetscCall(TaoSetHessian(tao,H_shell,H_shell,MatrixFreeHessian,(void*)&user));
137   } else {
138       PetscCall(TaoSetHessian(tao,user.H,user.H,FormHessian,(void*)&user));
139   }
140 
141   /* Set Variable bounds */
142   PetscCall(MSA_Plate(xl,xu,(void*)&user));
143   PetscCall(TaoSetVariableBounds(tao,xl,xu));
144 
145   /* Check for any tao command line options */
146   PetscCall(TaoSetFromOptions(tao));
147 
148   /* SOLVE THE APPLICATION */
149   PetscCall(TaoSolve(tao));
150 
151   PetscCall(TaoView(tao,PETSC_VIEWER_STDOUT_WORLD));
152 
153   /* Free TAO data structures */
154   PetscCall(TaoDestroy(&tao));
155 
156   /* Free PETSc data structures */
157   PetscCall(VecDestroy(&x));
158   PetscCall(VecDestroy(&xl));
159   PetscCall(VecDestroy(&xu));
160   PetscCall(MatDestroy(&user.H));
161   PetscCall(VecDestroy(&user.localX));
162   PetscCall(VecDestroy(&user.localV));
163   PetscCall(VecDestroy(&user.Bottom));
164   PetscCall(VecDestroy(&user.Top));
165   PetscCall(VecDestroy(&user.Left));
166   PetscCall(VecDestroy(&user.Right));
167   PetscCall(DMDestroy(&user.dm));
168   if (flg) {
169     PetscCall(MatDestroy(&H_shell));
170   }
171   PetscCall(PetscFinalize());
172   return 0;
173 }
174 
175 /*  FormFunctionGradient - Evaluates f(x) and gradient g(x).
176 
177     Input Parameters:
178 .   tao     - the Tao context
179 .   X      - input vector
180 .   userCtx - optional user-defined context, as set by TaoSetObjectiveAndGradient()
181 
182     Output Parameters:
183 .   fcn     - the function value
184 .   G      - vector containing the newly evaluated gradient
185 
186    Notes:
187    In this case, we discretize the domain and Create triangles. The
188    surface of each triangle is planar, whose surface area can be easily
189    computed.  The total surface area is found by sweeping through the grid
190    and computing the surface area of the two triangles that have their
191    right angle at the grid point.  The diagonal line segments on the
192    grid that define the triangles run from top left to lower right.
193    The numbering of points starts at the lower left and runs left to
194    right, then bottom to top.
195 */
196 PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *fcn, Vec G,void *userCtx)
197 {
198   AppCtx         *user = (AppCtx *) userCtx;
199   PetscInt       i,j,row;
200   PetscInt       mx=user->mx, my=user->my;
201   PetscInt       xs,xm,gxs,gxm,ys,ym,gys,gym;
202   PetscReal      ft=0;
203   PetscReal      zero=0.0;
204   PetscReal      hx=1.0/(mx+1),hy=1.0/(my+1), hydhx=hy/hx, hxdhy=hx/hy, area=0.5*hx*hy;
205   PetscReal      rhx=mx+1, rhy=my+1;
206   PetscReal      f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb;
207   PetscReal      df1dxc,df2dxc,df3dxc,df4dxc,df5dxc,df6dxc;
208   PetscReal      *g, *x,*left,*right,*bottom,*top;
209   Vec            localX = user->localX, localG = user->localV;
210 
211   /* Get local mesh boundaries */
212   PetscCall(DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL));
213   PetscCall(DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL));
214 
215   /* Scatter ghost points to local vector */
216   PetscCall(DMGlobalToLocalBegin(user->dm,X,INSERT_VALUES,localX));
217   PetscCall(DMGlobalToLocalEnd(user->dm,X,INSERT_VALUES,localX));
218 
219   /* Initialize vector to zero */
220   PetscCall(VecSet(localG, zero));
221 
222   /* Get pointers to vector data */
223   PetscCall(VecGetArray(localX,&x));
224   PetscCall(VecGetArray(localG,&g));
225   PetscCall(VecGetArray(user->Top,&top));
226   PetscCall(VecGetArray(user->Bottom,&bottom));
227   PetscCall(VecGetArray(user->Left,&left));
228   PetscCall(VecGetArray(user->Right,&right));
229 
230   /* Compute function over the locally owned part of the mesh */
231   for (j=ys; j<ys+ym; j++) {
232     for (i=xs; i< xs+xm; i++) {
233       row=(j-gys)*gxm + (i-gxs);
234 
235       xc = x[row];
236       xlt=xrb=xl=xr=xb=xt=xc;
237 
238       if (i==0) { /* left side */
239         xl= left[j-ys+1];
240         xlt = left[j-ys+2];
241       } else {
242         xl = x[row-1];
243       }
244 
245       if (j==0) { /* bottom side */
246         xb=bottom[i-xs+1];
247         xrb = bottom[i-xs+2];
248       } else {
249         xb = x[row-gxm];
250       }
251 
252       if (i+1 == gxs+gxm) { /* right side */
253         xr=right[j-ys+1];
254         xrb = right[j-ys];
255       } else {
256         xr = x[row+1];
257       }
258 
259       if (j+1==gys+gym) { /* top side */
260         xt=top[i-xs+1];
261         xlt = top[i-xs];
262       }else {
263         xt = x[row+gxm];
264       }
265 
266       if (i>gxs && j+1<gys+gym) {
267         xlt = x[row-1+gxm];
268       }
269       if (j>gys && i+1<gxs+gxm) {
270         xrb = x[row+1-gxm];
271       }
272 
273       d1 = (xc-xl);
274       d2 = (xc-xr);
275       d3 = (xc-xt);
276       d4 = (xc-xb);
277       d5 = (xr-xrb);
278       d6 = (xrb-xb);
279       d7 = (xlt-xl);
280       d8 = (xt-xlt);
281 
282       df1dxc = d1*hydhx;
283       df2dxc = (d1*hydhx + d4*hxdhy);
284       df3dxc = d3*hxdhy;
285       df4dxc = (d2*hydhx + d3*hxdhy);
286       df5dxc = d2*hydhx;
287       df6dxc = d4*hxdhy;
288 
289       d1 *= rhx;
290       d2 *= rhx;
291       d3 *= rhy;
292       d4 *= rhy;
293       d5 *= rhy;
294       d6 *= rhx;
295       d7 *= rhy;
296       d8 *= rhx;
297 
298       f1 = PetscSqrtScalar(1.0 + d1*d1 + d7*d7);
299       f2 = PetscSqrtScalar(1.0 + d1*d1 + d4*d4);
300       f3 = PetscSqrtScalar(1.0 + d3*d3 + d8*d8);
301       f4 = PetscSqrtScalar(1.0 + d3*d3 + d2*d2);
302       f5 = PetscSqrtScalar(1.0 + d2*d2 + d5*d5);
303       f6 = PetscSqrtScalar(1.0 + d4*d4 + d6*d6);
304 
305       ft = ft + (f2 + f4);
306 
307       df1dxc /= f1;
308       df2dxc /= f2;
309       df3dxc /= f3;
310       df4dxc /= f4;
311       df5dxc /= f5;
312       df6dxc /= f6;
313 
314       g[row] = (df1dxc+df2dxc+df3dxc+df4dxc+df5dxc+df6dxc) * 0.5;
315 
316     }
317   }
318 
319   /* Compute triangular areas along the border of the domain. */
320   if (xs==0) { /* left side */
321     for (j=ys; j<ys+ym; j++) {
322       d3=(left[j-ys+1] - left[j-ys+2])*rhy;
323       d2=(left[j-ys+1] - x[(j-gys)*gxm])*rhx;
324       ft = ft+PetscSqrtScalar(1.0 + d3*d3 + d2*d2);
325     }
326   }
327   if (ys==0) { /* bottom side */
328     for (i=xs; i<xs+xm; i++) {
329       d2=(bottom[i+1-xs]-bottom[i-xs+2])*rhx;
330       d3=(bottom[i-xs+1]-x[i-gxs])*rhy;
331       ft = ft+PetscSqrtScalar(1.0 + d3*d3 + d2*d2);
332     }
333   }
334 
335   if (xs+xm==mx) { /* right side */
336     for (j=ys; j< ys+ym; j++) {
337       d1=(x[(j+1-gys)*gxm-1]-right[j-ys+1])*rhx;
338       d4=(right[j-ys]-right[j-ys+1])*rhy;
339       ft = ft+PetscSqrtScalar(1.0 + d1*d1 + d4*d4);
340     }
341   }
342   if (ys+ym==my) { /* top side */
343     for (i=xs; i<xs+xm; i++) {
344       d1=(x[(gym-1)*gxm + i-gxs] - top[i-xs+1])*rhy;
345       d4=(top[i-xs+1] - top[i-xs])*rhx;
346       ft = ft+PetscSqrtScalar(1.0 + d1*d1 + d4*d4);
347     }
348   }
349 
350   if (ys==0 && xs==0) {
351     d1=(left[0]-left[1])*rhy;
352     d2=(bottom[0]-bottom[1])*rhx;
353     ft +=PetscSqrtScalar(1.0 + d1*d1 + d2*d2);
354   }
355   if (ys+ym == my && xs+xm == mx) {
356     d1=(right[ym+1] - right[ym])*rhy;
357     d2=(top[xm+1] - top[xm])*rhx;
358     ft +=PetscSqrtScalar(1.0 + d1*d1 + d2*d2);
359   }
360 
361   ft=ft*area;
362   PetscCallMPI(MPI_Allreduce(&ft,fcn,1,MPIU_REAL,MPIU_SUM,MPI_COMM_WORLD));
363 
364   /* Restore vectors */
365   PetscCall(VecRestoreArray(localX,&x));
366   PetscCall(VecRestoreArray(localG,&g));
367   PetscCall(VecRestoreArray(user->Left,&left));
368   PetscCall(VecRestoreArray(user->Top,&top));
369   PetscCall(VecRestoreArray(user->Bottom,&bottom));
370   PetscCall(VecRestoreArray(user->Right,&right));
371 
372   /* Scatter values to global vector */
373   PetscCall(DMLocalToGlobalBegin(user->dm,localG,INSERT_VALUES,G));
374   PetscCall(DMLocalToGlobalEnd(user->dm,localG,INSERT_VALUES,G));
375 
376   PetscCall(PetscLogFlops(70.0*xm*ym));
377 
378   return 0;
379 }
380 
381 /* ------------------------------------------------------------------- */
382 /*
383    FormHessian - Evaluates Hessian matrix.
384 
385    Input Parameters:
386 .  tao  - the Tao context
387 .  x    - input vector
388 .  ptr  - optional user-defined context, as set by TaoSetHessian()
389 
390    Output Parameters:
391 .  A    - Hessian matrix
392 .  B    - optionally different preconditioning matrix
393 
394    Notes:
395    Due to mesh point reordering with DMs, we must always work
396    with the local mesh points, and then transform them to the new
397    global numbering with the local-to-global mapping.  We cannot work
398    directly with the global numbers for the original uniprocessor mesh!
399 
400    Two methods are available for imposing this transformation
401    when setting matrix entries:
402      (A) MatSetValuesLocal(), using the local ordering (including
403          ghost points!)
404          - Do the following two steps once, before calling TaoSolve()
405            - Use DMGetISLocalToGlobalMapping() to extract the
406              local-to-global map from the DM
407            - Associate this map with the matrix by calling
408              MatSetLocalToGlobalMapping()
409          - Then set matrix entries using the local ordering
410            by calling MatSetValuesLocal()
411      (B) MatSetValues(), using the global ordering
412          - Use DMGetGlobalIndices() to extract the local-to-global map
413          - Then apply this map explicitly yourself
414          - Set matrix entries using the global ordering by calling
415            MatSetValues()
416    Option (A) seems cleaner/easier in many cases, and is the procedure
417    used in this example.
418 */
419 PetscErrorCode FormHessian(Tao tao,Vec X,Mat Hptr, Mat Hessian, void *ptr)
420 {
421   AppCtx         *user = (AppCtx *) ptr;
422   PetscInt       i,j,k,row;
423   PetscInt       mx=user->mx, my=user->my;
424   PetscInt       xs,xm,gxs,gxm,ys,ym,gys,gym,col[7];
425   PetscReal      hx=1.0/(mx+1), hy=1.0/(my+1), hydhx=hy/hx, hxdhy=hx/hy;
426   PetscReal      rhx=mx+1, rhy=my+1;
427   PetscReal      f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb;
428   PetscReal      hl,hr,ht,hb,hc,htl,hbr;
429   PetscReal      *x,*left,*right,*bottom,*top;
430   PetscReal      v[7];
431   Vec            localX = user->localX;
432   PetscBool      assembled;
433 
434   /* Set various matrix options */
435   PetscCall(MatSetOption(Hessian,MAT_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE));
436 
437   /* Initialize matrix entries to zero */
438   PetscCall(MatAssembled(Hessian,&assembled));
439   if (assembled) PetscCall(MatZeroEntries(Hessian));
440 
441   /* Get local mesh boundaries */
442   PetscCall(DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL));
443   PetscCall(DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL));
444 
445   /* Scatter ghost points to local vector */
446   PetscCall(DMGlobalToLocalBegin(user->dm,X,INSERT_VALUES,localX));
447   PetscCall(DMGlobalToLocalEnd(user->dm,X,INSERT_VALUES,localX));
448 
449   /* Get pointers to vector data */
450   PetscCall(VecGetArray(localX,&x));
451   PetscCall(VecGetArray(user->Top,&top));
452   PetscCall(VecGetArray(user->Bottom,&bottom));
453   PetscCall(VecGetArray(user->Left,&left));
454   PetscCall(VecGetArray(user->Right,&right));
455 
456   /* Compute Hessian over the locally owned part of the mesh */
457 
458   for (i=xs; i< xs+xm; i++) {
459 
460     for (j=ys; j<ys+ym; j++) {
461 
462       row=(j-gys)*gxm + (i-gxs);
463 
464       xc = x[row];
465       xlt=xrb=xl=xr=xb=xt=xc;
466 
467       /* Left side */
468       if (i==gxs) {
469         xl= left[j-ys+1];
470         xlt = left[j-ys+2];
471       } else {
472         xl = x[row-1];
473       }
474 
475       if (j==gys) {
476         xb=bottom[i-xs+1];
477         xrb = bottom[i-xs+2];
478       } else {
479         xb = x[row-gxm];
480       }
481 
482       if (i+1 == gxs+gxm) {
483         xr=right[j-ys+1];
484         xrb = right[j-ys];
485       } else {
486         xr = x[row+1];
487       }
488 
489       if (j+1==gys+gym) {
490         xt=top[i-xs+1];
491         xlt = top[i-xs];
492       }else {
493         xt = x[row+gxm];
494       }
495 
496       if (i>gxs && j+1<gys+gym) {
497         xlt = x[row-1+gxm];
498       }
499       if (j>gys && i+1<gxs+gxm) {
500         xrb = x[row+1-gxm];
501       }
502 
503       d1 = (xc-xl)*rhx;
504       d2 = (xc-xr)*rhx;
505       d3 = (xc-xt)*rhy;
506       d4 = (xc-xb)*rhy;
507       d5 = (xrb-xr)*rhy;
508       d6 = (xrb-xb)*rhx;
509       d7 = (xlt-xl)*rhy;
510       d8 = (xlt-xt)*rhx;
511 
512       f1 = PetscSqrtScalar(1.0 + d1*d1 + d7*d7);
513       f2 = PetscSqrtScalar(1.0 + d1*d1 + d4*d4);
514       f3 = PetscSqrtScalar(1.0 + d3*d3 + d8*d8);
515       f4 = PetscSqrtScalar(1.0 + d3*d3 + d2*d2);
516       f5 = PetscSqrtScalar(1.0 + d2*d2 + d5*d5);
517       f6 = PetscSqrtScalar(1.0 + d4*d4 + d6*d6);
518 
519       hl = (-hydhx*(1.0+d7*d7)+d1*d7)/(f1*f1*f1)+
520         (-hydhx*(1.0+d4*d4)+d1*d4)/(f2*f2*f2);
521       hr = (-hydhx*(1.0+d5*d5)+d2*d5)/(f5*f5*f5)+
522         (-hydhx*(1.0+d3*d3)+d2*d3)/(f4*f4*f4);
523       ht = (-hxdhy*(1.0+d8*d8)+d3*d8)/(f3*f3*f3)+
524         (-hxdhy*(1.0+d2*d2)+d2*d3)/(f4*f4*f4);
525       hb = (-hxdhy*(1.0+d6*d6)+d4*d6)/(f6*f6*f6)+
526         (-hxdhy*(1.0+d1*d1)+d1*d4)/(f2*f2*f2);
527 
528       hbr = -d2*d5/(f5*f5*f5) - d4*d6/(f6*f6*f6);
529       htl = -d1*d7/(f1*f1*f1) - d3*d8/(f3*f3*f3);
530 
531       hc = hydhx*(1.0+d7*d7)/(f1*f1*f1) + hxdhy*(1.0+d8*d8)/(f3*f3*f3) +
532         hydhx*(1.0+d5*d5)/(f5*f5*f5) + hxdhy*(1.0+d6*d6)/(f6*f6*f6) +
533         (hxdhy*(1.0+d1*d1)+hydhx*(1.0+d4*d4)-2*d1*d4)/(f2*f2*f2) +
534         (hxdhy*(1.0+d2*d2)+hydhx*(1.0+d3*d3)-2*d2*d3)/(f4*f4*f4);
535 
536       hl*=0.5; hr*=0.5; ht*=0.5; hb*=0.5; hbr*=0.5; htl*=0.5;  hc*=0.5;
537 
538       k=0;
539       if (j>0) {
540         v[k]=hb; col[k]=row - gxm; k++;
541       }
542 
543       if (j>0 && i < mx -1) {
544         v[k]=hbr; col[k]=row - gxm+1; k++;
545       }
546 
547       if (i>0) {
548         v[k]= hl; col[k]=row - 1; k++;
549       }
550 
551       v[k]= hc; col[k]=row; k++;
552 
553       if (i < mx-1) {
554         v[k]= hr; col[k]=row+1; k++;
555       }
556 
557       if (i>0 && j < my-1) {
558         v[k]= htl; col[k] = row+gxm-1; k++;
559       }
560 
561       if (j < my-1) {
562         v[k]= ht; col[k] = row+gxm; k++;
563       }
564 
565       /*
566          Set matrix values using local numbering, which was defined
567          earlier, in the main routine.
568       */
569       PetscCall(MatSetValuesLocal(Hessian,1,&row,k,col,v,INSERT_VALUES));
570 
571     }
572   }
573 
574   /* Restore vectors */
575   PetscCall(VecRestoreArray(localX,&x));
576   PetscCall(VecRestoreArray(user->Left,&left));
577   PetscCall(VecRestoreArray(user->Top,&top));
578   PetscCall(VecRestoreArray(user->Bottom,&bottom));
579   PetscCall(VecRestoreArray(user->Right,&right));
580 
581   /* Assemble the matrix */
582   PetscCall(MatAssemblyBegin(Hessian,MAT_FINAL_ASSEMBLY));
583   PetscCall(MatAssemblyEnd(Hessian,MAT_FINAL_ASSEMBLY));
584 
585   PetscCall(PetscLogFlops(199.0*xm*ym));
586   return 0;
587 }
588 
589 /* ------------------------------------------------------------------- */
590 /*
591    MSA_BoundaryConditions -  Calculates the boundary conditions for
592    the region.
593 
594    Input Parameter:
595 .  user - user-defined application context
596 
597    Output Parameter:
598 .  user - user-defined application context
599 */
600 static PetscErrorCode MSA_BoundaryConditions(AppCtx * user)
601 {
602   PetscInt   i,j,k,maxits=5,limit=0;
603   PetscInt   xs,ys,xm,ym,gxs,gys,gxm,gym;
604   PetscInt   mx=user->mx,my=user->my;
605   PetscInt   bsize=0, lsize=0, tsize=0, rsize=0;
606   PetscReal  one=1.0, two=2.0, three=3.0, scl=1.0, tol=1e-10;
607   PetscReal  fnorm,det,hx,hy,xt=0,yt=0;
608   PetscReal  u1,u2,nf1,nf2,njac11,njac12,njac21,njac22;
609   PetscReal  b=-0.5, t=0.5, l=-0.5, r=0.5;
610   PetscReal  *boundary;
611   PetscBool  flg;
612   Vec        Bottom,Top,Right,Left;
613 
614   /* Get local mesh boundaries */
615   PetscCall(DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL));
616   PetscCall(DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL));
617 
618   bsize=xm+2;
619   lsize=ym+2;
620   rsize=ym+2;
621   tsize=xm+2;
622 
623   PetscCall(VecCreateMPI(MPI_COMM_WORLD,bsize,PETSC_DECIDE,&Bottom));
624   PetscCall(VecCreateMPI(MPI_COMM_WORLD,tsize,PETSC_DECIDE,&Top));
625   PetscCall(VecCreateMPI(MPI_COMM_WORLD,lsize,PETSC_DECIDE,&Left));
626   PetscCall(VecCreateMPI(MPI_COMM_WORLD,rsize,PETSC_DECIDE,&Right));
627 
628   user->Top=Top;
629   user->Left=Left;
630   user->Bottom=Bottom;
631   user->Right=Right;
632 
633   hx= (r-l)/(mx+1); hy=(t-b)/(my+1);
634 
635   for (j=0; j<4; j++) {
636     if (j==0) {
637       yt=b;
638       xt=l+hx*xs;
639       limit=bsize;
640       VecGetArray(Bottom,&boundary);
641     } else if (j==1) {
642       yt=t;
643       xt=l+hx*xs;
644       limit=tsize;
645       VecGetArray(Top,&boundary);
646     } else if (j==2) {
647       yt=b+hy*ys;
648       xt=l;
649       limit=lsize;
650       VecGetArray(Left,&boundary);
651     } else if (j==3) {
652       yt=b+hy*ys;
653       xt=r;
654       limit=rsize;
655       VecGetArray(Right,&boundary);
656     }
657 
658     for (i=0; i<limit; i++) {
659       u1=xt;
660       u2=-yt;
661       for (k=0; k<maxits; k++) {
662         nf1=u1 + u1*u2*u2 - u1*u1*u1/three-xt;
663         nf2=-u2 - u1*u1*u2 + u2*u2*u2/three-yt;
664         fnorm=PetscSqrtScalar(nf1*nf1+nf2*nf2);
665         if (fnorm <= tol) break;
666         njac11=one+u2*u2-u1*u1;
667         njac12=two*u1*u2;
668         njac21=-two*u1*u2;
669         njac22=-one - u1*u1 + u2*u2;
670         det = njac11*njac22-njac21*njac12;
671         u1 = u1-(njac22*nf1-njac12*nf2)/det;
672         u2 = u2-(njac11*nf2-njac21*nf1)/det;
673       }
674 
675       boundary[i]=u1*u1-u2*u2;
676       if (j==0 || j==1) {
677         xt=xt+hx;
678       } else if (j==2 || j==3) {
679         yt=yt+hy;
680       }
681     }
682     if (j==0) {
683       PetscCall(VecRestoreArray(Bottom,&boundary));
684     } else if (j==1) {
685       PetscCall(VecRestoreArray(Top,&boundary));
686     } else if (j==2) {
687       PetscCall(VecRestoreArray(Left,&boundary));
688     } else if (j==3) {
689       PetscCall(VecRestoreArray(Right,&boundary));
690     }
691   }
692 
693   /* Scale the boundary if desired */
694 
695   PetscCall(PetscOptionsGetReal(NULL,NULL,"-bottom",&scl,&flg));
696   if (flg) PetscCall(VecScale(Bottom, scl));
697   PetscCall(PetscOptionsGetReal(NULL,NULL,"-top",&scl,&flg));
698   if (flg) PetscCall(VecScale(Top, scl));
699   PetscCall(PetscOptionsGetReal(NULL,NULL,"-right",&scl,&flg));
700   if (flg) PetscCall(VecScale(Right, scl));
701 
702   PetscCall(PetscOptionsGetReal(NULL,NULL,"-left",&scl,&flg));
703   if (flg) PetscCall(VecScale(Left, scl));
704   return 0;
705 }
706 
707 /* ------------------------------------------------------------------- */
708 /*
709    MSA_Plate -  Calculates an obstacle for surface to stretch over.
710 
711    Input Parameter:
712 .  user - user-defined application context
713 
714    Output Parameter:
715 .  user - user-defined application context
716 */
717 static PetscErrorCode MSA_Plate(Vec XL,Vec XU,void *ctx)
718 {
719   AppCtx         *user=(AppCtx *)ctx;
720   PetscInt       i,j,row;
721   PetscInt       xs,ys,xm,ym;
722   PetscInt       mx=user->mx, my=user->my, bmy, bmx;
723   PetscReal      t1,t2,t3;
724   PetscReal      *xl, lb=PETSC_NINFINITY, ub=PETSC_INFINITY;
725   PetscBool      cylinder;
726 
727   user->bmy = PetscMax(0,user->bmy);user->bmy = PetscMin(my,user->bmy);
728   user->bmx = PetscMax(0,user->bmx);user->bmx = PetscMin(mx,user->bmx);
729   bmy=user->bmy; bmx=user->bmx;
730 
731   PetscCall(DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL));
732 
733   PetscCall(VecSet(XL, lb));
734   PetscCall(VecSet(XU, ub));
735 
736   PetscCall(VecGetArray(XL,&xl));
737 
738   PetscCall(PetscOptionsHasName(NULL,NULL,"-cylinder",&cylinder));
739   /* Compute the optional lower box */
740   if (cylinder) {
741     for (i=xs; i< xs+xm; i++) {
742       for (j=ys; j<ys+ym; j++) {
743         row=(j-ys)*xm + (i-xs);
744         t1=(2.0*i-mx)*bmy;
745         t2=(2.0*j-my)*bmx;
746         t3=bmx*bmx*bmy*bmy;
747         if (t1*t1 + t2*t2 <= t3) {
748           xl[row] = user->bheight;
749         }
750       }
751     }
752   } else {
753     /* Compute the optional lower box */
754     for (i=xs; i< xs+xm; i++) {
755       for (j=ys; j<ys+ym; j++) {
756         row=(j-ys)*xm + (i-xs);
757         if (i>=(mx-bmx)/2 && i<mx-(mx-bmx)/2 &&
758             j>=(my-bmy)/2 && j<my-(my-bmy)/2) {
759           xl[row] = user->bheight;
760         }
761       }
762     }
763   }
764     PetscCall(VecRestoreArray(XL,&xl));
765 
766   return 0;
767 }
768 
769 /* ------------------------------------------------------------------- */
770 /*
771    MSA_InitialPoint - Calculates the initial guess in one of three ways.
772 
773    Input Parameters:
774 .  user - user-defined application context
775 .  X - vector for initial guess
776 
777    Output Parameters:
778 .  X - newly computed initial guess
779 */
780 static PetscErrorCode MSA_InitialPoint(AppCtx * user, Vec X)
781 {
782   PetscInt       start=-1,i,j;
783   PetscReal      zero=0.0;
784   PetscBool      flg;
785 
786   PetscCall(PetscOptionsGetInt(NULL,NULL,"-start",&start,&flg));
787   if (flg && start==0) { /* The zero vector is reasonable */
788     PetscCall(VecSet(X, zero));
789   } else if (flg && start>0) { /* Try a random start between -0.5 and 0.5 */
790     PetscRandom rctx;  PetscReal np5=-0.5;
791 
792     PetscCall(PetscRandomCreate(MPI_COMM_WORLD,&rctx));
793     for (i=0; i<start; i++) {
794       PetscCall(VecSetRandom(X, rctx));
795     }
796     PetscCall(PetscRandomDestroy(&rctx));
797     PetscCall(VecShift(X, np5));
798 
799   } else { /* Take an average of the boundary conditions */
800 
801     PetscInt row,xs,xm,gxs,gxm,ys,ym,gys,gym;
802     PetscInt mx=user->mx,my=user->my;
803     PetscReal *x,*left,*right,*bottom,*top;
804     Vec    localX = user->localX;
805 
806     /* Get local mesh boundaries */
807     PetscCall(DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL));
808     PetscCall(DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL));
809 
810     /* Get pointers to vector data */
811     PetscCall(VecGetArray(user->Top,&top));
812     PetscCall(VecGetArray(user->Bottom,&bottom));
813     PetscCall(VecGetArray(user->Left,&left));
814     PetscCall(VecGetArray(user->Right,&right));
815 
816     PetscCall(VecGetArray(localX,&x));
817     /* Perform local computations */
818     for (j=ys; j<ys+ym; j++) {
819       for (i=xs; i< xs+xm; i++) {
820         row=(j-gys)*gxm + (i-gxs);
821         x[row] = ((j+1)*bottom[i-xs+1]/my + (my-j+1)*top[i-xs+1]/(my+2)+(i+1)*left[j-ys+1]/mx + (mx-i+1)*right[j-ys+1]/(mx+2))/2.0;
822       }
823     }
824 
825     /* Restore vectors */
826     PetscCall(VecRestoreArray(localX,&x));
827 
828     PetscCall(VecRestoreArray(user->Left,&left));
829     PetscCall(VecRestoreArray(user->Top,&top));
830     PetscCall(VecRestoreArray(user->Bottom,&bottom));
831     PetscCall(VecRestoreArray(user->Right,&right));
832 
833     /* Scatter values into global vector */
834     PetscCall(DMLocalToGlobalBegin(user->dm,localX,INSERT_VALUES,X));
835     PetscCall(DMLocalToGlobalEnd(user->dm,localX,INSERT_VALUES,X));
836 
837   }
838   return 0;
839 }
840 
841 /* For testing matrix free submatrices */
842 PetscErrorCode MatrixFreeHessian(Tao tao, Vec x, Mat H, Mat Hpre, void *ptr)
843 {
844   AppCtx         *user = (AppCtx*)ptr;
845   PetscFunctionBegin;
846   PetscCall(FormHessian(tao,x,user->H,user->H,ptr));
847   PetscFunctionReturn(0);
848 }
849 PetscErrorCode MyMatMult(Mat H_shell, Vec X, Vec Y)
850 {
851   void           *ptr;
852   AppCtx         *user;
853   PetscFunctionBegin;
854   PetscCall(MatShellGetContext(H_shell,&ptr));
855   user = (AppCtx*)ptr;
856   PetscCall(MatMult(user->H,X,Y));
857   PetscFunctionReturn(0);
858 }
859 
860 /*TEST
861 
862    build:
863       requires: !complex
864 
865    test:
866       args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type tron -tao_gatol 1.e-5
867       requires: !single
868 
869    test:
870       suffix: 2
871       nsize: 2
872       args: -tao_smonitor -mx 8 -my 8 -bmx 2 -bmy 5 -bheight 0.3 -tao_type blmvm -tao_gatol 1.e-4
873       requires: !single
874 
875    test:
876       suffix: 3
877       nsize: 3
878       args: -tao_smonitor -mx 8 -my 12 -bmx 4 -bmy 10 -bheight 0.1 -tao_type tron -tao_gatol 1.e-5
879       requires: !single
880 
881    test:
882       suffix: 4
883       nsize: 3
884       args: -tao_smonitor -mx 8 -my 12 -bmx 4 -bmy 10 -bheight 0.1 -tao_subset_type mask -tao_type tron -tao_gatol 1.e-5
885       requires: !single
886 
887    test:
888       suffix: 5
889       nsize: 3
890       args: -tao_smonitor -mx 8 -my 12 -bmx 4 -bmy 10 -bheight 0.1 -tao_subset_type matrixfree -matrixfree -pc_type none -tao_type tron -tao_gatol 1.e-5
891       requires: !single
892 
893    test:
894       suffix: 6
895       nsize: 3
896       args: -tao_smonitor -mx 8 -my 12 -bmx 4 -bmy 10 -bheight 0.1 -tao_subset_type matrixfree -matrixfree -tao_type blmvm -tao_gatol 1.e-4
897       requires: !single
898 
899    test:
900       suffix: 7
901       nsize: 3
902       args: -tao_smonitor -mx 8 -my 12 -bmx 4 -bmy 10 -bheight 0.1 -tao_subset_type matrixfree -pc_type none -tao_type gpcg -tao_gatol 1.e-5
903       requires: !single
904 
905    test:
906       suffix: 8
907       args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bncg -tao_bncg_type gd -tao_gatol 1e-4
908       requires: !single
909 
910    test:
911       suffix: 9
912       args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bncg -tao_gatol 1e-4
913       requires: !single
914 
915    test:
916       suffix: 10
917       args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bnls -tao_gatol 1e-5
918       requires: !single
919 
920    test:
921       suffix: 11
922       args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bntr -tao_gatol 1e-5
923       requires: !single
924 
925    test:
926       suffix: 12
927       args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bntl -tao_gatol 1e-5
928       requires: !single
929 
930    test:
931       suffix: 13
932       args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bnls -tao_gatol 1e-5 -tao_bnk_max_cg_its 3
933       requires: !single
934 
935    test:
936       suffix: 14
937       args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bntr -tao_gatol 1e-5 -tao_bnk_max_cg_its 3
938       requires: !single
939 
940    test:
941       suffix: 15
942       args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bntl -tao_gatol 1e-5 -tao_bnk_max_cg_its 3
943       requires: !single
944 
945    test:
946      suffix: 16
947      args: -tao_smonitor -mx 8 -my 8 -bmx 2 -bmy 5 -bheight 0.3 -tao_gatol 1e-4 -tao_type bqnls
948      requires: !single
949 
950    test:
951      suffix: 17
952      args: -tao_smonitor -mx 8 -my 8 -bmx 2 -bmy 5 -bheight 0.3 -tao_gatol 1e-4 -tao_type bqnkls -tao_bqnk_mat_type lmvmbfgs
953      requires: !single
954 
955    test:
956      suffix: 18
957      args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bnls -tao_gatol 1e-5 -tao_mf_hessian
958      requires: !single
959 
960    test:
961      suffix: 19
962      args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bntr -tao_gatol 1e-5 -tao_mf_hessian
963      requires: !single
964 
965    test:
966      suffix: 20
967      args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bntl -tao_gatol 1e-5 -tao_mf_hessian
968 
969 TEST*/
970