1 /* Program usage: mpiexec -n 1 minsurf1 [-help] [all TAO options] */ 2 3 /* Include "petsctao.h" so we can use TAO solvers. */ 4 #include <petsctao.h> 5 6 static char help[] = 7 "This example demonstrates use of the TAO package to \n\ 8 solve an unconstrained minimization problem. This example is based on a \n\ 9 problem from the MINPACK-2 test suite. Given a rectangular 2-D domain and \n\ 10 boundary values along the edges of the domain, the objective is to find the\n\ 11 surface with the minimal area that satisfies the boundary conditions.\n\ 12 The command line options are:\n\ 13 -mx <xg>, where <xg> = number of grid points in the 1st coordinate direction\n\ 14 -my <yg>, where <yg> = number of grid points in the 2nd coordinate direction\n\ 15 -start <st>, where <st> =0 for zero vector, and an average of the boundary conditions otherwise \n\n"; 16 17 /*T 18 Concepts: TAO^Solving an unconstrained minimization problem 19 Routines: TaoCreate(); TaoSetType(); 20 Routines: TaoSetSolution(); 21 Routines: TaoSetObjectiveAndGradient(); 22 Routines: TaoSetHessian(); TaoSetFromOptions(); 23 Routines: TaoGetKSP(); TaoSolve(); 24 Routines: TaoDestroy(); 25 Processors: 1 26 T*/ 27 28 /* 29 User-defined application context - contains data needed by the 30 application-provided call-back routines, FormFunctionGradient() 31 and FormHessian(). 32 */ 33 typedef struct { 34 PetscInt mx, my; /* discretization in x, y directions */ 35 PetscReal *bottom, *top, *left, *right; /* boundary values */ 36 Mat H; 37 } AppCtx; 38 39 /* -------- User-defined Routines --------- */ 40 41 static PetscErrorCode MSA_BoundaryConditions(AppCtx*); 42 static PetscErrorCode MSA_InitialPoint(AppCtx*,Vec); 43 static PetscErrorCode QuadraticH(AppCtx*,Vec,Mat); 44 PetscErrorCode FormFunctionGradient(Tao,Vec,PetscReal*,Vec,void*); 45 PetscErrorCode FormHessian(Tao,Vec,Mat,Mat,void*); 46 47 int main(int argc, char **argv) 48 { 49 PetscInt N; /* Size of vector */ 50 PetscMPIInt size; /* Number of processors */ 51 Vec x; /* solution, gradient vectors */ 52 KSP ksp; /* PETSc Krylov subspace method */ 53 PetscBool flg; /* A return value when checking for user options */ 54 Tao tao; /* Tao solver context */ 55 AppCtx user; /* user-defined work context */ 56 57 /* Initialize TAO,PETSc */ 58 PetscCall(PetscInitialize(&argc, &argv,(char *)0,help)); 59 60 PetscCallMPI(MPI_Comm_size(MPI_COMM_WORLD,&size)); 61 PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Incorrect number of processors"); 62 63 /* Specify default dimension of the problem */ 64 user.mx = 4; user.my = 4; 65 66 /* Check for any command line arguments that override defaults */ 67 PetscCall(PetscOptionsGetInt(NULL,NULL,"-mx",&user.mx,&flg)); 68 PetscCall(PetscOptionsGetInt(NULL,NULL,"-my",&user.my,&flg)); 69 70 PetscCall(PetscPrintf(PETSC_COMM_SELF,"\n---- Minimum Surface Area Problem -----\n")); 71 PetscCall(PetscPrintf(PETSC_COMM_SELF,"mx: %D my: %D \n\n",user.mx,user.my)); 72 73 /* Calculate any derived values from parameters */ 74 N = user.mx*user.my; 75 76 /* Create TAO solver and set desired solution method */ 77 PetscCall(TaoCreate(PETSC_COMM_SELF,&tao)); 78 PetscCall(TaoSetType(tao,TAOLMVM)); 79 80 /* Initialize minsurf application data structure for use in the function evaluations */ 81 PetscCall(MSA_BoundaryConditions(&user)); /* Application specific routine */ 82 83 /* 84 Create a vector to hold the variables. Compute an initial solution. 85 Set this vector, which will be used by TAO. 86 */ 87 PetscCall(VecCreateSeq(PETSC_COMM_SELF,N,&x)); 88 PetscCall(MSA_InitialPoint(&user,x)); /* Application specific routine */ 89 PetscCall(TaoSetSolution(tao,x)); /* A TAO routine */ 90 91 /* Provide TAO routines for function, gradient, and Hessian evaluation */ 92 PetscCall(TaoSetObjectiveAndGradient(tao,NULL,FormFunctionGradient,(void *)&user)); 93 94 /* Create a matrix data structure to store the Hessian. This structure will be used by TAO */ 95 PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF,N,N,7,NULL,&(user.H))); 96 PetscCall(MatSetOption(user.H,MAT_SYMMETRIC,PETSC_TRUE)); 97 PetscCall(TaoSetHessian(tao,user.H,user.H,FormHessian,(void *)&user)); 98 99 /* Check for any TAO command line options */ 100 PetscCall(TaoSetFromOptions(tao)); 101 102 /* Limit the number of iterations in the KSP linear solver */ 103 PetscCall(TaoGetKSP(tao,&ksp)); 104 if (ksp) { 105 PetscCall(KSPSetTolerances(ksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,user.mx*user.my)); 106 } 107 108 /* SOLVE THE APPLICATION */ 109 PetscCall(TaoSolve(tao)); 110 111 PetscCall(TaoDestroy(&tao)); 112 PetscCall(VecDestroy(&x)); 113 PetscCall(MatDestroy(&user.H)); 114 PetscCall(PetscFree(user.bottom)); 115 PetscCall(PetscFree(user.top)); 116 PetscCall(PetscFree(user.left)); 117 PetscCall(PetscFree(user.right)); 118 119 PetscCall(PetscFinalize()); 120 return 0; 121 } 122 123 /* -------------------------------------------------------------------- */ 124 125 /* FormFunctionGradient - Evaluates function and corresponding gradient. 126 127 Input Parameters: 128 . tao - the Tao context 129 . X - input vector 130 . userCtx - optional user-defined context, as set by TaoSetFunctionGradient() 131 132 Output Parameters: 133 . fcn - the newly evaluated function 134 . G - vector containing the newly evaluated gradient 135 */ 136 PetscErrorCode FormFunctionGradient(Tao tao,Vec X,PetscReal *fcn,Vec G,void *userCtx) 137 { 138 AppCtx *user = (AppCtx *) userCtx; 139 PetscInt i,j,row; 140 PetscInt mx=user->mx, my=user->my; 141 PetscReal rhx=mx+1, rhy=my+1; 142 PetscReal hx=1.0/(mx+1),hy=1.0/(my+1), hydhx=hy/hx, hxdhy=hx/hy, area=0.5*hx*hy, ft=0; 143 PetscReal f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb; 144 PetscReal df1dxc,df2dxc,df3dxc,df4dxc,df5dxc,df6dxc; 145 PetscReal zero=0.0; 146 PetscScalar *g; 147 const PetscScalar *x; 148 149 PetscFunctionBeginUser; 150 PetscCall(VecSet(G, zero)); 151 152 PetscCall(VecGetArrayRead(X,&x)); 153 PetscCall(VecGetArray(G,&g)); 154 155 /* Compute function over the locally owned part of the mesh */ 156 for (j=0; j<my; j++) { 157 for (i=0; i< mx; i++) { 158 row=(j)*mx + (i); 159 xc = x[row]; 160 xlt=xrb=xl=xr=xb=xt=xc; 161 if (i==0) { /* left side */ 162 xl = user->left[j+1]; 163 xlt = user->left[j+2]; 164 } else { 165 xl = x[row-1]; 166 } 167 168 if (j==0) { /* bottom side */ 169 xb = user->bottom[i+1]; 170 xrb = user->bottom[i+2]; 171 } else { 172 xb = x[row-mx]; 173 } 174 175 if (i+1 == mx) { /* right side */ 176 xr = user->right[j+1]; 177 xrb = user->right[j]; 178 } else { 179 xr = x[row+1]; 180 } 181 182 if (j+1==0+my) { /* top side */ 183 xt = user->top[i+1]; 184 xlt = user->top[i]; 185 }else { 186 xt = x[row+mx]; 187 } 188 189 if (i>0 && j+1<my) { 190 xlt = x[row-1+mx]; 191 } 192 if (j>0 && i+1<mx) { 193 xrb = x[row+1-mx]; 194 } 195 196 d1 = (xc-xl); 197 d2 = (xc-xr); 198 d3 = (xc-xt); 199 d4 = (xc-xb); 200 d5 = (xr-xrb); 201 d6 = (xrb-xb); 202 d7 = (xlt-xl); 203 d8 = (xt-xlt); 204 205 df1dxc = d1*hydhx; 206 df2dxc = (d1*hydhx + d4*hxdhy); 207 df3dxc = d3*hxdhy; 208 df4dxc = (d2*hydhx + d3*hxdhy); 209 df5dxc = d2*hydhx; 210 df6dxc = d4*hxdhy; 211 212 d1 *= rhx; 213 d2 *= rhx; 214 d3 *= rhy; 215 d4 *= rhy; 216 d5 *= rhy; 217 d6 *= rhx; 218 d7 *= rhy; 219 d8 *= rhx; 220 221 f1 = PetscSqrtScalar(1.0 + d1*d1 + d7*d7); 222 f2 = PetscSqrtScalar(1.0 + d1*d1 + d4*d4); 223 f3 = PetscSqrtScalar(1.0 + d3*d3 + d8*d8); 224 f4 = PetscSqrtScalar(1.0 + d3*d3 + d2*d2); 225 f5 = PetscSqrtScalar(1.0 + d2*d2 + d5*d5); 226 f6 = PetscSqrtScalar(1.0 + d4*d4 + d6*d6); 227 228 ft = ft + (f2 + f4); 229 230 df1dxc /= f1; 231 df2dxc /= f2; 232 df3dxc /= f3; 233 df4dxc /= f4; 234 df5dxc /= f5; 235 df6dxc /= f6; 236 237 g[row] = (df1dxc+df2dxc+df3dxc+df4dxc+df5dxc+df6dxc)/2.0; 238 } 239 } 240 241 for (j=0; j<my; j++) { /* left side */ 242 d3 = (user->left[j+1] - user->left[j+2])*rhy; 243 d2 = (user->left[j+1] - x[j*mx])*rhx; 244 ft = ft+PetscSqrtScalar(1.0 + d3*d3 + d2*d2); 245 } 246 247 for (i=0; i<mx; i++) { /* bottom */ 248 d2 = (user->bottom[i+1]-user->bottom[i+2])*rhx; 249 d3 = (user->bottom[i+1]-x[i])*rhy; 250 ft = ft+PetscSqrtScalar(1.0 + d3*d3 + d2*d2); 251 } 252 253 for (j=0; j< my; j++) { /* right side */ 254 d1 = (x[(j+1)*mx-1]-user->right[j+1])*rhx; 255 d4 = (user->right[j]-user->right[j+1])*rhy; 256 ft = ft+PetscSqrtScalar(1.0 + d1*d1 + d4*d4); 257 } 258 259 for (i=0; i<mx; i++) { /* top side */ 260 d1 = (x[(my-1)*mx + i] - user->top[i+1])*rhy; 261 d4 = (user->top[i+1] - user->top[i])*rhx; 262 ft = ft+PetscSqrtScalar(1.0 + d1*d1 + d4*d4); 263 } 264 265 /* Bottom left corner */ 266 d1 = (user->left[0]-user->left[1])*rhy; 267 d2 = (user->bottom[0]-user->bottom[1])*rhx; 268 ft += PetscSqrtScalar(1.0 + d1*d1 + d2*d2); 269 270 /* Top right corner */ 271 d1 = (user->right[my+1] - user->right[my])*rhy; 272 d2 = (user->top[mx+1] - user->top[mx])*rhx; 273 ft += PetscSqrtScalar(1.0 + d1*d1 + d2*d2); 274 275 (*fcn)=ft*area; 276 277 /* Restore vectors */ 278 PetscCall(VecRestoreArrayRead(X,&x)); 279 PetscCall(VecRestoreArray(G,&g)); 280 PetscCall(PetscLogFlops(67.0*mx*my)); 281 PetscFunctionReturn(0); 282 } 283 284 /* ------------------------------------------------------------------- */ 285 /* 286 FormHessian - Evaluates the Hessian matrix. 287 288 Input Parameters: 289 . tao - the Tao context 290 . x - input vector 291 . ptr - optional user-defined context, as set by TaoSetHessian() 292 293 Output Parameters: 294 . H - Hessian matrix 295 . Hpre - optionally different preconditioning matrix 296 . flg - flag indicating matrix structure 297 298 */ 299 PetscErrorCode FormHessian(Tao tao,Vec X,Mat H, Mat Hpre, void *ptr) 300 { 301 AppCtx *user = (AppCtx *) ptr; 302 303 PetscFunctionBeginUser; 304 /* Evaluate the Hessian entries*/ 305 PetscCall(QuadraticH(user,X,H)); 306 PetscFunctionReturn(0); 307 } 308 309 /* ------------------------------------------------------------------- */ 310 /* 311 QuadraticH - Evaluates the Hessian matrix. 312 313 Input Parameters: 314 . user - user-defined context, as set by TaoSetHessian() 315 . X - input vector 316 317 Output Parameter: 318 . H - Hessian matrix 319 */ 320 PetscErrorCode QuadraticH(AppCtx *user, Vec X, Mat Hessian) 321 { 322 PetscInt i,j,k,row; 323 PetscInt mx=user->mx, my=user->my; 324 PetscInt col[7]; 325 PetscReal hx=1.0/(mx+1), hy=1.0/(my+1), hydhx=hy/hx, hxdhy=hx/hy; 326 PetscReal rhx=mx+1, rhy=my+1; 327 PetscReal f1,f2,f3,f4,f5,f6,d1,d2,d3,d4,d5,d6,d7,d8,xc,xl,xr,xt,xb,xlt,xrb; 328 PetscReal hl,hr,ht,hb,hc,htl,hbr; 329 const PetscScalar *x; 330 PetscReal v[7]; 331 332 PetscFunctionBeginUser; 333 /* Get pointers to vector data */ 334 PetscCall(VecGetArrayRead(X,&x)); 335 336 /* Initialize matrix entries to zero */ 337 PetscCall(MatZeroEntries(Hessian)); 338 339 /* Set various matrix options */ 340 PetscCall(MatSetOption(Hessian,MAT_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE)); 341 342 /* Compute Hessian over the locally owned part of the mesh */ 343 for (i=0; i< mx; i++) { 344 for (j=0; j<my; j++) { 345 346 row=(j)*mx + (i); 347 348 xc = x[row]; 349 xlt=xrb=xl=xr=xb=xt=xc; 350 351 /* Left side */ 352 if (i==0) { 353 xl= user->left[j+1]; 354 xlt = user->left[j+2]; 355 } else { 356 xl = x[row-1]; 357 } 358 359 if (j==0) { 360 xb=user->bottom[i+1]; 361 xrb = user->bottom[i+2]; 362 } else { 363 xb = x[row-mx]; 364 } 365 366 if (i+1 == mx) { 367 xr=user->right[j+1]; 368 xrb = user->right[j]; 369 } else { 370 xr = x[row+1]; 371 } 372 373 if (j+1==my) { 374 xt=user->top[i+1]; 375 xlt = user->top[i]; 376 }else { 377 xt = x[row+mx]; 378 } 379 380 if (i>0 && j+1<my) { 381 xlt = x[row-1+mx]; 382 } 383 if (j>0 && i+1<mx) { 384 xrb = x[row+1-mx]; 385 } 386 387 d1 = (xc-xl)*rhx; 388 d2 = (xc-xr)*rhx; 389 d3 = (xc-xt)*rhy; 390 d4 = (xc-xb)*rhy; 391 d5 = (xrb-xr)*rhy; 392 d6 = (xrb-xb)*rhx; 393 d7 = (xlt-xl)*rhy; 394 d8 = (xlt-xt)*rhx; 395 396 f1 = PetscSqrtScalar(1.0 + d1*d1 + d7*d7); 397 f2 = PetscSqrtScalar(1.0 + d1*d1 + d4*d4); 398 f3 = PetscSqrtScalar(1.0 + d3*d3 + d8*d8); 399 f4 = PetscSqrtScalar(1.0 + d3*d3 + d2*d2); 400 f5 = PetscSqrtScalar(1.0 + d2*d2 + d5*d5); 401 f6 = PetscSqrtScalar(1.0 + d4*d4 + d6*d6); 402 403 hl = (-hydhx*(1.0+d7*d7)+d1*d7)/(f1*f1*f1)+(-hydhx*(1.0+d4*d4)+d1*d4)/(f2*f2*f2); 404 hr = (-hydhx*(1.0+d5*d5)+d2*d5)/(f5*f5*f5)+(-hydhx*(1.0+d3*d3)+d2*d3)/(f4*f4*f4); 405 ht = (-hxdhy*(1.0+d8*d8)+d3*d8)/(f3*f3*f3)+(-hxdhy*(1.0+d2*d2)+d2*d3)/(f4*f4*f4); 406 hb = (-hxdhy*(1.0+d6*d6)+d4*d6)/(f6*f6*f6)+(-hxdhy*(1.0+d1*d1)+d1*d4)/(f2*f2*f2); 407 408 hbr = -d2*d5/(f5*f5*f5) - d4*d6/(f6*f6*f6); 409 htl = -d1*d7/(f1*f1*f1) - d3*d8/(f3*f3*f3); 410 411 hc = hydhx*(1.0+d7*d7)/(f1*f1*f1) + hxdhy*(1.0+d8*d8)/(f3*f3*f3) + hydhx*(1.0+d5*d5)/(f5*f5*f5) + hxdhy*(1.0+d6*d6)/(f6*f6*f6) + 412 (hxdhy*(1.0+d1*d1)+hydhx*(1.0+d4*d4)-2*d1*d4)/(f2*f2*f2) + (hxdhy*(1.0+d2*d2)+hydhx*(1.0+d3*d3)-2*d2*d3)/(f4*f4*f4); 413 414 hl*=0.5; hr*=0.5; ht*=0.5; hb*=0.5; hbr*=0.5; htl*=0.5; hc*=0.5; 415 416 k=0; 417 if (j>0) { 418 v[k]=hb; col[k]=row - mx; k++; 419 } 420 421 if (j>0 && i < mx -1) { 422 v[k]=hbr; col[k]=row - mx+1; k++; 423 } 424 425 if (i>0) { 426 v[k]= hl; col[k]=row - 1; k++; 427 } 428 429 v[k]= hc; col[k]=row; k++; 430 431 if (i < mx-1) { 432 v[k]= hr; col[k]=row+1; k++; 433 } 434 435 if (i>0 && j < my-1) { 436 v[k]= htl; col[k] = row+mx-1; k++; 437 } 438 439 if (j < my-1) { 440 v[k]= ht; col[k] = row+mx; k++; 441 } 442 443 /* 444 Set matrix values using local numbering, which was defined 445 earlier, in the main routine. 446 */ 447 PetscCall(MatSetValues(Hessian,1,&row,k,col,v,INSERT_VALUES)); 448 } 449 } 450 451 /* Restore vectors */ 452 PetscCall(VecRestoreArrayRead(X,&x)); 453 454 /* Assemble the matrix */ 455 PetscCall(MatAssemblyBegin(Hessian,MAT_FINAL_ASSEMBLY)); 456 PetscCall(MatAssemblyEnd(Hessian,MAT_FINAL_ASSEMBLY)); 457 458 PetscCall(PetscLogFlops(199.0*mx*my)); 459 PetscFunctionReturn(0); 460 } 461 462 /* ------------------------------------------------------------------- */ 463 /* 464 MSA_BoundaryConditions - Calculates the boundary conditions for 465 the region. 466 467 Input Parameter: 468 . user - user-defined application context 469 470 Output Parameter: 471 . user - user-defined application context 472 */ 473 static PetscErrorCode MSA_BoundaryConditions(AppCtx * user) 474 { 475 PetscInt i,j,k,limit=0; 476 PetscInt maxits=5; 477 PetscInt mx=user->mx,my=user->my; 478 PetscInt bsize=0, lsize=0, tsize=0, rsize=0; 479 PetscReal one=1.0, two=2.0, three=3.0, tol=1e-10; 480 PetscReal fnorm,det,hx,hy,xt=0,yt=0; 481 PetscReal u1,u2,nf1,nf2,njac11,njac12,njac21,njac22; 482 PetscReal b=-0.5, t=0.5, l=-0.5, r=0.5; 483 PetscReal *boundary; 484 485 PetscFunctionBeginUser; 486 bsize=mx+2; lsize=my+2; rsize=my+2; tsize=mx+2; 487 488 PetscCall(PetscMalloc1(bsize,&user->bottom)); 489 PetscCall(PetscMalloc1(tsize,&user->top)); 490 PetscCall(PetscMalloc1(lsize,&user->left)); 491 PetscCall(PetscMalloc1(rsize,&user->right)); 492 493 hx= (r-l)/(mx+1); hy=(t-b)/(my+1); 494 495 for (j=0; j<4; j++) { 496 if (j==0) { 497 yt=b; 498 xt=l; 499 limit=bsize; 500 boundary=user->bottom; 501 } else if (j==1) { 502 yt=t; 503 xt=l; 504 limit=tsize; 505 boundary=user->top; 506 } else if (j==2) { 507 yt=b; 508 xt=l; 509 limit=lsize; 510 boundary=user->left; 511 } else { /* if (j==3) */ 512 yt=b; 513 xt=r; 514 limit=rsize; 515 boundary=user->right; 516 } 517 518 for (i=0; i<limit; i++) { 519 u1=xt; 520 u2=-yt; 521 for (k=0; k<maxits; k++) { 522 nf1=u1 + u1*u2*u2 - u1*u1*u1/three-xt; 523 nf2=-u2 - u1*u1*u2 + u2*u2*u2/three-yt; 524 fnorm=PetscSqrtScalar(nf1*nf1+nf2*nf2); 525 if (fnorm <= tol) break; 526 njac11=one+u2*u2-u1*u1; 527 njac12=two*u1*u2; 528 njac21=-two*u1*u2; 529 njac22=-one - u1*u1 + u2*u2; 530 det = njac11*njac22-njac21*njac12; 531 u1 = u1-(njac22*nf1-njac12*nf2)/det; 532 u2 = u2-(njac11*nf2-njac21*nf1)/det; 533 } 534 535 boundary[i]=u1*u1-u2*u2; 536 if (j==0 || j==1) { 537 xt=xt+hx; 538 } else { /* if (j==2 || j==3) */ 539 yt=yt+hy; 540 } 541 } 542 } 543 PetscFunctionReturn(0); 544 } 545 546 /* ------------------------------------------------------------------- */ 547 /* 548 MSA_InitialPoint - Calculates the initial guess in one of three ways. 549 550 Input Parameters: 551 . user - user-defined application context 552 . X - vector for initial guess 553 554 Output Parameters: 555 . X - newly computed initial guess 556 */ 557 static PetscErrorCode MSA_InitialPoint(AppCtx * user, Vec X) 558 { 559 PetscInt start=-1,i,j; 560 PetscReal zero=0.0; 561 PetscBool flg; 562 563 PetscFunctionBeginUser; 564 PetscCall(VecSet(X, zero)); 565 PetscCall(PetscOptionsGetInt(NULL,NULL,"-start",&start,&flg)); 566 567 if (flg && start==0) { /* The zero vector is reasonable */ 568 PetscCall(VecSet(X, zero)); 569 } else { /* Take an average of the boundary conditions */ 570 PetscInt row; 571 PetscInt mx=user->mx,my=user->my; 572 PetscReal *x; 573 574 /* Get pointers to vector data */ 575 PetscCall(VecGetArray(X,&x)); 576 /* Perform local computations */ 577 for (j=0; j<my; j++) { 578 for (i=0; i< mx; i++) { 579 row=(j)*mx + (i); 580 x[row] = (((j+1)*user->bottom[i+1]+(my-j+1)*user->top[i+1])/(my+2)+ ((i+1)*user->left[j+1]+(mx-i+1)*user->right[j+1])/(mx+2))/2.0; 581 } 582 } 583 /* Restore vectors */ 584 PetscCall(VecRestoreArray(X,&x)); 585 } 586 PetscFunctionReturn(0); 587 } 588 589 /*TEST 590 591 build: 592 requires: !complex 593 594 test: 595 args: -tao_smonitor -tao_type nls -mx 10 -my 8 -tao_gatol 1.e-4 596 requires: !single 597 598 test: 599 suffix: 2 600 args: -tao_smonitor -tao_type bmrm -mx 10 -my 8 -tao_gatol 1.e-3 601 requires: !single 602 603 test: 604 suffix: 3 605 args: -tao_smonitor -tao_type lmvm -mx 10 -my 8 -tao_gatol 1.e-3 606 requires: !single 607 608 test: 609 suffix: 4 610 args: -tao_smonitor -tao_type bntr -mx 10 -my 8 -tao_gatol 1.e-4 611 requires: !single 612 613 test: 614 suffix: 5 615 args: -tao_smonitor -tao_type bntl -mx 10 -my 8 -tao_gatol 1.e-4 616 requires: !single 617 618 test: 619 suffix: 6 620 args: -tao_smonitor -tao_type bnls -mx 10 -my 8 -tao_gatol 1.e-4 621 requires: !single 622 623 test: 624 suffix: 7 625 args: -tao_smonitor -tao_type bntr -mx 10 -my 8 -tao_bnk_max_cg_its 3 -tao_gatol 1.e-4 626 requires: !single 627 628 test: 629 suffix: 8 630 args: -tao_smonitor -tao_type bntl -mx 10 -my 8 -tao_bnk_max_cg_its 3 -tao_gatol 1.e-4 631 requires: !single 632 633 test: 634 suffix: 9 635 args: -tao_smonitor -tao_type bnls -mx 10 -my 8 -tao_bnk_max_cg_its 3 -tao_gatol 1.e-4 636 requires: !single 637 638 test: 639 suffix: 10 640 args: -tao_smonitor -tao_type bnls -mx 10 -my 8 -tao_bnk_max_cg_its 3 -tao_gatol 1.e-4 -tao_mf_hessian 641 642 test: 643 suffix: 11 644 args: -tao_smonitor -tao_type bntr -mx 10 -my 8 -tao_gatol 1.e-4 -tao_mf_hessian 645 requires: !single 646 647 test: 648 suffix: 12 649 args: -tao_smonitor -tao_type bntl -mx 10 -my 8 -tao_gatol 1.e-4 -tao_mf_hessian 650 requires: !single 651 652 TEST*/ 653