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