xref: /petsc/src/tao/bound/tutorials/jbearing2.c (revision 40badf4fbc550ac1f60bd080eaff6de6d55b946d)
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(); TaoSetObjectiveAndGradient();
28    Routines: TaoSetHessian();
29    Routines: TaoSetVariableBounds();
30    Routines: TaoSetMonitor(); TaoSetConvergenceTest();
31    Routines: TaoSetSolution();
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   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-mx",&user.nx,&flg));
86   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-my",&user.ny,&flg));
87   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-ecc",&user.ecc,&flg));
88   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-b",&user.b,&flg));
89   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-test_getdiagonal",&testgetdiag,NULL));
90 
91   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\n---- Journal Bearing Problem SHB-----\n"));
92   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"mx: %D,  my: %D,  ecc: %g \n\n",user.nx,user.ny,(double)user.ecc));
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   CHKERRQ(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));
103   CHKERRQ(DMSetFromOptions(user.dm));
104   CHKERRQ(DMSetUp(user.dm));
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   CHKERRQ(DMCreateGlobalVector(user.dm,&x)); /* Solution */
113   CHKERRQ(VecDuplicate(x,&user.B)); /* Linear objective */
114 
115   /*  Create matrix user.A to store quadratic, Create a local ordering scheme. */
116   CHKERRQ(VecGetLocalSize(x,&m));
117   CHKERRQ(DMCreateMatrix(user.dm,&user.A));
118 
119   if (testgetdiag) {
120     CHKERRQ(MatSetOperation(user.A,MATOP_GET_DIAGONAL,NULL));
121   }
122 
123   /* User defined function -- compute linear term of quadratic */
124   CHKERRQ(ComputeB(&user));
125 
126   /* The TAO code begins here */
127 
128   /*
129      Create the optimization solver
130      Suitable methods: TAOGPCG, TAOBQPIP, TAOTRON, TAOBLMVM
131   */
132   CHKERRQ(TaoCreate(PETSC_COMM_WORLD,&tao));
133   CHKERRQ(TaoSetType(tao,TAOBLMVM));
134 
135   /* Set the initial vector */
136   CHKERRQ(VecSet(x, zero));
137   CHKERRQ(TaoSetSolution(tao,x));
138 
139   /* Set the user function, gradient, hessian evaluation routines and data structures */
140   CHKERRQ(TaoSetObjectiveAndGradient(tao,NULL,FormFunctionGradient,(void*) &user));
141 
142   CHKERRQ(TaoSetHessian(tao,user.A,user.A,FormHessian,(void*)&user));
143 
144   /* Set a routine that defines the bounds */
145   CHKERRQ(VecDuplicate(x,&xl));
146   CHKERRQ(VecDuplicate(x,&xu));
147   CHKERRQ(VecSet(xl, zero));
148   CHKERRQ(VecSet(xu, d1000));
149   CHKERRQ(TaoSetVariableBounds(tao,xl,xu));
150 
151   CHKERRQ(TaoGetKSP(tao,&ksp));
152   if (ksp) {
153     CHKERRQ(KSPSetType(ksp,KSPCG));
154   }
155 
156   CHKERRQ(PetscOptionsHasName(NULL,NULL,"-testmonitor",&flg));
157   if (flg) {
158     CHKERRQ(TaoSetMonitor(tao,Monitor,&user,NULL));
159   }
160   CHKERRQ(PetscOptionsHasName(NULL,NULL,"-testconvergence",&flg));
161   if (flg) {
162     CHKERRQ(TaoSetConvergenceTest(tao,ConvergenceTest,&user));
163   }
164 
165   /* Check for any tao command line options */
166   CHKERRQ(TaoSetFromOptions(tao));
167 
168   /* Solve the bound constrained problem */
169   CHKERRQ(TaoSolve(tao));
170 
171   /* Free PETSc data structures */
172   CHKERRQ(VecDestroy(&x));
173   CHKERRQ(VecDestroy(&xl));
174   CHKERRQ(VecDestroy(&xu));
175   CHKERRQ(MatDestroy(&user.A));
176   CHKERRQ(VecDestroy(&user.B));
177 
178   /* Free TAO data structures */
179   CHKERRQ(TaoDestroy(&tao));
180   CHKERRQ(DMDestroy(&user.dm));
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   PetscInt       i,j,k;
194   PetscInt       nx,ny,xs,xm,gxs,gxm,ys,ym,gys,gym;
195   PetscReal      two=2.0, pi=4.0*atan(1.0);
196   PetscReal      hx,hy,ehxhy;
197   PetscReal      temp,*b;
198   PetscReal      ecc=user->ecc;
199 
200   PetscFunctionBegin;
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   CHKERRQ(DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL));
211   CHKERRQ(DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL));
212 
213   /* Compute the linear term in the objective function */
214   CHKERRQ(VecGetArray(user->B,&b));
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   CHKERRQ(VecRestoreArray(user->B,&b));
223   CHKERRQ(PetscLogFlops(5.0*xm*ym+3.0*xm));
224   PetscFunctionReturn(0);
225 }
226 
227 PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *fcn,Vec G,void *ptr)
228 {
229   AppCtx*        user=(AppCtx*)ptr;
230   PetscInt       i,j,k,kk;
231   PetscInt       col[5],row,nx,ny,xs,xm,gxs,gxm,ys,ym,gys,gym;
232   PetscReal      one=1.0, two=2.0, six=6.0,pi=4.0*atan(1.0);
233   PetscReal      hx,hy,hxhy,hxhx,hyhy;
234   PetscReal      xi,v[5];
235   PetscReal      ecc=user->ecc, trule1,trule2,trule3,trule4,trule5,trule6;
236   PetscReal      vmiddle, vup, vdown, vleft, vright;
237   PetscReal      tt,f1,f2;
238   PetscReal      *x,*g,zero=0.0;
239   Vec            localX;
240 
241   PetscFunctionBegin;
242   nx=user->nx;
243   ny=user->ny;
244   hx=two*pi/(nx+1.0);
245   hy=two*user->b/(ny+1.0);
246   hxhy=hx*hy;
247   hxhx=one/(hx*hx);
248   hyhy=one/(hy*hy);
249 
250   CHKERRQ(DMGetLocalVector(user->dm,&localX));
251 
252   CHKERRQ(DMGlobalToLocalBegin(user->dm,X,INSERT_VALUES,localX));
253   CHKERRQ(DMGlobalToLocalEnd(user->dm,X,INSERT_VALUES,localX));
254 
255   CHKERRQ(VecSet(G, zero));
256   /*
257     Get local grid boundaries
258   */
259   CHKERRQ(DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL));
260   CHKERRQ(DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL));
261 
262   CHKERRQ(VecGetArray(localX,&x));
263   CHKERRQ(VecGetArray(G,&g));
264 
265   for (i=xs; i< xs+xm; i++) {
266     xi=(i+1)*hx;
267     trule1=hxhy*(p(xi,ecc) + p(xi+hx,ecc) + p(xi,ecc)) / six; /* L(i,j) */
268     trule2=hxhy*(p(xi,ecc) + p(xi-hx,ecc) + p(xi,ecc)) / six; /* U(i,j) */
269     trule3=hxhy*(p(xi,ecc) + p(xi+hx,ecc) + p(xi+hx,ecc)) / six; /* U(i+1,j) */
270     trule4=hxhy*(p(xi,ecc) + p(xi-hx,ecc) + p(xi-hx,ecc)) / six; /* L(i-1,j) */
271     trule5=trule1; /* L(i,j-1) */
272     trule6=trule2; /* U(i,j+1) */
273 
274     vdown=-(trule5+trule2)*hyhy;
275     vleft=-hxhx*(trule2+trule4);
276     vright= -hxhx*(trule1+trule3);
277     vup=-hyhy*(trule1+trule6);
278     vmiddle=(hxhx)*(trule1+trule2+trule3+trule4)+hyhy*(trule1+trule2+trule5+trule6);
279 
280     for (j=ys; j<ys+ym; j++) {
281 
282       row=(j-gys)*gxm + (i-gxs);
283        v[0]=0; v[1]=0; v[2]=0; v[3]=0; v[4]=0;
284 
285        k=0;
286        if (j>gys) {
287          v[k]=vdown; col[k]=row - gxm; k++;
288        }
289 
290        if (i>gxs) {
291          v[k]= vleft; col[k]=row - 1; k++;
292        }
293 
294        v[k]= vmiddle; col[k]=row; k++;
295 
296        if (i+1 < gxs+gxm) {
297          v[k]= vright; col[k]=row+1; k++;
298        }
299 
300        if (j+1 <gys+gym) {
301          v[k]= vup; col[k] = row+gxm; k++;
302        }
303        tt=0;
304        for (kk=0;kk<k;kk++) {
305          tt+=v[kk]*x[col[kk]];
306        }
307        row=(j-ys)*xm + (i-xs);
308        g[row]=tt;
309 
310      }
311 
312   }
313 
314   CHKERRQ(VecRestoreArray(localX,&x));
315   CHKERRQ(VecRestoreArray(G,&g));
316 
317   CHKERRQ(DMRestoreLocalVector(user->dm,&localX));
318 
319   CHKERRQ(VecDot(X,G,&f1));
320   CHKERRQ(VecDot(user->B,X,&f2));
321   CHKERRQ(VecAXPY(G, one, user->B));
322   *fcn = f1/2.0 + f2;
323 
324   CHKERRQ(PetscLogFlops((91 + 10.0*ym) * xm));
325   PetscFunctionReturn(0);
326 
327 }
328 
329 /*
330    FormHessian computes the quadratic term in the quadratic objective function
331    Notice that the objective function in this problem is quadratic (therefore a constant
332    hessian).  If using a nonquadratic solver, then you might want to reconsider this function
333 */
334 PetscErrorCode FormHessian(Tao tao,Vec X,Mat hes, Mat Hpre, void *ptr)
335 {
336   AppCtx*        user=(AppCtx*)ptr;
337   PetscInt       i,j,k;
338   PetscInt       col[5],row,nx,ny,xs,xm,gxs,gxm,ys,ym,gys,gym;
339   PetscReal      one=1.0, two=2.0, six=6.0,pi=4.0*atan(1.0);
340   PetscReal      hx,hy,hxhy,hxhx,hyhy;
341   PetscReal      xi,v[5];
342   PetscReal      ecc=user->ecc, trule1,trule2,trule3,trule4,trule5,trule6;
343   PetscReal      vmiddle, vup, vdown, vleft, vright;
344   PetscBool      assembled;
345 
346   PetscFunctionBegin;
347   nx=user->nx;
348   ny=user->ny;
349   hx=two*pi/(nx+1.0);
350   hy=two*user->b/(ny+1.0);
351   hxhy=hx*hy;
352   hxhx=one/(hx*hx);
353   hyhy=one/(hy*hy);
354 
355   /*
356     Get local grid boundaries
357   */
358   CHKERRQ(DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL));
359   CHKERRQ(DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL));
360   CHKERRQ(MatAssembled(hes,&assembled));
361   if (assembled) CHKERRQ(MatZeroEntries(hes));
362 
363   for (i=xs; i< xs+xm; i++) {
364     xi=(i+1)*hx;
365     trule1=hxhy*(p(xi,ecc) + p(xi+hx,ecc) + p(xi,ecc)) / six; /* L(i,j) */
366     trule2=hxhy*(p(xi,ecc) + p(xi-hx,ecc) + p(xi,ecc)) / six; /* U(i,j) */
367     trule3=hxhy*(p(xi,ecc) + p(xi+hx,ecc) + p(xi+hx,ecc)) / six; /* U(i+1,j) */
368     trule4=hxhy*(p(xi,ecc) + p(xi-hx,ecc) + p(xi-hx,ecc)) / six; /* L(i-1,j) */
369     trule5=trule1; /* L(i,j-1) */
370     trule6=trule2; /* U(i,j+1) */
371 
372     vdown=-(trule5+trule2)*hyhy;
373     vleft=-hxhx*(trule2+trule4);
374     vright= -hxhx*(trule1+trule3);
375     vup=-hyhy*(trule1+trule6);
376     vmiddle=(hxhx)*(trule1+trule2+trule3+trule4)+hyhy*(trule1+trule2+trule5+trule6);
377     v[0]=0; v[1]=0; v[2]=0; v[3]=0; v[4]=0;
378 
379     for (j=ys; j<ys+ym; j++) {
380       row=(j-gys)*gxm + (i-gxs);
381 
382       k=0;
383       if (j>gys) {
384         v[k]=vdown; col[k]=row - gxm; k++;
385       }
386 
387       if (i>gxs) {
388         v[k]= vleft; col[k]=row - 1; k++;
389       }
390 
391       v[k]= vmiddle; col[k]=row; k++;
392 
393       if (i+1 < gxs+gxm) {
394         v[k]= vright; col[k]=row+1; k++;
395       }
396 
397       if (j+1 <gys+gym) {
398         v[k]= vup; col[k] = row+gxm; k++;
399       }
400       CHKERRQ(MatSetValuesLocal(hes,1,&row,k,col,v,INSERT_VALUES));
401 
402     }
403 
404   }
405 
406   /*
407      Assemble matrix, using the 2-step process:
408      MatAssemblyBegin(), MatAssemblyEnd().
409      By placing code between these two statements, computations can be
410      done while messages are in transition.
411   */
412   CHKERRQ(MatAssemblyBegin(hes,MAT_FINAL_ASSEMBLY));
413   CHKERRQ(MatAssemblyEnd(hes,MAT_FINAL_ASSEMBLY));
414 
415   /*
416     Tell the matrix we will never add a new nonzero location to the
417     matrix. If we do it will generate an error.
418   */
419   CHKERRQ(MatSetOption(hes,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE));
420   CHKERRQ(MatSetOption(hes,MAT_SYMMETRIC,PETSC_TRUE));
421 
422   CHKERRQ(PetscLogFlops(9.0*xm*ym+49.0*xm));
423   PetscFunctionReturn(0);
424 }
425 
426 PetscErrorCode Monitor(Tao tao, void *ctx)
427 {
428   PetscInt           its;
429   PetscReal          f,gnorm,cnorm,xdiff;
430   TaoConvergedReason reason;
431 
432   PetscFunctionBegin;
433   CHKERRQ(TaoGetSolutionStatus(tao, &its, &f, &gnorm, &cnorm, &xdiff, &reason));
434   if (!(its%5)) {
435     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"iteration=%D\tf=%g\n",its,(double)f));
436   }
437   PetscFunctionReturn(0);
438 }
439 
440 PetscErrorCode ConvergenceTest(Tao tao, void *ctx)
441 {
442   PetscInt           its;
443   PetscReal          f,gnorm,cnorm,xdiff;
444   TaoConvergedReason reason;
445 
446   PetscFunctionBegin;
447   CHKERRQ(TaoGetSolutionStatus(tao, &its, &f, &gnorm, &cnorm, &xdiff, &reason));
448   if (its == 100) {
449     CHKERRQ(TaoSetConvergedReason(tao,TAO_DIVERGED_MAXITS));
450   }
451   PetscFunctionReturn(0);
452 
453 }
454 
455 /*TEST
456 
457    build:
458       requires: !complex
459 
460    test:
461       args: -tao_smonitor -mx 8 -my 12 -tao_type tron -tao_gatol 1.e-5
462       requires: !single
463 
464    test:
465       suffix: 2
466       nsize: 2
467       args: -tao_smonitor -mx 50 -my 50 -ecc 0.99 -tao_type gpcg -tao_gatol 1.e-5
468       requires: !single
469 
470    test:
471       suffix: 3
472       nsize: 2
473       args: -tao_smonitor -mx 10 -my 16 -ecc 0.9 -tao_type bqpip -tao_gatol 1.e-4
474       requires: !single
475 
476    test:
477       suffix: 4
478       nsize: 2
479       args: -tao_smonitor -mx 10 -my 16 -ecc 0.9 -tao_type bqpip -tao_gatol 1.e-4 -test_getdiagonal
480       output_file: output/jbearing2_3.out
481       requires: !single
482 
483    test:
484       suffix: 5
485       args: -tao_smonitor -mx 8 -my 12 -tao_type bncg -tao_bncg_type gd -tao_gatol 1e-4
486       requires: !single
487 
488    test:
489       suffix: 6
490       args: -tao_smonitor -mx 8 -my 12 -tao_type bncg -tao_gatol 1e-4
491       requires: !single
492 
493    test:
494       suffix: 7
495       args: -tao_smonitor -mx 8 -my 12 -tao_type bnls -tao_gatol 1e-5
496       requires: !single
497 
498    test:
499       suffix: 8
500       args: -tao_smonitor -mx 8 -my 12 -tao_type bntr -tao_gatol 1e-5
501       requires: !single
502 
503    test:
504       suffix: 9
505       args: -tao_smonitor -mx 8 -my 12 -tao_type bntl -tao_gatol 1e-5
506       requires: !single
507 
508    test:
509       suffix: 10
510       args: -tao_smonitor -mx 8 -my 12 -tao_type bnls -tao_gatol 1e-5 -tao_bnk_max_cg_its 3
511       requires: !single
512 
513    test:
514       suffix: 11
515       args: -tao_smonitor -mx 8 -my 12 -tao_type bntr -tao_gatol 1e-5 -tao_bnk_max_cg_its 3
516       requires: !single
517 
518    test:
519       suffix: 12
520       args: -tao_smonitor -mx 8 -my 12 -tao_type bntl -tao_gatol 1e-5 -tao_bnk_max_cg_its 3
521       requires: !single
522 
523    test:
524      suffix: 13
525      args: -tao_smonitor -mx 8 -my 12 -tao_gatol 1e-4 -tao_type bqnls
526      requires: !single
527 
528    test:
529      suffix: 14
530      args: -tao_smonitor -mx 8 -my 12 -tao_gatol 1e-4 -tao_type blmvm
531      requires: !single
532 
533    test:
534      suffix: 15
535      args: -tao_smonitor -mx 8 -my 12 -tao_gatol 1e-4 -tao_type bqnkls -tao_bqnk_mat_type lmvmbfgs
536      requires: !single
537 
538    test:
539      suffix: 16
540      args: -tao_smonitor -mx 8 -my 12 -tao_gatol 1e-4 -tao_type bqnktr -tao_bqnk_mat_type lmvmsr1
541      requires: !single
542 
543    test:
544      suffix: 17
545      args: -tao_smonitor -mx 8 -my 12 -tao_gatol 1e-4 -tao_type bqnls -tao_bqnls_mat_lmvm_scale_type scalar -tao_view
546      requires: !single
547 
548    test:
549      suffix: 18
550      args: -tao_smonitor -mx 8 -my 12 -tao_gatol 1e-4 -tao_type bqnls -tao_bqnls_mat_lmvm_scale_type none -tao_view
551      requires: !single
552 
553    test:
554      suffix: 19
555      args: -tao_smonitor -mx 8 -my 12 -tao_type bnls -tao_gatol 1e-5 -tao_mf_hessian
556      requires: !single
557 
558    test:
559       suffix: 20
560       args: -tao_smonitor -mx 8 -my 12 -tao_type bntr -tao_gatol 1e-5 -tao_mf_hessian
561       requires: !single
562 
563    test:
564       suffix: 21
565       args: -tao_smonitor -mx 8 -my 12 -tao_type bntl -tao_gatol 1e-5 -tao_mf_hessian
566       requires: !single
567 TEST*/
568