1 2 static char help[] = "Solves the nonlinear system, the Bratu (SFI - solid fuel ignition) problem in a 2D rectangular domain.\n\ 3 This example also illustrates the use of matrix coloring. Runtime options include:\n\ 4 -par <parameter>, where <parameter> indicates the problem's nonlinearity\n\ 5 problem SFI: <parameter> = Bratu parameter (0 <= par <= 6.81)\n\ 6 -mx <xg>, where <xg> = number of grid points in the x-direction\n\ 7 -my <yg>, where <yg> = number of grid points in the y-direction\n\n"; 8 9 /* 10 11 Solid Fuel Ignition (SFI) problem. This problem is modeled by 12 the partial differential equation 13 14 -Laplacian u - lambda*exp(u) = 0, 0 < x,y < 1, 15 16 with boundary conditions 17 18 u = 0 for x = 0, x = 1, y = 0, y = 1. 19 20 A finite difference approximation with the usual 5-point stencil 21 is used to discretize the boundary value problem to obtain a nonlinear 22 system of equations. 23 24 The parallel version of this code is snes/tutorials/ex5.c 25 26 */ 27 28 /* 29 Include "petscsnes.h" so that we can use SNES solvers. Note that 30 this file automatically includes: 31 petscsys.h - base PETSc routines petscvec.h - vectors 32 petscmat.h - matrices 33 petscis.h - index sets petscksp.h - Krylov subspace methods 34 petscviewer.h - viewers petscpc.h - preconditioners 35 petscksp.h - linear solvers 36 */ 37 38 #include <petscsnes.h> 39 40 /* 41 User-defined application context - contains data needed by the 42 application-provided call-back routines, FormJacobian() and 43 FormFunction(). 44 */ 45 typedef struct { 46 PetscReal param; /* test problem parameter */ 47 PetscInt mx; /* Discretization in x-direction */ 48 PetscInt my; /* Discretization in y-direction */ 49 } AppCtx; 50 51 /* 52 User-defined routines 53 */ 54 extern PetscErrorCode FormJacobian(SNES, Vec, Mat, Mat, void *); 55 extern PetscErrorCode FormFunction(SNES, Vec, Vec, void *); 56 extern PetscErrorCode FormInitialGuess(AppCtx *, Vec); 57 extern PetscErrorCode ConvergenceTest(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *); 58 extern PetscErrorCode ConvergenceDestroy(void *); 59 extern PetscErrorCode postcheck(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *, void *); 60 61 int main(int argc, char **argv) 62 { 63 SNES snes; /* nonlinear solver context */ 64 Vec x, r; /* solution, residual vectors */ 65 Mat J; /* Jacobian matrix */ 66 AppCtx user; /* user-defined application context */ 67 PetscInt i, its, N, hist_its[50]; 68 PetscMPIInt size; 69 PetscReal bratu_lambda_max = 6.81, bratu_lambda_min = 0., history[50]; 70 MatFDColoring fdcoloring; 71 PetscBool matrix_free = PETSC_FALSE, flg, fd_coloring = PETSC_FALSE, use_convergence_test = PETSC_FALSE, pc = PETSC_FALSE, prunejacobian = PETSC_FALSE; 72 KSP ksp; 73 PetscInt *testarray; 74 75 PetscFunctionBeginUser; 76 PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 77 PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 78 PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!"); 79 80 /* 81 Initialize problem parameters 82 */ 83 user.mx = 4; 84 user.my = 4; 85 user.param = 6.0; 86 PetscCall(PetscOptionsGetInt(NULL, NULL, "-mx", &user.mx, NULL)); 87 PetscCall(PetscOptionsGetInt(NULL, NULL, "-my", &user.my, NULL)); 88 PetscCall(PetscOptionsGetReal(NULL, NULL, "-par", &user.param, NULL)); 89 PetscCall(PetscOptionsGetBool(NULL, NULL, "-pc", &pc, NULL)); 90 PetscCheck(user.param < bratu_lambda_max && user.param > bratu_lambda_min, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Lambda is out of range"); 91 N = user.mx * user.my; 92 PetscCall(PetscOptionsGetBool(NULL, NULL, "-use_convergence_test", &use_convergence_test, NULL)); 93 PetscCall(PetscOptionsGetBool(NULL, NULL, "-prune_jacobian", &prunejacobian, NULL)); 94 95 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 96 Create nonlinear solver context 97 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 98 99 PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 100 101 if (pc) { 102 PetscCall(SNESSetType(snes, SNESNEWTONTR)); 103 PetscCall(SNESNewtonTRSetPostCheck(snes, postcheck, NULL)); 104 } 105 106 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 107 Create vector data structures; set function evaluation routine 108 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 109 110 PetscCall(VecCreate(PETSC_COMM_WORLD, &x)); 111 PetscCall(VecSetSizes(x, PETSC_DECIDE, N)); 112 PetscCall(VecSetFromOptions(x)); 113 PetscCall(VecDuplicate(x, &r)); 114 115 /* 116 Set function evaluation routine and vector. Whenever the nonlinear 117 solver needs to evaluate the nonlinear function, it will call this 118 routine. 119 - Note that the final routine argument is the user-defined 120 context that provides application-specific data for the 121 function evaluation routine. 122 */ 123 PetscCall(SNESSetFunction(snes, r, FormFunction, (void *)&user)); 124 125 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 126 Create matrix data structure; set Jacobian evaluation routine 127 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 128 129 /* 130 Create matrix. Here we only approximately preallocate storage space 131 for the Jacobian. See the users manual for a discussion of better 132 techniques for preallocating matrix memory. 133 */ 134 PetscCall(PetscOptionsGetBool(NULL, NULL, "-snes_mf", &matrix_free, NULL)); 135 if (!matrix_free) { 136 PetscBool matrix_free_operator = PETSC_FALSE; 137 PetscCall(PetscOptionsGetBool(NULL, NULL, "-snes_mf_operator", &matrix_free_operator, NULL)); 138 if (matrix_free_operator) matrix_free = PETSC_FALSE; 139 } 140 if (!matrix_free) PetscCall(MatCreateSeqAIJ(PETSC_COMM_WORLD, N, N, 5, NULL, &J)); 141 142 /* 143 This option will cause the Jacobian to be computed via finite differences 144 efficiently using a coloring of the columns of the matrix. 145 */ 146 PetscCall(PetscOptionsGetBool(NULL, NULL, "-snes_fd_coloring", &fd_coloring, NULL)); 147 PetscCheck(!matrix_free || !fd_coloring, PETSC_COMM_WORLD, PETSC_ERR_ARG_INCOMP, "Use only one of -snes_mf, -snes_fd_coloring options!\nYou can do -snes_mf_operator -snes_fd_coloring"); 148 149 if (fd_coloring) { 150 ISColoring iscoloring; 151 MatColoring mc; 152 if (prunejacobian) { 153 /* Initialize x with random nonzero values so that the nonzeros in the Jacobian 154 can better reflect the sparsity structure of the Jacobian. */ 155 PetscRandom rctx; 156 PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rctx)); 157 PetscCall(PetscRandomSetInterval(rctx, 1.0, 2.0)); 158 PetscCall(VecSetRandom(x, rctx)); 159 PetscCall(PetscRandomDestroy(&rctx)); 160 } 161 162 /* 163 This initializes the nonzero structure of the Jacobian. This is artificial 164 because clearly if we had a routine to compute the Jacobian we won't need 165 to use finite differences. 166 */ 167 PetscCall(FormJacobian(snes, x, J, J, &user)); 168 169 /* 170 Color the matrix, i.e. determine groups of columns that share no common 171 rows. These columns in the Jacobian can all be computed simultaneously. 172 */ 173 PetscCall(MatColoringCreate(J, &mc)); 174 PetscCall(MatColoringSetType(mc, MATCOLORINGSL)); 175 PetscCall(MatColoringSetFromOptions(mc)); 176 PetscCall(MatColoringApply(mc, &iscoloring)); 177 PetscCall(MatColoringDestroy(&mc)); 178 /* 179 Create the data structure that SNESComputeJacobianDefaultColor() uses 180 to compute the actual Jacobians via finite differences. 181 */ 182 PetscCall(MatFDColoringCreate(J, iscoloring, &fdcoloring)); 183 PetscCall(MatFDColoringSetFunction(fdcoloring, (PetscErrorCode(*)(void))FormFunction, &user)); 184 PetscCall(MatFDColoringSetFromOptions(fdcoloring)); 185 PetscCall(MatFDColoringSetUp(J, iscoloring, fdcoloring)); 186 /* 187 Tell SNES to use the routine SNESComputeJacobianDefaultColor() 188 to compute Jacobians. 189 */ 190 PetscCall(SNESSetJacobian(snes, J, J, SNESComputeJacobianDefaultColor, fdcoloring)); 191 PetscCall(ISColoringDestroy(&iscoloring)); 192 if (prunejacobian) PetscCall(SNESPruneJacobianColor(snes, J, J)); 193 } 194 /* 195 Set Jacobian matrix data structure and default Jacobian evaluation 196 routine. Whenever the nonlinear solver needs to compute the 197 Jacobian matrix, it will call this routine. 198 - Note that the final routine argument is the user-defined 199 context that provides application-specific data for the 200 Jacobian evaluation routine. 201 - The user can override with: 202 -snes_fd : default finite differencing approximation of Jacobian 203 -snes_mf : matrix-free Newton-Krylov method with no preconditioning 204 (unless user explicitly sets preconditioner) 205 -snes_mf_operator : form preconditioning matrix as set by the user, 206 but use matrix-free approx for Jacobian-vector 207 products within Newton-Krylov method 208 */ 209 else if (!matrix_free) { 210 PetscCall(SNESSetJacobian(snes, J, J, FormJacobian, (void *)&user)); 211 } 212 213 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 214 Customize nonlinear solver; set runtime options 215 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 216 217 /* 218 Set runtime options (e.g., -snes_monitor -snes_rtol <rtol> -ksp_type <type>) 219 */ 220 PetscCall(SNESSetFromOptions(snes)); 221 222 /* 223 Set array that saves the function norms. This array is intended 224 when the user wants to save the convergence history for later use 225 rather than just to view the function norms via -snes_monitor. 226 */ 227 PetscCall(SNESSetConvergenceHistory(snes, history, hist_its, 50, PETSC_TRUE)); 228 229 /* 230 Add a user provided convergence test; this is to test that SNESNEWTONTR properly calls the 231 user provided test before the specialized test. The convergence context is just an array to 232 test that it gets properly freed at the end 233 */ 234 if (use_convergence_test) { 235 PetscCall(SNESGetKSP(snes, &ksp)); 236 PetscCall(PetscMalloc1(5, &testarray)); 237 PetscCall(KSPSetConvergenceTest(ksp, ConvergenceTest, testarray, ConvergenceDestroy)); 238 } 239 240 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 241 Evaluate initial guess; then solve nonlinear system 242 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 243 /* 244 Note: The user should initialize the vector, x, with the initial guess 245 for the nonlinear solver prior to calling SNESSolve(). In particular, 246 to employ an initial guess of zero, the user should explicitly set 247 this vector to zero by calling VecSet(). 248 */ 249 PetscCall(FormInitialGuess(&user, x)); 250 PetscCall(SNESSolve(snes, NULL, x)); 251 PetscCall(SNESGetIterationNumber(snes, &its)); 252 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Number of SNES iterations = %" PetscInt_FMT "\n", its)); 253 254 /* 255 Print the convergence history. This is intended just to demonstrate 256 use of the data attained via SNESSetConvergenceHistory(). 257 */ 258 PetscCall(PetscOptionsHasName(NULL, NULL, "-print_history", &flg)); 259 if (flg) { 260 for (i = 0; i < its + 1; i++) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "iteration %" PetscInt_FMT ": Linear iterations %" PetscInt_FMT " Function norm = %g\n", i, hist_its[i], (double)history[i])); 261 } 262 263 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 264 Free work space. All PETSc objects should be destroyed when they 265 are no longer needed. 266 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 267 268 if (!matrix_free) PetscCall(MatDestroy(&J)); 269 if (fd_coloring) PetscCall(MatFDColoringDestroy(&fdcoloring)); 270 PetscCall(VecDestroy(&x)); 271 PetscCall(VecDestroy(&r)); 272 PetscCall(SNESDestroy(&snes)); 273 PetscCall(PetscFinalize()); 274 return 0; 275 } 276 277 /* 278 FormInitialGuess - Forms initial approximation. 279 280 Input Parameters: 281 user - user-defined application context 282 X - vector 283 284 Output Parameter: 285 X - vector 286 */ 287 PetscErrorCode FormInitialGuess(AppCtx *user, Vec X) 288 { 289 PetscInt i, j, row, mx, my; 290 PetscReal lambda, temp1, temp, hx, hy; 291 PetscScalar *x; 292 293 mx = user->mx; 294 my = user->my; 295 lambda = user->param; 296 297 hx = 1.0 / (PetscReal)(mx - 1); 298 hy = 1.0 / (PetscReal)(my - 1); 299 300 /* 301 Get a pointer to vector data. 302 - For default PETSc vectors, VecGetArray() returns a pointer to 303 the data array. Otherwise, the routine is implementation dependent. 304 - You MUST call VecRestoreArray() when you no longer need access to 305 the array. 306 */ 307 PetscCall(VecGetArray(X, &x)); 308 temp1 = lambda / (lambda + 1.0); 309 for (j = 0; j < my; j++) { 310 temp = (PetscReal)(PetscMin(j, my - j - 1)) * hy; 311 for (i = 0; i < mx; i++) { 312 row = i + j * mx; 313 if (i == 0 || j == 0 || i == mx - 1 || j == my - 1) { 314 x[row] = 0.0; 315 continue; 316 } 317 x[row] = temp1 * PetscSqrtReal(PetscMin((PetscReal)(PetscMin(i, mx - i - 1)) * hx, temp)); 318 } 319 } 320 321 /* 322 Restore vector 323 */ 324 PetscCall(VecRestoreArray(X, &x)); 325 return 0; 326 } 327 328 /* 329 FormFunction - Evaluates nonlinear function, F(x). 330 331 Input Parameters: 332 . snes - the SNES context 333 . X - input vector 334 . ptr - optional user-defined context, as set by SNESSetFunction() 335 336 Output Parameter: 337 . F - function vector 338 */ 339 PetscErrorCode FormFunction(SNES snes, Vec X, Vec F, void *ptr) 340 { 341 AppCtx *user = (AppCtx *)ptr; 342 PetscInt i, j, row, mx, my; 343 PetscReal two = 2.0, one = 1.0, lambda, hx, hy, hxdhy, hydhx; 344 PetscScalar ut, ub, ul, ur, u, uxx, uyy, sc, *f; 345 const PetscScalar *x; 346 347 mx = user->mx; 348 my = user->my; 349 lambda = user->param; 350 hx = one / (PetscReal)(mx - 1); 351 hy = one / (PetscReal)(my - 1); 352 sc = hx * hy; 353 hxdhy = hx / hy; 354 hydhx = hy / hx; 355 356 /* 357 Get pointers to vector data 358 */ 359 PetscCall(VecGetArrayRead(X, &x)); 360 PetscCall(VecGetArray(F, &f)); 361 362 /* 363 Compute function 364 */ 365 for (j = 0; j < my; j++) { 366 for (i = 0; i < mx; i++) { 367 row = i + j * mx; 368 if (i == 0 || j == 0 || i == mx - 1 || j == my - 1) { 369 f[row] = x[row]; 370 continue; 371 } 372 u = x[row]; 373 ub = x[row - mx]; 374 ul = x[row - 1]; 375 ut = x[row + mx]; 376 ur = x[row + 1]; 377 uxx = (-ur + two * u - ul) * hydhx; 378 uyy = (-ut + two * u - ub) * hxdhy; 379 f[row] = uxx + uyy - sc * lambda * PetscExpScalar(u); 380 } 381 } 382 383 /* 384 Restore vectors 385 */ 386 PetscCall(VecRestoreArrayRead(X, &x)); 387 PetscCall(VecRestoreArray(F, &f)); 388 return 0; 389 } 390 391 /* 392 FormJacobian - Evaluates Jacobian matrix. 393 394 Input Parameters: 395 . snes - the SNES context 396 . x - input vector 397 . ptr - optional user-defined context, as set by SNESSetJacobian() 398 399 Output Parameters: 400 . A - Jacobian matrix 401 . B - optionally different preconditioning matrix 402 . flag - flag indicating matrix structure 403 */ 404 PetscErrorCode FormJacobian(SNES snes, Vec X, Mat J, Mat jac, void *ptr) 405 { 406 AppCtx *user = (AppCtx *)ptr; /* user-defined applicatin context */ 407 PetscInt i, j, row, mx, my, col[5]; 408 PetscScalar two = 2.0, one = 1.0, lambda, v[5], sc; 409 const PetscScalar *x; 410 PetscReal hx, hy, hxdhy, hydhx; 411 412 mx = user->mx; 413 my = user->my; 414 lambda = user->param; 415 hx = 1.0 / (PetscReal)(mx - 1); 416 hy = 1.0 / (PetscReal)(my - 1); 417 sc = hx * hy; 418 hxdhy = hx / hy; 419 hydhx = hy / hx; 420 421 /* 422 Get pointer to vector data 423 */ 424 PetscCall(VecGetArrayRead(X, &x)); 425 426 /* 427 Compute entries of the Jacobian 428 */ 429 for (j = 0; j < my; j++) { 430 for (i = 0; i < mx; i++) { 431 row = i + j * mx; 432 if (i == 0 || j == 0 || i == mx - 1 || j == my - 1) { 433 PetscCall(MatSetValues(jac, 1, &row, 1, &row, &one, INSERT_VALUES)); 434 continue; 435 } 436 v[0] = -hxdhy; 437 col[0] = row - mx; 438 v[1] = -hydhx; 439 col[1] = row - 1; 440 v[2] = two * (hydhx + hxdhy) - sc * lambda * PetscExpScalar(x[row]); 441 col[2] = row; 442 v[3] = -hydhx; 443 col[3] = row + 1; 444 v[4] = -hxdhy; 445 col[4] = row + mx; 446 PetscCall(MatSetValues(jac, 1, &row, 5, col, v, INSERT_VALUES)); 447 } 448 } 449 450 /* 451 Restore vector 452 */ 453 PetscCall(VecRestoreArrayRead(X, &x)); 454 455 /* 456 Assemble matrix 457 */ 458 PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY)); 459 PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY)); 460 461 if (jac != J) { 462 PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 463 PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 464 } 465 466 return 0; 467 } 468 469 PetscErrorCode ConvergenceTest(KSP ksp, PetscInt it, PetscReal nrm, KSPConvergedReason *reason, void *ctx) 470 { 471 PetscFunctionBegin; 472 *reason = KSP_CONVERGED_ITERATING; 473 if (it > 1) { 474 *reason = KSP_CONVERGED_ITS; 475 PetscCall(PetscInfo(NULL, "User provided convergence test returning after 2 iterations\n")); 476 } 477 PetscFunctionReturn(0); 478 } 479 480 PetscErrorCode ConvergenceDestroy(void *ctx) 481 { 482 PetscFunctionBegin; 483 PetscCall(PetscInfo(NULL, "User provided convergence destroy called\n")); 484 PetscCall(PetscFree(ctx)); 485 PetscFunctionReturn(0); 486 } 487 488 PetscErrorCode postcheck(SNES snes, Vec x, Vec y, Vec w, PetscBool *changed_y, PetscBool *changed_w, void *ctx) 489 { 490 PetscReal norm; 491 Vec tmp; 492 493 PetscFunctionBegin; 494 PetscCall(VecDuplicate(x, &tmp)); 495 PetscCall(VecWAXPY(tmp, -1.0, x, w)); 496 PetscCall(VecNorm(tmp, NORM_2, &norm)); 497 PetscCall(VecDestroy(&tmp)); 498 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Norm of search step %g\n", (double)norm)); 499 PetscFunctionReturn(0); 500 } 501 502 /*TEST 503 504 build: 505 requires: !single 506 507 test: 508 args: -ksp_gmres_cgs_refinement_type refine_always 509 510 test: 511 suffix: 2 512 args: -snes_monitor_short -snes_type newtontr -ksp_gmres_cgs_refinement_type refine_always 513 514 test: 515 suffix: 2a 516 filter: grep -i KSPConvergedDefault > /dev/null && echo "Found KSPConvergedDefault" 517 args: -snes_monitor_short -snes_type newtontr -ksp_gmres_cgs_refinement_type refine_always -info 518 requires: defined(PETSC_USE_INFO) 519 520 test: 521 suffix: 2b 522 filter: grep -i "User provided convergence test" > /dev/null && echo "Found User provided convergence test" 523 args: -snes_monitor_short -snes_type newtontr -ksp_gmres_cgs_refinement_type refine_always -use_convergence_test -info 524 requires: defined(PETSC_USE_INFO) 525 526 test: 527 suffix: 3 528 args: -snes_monitor_short -mat_coloring_type sl -snes_fd_coloring -mx 8 -my 11 -ksp_gmres_cgs_refinement_type refine_always 529 530 test: 531 suffix: 4 532 args: -pc -par 6.807 -snes_monitor -snes_converged_reason 533 534 test: 535 suffix: 5 536 args: -snes_monitor_short -mat_coloring_type sl -snes_fd_coloring -mx 8 -my 11 -ksp_gmres_cgs_refinement_type refine_always -prune_jacobian 537 output_file: output/ex1_3.out 538 TEST*/ 539