xref: /petsc/src/tao/unconstrained/tutorials/rosenbrock1.c (revision a0254a939f1187a8a30e788ec1e80a5d3aab8d9e)
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 /*
14c4762a1bSJed Brown    User-defined application context - contains data needed by the
15c4762a1bSJed Brown    application-provided call-back routines that evaluate the function,
16c4762a1bSJed Brown    gradient, and hessian.
17c4762a1bSJed Brown */
18c4762a1bSJed Brown typedef struct {
19c4762a1bSJed Brown   PetscInt  n;     /* dimension */
20c4762a1bSJed Brown   PetscReal alpha; /* condition parameter */
21c4762a1bSJed Brown   PetscBool chained;
22c4762a1bSJed Brown } AppCtx;
23c4762a1bSJed Brown 
24c4762a1bSJed Brown /* -------------- User-defined routines ---------- */
25c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao, Vec, PetscReal *, Vec, void *);
26c4762a1bSJed Brown PetscErrorCode FormHessian(Tao, Vec, Mat, Mat, void *);
27c4762a1bSJed Brown 
28d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
29d71ae5a4SJacob Faibussowitsch {
30c4762a1bSJed Brown   PetscReal   zero = 0.0;
31c4762a1bSJed Brown   Vec         x; /* solution vector */
32c4762a1bSJed Brown   Mat         H;
33c4762a1bSJed Brown   Tao         tao; /* Tao solver context */
34c4762a1bSJed Brown   PetscBool   flg, test_lmvm = PETSC_FALSE;
35c4762a1bSJed Brown   PetscMPIInt size; /* number of processes running */
36c4762a1bSJed Brown   AppCtx      user; /* user-defined application context */
37c4762a1bSJed Brown   KSP         ksp;
38c4762a1bSJed Brown   PC          pc;
39c4762a1bSJed Brown   Mat         M;
40c4762a1bSJed Brown   Vec         in, out, out2;
41c4762a1bSJed Brown   PetscReal   mult_solve_dist;
42c4762a1bSJed Brown 
43c4762a1bSJed Brown   /* Initialize TAO and PETSc */
44327415f7SBarry Smith   PetscFunctionBeginUser;
459566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
469566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
473c859ba3SBarry Smith   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Incorrect number of processors");
48c4762a1bSJed Brown 
49c4762a1bSJed Brown   /* Initialize problem parameters */
509371c9d4SSatish Balay   user.n       = 2;
519371c9d4SSatish Balay   user.alpha   = 99.0;
529371c9d4SSatish Balay   user.chained = PETSC_FALSE;
53c4762a1bSJed Brown   /* Check for command line arguments to override defaults */
549566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &user.n, &flg));
559566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-alpha", &user.alpha, &flg));
569566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-chained", &user.chained, &flg));
579566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_lmvm", &test_lmvm, &flg));
58c4762a1bSJed Brown 
59c4762a1bSJed Brown   /* Allocate vectors for the solution and gradient */
609566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF, user.n, &x));
619566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqBAIJ(PETSC_COMM_SELF, 2, user.n, user.n, 1, NULL, &H));
62c4762a1bSJed Brown 
63c4762a1bSJed Brown   /* The TAO code begins here */
64c4762a1bSJed Brown 
65c4762a1bSJed Brown   /* Create TAO solver with desired solution method */
669566063dSJacob Faibussowitsch   PetscCall(TaoCreate(PETSC_COMM_SELF, &tao));
679566063dSJacob Faibussowitsch   PetscCall(TaoSetType(tao, TAOLMVM));
68c4762a1bSJed Brown 
69c4762a1bSJed Brown   /* Set solution vec and an initial guess */
709566063dSJacob Faibussowitsch   PetscCall(VecSet(x, zero));
719566063dSJacob Faibussowitsch   PetscCall(TaoSetSolution(tao, x));
72c4762a1bSJed Brown 
73c4762a1bSJed Brown   /* Set routines for function, gradient, hessian evaluation */
749566063dSJacob Faibussowitsch   PetscCall(TaoSetObjectiveAndGradient(tao, NULL, FormFunctionGradient, &user));
759566063dSJacob Faibussowitsch   PetscCall(TaoSetHessian(tao, H, H, FormHessian, &user));
76c4762a1bSJed Brown 
77c4762a1bSJed Brown   /* Test the LMVM matrix */
7848a46eb9SPierre Jolivet   if (test_lmvm) PetscCall(PetscOptionsSetValue(NULL, "-tao_type", "bqnktr"));
79c4762a1bSJed Brown 
80c4762a1bSJed Brown   /* Check for TAO command line options */
819566063dSJacob Faibussowitsch   PetscCall(TaoSetFromOptions(tao));
82c4762a1bSJed Brown 
83c4762a1bSJed Brown   /* SOLVE THE APPLICATION */
849566063dSJacob Faibussowitsch   PetscCall(TaoSolve(tao));
85c4762a1bSJed Brown 
86c4762a1bSJed Brown   /* Test the LMVM matrix */
87c4762a1bSJed Brown   if (test_lmvm) {
889566063dSJacob Faibussowitsch     PetscCall(TaoGetKSP(tao, &ksp));
899566063dSJacob Faibussowitsch     PetscCall(KSPGetPC(ksp, &pc));
909566063dSJacob Faibussowitsch     PetscCall(PCLMVMGetMatLMVM(pc, &M));
919566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(x, &in));
929566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(x, &out));
939566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(x, &out2));
949566063dSJacob Faibussowitsch     PetscCall(VecSet(in, 1.0));
959566063dSJacob Faibussowitsch     PetscCall(MatMult(M, in, out));
969566063dSJacob Faibussowitsch     PetscCall(MatSolve(M, out, out2));
979566063dSJacob Faibussowitsch     PetscCall(VecAXPY(out2, -1.0, in));
989566063dSJacob Faibussowitsch     PetscCall(VecNorm(out2, NORM_2, &mult_solve_dist));
99c4762a1bSJed Brown     if (mult_solve_dist < 1.e-11) {
1009566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PetscObjectComm((PetscObject)tao), "error between LMVM MatMult and MatSolve: < 1.e-11\n"));
101c4762a1bSJed Brown     } else if (mult_solve_dist < 1.e-6) {
1029566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PetscObjectComm((PetscObject)tao), "error between LMVM MatMult and MatSolve: < 1.e-6\n"));
103c4762a1bSJed Brown     } else {
1049566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PetscObjectComm((PetscObject)tao), "error between LMVM MatMult and MatSolve: %e\n", (double)mult_solve_dist));
105c4762a1bSJed Brown     }
1069566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&in));
1079566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&out));
1089566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&out2));
109c4762a1bSJed Brown   }
110c4762a1bSJed Brown 
1119566063dSJacob Faibussowitsch   PetscCall(TaoDestroy(&tao));
1129566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
1139566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&H));
114c4762a1bSJed Brown 
1159566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
116b122ec5aSJacob Faibussowitsch   return 0;
117c4762a1bSJed Brown }
118c4762a1bSJed Brown 
119c4762a1bSJed Brown /* -------------------------------------------------------------------- */
120c4762a1bSJed Brown /*
121c4762a1bSJed Brown     FormFunctionGradient - Evaluates the function, f(X), and gradient, G(X).
122c4762a1bSJed Brown 
123c4762a1bSJed Brown     Input Parameters:
124c4762a1bSJed Brown .   tao  - the Tao context
125c4762a1bSJed Brown .   X    - input vector
126c4762a1bSJed Brown .   ptr  - optional user-defined context, as set by TaoSetFunctionGradient()
127c4762a1bSJed Brown 
128c4762a1bSJed Brown     Output Parameters:
129c4762a1bSJed Brown .   G - vector containing the newly evaluated gradient
130c4762a1bSJed Brown .   f - function value
131c4762a1bSJed Brown 
132c4762a1bSJed Brown     Note:
133c4762a1bSJed Brown     Some optimization methods ask for the function and the gradient evaluation
134c4762a1bSJed Brown     at the same time.  Evaluating both at once may be more efficient that
135c4762a1bSJed Brown     evaluating each separately.
136c4762a1bSJed Brown */
137d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *f, Vec G, void *ptr)
138d71ae5a4SJacob Faibussowitsch {
139c4762a1bSJed Brown   AppCtx            *user = (AppCtx *)ptr;
140c4762a1bSJed Brown   PetscInt           i, nn = user->n / 2;
141c4762a1bSJed Brown   PetscReal          ff = 0, t1, t2, alpha = user->alpha;
142c4762a1bSJed Brown   PetscScalar       *g;
143c4762a1bSJed Brown   const PetscScalar *x;
144c4762a1bSJed Brown 
145c4762a1bSJed Brown   PetscFunctionBeginUser;
146c4762a1bSJed Brown   /* Get pointers to vector data */
1479566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
1489566063dSJacob Faibussowitsch   PetscCall(VecGetArray(G, &g));
149c4762a1bSJed Brown 
150c4762a1bSJed Brown   /* Compute G(X) */
151c4762a1bSJed Brown   if (user->chained) {
152c4762a1bSJed Brown     g[0] = 0;
153c4762a1bSJed Brown     for (i = 0; i < user->n - 1; i++) {
154c4762a1bSJed Brown       t1 = x[i + 1] - x[i] * x[i];
155c4762a1bSJed Brown       ff += PetscSqr(1 - x[i]) + alpha * t1 * t1;
156c4762a1bSJed Brown       g[i] += -2 * (1 - x[i]) + 2 * alpha * t1 * (-2 * x[i]);
157c4762a1bSJed Brown       g[i + 1] = 2 * alpha * t1;
158c4762a1bSJed Brown     }
159c4762a1bSJed Brown   } else {
160c4762a1bSJed Brown     for (i = 0; i < nn; i++) {
1619371c9d4SSatish Balay       t1 = x[2 * i + 1] - x[2 * i] * x[2 * i];
1629371c9d4SSatish Balay       t2 = 1 - x[2 * i];
163c4762a1bSJed Brown       ff += alpha * t1 * t1 + t2 * t2;
164c4762a1bSJed Brown       g[2 * i]     = -4 * alpha * t1 * x[2 * i] - 2.0 * t2;
165c4762a1bSJed Brown       g[2 * i + 1] = 2 * alpha * t1;
166c4762a1bSJed Brown     }
167c4762a1bSJed Brown   }
168c4762a1bSJed Brown 
169c4762a1bSJed Brown   /* Restore vectors */
1709566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
1719566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(G, &g));
172c4762a1bSJed Brown   *f = ff;
173c4762a1bSJed Brown 
1749566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(15.0 * nn));
1753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
176c4762a1bSJed Brown }
177c4762a1bSJed Brown 
178c4762a1bSJed Brown /* ------------------------------------------------------------------- */
179c4762a1bSJed Brown /*
180c4762a1bSJed Brown    FormHessian - Evaluates Hessian matrix.
181c4762a1bSJed Brown 
182c4762a1bSJed Brown    Input Parameters:
183c4762a1bSJed Brown .  tao   - the Tao context
184c4762a1bSJed Brown .  x     - input vector
185c4762a1bSJed Brown .  ptr   - optional user-defined context, as set by TaoSetHessian()
186c4762a1bSJed Brown 
187c4762a1bSJed Brown    Output Parameters:
188c4762a1bSJed Brown .  H     - Hessian matrix
189c4762a1bSJed Brown 
190c4762a1bSJed Brown    Note:  Providing the Hessian may not be necessary.  Only some solvers
191c4762a1bSJed Brown    require this matrix.
192c4762a1bSJed Brown */
193d71ae5a4SJacob Faibussowitsch PetscErrorCode FormHessian(Tao tao, Vec X, Mat H, Mat Hpre, void *ptr)
194d71ae5a4SJacob Faibussowitsch {
195c4762a1bSJed Brown   AppCtx            *user = (AppCtx *)ptr;
196c4762a1bSJed Brown   PetscInt           i, ind[2];
197c4762a1bSJed Brown   PetscReal          alpha = user->alpha;
198c4762a1bSJed Brown   PetscReal          v[2][2];
199c4762a1bSJed Brown   const PetscScalar *x;
200c4762a1bSJed Brown   PetscBool          assembled;
201c4762a1bSJed Brown 
202c4762a1bSJed Brown   PetscFunctionBeginUser;
203c4762a1bSJed Brown   /* Zero existing matrix entries */
2049566063dSJacob Faibussowitsch   PetscCall(MatAssembled(H, &assembled));
2059566063dSJacob Faibussowitsch   if (assembled) PetscCall(MatZeroEntries(H));
206c4762a1bSJed Brown 
207c4762a1bSJed Brown   /* Get a pointer to vector data */
2089566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
209c4762a1bSJed Brown 
210c4762a1bSJed Brown   /* Compute H(X) entries */
211c4762a1bSJed Brown   if (user->chained) {
2129566063dSJacob Faibussowitsch     PetscCall(MatZeroEntries(H));
213c4762a1bSJed Brown     for (i = 0; i < user->n - 1; i++) {
214c4762a1bSJed Brown       PetscScalar t1 = x[i + 1] - x[i] * x[i];
215c4762a1bSJed Brown       v[0][0]        = 2 + 2 * alpha * (t1 * (-2) - 2 * x[i]);
216c4762a1bSJed Brown       v[0][1]        = 2 * alpha * (-2 * x[i]);
217c4762a1bSJed Brown       v[1][0]        = 2 * alpha * (-2 * x[i]);
218c4762a1bSJed Brown       v[1][1]        = 2 * alpha * t1;
2199371c9d4SSatish Balay       ind[0]         = i;
2209371c9d4SSatish Balay       ind[1]         = i + 1;
2219566063dSJacob Faibussowitsch       PetscCall(MatSetValues(H, 2, ind, 2, ind, v[0], ADD_VALUES));
222c4762a1bSJed Brown     }
223c4762a1bSJed Brown   } else {
224c4762a1bSJed Brown     for (i = 0; i < user->n / 2; i++) {
225c4762a1bSJed Brown       v[1][1] = 2 * alpha;
226c4762a1bSJed Brown       v[0][0] = -4 * alpha * (x[2 * i + 1] - 3 * x[2 * i] * x[2 * i]) + 2;
227c4762a1bSJed Brown       v[1][0] = v[0][1] = -4.0 * alpha * x[2 * i];
2289371c9d4SSatish Balay       ind[0]            = 2 * i;
2299371c9d4SSatish Balay       ind[1]            = 2 * i + 1;
2309566063dSJacob Faibussowitsch       PetscCall(MatSetValues(H, 2, ind, 2, ind, v[0], INSERT_VALUES));
231c4762a1bSJed Brown     }
232c4762a1bSJed Brown   }
2339566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
234c4762a1bSJed Brown 
235c4762a1bSJed Brown   /* Assemble matrix */
2369566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(H, MAT_FINAL_ASSEMBLY));
2379566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(H, MAT_FINAL_ASSEMBLY));
2389566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(9.0 * user->n / 2.0));
2393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
240c4762a1bSJed Brown }
241c4762a1bSJed Brown 
242c4762a1bSJed Brown /*TEST
243c4762a1bSJed Brown 
244c4762a1bSJed Brown    build:
245c4762a1bSJed Brown       requires: !complex
246c4762a1bSJed Brown 
247c4762a1bSJed Brown    test:
248c4762a1bSJed Brown       args: -tao_smonitor -tao_type nls -tao_gatol 1.e-4
249c4762a1bSJed Brown       requires: !single
250c4762a1bSJed Brown 
251c4762a1bSJed Brown    test:
252c4762a1bSJed Brown       suffix: 2
253c4762a1bSJed Brown       args: -tao_smonitor -tao_type lmvm -tao_gatol 1.e-3
254c4762a1bSJed Brown 
255c4762a1bSJed Brown    test:
256c4762a1bSJed Brown       suffix: 3
257c4762a1bSJed Brown       args: -tao_smonitor -tao_type ntr -tao_gatol 1.e-4
258c4762a1bSJed Brown       requires: !single
259c4762a1bSJed Brown 
260c4762a1bSJed Brown    test:
261c4762a1bSJed Brown       suffix: 4
262c4762a1bSJed Brown       args: -tao_smonitor -tao_type ntr -tao_mf_hessian -tao_ntr_pc_type none -tao_gatol 1.e-4
263c4762a1bSJed Brown 
264c4762a1bSJed Brown    test:
265c4762a1bSJed Brown       suffix: 5
266c4762a1bSJed Brown       args: -tao_smonitor -tao_type bntr -tao_gatol 1.e-4
267c4762a1bSJed Brown 
268c4762a1bSJed Brown    test:
269c4762a1bSJed Brown       suffix: 6
270c4762a1bSJed Brown       args: -tao_smonitor -tao_type bntl -tao_gatol 1.e-4
271c4762a1bSJed Brown 
272c4762a1bSJed Brown    test:
273c4762a1bSJed Brown       suffix: 7
274c4762a1bSJed Brown       args: -tao_smonitor -tao_type bnls -tao_gatol 1.e-4
275c4762a1bSJed Brown 
276c4762a1bSJed Brown    test:
277c4762a1bSJed Brown       suffix: 8
278c4762a1bSJed Brown       args: -tao_smonitor -tao_type bntr -tao_bnk_max_cg_its 3 -tao_gatol 1.e-4
279c4762a1bSJed Brown 
280c4762a1bSJed Brown    test:
281c4762a1bSJed Brown       suffix: 9
282c4762a1bSJed Brown       args: -tao_smonitor -tao_type bntl -tao_bnk_max_cg_its 3 -tao_gatol 1.e-4
283c4762a1bSJed Brown 
284c4762a1bSJed Brown    test:
285c4762a1bSJed Brown       suffix: 10
286c4762a1bSJed Brown       args: -tao_smonitor -tao_type bnls -tao_bnk_max_cg_its 3 -tao_gatol 1.e-4
287c4762a1bSJed Brown 
288c4762a1bSJed Brown    test:
289c4762a1bSJed Brown       suffix: 11
290864588a7SAlp Dener       args: -test_lmvm -tao_max_it 10 -tao_bqnk_mat_type lmvmbroyden
291c4762a1bSJed Brown 
292c4762a1bSJed Brown    test:
293c4762a1bSJed Brown       suffix: 12
294864588a7SAlp Dener       args: -test_lmvm -tao_max_it 10 -tao_bqnk_mat_type lmvmbadbroyden
295c4762a1bSJed Brown 
296c4762a1bSJed Brown    test:
297c4762a1bSJed Brown      suffix: 13
298864588a7SAlp Dener      args: -test_lmvm -tao_max_it 10 -tao_bqnk_mat_type lmvmsymbroyden
299c4762a1bSJed Brown 
300c4762a1bSJed Brown    test:
301c4762a1bSJed Brown      suffix: 14
302c4762a1bSJed Brown      args: -test_lmvm -tao_max_it 10 -tao_bqnk_mat_type lmvmbfgs
303c4762a1bSJed Brown 
304c4762a1bSJed Brown    test:
305c4762a1bSJed Brown      suffix: 15
306c4762a1bSJed Brown      args: -test_lmvm -tao_max_it 10 -tao_bqnk_mat_type lmvmdfp
307c4762a1bSJed Brown 
308c4762a1bSJed Brown    test:
309c4762a1bSJed Brown      suffix: 16
310c4762a1bSJed Brown      args: -test_lmvm -tao_max_it 10 -tao_bqnk_mat_type lmvmsr1
311c4762a1bSJed Brown 
312c4762a1bSJed Brown    test:
313c4762a1bSJed Brown      suffix: 17
314c4762a1bSJed Brown      args: -tao_smonitor -tao_gatol 1e-4 -tao_type bqnls
315c4762a1bSJed Brown 
316c4762a1bSJed Brown    test:
317c4762a1bSJed Brown      suffix: 18
318c4762a1bSJed Brown      args: -tao_smonitor -tao_gatol 1e-4 -tao_type blmvm
319c4762a1bSJed Brown 
320c4762a1bSJed Brown    test:
321c4762a1bSJed Brown      suffix: 19
322c4762a1bSJed Brown      args: -tao_smonitor -tao_gatol 1e-4 -tao_type bqnktr -tao_bqnk_mat_type lmvmsr1
323c4762a1bSJed Brown 
324c4762a1bSJed Brown    test:
325c4762a1bSJed Brown      suffix: 20
326c4762a1bSJed Brown      args: -tao_monitor -tao_gatol 1e-4 -tao_type blmvm -tao_ls_monitor
327c4762a1bSJed Brown 
328c4762a1bSJed Brown    test:
329c4762a1bSJed Brown      suffix: 21
330864588a7SAlp Dener      args: -test_lmvm -tao_max_it 10 -tao_bqnk_mat_type lmvmsymbadbroyden
331c4762a1bSJed Brown 
33205579b36STristan Konolige    test:
33305579b36STristan Konolige      suffix: 22
33405579b36STristan Konolige      args: -tao_max_it 1 -tao_converged_reason
33505579b36STristan Konolige 
33605579b36STristan Konolige    test:
33705579b36STristan Konolige      suffix: 23
33805579b36STristan Konolige      args: -tao_max_funcs 0 -tao_converged_reason
33905579b36STristan Konolige 
34005579b36STristan Konolige    test:
34105579b36STristan Konolige      suffix: 24
34205579b36STristan Konolige      args: -tao_gatol 10 -tao_converged_reason
34305579b36STristan Konolige 
34405579b36STristan Konolige    test:
34505579b36STristan Konolige      suffix: 25
34605579b36STristan Konolige      args: -tao_grtol 10 -tao_converged_reason
34705579b36STristan Konolige 
34805579b36STristan Konolige    test:
34905579b36STristan Konolige      suffix: 26
35005579b36STristan Konolige      args: -tao_gttol 10 -tao_converged_reason
35105579b36STristan Konolige 
35205579b36STristan Konolige    test:
35305579b36STristan Konolige      suffix: 27
35405579b36STristan Konolige      args: -tao_steptol 10 -tao_converged_reason
35505579b36STristan Konolige 
35605579b36STristan Konolige    test:
35705579b36STristan Konolige      suffix: 28
35805579b36STristan Konolige      args: -tao_fmin 10 -tao_converged_reason
35905579b36STristan Konolige 
360f4f59681SStefano Zampini    test:
361f4f59681SStefano Zampini      suffix: snes
362f4f59681SStefano Zampini      args: -snes_monitor ::ascii_info_detail -tao_type snes -snes_type newtontr -snes_atol 1.e-4 -pc_type none -tao_mf_hessian -ksp_type cg
363f4f59681SStefano Zampini 
364f4f59681SStefano Zampini    test:
365f4f59681SStefano Zampini      suffix: snes_ls_armijo
366f4f59681SStefano Zampini      args: -snes_monitor ::ascii_info_detail -tao_type snes -snes_type newtonls -snes_atol 1.e-4 -pc_type none -tao_mf_hessian -snes_linesearch_monitor -snes_linesearch_order 1
367f4f59681SStefano Zampini 
3684b0a5c37SStefano Zampini    test:
3694b0a5c37SStefano Zampini      suffix: snes_tr_cgnegcurve_kmdc
3704b0a5c37SStefano Zampini      args: -snes_monitor ::ascii_info_detail -tao_type snes -snes_type newtontr -snes_atol 1.e-4 -pc_type none  -ksp_type cg -snes_tr_kmdc 0.01 -ksp_converged_neg_curve -ksp_converged_reason
3714b0a5c37SStefano Zampini      requires: !single
3724b0a5c37SStefano Zampini 
373*a0254a93SStefano Zampini    test:
374*a0254a93SStefano Zampini      suffix: snes_ls_lmvm
375*a0254a93SStefano Zampini      args: -snes_monitor ::ascii_info_detail -tao_type snes -snes_type newtonls -snes_atol 1.e-4 -pc_type lmvm -tao_mf_hessian
376*a0254a93SStefano Zampini 
377c4762a1bSJed Brown TEST*/
378