1 /* Program usage: mpiexec -n 1 eptorsion1 [-help] [all TAO options] */ 2 3 /* ---------------------------------------------------------------------- 4 5 Elastic-plastic torsion problem. 6 7 The elastic plastic torsion problem arises from the determination 8 of the stress field on an infinitely long cylindrical bar, which is 9 equivalent to the solution of the following problem: 10 11 min{ .5 * integral(||gradient(v(x))||^2 dx) - C * integral(v(x) dx)} 12 13 where C is the torsion angle per unit length. 14 15 The multiprocessor version of this code is eptorsion2.c. 16 17 ---------------------------------------------------------------------- */ 18 19 /* 20 Include "petsctao.h" so that we can use TAO solvers. Note that this 21 file automatically includes files for lower-level support, such as those 22 provided by the PETSc library: 23 petsc.h - base PETSc routines petscvec.h - vectors 24 petscsys.h - system routines petscmat.h - matrices 25 petscis.h - index sets petscksp.h - Krylov subspace methods 26 petscviewer.h - viewers petscpc.h - preconditioners 27 */ 28 29 #include <petsctao.h> 30 31 static char help[] = "Demonstrates use of the TAO package to solve \n\ 32 unconstrained minimization problems on a single processor. This example \n\ 33 is based on the Elastic-Plastic Torsion (dept) problem from the MINPACK-2 \n\ 34 test suite.\n\ 35 The command line options are:\n\ 36 -mx <xg>, where <xg> = number of grid points in the 1st coordinate direction\n\ 37 -my <yg>, where <yg> = number of grid points in the 2nd coordinate direction\n\ 38 -par <param>, where <param> = angle of twist per unit length\n\n"; 39 40 /* 41 User-defined application context - contains data needed by the 42 application-provided call-back routines, FormFunction(), 43 FormGradient(), and FormHessian(). 44 */ 45 46 typedef struct { 47 PetscReal param; /* nonlinearity parameter */ 48 PetscInt mx, my; /* discretization in x- and y-directions */ 49 PetscInt ndim; /* problem dimension */ 50 Vec s, y, xvec; /* work space for computing Hessian */ 51 PetscReal hx, hy; /* mesh spacing in x- and y-directions */ 52 } AppCtx; 53 54 /* -------- User-defined Routines --------- */ 55 56 PetscErrorCode FormInitialGuess(AppCtx *, Vec); 57 PetscErrorCode FormFunction(Tao, Vec, PetscReal *, void *); 58 PetscErrorCode FormGradient(Tao, Vec, Vec, void *); 59 PetscErrorCode FormHessian(Tao, Vec, Mat, Mat, void *); 60 PetscErrorCode HessianProductMat(Mat, Vec, Vec); 61 PetscErrorCode HessianProduct(void *, Vec, Vec); 62 PetscErrorCode MatrixFreeHessian(Tao, Vec, Mat, Mat, void *); 63 PetscErrorCode FormFunctionGradient(Tao, Vec, PetscReal *, Vec, void *); 64 65 PetscErrorCode main(int argc, char **argv) { 66 PetscInt mx = 10; /* discretization in x-direction */ 67 PetscInt my = 10; /* discretization in y-direction */ 68 Vec x; /* solution, gradient vectors */ 69 PetscBool flg; /* A return value when checking for use options */ 70 Tao tao; /* Tao solver context */ 71 Mat H; /* Hessian matrix */ 72 AppCtx user; /* application context */ 73 PetscMPIInt size; /* number of processes */ 74 PetscReal one = 1.0; 75 76 PetscBool test_lmvm = PETSC_FALSE; 77 KSP ksp; 78 PC pc; 79 Mat M; 80 Vec in, out, out2; 81 PetscReal mult_solve_dist; 82 83 /* Initialize TAO,PETSc */ 84 PetscFunctionBeginUser; 85 PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 86 PetscCallMPI(MPI_Comm_size(MPI_COMM_WORLD, &size)); 87 PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Incorrect number of processors"); 88 89 /* Specify default parameters for the problem, check for command-line overrides */ 90 user.param = 5.0; 91 PetscCall(PetscOptionsGetInt(NULL, NULL, "-my", &my, &flg)); 92 PetscCall(PetscOptionsGetInt(NULL, NULL, "-mx", &mx, &flg)); 93 PetscCall(PetscOptionsGetReal(NULL, NULL, "-par", &user.param, &flg)); 94 PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_lmvm", &test_lmvm, &flg)); 95 96 PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n---- Elastic-Plastic Torsion Problem -----\n")); 97 PetscCall(PetscPrintf(PETSC_COMM_SELF, "mx: %" PetscInt_FMT " my: %" PetscInt_FMT " \n\n", mx, my)); 98 user.ndim = mx * my; 99 user.mx = mx; 100 user.my = my; 101 user.hx = one / (mx + 1); 102 user.hy = one / (my + 1); 103 104 /* Allocate vectors */ 105 PetscCall(VecCreateSeq(PETSC_COMM_SELF, user.ndim, &user.y)); 106 PetscCall(VecDuplicate(user.y, &user.s)); 107 PetscCall(VecDuplicate(user.y, &x)); 108 109 /* Create TAO solver and set desired solution method */ 110 PetscCall(TaoCreate(PETSC_COMM_SELF, &tao)); 111 PetscCall(TaoSetType(tao, TAOLMVM)); 112 113 /* Set solution vector with an initial guess */ 114 PetscCall(FormInitialGuess(&user, x)); 115 PetscCall(TaoSetSolution(tao, x)); 116 117 /* Set routine for function and gradient evaluation */ 118 PetscCall(TaoSetObjectiveAndGradient(tao, NULL, FormFunctionGradient, (void *)&user)); 119 120 /* From command line options, determine if using matrix-free hessian */ 121 PetscCall(PetscOptionsHasName(NULL, NULL, "-my_tao_mf", &flg)); 122 if (flg) { 123 PetscCall(MatCreateShell(PETSC_COMM_SELF, user.ndim, user.ndim, user.ndim, user.ndim, (void *)&user, &H)); 124 PetscCall(MatShellSetOperation(H, MATOP_MULT, (void (*)(void))HessianProductMat)); 125 PetscCall(MatSetOption(H, MAT_SYMMETRIC, PETSC_TRUE)); 126 127 PetscCall(TaoSetHessian(tao, H, H, MatrixFreeHessian, (void *)&user)); 128 } else { 129 PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, user.ndim, user.ndim, 5, NULL, &H)); 130 PetscCall(MatSetOption(H, MAT_SYMMETRIC, PETSC_TRUE)); 131 PetscCall(TaoSetHessian(tao, H, H, FormHessian, (void *)&user)); 132 } 133 134 /* Test the LMVM matrix */ 135 if (test_lmvm) { 136 PetscCall(PetscOptionsSetValue(NULL, "-tao_type", "bntr")); 137 PetscCall(PetscOptionsSetValue(NULL, "-tao_bnk_pc_type", "lmvm")); 138 } 139 140 /* Check for any TAO command line options */ 141 PetscCall(TaoSetFromOptions(tao)); 142 143 /* SOLVE THE APPLICATION */ 144 PetscCall(TaoSolve(tao)); 145 146 /* Test the LMVM matrix */ 147 if (test_lmvm) { 148 PetscCall(TaoGetKSP(tao, &ksp)); 149 PetscCall(KSPGetPC(ksp, &pc)); 150 PetscCall(PCLMVMGetMatLMVM(pc, &M)); 151 PetscCall(VecDuplicate(x, &in)); 152 PetscCall(VecDuplicate(x, &out)); 153 PetscCall(VecDuplicate(x, &out2)); 154 PetscCall(VecSet(in, 5.0)); 155 PetscCall(MatMult(M, in, out)); 156 PetscCall(MatSolve(M, out, out2)); 157 PetscCall(VecAXPY(out2, -1.0, in)); 158 PetscCall(VecNorm(out2, NORM_2, &mult_solve_dist)); 159 PetscCall(PetscPrintf(PetscObjectComm((PetscObject)tao), "error between MatMult and MatSolve: %e\n", (double)mult_solve_dist)); 160 PetscCall(VecDestroy(&in)); 161 PetscCall(VecDestroy(&out)); 162 PetscCall(VecDestroy(&out2)); 163 } 164 165 PetscCall(TaoDestroy(&tao)); 166 PetscCall(VecDestroy(&user.s)); 167 PetscCall(VecDestroy(&user.y)); 168 PetscCall(VecDestroy(&x)); 169 PetscCall(MatDestroy(&H)); 170 171 PetscCall(PetscFinalize()); 172 return 0; 173 } 174 175 /* ------------------------------------------------------------------- */ 176 /* 177 FormInitialGuess - Computes an initial approximation to the solution. 178 179 Input Parameters: 180 . user - user-defined application context 181 . X - vector 182 183 Output Parameters: 184 . X - vector 185 */ 186 PetscErrorCode FormInitialGuess(AppCtx *user, Vec X) { 187 PetscReal hx = user->hx, hy = user->hy, temp; 188 PetscReal val; 189 PetscInt i, j, k, nx = user->mx, ny = user->my; 190 191 /* Compute initial guess */ 192 PetscFunctionBeginUser; 193 for (j = 0; j < ny; j++) { 194 temp = PetscMin(j + 1, ny - j) * hy; 195 for (i = 0; i < nx; i++) { 196 k = nx * j + i; 197 val = PetscMin((PetscMin(i + 1, nx - i)) * hx, temp); 198 PetscCall(VecSetValues(X, 1, &k, &val, ADD_VALUES)); 199 } 200 } 201 PetscCall(VecAssemblyBegin(X)); 202 PetscCall(VecAssemblyEnd(X)); 203 PetscFunctionReturn(0); 204 } 205 206 /* ------------------------------------------------------------------- */ 207 /* 208 FormFunctionGradient - Evaluates the function and corresponding gradient. 209 210 Input Parameters: 211 tao - the Tao context 212 X - the input vector 213 ptr - optional user-defined context, as set by TaoSetFunction() 214 215 Output Parameters: 216 f - the newly evaluated function 217 G - the newly evaluated gradient 218 */ 219 PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *f, Vec G, void *ptr) { 220 PetscFunctionBeginUser; 221 PetscCall(FormFunction(tao, X, f, ptr)); 222 PetscCall(FormGradient(tao, X, G, ptr)); 223 PetscFunctionReturn(0); 224 } 225 226 /* ------------------------------------------------------------------- */ 227 /* 228 FormFunction - Evaluates the function, f(X). 229 230 Input Parameters: 231 . tao - the Tao context 232 . X - the input vector 233 . ptr - optional user-defined context, as set by TaoSetFunction() 234 235 Output Parameters: 236 . f - the newly evaluated function 237 */ 238 PetscErrorCode FormFunction(Tao tao, Vec X, PetscReal *f, void *ptr) { 239 AppCtx *user = (AppCtx *)ptr; 240 PetscReal hx = user->hx, hy = user->hy, area, three = 3.0, p5 = 0.5; 241 PetscReal zero = 0.0, vb, vl, vr, vt, dvdx, dvdy, flin = 0.0, fquad = 0.0; 242 PetscReal v, cdiv3 = user->param / three; 243 const PetscScalar *x; 244 PetscInt nx = user->mx, ny = user->my, i, j, k; 245 246 PetscFunctionBeginUser; 247 /* Get pointer to vector data */ 248 PetscCall(VecGetArrayRead(X, &x)); 249 250 /* Compute function contributions over the lower triangular elements */ 251 for (j = -1; j < ny; j++) { 252 for (i = -1; i < nx; i++) { 253 k = nx * j + i; 254 v = zero; 255 vr = zero; 256 vt = zero; 257 if (i >= 0 && j >= 0) v = x[k]; 258 if (i < nx - 1 && j > -1) vr = x[k + 1]; 259 if (i > -1 && j < ny - 1) vt = x[k + nx]; 260 dvdx = (vr - v) / hx; 261 dvdy = (vt - v) / hy; 262 fquad += dvdx * dvdx + dvdy * dvdy; 263 flin -= cdiv3 * (v + vr + vt); 264 } 265 } 266 267 /* Compute function contributions over the upper triangular elements */ 268 for (j = 0; j <= ny; j++) { 269 for (i = 0; i <= nx; i++) { 270 k = nx * j + i; 271 vb = zero; 272 vl = zero; 273 v = zero; 274 if (i < nx && j > 0) vb = x[k - nx]; 275 if (i > 0 && j < ny) vl = x[k - 1]; 276 if (i < nx && j < ny) v = x[k]; 277 dvdx = (v - vl) / hx; 278 dvdy = (v - vb) / hy; 279 fquad = fquad + dvdx * dvdx + dvdy * dvdy; 280 flin = flin - cdiv3 * (vb + vl + v); 281 } 282 } 283 284 /* Restore vector */ 285 PetscCall(VecRestoreArrayRead(X, &x)); 286 287 /* Scale the function */ 288 area = p5 * hx * hy; 289 *f = area * (p5 * fquad + flin); 290 291 PetscCall(PetscLogFlops(24.0 * nx * ny)); 292 PetscFunctionReturn(0); 293 } 294 295 /* ------------------------------------------------------------------- */ 296 /* 297 FormGradient - Evaluates the gradient, G(X). 298 299 Input Parameters: 300 . tao - the Tao context 301 . X - input vector 302 . ptr - optional user-defined context 303 304 Output Parameters: 305 . G - vector containing the newly evaluated gradient 306 */ 307 PetscErrorCode FormGradient(Tao tao, Vec X, Vec G, void *ptr) { 308 AppCtx *user = (AppCtx *)ptr; 309 PetscReal zero = 0.0, p5 = 0.5, three = 3.0, area, val; 310 PetscInt nx = user->mx, ny = user->my, ind, i, j, k; 311 PetscReal hx = user->hx, hy = user->hy; 312 PetscReal vb, vl, vr, vt, dvdx, dvdy; 313 PetscReal v, cdiv3 = user->param / three; 314 const PetscScalar *x; 315 316 PetscFunctionBeginUser; 317 /* Initialize gradient to zero */ 318 PetscCall(VecSet(G, zero)); 319 320 /* Get pointer to vector data */ 321 PetscCall(VecGetArrayRead(X, &x)); 322 323 /* Compute gradient contributions over the lower triangular elements */ 324 for (j = -1; j < ny; j++) { 325 for (i = -1; i < nx; i++) { 326 k = nx * j + i; 327 v = zero; 328 vr = zero; 329 vt = zero; 330 if (i >= 0 && j >= 0) v = x[k]; 331 if (i < nx - 1 && j > -1) vr = x[k + 1]; 332 if (i > -1 && j < ny - 1) vt = x[k + nx]; 333 dvdx = (vr - v) / hx; 334 dvdy = (vt - v) / hy; 335 if (i != -1 && j != -1) { 336 ind = k; 337 val = -dvdx / hx - dvdy / hy - cdiv3; 338 PetscCall(VecSetValues(G, 1, &ind, &val, ADD_VALUES)); 339 } 340 if (i != nx - 1 && j != -1) { 341 ind = k + 1; 342 val = dvdx / hx - cdiv3; 343 PetscCall(VecSetValues(G, 1, &ind, &val, ADD_VALUES)); 344 } 345 if (i != -1 && j != ny - 1) { 346 ind = k + nx; 347 val = dvdy / hy - cdiv3; 348 PetscCall(VecSetValues(G, 1, &ind, &val, ADD_VALUES)); 349 } 350 } 351 } 352 353 /* Compute gradient contributions over the upper triangular elements */ 354 for (j = 0; j <= ny; j++) { 355 for (i = 0; i <= nx; i++) { 356 k = nx * j + i; 357 vb = zero; 358 vl = zero; 359 v = zero; 360 if (i < nx && j > 0) vb = x[k - nx]; 361 if (i > 0 && j < ny) vl = x[k - 1]; 362 if (i < nx && j < ny) v = x[k]; 363 dvdx = (v - vl) / hx; 364 dvdy = (v - vb) / hy; 365 if (i != nx && j != 0) { 366 ind = k - nx; 367 val = -dvdy / hy - cdiv3; 368 PetscCall(VecSetValues(G, 1, &ind, &val, ADD_VALUES)); 369 } 370 if (i != 0 && j != ny) { 371 ind = k - 1; 372 val = -dvdx / hx - cdiv3; 373 PetscCall(VecSetValues(G, 1, &ind, &val, ADD_VALUES)); 374 } 375 if (i != nx && j != ny) { 376 ind = k; 377 val = dvdx / hx + dvdy / hy - cdiv3; 378 PetscCall(VecSetValues(G, 1, &ind, &val, ADD_VALUES)); 379 } 380 } 381 } 382 PetscCall(VecRestoreArrayRead(X, &x)); 383 384 /* Assemble gradient vector */ 385 PetscCall(VecAssemblyBegin(G)); 386 PetscCall(VecAssemblyEnd(G)); 387 388 /* Scale the gradient */ 389 area = p5 * hx * hy; 390 PetscCall(VecScale(G, area)); 391 PetscCall(PetscLogFlops(24.0 * nx * ny)); 392 PetscFunctionReturn(0); 393 } 394 395 /* ------------------------------------------------------------------- */ 396 /* 397 FormHessian - Forms the Hessian matrix. 398 399 Input Parameters: 400 . tao - the Tao context 401 . X - the input vector 402 . ptr - optional user-defined context, as set by TaoSetHessian() 403 404 Output Parameters: 405 . H - Hessian matrix 406 . PrecH - optionally different preconditioning Hessian 407 . flag - flag indicating matrix structure 408 409 Notes: 410 This routine is intended simply as an example of the interface 411 to a Hessian evaluation routine. Since this example compute the 412 Hessian a column at a time, it is not particularly efficient and 413 is not recommended. 414 */ 415 PetscErrorCode FormHessian(Tao tao, Vec X, Mat H, Mat Hpre, void *ptr) { 416 AppCtx *user = (AppCtx *)ptr; 417 PetscInt i, j, ndim = user->ndim; 418 PetscReal *y, zero = 0.0, one = 1.0; 419 PetscBool assembled; 420 421 PetscFunctionBeginUser; 422 user->xvec = X; 423 424 /* Initialize Hessian entries and work vector to zero */ 425 PetscCall(MatAssembled(H, &assembled)); 426 if (assembled) PetscCall(MatZeroEntries(H)); 427 428 PetscCall(VecSet(user->s, zero)); 429 430 /* Loop over matrix columns to compute entries of the Hessian */ 431 for (j = 0; j < ndim; j++) { 432 PetscCall(VecSetValues(user->s, 1, &j, &one, INSERT_VALUES)); 433 PetscCall(VecAssemblyBegin(user->s)); 434 PetscCall(VecAssemblyEnd(user->s)); 435 436 PetscCall(HessianProduct(ptr, user->s, user->y)); 437 438 PetscCall(VecSetValues(user->s, 1, &j, &zero, INSERT_VALUES)); 439 PetscCall(VecAssemblyBegin(user->s)); 440 PetscCall(VecAssemblyEnd(user->s)); 441 442 PetscCall(VecGetArray(user->y, &y)); 443 for (i = 0; i < ndim; i++) { 444 if (y[i] != zero) { PetscCall(MatSetValues(H, 1, &i, 1, &j, &y[i], ADD_VALUES)); } 445 } 446 PetscCall(VecRestoreArray(user->y, &y)); 447 } 448 PetscCall(MatAssemblyBegin(H, MAT_FINAL_ASSEMBLY)); 449 PetscCall(MatAssemblyEnd(H, MAT_FINAL_ASSEMBLY)); 450 PetscFunctionReturn(0); 451 } 452 453 /* ------------------------------------------------------------------- */ 454 /* 455 MatrixFreeHessian - Sets a pointer for use in computing Hessian-vector 456 products. 457 458 Input Parameters: 459 . tao - the Tao context 460 . X - the input vector 461 . ptr - optional user-defined context, as set by TaoSetHessian() 462 463 Output Parameters: 464 . H - Hessian matrix 465 . PrecH - optionally different preconditioning Hessian 466 . flag - flag indicating matrix structure 467 */ 468 PetscErrorCode MatrixFreeHessian(Tao tao, Vec X, Mat H, Mat PrecH, void *ptr) { 469 AppCtx *user = (AppCtx *)ptr; 470 471 /* Sets location of vector for use in computing matrix-vector products of the form H(X)*y */ 472 PetscFunctionBeginUser; 473 user->xvec = X; 474 PetscFunctionReturn(0); 475 } 476 477 /* ------------------------------------------------------------------- */ 478 /* 479 HessianProductMat - Computes the matrix-vector product 480 y = mat*svec. 481 482 Input Parameters: 483 . mat - input matrix 484 . svec - input vector 485 486 Output Parameters: 487 . y - solution vector 488 */ 489 PetscErrorCode HessianProductMat(Mat mat, Vec svec, Vec y) { 490 void *ptr; 491 492 PetscFunctionBeginUser; 493 PetscCall(MatShellGetContext(mat, &ptr)); 494 PetscCall(HessianProduct(ptr, svec, y)); 495 PetscFunctionReturn(0); 496 } 497 498 /* ------------------------------------------------------------------- */ 499 /* 500 Hessian Product - Computes the matrix-vector product: 501 y = f''(x)*svec. 502 503 Input Parameters: 504 . ptr - pointer to the user-defined context 505 . svec - input vector 506 507 Output Parameters: 508 . y - product vector 509 */ 510 PetscErrorCode HessianProduct(void *ptr, Vec svec, Vec y) { 511 AppCtx *user = (AppCtx *)ptr; 512 PetscReal p5 = 0.5, zero = 0.0, one = 1.0, hx, hy, val, area; 513 const PetscScalar *x, *s; 514 PetscReal v, vb, vl, vr, vt, hxhx, hyhy; 515 PetscInt nx, ny, i, j, k, ind; 516 517 PetscFunctionBeginUser; 518 nx = user->mx; 519 ny = user->my; 520 hx = user->hx; 521 hy = user->hy; 522 hxhx = one / (hx * hx); 523 hyhy = one / (hy * hy); 524 525 /* Get pointers to vector data */ 526 PetscCall(VecGetArrayRead(user->xvec, &x)); 527 PetscCall(VecGetArrayRead(svec, &s)); 528 529 /* Initialize product vector to zero */ 530 PetscCall(VecSet(y, zero)); 531 532 /* Compute f''(x)*s over the lower triangular elements */ 533 for (j = -1; j < ny; j++) { 534 for (i = -1; i < nx; i++) { 535 k = nx * j + i; 536 v = zero; 537 vr = zero; 538 vt = zero; 539 if (i != -1 && j != -1) v = s[k]; 540 if (i != nx - 1 && j != -1) { 541 vr = s[k + 1]; 542 ind = k + 1; 543 val = hxhx * (vr - v); 544 PetscCall(VecSetValues(y, 1, &ind, &val, ADD_VALUES)); 545 } 546 if (i != -1 && j != ny - 1) { 547 vt = s[k + nx]; 548 ind = k + nx; 549 val = hyhy * (vt - v); 550 PetscCall(VecSetValues(y, 1, &ind, &val, ADD_VALUES)); 551 } 552 if (i != -1 && j != -1) { 553 ind = k; 554 val = hxhx * (v - vr) + hyhy * (v - vt); 555 PetscCall(VecSetValues(y, 1, &ind, &val, ADD_VALUES)); 556 } 557 } 558 } 559 560 /* Compute f''(x)*s over the upper triangular elements */ 561 for (j = 0; j <= ny; j++) { 562 for (i = 0; i <= nx; i++) { 563 k = nx * j + i; 564 v = zero; 565 vl = zero; 566 vb = zero; 567 if (i != nx && j != ny) v = s[k]; 568 if (i != nx && j != 0) { 569 vb = s[k - nx]; 570 ind = k - nx; 571 val = hyhy * (vb - v); 572 PetscCall(VecSetValues(y, 1, &ind, &val, ADD_VALUES)); 573 } 574 if (i != 0 && j != ny) { 575 vl = s[k - 1]; 576 ind = k - 1; 577 val = hxhx * (vl - v); 578 PetscCall(VecSetValues(y, 1, &ind, &val, ADD_VALUES)); 579 } 580 if (i != nx && j != ny) { 581 ind = k; 582 val = hxhx * (v - vl) + hyhy * (v - vb); 583 PetscCall(VecSetValues(y, 1, &ind, &val, ADD_VALUES)); 584 } 585 } 586 } 587 /* Restore vector data */ 588 PetscCall(VecRestoreArrayRead(svec, &s)); 589 PetscCall(VecRestoreArrayRead(user->xvec, &x)); 590 591 /* Assemble vector */ 592 PetscCall(VecAssemblyBegin(y)); 593 PetscCall(VecAssemblyEnd(y)); 594 595 /* Scale resulting vector by area */ 596 area = p5 * hx * hy; 597 PetscCall(VecScale(y, area)); 598 PetscCall(PetscLogFlops(18.0 * nx * ny)); 599 PetscFunctionReturn(0); 600 } 601 602 /*TEST 603 604 build: 605 requires: !complex 606 607 test: 608 suffix: 1 609 args: -tao_smonitor -tao_type ntl -tao_gatol 1.e-4 610 611 test: 612 suffix: 2 613 args: -tao_smonitor -tao_type ntr -tao_gatol 1.e-4 614 615 test: 616 suffix: 3 617 args: -tao_smonitor -tao_type bntr -tao_gatol 1.e-4 -my_tao_mf -tao_test_hessian 618 619 test: 620 suffix: 4 621 args: -tao_smonitor -tao_gatol 1e-3 -tao_type bqnls 622 623 test: 624 suffix: 5 625 args: -tao_smonitor -tao_gatol 1e-3 -tao_type blmvm 626 627 test: 628 suffix: 6 629 args: -tao_smonitor -tao_gatol 1e-3 -tao_type bqnktr -tao_bqnk_mat_type lmvmsr1 630 631 TEST*/ 632