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; 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 94 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 95 Create nonlinear solver context 96 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 97 98 PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 99 100 if (pc) { 101 PetscCall(SNESSetType(snes, SNESNEWTONTR)); 102 PetscCall(SNESNewtonTRSetPostCheck(snes, postcheck, NULL)); 103 } 104 105 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 106 Create vector data structures; set function evaluation routine 107 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 108 109 PetscCall(VecCreate(PETSC_COMM_WORLD, &x)); 110 PetscCall(VecSetSizes(x, PETSC_DECIDE, N)); 111 PetscCall(VecSetFromOptions(x)); 112 PetscCall(VecDuplicate(x, &r)); 113 114 /* 115 Set function evaluation routine and vector. Whenever the nonlinear 116 solver needs to evaluate the nonlinear function, it will call this 117 routine. 118 - Note that the final routine argument is the user-defined 119 context that provides application-specific data for the 120 function evaluation routine. 121 */ 122 PetscCall(SNESSetFunction(snes, r, FormFunction, (void *)&user)); 123 124 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 125 Create matrix data structure; set Jacobian evaluation routine 126 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 127 128 /* 129 Create matrix. Here we only approximately preallocate storage space 130 for the Jacobian. See the users manual for a discussion of better 131 techniques for preallocating matrix memory. 132 */ 133 PetscCall(PetscOptionsGetBool(NULL, NULL, "-snes_mf", &matrix_free, NULL)); 134 if (!matrix_free) { 135 PetscBool matrix_free_operator = PETSC_FALSE; 136 PetscCall(PetscOptionsGetBool(NULL, NULL, "-snes_mf_operator", &matrix_free_operator, NULL)); 137 if (matrix_free_operator) matrix_free = PETSC_FALSE; 138 } 139 if (!matrix_free) PetscCall(MatCreateSeqAIJ(PETSC_COMM_WORLD, N, N, 5, NULL, &J)); 140 141 /* 142 This option will cause the Jacobian to be computed via finite differences 143 efficiently using a coloring of the columns of the matrix. 144 */ 145 PetscCall(PetscOptionsGetBool(NULL, NULL, "-snes_fd_coloring", &fd_coloring, NULL)); 146 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"); 147 148 if (fd_coloring) { 149 ISColoring iscoloring; 150 MatColoring mc; 151 152 /* 153 This initializes the nonzero structure of the Jacobian. This is artificial 154 because clearly if we had a routine to compute the Jacobian we won't need 155 to use finite differences. 156 */ 157 PetscCall(FormJacobian(snes, x, J, J, &user)); 158 159 /* 160 Color the matrix, i.e. determine groups of columns that share no common 161 rows. These columns in the Jacobian can all be computed simultaneously. 162 */ 163 PetscCall(MatColoringCreate(J, &mc)); 164 PetscCall(MatColoringSetType(mc, MATCOLORINGSL)); 165 PetscCall(MatColoringSetFromOptions(mc)); 166 PetscCall(MatColoringApply(mc, &iscoloring)); 167 PetscCall(MatColoringDestroy(&mc)); 168 /* 169 Create the data structure that SNESComputeJacobianDefaultColor() uses 170 to compute the actual Jacobians via finite differences. 171 */ 172 PetscCall(MatFDColoringCreate(J, iscoloring, &fdcoloring)); 173 PetscCall(MatFDColoringSetFunction(fdcoloring, (PetscErrorCode(*)(void))FormFunction, &user)); 174 PetscCall(MatFDColoringSetFromOptions(fdcoloring)); 175 PetscCall(MatFDColoringSetUp(J, iscoloring, fdcoloring)); 176 /* 177 Tell SNES to use the routine SNESComputeJacobianDefaultColor() 178 to compute Jacobians. 179 */ 180 PetscCall(SNESSetJacobian(snes, J, J, SNESComputeJacobianDefaultColor, fdcoloring)); 181 PetscCall(ISColoringDestroy(&iscoloring)); 182 } 183 /* 184 Set Jacobian matrix data structure and default Jacobian evaluation 185 routine. Whenever the nonlinear solver needs to compute the 186 Jacobian matrix, it will call this routine. 187 - Note that the final routine argument is the user-defined 188 context that provides application-specific data for the 189 Jacobian evaluation routine. 190 - The user can override with: 191 -snes_fd : default finite differencing approximation of Jacobian 192 -snes_mf : matrix-free Newton-Krylov method with no preconditioning 193 (unless user explicitly sets preconditioner) 194 -snes_mf_operator : form preconditioning matrix as set by the user, 195 but use matrix-free approx for Jacobian-vector 196 products within Newton-Krylov method 197 */ 198 else if (!matrix_free) { 199 PetscCall(SNESSetJacobian(snes, J, J, FormJacobian, (void *)&user)); 200 } 201 202 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 203 Customize nonlinear solver; set runtime options 204 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 205 206 /* 207 Set runtime options (e.g., -snes_monitor -snes_rtol <rtol> -ksp_type <type>) 208 */ 209 PetscCall(SNESSetFromOptions(snes)); 210 211 /* 212 Set array that saves the function norms. This array is intended 213 when the user wants to save the convergence history for later use 214 rather than just to view the function norms via -snes_monitor. 215 */ 216 PetscCall(SNESSetConvergenceHistory(snes, history, hist_its, 50, PETSC_TRUE)); 217 218 /* 219 Add a user provided convergence test; this is to test that SNESNEWTONTR properly calls the 220 user provided test before the specialized test. The convergence context is just an array to 221 test that it gets properly freed at the end 222 */ 223 if (use_convergence_test) { 224 PetscCall(SNESGetKSP(snes, &ksp)); 225 PetscCall(PetscMalloc1(5, &testarray)); 226 PetscCall(KSPSetConvergenceTest(ksp, ConvergenceTest, testarray, ConvergenceDestroy)); 227 } 228 229 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 230 Evaluate initial guess; then solve nonlinear system 231 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 232 /* 233 Note: The user should initialize the vector, x, with the initial guess 234 for the nonlinear solver prior to calling SNESSolve(). In particular, 235 to employ an initial guess of zero, the user should explicitly set 236 this vector to zero by calling VecSet(). 237 */ 238 PetscCall(FormInitialGuess(&user, x)); 239 PetscCall(SNESSolve(snes, NULL, x)); 240 PetscCall(SNESGetIterationNumber(snes, &its)); 241 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Number of SNES iterations = %" PetscInt_FMT "\n", its)); 242 243 /* 244 Print the convergence history. This is intended just to demonstrate 245 use of the data attained via SNESSetConvergenceHistory(). 246 */ 247 PetscCall(PetscOptionsHasName(NULL, NULL, "-print_history", &flg)); 248 if (flg) { 249 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])); 250 } 251 252 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 253 Free work space. All PETSc objects should be destroyed when they 254 are no longer needed. 255 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 256 257 if (!matrix_free) PetscCall(MatDestroy(&J)); 258 if (fd_coloring) PetscCall(MatFDColoringDestroy(&fdcoloring)); 259 PetscCall(VecDestroy(&x)); 260 PetscCall(VecDestroy(&r)); 261 PetscCall(SNESDestroy(&snes)); 262 PetscCall(PetscFinalize()); 263 return 0; 264 } 265 266 /* 267 FormInitialGuess - Forms initial approximation. 268 269 Input Parameters: 270 user - user-defined application context 271 X - vector 272 273 Output Parameter: 274 X - vector 275 */ 276 PetscErrorCode FormInitialGuess(AppCtx *user, Vec X) 277 { 278 PetscInt i, j, row, mx, my; 279 PetscReal lambda, temp1, temp, hx, hy; 280 PetscScalar *x; 281 282 mx = user->mx; 283 my = user->my; 284 lambda = user->param; 285 286 hx = 1.0 / (PetscReal)(mx - 1); 287 hy = 1.0 / (PetscReal)(my - 1); 288 289 /* 290 Get a pointer to vector data. 291 - For default PETSc vectors, VecGetArray() returns a pointer to 292 the data array. Otherwise, the routine is implementation dependent. 293 - You MUST call VecRestoreArray() when you no longer need access to 294 the array. 295 */ 296 PetscCall(VecGetArray(X, &x)); 297 temp1 = lambda / (lambda + 1.0); 298 for (j = 0; j < my; j++) { 299 temp = (PetscReal)(PetscMin(j, my - j - 1)) * hy; 300 for (i = 0; i < mx; i++) { 301 row = i + j * mx; 302 if (i == 0 || j == 0 || i == mx - 1 || j == my - 1) { 303 x[row] = 0.0; 304 continue; 305 } 306 x[row] = temp1 * PetscSqrtReal(PetscMin((PetscReal)(PetscMin(i, mx - i - 1)) * hx, temp)); 307 } 308 } 309 310 /* 311 Restore vector 312 */ 313 PetscCall(VecRestoreArray(X, &x)); 314 return 0; 315 } 316 317 /* 318 FormFunction - Evaluates nonlinear function, F(x). 319 320 Input Parameters: 321 . snes - the SNES context 322 . X - input vector 323 . ptr - optional user-defined context, as set by SNESSetFunction() 324 325 Output Parameter: 326 . F - function vector 327 */ 328 PetscErrorCode FormFunction(SNES snes, Vec X, Vec F, void *ptr) 329 { 330 AppCtx *user = (AppCtx *)ptr; 331 PetscInt i, j, row, mx, my; 332 PetscReal two = 2.0, one = 1.0, lambda, hx, hy, hxdhy, hydhx; 333 PetscScalar ut, ub, ul, ur, u, uxx, uyy, sc, *f; 334 const PetscScalar *x; 335 336 mx = user->mx; 337 my = user->my; 338 lambda = user->param; 339 hx = one / (PetscReal)(mx - 1); 340 hy = one / (PetscReal)(my - 1); 341 sc = hx * hy; 342 hxdhy = hx / hy; 343 hydhx = hy / hx; 344 345 /* 346 Get pointers to vector data 347 */ 348 PetscCall(VecGetArrayRead(X, &x)); 349 PetscCall(VecGetArray(F, &f)); 350 351 /* 352 Compute function 353 */ 354 for (j = 0; j < my; j++) { 355 for (i = 0; i < mx; i++) { 356 row = i + j * mx; 357 if (i == 0 || j == 0 || i == mx - 1 || j == my - 1) { 358 f[row] = x[row]; 359 continue; 360 } 361 u = x[row]; 362 ub = x[row - mx]; 363 ul = x[row - 1]; 364 ut = x[row + mx]; 365 ur = x[row + 1]; 366 uxx = (-ur + two * u - ul) * hydhx; 367 uyy = (-ut + two * u - ub) * hxdhy; 368 f[row] = uxx + uyy - sc * lambda * PetscExpScalar(u); 369 } 370 } 371 372 /* 373 Restore vectors 374 */ 375 PetscCall(VecRestoreArrayRead(X, &x)); 376 PetscCall(VecRestoreArray(F, &f)); 377 return 0; 378 } 379 380 /* 381 FormJacobian - Evaluates Jacobian matrix. 382 383 Input Parameters: 384 . snes - the SNES context 385 . x - input vector 386 . ptr - optional user-defined context, as set by SNESSetJacobian() 387 388 Output Parameters: 389 . A - Jacobian matrix 390 . B - optionally different preconditioning matrix 391 . flag - flag indicating matrix structure 392 */ 393 PetscErrorCode FormJacobian(SNES snes, Vec X, Mat J, Mat jac, void *ptr) 394 { 395 AppCtx *user = (AppCtx *)ptr; /* user-defined applicatin context */ 396 PetscInt i, j, row, mx, my, col[5]; 397 PetscScalar two = 2.0, one = 1.0, lambda, v[5], sc; 398 const PetscScalar *x; 399 PetscReal hx, hy, hxdhy, hydhx; 400 401 mx = user->mx; 402 my = user->my; 403 lambda = user->param; 404 hx = 1.0 / (PetscReal)(mx - 1); 405 hy = 1.0 / (PetscReal)(my - 1); 406 sc = hx * hy; 407 hxdhy = hx / hy; 408 hydhx = hy / hx; 409 410 /* 411 Get pointer to vector data 412 */ 413 PetscCall(VecGetArrayRead(X, &x)); 414 415 /* 416 Compute entries of the Jacobian 417 */ 418 for (j = 0; j < my; j++) { 419 for (i = 0; i < mx; i++) { 420 row = i + j * mx; 421 if (i == 0 || j == 0 || i == mx - 1 || j == my - 1) { 422 PetscCall(MatSetValues(jac, 1, &row, 1, &row, &one, INSERT_VALUES)); 423 continue; 424 } 425 v[0] = -hxdhy; 426 col[0] = row - mx; 427 v[1] = -hydhx; 428 col[1] = row - 1; 429 v[2] = two * (hydhx + hxdhy) - sc * lambda * PetscExpScalar(x[row]); 430 col[2] = row; 431 v[3] = -hydhx; 432 col[3] = row + 1; 433 v[4] = -hxdhy; 434 col[4] = row + mx; 435 PetscCall(MatSetValues(jac, 1, &row, 5, col, v, INSERT_VALUES)); 436 } 437 } 438 439 /* 440 Restore vector 441 */ 442 PetscCall(VecRestoreArrayRead(X, &x)); 443 444 /* 445 Assemble matrix 446 */ 447 PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY)); 448 PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY)); 449 450 if (jac != J) { 451 PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 452 PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 453 } 454 455 return 0; 456 } 457 458 PetscErrorCode ConvergenceTest(KSP ksp, PetscInt it, PetscReal nrm, KSPConvergedReason *reason, void *ctx) 459 { 460 PetscFunctionBegin; 461 *reason = KSP_CONVERGED_ITERATING; 462 if (it > 1) { 463 *reason = KSP_CONVERGED_ITS; 464 PetscCall(PetscInfo(NULL, "User provided convergence test returning after 2 iterations\n")); 465 } 466 PetscFunctionReturn(0); 467 } 468 469 PetscErrorCode ConvergenceDestroy(void *ctx) 470 { 471 PetscFunctionBegin; 472 PetscCall(PetscInfo(NULL, "User provided convergence destroy called\n")); 473 PetscCall(PetscFree(ctx)); 474 PetscFunctionReturn(0); 475 } 476 477 PetscErrorCode postcheck(SNES snes, Vec x, Vec y, Vec w, PetscBool *changed_y, PetscBool *changed_w, void *ctx) 478 { 479 PetscReal norm; 480 Vec tmp; 481 482 PetscFunctionBegin; 483 PetscCall(VecDuplicate(x, &tmp)); 484 PetscCall(VecWAXPY(tmp, -1.0, x, w)); 485 PetscCall(VecNorm(tmp, NORM_2, &norm)); 486 PetscCall(VecDestroy(&tmp)); 487 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Norm of search step %g\n", (double)norm)); 488 PetscFunctionReturn(0); 489 } 490 491 /*TEST 492 493 build: 494 requires: !single 495 496 test: 497 args: -ksp_gmres_cgs_refinement_type refine_always 498 499 test: 500 suffix: 2 501 args: -snes_monitor_short -snes_type newtontr -ksp_gmres_cgs_refinement_type refine_always 502 503 test: 504 suffix: 2a 505 filter: grep -i KSPConvergedDefault > /dev/null && echo "Found KSPConvergedDefault" 506 args: -snes_monitor_short -snes_type newtontr -ksp_gmres_cgs_refinement_type refine_always -info 507 requires: defined(PETSC_USE_INFO) 508 509 test: 510 suffix: 2b 511 filter: grep -i "User provided convergence test" > /dev/null && echo "Found User provided convergence test" 512 args: -snes_monitor_short -snes_type newtontr -ksp_gmres_cgs_refinement_type refine_always -use_convergence_test -info 513 requires: defined(PETSC_USE_INFO) 514 515 test: 516 suffix: 3 517 args: -snes_monitor_short -mat_coloring_type sl -snes_fd_coloring -mx 8 -my 11 -ksp_gmres_cgs_refinement_type refine_always 518 519 test: 520 suffix: 4 521 args: -pc -par 6.807 -snes_monitor -snes_converged_reason 522 523 TEST*/ 524