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