xref: /petsc/src/tao/unconstrained/tutorials/rosenbrock1.c (revision ca0c957dd390e0254b8e0ddfb1485a321c0a7f54)
1c4762a1bSJed Brown /* Program usage: mpiexec -n 1 rosenbrock1 [-help] [all TAO options] */
2c4762a1bSJed Brown 
3c4762a1bSJed Brown /*  Include "petsctao.h" so we can use TAO solvers.  */
4c4762a1bSJed Brown #include <petsctao.h>
5c4762a1bSJed Brown 
6c4762a1bSJed Brown static  char help[] = "This example demonstrates use of the TAO package to \n\
7c4762a1bSJed Brown solve an unconstrained minimization problem on a single processor.  We \n\
8c4762a1bSJed Brown minimize the extended Rosenbrock function: \n\
9c4762a1bSJed Brown    sum_{i=0}^{n/2-1} ( alpha*(x_{2i+1}-x_{2i}^2)^2 + (1-x_{2i})^2 ) \n\
10c4762a1bSJed Brown or the chained Rosenbrock function:\n\
11c4762a1bSJed Brown    sum_{i=0}^{n-1} alpha*(x_{i+1} - x_i^2)^2 + (1 - x_i)^2\n";
12c4762a1bSJed Brown 
13c4762a1bSJed Brown /*T
14c4762a1bSJed Brown    Concepts: TAO^Solving an unconstrained minimization problem
15c4762a1bSJed Brown    Routines: TaoCreate();
16c4762a1bSJed Brown    Routines: TaoSetType(); TaoSetObjectiveAndGradientRoutine();
17c4762a1bSJed Brown    Routines: TaoSetHessianRoutine();
18c4762a1bSJed Brown    Routines: TaoSetInitialVector();
19c4762a1bSJed Brown    Routines: TaoSetFromOptions();
20c4762a1bSJed Brown    Routines: TaoSolve();
21c4762a1bSJed Brown    Routines: TaoDestroy();
22c4762a1bSJed Brown    Processors: 1
23c4762a1bSJed Brown T*/
24c4762a1bSJed Brown 
25c4762a1bSJed Brown 
26c4762a1bSJed Brown 
27c4762a1bSJed Brown 
28c4762a1bSJed Brown /*
29c4762a1bSJed Brown    User-defined application context - contains data needed by the
30c4762a1bSJed Brown    application-provided call-back routines that evaluate the function,
31c4762a1bSJed Brown    gradient, and hessian.
32c4762a1bSJed Brown */
33c4762a1bSJed Brown typedef struct {
34c4762a1bSJed Brown   PetscInt  n;          /* dimension */
35c4762a1bSJed Brown   PetscReal alpha;   /* condition parameter */
36c4762a1bSJed Brown   PetscBool chained;
37c4762a1bSJed Brown } AppCtx;
38c4762a1bSJed Brown 
39c4762a1bSJed Brown /* -------------- User-defined routines ---------- */
40c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao,Vec,PetscReal*,Vec,void*);
41c4762a1bSJed Brown PetscErrorCode FormHessian(Tao,Vec,Mat,Mat,void*);
42c4762a1bSJed Brown 
43c4762a1bSJed Brown int main(int argc,char **argv)
44c4762a1bSJed Brown {
45c4762a1bSJed Brown   PetscErrorCode     ierr;                  /* used to check for functions returning nonzeros */
46c4762a1bSJed Brown   PetscReal          zero=0.0;
47c4762a1bSJed Brown   Vec                x;                     /* solution vector */
48c4762a1bSJed Brown   Mat                H;
49c4762a1bSJed Brown   Tao                tao;                   /* Tao solver context */
50c4762a1bSJed Brown   PetscBool          flg, test_lmvm = PETSC_FALSE;
51c4762a1bSJed Brown   PetscMPIInt        size;                  /* number of processes running */
52c4762a1bSJed Brown   AppCtx             user;                  /* user-defined application context */
53c4762a1bSJed Brown   KSP                ksp;
54c4762a1bSJed Brown   PC                 pc;
55c4762a1bSJed Brown   Mat                M;
56c4762a1bSJed Brown   Vec                in, out, out2;
57c4762a1bSJed Brown   PetscReal          mult_solve_dist;
58c4762a1bSJed Brown 
59c4762a1bSJed Brown   /* Initialize TAO and PETSc */
60c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
61c4762a1bSJed Brown   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
62c4762a1bSJed Brown   if (size >1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Incorrect number of processors");
63c4762a1bSJed Brown 
64c4762a1bSJed Brown   /* Initialize problem parameters */
65c4762a1bSJed Brown   user.n = 2; user.alpha = 99.0; user.chained = PETSC_FALSE;
66c4762a1bSJed Brown   /* Check for command line arguments to override defaults */
67c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-n",&user.n,&flg);CHKERRQ(ierr);
68c4762a1bSJed Brown   ierr = PetscOptionsGetReal(NULL,NULL,"-alpha",&user.alpha,&flg);CHKERRQ(ierr);
69c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-chained",&user.chained,&flg);CHKERRQ(ierr);
70c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-test_lmvm",&test_lmvm,&flg);CHKERRQ(ierr);
71c4762a1bSJed Brown 
72c4762a1bSJed Brown   /* Allocate vectors for the solution and gradient */
73c4762a1bSJed Brown   ierr = VecCreateSeq(PETSC_COMM_SELF,user.n,&x);CHKERRQ(ierr);
74c4762a1bSJed Brown   ierr = MatCreateSeqBAIJ(PETSC_COMM_SELF,2,user.n,user.n,1,NULL,&H);CHKERRQ(ierr);
75c4762a1bSJed Brown 
76c4762a1bSJed Brown   /* The TAO code begins here */
77c4762a1bSJed Brown 
78c4762a1bSJed Brown   /* Create TAO solver with desired solution method */
79c4762a1bSJed Brown   ierr = TaoCreate(PETSC_COMM_SELF,&tao);CHKERRQ(ierr);
80c4762a1bSJed Brown   ierr = TaoSetType(tao,TAOLMVM);CHKERRQ(ierr);
81c4762a1bSJed Brown 
82c4762a1bSJed Brown   /* Set solution vec and an initial guess */
83c4762a1bSJed Brown   ierr = VecSet(x, zero);CHKERRQ(ierr);
84c4762a1bSJed Brown   ierr = TaoSetInitialVector(tao,x);CHKERRQ(ierr);
85c4762a1bSJed Brown 
86c4762a1bSJed Brown   /* Set routines for function, gradient, hessian evaluation */
87c4762a1bSJed Brown   ierr = TaoSetObjectiveAndGradientRoutine(tao,FormFunctionGradient,&user);CHKERRQ(ierr);
88c4762a1bSJed Brown   ierr = TaoSetHessianRoutine(tao,H,H,FormHessian,&user);CHKERRQ(ierr);
89c4762a1bSJed Brown 
90c4762a1bSJed Brown   /* Test the LMVM matrix */
91c4762a1bSJed Brown   if (test_lmvm) {
92c4762a1bSJed Brown     ierr = PetscOptionsSetValue(NULL, "-tao_type", "bqnktr");CHKERRQ(ierr);
93c4762a1bSJed Brown   }
94c4762a1bSJed Brown 
95c4762a1bSJed Brown   /* Check for TAO command line options */
96c4762a1bSJed Brown   ierr = TaoSetFromOptions(tao);CHKERRQ(ierr);
97c4762a1bSJed Brown 
98c4762a1bSJed Brown   /* SOLVE THE APPLICATION */
99c4762a1bSJed Brown   ierr = TaoSolve(tao);CHKERRQ(ierr);
100c4762a1bSJed Brown 
101c4762a1bSJed Brown   /* Test the LMVM matrix */
102c4762a1bSJed Brown   if (test_lmvm) {
103c4762a1bSJed Brown     ierr = TaoGetKSP(tao, &ksp);CHKERRQ(ierr);
104c4762a1bSJed Brown     ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
105c4762a1bSJed Brown     ierr = PCLMVMGetMatLMVM(pc, &M);CHKERRQ(ierr);
106c4762a1bSJed Brown     ierr = VecDuplicate(x, &in);CHKERRQ(ierr);
107c4762a1bSJed Brown     ierr = VecDuplicate(x, &out);CHKERRQ(ierr);
108c4762a1bSJed Brown     ierr = VecDuplicate(x, &out2);CHKERRQ(ierr);
109c4762a1bSJed Brown     ierr = VecSet(in, 1.0);CHKERRQ(ierr);
110c4762a1bSJed Brown     ierr = MatMult(M, in, out);CHKERRQ(ierr);
111c4762a1bSJed Brown     ierr = MatSolve(M, out, out2);CHKERRQ(ierr);
112c4762a1bSJed Brown     ierr = VecAXPY(out2, -1.0, in);CHKERRQ(ierr);
113c4762a1bSJed Brown     ierr = VecNorm(out2, NORM_2, &mult_solve_dist);CHKERRQ(ierr);
114c4762a1bSJed Brown     if (mult_solve_dist < 1.e-11) {
115c4762a1bSJed Brown       ierr = PetscPrintf(PetscObjectComm((PetscObject)tao), "error between LMVM MatMult and MatSolve: < 1.e-11\n");CHKERRQ(ierr);
116c4762a1bSJed Brown     } else if(mult_solve_dist < 1.e-6) {
117c4762a1bSJed Brown       ierr = PetscPrintf(PetscObjectComm((PetscObject)tao), "error between LMVM MatMult and MatSolve: < 1.e-6\n");CHKERRQ(ierr);
118c4762a1bSJed Brown     } else {
119c4762a1bSJed Brown       ierr = PetscPrintf(PetscObjectComm((PetscObject)tao), "error between LMVM MatMult and MatSolve: %e\n", (double)mult_solve_dist);CHKERRQ(ierr);
120c4762a1bSJed Brown     }
121c4762a1bSJed Brown     ierr = VecDestroy(&in);CHKERRQ(ierr);
122c4762a1bSJed Brown     ierr = VecDestroy(&out);CHKERRQ(ierr);
123c4762a1bSJed Brown     ierr = VecDestroy(&out2);CHKERRQ(ierr);
124c4762a1bSJed Brown   }
125c4762a1bSJed Brown 
126c4762a1bSJed Brown   ierr = TaoDestroy(&tao);CHKERRQ(ierr);
127c4762a1bSJed Brown   ierr = VecDestroy(&x);CHKERRQ(ierr);
128c4762a1bSJed Brown   ierr = MatDestroy(&H);CHKERRQ(ierr);
129c4762a1bSJed Brown 
130c4762a1bSJed Brown   ierr = PetscFinalize();
131c4762a1bSJed Brown   return ierr;
132c4762a1bSJed Brown }
133c4762a1bSJed Brown 
134c4762a1bSJed Brown /* -------------------------------------------------------------------- */
135c4762a1bSJed Brown /*
136c4762a1bSJed Brown     FormFunctionGradient - Evaluates the function, f(X), and gradient, G(X).
137c4762a1bSJed Brown 
138c4762a1bSJed Brown     Input Parameters:
139c4762a1bSJed Brown .   tao  - the Tao context
140c4762a1bSJed Brown .   X    - input vector
141c4762a1bSJed Brown .   ptr  - optional user-defined context, as set by TaoSetFunctionGradient()
142c4762a1bSJed Brown 
143c4762a1bSJed Brown     Output Parameters:
144c4762a1bSJed Brown .   G - vector containing the newly evaluated gradient
145c4762a1bSJed Brown .   f - function value
146c4762a1bSJed Brown 
147c4762a1bSJed Brown     Note:
148c4762a1bSJed Brown     Some optimization methods ask for the function and the gradient evaluation
149c4762a1bSJed Brown     at the same time.  Evaluating both at once may be more efficient that
150c4762a1bSJed Brown     evaluating each separately.
151c4762a1bSJed Brown */
152c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao tao,Vec X,PetscReal *f, Vec G,void *ptr)
153c4762a1bSJed Brown {
154c4762a1bSJed Brown   AppCtx            *user = (AppCtx *) ptr;
155c4762a1bSJed Brown   PetscInt          i,nn=user->n/2;
156c4762a1bSJed Brown   PetscErrorCode    ierr;
157c4762a1bSJed Brown   PetscReal         ff=0,t1,t2,alpha=user->alpha;
158c4762a1bSJed Brown   PetscScalar       *g;
159c4762a1bSJed Brown   const PetscScalar *x;
160c4762a1bSJed Brown 
161c4762a1bSJed Brown   PetscFunctionBeginUser;
162c4762a1bSJed Brown   /* Get pointers to vector data */
163c4762a1bSJed Brown   ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
164c4762a1bSJed Brown   ierr = VecGetArray(G,&g);CHKERRQ(ierr);
165c4762a1bSJed Brown 
166c4762a1bSJed Brown   /* Compute G(X) */
167c4762a1bSJed Brown   if (user->chained) {
168c4762a1bSJed Brown     g[0] = 0;
169c4762a1bSJed Brown     for (i=0; i<user->n-1; i++) {
170c4762a1bSJed Brown       t1 = x[i+1] - x[i]*x[i];
171c4762a1bSJed Brown       ff += PetscSqr(1 - x[i]) + alpha*t1*t1;
172c4762a1bSJed Brown       g[i] += -2*(1 - x[i]) + 2*alpha*t1*(-2*x[i]);
173c4762a1bSJed Brown       g[i+1] = 2*alpha*t1;
174c4762a1bSJed Brown     }
175c4762a1bSJed Brown   } else {
176c4762a1bSJed Brown     for (i=0; i<nn; i++){
177c4762a1bSJed Brown       t1 = x[2*i+1]-x[2*i]*x[2*i]; t2= 1-x[2*i];
178c4762a1bSJed Brown       ff += alpha*t1*t1 + t2*t2;
179c4762a1bSJed Brown       g[2*i] = -4*alpha*t1*x[2*i]-2.0*t2;
180c4762a1bSJed Brown       g[2*i+1] = 2*alpha*t1;
181c4762a1bSJed Brown     }
182c4762a1bSJed Brown   }
183c4762a1bSJed Brown 
184c4762a1bSJed Brown   /* Restore vectors */
185c4762a1bSJed Brown   ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
186c4762a1bSJed Brown   ierr = VecRestoreArray(G,&g);CHKERRQ(ierr);
187c4762a1bSJed Brown   *f   = ff;
188c4762a1bSJed Brown 
189*ca0c957dSBarry Smith   ierr = PetscLogFlops(15.0*nn);CHKERRQ(ierr);
190c4762a1bSJed Brown   PetscFunctionReturn(0);
191c4762a1bSJed Brown }
192c4762a1bSJed Brown 
193c4762a1bSJed Brown /* ------------------------------------------------------------------- */
194c4762a1bSJed Brown /*
195c4762a1bSJed Brown    FormHessian - Evaluates Hessian matrix.
196c4762a1bSJed Brown 
197c4762a1bSJed Brown    Input Parameters:
198c4762a1bSJed Brown .  tao   - the Tao context
199c4762a1bSJed Brown .  x     - input vector
200c4762a1bSJed Brown .  ptr   - optional user-defined context, as set by TaoSetHessian()
201c4762a1bSJed Brown 
202c4762a1bSJed Brown    Output Parameters:
203c4762a1bSJed Brown .  H     - Hessian matrix
204c4762a1bSJed Brown 
205c4762a1bSJed Brown    Note:  Providing the Hessian may not be necessary.  Only some solvers
206c4762a1bSJed Brown    require this matrix.
207c4762a1bSJed Brown */
208c4762a1bSJed Brown PetscErrorCode FormHessian(Tao tao,Vec X,Mat H, Mat Hpre, void *ptr)
209c4762a1bSJed Brown {
210c4762a1bSJed Brown   AppCtx            *user = (AppCtx*)ptr;
211c4762a1bSJed Brown   PetscErrorCode    ierr;
212c4762a1bSJed Brown   PetscInt          i, ind[2];
213c4762a1bSJed Brown   PetscReal         alpha=user->alpha;
214c4762a1bSJed Brown   PetscReal         v[2][2];
215c4762a1bSJed Brown   const PetscScalar *x;
216c4762a1bSJed Brown   PetscBool         assembled;
217c4762a1bSJed Brown 
218c4762a1bSJed Brown   PetscFunctionBeginUser;
219c4762a1bSJed Brown   /* Zero existing matrix entries */
220c4762a1bSJed Brown   ierr = MatAssembled(H,&assembled);CHKERRQ(ierr);
221c4762a1bSJed Brown   if (assembled){ierr = MatZeroEntries(H); CHKERRQ(ierr);}
222c4762a1bSJed Brown 
223c4762a1bSJed Brown   /* Get a pointer to vector data */
224c4762a1bSJed Brown   ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
225c4762a1bSJed Brown 
226c4762a1bSJed Brown   /* Compute H(X) entries */
227c4762a1bSJed Brown   if (user->chained) {
228c4762a1bSJed Brown     ierr = MatZeroEntries(H);CHKERRQ(ierr);
229c4762a1bSJed Brown     for (i=0; i<user->n-1; i++) {
230c4762a1bSJed Brown       PetscScalar t1 = x[i+1] - x[i]*x[i];
231c4762a1bSJed Brown       v[0][0] = 2 + 2*alpha*(t1*(-2) - 2*x[i]);
232c4762a1bSJed Brown       v[0][1] = 2*alpha*(-2*x[i]);
233c4762a1bSJed Brown       v[1][0] = 2*alpha*(-2*x[i]);
234c4762a1bSJed Brown       v[1][1] = 2*alpha*t1;
235c4762a1bSJed Brown       ind[0] = i; ind[1] = i+1;
236c4762a1bSJed Brown       ierr = MatSetValues(H,2,ind,2,ind,v[0],ADD_VALUES);CHKERRQ(ierr);
237c4762a1bSJed Brown     }
238c4762a1bSJed Brown   } else {
239c4762a1bSJed Brown     for (i=0; i<user->n/2; i++){
240c4762a1bSJed Brown       v[1][1] = 2*alpha;
241c4762a1bSJed Brown       v[0][0] = -4*alpha*(x[2*i+1]-3*x[2*i]*x[2*i]) + 2;
242c4762a1bSJed Brown       v[1][0] = v[0][1] = -4.0*alpha*x[2*i];
243c4762a1bSJed Brown       ind[0]=2*i; ind[1]=2*i+1;
244c4762a1bSJed Brown       ierr = MatSetValues(H,2,ind,2,ind,v[0],INSERT_VALUES);CHKERRQ(ierr);
245c4762a1bSJed Brown     }
246c4762a1bSJed Brown   }
247c4762a1bSJed Brown   ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
248c4762a1bSJed Brown 
249c4762a1bSJed Brown   /* Assemble matrix */
250c4762a1bSJed Brown   ierr = MatAssemblyBegin(H,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
251c4762a1bSJed Brown   ierr = MatAssemblyEnd(H,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
252c4762a1bSJed Brown   ierr = PetscLogFlops(9.0*user->n/2.0);CHKERRQ(ierr);
253c4762a1bSJed Brown   PetscFunctionReturn(0);
254c4762a1bSJed Brown }
255c4762a1bSJed Brown 
256c4762a1bSJed Brown 
257c4762a1bSJed Brown /*TEST
258c4762a1bSJed Brown 
259c4762a1bSJed Brown    build:
260c4762a1bSJed Brown       requires: !complex
261c4762a1bSJed Brown 
262c4762a1bSJed Brown    test:
263c4762a1bSJed Brown       args: -tao_smonitor -tao_type nls -tao_gatol 1.e-4
264c4762a1bSJed Brown       requires: !single
265c4762a1bSJed Brown 
266c4762a1bSJed Brown    test:
267c4762a1bSJed Brown       suffix: 2
268c4762a1bSJed Brown       args: -tao_smonitor -tao_type lmvm -tao_gatol 1.e-3
269c4762a1bSJed Brown 
270c4762a1bSJed Brown    test:
271c4762a1bSJed Brown       suffix: 3
272c4762a1bSJed Brown       args: -tao_smonitor -tao_type ntr -tao_gatol 1.e-4
273c4762a1bSJed Brown       requires: !single
274c4762a1bSJed Brown 
275c4762a1bSJed Brown    test:
276c4762a1bSJed Brown       suffix: 4
277c4762a1bSJed Brown       args: -tao_smonitor -tao_type ntr -tao_mf_hessian -tao_ntr_pc_type none -tao_gatol 1.e-4
278c4762a1bSJed Brown 
279c4762a1bSJed Brown    test:
280c4762a1bSJed Brown       suffix: 5
281c4762a1bSJed Brown       args: -tao_smonitor -tao_type bntr -tao_gatol 1.e-4
282c4762a1bSJed Brown 
283c4762a1bSJed Brown    test:
284c4762a1bSJed Brown       suffix: 6
285c4762a1bSJed Brown       args: -tao_smonitor -tao_type bntl -tao_gatol 1.e-4
286c4762a1bSJed Brown 
287c4762a1bSJed Brown    test:
288c4762a1bSJed Brown       suffix: 7
289c4762a1bSJed Brown       args: -tao_smonitor -tao_type bnls -tao_gatol 1.e-4
290c4762a1bSJed Brown 
291c4762a1bSJed Brown    test:
292c4762a1bSJed Brown       suffix: 8
293c4762a1bSJed Brown       args: -tao_smonitor -tao_type bntr -tao_bnk_max_cg_its 3 -tao_gatol 1.e-4
294c4762a1bSJed Brown 
295c4762a1bSJed Brown    test:
296c4762a1bSJed Brown       suffix: 9
297c4762a1bSJed Brown       args: -tao_smonitor -tao_type bntl -tao_bnk_max_cg_its 3 -tao_gatol 1.e-4
298c4762a1bSJed Brown 
299c4762a1bSJed Brown    test:
300c4762a1bSJed Brown       suffix: 10
301c4762a1bSJed Brown       args: -tao_smonitor -tao_type bnls -tao_bnk_max_cg_its 3 -tao_gatol 1.e-4
302c4762a1bSJed Brown 
303c4762a1bSJed Brown    test:
304c4762a1bSJed Brown       suffix: 11
305864588a7SAlp Dener       args: -test_lmvm -tao_max_it 10 -tao_bqnk_mat_type lmvmbroyden
306c4762a1bSJed Brown 
307c4762a1bSJed Brown    test:
308c4762a1bSJed Brown       suffix: 12
309864588a7SAlp Dener       args: -test_lmvm -tao_max_it 10 -tao_bqnk_mat_type lmvmbadbroyden
310c4762a1bSJed Brown 
311c4762a1bSJed Brown    test:
312c4762a1bSJed Brown      suffix: 13
313864588a7SAlp Dener      args: -test_lmvm -tao_max_it 10 -tao_bqnk_mat_type lmvmsymbroyden
314c4762a1bSJed Brown 
315c4762a1bSJed Brown    test:
316c4762a1bSJed Brown      suffix: 14
317c4762a1bSJed Brown      args: -test_lmvm -tao_max_it 10 -tao_bqnk_mat_type lmvmbfgs
318c4762a1bSJed Brown 
319c4762a1bSJed Brown    test:
320c4762a1bSJed Brown      suffix: 15
321c4762a1bSJed Brown      args: -test_lmvm -tao_max_it 10 -tao_bqnk_mat_type lmvmdfp
322c4762a1bSJed Brown 
323c4762a1bSJed Brown    test:
324c4762a1bSJed Brown      suffix: 16
325c4762a1bSJed Brown      args: -test_lmvm -tao_max_it 10 -tao_bqnk_mat_type lmvmsr1
326c4762a1bSJed Brown 
327c4762a1bSJed Brown    test:
328c4762a1bSJed Brown      suffix: 17
329c4762a1bSJed Brown      args: -tao_smonitor -tao_gatol 1e-4 -tao_type bqnls
330c4762a1bSJed Brown 
331c4762a1bSJed Brown    test:
332c4762a1bSJed Brown      suffix: 18
333c4762a1bSJed Brown      args: -tao_smonitor -tao_gatol 1e-4 -tao_type blmvm
334c4762a1bSJed Brown 
335c4762a1bSJed Brown    test:
336c4762a1bSJed Brown      suffix: 19
337c4762a1bSJed Brown      args: -tao_smonitor -tao_gatol 1e-4 -tao_type bqnktr -tao_bqnk_mat_type lmvmsr1
338c4762a1bSJed Brown 
339c4762a1bSJed Brown    test:
340c4762a1bSJed Brown      suffix: 20
341c4762a1bSJed Brown      args: -tao_monitor -tao_gatol 1e-4 -tao_type blmvm -tao_ls_monitor
342c4762a1bSJed Brown 
343c4762a1bSJed Brown    test:
344c4762a1bSJed Brown      suffix: 21
345864588a7SAlp Dener      args: -test_lmvm -tao_max_it 10 -tao_bqnk_mat_type lmvmsymbadbroyden
346c4762a1bSJed Brown 
34705579b36STristan Konolige    test:
34805579b36STristan Konolige      suffix: 22
34905579b36STristan Konolige      args: -tao_max_it 1 -tao_converged_reason
35005579b36STristan Konolige 
35105579b36STristan Konolige    test:
35205579b36STristan Konolige      suffix: 23
35305579b36STristan Konolige      args: -tao_max_funcs 0 -tao_converged_reason
35405579b36STristan Konolige 
35505579b36STristan Konolige    test:
35605579b36STristan Konolige      suffix: 24
35705579b36STristan Konolige      args: -tao_gatol 10 -tao_converged_reason
35805579b36STristan Konolige 
35905579b36STristan Konolige    test:
36005579b36STristan Konolige      suffix: 25
36105579b36STristan Konolige      args: -tao_grtol 10 -tao_converged_reason
36205579b36STristan Konolige 
36305579b36STristan Konolige    test:
36405579b36STristan Konolige      suffix: 26
36505579b36STristan Konolige      args: -tao_gttol 10 -tao_converged_reason
36605579b36STristan Konolige 
36705579b36STristan Konolige    test:
36805579b36STristan Konolige      suffix: 27
36905579b36STristan Konolige      args: -tao_steptol 10 -tao_converged_reason
37005579b36STristan Konolige 
37105579b36STristan Konolige    test:
37205579b36STristan Konolige      suffix: 28
37305579b36STristan Konolige      args: -tao_fmin 10 -tao_converged_reason
37405579b36STristan Konolige 
375c4762a1bSJed Brown TEST*/
376