1 /* 2 Include "petsctao.h" so we can use TAO solvers 3 Include "petscdmda.h" so that we can use distributed arrays (DMs) for managing 4 Include "petscksp.h" so we can set KSP type 5 the parallel mesh. 6 */ 7 8 #include <petsctao.h> 9 #include <petscdmda.h> 10 11 static char help[] = "This example demonstrates use of the TAO package to \n\ 12 solve a bound constrained minimization problem. This example is based on \n\ 13 the problem DPJB from the MINPACK-2 test suite. This pressure journal \n\ 14 bearing problem is an example of elliptic variational problem defined over \n\ 15 a two dimensional rectangle. By discretizing the domain into triangular \n\ 16 elements, the pressure surrounding the journal bearing is defined as the \n\ 17 minimum of a quadratic function whose variables are bounded below by zero.\n\ 18 The command line options are:\n\ 19 -mx <xg>, where <xg> = number of grid points in the 1st coordinate direction\n\ 20 -my <yg>, where <yg> = number of grid points in the 2nd coordinate direction\n\ 21 \n"; 22 23 /* 24 User-defined application context - contains data needed by the 25 application-provided call-back routines, FormFunctionGradient(), 26 FormHessian(). 27 */ 28 typedef struct { 29 /* problem parameters */ 30 PetscReal ecc; /* test problem parameter */ 31 PetscReal b; /* A dimension of journal bearing */ 32 PetscInt nx, ny; /* discretization in x, y directions */ 33 34 /* Working space */ 35 DM dm; /* distributed array data structure */ 36 Mat A; /* Quadratic Objective term */ 37 Vec B; /* Linear Objective term */ 38 } AppCtx; 39 40 /* User-defined routines */ 41 static PetscReal p(PetscReal xi, PetscReal ecc); 42 static PetscErrorCode FormFunctionGradient(Tao, Vec, PetscReal *, Vec, void *); 43 static PetscErrorCode FormHessian(Tao, Vec, Mat, Mat, void *); 44 static PetscErrorCode ComputeB(AppCtx *); 45 static PetscErrorCode Monitor(Tao, void *); 46 static PetscErrorCode ConvergenceTest(Tao, void *); 47 48 int main(int argc, char **argv) { 49 PetscInt Nx, Ny; /* number of processors in x- and y- directions */ 50 PetscInt m; /* number of local elements in vectors */ 51 Vec x; /* variables vector */ 52 Vec xl, xu; /* bounds vectors */ 53 PetscReal d1000 = 1000; 54 PetscBool flg, testgetdiag; /* A return variable when checking for user options */ 55 Tao tao; /* Tao solver context */ 56 KSP ksp; 57 AppCtx user; /* user-defined work context */ 58 PetscReal zero = 0.0; /* lower bound on all variables */ 59 60 /* Initialize PETSC and TAO */ 61 PetscFunctionBeginUser; 62 PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 63 64 /* Set the default values for the problem parameters */ 65 user.nx = 50; 66 user.ny = 50; 67 user.ecc = 0.1; 68 user.b = 10.0; 69 testgetdiag = PETSC_FALSE; 70 71 /* Check for any command line arguments that override defaults */ 72 PetscCall(PetscOptionsGetInt(NULL, NULL, "-mx", &user.nx, &flg)); 73 PetscCall(PetscOptionsGetInt(NULL, NULL, "-my", &user.ny, &flg)); 74 PetscCall(PetscOptionsGetReal(NULL, NULL, "-ecc", &user.ecc, &flg)); 75 PetscCall(PetscOptionsGetReal(NULL, NULL, "-b", &user.b, &flg)); 76 PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_getdiagonal", &testgetdiag, NULL)); 77 78 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n---- Journal Bearing Problem SHB-----\n")); 79 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "mx: %" PetscInt_FMT ", my: %" PetscInt_FMT ", ecc: %g \n\n", user.nx, user.ny, (double)user.ecc)); 80 81 /* Let Petsc determine the grid division */ 82 Nx = PETSC_DECIDE; 83 Ny = PETSC_DECIDE; 84 85 /* 86 A two dimensional distributed array will help define this problem, 87 which derives from an elliptic PDE on two dimensional domain. From 88 the distributed array, Create the vectors. 89 */ 90 PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_STAR, user.nx, user.ny, Nx, Ny, 1, 1, NULL, NULL, &user.dm)); 91 PetscCall(DMSetFromOptions(user.dm)); 92 PetscCall(DMSetUp(user.dm)); 93 94 /* 95 Extract global and local vectors from DM; the vector user.B is 96 used solely as work space for the evaluation of the function, 97 gradient, and Hessian. Duplicate for remaining vectors that are 98 the same types. 99 */ 100 PetscCall(DMCreateGlobalVector(user.dm, &x)); /* Solution */ 101 PetscCall(VecDuplicate(x, &user.B)); /* Linear objective */ 102 103 /* Create matrix user.A to store quadratic, Create a local ordering scheme. */ 104 PetscCall(VecGetLocalSize(x, &m)); 105 PetscCall(DMCreateMatrix(user.dm, &user.A)); 106 107 if (testgetdiag) PetscCall(MatSetOperation(user.A, MATOP_GET_DIAGONAL, NULL)); 108 109 /* User defined function -- compute linear term of quadratic */ 110 PetscCall(ComputeB(&user)); 111 112 /* The TAO code begins here */ 113 114 /* 115 Create the optimization solver 116 Suitable methods: TAOGPCG, TAOBQPIP, TAOTRON, TAOBLMVM 117 */ 118 PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao)); 119 PetscCall(TaoSetType(tao, TAOBLMVM)); 120 121 /* Set the initial vector */ 122 PetscCall(VecSet(x, zero)); 123 PetscCall(TaoSetSolution(tao, x)); 124 125 /* Set the user function, gradient, hessian evaluation routines and data structures */ 126 PetscCall(TaoSetObjectiveAndGradient(tao, NULL, FormFunctionGradient, (void *)&user)); 127 128 PetscCall(TaoSetHessian(tao, user.A, user.A, FormHessian, (void *)&user)); 129 130 /* Set a routine that defines the bounds */ 131 PetscCall(VecDuplicate(x, &xl)); 132 PetscCall(VecDuplicate(x, &xu)); 133 PetscCall(VecSet(xl, zero)); 134 PetscCall(VecSet(xu, d1000)); 135 PetscCall(TaoSetVariableBounds(tao, xl, xu)); 136 137 PetscCall(TaoGetKSP(tao, &ksp)); 138 if (ksp) PetscCall(KSPSetType(ksp, KSPCG)); 139 140 PetscCall(PetscOptionsHasName(NULL, NULL, "-testmonitor", &flg)); 141 if (flg) { PetscCall(TaoSetMonitor(tao, Monitor, &user, NULL)); } 142 PetscCall(PetscOptionsHasName(NULL, NULL, "-testconvergence", &flg)); 143 if (flg) { PetscCall(TaoSetConvergenceTest(tao, ConvergenceTest, &user)); } 144 145 /* Check for any tao command line options */ 146 PetscCall(TaoSetFromOptions(tao)); 147 148 /* Solve the bound constrained problem */ 149 PetscCall(TaoSolve(tao)); 150 151 /* Free PETSc data structures */ 152 PetscCall(VecDestroy(&x)); 153 PetscCall(VecDestroy(&xl)); 154 PetscCall(VecDestroy(&xu)); 155 PetscCall(MatDestroy(&user.A)); 156 PetscCall(VecDestroy(&user.B)); 157 158 /* Free TAO data structures */ 159 PetscCall(TaoDestroy(&tao)); 160 PetscCall(DMDestroy(&user.dm)); 161 PetscCall(PetscFinalize()); 162 return 0; 163 } 164 165 static PetscReal p(PetscReal xi, PetscReal ecc) { 166 PetscReal t = 1.0 + ecc * PetscCosScalar(xi); 167 return (t * t * t); 168 } 169 170 PetscErrorCode ComputeB(AppCtx *user) { 171 PetscInt i, j, k; 172 PetscInt nx, ny, xs, xm, gxs, gxm, ys, ym, gys, gym; 173 PetscReal two = 2.0, pi = 4.0 * atan(1.0); 174 PetscReal hx, hy, ehxhy; 175 PetscReal temp, *b; 176 PetscReal ecc = user->ecc; 177 178 PetscFunctionBegin; 179 nx = user->nx; 180 ny = user->ny; 181 hx = two * pi / (nx + 1.0); 182 hy = two * user->b / (ny + 1.0); 183 ehxhy = ecc * hx * hy; 184 185 /* 186 Get local grid boundaries 187 */ 188 PetscCall(DMDAGetCorners(user->dm, &xs, &ys, NULL, &xm, &ym, NULL)); 189 PetscCall(DMDAGetGhostCorners(user->dm, &gxs, &gys, NULL, &gxm, &gym, NULL)); 190 191 /* Compute the linear term in the objective function */ 192 PetscCall(VecGetArray(user->B, &b)); 193 for (i = xs; i < xs + xm; i++) { 194 temp = PetscSinScalar((i + 1) * hx); 195 for (j = ys; j < ys + ym; j++) { 196 k = xm * (j - ys) + (i - xs); 197 b[k] = -ehxhy * temp; 198 } 199 } 200 PetscCall(VecRestoreArray(user->B, &b)); 201 PetscCall(PetscLogFlops(5.0 * xm * ym + 3.0 * xm)); 202 PetscFunctionReturn(0); 203 } 204 205 PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *fcn, Vec G, void *ptr) { 206 AppCtx *user = (AppCtx *)ptr; 207 PetscInt i, j, k, kk; 208 PetscInt col[5], row, nx, ny, xs, xm, gxs, gxm, ys, ym, gys, gym; 209 PetscReal one = 1.0, two = 2.0, six = 6.0, pi = 4.0 * atan(1.0); 210 PetscReal hx, hy, hxhy, hxhx, hyhy; 211 PetscReal xi, v[5]; 212 PetscReal ecc = user->ecc, trule1, trule2, trule3, trule4, trule5, trule6; 213 PetscReal vmiddle, vup, vdown, vleft, vright; 214 PetscReal tt, f1, f2; 215 PetscReal *x, *g, zero = 0.0; 216 Vec localX; 217 218 PetscFunctionBegin; 219 nx = user->nx; 220 ny = user->ny; 221 hx = two * pi / (nx + 1.0); 222 hy = two * user->b / (ny + 1.0); 223 hxhy = hx * hy; 224 hxhx = one / (hx * hx); 225 hyhy = one / (hy * hy); 226 227 PetscCall(DMGetLocalVector(user->dm, &localX)); 228 229 PetscCall(DMGlobalToLocalBegin(user->dm, X, INSERT_VALUES, localX)); 230 PetscCall(DMGlobalToLocalEnd(user->dm, X, INSERT_VALUES, localX)); 231 232 PetscCall(VecSet(G, zero)); 233 /* 234 Get local grid boundaries 235 */ 236 PetscCall(DMDAGetCorners(user->dm, &xs, &ys, NULL, &xm, &ym, NULL)); 237 PetscCall(DMDAGetGhostCorners(user->dm, &gxs, &gys, NULL, &gxm, &gym, NULL)); 238 239 PetscCall(VecGetArray(localX, &x)); 240 PetscCall(VecGetArray(G, &g)); 241 242 for (i = xs; i < xs + xm; i++) { 243 xi = (i + 1) * hx; 244 trule1 = hxhy * (p(xi, ecc) + p(xi + hx, ecc) + p(xi, ecc)) / six; /* L(i,j) */ 245 trule2 = hxhy * (p(xi, ecc) + p(xi - hx, ecc) + p(xi, ecc)) / six; /* U(i,j) */ 246 trule3 = hxhy * (p(xi, ecc) + p(xi + hx, ecc) + p(xi + hx, ecc)) / six; /* U(i+1,j) */ 247 trule4 = hxhy * (p(xi, ecc) + p(xi - hx, ecc) + p(xi - hx, ecc)) / six; /* L(i-1,j) */ 248 trule5 = trule1; /* L(i,j-1) */ 249 trule6 = trule2; /* U(i,j+1) */ 250 251 vdown = -(trule5 + trule2) * hyhy; 252 vleft = -hxhx * (trule2 + trule4); 253 vright = -hxhx * (trule1 + trule3); 254 vup = -hyhy * (trule1 + trule6); 255 vmiddle = (hxhx) * (trule1 + trule2 + trule3 + trule4) + hyhy * (trule1 + trule2 + trule5 + trule6); 256 257 for (j = ys; j < ys + ym; j++) { 258 row = (j - gys) * gxm + (i - gxs); 259 v[0] = 0; 260 v[1] = 0; 261 v[2] = 0; 262 v[3] = 0; 263 v[4] = 0; 264 265 k = 0; 266 if (j > gys) { 267 v[k] = vdown; 268 col[k] = row - gxm; 269 k++; 270 } 271 272 if (i > gxs) { 273 v[k] = vleft; 274 col[k] = row - 1; 275 k++; 276 } 277 278 v[k] = vmiddle; 279 col[k] = row; 280 k++; 281 282 if (i + 1 < gxs + gxm) { 283 v[k] = vright; 284 col[k] = row + 1; 285 k++; 286 } 287 288 if (j + 1 < gys + gym) { 289 v[k] = vup; 290 col[k] = row + gxm; 291 k++; 292 } 293 tt = 0; 294 for (kk = 0; kk < k; kk++) { tt += v[kk] * x[col[kk]]; } 295 row = (j - ys) * xm + (i - xs); 296 g[row] = tt; 297 } 298 } 299 300 PetscCall(VecRestoreArray(localX, &x)); 301 PetscCall(VecRestoreArray(G, &g)); 302 303 PetscCall(DMRestoreLocalVector(user->dm, &localX)); 304 305 PetscCall(VecDot(X, G, &f1)); 306 PetscCall(VecDot(user->B, X, &f2)); 307 PetscCall(VecAXPY(G, one, user->B)); 308 *fcn = f1 / 2.0 + f2; 309 310 PetscCall(PetscLogFlops((91 + 10.0 * ym) * xm)); 311 PetscFunctionReturn(0); 312 } 313 314 /* 315 FormHessian computes the quadratic term in the quadratic objective function 316 Notice that the objective function in this problem is quadratic (therefore a constant 317 hessian). If using a nonquadratic solver, then you might want to reconsider this function 318 */ 319 PetscErrorCode FormHessian(Tao tao, Vec X, Mat hes, Mat Hpre, void *ptr) { 320 AppCtx *user = (AppCtx *)ptr; 321 PetscInt i, j, k; 322 PetscInt col[5], row, nx, ny, xs, xm, gxs, gxm, ys, ym, gys, gym; 323 PetscReal one = 1.0, two = 2.0, six = 6.0, pi = 4.0 * atan(1.0); 324 PetscReal hx, hy, hxhy, hxhx, hyhy; 325 PetscReal xi, v[5]; 326 PetscReal ecc = user->ecc, trule1, trule2, trule3, trule4, trule5, trule6; 327 PetscReal vmiddle, vup, vdown, vleft, vright; 328 PetscBool assembled; 329 330 PetscFunctionBegin; 331 nx = user->nx; 332 ny = user->ny; 333 hx = two * pi / (nx + 1.0); 334 hy = two * user->b / (ny + 1.0); 335 hxhy = hx * hy; 336 hxhx = one / (hx * hx); 337 hyhy = one / (hy * hy); 338 339 /* 340 Get local grid boundaries 341 */ 342 PetscCall(DMDAGetCorners(user->dm, &xs, &ys, NULL, &xm, &ym, NULL)); 343 PetscCall(DMDAGetGhostCorners(user->dm, &gxs, &gys, NULL, &gxm, &gym, NULL)); 344 PetscCall(MatAssembled(hes, &assembled)); 345 if (assembled) PetscCall(MatZeroEntries(hes)); 346 347 for (i = xs; i < xs + xm; i++) { 348 xi = (i + 1) * hx; 349 trule1 = hxhy * (p(xi, ecc) + p(xi + hx, ecc) + p(xi, ecc)) / six; /* L(i,j) */ 350 trule2 = hxhy * (p(xi, ecc) + p(xi - hx, ecc) + p(xi, ecc)) / six; /* U(i,j) */ 351 trule3 = hxhy * (p(xi, ecc) + p(xi + hx, ecc) + p(xi + hx, ecc)) / six; /* U(i+1,j) */ 352 trule4 = hxhy * (p(xi, ecc) + p(xi - hx, ecc) + p(xi - hx, ecc)) / six; /* L(i-1,j) */ 353 trule5 = trule1; /* L(i,j-1) */ 354 trule6 = trule2; /* U(i,j+1) */ 355 356 vdown = -(trule5 + trule2) * hyhy; 357 vleft = -hxhx * (trule2 + trule4); 358 vright = -hxhx * (trule1 + trule3); 359 vup = -hyhy * (trule1 + trule6); 360 vmiddle = (hxhx) * (trule1 + trule2 + trule3 + trule4) + hyhy * (trule1 + trule2 + trule5 + trule6); 361 v[0] = 0; 362 v[1] = 0; 363 v[2] = 0; 364 v[3] = 0; 365 v[4] = 0; 366 367 for (j = ys; j < ys + ym; j++) { 368 row = (j - gys) * gxm + (i - gxs); 369 370 k = 0; 371 if (j > gys) { 372 v[k] = vdown; 373 col[k] = row - gxm; 374 k++; 375 } 376 377 if (i > gxs) { 378 v[k] = vleft; 379 col[k] = row - 1; 380 k++; 381 } 382 383 v[k] = vmiddle; 384 col[k] = row; 385 k++; 386 387 if (i + 1 < gxs + gxm) { 388 v[k] = vright; 389 col[k] = row + 1; 390 k++; 391 } 392 393 if (j + 1 < gys + gym) { 394 v[k] = vup; 395 col[k] = row + gxm; 396 k++; 397 } 398 PetscCall(MatSetValuesLocal(hes, 1, &row, k, col, v, INSERT_VALUES)); 399 } 400 } 401 402 /* 403 Assemble matrix, using the 2-step process: 404 MatAssemblyBegin(), MatAssemblyEnd(). 405 By placing code between these two statements, computations can be 406 done while messages are in transition. 407 */ 408 PetscCall(MatAssemblyBegin(hes, MAT_FINAL_ASSEMBLY)); 409 PetscCall(MatAssemblyEnd(hes, MAT_FINAL_ASSEMBLY)); 410 411 /* 412 Tell the matrix we will never add a new nonzero location to the 413 matrix. If we do it will generate an error. 414 */ 415 PetscCall(MatSetOption(hes, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE)); 416 PetscCall(MatSetOption(hes, MAT_SYMMETRIC, PETSC_TRUE)); 417 418 PetscCall(PetscLogFlops(9.0 * xm * ym + 49.0 * xm)); 419 PetscFunctionReturn(0); 420 } 421 422 PetscErrorCode Monitor(Tao tao, void *ctx) { 423 PetscInt its; 424 PetscReal f, gnorm, cnorm, xdiff; 425 TaoConvergedReason reason; 426 427 PetscFunctionBegin; 428 PetscCall(TaoGetSolutionStatus(tao, &its, &f, &gnorm, &cnorm, &xdiff, &reason)); 429 if (!(its % 5)) { PetscCall(PetscPrintf(PETSC_COMM_WORLD, "iteration=%" PetscInt_FMT "\tf=%g\n", its, (double)f)); } 430 PetscFunctionReturn(0); 431 } 432 433 PetscErrorCode ConvergenceTest(Tao tao, void *ctx) { 434 PetscInt its; 435 PetscReal f, gnorm, cnorm, xdiff; 436 TaoConvergedReason reason; 437 438 PetscFunctionBegin; 439 PetscCall(TaoGetSolutionStatus(tao, &its, &f, &gnorm, &cnorm, &xdiff, &reason)); 440 if (its == 100) { PetscCall(TaoSetConvergedReason(tao, TAO_DIVERGED_MAXITS)); } 441 PetscFunctionReturn(0); 442 } 443 444 /*TEST 445 446 build: 447 requires: !complex 448 449 test: 450 args: -tao_smonitor -mx 8 -my 12 -tao_type tron -tao_gatol 1.e-5 451 requires: !single 452 453 test: 454 suffix: 2 455 nsize: 2 456 args: -tao_smonitor -mx 50 -my 50 -ecc 0.99 -tao_type gpcg -tao_gatol 1.e-5 457 requires: !single 458 459 test: 460 suffix: 3 461 nsize: 2 462 args: -tao_smonitor -mx 10 -my 16 -ecc 0.9 -tao_type bqpip -tao_gatol 1.e-4 463 requires: !single 464 465 test: 466 suffix: 4 467 nsize: 2 468 args: -tao_smonitor -mx 10 -my 16 -ecc 0.9 -tao_type bqpip -tao_gatol 1.e-4 -test_getdiagonal 469 output_file: output/jbearing2_3.out 470 requires: !single 471 472 test: 473 suffix: 5 474 args: -tao_smonitor -mx 8 -my 12 -tao_type bncg -tao_bncg_type gd -tao_gatol 1e-4 475 requires: !single 476 477 test: 478 suffix: 6 479 args: -tao_smonitor -mx 8 -my 12 -tao_type bncg -tao_gatol 1e-4 480 requires: !single 481 482 test: 483 suffix: 7 484 args: -tao_smonitor -mx 8 -my 12 -tao_type bnls -tao_gatol 1e-5 485 requires: !single 486 487 test: 488 suffix: 8 489 args: -tao_smonitor -mx 8 -my 12 -tao_type bntr -tao_gatol 1e-5 490 requires: !single 491 492 test: 493 suffix: 9 494 args: -tao_smonitor -mx 8 -my 12 -tao_type bntl -tao_gatol 1e-5 495 requires: !single 496 497 test: 498 suffix: 10 499 args: -tao_smonitor -mx 8 -my 12 -tao_type bnls -tao_gatol 1e-5 -tao_bnk_max_cg_its 3 500 requires: !single 501 502 test: 503 suffix: 11 504 args: -tao_smonitor -mx 8 -my 12 -tao_type bntr -tao_gatol 1e-5 -tao_bnk_max_cg_its 3 505 requires: !single 506 507 test: 508 suffix: 12 509 args: -tao_smonitor -mx 8 -my 12 -tao_type bntl -tao_gatol 1e-5 -tao_bnk_max_cg_its 3 510 requires: !single 511 512 test: 513 suffix: 13 514 args: -tao_smonitor -mx 8 -my 12 -tao_gatol 1e-4 -tao_type bqnls 515 requires: !single 516 517 test: 518 suffix: 14 519 args: -tao_smonitor -mx 8 -my 12 -tao_gatol 1e-4 -tao_type blmvm 520 requires: !single 521 522 test: 523 suffix: 15 524 args: -tao_smonitor -mx 8 -my 12 -tao_gatol 1e-4 -tao_type bqnkls -tao_bqnk_mat_type lmvmbfgs 525 requires: !single 526 527 test: 528 suffix: 16 529 args: -tao_smonitor -mx 8 -my 12 -tao_gatol 1e-4 -tao_type bqnktr -tao_bqnk_mat_type lmvmsr1 530 requires: !single 531 532 test: 533 suffix: 17 534 args: -tao_smonitor -mx 8 -my 12 -tao_gatol 1e-4 -tao_type bqnls -tao_bqnls_mat_lmvm_scale_type scalar -tao_view 535 requires: !single 536 537 test: 538 suffix: 18 539 args: -tao_smonitor -mx 8 -my 12 -tao_gatol 1e-4 -tao_type bqnls -tao_bqnls_mat_lmvm_scale_type none -tao_view 540 requires: !single 541 542 test: 543 suffix: 19 544 args: -tao_smonitor -mx 8 -my 12 -tao_type bnls -tao_gatol 1e-5 -tao_mf_hessian 545 requires: !single 546 547 test: 548 suffix: 20 549 args: -tao_smonitor -mx 8 -my 12 -tao_type bntr -tao_gatol 1e-5 -tao_mf_hessian 550 requires: !single 551 552 test: 553 suffix: 21 554 args: -tao_smonitor -mx 8 -my 12 -tao_type bntl -tao_gatol 1e-5 -tao_mf_hessian 555 requires: !single 556 TEST*/ 557