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