xref: /petsc/src/tao/bound/tutorials/jbearing2.c (revision 030f984af8d8bb4c203755d35bded3c05b3d83ce)
1 /*
2   Include "petsctao.h" so we can use TAO solvers
3   Include "petscdmda.h" so that we can use distributed arrays (DMs) for managing
4   Include "petscksp.h" so we can set KSP type
5   the parallel mesh.
6 */
7 
8 #include <petsctao.h>
9 #include <petscdmda.h>
10 
11 static  char help[]=
12 "This example demonstrates use of the TAO package to \n\
13 solve a bound constrained minimization problem.  This example is based on \n\
14 the problem DPJB from the MINPACK-2 test suite.  This pressure journal \n\
15 bearing problem is an example of elliptic variational problem defined over \n\
16 a two dimensional rectangle.  By discretizing the domain into triangular \n\
17 elements, the pressure surrounding the journal bearing is defined as the \n\
18 minimum of a quadratic function whose variables are bounded below by zero.\n\
19 The command line options are:\n\
20   -mx <xg>, where <xg> = number of grid points in the 1st coordinate direction\n\
21   -my <yg>, where <yg> = number of grid points in the 2nd coordinate direction\n\
22  \n";
23 
24 /*T
25    Concepts: TAO^Solving a bound constrained minimization problem
26    Routines: TaoCreate();
27    Routines: TaoSetType(); TaoSetObjectiveAndGradientRoutine();
28    Routines: TaoSetHessianRoutine();
29    Routines: TaoSetVariableBounds();
30    Routines: TaoSetMonitor(); TaoSetConvergenceTest();
31    Routines: TaoSetInitialVector();
32    Routines: TaoSetFromOptions();
33    Routines: TaoSolve();
34    Routines: TaoDestroy();
35    Processors: n
36 T*/
37 
38 /*
39    User-defined application context - contains data needed by the
40    application-provided call-back routines, FormFunctionGradient(),
41    FormHessian().
42 */
43 typedef struct {
44   /* problem parameters */
45   PetscReal      ecc;          /* test problem parameter */
46   PetscReal      b;            /* A dimension of journal bearing */
47   PetscInt       nx,ny;        /* discretization in x, y directions */
48 
49   /* Working space */
50   DM          dm;           /* distributed array data structure */
51   Mat         A;            /* Quadratic Objective term */
52   Vec         B;            /* Linear Objective term */
53 } AppCtx;
54 
55 /* User-defined routines */
56 static PetscReal p(PetscReal xi, PetscReal ecc);
57 static PetscErrorCode FormFunctionGradient(Tao, Vec, PetscReal *,Vec,void *);
58 static PetscErrorCode FormHessian(Tao,Vec,Mat, Mat, void *);
59 static PetscErrorCode ComputeB(AppCtx*);
60 static PetscErrorCode Monitor(Tao, void*);
61 static PetscErrorCode ConvergenceTest(Tao, void*);
62 
63 int main(int argc, char **argv)
64 {
65   PetscErrorCode     ierr;            /* used to check for functions returning nonzeros */
66   PetscInt           Nx, Ny;          /* number of processors in x- and y- directions */
67   PetscInt           m;               /* number of local elements in vectors */
68   Vec                x;               /* variables vector */
69   Vec                xl,xu;           /* bounds vectors */
70   PetscReal          d1000 = 1000;
71   PetscBool          flg,testgetdiag; /* A return variable when checking for user options */
72   Tao                tao;             /* Tao solver context */
73   KSP                ksp;
74   AppCtx             user;            /* user-defined work context */
75   PetscReal          zero = 0.0;      /* lower bound on all variables */
76 
77   /* Initialize PETSC and TAO */
78   ierr = PetscInitialize(&argc, &argv,(char *)0,help);if (ierr) return ierr;
79 
80   /* Set the default values for the problem parameters */
81   user.nx = 50; user.ny = 50; user.ecc = 0.1; user.b = 10.0;
82   testgetdiag = PETSC_FALSE;
83 
84   /* Check for any command line arguments that override defaults */
85   ierr = PetscOptionsGetInt(NULL,NULL,"-mx",&user.nx,&flg);CHKERRQ(ierr);
86   ierr = PetscOptionsGetInt(NULL,NULL,"-my",&user.ny,&flg);CHKERRQ(ierr);
87   ierr = PetscOptionsGetReal(NULL,NULL,"-ecc",&user.ecc,&flg);CHKERRQ(ierr);
88   ierr = PetscOptionsGetReal(NULL,NULL,"-b",&user.b,&flg);CHKERRQ(ierr);
89   ierr = PetscOptionsGetBool(NULL,NULL,"-test_getdiagonal",&testgetdiag,NULL);CHKERRQ(ierr);
90 
91   ierr = PetscPrintf(PETSC_COMM_WORLD,"\n---- Journal Bearing Problem SHB-----\n");CHKERRQ(ierr);
92   ierr = PetscPrintf(PETSC_COMM_WORLD,"mx: %D,  my: %D,  ecc: %g \n\n",user.nx,user.ny,(double)user.ecc);CHKERRQ(ierr);
93 
94   /* Let Petsc determine the grid division */
95   Nx = PETSC_DECIDE; Ny = PETSC_DECIDE;
96 
97   /*
98      A two dimensional distributed array will help define this problem,
99      which derives from an elliptic PDE on two dimensional domain.  From
100      the distributed array, Create the vectors.
101   */
102   ierr = DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,user.nx,user.ny,Nx,Ny,1,1,NULL,NULL,&user.dm);CHKERRQ(ierr);
103   ierr = DMSetFromOptions(user.dm);CHKERRQ(ierr);
104   ierr = DMSetUp(user.dm);CHKERRQ(ierr);
105 
106   /*
107      Extract global and local vectors from DM; the vector user.B is
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   ierr = DMCreateGlobalVector(user.dm,&x);CHKERRQ(ierr); /* Solution */
113   ierr = VecDuplicate(x,&user.B);CHKERRQ(ierr); /* Linear objective */
114 
115   /*  Create matrix user.A to store quadratic, Create a local ordering scheme. */
116   ierr = VecGetLocalSize(x,&m);CHKERRQ(ierr);
117   ierr = DMCreateMatrix(user.dm,&user.A);CHKERRQ(ierr);
118 
119   if (testgetdiag) {
120     ierr = MatSetOperation(user.A,MATOP_GET_DIAGONAL,NULL);CHKERRQ(ierr);
121   }
122 
123   /* User defined function -- compute linear term of quadratic */
124   ierr = ComputeB(&user);CHKERRQ(ierr);
125 
126   /* The TAO code begins here */
127 
128   /*
129      Create the optimization solver
130      Suitable methods: TAOGPCG, TAOBQPIP, TAOTRON, TAOBLMVM
131   */
132   ierr = TaoCreate(PETSC_COMM_WORLD,&tao);CHKERRQ(ierr);
133   ierr = TaoSetType(tao,TAOBLMVM);CHKERRQ(ierr);
134 
135   /* Set the initial vector */
136   ierr = VecSet(x, zero);CHKERRQ(ierr);
137   ierr = TaoSetInitialVector(tao,x);CHKERRQ(ierr);
138 
139   /* Set the user function, gradient, hessian evaluation routines and data structures */
140   ierr = TaoSetObjectiveAndGradientRoutine(tao,FormFunctionGradient,(void*) &user);CHKERRQ(ierr);
141 
142   ierr = TaoSetHessianRoutine(tao,user.A,user.A,FormHessian,(void*)&user);CHKERRQ(ierr);
143 
144   /* Set a routine that defines the bounds */
145   ierr = VecDuplicate(x,&xl);CHKERRQ(ierr);
146   ierr = VecDuplicate(x,&xu);CHKERRQ(ierr);
147   ierr = VecSet(xl, zero);CHKERRQ(ierr);
148   ierr = VecSet(xu, d1000);CHKERRQ(ierr);
149   ierr = TaoSetVariableBounds(tao,xl,xu);CHKERRQ(ierr);
150 
151   ierr = TaoGetKSP(tao,&ksp);CHKERRQ(ierr);
152   if (ksp) {
153     ierr = KSPSetType(ksp,KSPCG);CHKERRQ(ierr);
154   }
155 
156   ierr = PetscOptionsHasName(NULL,NULL,"-testmonitor",&flg);CHKERRQ(ierr);
157   if (flg) {
158     ierr = TaoSetMonitor(tao,Monitor,&user,NULL);CHKERRQ(ierr);
159   }
160   ierr = PetscOptionsHasName(NULL,NULL,"-testconvergence",&flg);CHKERRQ(ierr);
161   if (flg) {
162     ierr = TaoSetConvergenceTest(tao,ConvergenceTest,&user);CHKERRQ(ierr);
163   }
164 
165   /* Check for any tao command line options */
166   ierr = TaoSetFromOptions(tao);CHKERRQ(ierr);
167 
168   /* Solve the bound constrained problem */
169   ierr = TaoSolve(tao);CHKERRQ(ierr);
170 
171   /* Free PETSc data structures */
172   ierr = VecDestroy(&x);CHKERRQ(ierr);
173   ierr = VecDestroy(&xl);CHKERRQ(ierr);
174   ierr = VecDestroy(&xu);CHKERRQ(ierr);
175   ierr = MatDestroy(&user.A);CHKERRQ(ierr);
176   ierr = VecDestroy(&user.B);CHKERRQ(ierr);
177 
178   /* Free TAO data structures */
179   ierr = TaoDestroy(&tao);CHKERRQ(ierr);
180   ierr = DMDestroy(&user.dm);CHKERRQ(ierr);
181   ierr = PetscFinalize();
182   return ierr;
183 }
184 
185 static PetscReal p(PetscReal xi, PetscReal ecc)
186 {
187   PetscReal t=1.0+ecc*PetscCosScalar(xi);
188   return (t*t*t);
189 }
190 
191 PetscErrorCode ComputeB(AppCtx* user)
192 {
193   PetscErrorCode ierr;
194   PetscInt       i,j,k;
195   PetscInt       nx,ny,xs,xm,gxs,gxm,ys,ym,gys,gym;
196   PetscReal      two=2.0, pi=4.0*atan(1.0);
197   PetscReal      hx,hy,ehxhy;
198   PetscReal      temp,*b;
199   PetscReal      ecc=user->ecc;
200 
201   nx=user->nx;
202   ny=user->ny;
203   hx=two*pi/(nx+1.0);
204   hy=two*user->b/(ny+1.0);
205   ehxhy = ecc*hx*hy;
206 
207   /*
208      Get local grid boundaries
209   */
210   ierr = DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr);
211   ierr = DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL);CHKERRQ(ierr);
212 
213   /* Compute the linear term in the objective function */
214   ierr = VecGetArray(user->B,&b);CHKERRQ(ierr);
215   for (i=xs; i<xs+xm; i++) {
216     temp=PetscSinScalar((i+1)*hx);
217     for (j=ys; j<ys+ym; j++) {
218       k=xm*(j-ys)+(i-xs);
219       b[k]=  - ehxhy*temp;
220     }
221   }
222   ierr = VecRestoreArray(user->B,&b);CHKERRQ(ierr);
223   ierr = PetscLogFlops(5.0*xm*ym+3.0*xm);CHKERRQ(ierr);
224 
225   return 0;
226 }
227 
228 PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *fcn,Vec G,void *ptr)
229 {
230   AppCtx*        user=(AppCtx*)ptr;
231   PetscErrorCode ierr;
232   PetscInt       i,j,k,kk;
233   PetscInt       col[5],row,nx,ny,xs,xm,gxs,gxm,ys,ym,gys,gym;
234   PetscReal      one=1.0, two=2.0, six=6.0,pi=4.0*atan(1.0);
235   PetscReal      hx,hy,hxhy,hxhx,hyhy;
236   PetscReal      xi,v[5];
237   PetscReal      ecc=user->ecc, trule1,trule2,trule3,trule4,trule5,trule6;
238   PetscReal      vmiddle, vup, vdown, vleft, vright;
239   PetscReal      tt,f1,f2;
240   PetscReal      *x,*g,zero=0.0;
241   Vec            localX;
242 
243   nx=user->nx;
244   ny=user->ny;
245   hx=two*pi/(nx+1.0);
246   hy=two*user->b/(ny+1.0);
247   hxhy=hx*hy;
248   hxhx=one/(hx*hx);
249   hyhy=one/(hy*hy);
250 
251   ierr = DMGetLocalVector(user->dm,&localX);CHKERRQ(ierr);
252 
253   ierr = DMGlobalToLocalBegin(user->dm,X,INSERT_VALUES,localX);CHKERRQ(ierr);
254   ierr = DMGlobalToLocalEnd(user->dm,X,INSERT_VALUES,localX);CHKERRQ(ierr);
255 
256   ierr = VecSet(G, zero);CHKERRQ(ierr);
257   /*
258     Get local grid boundaries
259   */
260   ierr = DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr);
261   ierr = DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL);CHKERRQ(ierr);
262 
263   ierr = VecGetArray(localX,&x);CHKERRQ(ierr);
264   ierr = VecGetArray(G,&g);CHKERRQ(ierr);
265 
266   for (i=xs; i< xs+xm; i++) {
267     xi=(i+1)*hx;
268     trule1=hxhy*(p(xi,ecc) + p(xi+hx,ecc) + p(xi,ecc)) / six; /* L(i,j) */
269     trule2=hxhy*(p(xi,ecc) + p(xi-hx,ecc) + p(xi,ecc)) / six; /* U(i,j) */
270     trule3=hxhy*(p(xi,ecc) + p(xi+hx,ecc) + p(xi+hx,ecc)) / six; /* U(i+1,j) */
271     trule4=hxhy*(p(xi,ecc) + p(xi-hx,ecc) + p(xi-hx,ecc)) / six; /* L(i-1,j) */
272     trule5=trule1; /* L(i,j-1) */
273     trule6=trule2; /* U(i,j+1) */
274 
275     vdown=-(trule5+trule2)*hyhy;
276     vleft=-hxhx*(trule2+trule4);
277     vright= -hxhx*(trule1+trule3);
278     vup=-hyhy*(trule1+trule6);
279     vmiddle=(hxhx)*(trule1+trule2+trule3+trule4)+hyhy*(trule1+trule2+trule5+trule6);
280 
281     for (j=ys; j<ys+ym; j++) {
282 
283       row=(j-gys)*gxm + (i-gxs);
284        v[0]=0; v[1]=0; v[2]=0; v[3]=0; v[4]=0;
285 
286        k=0;
287        if (j>gys) {
288          v[k]=vdown; col[k]=row - gxm; k++;
289        }
290 
291        if (i>gxs) {
292          v[k]= vleft; col[k]=row - 1; k++;
293        }
294 
295        v[k]= vmiddle; col[k]=row; k++;
296 
297        if (i+1 < gxs+gxm) {
298          v[k]= vright; col[k]=row+1; k++;
299        }
300 
301        if (j+1 <gys+gym) {
302          v[k]= vup; col[k] = row+gxm; k++;
303        }
304        tt=0;
305        for (kk=0;kk<k;kk++) {
306          tt+=v[kk]*x[col[kk]];
307        }
308        row=(j-ys)*xm + (i-xs);
309        g[row]=tt;
310 
311      }
312 
313   }
314 
315   ierr = VecRestoreArray(localX,&x);CHKERRQ(ierr);
316   ierr = VecRestoreArray(G,&g);CHKERRQ(ierr);
317 
318   ierr = DMRestoreLocalVector(user->dm,&localX);CHKERRQ(ierr);
319 
320   ierr = VecDot(X,G,&f1);CHKERRQ(ierr);
321   ierr = VecDot(user->B,X,&f2);CHKERRQ(ierr);
322   ierr = VecAXPY(G, one, user->B);CHKERRQ(ierr);
323   *fcn = f1/2.0 + f2;
324 
325   ierr = PetscLogFlops((91 + 10.0*ym) * xm);CHKERRQ(ierr);
326   return 0;
327 
328 }
329 
330 /*
331    FormHessian computes the quadratic term in the quadratic objective function
332    Notice that the objective function in this problem is quadratic (therefore a constant
333    hessian).  If using a nonquadratic solver, then you might want to reconsider this function
334 */
335 PetscErrorCode FormHessian(Tao tao,Vec X,Mat hes, Mat Hpre, void *ptr)
336 {
337   AppCtx*        user=(AppCtx*)ptr;
338   PetscErrorCode ierr;
339   PetscInt       i,j,k;
340   PetscInt       col[5],row,nx,ny,xs,xm,gxs,gxm,ys,ym,gys,gym;
341   PetscReal      one=1.0, two=2.0, six=6.0,pi=4.0*atan(1.0);
342   PetscReal      hx,hy,hxhy,hxhx,hyhy;
343   PetscReal      xi,v[5];
344   PetscReal      ecc=user->ecc, trule1,trule2,trule3,trule4,trule5,trule6;
345   PetscReal      vmiddle, vup, vdown, vleft, vright;
346   PetscBool      assembled;
347 
348   nx=user->nx;
349   ny=user->ny;
350   hx=two*pi/(nx+1.0);
351   hy=two*user->b/(ny+1.0);
352   hxhy=hx*hy;
353   hxhx=one/(hx*hx);
354   hyhy=one/(hy*hy);
355 
356   /*
357     Get local grid boundaries
358   */
359   ierr = DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr);
360   ierr = DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL);CHKERRQ(ierr);
361   ierr = MatAssembled(hes,&assembled);CHKERRQ(ierr);
362   if (assembled) {ierr = MatZeroEntries(hes);CHKERRQ(ierr);}
363 
364   for (i=xs; i< xs+xm; i++) {
365     xi=(i+1)*hx;
366     trule1=hxhy*(p(xi,ecc) + p(xi+hx,ecc) + p(xi,ecc)) / six; /* L(i,j) */
367     trule2=hxhy*(p(xi,ecc) + p(xi-hx,ecc) + p(xi,ecc)) / six; /* U(i,j) */
368     trule3=hxhy*(p(xi,ecc) + p(xi+hx,ecc) + p(xi+hx,ecc)) / six; /* U(i+1,j) */
369     trule4=hxhy*(p(xi,ecc) + p(xi-hx,ecc) + p(xi-hx,ecc)) / six; /* L(i-1,j) */
370     trule5=trule1; /* L(i,j-1) */
371     trule6=trule2; /* U(i,j+1) */
372 
373     vdown=-(trule5+trule2)*hyhy;
374     vleft=-hxhx*(trule2+trule4);
375     vright= -hxhx*(trule1+trule3);
376     vup=-hyhy*(trule1+trule6);
377     vmiddle=(hxhx)*(trule1+trule2+trule3+trule4)+hyhy*(trule1+trule2+trule5+trule6);
378     v[0]=0; v[1]=0; v[2]=0; v[3]=0; v[4]=0;
379 
380     for (j=ys; j<ys+ym; j++) {
381       row=(j-gys)*gxm + (i-gxs);
382 
383       k=0;
384       if (j>gys) {
385         v[k]=vdown; col[k]=row - gxm; k++;
386       }
387 
388       if (i>gxs) {
389         v[k]= vleft; col[k]=row - 1; k++;
390       }
391 
392       v[k]= vmiddle; col[k]=row; k++;
393 
394       if (i+1 < gxs+gxm) {
395         v[k]= vright; col[k]=row+1; k++;
396       }
397 
398       if (j+1 <gys+gym) {
399         v[k]= vup; col[k] = row+gxm; k++;
400       }
401       ierr = MatSetValuesLocal(hes,1,&row,k,col,v,INSERT_VALUES);CHKERRQ(ierr);
402 
403     }
404 
405   }
406 
407   /*
408      Assemble matrix, using the 2-step process:
409      MatAssemblyBegin(), MatAssemblyEnd().
410      By placing code between these two statements, computations can be
411      done while messages are in transition.
412   */
413   ierr = MatAssemblyBegin(hes,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
414   ierr = MatAssemblyEnd(hes,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
415 
416   /*
417     Tell the matrix we will never add a new nonzero location to the
418     matrix. If we do it will generate an error.
419   */
420   ierr = MatSetOption(hes,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
421   ierr = MatSetOption(hes,MAT_SYMMETRIC,PETSC_TRUE);CHKERRQ(ierr);
422 
423   ierr = PetscLogFlops(9.0*xm*ym+49.0*xm);CHKERRQ(ierr);
424   ierr = MatNorm(hes,NORM_1,&hx);CHKERRQ(ierr);
425   return 0;
426 }
427 
428 PetscErrorCode Monitor(Tao tao, void *ctx)
429 {
430   PetscErrorCode     ierr;
431   PetscInt           its;
432   PetscReal          f,gnorm,cnorm,xdiff;
433   TaoConvergedReason reason;
434 
435   PetscFunctionBegin;
436   ierr = TaoGetSolutionStatus(tao, &its, &f, &gnorm, &cnorm, &xdiff, &reason);CHKERRQ(ierr);
437   if (!(its%5)) {
438     ierr = PetscPrintf(PETSC_COMM_WORLD,"iteration=%D\tf=%g\n",its,(double)f);CHKERRQ(ierr);
439   }
440   PetscFunctionReturn(0);
441 }
442 
443 PetscErrorCode ConvergenceTest(Tao tao, void *ctx)
444 {
445   PetscErrorCode     ierr;
446   PetscInt           its;
447   PetscReal          f,gnorm,cnorm,xdiff;
448   TaoConvergedReason reason;
449 
450   PetscFunctionBegin;
451   ierr = TaoGetSolutionStatus(tao, &its, &f, &gnorm, &cnorm, &xdiff, &reason);CHKERRQ(ierr);
452   if (its == 100) {
453     ierr = TaoSetConvergedReason(tao,TAO_DIVERGED_MAXITS);CHKERRQ(ierr);
454   }
455   PetscFunctionReturn(0);
456 
457 }
458 
459 /*TEST
460 
461    build:
462       requires: !complex
463 
464    test:
465       args: -tao_smonitor -mx 8 -my 12 -tao_type tron -tao_gatol 1.e-5
466       requires: !single
467 
468    test:
469       suffix: 2
470       nsize: 2
471       args: -tao_smonitor -mx 50 -my 50 -ecc 0.99 -tao_type gpcg -tao_gatol 1.e-5
472       requires: !single
473 
474    test:
475       suffix: 3
476       nsize: 2
477       args: -tao_smonitor -mx 10 -my 16 -ecc 0.9 -tao_type bqpip -tao_gatol 1.e-4
478       requires: !single
479 
480    test:
481       suffix: 4
482       nsize: 2
483       args: -tao_smonitor -mx 10 -my 16 -ecc 0.9 -tao_type bqpip -tao_gatol 1.e-4 -test_getdiagonal
484       output_file: output/jbearing2_3.out
485       requires: !single
486 
487    test:
488       suffix: 5
489       args: -tao_smonitor -mx 8 -my 12 -tao_type bncg -tao_bncg_type gd -tao_gatol 1e-4
490       requires: !single
491 
492    test:
493       suffix: 6
494       args: -tao_smonitor -mx 8 -my 12 -tao_type bncg -tao_gatol 1e-4
495       requires: !single
496 
497    test:
498       suffix: 7
499       args: -tao_smonitor -mx 8 -my 12 -tao_type bnls -tao_gatol 1e-5
500       requires: !single
501 
502    test:
503       suffix: 8
504       args: -tao_smonitor -mx 8 -my 12 -tao_type bntr -tao_gatol 1e-5
505       requires: !single
506 
507    test:
508       suffix: 9
509       args: -tao_smonitor -mx 8 -my 12 -tao_type bntl -tao_gatol 1e-5
510       requires: !single
511 
512    test:
513       suffix: 10
514       args: -tao_smonitor -mx 8 -my 12 -tao_type bnls -tao_gatol 1e-5 -tao_bnk_max_cg_its 3
515       requires: !single
516 
517    test:
518       suffix: 11
519       args: -tao_smonitor -mx 8 -my 12 -tao_type bntr -tao_gatol 1e-5 -tao_bnk_max_cg_its 3
520       requires: !single
521 
522    test:
523       suffix: 12
524       args: -tao_smonitor -mx 8 -my 12 -tao_type bntl -tao_gatol 1e-5 -tao_bnk_max_cg_its 3
525       requires: !single
526 
527    test:
528      suffix: 13
529      args: -tao_smonitor -mx 8 -my 12 -tao_gatol 1e-4 -tao_type bqnls
530      requires: !single
531 
532    test:
533      suffix: 14
534      args: -tao_smonitor -mx 8 -my 12 -tao_gatol 1e-4 -tao_type blmvm
535      requires: !single
536 
537    test:
538      suffix: 15
539      args: -tao_smonitor -mx 8 -my 12 -tao_gatol 1e-4 -tao_type bqnkls -tao_bqnk_mat_type lmvmbfgs
540      requires: !single
541 
542    test:
543      suffix: 16
544      args: -tao_smonitor -mx 8 -my 12 -tao_gatol 1e-4 -tao_type bqnktr -tao_bqnk_mat_type lmvmsr1
545      requires: !single
546 
547    test:
548      suffix: 17
549      args: -tao_smonitor -mx 8 -my 12 -tao_gatol 1e-4 -tao_type bqnls -tao_bqnls_mat_lmvm_scale_type scalar -tao_view
550      requires: !single
551 
552    test:
553      suffix: 18
554      args: -tao_smonitor -mx 8 -my 12 -tao_gatol 1e-4 -tao_type bqnls -tao_bqnls_mat_lmvm_scale_type none -tao_view
555      requires: !single
556 
557 TEST*/
558