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 PetscReal zero=0.0; 43c4762a1bSJed Brown Vec x; /* solution vector */ 44c4762a1bSJed Brown Mat H; 45c4762a1bSJed Brown Tao tao; /* Tao solver context */ 46c4762a1bSJed Brown PetscBool flg, test_lmvm = PETSC_FALSE; 47c4762a1bSJed Brown PetscMPIInt size; /* number of processes running */ 48c4762a1bSJed Brown AppCtx user; /* user-defined application context */ 49c4762a1bSJed Brown KSP ksp; 50c4762a1bSJed Brown PC pc; 51c4762a1bSJed Brown Mat M; 52c4762a1bSJed Brown Vec in, out, out2; 53c4762a1bSJed Brown PetscReal mult_solve_dist; 54c4762a1bSJed Brown 55c4762a1bSJed Brown /* Initialize TAO and PETSc */ 56*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help)); 575f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 583c859ba3SBarry Smith PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Incorrect number of processors"); 59c4762a1bSJed Brown 60c4762a1bSJed Brown /* Initialize problem parameters */ 61c4762a1bSJed Brown user.n = 2; user.alpha = 99.0; user.chained = PETSC_FALSE; 62c4762a1bSJed Brown /* Check for command line arguments to override defaults */ 635f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n",&user.n,&flg)); 645f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-alpha",&user.alpha,&flg)); 655f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-chained",&user.chained,&flg)); 665f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-test_lmvm",&test_lmvm,&flg)); 67c4762a1bSJed Brown 68c4762a1bSJed Brown /* Allocate vectors for the solution and gradient */ 695f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,user.n,&x)); 705f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSeqBAIJ(PETSC_COMM_SELF,2,user.n,user.n,1,NULL,&H)); 71c4762a1bSJed Brown 72c4762a1bSJed Brown /* The TAO code begins here */ 73c4762a1bSJed Brown 74c4762a1bSJed Brown /* Create TAO solver with desired solution method */ 755f80ce2aSJacob Faibussowitsch CHKERRQ(TaoCreate(PETSC_COMM_SELF,&tao)); 765f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSetType(tao,TAOLMVM)); 77c4762a1bSJed Brown 78c4762a1bSJed Brown /* Set solution vec and an initial guess */ 795f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(x, zero)); 805f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSetSolution(tao,x)); 81c4762a1bSJed Brown 82c4762a1bSJed Brown /* Set routines for function, gradient, hessian evaluation */ 835f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSetObjectiveAndGradient(tao,NULL,FormFunctionGradient,&user)); 845f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSetHessian(tao,H,H,FormHessian,&user)); 85c4762a1bSJed Brown 86c4762a1bSJed Brown /* Test the LMVM matrix */ 87c4762a1bSJed Brown if (test_lmvm) { 885f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsSetValue(NULL, "-tao_type", "bqnktr")); 89c4762a1bSJed Brown } 90c4762a1bSJed Brown 91c4762a1bSJed Brown /* Check for TAO command line options */ 925f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSetFromOptions(tao)); 93c4762a1bSJed Brown 94c4762a1bSJed Brown /* SOLVE THE APPLICATION */ 955f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSolve(tao)); 96c4762a1bSJed Brown 97c4762a1bSJed Brown /* Test the LMVM matrix */ 98c4762a1bSJed Brown if (test_lmvm) { 995f80ce2aSJacob Faibussowitsch CHKERRQ(TaoGetKSP(tao, &ksp)); 1005f80ce2aSJacob Faibussowitsch CHKERRQ(KSPGetPC(ksp, &pc)); 1015f80ce2aSJacob Faibussowitsch CHKERRQ(PCLMVMGetMatLMVM(pc, &M)); 1025f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(x, &in)); 1035f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(x, &out)); 1045f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(x, &out2)); 1055f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(in, 1.0)); 1065f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(M, in, out)); 1075f80ce2aSJacob Faibussowitsch CHKERRQ(MatSolve(M, out, out2)); 1085f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(out2, -1.0, in)); 1095f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(out2, NORM_2, &mult_solve_dist)); 110c4762a1bSJed Brown if (mult_solve_dist < 1.e-11) { 1115f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PetscObjectComm((PetscObject)tao), "error between LMVM MatMult and MatSolve: < 1.e-11\n")); 112c4762a1bSJed Brown } else if (mult_solve_dist < 1.e-6) { 1135f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PetscObjectComm((PetscObject)tao), "error between LMVM MatMult and MatSolve: < 1.e-6\n")); 114c4762a1bSJed Brown } else { 1155f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PetscObjectComm((PetscObject)tao), "error between LMVM MatMult and MatSolve: %e\n", (double)mult_solve_dist)); 116c4762a1bSJed Brown } 1175f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&in)); 1185f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&out)); 1195f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&out2)); 120c4762a1bSJed Brown } 121c4762a1bSJed Brown 1225f80ce2aSJacob Faibussowitsch CHKERRQ(TaoDestroy(&tao)); 1235f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&x)); 1245f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&H)); 125c4762a1bSJed Brown 126*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 127*b122ec5aSJacob Faibussowitsch return 0; 128c4762a1bSJed Brown } 129c4762a1bSJed Brown 130c4762a1bSJed Brown /* -------------------------------------------------------------------- */ 131c4762a1bSJed Brown /* 132c4762a1bSJed Brown FormFunctionGradient - Evaluates the function, f(X), and gradient, G(X). 133c4762a1bSJed Brown 134c4762a1bSJed Brown Input Parameters: 135c4762a1bSJed Brown . tao - the Tao context 136c4762a1bSJed Brown . X - input vector 137c4762a1bSJed Brown . ptr - optional user-defined context, as set by TaoSetFunctionGradient() 138c4762a1bSJed Brown 139c4762a1bSJed Brown Output Parameters: 140c4762a1bSJed Brown . G - vector containing the newly evaluated gradient 141c4762a1bSJed Brown . f - function value 142c4762a1bSJed Brown 143c4762a1bSJed Brown Note: 144c4762a1bSJed Brown Some optimization methods ask for the function and the gradient evaluation 145c4762a1bSJed Brown at the same time. Evaluating both at once may be more efficient that 146c4762a1bSJed Brown evaluating each separately. 147c4762a1bSJed Brown */ 148c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao tao,Vec X,PetscReal *f, Vec G,void *ptr) 149c4762a1bSJed Brown { 150c4762a1bSJed Brown AppCtx *user = (AppCtx *) ptr; 151c4762a1bSJed Brown PetscInt i,nn=user->n/2; 152c4762a1bSJed Brown PetscReal ff=0,t1,t2,alpha=user->alpha; 153c4762a1bSJed Brown PetscScalar *g; 154c4762a1bSJed Brown const PetscScalar *x; 155c4762a1bSJed Brown 156c4762a1bSJed Brown PetscFunctionBeginUser; 157c4762a1bSJed Brown /* Get pointers to vector data */ 1585f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(X,&x)); 1595f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(G,&g)); 160c4762a1bSJed Brown 161c4762a1bSJed Brown /* Compute G(X) */ 162c4762a1bSJed Brown if (user->chained) { 163c4762a1bSJed Brown g[0] = 0; 164c4762a1bSJed Brown for (i=0; i<user->n-1; i++) { 165c4762a1bSJed Brown t1 = x[i+1] - x[i]*x[i]; 166c4762a1bSJed Brown ff += PetscSqr(1 - x[i]) + alpha*t1*t1; 167c4762a1bSJed Brown g[i] += -2*(1 - x[i]) + 2*alpha*t1*(-2*x[i]); 168c4762a1bSJed Brown g[i+1] = 2*alpha*t1; 169c4762a1bSJed Brown } 170c4762a1bSJed Brown } else { 171c4762a1bSJed Brown for (i=0; i<nn; i++) { 172c4762a1bSJed Brown t1 = x[2*i+1]-x[2*i]*x[2*i]; t2= 1-x[2*i]; 173c4762a1bSJed Brown ff += alpha*t1*t1 + t2*t2; 174c4762a1bSJed Brown g[2*i] = -4*alpha*t1*x[2*i]-2.0*t2; 175c4762a1bSJed Brown g[2*i+1] = 2*alpha*t1; 176c4762a1bSJed Brown } 177c4762a1bSJed Brown } 178c4762a1bSJed Brown 179c4762a1bSJed Brown /* Restore vectors */ 1805f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(X,&x)); 1815f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(G,&g)); 182c4762a1bSJed Brown *f = ff; 183c4762a1bSJed Brown 1845f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogFlops(15.0*nn)); 185c4762a1bSJed Brown PetscFunctionReturn(0); 186c4762a1bSJed Brown } 187c4762a1bSJed Brown 188c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 189c4762a1bSJed Brown /* 190c4762a1bSJed Brown FormHessian - Evaluates Hessian matrix. 191c4762a1bSJed Brown 192c4762a1bSJed Brown Input Parameters: 193c4762a1bSJed Brown . tao - the Tao context 194c4762a1bSJed Brown . x - input vector 195c4762a1bSJed Brown . ptr - optional user-defined context, as set by TaoSetHessian() 196c4762a1bSJed Brown 197c4762a1bSJed Brown Output Parameters: 198c4762a1bSJed Brown . H - Hessian matrix 199c4762a1bSJed Brown 200c4762a1bSJed Brown Note: Providing the Hessian may not be necessary. Only some solvers 201c4762a1bSJed Brown require this matrix. 202c4762a1bSJed Brown */ 203c4762a1bSJed Brown PetscErrorCode FormHessian(Tao tao,Vec X,Mat H, Mat Hpre, void *ptr) 204c4762a1bSJed Brown { 205c4762a1bSJed Brown AppCtx *user = (AppCtx*)ptr; 206c4762a1bSJed Brown PetscInt i, ind[2]; 207c4762a1bSJed Brown PetscReal alpha=user->alpha; 208c4762a1bSJed Brown PetscReal v[2][2]; 209c4762a1bSJed Brown const PetscScalar *x; 210c4762a1bSJed Brown PetscBool assembled; 211c4762a1bSJed Brown 212c4762a1bSJed Brown PetscFunctionBeginUser; 213c4762a1bSJed Brown /* Zero existing matrix entries */ 2145f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssembled(H,&assembled)); 2155f80ce2aSJacob Faibussowitsch if (assembled) CHKERRQ(MatZeroEntries(H)); 216c4762a1bSJed Brown 217c4762a1bSJed Brown /* Get a pointer to vector data */ 2185f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(X,&x)); 219c4762a1bSJed Brown 220c4762a1bSJed Brown /* Compute H(X) entries */ 221c4762a1bSJed Brown if (user->chained) { 2225f80ce2aSJacob Faibussowitsch CHKERRQ(MatZeroEntries(H)); 223c4762a1bSJed Brown for (i=0; i<user->n-1; i++) { 224c4762a1bSJed Brown PetscScalar t1 = x[i+1] - x[i]*x[i]; 225c4762a1bSJed Brown v[0][0] = 2 + 2*alpha*(t1*(-2) - 2*x[i]); 226c4762a1bSJed Brown v[0][1] = 2*alpha*(-2*x[i]); 227c4762a1bSJed Brown v[1][0] = 2*alpha*(-2*x[i]); 228c4762a1bSJed Brown v[1][1] = 2*alpha*t1; 229c4762a1bSJed Brown ind[0] = i; ind[1] = i+1; 2305f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(H,2,ind,2,ind,v[0],ADD_VALUES)); 231c4762a1bSJed Brown } 232c4762a1bSJed Brown } else { 233c4762a1bSJed Brown for (i=0; i<user->n/2; i++) { 234c4762a1bSJed Brown v[1][1] = 2*alpha; 235c4762a1bSJed Brown v[0][0] = -4*alpha*(x[2*i+1]-3*x[2*i]*x[2*i]) + 2; 236c4762a1bSJed Brown v[1][0] = v[0][1] = -4.0*alpha*x[2*i]; 237c4762a1bSJed Brown ind[0]=2*i; ind[1]=2*i+1; 2385f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(H,2,ind,2,ind,v[0],INSERT_VALUES)); 239c4762a1bSJed Brown } 240c4762a1bSJed Brown } 2415f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(X,&x)); 242c4762a1bSJed Brown 243c4762a1bSJed Brown /* Assemble matrix */ 2445f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(H,MAT_FINAL_ASSEMBLY)); 2455f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(H,MAT_FINAL_ASSEMBLY)); 2465f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogFlops(9.0*user->n/2.0)); 247c4762a1bSJed Brown PetscFunctionReturn(0); 248c4762a1bSJed Brown } 249c4762a1bSJed Brown 250c4762a1bSJed Brown /*TEST 251c4762a1bSJed Brown 252c4762a1bSJed Brown build: 253c4762a1bSJed Brown requires: !complex 254c4762a1bSJed Brown 255c4762a1bSJed Brown test: 256c4762a1bSJed Brown args: -tao_smonitor -tao_type nls -tao_gatol 1.e-4 257c4762a1bSJed Brown requires: !single 258c4762a1bSJed Brown 259c4762a1bSJed Brown test: 260c4762a1bSJed Brown suffix: 2 261c4762a1bSJed Brown args: -tao_smonitor -tao_type lmvm -tao_gatol 1.e-3 262c4762a1bSJed Brown 263c4762a1bSJed Brown test: 264c4762a1bSJed Brown suffix: 3 265c4762a1bSJed Brown args: -tao_smonitor -tao_type ntr -tao_gatol 1.e-4 266c4762a1bSJed Brown requires: !single 267c4762a1bSJed Brown 268c4762a1bSJed Brown test: 269c4762a1bSJed Brown suffix: 4 270c4762a1bSJed Brown args: -tao_smonitor -tao_type ntr -tao_mf_hessian -tao_ntr_pc_type none -tao_gatol 1.e-4 271c4762a1bSJed Brown 272c4762a1bSJed Brown test: 273c4762a1bSJed Brown suffix: 5 274c4762a1bSJed Brown args: -tao_smonitor -tao_type bntr -tao_gatol 1.e-4 275c4762a1bSJed Brown 276c4762a1bSJed Brown test: 277c4762a1bSJed Brown suffix: 6 278c4762a1bSJed Brown args: -tao_smonitor -tao_type bntl -tao_gatol 1.e-4 279c4762a1bSJed Brown 280c4762a1bSJed Brown test: 281c4762a1bSJed Brown suffix: 7 282c4762a1bSJed Brown args: -tao_smonitor -tao_type bnls -tao_gatol 1.e-4 283c4762a1bSJed Brown 284c4762a1bSJed Brown test: 285c4762a1bSJed Brown suffix: 8 286c4762a1bSJed Brown args: -tao_smonitor -tao_type bntr -tao_bnk_max_cg_its 3 -tao_gatol 1.e-4 287c4762a1bSJed Brown 288c4762a1bSJed Brown test: 289c4762a1bSJed Brown suffix: 9 290c4762a1bSJed Brown args: -tao_smonitor -tao_type bntl -tao_bnk_max_cg_its 3 -tao_gatol 1.e-4 291c4762a1bSJed Brown 292c4762a1bSJed Brown test: 293c4762a1bSJed Brown suffix: 10 294c4762a1bSJed Brown args: -tao_smonitor -tao_type bnls -tao_bnk_max_cg_its 3 -tao_gatol 1.e-4 295c4762a1bSJed Brown 296c4762a1bSJed Brown test: 297c4762a1bSJed Brown suffix: 11 298864588a7SAlp Dener args: -test_lmvm -tao_max_it 10 -tao_bqnk_mat_type lmvmbroyden 299c4762a1bSJed Brown 300c4762a1bSJed Brown test: 301c4762a1bSJed Brown suffix: 12 302864588a7SAlp Dener args: -test_lmvm -tao_max_it 10 -tao_bqnk_mat_type lmvmbadbroyden 303c4762a1bSJed Brown 304c4762a1bSJed Brown test: 305c4762a1bSJed Brown suffix: 13 306864588a7SAlp Dener args: -test_lmvm -tao_max_it 10 -tao_bqnk_mat_type lmvmsymbroyden 307c4762a1bSJed Brown 308c4762a1bSJed Brown test: 309c4762a1bSJed Brown suffix: 14 310c4762a1bSJed Brown args: -test_lmvm -tao_max_it 10 -tao_bqnk_mat_type lmvmbfgs 311c4762a1bSJed Brown 312c4762a1bSJed Brown test: 313c4762a1bSJed Brown suffix: 15 314c4762a1bSJed Brown args: -test_lmvm -tao_max_it 10 -tao_bqnk_mat_type lmvmdfp 315c4762a1bSJed Brown 316c4762a1bSJed Brown test: 317c4762a1bSJed Brown suffix: 16 318c4762a1bSJed Brown args: -test_lmvm -tao_max_it 10 -tao_bqnk_mat_type lmvmsr1 319c4762a1bSJed Brown 320c4762a1bSJed Brown test: 321c4762a1bSJed Brown suffix: 17 322c4762a1bSJed Brown args: -tao_smonitor -tao_gatol 1e-4 -tao_type bqnls 323c4762a1bSJed Brown 324c4762a1bSJed Brown test: 325c4762a1bSJed Brown suffix: 18 326c4762a1bSJed Brown args: -tao_smonitor -tao_gatol 1e-4 -tao_type blmvm 327c4762a1bSJed Brown 328c4762a1bSJed Brown test: 329c4762a1bSJed Brown suffix: 19 330c4762a1bSJed Brown args: -tao_smonitor -tao_gatol 1e-4 -tao_type bqnktr -tao_bqnk_mat_type lmvmsr1 331c4762a1bSJed Brown 332c4762a1bSJed Brown test: 333c4762a1bSJed Brown suffix: 20 334c4762a1bSJed Brown args: -tao_monitor -tao_gatol 1e-4 -tao_type blmvm -tao_ls_monitor 335c4762a1bSJed Brown 336c4762a1bSJed Brown test: 337c4762a1bSJed Brown suffix: 21 338864588a7SAlp Dener args: -test_lmvm -tao_max_it 10 -tao_bqnk_mat_type lmvmsymbadbroyden 339c4762a1bSJed Brown 34005579b36STristan Konolige test: 34105579b36STristan Konolige suffix: 22 34205579b36STristan Konolige args: -tao_max_it 1 -tao_converged_reason 34305579b36STristan Konolige 34405579b36STristan Konolige test: 34505579b36STristan Konolige suffix: 23 34605579b36STristan Konolige args: -tao_max_funcs 0 -tao_converged_reason 34705579b36STristan Konolige 34805579b36STristan Konolige test: 34905579b36STristan Konolige suffix: 24 35005579b36STristan Konolige args: -tao_gatol 10 -tao_converged_reason 35105579b36STristan Konolige 35205579b36STristan Konolige test: 35305579b36STristan Konolige suffix: 25 35405579b36STristan Konolige args: -tao_grtol 10 -tao_converged_reason 35505579b36STristan Konolige 35605579b36STristan Konolige test: 35705579b36STristan Konolige suffix: 26 35805579b36STristan Konolige args: -tao_gttol 10 -tao_converged_reason 35905579b36STristan Konolige 36005579b36STristan Konolige test: 36105579b36STristan Konolige suffix: 27 36205579b36STristan Konolige args: -tao_steptol 10 -tao_converged_reason 36305579b36STristan Konolige 36405579b36STristan Konolige test: 36505579b36STristan Konolige suffix: 28 36605579b36STristan Konolige args: -tao_fmin 10 -tao_converged_reason 36705579b36STristan Konolige 368c4762a1bSJed Brown TEST*/ 369