xref: /petsc/src/tao/unconstrained/tutorials/minsurf1.c (revision ebead697dbf761eb322f829370bbe90b3bd93fa3)
1 /* Program usage: mpiexec -n 1 minsurf1 [-help] [all TAO options] */
2 
3 /*  Include "petsctao.h" so we can use TAO solvers.  */
4 #include <petsctao.h>
5 
6 static char  help[] =
7 "This example demonstrates use of the TAO package to \n\
8 solve an unconstrained minimization problem.  This example is based on a \n\
9 problem from the MINPACK-2 test suite.  Given a rectangular 2-D domain and \n\
10 boundary values along the edges of the domain, the objective is to find the\n\
11 surface with the minimal area that satisfies the boundary conditions.\n\
12 The command line options are:\n\
13   -mx <xg>, where <xg> = number of grid points in the 1st coordinate direction\n\
14   -my <yg>, where <yg> = number of grid points in the 2nd coordinate direction\n\
15   -start <st>, where <st> =0 for zero vector, and an average of the boundary conditions otherwise \n\n";
16 
17 /*
18    User-defined application context - contains data needed by the
19    application-provided call-back routines, FormFunctionGradient()
20    and FormHessian().
21 */
22 typedef struct {
23   PetscInt    mx, my;                 /* discretization in x, y directions */
24   PetscReal   *bottom, *top, *left, *right;             /* boundary values */
25   Mat         H;
26 } AppCtx;
27 
28 /* -------- User-defined Routines --------- */
29 
30 static PetscErrorCode MSA_BoundaryConditions(AppCtx*);
31 static PetscErrorCode MSA_InitialPoint(AppCtx*,Vec);
32 static PetscErrorCode QuadraticH(AppCtx*,Vec,Mat);
33 PetscErrorCode FormFunctionGradient(Tao,Vec,PetscReal*,Vec,void*);
34 PetscErrorCode FormHessian(Tao,Vec,Mat,Mat,void*);
35 
36 int main(int argc, char **argv)
37 {
38   PetscInt           N;                 /* Size of vector */
39   PetscMPIInt        size;              /* Number of processors */
40   Vec                x;                 /* solution, gradient vectors */
41   KSP                ksp;               /*  PETSc Krylov subspace method */
42   PetscBool          flg;               /* A return value when checking for user options */
43   Tao                tao;               /* Tao solver context */
44   AppCtx             user;              /* user-defined work context */
45 
46   /* Initialize TAO,PETSc */
47   PetscFunctionBeginUser;
48   PetscCall(PetscInitialize(&argc, &argv,(char *)0,help));
49 
50   PetscCallMPI(MPI_Comm_size(MPI_COMM_WORLD,&size));
51   PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Incorrect number of processors");
52 
53   /* Specify default dimension of the problem */
54   user.mx = 4; user.my = 4;
55 
56   /* Check for any command line arguments that override defaults */
57   PetscCall(PetscOptionsGetInt(NULL,NULL,"-mx",&user.mx,&flg));
58   PetscCall(PetscOptionsGetInt(NULL,NULL,"-my",&user.my,&flg));
59 
60   PetscCall(PetscPrintf(PETSC_COMM_SELF,"\n---- Minimum Surface Area Problem -----\n"));
61   PetscCall(PetscPrintf(PETSC_COMM_SELF,"mx: %" PetscInt_FMT "     my: %" PetscInt_FMT "   \n\n",user.mx,user.my));
62 
63   /* Calculate any derived values from parameters */
64   N    = user.mx*user.my;
65 
66   /* Create TAO solver and set desired solution method  */
67   PetscCall(TaoCreate(PETSC_COMM_SELF,&tao));
68   PetscCall(TaoSetType(tao,TAOLMVM));
69 
70   /* Initialize minsurf application data structure for use in the function evaluations  */
71   PetscCall(MSA_BoundaryConditions(&user));            /* Application specific routine */
72 
73   /*
74      Create a vector to hold the variables.  Compute an initial solution.
75      Set this vector, which will be used by TAO.
76   */
77   PetscCall(VecCreateSeq(PETSC_COMM_SELF,N,&x));
78   PetscCall(MSA_InitialPoint(&user,x));                /* Application specific routine */
79   PetscCall(TaoSetSolution(tao,x));   /* A TAO routine                */
80 
81   /* Provide TAO routines for function, gradient, and Hessian evaluation */
82   PetscCall(TaoSetObjectiveAndGradient(tao,NULL,FormFunctionGradient,(void *)&user));
83 
84   /* Create a matrix data structure to store the Hessian.  This structure will be used by TAO */
85   PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF,N,N,7,NULL,&(user.H)));
86   PetscCall(MatSetOption(user.H,MAT_SYMMETRIC,PETSC_TRUE));
87   PetscCall(TaoSetHessian(tao,user.H,user.H,FormHessian,(void *)&user));
88 
89   /* Check for any TAO command line options */
90   PetscCall(TaoSetFromOptions(tao));
91 
92   /* Limit the number of iterations in the KSP linear solver */
93   PetscCall(TaoGetKSP(tao,&ksp));
94   if (ksp) PetscCall(KSPSetTolerances(ksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,user.mx*user.my));
95 
96   /* SOLVE THE APPLICATION */
97   PetscCall(TaoSolve(tao));
98 
99   PetscCall(TaoDestroy(&tao));
100   PetscCall(VecDestroy(&x));
101   PetscCall(MatDestroy(&user.H));
102   PetscCall(PetscFree(user.bottom));
103   PetscCall(PetscFree(user.top));
104   PetscCall(PetscFree(user.left));
105   PetscCall(PetscFree(user.right));
106 
107   PetscCall(PetscFinalize());
108   return 0;
109 }
110 
111 /* -------------------------------------------------------------------- */
112 
113 /*  FormFunctionGradient - Evaluates function and corresponding gradient.
114 
115     Input Parameters:
116 .   tao     - the Tao context
117 .   X       - input vector
118 .   userCtx - optional user-defined context, as set by TaoSetFunctionGradient()
119 
120     Output Parameters:
121 .   fcn     - the newly evaluated function
122 .   G       - vector containing the newly evaluated gradient
123 */
124 PetscErrorCode FormFunctionGradient(Tao tao,Vec X,PetscReal *fcn,Vec G,void *userCtx)
125 {
126   AppCtx            *user = (AppCtx *) userCtx;
127   PetscInt          i,j,row;
128   PetscInt          mx=user->mx, my=user->my;
129   PetscReal         rhx=mx+1, rhy=my+1;
130   PetscReal         hx=1.0/(mx+1),hy=1.0/(my+1), hydhx=hy/hx, hxdhy=hx/hy, area=0.5*hx*hy, ft=0;
131   PetscReal         f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb;
132   PetscReal         df1dxc,df2dxc,df3dxc,df4dxc,df5dxc,df6dxc;
133   PetscReal         zero=0.0;
134   PetscScalar       *g;
135   const PetscScalar *x;
136 
137   PetscFunctionBeginUser;
138   PetscCall(VecSet(G, zero));
139 
140   PetscCall(VecGetArrayRead(X,&x));
141   PetscCall(VecGetArray(G,&g));
142 
143   /* Compute function over the locally owned part of the mesh */
144   for (j=0; j<my; j++) {
145     for (i=0; i< mx; i++) {
146       row=(j)*mx + (i);
147       xc = x[row];
148       xlt=xrb=xl=xr=xb=xt=xc;
149       if (i==0) { /* left side */
150         xl  = user->left[j+1];
151         xlt = user->left[j+2];
152       } else {
153         xl = x[row-1];
154       }
155 
156       if (j==0) { /* bottom side */
157         xb  = user->bottom[i+1];
158         xrb = user->bottom[i+2];
159       } else {
160         xb = x[row-mx];
161       }
162 
163       if (i+1 == mx) { /* right side */
164         xr  = user->right[j+1];
165         xrb = user->right[j];
166       } else {
167         xr = x[row+1];
168       }
169 
170       if (j+1==0+my) { /* top side */
171         xt  = user->top[i+1];
172         xlt = user->top[i];
173       }else {
174         xt = x[row+mx];
175       }
176 
177       if (i>0 && j+1<my) {
178         xlt = x[row-1+mx];
179       }
180       if (j>0 && i+1<mx) {
181         xrb = x[row+1-mx];
182       }
183 
184       d1 = (xc-xl);
185       d2 = (xc-xr);
186       d3 = (xc-xt);
187       d4 = (xc-xb);
188       d5 = (xr-xrb);
189       d6 = (xrb-xb);
190       d7 = (xlt-xl);
191       d8 = (xt-xlt);
192 
193       df1dxc = d1*hydhx;
194       df2dxc = (d1*hydhx + d4*hxdhy);
195       df3dxc = d3*hxdhy;
196       df4dxc = (d2*hydhx + d3*hxdhy);
197       df5dxc = d2*hydhx;
198       df6dxc = d4*hxdhy;
199 
200       d1 *= rhx;
201       d2 *= rhx;
202       d3 *= rhy;
203       d4 *= rhy;
204       d5 *= rhy;
205       d6 *= rhx;
206       d7 *= rhy;
207       d8 *= rhx;
208 
209       f1 = PetscSqrtScalar(1.0 + d1*d1 + d7*d7);
210       f2 = PetscSqrtScalar(1.0 + d1*d1 + d4*d4);
211       f3 = PetscSqrtScalar(1.0 + d3*d3 + d8*d8);
212       f4 = PetscSqrtScalar(1.0 + d3*d3 + d2*d2);
213       f5 = PetscSqrtScalar(1.0 + d2*d2 + d5*d5);
214       f6 = PetscSqrtScalar(1.0 + d4*d4 + d6*d6);
215 
216       ft = ft + (f2 + f4);
217 
218       df1dxc /= f1;
219       df2dxc /= f2;
220       df3dxc /= f3;
221       df4dxc /= f4;
222       df5dxc /= f5;
223       df6dxc /= f6;
224 
225       g[row] = (df1dxc+df2dxc+df3dxc+df4dxc+df5dxc+df6dxc)/2.0;
226     }
227   }
228 
229   for (j=0; j<my; j++) {   /* left side */
230     d3 = (user->left[j+1] - user->left[j+2])*rhy;
231     d2 = (user->left[j+1] - x[j*mx])*rhx;
232     ft = ft+PetscSqrtScalar(1.0 + d3*d3 + d2*d2);
233   }
234 
235   for (i=0; i<mx; i++) { /* bottom */
236     d2 = (user->bottom[i+1]-user->bottom[i+2])*rhx;
237     d3 = (user->bottom[i+1]-x[i])*rhy;
238     ft = ft+PetscSqrtScalar(1.0 + d3*d3 + d2*d2);
239   }
240 
241   for (j=0; j< my; j++) { /* right side */
242     d1 = (x[(j+1)*mx-1]-user->right[j+1])*rhx;
243     d4 = (user->right[j]-user->right[j+1])*rhy;
244     ft = ft+PetscSqrtScalar(1.0 + d1*d1 + d4*d4);
245   }
246 
247   for (i=0; i<mx; i++) { /* top side */
248     d1 = (x[(my-1)*mx + i] - user->top[i+1])*rhy;
249     d4 = (user->top[i+1] - user->top[i])*rhx;
250     ft = ft+PetscSqrtScalar(1.0 + d1*d1 + d4*d4);
251   }
252 
253   /* Bottom left corner */
254   d1  = (user->left[0]-user->left[1])*rhy;
255   d2  = (user->bottom[0]-user->bottom[1])*rhx;
256   ft += PetscSqrtScalar(1.0 + d1*d1 + d2*d2);
257 
258   /* Top right corner */
259   d1  = (user->right[my+1] - user->right[my])*rhy;
260   d2  = (user->top[mx+1] - user->top[mx])*rhx;
261   ft += PetscSqrtScalar(1.0 + d1*d1 + d2*d2);
262 
263   (*fcn)=ft*area;
264 
265   /* Restore vectors */
266   PetscCall(VecRestoreArrayRead(X,&x));
267   PetscCall(VecRestoreArray(G,&g));
268   PetscCall(PetscLogFlops(67.0*mx*my));
269   PetscFunctionReturn(0);
270 }
271 
272 /* ------------------------------------------------------------------- */
273 /*
274    FormHessian - Evaluates the Hessian matrix.
275 
276    Input Parameters:
277 .  tao  - the Tao context
278 .  x    - input vector
279 .  ptr  - optional user-defined context, as set by TaoSetHessian()
280 
281    Output Parameters:
282 .  H    - Hessian matrix
283 .  Hpre - optionally different preconditioning matrix
284 .  flg  - flag indicating matrix structure
285 
286 */
287 PetscErrorCode FormHessian(Tao tao,Vec X,Mat H, Mat Hpre, void *ptr)
288 {
289   AppCtx         *user = (AppCtx *) ptr;
290 
291   PetscFunctionBeginUser;
292   /* Evaluate the Hessian entries*/
293   PetscCall(QuadraticH(user,X,H));
294   PetscFunctionReturn(0);
295 }
296 
297 /* ------------------------------------------------------------------- */
298 /*
299    QuadraticH - Evaluates the Hessian matrix.
300 
301    Input Parameters:
302 .  user - user-defined context, as set by TaoSetHessian()
303 .  X    - input vector
304 
305    Output Parameter:
306 .  H    - Hessian matrix
307 */
308 PetscErrorCode QuadraticH(AppCtx *user, Vec X, Mat Hessian)
309 {
310   PetscInt          i,j,k,row;
311   PetscInt          mx=user->mx, my=user->my;
312   PetscInt          col[7];
313   PetscReal         hx=1.0/(mx+1), hy=1.0/(my+1), hydhx=hy/hx, hxdhy=hx/hy;
314   PetscReal         rhx=mx+1, rhy=my+1;
315   PetscReal         f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb;
316   PetscReal         hl,hr,ht,hb,hc,htl,hbr;
317   const PetscScalar *x;
318   PetscReal         v[7];
319 
320   PetscFunctionBeginUser;
321   /* Get pointers to vector data */
322   PetscCall(VecGetArrayRead(X,&x));
323 
324   /* Initialize matrix entries to zero */
325   PetscCall(MatZeroEntries(Hessian));
326 
327   /* Set various matrix options */
328   PetscCall(MatSetOption(Hessian,MAT_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE));
329 
330   /* Compute Hessian over the locally owned part of the mesh */
331   for (i=0; i< mx; i++) {
332     for (j=0; j<my; j++) {
333 
334       row=(j)*mx + (i);
335 
336       xc = x[row];
337       xlt=xrb=xl=xr=xb=xt=xc;
338 
339       /* Left side */
340       if (i==0) {
341         xl= user->left[j+1];
342         xlt = user->left[j+2];
343       } else {
344         xl = x[row-1];
345       }
346 
347       if (j==0) {
348         xb=user->bottom[i+1];
349         xrb = user->bottom[i+2];
350       } else {
351         xb = x[row-mx];
352       }
353 
354       if (i+1 == mx) {
355         xr=user->right[j+1];
356         xrb = user->right[j];
357       } else {
358         xr = x[row+1];
359       }
360 
361       if (j+1==my) {
362         xt=user->top[i+1];
363         xlt = user->top[i];
364       }else {
365         xt = x[row+mx];
366       }
367 
368       if (i>0 && j+1<my) {
369         xlt = x[row-1+mx];
370       }
371       if (j>0 && i+1<mx) {
372         xrb = x[row+1-mx];
373       }
374 
375       d1 = (xc-xl)*rhx;
376       d2 = (xc-xr)*rhx;
377       d3 = (xc-xt)*rhy;
378       d4 = (xc-xb)*rhy;
379       d5 = (xrb-xr)*rhy;
380       d6 = (xrb-xb)*rhx;
381       d7 = (xlt-xl)*rhy;
382       d8 = (xlt-xt)*rhx;
383 
384       f1 = PetscSqrtScalar(1.0 + d1*d1 + d7*d7);
385       f2 = PetscSqrtScalar(1.0 + d1*d1 + d4*d4);
386       f3 = PetscSqrtScalar(1.0 + d3*d3 + d8*d8);
387       f4 = PetscSqrtScalar(1.0 + d3*d3 + d2*d2);
388       f5 = PetscSqrtScalar(1.0 + d2*d2 + d5*d5);
389       f6 = PetscSqrtScalar(1.0 + d4*d4 + d6*d6);
390 
391       hl = (-hydhx*(1.0+d7*d7)+d1*d7)/(f1*f1*f1)+(-hydhx*(1.0+d4*d4)+d1*d4)/(f2*f2*f2);
392       hr = (-hydhx*(1.0+d5*d5)+d2*d5)/(f5*f5*f5)+(-hydhx*(1.0+d3*d3)+d2*d3)/(f4*f4*f4);
393       ht = (-hxdhy*(1.0+d8*d8)+d3*d8)/(f3*f3*f3)+(-hxdhy*(1.0+d2*d2)+d2*d3)/(f4*f4*f4);
394       hb = (-hxdhy*(1.0+d6*d6)+d4*d6)/(f6*f6*f6)+(-hxdhy*(1.0+d1*d1)+d1*d4)/(f2*f2*f2);
395 
396       hbr = -d2*d5/(f5*f5*f5) - d4*d6/(f6*f6*f6);
397       htl = -d1*d7/(f1*f1*f1) - d3*d8/(f3*f3*f3);
398 
399       hc = hydhx*(1.0+d7*d7)/(f1*f1*f1) + hxdhy*(1.0+d8*d8)/(f3*f3*f3) + hydhx*(1.0+d5*d5)/(f5*f5*f5) + hxdhy*(1.0+d6*d6)/(f6*f6*f6) +
400            (hxdhy*(1.0+d1*d1)+hydhx*(1.0+d4*d4)-2*d1*d4)/(f2*f2*f2) +  (hxdhy*(1.0+d2*d2)+hydhx*(1.0+d3*d3)-2*d2*d3)/(f4*f4*f4);
401 
402       hl*=0.5; hr*=0.5; ht*=0.5; hb*=0.5; hbr*=0.5; htl*=0.5;  hc*=0.5;
403 
404       k=0;
405       if (j>0) {
406         v[k]=hb; col[k]=row - mx; k++;
407       }
408 
409       if (j>0 && i < mx -1) {
410         v[k]=hbr; col[k]=row - mx+1; k++;
411       }
412 
413       if (i>0) {
414         v[k]= hl; col[k]=row - 1; k++;
415       }
416 
417       v[k]= hc; col[k]=row; k++;
418 
419       if (i < mx-1) {
420         v[k]= hr; col[k]=row+1; k++;
421       }
422 
423       if (i>0 && j < my-1) {
424         v[k]= htl; col[k] = row+mx-1; k++;
425       }
426 
427       if (j < my-1) {
428         v[k]= ht; col[k] = row+mx; k++;
429       }
430 
431       /*
432          Set matrix values using local numbering, which was defined
433          earlier, in the main routine.
434       */
435       PetscCall(MatSetValues(Hessian,1,&row,k,col,v,INSERT_VALUES));
436     }
437   }
438 
439   /* Restore vectors */
440   PetscCall(VecRestoreArrayRead(X,&x));
441 
442   /* Assemble the matrix */
443   PetscCall(MatAssemblyBegin(Hessian,MAT_FINAL_ASSEMBLY));
444   PetscCall(MatAssemblyEnd(Hessian,MAT_FINAL_ASSEMBLY));
445 
446   PetscCall(PetscLogFlops(199.0*mx*my));
447   PetscFunctionReturn(0);
448 }
449 
450 /* ------------------------------------------------------------------- */
451 /*
452    MSA_BoundaryConditions -  Calculates the boundary conditions for
453    the region.
454 
455    Input Parameter:
456 .  user - user-defined application context
457 
458    Output Parameter:
459 .  user - user-defined application context
460 */
461 static PetscErrorCode MSA_BoundaryConditions(AppCtx * user)
462 {
463   PetscInt       i,j,k,limit=0;
464   PetscInt       maxits=5;
465   PetscInt       mx=user->mx,my=user->my;
466   PetscInt       bsize=0, lsize=0, tsize=0, rsize=0;
467   PetscReal      one=1.0, two=2.0, three=3.0, tol=1e-10;
468   PetscReal      fnorm,det,hx,hy,xt=0,yt=0;
469   PetscReal      u1,u2,nf1,nf2,njac11,njac12,njac21,njac22;
470   PetscReal      b=-0.5, t=0.5, l=-0.5, r=0.5;
471   PetscReal      *boundary;
472 
473   PetscFunctionBeginUser;
474   bsize=mx+2; lsize=my+2; rsize=my+2; tsize=mx+2;
475 
476   PetscCall(PetscMalloc1(bsize,&user->bottom));
477   PetscCall(PetscMalloc1(tsize,&user->top));
478   PetscCall(PetscMalloc1(lsize,&user->left));
479   PetscCall(PetscMalloc1(rsize,&user->right));
480 
481   hx= (r-l)/(mx+1); hy=(t-b)/(my+1);
482 
483   for (j=0; j<4; j++) {
484     if (j==0) {
485       yt=b;
486       xt=l;
487       limit=bsize;
488       boundary=user->bottom;
489     } else if (j==1) {
490       yt=t;
491       xt=l;
492       limit=tsize;
493       boundary=user->top;
494     } else if (j==2) {
495       yt=b;
496       xt=l;
497       limit=lsize;
498       boundary=user->left;
499     } else {  /*  if (j==3) */
500       yt=b;
501       xt=r;
502       limit=rsize;
503       boundary=user->right;
504     }
505 
506     for (i=0; i<limit; i++) {
507       u1=xt;
508       u2=-yt;
509       for (k=0; k<maxits; k++) {
510         nf1=u1 + u1*u2*u2 - u1*u1*u1/three-xt;
511         nf2=-u2 - u1*u1*u2 + u2*u2*u2/three-yt;
512         fnorm=PetscSqrtScalar(nf1*nf1+nf2*nf2);
513         if (fnorm <= tol) break;
514         njac11=one+u2*u2-u1*u1;
515         njac12=two*u1*u2;
516         njac21=-two*u1*u2;
517         njac22=-one - u1*u1 + u2*u2;
518         det = njac11*njac22-njac21*njac12;
519         u1 = u1-(njac22*nf1-njac12*nf2)/det;
520         u2 = u2-(njac11*nf2-njac21*nf1)/det;
521       }
522 
523       boundary[i]=u1*u1-u2*u2;
524       if (j==0 || j==1) {
525         xt=xt+hx;
526       } else { /*  if (j==2 || j==3) */
527         yt=yt+hy;
528       }
529     }
530   }
531   PetscFunctionReturn(0);
532 }
533 
534 /* ------------------------------------------------------------------- */
535 /*
536    MSA_InitialPoint - Calculates the initial guess in one of three ways.
537 
538    Input Parameters:
539 .  user - user-defined application context
540 .  X - vector for initial guess
541 
542    Output Parameters:
543 .  X - newly computed initial guess
544 */
545 static PetscErrorCode MSA_InitialPoint(AppCtx * user, Vec X)
546 {
547   PetscInt       start=-1,i,j;
548   PetscReal      zero=0.0;
549   PetscBool      flg;
550 
551   PetscFunctionBeginUser;
552   PetscCall(VecSet(X, zero));
553   PetscCall(PetscOptionsGetInt(NULL,NULL,"-start",&start,&flg));
554 
555   if (flg && start==0) { /* The zero vector is reasonable */
556      PetscCall(VecSet(X, zero));
557    } else { /* Take an average of the boundary conditions */
558     PetscInt    row;
559     PetscInt    mx=user->mx,my=user->my;
560     PetscReal *x;
561 
562     /* Get pointers to vector data */
563     PetscCall(VecGetArray(X,&x));
564     /* Perform local computations */
565     for (j=0; j<my; j++) {
566       for (i=0; i< mx; i++) {
567         row=(j)*mx + (i);
568         x[row] = (((j+1)*user->bottom[i+1]+(my-j+1)*user->top[i+1])/(my+2)+ ((i+1)*user->left[j+1]+(mx-i+1)*user->right[j+1])/(mx+2))/2.0;
569       }
570     }
571     /* Restore vectors */
572     PetscCall(VecRestoreArray(X,&x));
573   }
574   PetscFunctionReturn(0);
575 }
576 
577 /*TEST
578 
579    build:
580       requires: !complex
581 
582    test:
583       args: -tao_smonitor -tao_type nls -mx 10 -my 8 -tao_gatol 1.e-4
584       requires: !single
585 
586    test:
587       suffix: 2
588       args: -tao_smonitor -tao_type bmrm -mx 10 -my 8 -tao_gatol 1.e-3
589       requires: !single
590 
591    test:
592       suffix: 3
593       args: -tao_smonitor -tao_type lmvm -mx 10 -my 8 -tao_gatol 1.e-3
594       requires: !single
595 
596    test:
597       suffix: 4
598       args: -tao_smonitor -tao_type bntr -mx 10 -my 8 -tao_gatol 1.e-4
599       requires: !single
600 
601    test:
602       suffix: 4_ew
603       args: -tao_smonitor -tao_type bntr -tao_ksp_ew -mx 10 -my 8 -tao_gatol 1.e-4
604       requires: !single
605 
606    test:
607       suffix: 5
608       args: -tao_smonitor -tao_type bntl -mx 10 -my 8 -tao_gatol 1.e-4
609       requires: !single
610 
611    test:
612       suffix: 5_ew
613       args: -tao_smonitor -tao_type bntl -tao_ksp_ew -mx 10 -my 8 -tao_gatol 1.e-4
614       requires: !single
615 
616    test:
617       suffix: 6
618       args: -tao_smonitor -tao_type bnls -mx 10 -my 8 -tao_gatol 1.e-4
619       requires: !single
620 
621    test:
622       suffix: 6_ew
623       args: -tao_smonitor -tao_type bnls -tao_ksp_ew -tao_bnk_ksp_ew_version 3 -mx 10 -my 8 -tao_gatol 1.e-4
624       requires: !single
625 
626    test:
627       suffix: 7
628       args: -tao_smonitor -tao_type bntr -mx 10 -my 8 -tao_bnk_max_cg_its 3 -tao_gatol 1.e-4
629       requires: !single
630 
631    test:
632       suffix: 8
633       args: -tao_smonitor -tao_type bntl -mx 10 -my 8 -tao_bnk_max_cg_its 3 -tao_gatol 1.e-4
634       requires: !single
635 
636    test:
637       suffix: 9
638       args: -tao_smonitor -tao_type bnls -mx 10 -my 8 -tao_bnk_max_cg_its 3 -tao_gatol 1.e-4
639       requires: !single
640 
641    test:
642       suffix: 10
643       args: -tao_smonitor -tao_type bnls -mx 10 -my 8 -tao_bnk_max_cg_its 3 -tao_gatol 1.e-4 -tao_mf_hessian
644 
645    test:
646       suffix: 11
647       args: -tao_smonitor -tao_type bntr -mx 10 -my 8 -tao_gatol 1.e-4 -tao_mf_hessian
648       requires: !single
649 
650    test:
651       suffix: 12
652       args: -tao_smonitor -tao_type bntl -mx 10 -my 8 -tao_gatol 1.e-4 -tao_mf_hessian
653       requires: !single
654 
655 TEST*/
656