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