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