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