xref: /petsc/src/tao/bound/tutorials/jbearing2.c (revision d21efd2e5d911db017a545648c4fa4838359bb2d)
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   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 = TaoSetSolution(tao,x);CHKERRQ(ierr);
138 
139   /* Set the user function, gradient, hessian evaluation routines and data structures */
140   ierr = TaoSetObjectiveAndGradient(tao,NULL,FormFunctionGradient,(void*) &user);CHKERRQ(ierr);
141 
142   ierr = TaoSetHessian(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   PetscFunctionBegin;
202   nx=user->nx;
203   ny=user->ny;
204   hx=two*pi/(nx+1.0);
205   hy=two*user->b/(ny+1.0);
206   ehxhy = ecc*hx*hy;
207 
208   /*
209      Get local grid boundaries
210   */
211   ierr = DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr);
212   ierr = DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL);CHKERRQ(ierr);
213 
214   /* Compute the linear term in the objective function */
215   ierr = VecGetArray(user->B,&b);CHKERRQ(ierr);
216   for (i=xs; i<xs+xm; i++) {
217     temp=PetscSinScalar((i+1)*hx);
218     for (j=ys; j<ys+ym; j++) {
219       k=xm*(j-ys)+(i-xs);
220       b[k]=  - ehxhy*temp;
221     }
222   }
223   ierr = VecRestoreArray(user->B,&b);CHKERRQ(ierr);
224   ierr = PetscLogFlops(5.0*xm*ym+3.0*xm);CHKERRQ(ierr);
225   PetscFunctionReturn(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   PetscFunctionBegin;
244   nx=user->nx;
245   ny=user->ny;
246   hx=two*pi/(nx+1.0);
247   hy=two*user->b/(ny+1.0);
248   hxhy=hx*hy;
249   hxhx=one/(hx*hx);
250   hyhy=one/(hy*hy);
251 
252   ierr = DMGetLocalVector(user->dm,&localX);CHKERRQ(ierr);
253 
254   ierr = DMGlobalToLocalBegin(user->dm,X,INSERT_VALUES,localX);CHKERRQ(ierr);
255   ierr = DMGlobalToLocalEnd(user->dm,X,INSERT_VALUES,localX);CHKERRQ(ierr);
256 
257   ierr = VecSet(G, zero);CHKERRQ(ierr);
258   /*
259     Get local grid boundaries
260   */
261   ierr = DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr);
262   ierr = DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL);CHKERRQ(ierr);
263 
264   ierr = VecGetArray(localX,&x);CHKERRQ(ierr);
265   ierr = VecGetArray(G,&g);CHKERRQ(ierr);
266 
267   for (i=xs; i< xs+xm; i++) {
268     xi=(i+1)*hx;
269     trule1=hxhy*(p(xi,ecc) + p(xi+hx,ecc) + p(xi,ecc)) / six; /* L(i,j) */
270     trule2=hxhy*(p(xi,ecc) + p(xi-hx,ecc) + p(xi,ecc)) / six; /* U(i,j) */
271     trule3=hxhy*(p(xi,ecc) + p(xi+hx,ecc) + p(xi+hx,ecc)) / six; /* U(i+1,j) */
272     trule4=hxhy*(p(xi,ecc) + p(xi-hx,ecc) + p(xi-hx,ecc)) / six; /* L(i-1,j) */
273     trule5=trule1; /* L(i,j-1) */
274     trule6=trule2; /* U(i,j+1) */
275 
276     vdown=-(trule5+trule2)*hyhy;
277     vleft=-hxhx*(trule2+trule4);
278     vright= -hxhx*(trule1+trule3);
279     vup=-hyhy*(trule1+trule6);
280     vmiddle=(hxhx)*(trule1+trule2+trule3+trule4)+hyhy*(trule1+trule2+trule5+trule6);
281 
282     for (j=ys; j<ys+ym; j++) {
283 
284       row=(j-gys)*gxm + (i-gxs);
285        v[0]=0; v[1]=0; v[2]=0; v[3]=0; v[4]=0;
286 
287        k=0;
288        if (j>gys) {
289          v[k]=vdown; col[k]=row - gxm; k++;
290        }
291 
292        if (i>gxs) {
293          v[k]= vleft; col[k]=row - 1; k++;
294        }
295 
296        v[k]= vmiddle; col[k]=row; k++;
297 
298        if (i+1 < gxs+gxm) {
299          v[k]= vright; col[k]=row+1; k++;
300        }
301 
302        if (j+1 <gys+gym) {
303          v[k]= vup; col[k] = row+gxm; k++;
304        }
305        tt=0;
306        for (kk=0;kk<k;kk++) {
307          tt+=v[kk]*x[col[kk]];
308        }
309        row=(j-ys)*xm + (i-xs);
310        g[row]=tt;
311 
312      }
313 
314   }
315 
316   ierr = VecRestoreArray(localX,&x);CHKERRQ(ierr);
317   ierr = VecRestoreArray(G,&g);CHKERRQ(ierr);
318 
319   ierr = DMRestoreLocalVector(user->dm,&localX);CHKERRQ(ierr);
320 
321   ierr = VecDot(X,G,&f1);CHKERRQ(ierr);
322   ierr = VecDot(user->B,X,&f2);CHKERRQ(ierr);
323   ierr = VecAXPY(G, one, user->B);CHKERRQ(ierr);
324   *fcn = f1/2.0 + f2;
325 
326   ierr = PetscLogFlops((91 + 10.0*ym) * xm);CHKERRQ(ierr);
327   PetscFunctionReturn(0);
328 
329 }
330 
331 /*
332    FormHessian computes the quadratic term in the quadratic objective function
333    Notice that the objective function in this problem is quadratic (therefore a constant
334    hessian).  If using a nonquadratic solver, then you might want to reconsider this function
335 */
336 PetscErrorCode FormHessian(Tao tao,Vec X,Mat hes, Mat Hpre, void *ptr)
337 {
338   AppCtx*        user=(AppCtx*)ptr;
339   PetscErrorCode ierr;
340   PetscInt       i,j,k;
341   PetscInt       col[5],row,nx,ny,xs,xm,gxs,gxm,ys,ym,gys,gym;
342   PetscReal      one=1.0, two=2.0, six=6.0,pi=4.0*atan(1.0);
343   PetscReal      hx,hy,hxhy,hxhx,hyhy;
344   PetscReal      xi,v[5];
345   PetscReal      ecc=user->ecc, trule1,trule2,trule3,trule4,trule5,trule6;
346   PetscReal      vmiddle, vup, vdown, vleft, vright;
347   PetscBool      assembled;
348 
349   PetscFunctionBegin;
350   nx=user->nx;
351   ny=user->ny;
352   hx=two*pi/(nx+1.0);
353   hy=two*user->b/(ny+1.0);
354   hxhy=hx*hy;
355   hxhx=one/(hx*hx);
356   hyhy=one/(hy*hy);
357 
358   /*
359     Get local grid boundaries
360   */
361   ierr = DMDAGetCorners(user->dm,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr);
362   ierr = DMDAGetGhostCorners(user->dm,&gxs,&gys,NULL,&gxm,&gym,NULL);CHKERRQ(ierr);
363   ierr = MatAssembled(hes,&assembled);CHKERRQ(ierr);
364   if (assembled) {ierr = MatZeroEntries(hes);CHKERRQ(ierr);}
365 
366   for (i=xs; i< xs+xm; i++) {
367     xi=(i+1)*hx;
368     trule1=hxhy*(p(xi,ecc) + p(xi+hx,ecc) + p(xi,ecc)) / six; /* L(i,j) */
369     trule2=hxhy*(p(xi,ecc) + p(xi-hx,ecc) + p(xi,ecc)) / six; /* U(i,j) */
370     trule3=hxhy*(p(xi,ecc) + p(xi+hx,ecc) + p(xi+hx,ecc)) / six; /* U(i+1,j) */
371     trule4=hxhy*(p(xi,ecc) + p(xi-hx,ecc) + p(xi-hx,ecc)) / six; /* L(i-1,j) */
372     trule5=trule1; /* L(i,j-1) */
373     trule6=trule2; /* U(i,j+1) */
374 
375     vdown=-(trule5+trule2)*hyhy;
376     vleft=-hxhx*(trule2+trule4);
377     vright= -hxhx*(trule1+trule3);
378     vup=-hyhy*(trule1+trule6);
379     vmiddle=(hxhx)*(trule1+trule2+trule3+trule4)+hyhy*(trule1+trule2+trule5+trule6);
380     v[0]=0; v[1]=0; v[2]=0; v[3]=0; v[4]=0;
381 
382     for (j=ys; j<ys+ym; j++) {
383       row=(j-gys)*gxm + (i-gxs);
384 
385       k=0;
386       if (j>gys) {
387         v[k]=vdown; col[k]=row - gxm; k++;
388       }
389 
390       if (i>gxs) {
391         v[k]= vleft; col[k]=row - 1; k++;
392       }
393 
394       v[k]= vmiddle; col[k]=row; k++;
395 
396       if (i+1 < gxs+gxm) {
397         v[k]= vright; col[k]=row+1; k++;
398       }
399 
400       if (j+1 <gys+gym) {
401         v[k]= vup; col[k] = row+gxm; k++;
402       }
403       ierr = MatSetValuesLocal(hes,1,&row,k,col,v,INSERT_VALUES);CHKERRQ(ierr);
404 
405     }
406 
407   }
408 
409   /*
410      Assemble matrix, using the 2-step process:
411      MatAssemblyBegin(), MatAssemblyEnd().
412      By placing code between these two statements, computations can be
413      done while messages are in transition.
414   */
415   ierr = MatAssemblyBegin(hes,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
416   ierr = MatAssemblyEnd(hes,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
417 
418   /*
419     Tell the matrix we will never add a new nonzero location to the
420     matrix. If we do it will generate an error.
421   */
422   ierr = MatSetOption(hes,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
423   ierr = MatSetOption(hes,MAT_SYMMETRIC,PETSC_TRUE);CHKERRQ(ierr);
424 
425   ierr = PetscLogFlops(9.0*xm*ym+49.0*xm);CHKERRQ(ierr);
426   PetscFunctionReturn(0);
427 }
428 
429 PetscErrorCode Monitor(Tao tao, void *ctx)
430 {
431   PetscErrorCode     ierr;
432   PetscInt           its;
433   PetscReal          f,gnorm,cnorm,xdiff;
434   TaoConvergedReason reason;
435 
436   PetscFunctionBegin;
437   ierr = TaoGetSolutionStatus(tao, &its, &f, &gnorm, &cnorm, &xdiff, &reason);CHKERRQ(ierr);
438   if (!(its%5)) {
439     ierr = PetscPrintf(PETSC_COMM_WORLD,"iteration=%D\tf=%g\n",its,(double)f);CHKERRQ(ierr);
440   }
441   PetscFunctionReturn(0);
442 }
443 
444 PetscErrorCode ConvergenceTest(Tao tao, void *ctx)
445 {
446   PetscErrorCode     ierr;
447   PetscInt           its;
448   PetscReal          f,gnorm,cnorm,xdiff;
449   TaoConvergedReason reason;
450 
451   PetscFunctionBegin;
452   ierr = TaoGetSolutionStatus(tao, &its, &f, &gnorm, &cnorm, &xdiff, &reason);CHKERRQ(ierr);
453   if (its == 100) {
454     ierr = TaoSetConvergedReason(tao,TAO_DIVERGED_MAXITS);CHKERRQ(ierr);
455   }
456   PetscFunctionReturn(0);
457 
458 }
459 
460 /*TEST
461 
462    build:
463       requires: !complex
464 
465    test:
466       args: -tao_smonitor -mx 8 -my 12 -tao_type tron -tao_gatol 1.e-5
467       requires: !single
468 
469    test:
470       suffix: 2
471       nsize: 2
472       args: -tao_smonitor -mx 50 -my 50 -ecc 0.99 -tao_type gpcg -tao_gatol 1.e-5
473       requires: !single
474 
475    test:
476       suffix: 3
477       nsize: 2
478       args: -tao_smonitor -mx 10 -my 16 -ecc 0.9 -tao_type bqpip -tao_gatol 1.e-4
479       requires: !single
480 
481    test:
482       suffix: 4
483       nsize: 2
484       args: -tao_smonitor -mx 10 -my 16 -ecc 0.9 -tao_type bqpip -tao_gatol 1.e-4 -test_getdiagonal
485       output_file: output/jbearing2_3.out
486       requires: !single
487 
488    test:
489       suffix: 5
490       args: -tao_smonitor -mx 8 -my 12 -tao_type bncg -tao_bncg_type gd -tao_gatol 1e-4
491       requires: !single
492 
493    test:
494       suffix: 6
495       args: -tao_smonitor -mx 8 -my 12 -tao_type bncg -tao_gatol 1e-4
496       requires: !single
497 
498    test:
499       suffix: 7
500       args: -tao_smonitor -mx 8 -my 12 -tao_type bnls -tao_gatol 1e-5
501       requires: !single
502 
503    test:
504       suffix: 8
505       args: -tao_smonitor -mx 8 -my 12 -tao_type bntr -tao_gatol 1e-5
506       requires: !single
507 
508    test:
509       suffix: 9
510       args: -tao_smonitor -mx 8 -my 12 -tao_type bntl -tao_gatol 1e-5
511       requires: !single
512 
513    test:
514       suffix: 10
515       args: -tao_smonitor -mx 8 -my 12 -tao_type bnls -tao_gatol 1e-5 -tao_bnk_max_cg_its 3
516       requires: !single
517 
518    test:
519       suffix: 11
520       args: -tao_smonitor -mx 8 -my 12 -tao_type bntr -tao_gatol 1e-5 -tao_bnk_max_cg_its 3
521       requires: !single
522 
523    test:
524       suffix: 12
525       args: -tao_smonitor -mx 8 -my 12 -tao_type bntl -tao_gatol 1e-5 -tao_bnk_max_cg_its 3
526       requires: !single
527 
528    test:
529      suffix: 13
530      args: -tao_smonitor -mx 8 -my 12 -tao_gatol 1e-4 -tao_type bqnls
531      requires: !single
532 
533    test:
534      suffix: 14
535      args: -tao_smonitor -mx 8 -my 12 -tao_gatol 1e-4 -tao_type blmvm
536      requires: !single
537 
538    test:
539      suffix: 15
540      args: -tao_smonitor -mx 8 -my 12 -tao_gatol 1e-4 -tao_type bqnkls -tao_bqnk_mat_type lmvmbfgs
541      requires: !single
542 
543    test:
544      suffix: 16
545      args: -tao_smonitor -mx 8 -my 12 -tao_gatol 1e-4 -tao_type bqnktr -tao_bqnk_mat_type lmvmsr1
546      requires: !single
547 
548    test:
549      suffix: 17
550      args: -tao_smonitor -mx 8 -my 12 -tao_gatol 1e-4 -tao_type bqnls -tao_bqnls_mat_lmvm_scale_type scalar -tao_view
551      requires: !single
552 
553    test:
554      suffix: 18
555      args: -tao_smonitor -mx 8 -my 12 -tao_gatol 1e-4 -tao_type bqnls -tao_bqnls_mat_lmvm_scale_type none -tao_view
556      requires: !single
557 
558    test:
559      suffix: 19
560      args: -tao_smonitor -mx 8 -my 12 -tao_type bnls -tao_gatol 1e-5 -tao_mf_hessian
561      requires: !single
562 
563    test:
564       suffix: 20
565       args: -tao_smonitor -mx 8 -my 12 -tao_type bntr -tao_gatol 1e-5 -tao_mf_hessian
566       requires: !single
567 
568    test:
569       suffix: 21
570       args: -tao_smonitor -mx 8 -my 12 -tao_type bntl -tao_gatol 1e-5 -tao_mf_hessian
571       requires: !single
572 TEST*/
573