1 #include <petscdmda.h> 2 #include <petsctao.h> 3 4 static char help[] = "This example demonstrates use of the TAO package to \n\ 5 solve an unconstrained minimization problem. This example is based on a \n\ 6 problem from the MINPACK-2 test suite. Given a rectangular 2-D domain, \n\ 7 boundary values along the edges of the domain, and a plate represented by \n\ 8 lower boundary conditions, the objective is to find the\n\ 9 surface with the minimal area that satisfies the boundary conditions.\n\ 10 The command line options are:\n\ 11 -mx <xg>, where <xg> = number of grid points in the 1st coordinate direction\n\ 12 -my <yg>, where <yg> = number of grid points in the 2nd coordinate direction\n\ 13 -bmx <bxg>, where <bxg> = number of grid points under plate in 1st direction\n\ 14 -bmy <byg>, where <byg> = number of grid points under plate in 2nd direction\n\ 15 -bheight <ht>, where <ht> = height of the plate\n\ 16 -start <st>, where <st> =0 for zero vector, <st> >0 for random start, and <st> <0 \n\ 17 for an average of the boundary conditions\n\n"; 18 19 /* 20 User-defined application context - contains data needed by the 21 application-provided call-back routines, FormFunctionGradient(), 22 FormHessian(). 23 */ 24 typedef struct { 25 /* problem parameters */ 26 PetscReal bheight; /* Height of plate under the surface */ 27 PetscInt mx, my; /* discretization in x, y directions */ 28 PetscInt bmx, bmy; /* Size of plate under the surface */ 29 Vec Bottom, Top, Left, Right; /* boundary values */ 30 31 /* Working space */ 32 Vec localX, localV; /* ghosted local vector */ 33 DM dm; /* distributed array data structure */ 34 Mat H; 35 } AppCtx; 36 37 /* -------- User-defined Routines --------- */ 38 39 static PetscErrorCode MSA_BoundaryConditions(AppCtx *); 40 static PetscErrorCode MSA_InitialPoint(AppCtx *, Vec); 41 static PetscErrorCode MSA_Plate(Vec, Vec, void *); 42 PetscErrorCode FormFunctionGradient(Tao, Vec, PetscReal *, Vec, void *); 43 PetscErrorCode FormHessian(Tao, Vec, Mat, Mat, void *); 44 45 /* For testing matrix free submatrices */ 46 PetscErrorCode MatrixFreeHessian(Tao, Vec, Mat, Mat, void *); 47 PetscErrorCode MyMatMult(Mat, Vec, Vec); 48 49 int main(int argc, char **argv) 50 { 51 PetscInt Nx, Ny; /* number of processors in x- and y- directions */ 52 PetscInt m, N; /* number of local and global elements in vectors */ 53 Vec x, xl, xu; /* solution vector and bounds*/ 54 PetscBool flg; /* A return variable when checking for user options */ 55 Tao tao; /* Tao solver context */ 56 ISLocalToGlobalMapping isltog; /* local-to-global mapping object */ 57 Mat H_shell; /* to test matrix-free submatrices */ 58 AppCtx user; /* user-defined work context */ 59 60 /* Initialize PETSc, TAO */ 61 PetscFunctionBeginUser; 62 PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 63 64 /* Specify default dimension of the problem */ 65 user.mx = 10; 66 user.my = 10; 67 user.bheight = 0.1; 68 69 /* Check for any command line arguments that override defaults */ 70 PetscCall(PetscOptionsGetInt(NULL, NULL, "-mx", &user.mx, &flg)); 71 PetscCall(PetscOptionsGetInt(NULL, NULL, "-my", &user.my, &flg)); 72 PetscCall(PetscOptionsGetReal(NULL, NULL, "-bheight", &user.bheight, &flg)); 73 74 user.bmx = user.mx / 2; 75 user.bmy = user.my / 2; 76 PetscCall(PetscOptionsGetInt(NULL, NULL, "-bmx", &user.bmx, &flg)); 77 PetscCall(PetscOptionsGetInt(NULL, NULL, "-bmy", &user.bmy, &flg)); 78 79 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n---- Minimum Surface Area With Plate Problem -----\n")); 80 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "mx:%" PetscInt_FMT ", my:%" PetscInt_FMT ", bmx:%" PetscInt_FMT ", bmy:%" PetscInt_FMT ", height:%g\n", user.mx, user.my, user.bmx, user.bmy, (double)user.bheight)); 81 82 /* Calculate any derived values from parameters */ 83 N = user.mx * user.my; 84 85 /* Let Petsc determine the dimensions of the local vectors */ 86 Nx = PETSC_DECIDE; 87 Ny = PETSC_DECIDE; 88 89 /* 90 A two dimensional distributed array will help define this problem, 91 which derives from an elliptic PDE on two dimensional domain. From 92 the distributed array, Create the vectors. 93 */ 94 PetscCall(DMDACreate2d(MPI_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, user.mx, user.my, Nx, Ny, 1, 1, NULL, NULL, &user.dm)); 95 PetscCall(DMSetFromOptions(user.dm)); 96 PetscCall(DMSetUp(user.dm)); 97 /* 98 Extract global and local vectors from DM; The local vectors are 99 used solely as work space for the evaluation of the function, 100 gradient, and Hessian. Duplicate for remaining vectors that are 101 the same types. 102 */ 103 PetscCall(DMCreateGlobalVector(user.dm, &x)); /* Solution */ 104 PetscCall(DMCreateLocalVector(user.dm, &user.localX)); 105 PetscCall(VecDuplicate(user.localX, &user.localV)); 106 107 PetscCall(VecDuplicate(x, &xl)); 108 PetscCall(VecDuplicate(x, &xu)); 109 110 /* The TAO code begins here */ 111 112 /* 113 Create TAO solver and set desired solution method 114 The method must either be TAOTRON or TAOBLMVM 115 If TAOBLMVM is used, then hessian function is not called. 116 */ 117 PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao)); 118 PetscCall(TaoSetType(tao, TAOBLMVM)); 119 120 /* Set initial solution guess; */ 121 PetscCall(MSA_BoundaryConditions(&user)); 122 PetscCall(MSA_InitialPoint(&user, x)); 123 PetscCall(TaoSetSolution(tao, x)); 124 125 /* Set routines for function, gradient and hessian evaluation */ 126 PetscCall(TaoSetObjectiveAndGradient(tao, NULL, FormFunctionGradient, (void *)&user)); 127 128 PetscCall(VecGetLocalSize(x, &m)); 129 PetscCall(MatCreateAIJ(MPI_COMM_WORLD, m, m, N, N, 7, NULL, 3, NULL, &(user.H))); 130 PetscCall(MatSetOption(user.H, MAT_SYMMETRIC, PETSC_TRUE)); 131 132 PetscCall(DMGetLocalToGlobalMapping(user.dm, &isltog)); 133 PetscCall(MatSetLocalToGlobalMapping(user.H, isltog, isltog)); 134 PetscCall(PetscOptionsHasName(NULL, NULL, "-matrixfree", &flg)); 135 if (flg) { 136 PetscCall(MatCreateShell(PETSC_COMM_WORLD, m, m, N, N, (void *)&user, &H_shell)); 137 PetscCall(MatShellSetOperation(H_shell, MATOP_MULT, (void (*)(void))MyMatMult)); 138 PetscCall(MatSetOption(H_shell, MAT_SYMMETRIC, PETSC_TRUE)); 139 PetscCall(TaoSetHessian(tao, H_shell, H_shell, MatrixFreeHessian, (void *)&user)); 140 } else { 141 PetscCall(TaoSetHessian(tao, user.H, user.H, FormHessian, (void *)&user)); 142 } 143 144 /* Set Variable bounds */ 145 PetscCall(MSA_Plate(xl, xu, (void *)&user)); 146 PetscCall(TaoSetVariableBounds(tao, xl, xu)); 147 148 /* Check for any tao command line options */ 149 PetscCall(TaoSetFromOptions(tao)); 150 151 /* SOLVE THE APPLICATION */ 152 PetscCall(TaoSolve(tao)); 153 154 PetscCall(TaoView(tao, PETSC_VIEWER_STDOUT_WORLD)); 155 156 /* Free TAO data structures */ 157 PetscCall(TaoDestroy(&tao)); 158 159 /* Free PETSc data structures */ 160 PetscCall(VecDestroy(&x)); 161 PetscCall(VecDestroy(&xl)); 162 PetscCall(VecDestroy(&xu)); 163 PetscCall(MatDestroy(&user.H)); 164 PetscCall(VecDestroy(&user.localX)); 165 PetscCall(VecDestroy(&user.localV)); 166 PetscCall(VecDestroy(&user.Bottom)); 167 PetscCall(VecDestroy(&user.Top)); 168 PetscCall(VecDestroy(&user.Left)); 169 PetscCall(VecDestroy(&user.Right)); 170 PetscCall(DMDestroy(&user.dm)); 171 if (flg) PetscCall(MatDestroy(&H_shell)); 172 PetscCall(PetscFinalize()); 173 return 0; 174 } 175 176 /* FormFunctionGradient - Evaluates f(x) and gradient g(x). 177 178 Input Parameters: 179 . tao - the Tao context 180 . X - input vector 181 . userCtx - optional user-defined context, as set by TaoSetObjectiveAndGradient() 182 183 Output Parameters: 184 . fcn - the function value 185 . G - vector containing the newly evaluated gradient 186 187 Notes: 188 In this case, we discretize the domain and Create triangles. The 189 surface of each triangle is planar, whose surface area can be easily 190 computed. The total surface area is found by sweeping through the grid 191 and computing the surface area of the two triangles that have their 192 right angle at the grid point. The diagonal line segments on the 193 grid that define the triangles run from top left to lower right. 194 The numbering of points starts at the lower left and runs left to 195 right, then bottom to top. 196 */ 197 PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *fcn, Vec G, void *userCtx) 198 { 199 AppCtx *user = (AppCtx *)userCtx; 200 PetscInt i, j, row; 201 PetscInt mx = user->mx, my = user->my; 202 PetscInt xs, xm, gxs, gxm, ys, ym, gys, gym; 203 PetscReal ft = 0; 204 PetscReal zero = 0.0; 205 PetscReal hx = 1.0 / (mx + 1), hy = 1.0 / (my + 1), hydhx = hy / hx, hxdhy = hx / hy, area = 0.5 * hx * hy; 206 PetscReal rhx = mx + 1, rhy = my + 1; 207 PetscReal f1, f2, f3, f4, f5, f6, d1, d2, d3, d4, d5, d6, d7, d8, xc, xl, xr, xt, xb, xlt, xrb; 208 PetscReal df1dxc, df2dxc, df3dxc, df4dxc, df5dxc, df6dxc; 209 PetscReal *g, *x, *left, *right, *bottom, *top; 210 Vec localX = user->localX, localG = user->localV; 211 212 /* Get local mesh boundaries */ 213 PetscCall(DMDAGetCorners(user->dm, &xs, &ys, NULL, &xm, &ym, NULL)); 214 PetscCall(DMDAGetGhostCorners(user->dm, &gxs, &gys, NULL, &gxm, &gym, NULL)); 215 216 /* Scatter ghost points to local vector */ 217 PetscCall(DMGlobalToLocalBegin(user->dm, X, INSERT_VALUES, localX)); 218 PetscCall(DMGlobalToLocalEnd(user->dm, X, INSERT_VALUES, localX)); 219 220 /* Initialize vector to zero */ 221 PetscCall(VecSet(localG, zero)); 222 223 /* Get pointers to vector data */ 224 PetscCall(VecGetArray(localX, &x)); 225 PetscCall(VecGetArray(localG, &g)); 226 PetscCall(VecGetArray(user->Top, &top)); 227 PetscCall(VecGetArray(user->Bottom, &bottom)); 228 PetscCall(VecGetArray(user->Left, &left)); 229 PetscCall(VecGetArray(user->Right, &right)); 230 231 /* Compute function over the locally owned part of the mesh */ 232 for (j = ys; j < ys + ym; j++) { 233 for (i = xs; i < xs + xm; i++) { 234 row = (j - gys) * gxm + (i - gxs); 235 236 xc = x[row]; 237 xlt = xrb = xl = xr = xb = xt = xc; 238 239 if (i == 0) { /* left side */ 240 xl = left[j - ys + 1]; 241 xlt = left[j - ys + 2]; 242 } else { 243 xl = x[row - 1]; 244 } 245 246 if (j == 0) { /* bottom side */ 247 xb = bottom[i - xs + 1]; 248 xrb = bottom[i - xs + 2]; 249 } else { 250 xb = x[row - gxm]; 251 } 252 253 if (i + 1 == gxs + gxm) { /* right side */ 254 xr = right[j - ys + 1]; 255 xrb = right[j - ys]; 256 } else { 257 xr = x[row + 1]; 258 } 259 260 if (j + 1 == gys + gym) { /* top side */ 261 xt = top[i - xs + 1]; 262 xlt = top[i - xs]; 263 } else { 264 xt = x[row + gxm]; 265 } 266 267 if (i > gxs && j + 1 < gys + gym) xlt = x[row - 1 + gxm]; 268 if (j > gys && i + 1 < gxs + gxm) xrb = x[row + 1 - gxm]; 269 270 d1 = (xc - xl); 271 d2 = (xc - xr); 272 d3 = (xc - xt); 273 d4 = (xc - xb); 274 d5 = (xr - xrb); 275 d6 = (xrb - xb); 276 d7 = (xlt - xl); 277 d8 = (xt - xlt); 278 279 df1dxc = d1 * hydhx; 280 df2dxc = (d1 * hydhx + d4 * hxdhy); 281 df3dxc = d3 * hxdhy; 282 df4dxc = (d2 * hydhx + d3 * hxdhy); 283 df5dxc = d2 * hydhx; 284 df6dxc = d4 * hxdhy; 285 286 d1 *= rhx; 287 d2 *= rhx; 288 d3 *= rhy; 289 d4 *= rhy; 290 d5 *= rhy; 291 d6 *= rhx; 292 d7 *= rhy; 293 d8 *= rhx; 294 295 f1 = PetscSqrtScalar(1.0 + d1 * d1 + d7 * d7); 296 f2 = PetscSqrtScalar(1.0 + d1 * d1 + d4 * d4); 297 f3 = PetscSqrtScalar(1.0 + d3 * d3 + d8 * d8); 298 f4 = PetscSqrtScalar(1.0 + d3 * d3 + d2 * d2); 299 f5 = PetscSqrtScalar(1.0 + d2 * d2 + d5 * d5); 300 f6 = PetscSqrtScalar(1.0 + d4 * d4 + d6 * d6); 301 302 ft = ft + (f2 + f4); 303 304 df1dxc /= f1; 305 df2dxc /= f2; 306 df3dxc /= f3; 307 df4dxc /= f4; 308 df5dxc /= f5; 309 df6dxc /= f6; 310 311 g[row] = (df1dxc + df2dxc + df3dxc + df4dxc + df5dxc + df6dxc) * 0.5; 312 } 313 } 314 315 /* Compute triangular areas along the border of the domain. */ 316 if (xs == 0) { /* left side */ 317 for (j = ys; j < ys + ym; j++) { 318 d3 = (left[j - ys + 1] - left[j - ys + 2]) * rhy; 319 d2 = (left[j - ys + 1] - x[(j - gys) * gxm]) * rhx; 320 ft = ft + PetscSqrtScalar(1.0 + d3 * d3 + d2 * d2); 321 } 322 } 323 if (ys == 0) { /* bottom side */ 324 for (i = xs; i < xs + xm; i++) { 325 d2 = (bottom[i + 1 - xs] - bottom[i - xs + 2]) * rhx; 326 d3 = (bottom[i - xs + 1] - x[i - gxs]) * rhy; 327 ft = ft + PetscSqrtScalar(1.0 + d3 * d3 + d2 * d2); 328 } 329 } 330 331 if (xs + xm == mx) { /* right side */ 332 for (j = ys; j < ys + ym; j++) { 333 d1 = (x[(j + 1 - gys) * gxm - 1] - right[j - ys + 1]) * rhx; 334 d4 = (right[j - ys] - right[j - ys + 1]) * rhy; 335 ft = ft + PetscSqrtScalar(1.0 + d1 * d1 + d4 * d4); 336 } 337 } 338 if (ys + ym == my) { /* top side */ 339 for (i = xs; i < xs + xm; i++) { 340 d1 = (x[(gym - 1) * gxm + i - gxs] - top[i - xs + 1]) * rhy; 341 d4 = (top[i - xs + 1] - top[i - xs]) * rhx; 342 ft = ft + PetscSqrtScalar(1.0 + d1 * d1 + d4 * d4); 343 } 344 } 345 346 if (ys == 0 && xs == 0) { 347 d1 = (left[0] - left[1]) * rhy; 348 d2 = (bottom[0] - bottom[1]) * rhx; 349 ft += PetscSqrtScalar(1.0 + d1 * d1 + d2 * d2); 350 } 351 if (ys + ym == my && xs + xm == mx) { 352 d1 = (right[ym + 1] - right[ym]) * rhy; 353 d2 = (top[xm + 1] - top[xm]) * rhx; 354 ft += PetscSqrtScalar(1.0 + d1 * d1 + d2 * d2); 355 } 356 357 ft = ft * area; 358 PetscCallMPI(MPI_Allreduce(&ft, fcn, 1, MPIU_REAL, MPIU_SUM, MPI_COMM_WORLD)); 359 360 /* Restore vectors */ 361 PetscCall(VecRestoreArray(localX, &x)); 362 PetscCall(VecRestoreArray(localG, &g)); 363 PetscCall(VecRestoreArray(user->Left, &left)); 364 PetscCall(VecRestoreArray(user->Top, &top)); 365 PetscCall(VecRestoreArray(user->Bottom, &bottom)); 366 PetscCall(VecRestoreArray(user->Right, &right)); 367 368 /* Scatter values to global vector */ 369 PetscCall(DMLocalToGlobalBegin(user->dm, localG, INSERT_VALUES, G)); 370 PetscCall(DMLocalToGlobalEnd(user->dm, localG, INSERT_VALUES, G)); 371 372 PetscCall(PetscLogFlops(70.0 * xm * ym)); 373 374 return 0; 375 } 376 377 /* ------------------------------------------------------------------- */ 378 /* 379 FormHessian - Evaluates Hessian matrix. 380 381 Input Parameters: 382 . tao - the Tao context 383 . x - input vector 384 . ptr - optional user-defined context, as set by TaoSetHessian() 385 386 Output Parameters: 387 . A - Hessian matrix 388 . B - optionally different preconditioning matrix 389 390 Notes: 391 Due to mesh point reordering with DMs, we must always work 392 with the local mesh points, and then transform them to the new 393 global numbering with the local-to-global mapping. We cannot work 394 directly with the global numbers for the original uniprocessor mesh! 395 396 Two methods are available for imposing this transformation 397 when setting matrix entries: 398 (A) MatSetValuesLocal(), using the local ordering (including 399 ghost points!) 400 - Do the following two steps once, before calling TaoSolve() 401 - Use DMGetISLocalToGlobalMapping() to extract the 402 local-to-global map from the DM 403 - Associate this map with the matrix by calling 404 MatSetLocalToGlobalMapping() 405 - Then set matrix entries using the local ordering 406 by calling MatSetValuesLocal() 407 (B) MatSetValues(), using the global ordering 408 - Use DMGetGlobalIndices() to extract the local-to-global map 409 - Then apply this map explicitly yourself 410 - Set matrix entries using the global ordering by calling 411 MatSetValues() 412 Option (A) seems cleaner/easier in many cases, and is the procedure 413 used in this example. 414 */ 415 PetscErrorCode FormHessian(Tao tao, Vec X, Mat Hptr, Mat Hessian, void *ptr) 416 { 417 AppCtx *user = (AppCtx *)ptr; 418 PetscInt i, j, k, row; 419 PetscInt mx = user->mx, my = user->my; 420 PetscInt xs, xm, gxs, gxm, ys, ym, gys, gym, col[7]; 421 PetscReal hx = 1.0 / (mx + 1), hy = 1.0 / (my + 1), hydhx = hy / hx, hxdhy = hx / hy; 422 PetscReal rhx = mx + 1, rhy = my + 1; 423 PetscReal f1, f2, f3, f4, f5, f6, d1, d2, d3, d4, d5, d6, d7, d8, xc, xl, xr, xt, xb, xlt, xrb; 424 PetscReal hl, hr, ht, hb, hc, htl, hbr; 425 PetscReal *x, *left, *right, *bottom, *top; 426 PetscReal v[7]; 427 Vec localX = user->localX; 428 PetscBool assembled; 429 430 /* Set various matrix options */ 431 PetscCall(MatSetOption(Hessian, MAT_IGNORE_OFF_PROC_ENTRIES, PETSC_TRUE)); 432 433 /* Initialize matrix entries to zero */ 434 PetscCall(MatAssembled(Hessian, &assembled)); 435 if (assembled) PetscCall(MatZeroEntries(Hessian)); 436 437 /* Get local mesh boundaries */ 438 PetscCall(DMDAGetCorners(user->dm, &xs, &ys, NULL, &xm, &ym, NULL)); 439 PetscCall(DMDAGetGhostCorners(user->dm, &gxs, &gys, NULL, &gxm, &gym, NULL)); 440 441 /* Scatter ghost points to local vector */ 442 PetscCall(DMGlobalToLocalBegin(user->dm, X, INSERT_VALUES, localX)); 443 PetscCall(DMGlobalToLocalEnd(user->dm, X, INSERT_VALUES, localX)); 444 445 /* Get pointers to vector data */ 446 PetscCall(VecGetArray(localX, &x)); 447 PetscCall(VecGetArray(user->Top, &top)); 448 PetscCall(VecGetArray(user->Bottom, &bottom)); 449 PetscCall(VecGetArray(user->Left, &left)); 450 PetscCall(VecGetArray(user->Right, &right)); 451 452 /* Compute Hessian over the locally owned part of the mesh */ 453 454 for (i = xs; i < xs + xm; i++) { 455 for (j = ys; j < ys + ym; j++) { 456 row = (j - gys) * gxm + (i - gxs); 457 458 xc = x[row]; 459 xlt = xrb = xl = xr = xb = xt = xc; 460 461 /* Left side */ 462 if (i == gxs) { 463 xl = left[j - ys + 1]; 464 xlt = left[j - ys + 2]; 465 } else { 466 xl = x[row - 1]; 467 } 468 469 if (j == gys) { 470 xb = bottom[i - xs + 1]; 471 xrb = bottom[i - xs + 2]; 472 } else { 473 xb = x[row - gxm]; 474 } 475 476 if (i + 1 == gxs + gxm) { 477 xr = right[j - ys + 1]; 478 xrb = right[j - ys]; 479 } else { 480 xr = x[row + 1]; 481 } 482 483 if (j + 1 == gys + gym) { 484 xt = top[i - xs + 1]; 485 xlt = top[i - xs]; 486 } else { 487 xt = x[row + gxm]; 488 } 489 490 if (i > gxs && j + 1 < gys + gym) xlt = x[row - 1 + gxm]; 491 if (j > gys && i + 1 < gxs + gxm) xrb = x[row + 1 - gxm]; 492 493 d1 = (xc - xl) * rhx; 494 d2 = (xc - xr) * rhx; 495 d3 = (xc - xt) * rhy; 496 d4 = (xc - xb) * rhy; 497 d5 = (xrb - xr) * rhy; 498 d6 = (xrb - xb) * rhx; 499 d7 = (xlt - xl) * rhy; 500 d8 = (xlt - xt) * rhx; 501 502 f1 = PetscSqrtScalar(1.0 + d1 * d1 + d7 * d7); 503 f2 = PetscSqrtScalar(1.0 + d1 * d1 + d4 * d4); 504 f3 = PetscSqrtScalar(1.0 + d3 * d3 + d8 * d8); 505 f4 = PetscSqrtScalar(1.0 + d3 * d3 + d2 * d2); 506 f5 = PetscSqrtScalar(1.0 + d2 * d2 + d5 * d5); 507 f6 = PetscSqrtScalar(1.0 + d4 * d4 + d6 * d6); 508 509 hl = (-hydhx * (1.0 + d7 * d7) + d1 * d7) / (f1 * f1 * f1) + (-hydhx * (1.0 + d4 * d4) + d1 * d4) / (f2 * f2 * f2); 510 hr = (-hydhx * (1.0 + d5 * d5) + d2 * d5) / (f5 * f5 * f5) + (-hydhx * (1.0 + d3 * d3) + d2 * d3) / (f4 * f4 * f4); 511 ht = (-hxdhy * (1.0 + d8 * d8) + d3 * d8) / (f3 * f3 * f3) + (-hxdhy * (1.0 + d2 * d2) + d2 * d3) / (f4 * f4 * f4); 512 hb = (-hxdhy * (1.0 + d6 * d6) + d4 * d6) / (f6 * f6 * f6) + (-hxdhy * (1.0 + d1 * d1) + d1 * d4) / (f2 * f2 * f2); 513 514 hbr = -d2 * d5 / (f5 * f5 * f5) - d4 * d6 / (f6 * f6 * f6); 515 htl = -d1 * d7 / (f1 * f1 * f1) - d3 * d8 / (f3 * f3 * f3); 516 517 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) + (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); 518 519 hl *= 0.5; 520 hr *= 0.5; 521 ht *= 0.5; 522 hb *= 0.5; 523 hbr *= 0.5; 524 htl *= 0.5; 525 hc *= 0.5; 526 527 k = 0; 528 if (j > 0) { 529 v[k] = hb; 530 col[k] = row - gxm; 531 k++; 532 } 533 534 if (j > 0 && i < mx - 1) { 535 v[k] = hbr; 536 col[k] = row - gxm + 1; 537 k++; 538 } 539 540 if (i > 0) { 541 v[k] = hl; 542 col[k] = row - 1; 543 k++; 544 } 545 546 v[k] = hc; 547 col[k] = row; 548 k++; 549 550 if (i < mx - 1) { 551 v[k] = hr; 552 col[k] = row + 1; 553 k++; 554 } 555 556 if (i > 0 && j < my - 1) { 557 v[k] = htl; 558 col[k] = row + gxm - 1; 559 k++; 560 } 561 562 if (j < my - 1) { 563 v[k] = ht; 564 col[k] = row + gxm; 565 k++; 566 } 567 568 /* 569 Set matrix values using local numbering, which was defined 570 earlier, in the main routine. 571 */ 572 PetscCall(MatSetValuesLocal(Hessian, 1, &row, k, col, v, INSERT_VALUES)); 573 } 574 } 575 576 /* Restore vectors */ 577 PetscCall(VecRestoreArray(localX, &x)); 578 PetscCall(VecRestoreArray(user->Left, &left)); 579 PetscCall(VecRestoreArray(user->Top, &top)); 580 PetscCall(VecRestoreArray(user->Bottom, &bottom)); 581 PetscCall(VecRestoreArray(user->Right, &right)); 582 583 /* Assemble the matrix */ 584 PetscCall(MatAssemblyBegin(Hessian, MAT_FINAL_ASSEMBLY)); 585 PetscCall(MatAssemblyEnd(Hessian, MAT_FINAL_ASSEMBLY)); 586 587 PetscCall(PetscLogFlops(199.0 * xm * ym)); 588 return 0; 589 } 590 591 /* ------------------------------------------------------------------- */ 592 /* 593 MSA_BoundaryConditions - Calculates the boundary conditions for 594 the region. 595 596 Input Parameter: 597 . user - user-defined application context 598 599 Output Parameter: 600 . user - user-defined application context 601 */ 602 static PetscErrorCode MSA_BoundaryConditions(AppCtx *user) 603 { 604 PetscInt i, j, k, maxits = 5, limit = 0; 605 PetscInt xs, ys, xm, ym, gxs, gys, gxm, gym; 606 PetscInt mx = user->mx, my = user->my; 607 PetscInt bsize = 0, lsize = 0, tsize = 0, rsize = 0; 608 PetscReal one = 1.0, two = 2.0, three = 3.0, scl = 1.0, tol = 1e-10; 609 PetscReal fnorm, det, hx, hy, xt = 0, yt = 0; 610 PetscReal u1, u2, nf1, nf2, njac11, njac12, njac21, njac22; 611 PetscReal b = -0.5, t = 0.5, l = -0.5, r = 0.5; 612 PetscReal *boundary; 613 PetscBool flg; 614 Vec Bottom, Top, Right, Left; 615 616 /* Get local mesh boundaries */ 617 PetscCall(DMDAGetCorners(user->dm, &xs, &ys, NULL, &xm, &ym, NULL)); 618 PetscCall(DMDAGetGhostCorners(user->dm, &gxs, &gys, NULL, &gxm, &gym, NULL)); 619 620 bsize = xm + 2; 621 lsize = ym + 2; 622 rsize = ym + 2; 623 tsize = xm + 2; 624 625 PetscCall(VecCreateMPI(MPI_COMM_WORLD, bsize, PETSC_DECIDE, &Bottom)); 626 PetscCall(VecCreateMPI(MPI_COMM_WORLD, tsize, PETSC_DECIDE, &Top)); 627 PetscCall(VecCreateMPI(MPI_COMM_WORLD, lsize, PETSC_DECIDE, &Left)); 628 PetscCall(VecCreateMPI(MPI_COMM_WORLD, rsize, PETSC_DECIDE, &Right)); 629 630 user->Top = Top; 631 user->Left = Left; 632 user->Bottom = Bottom; 633 user->Right = Right; 634 635 hx = (r - l) / (mx + 1); 636 hy = (t - b) / (my + 1); 637 638 for (j = 0; j < 4; j++) { 639 if (j == 0) { 640 yt = b; 641 xt = l + hx * xs; 642 limit = bsize; 643 VecGetArray(Bottom, &boundary); 644 } else if (j == 1) { 645 yt = t; 646 xt = l + hx * xs; 647 limit = tsize; 648 VecGetArray(Top, &boundary); 649 } else if (j == 2) { 650 yt = b + hy * ys; 651 xt = l; 652 limit = lsize; 653 VecGetArray(Left, &boundary); 654 } else if (j == 3) { 655 yt = b + hy * ys; 656 xt = r; 657 limit = rsize; 658 VecGetArray(Right, &boundary); 659 } 660 661 for (i = 0; i < limit; i++) { 662 u1 = xt; 663 u2 = -yt; 664 for (k = 0; k < maxits; k++) { 665 nf1 = u1 + u1 * u2 * u2 - u1 * u1 * u1 / three - xt; 666 nf2 = -u2 - u1 * u1 * u2 + u2 * u2 * u2 / three - yt; 667 fnorm = PetscSqrtScalar(nf1 * nf1 + nf2 * nf2); 668 if (fnorm <= tol) break; 669 njac11 = one + u2 * u2 - u1 * u1; 670 njac12 = two * u1 * u2; 671 njac21 = -two * u1 * u2; 672 njac22 = -one - u1 * u1 + u2 * u2; 673 det = njac11 * njac22 - njac21 * njac12; 674 u1 = u1 - (njac22 * nf1 - njac12 * nf2) / det; 675 u2 = u2 - (njac11 * nf2 - njac21 * nf1) / det; 676 } 677 678 boundary[i] = u1 * u1 - u2 * u2; 679 if (j == 0 || j == 1) { 680 xt = xt + hx; 681 } else if (j == 2 || j == 3) { 682 yt = yt + hy; 683 } 684 } 685 if (j == 0) { 686 PetscCall(VecRestoreArray(Bottom, &boundary)); 687 } else if (j == 1) { 688 PetscCall(VecRestoreArray(Top, &boundary)); 689 } else if (j == 2) { 690 PetscCall(VecRestoreArray(Left, &boundary)); 691 } else if (j == 3) { 692 PetscCall(VecRestoreArray(Right, &boundary)); 693 } 694 } 695 696 /* Scale the boundary if desired */ 697 698 PetscCall(PetscOptionsGetReal(NULL, NULL, "-bottom", &scl, &flg)); 699 if (flg) PetscCall(VecScale(Bottom, scl)); 700 PetscCall(PetscOptionsGetReal(NULL, NULL, "-top", &scl, &flg)); 701 if (flg) PetscCall(VecScale(Top, scl)); 702 PetscCall(PetscOptionsGetReal(NULL, NULL, "-right", &scl, &flg)); 703 if (flg) PetscCall(VecScale(Right, scl)); 704 705 PetscCall(PetscOptionsGetReal(NULL, NULL, "-left", &scl, &flg)); 706 if (flg) PetscCall(VecScale(Left, scl)); 707 return 0; 708 } 709 710 /* ------------------------------------------------------------------- */ 711 /* 712 MSA_Plate - Calculates an obstacle for surface to stretch over. 713 714 Input Parameter: 715 . user - user-defined application context 716 717 Output Parameter: 718 . user - user-defined application context 719 */ 720 static PetscErrorCode MSA_Plate(Vec XL, Vec XU, void *ctx) 721 { 722 AppCtx *user = (AppCtx *)ctx; 723 PetscInt i, j, row; 724 PetscInt xs, ys, xm, ym; 725 PetscInt mx = user->mx, my = user->my, bmy, bmx; 726 PetscReal t1, t2, t3; 727 PetscReal *xl, lb = PETSC_NINFINITY, ub = PETSC_INFINITY; 728 PetscBool cylinder; 729 730 user->bmy = PetscMax(0, user->bmy); 731 user->bmy = PetscMin(my, user->bmy); 732 user->bmx = PetscMax(0, user->bmx); 733 user->bmx = PetscMin(mx, user->bmx); 734 bmy = user->bmy; 735 bmx = user->bmx; 736 737 PetscCall(DMDAGetCorners(user->dm, &xs, &ys, NULL, &xm, &ym, NULL)); 738 739 PetscCall(VecSet(XL, lb)); 740 PetscCall(VecSet(XU, ub)); 741 742 PetscCall(VecGetArray(XL, &xl)); 743 744 PetscCall(PetscOptionsHasName(NULL, NULL, "-cylinder", &cylinder)); 745 /* Compute the optional lower box */ 746 if (cylinder) { 747 for (i = xs; i < xs + xm; i++) { 748 for (j = ys; j < ys + ym; j++) { 749 row = (j - ys) * xm + (i - xs); 750 t1 = (2.0 * i - mx) * bmy; 751 t2 = (2.0 * j - my) * bmx; 752 t3 = bmx * bmx * bmy * bmy; 753 if (t1 * t1 + t2 * t2 <= t3) xl[row] = user->bheight; 754 } 755 } 756 } else { 757 /* Compute the optional lower box */ 758 for (i = xs; i < xs + xm; i++) { 759 for (j = ys; j < ys + ym; j++) { 760 row = (j - ys) * xm + (i - xs); 761 if (i >= (mx - bmx) / 2 && i < mx - (mx - bmx) / 2 && j >= (my - bmy) / 2 && j < my - (my - bmy) / 2) xl[row] = user->bheight; 762 } 763 } 764 } 765 PetscCall(VecRestoreArray(XL, &xl)); 766 767 return 0; 768 } 769 770 /* ------------------------------------------------------------------- */ 771 /* 772 MSA_InitialPoint - Calculates the initial guess in one of three ways. 773 774 Input Parameters: 775 . user - user-defined application context 776 . X - vector for initial guess 777 778 Output Parameters: 779 . X - newly computed initial guess 780 */ 781 static PetscErrorCode MSA_InitialPoint(AppCtx *user, Vec X) 782 { 783 PetscInt start = -1, i, j; 784 PetscReal zero = 0.0; 785 PetscBool flg; 786 787 PetscCall(PetscOptionsGetInt(NULL, NULL, "-start", &start, &flg)); 788 if (flg && start == 0) { /* The zero vector is reasonable */ 789 PetscCall(VecSet(X, zero)); 790 } else if (flg && start > 0) { /* Try a random start between -0.5 and 0.5 */ 791 PetscRandom rctx; 792 PetscReal np5 = -0.5; 793 794 PetscCall(PetscRandomCreate(MPI_COMM_WORLD, &rctx)); 795 for (i = 0; i < start; i++) PetscCall(VecSetRandom(X, rctx)); 796 PetscCall(PetscRandomDestroy(&rctx)); 797 PetscCall(VecShift(X, np5)); 798 799 } else { /* Take an average of the boundary conditions */ 800 801 PetscInt row, xs, xm, gxs, gxm, ys, ym, gys, gym; 802 PetscInt mx = user->mx, my = user->my; 803 PetscReal *x, *left, *right, *bottom, *top; 804 Vec localX = user->localX; 805 806 /* Get local mesh boundaries */ 807 PetscCall(DMDAGetCorners(user->dm, &xs, &ys, NULL, &xm, &ym, NULL)); 808 PetscCall(DMDAGetGhostCorners(user->dm, &gxs, &gys, NULL, &gxm, &gym, NULL)); 809 810 /* Get pointers to vector data */ 811 PetscCall(VecGetArray(user->Top, &top)); 812 PetscCall(VecGetArray(user->Bottom, &bottom)); 813 PetscCall(VecGetArray(user->Left, &left)); 814 PetscCall(VecGetArray(user->Right, &right)); 815 816 PetscCall(VecGetArray(localX, &x)); 817 /* Perform local computations */ 818 for (j = ys; j < ys + ym; j++) { 819 for (i = xs; i < xs + xm; i++) { 820 row = (j - gys) * gxm + (i - gxs); 821 x[row] = ((j + 1) * bottom[i - xs + 1] / my + (my - j + 1) * top[i - xs + 1] / (my + 2) + (i + 1) * left[j - ys + 1] / mx + (mx - i + 1) * right[j - ys + 1] / (mx + 2)) / 2.0; 822 } 823 } 824 825 /* Restore vectors */ 826 PetscCall(VecRestoreArray(localX, &x)); 827 828 PetscCall(VecRestoreArray(user->Left, &left)); 829 PetscCall(VecRestoreArray(user->Top, &top)); 830 PetscCall(VecRestoreArray(user->Bottom, &bottom)); 831 PetscCall(VecRestoreArray(user->Right, &right)); 832 833 /* Scatter values into global vector */ 834 PetscCall(DMLocalToGlobalBegin(user->dm, localX, INSERT_VALUES, X)); 835 PetscCall(DMLocalToGlobalEnd(user->dm, localX, INSERT_VALUES, X)); 836 } 837 return 0; 838 } 839 840 /* For testing matrix free submatrices */ 841 PetscErrorCode MatrixFreeHessian(Tao tao, Vec x, Mat H, Mat Hpre, void *ptr) 842 { 843 AppCtx *user = (AppCtx *)ptr; 844 PetscFunctionBegin; 845 PetscCall(FormHessian(tao, x, user->H, user->H, ptr)); 846 PetscFunctionReturn(0); 847 } 848 PetscErrorCode MyMatMult(Mat H_shell, Vec X, Vec Y) 849 { 850 void *ptr; 851 AppCtx *user; 852 PetscFunctionBegin; 853 PetscCall(MatShellGetContext(H_shell, &ptr)); 854 user = (AppCtx *)ptr; 855 PetscCall(MatMult(user->H, X, Y)); 856 PetscFunctionReturn(0); 857 } 858 859 /*TEST 860 861 build: 862 requires: !complex 863 864 test: 865 args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type tron -tao_gatol 1.e-5 866 requires: !single 867 868 test: 869 suffix: 2 870 nsize: 2 871 args: -tao_smonitor -mx 8 -my 8 -bmx 2 -bmy 5 -bheight 0.3 -tao_type blmvm -tao_gatol 1.e-4 872 requires: !single 873 874 test: 875 suffix: 3 876 nsize: 3 877 args: -tao_smonitor -mx 8 -my 12 -bmx 4 -bmy 10 -bheight 0.1 -tao_type tron -tao_gatol 1.e-5 878 requires: !single 879 880 test: 881 suffix: 4 882 nsize: 3 883 args: -tao_smonitor -mx 8 -my 12 -bmx 4 -bmy 10 -bheight 0.1 -tao_subset_type mask -tao_type tron -tao_gatol 1.e-5 884 requires: !single 885 886 test: 887 suffix: 5 888 nsize: 3 889 args: -tao_smonitor -mx 8 -my 12 -bmx 4 -bmy 10 -bheight 0.1 -tao_subset_type matrixfree -matrixfree -pc_type none -tao_type tron -tao_gatol 1.e-5 890 requires: !single 891 892 test: 893 suffix: 6 894 nsize: 3 895 args: -tao_smonitor -mx 8 -my 12 -bmx 4 -bmy 10 -bheight 0.1 -tao_subset_type matrixfree -matrixfree -tao_type blmvm -tao_gatol 1.e-4 896 requires: !single 897 898 test: 899 suffix: 7 900 nsize: 3 901 args: -tao_smonitor -mx 8 -my 12 -bmx 4 -bmy 10 -bheight 0.1 -tao_subset_type matrixfree -pc_type none -tao_type gpcg -tao_gatol 1.e-5 902 requires: !single 903 904 test: 905 suffix: 8 906 args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bncg -tao_bncg_type gd -tao_gatol 1e-4 907 requires: !single 908 909 test: 910 suffix: 9 911 args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bncg -tao_gatol 1e-4 912 requires: !single 913 914 test: 915 suffix: 10 916 args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bnls -tao_gatol 1e-5 917 requires: !single 918 919 test: 920 suffix: 11 921 args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bntr -tao_gatol 1e-5 922 requires: !single 923 924 test: 925 suffix: 12 926 args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bntl -tao_gatol 1e-5 927 requires: !single 928 929 test: 930 suffix: 13 931 args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bnls -tao_gatol 1e-5 -tao_bnk_max_cg_its 3 932 requires: !single 933 934 test: 935 suffix: 14 936 args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bntr -tao_gatol 1e-5 -tao_bnk_max_cg_its 3 937 requires: !single 938 939 test: 940 suffix: 15 941 args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bntl -tao_gatol 1e-5 -tao_bnk_max_cg_its 3 942 requires: !single 943 944 test: 945 suffix: 16 946 args: -tao_smonitor -mx 8 -my 8 -bmx 2 -bmy 5 -bheight 0.3 -tao_gatol 1e-4 -tao_type bqnls 947 requires: !single 948 949 test: 950 suffix: 17 951 args: -tao_smonitor -mx 8 -my 8 -bmx 2 -bmy 5 -bheight 0.3 -tao_gatol 1e-4 -tao_type bqnkls -tao_bqnk_mat_type lmvmbfgs 952 requires: !single 953 954 test: 955 suffix: 18 956 args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bnls -tao_gatol 1e-5 -tao_mf_hessian 957 requires: !single 958 959 test: 960 suffix: 19 961 args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bntr -tao_gatol 1e-5 -tao_mf_hessian 962 requires: !single 963 964 test: 965 suffix: 20 966 args: -tao_smonitor -mx 8 -my 6 -bmx 3 -bmy 3 -bheight 0.2 -tao_type bntl -tao_gatol 1e-5 -tao_mf_hessian 967 968 TEST*/ 969