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; user.my = 4; user.param = 6.0; 84 PetscCall(PetscOptionsGetInt(NULL,NULL,"-mx",&user.mx,NULL)); 85 PetscCall(PetscOptionsGetInt(NULL,NULL,"-my",&user.my,NULL)); 86 PetscCall(PetscOptionsGetReal(NULL,NULL,"-par",&user.param,NULL)); 87 PetscCall(PetscOptionsGetBool(NULL,NULL,"-pc",&pc,NULL)); 88 PetscCheck(user.param < bratu_lambda_max && user.param > bratu_lambda_min,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lambda is out of range"); 89 N = user.mx*user.my; 90 PetscCall(PetscOptionsGetBool(NULL,NULL,"-use_convergence_test",&use_convergence_test,NULL)); 91 92 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 93 Create nonlinear solver context 94 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 95 96 PetscCall(SNESCreate(PETSC_COMM_WORLD,&snes)); 97 98 if (pc) { 99 PetscCall(SNESSetType(snes,SNESNEWTONTR)); 100 PetscCall(SNESNewtonTRSetPostCheck(snes, postcheck,NULL)); 101 } 102 103 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 104 Create vector data structures; set function evaluation routine 105 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 106 107 PetscCall(VecCreate(PETSC_COMM_WORLD,&x)); 108 PetscCall(VecSetSizes(x,PETSC_DECIDE,N)); 109 PetscCall(VecSetFromOptions(x)); 110 PetscCall(VecDuplicate(x,&r)); 111 112 /* 113 Set function evaluation routine and vector. Whenever the nonlinear 114 solver needs to evaluate the nonlinear function, it will call this 115 routine. 116 - Note that the final routine argument is the user-defined 117 context that provides application-specific data for the 118 function evaluation routine. 119 */ 120 PetscCall(SNESSetFunction(snes,r,FormFunction,(void*)&user)); 121 122 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 123 Create matrix data structure; set Jacobian evaluation routine 124 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 125 126 /* 127 Create matrix. Here we only approximately preallocate storage space 128 for the Jacobian. See the users manual for a discussion of better 129 techniques for preallocating matrix memory. 130 */ 131 PetscCall(PetscOptionsGetBool(NULL,NULL,"-snes_mf",&matrix_free,NULL)); 132 if (!matrix_free) { 133 PetscBool matrix_free_operator = PETSC_FALSE; 134 PetscCall(PetscOptionsGetBool(NULL,NULL,"-snes_mf_operator",&matrix_free_operator,NULL)); 135 if (matrix_free_operator) matrix_free = PETSC_FALSE; 136 } 137 if (!matrix_free) { 138 PetscCall(MatCreateSeqAIJ(PETSC_COMM_WORLD,N,N,5,NULL,&J)); 139 } 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++) { 250 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"iteration %" PetscInt_FMT ": Linear iterations %" PetscInt_FMT " Function norm = %g\n",i,hist_its[i],(double)history[i])); 251 } 252 } 253 254 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 255 Free work space. All PETSc objects should be destroyed when they 256 are no longer needed. 257 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 258 259 if (!matrix_free) { 260 PetscCall(MatDestroy(&J)); 261 } 262 if (fd_coloring) { 263 PetscCall(MatFDColoringDestroy(&fdcoloring)); 264 } 265 PetscCall(VecDestroy(&x)); 266 PetscCall(VecDestroy(&r)); 267 PetscCall(SNESDestroy(&snes)); 268 PetscCall(PetscFinalize()); 269 return 0; 270 } 271 /* ------------------------------------------------------------------- */ 272 /* 273 FormInitialGuess - Forms initial approximation. 274 275 Input Parameters: 276 user - user-defined application context 277 X - vector 278 279 Output Parameter: 280 X - vector 281 */ 282 PetscErrorCode FormInitialGuess(AppCtx *user,Vec X) 283 { 284 PetscInt i,j,row,mx,my; 285 PetscReal lambda,temp1,temp,hx,hy; 286 PetscScalar *x; 287 288 mx = user->mx; 289 my = user->my; 290 lambda = user->param; 291 292 hx = 1.0 / (PetscReal)(mx-1); 293 hy = 1.0 / (PetscReal)(my-1); 294 295 /* 296 Get a pointer to vector data. 297 - For default PETSc vectors, VecGetArray() returns a pointer to 298 the data array. Otherwise, the routine is implementation dependent. 299 - You MUST call VecRestoreArray() when you no longer need access to 300 the array. 301 */ 302 PetscCall(VecGetArray(X,&x)); 303 temp1 = lambda/(lambda + 1.0); 304 for (j=0; j<my; j++) { 305 temp = (PetscReal)(PetscMin(j,my-j-1))*hy; 306 for (i=0; i<mx; i++) { 307 row = i + j*mx; 308 if (i == 0 || j == 0 || i == mx-1 || j == my-1) { 309 x[row] = 0.0; 310 continue; 311 } 312 x[row] = temp1*PetscSqrtReal(PetscMin((PetscReal)(PetscMin(i,mx-i-1))*hx,temp)); 313 } 314 } 315 316 /* 317 Restore vector 318 */ 319 PetscCall(VecRestoreArray(X,&x)); 320 return 0; 321 } 322 /* ------------------------------------------------------------------- */ 323 /* 324 FormFunction - Evaluates nonlinear function, F(x). 325 326 Input Parameters: 327 . snes - the SNES context 328 . X - input vector 329 . ptr - optional user-defined context, as set by SNESSetFunction() 330 331 Output Parameter: 332 . F - function vector 333 */ 334 PetscErrorCode FormFunction(SNES snes,Vec X,Vec F,void *ptr) 335 { 336 AppCtx *user = (AppCtx*)ptr; 337 PetscInt i,j,row,mx,my; 338 PetscReal two = 2.0,one = 1.0,lambda,hx,hy,hxdhy,hydhx; 339 PetscScalar ut,ub,ul,ur,u,uxx,uyy,sc,*f; 340 const PetscScalar *x; 341 342 mx = user->mx; 343 my = user->my; 344 lambda = user->param; 345 hx = one / (PetscReal)(mx-1); 346 hy = one / (PetscReal)(my-1); 347 sc = hx*hy; 348 hxdhy = hx/hy; 349 hydhx = hy/hx; 350 351 /* 352 Get pointers to vector data 353 */ 354 PetscCall(VecGetArrayRead(X,&x)); 355 PetscCall(VecGetArray(F,&f)); 356 357 /* 358 Compute function 359 */ 360 for (j=0; j<my; j++) { 361 for (i=0; i<mx; i++) { 362 row = i + j*mx; 363 if (i == 0 || j == 0 || i == mx-1 || j == my-1) { 364 f[row] = x[row]; 365 continue; 366 } 367 u = x[row]; 368 ub = x[row - mx]; 369 ul = x[row - 1]; 370 ut = x[row + mx]; 371 ur = x[row + 1]; 372 uxx = (-ur + two*u - ul)*hydhx; 373 uyy = (-ut + two*u - ub)*hxdhy; 374 f[row] = uxx + uyy - sc*lambda*PetscExpScalar(u); 375 } 376 } 377 378 /* 379 Restore vectors 380 */ 381 PetscCall(VecRestoreArrayRead(X,&x)); 382 PetscCall(VecRestoreArray(F,&f)); 383 return 0; 384 } 385 /* ------------------------------------------------------------------- */ 386 /* 387 FormJacobian - Evaluates Jacobian matrix. 388 389 Input Parameters: 390 . snes - the SNES context 391 . x - input vector 392 . ptr - optional user-defined context, as set by SNESSetJacobian() 393 394 Output Parameters: 395 . A - Jacobian matrix 396 . B - optionally different preconditioning matrix 397 . flag - flag indicating matrix structure 398 */ 399 PetscErrorCode FormJacobian(SNES snes,Vec X,Mat J,Mat jac,void *ptr) 400 { 401 AppCtx *user = (AppCtx*)ptr; /* user-defined applicatin context */ 402 PetscInt i,j,row,mx,my,col[5]; 403 PetscScalar two = 2.0,one = 1.0,lambda,v[5],sc; 404 const PetscScalar *x; 405 PetscReal hx,hy,hxdhy,hydhx; 406 407 mx = user->mx; 408 my = user->my; 409 lambda = user->param; 410 hx = 1.0 / (PetscReal)(mx-1); 411 hy = 1.0 / (PetscReal)(my-1); 412 sc = hx*hy; 413 hxdhy = hx/hy; 414 hydhx = hy/hx; 415 416 /* 417 Get pointer to vector data 418 */ 419 PetscCall(VecGetArrayRead(X,&x)); 420 421 /* 422 Compute entries of the Jacobian 423 */ 424 for (j=0; j<my; j++) { 425 for (i=0; i<mx; i++) { 426 row = i + j*mx; 427 if (i == 0 || j == 0 || i == mx-1 || j == my-1) { 428 PetscCall(MatSetValues(jac,1,&row,1,&row,&one,INSERT_VALUES)); 429 continue; 430 } 431 v[0] = -hxdhy; col[0] = row - mx; 432 v[1] = -hydhx; col[1] = row - 1; 433 v[2] = two*(hydhx + hxdhy) - sc*lambda*PetscExpScalar(x[row]); col[2] = row; 434 v[3] = -hydhx; col[3] = row + 1; 435 v[4] = -hxdhy; col[4] = row + mx; 436 PetscCall(MatSetValues(jac,1,&row,5,col,v,INSERT_VALUES)); 437 } 438 } 439 440 /* 441 Restore vector 442 */ 443 PetscCall(VecRestoreArrayRead(X,&x)); 444 445 /* 446 Assemble matrix 447 */ 448 PetscCall(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY)); 449 PetscCall(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY)); 450 451 if (jac != J) { 452 PetscCall(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY)); 453 PetscCall(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY)); 454 } 455 456 return 0; 457 } 458 459 PetscErrorCode ConvergenceTest(KSP ksp,PetscInt it,PetscReal nrm,KSPConvergedReason *reason,void *ctx) 460 { 461 PetscFunctionBegin; 462 *reason = KSP_CONVERGED_ITERATING; 463 if (it > 1) { 464 *reason = KSP_CONVERGED_ITS; 465 PetscCall(PetscInfo(NULL,"User provided convergence test returning after 2 iterations\n")); 466 } 467 PetscFunctionReturn(0); 468 } 469 470 PetscErrorCode ConvergenceDestroy(void* ctx) 471 { 472 PetscFunctionBegin; 473 PetscCall(PetscInfo(NULL,"User provided convergence destroy called\n")); 474 PetscCall(PetscFree(ctx)); 475 PetscFunctionReturn(0); 476 } 477 478 PetscErrorCode postcheck(SNES snes,Vec x,Vec y,Vec w,PetscBool *changed_y,PetscBool *changed_w,void *ctx) 479 { 480 PetscReal norm; 481 Vec tmp; 482 483 PetscFunctionBegin; 484 PetscCall(VecDuplicate(x,&tmp)); 485 PetscCall(VecWAXPY(tmp,-1.0,x,w)); 486 PetscCall(VecNorm(tmp,NORM_2,&norm)); 487 PetscCall(VecDestroy(&tmp)); 488 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Norm of search step %g\n",(double)norm)); 489 PetscFunctionReturn(0); 490 } 491 492 /*TEST 493 494 build: 495 requires: !single 496 497 test: 498 args: -ksp_gmres_cgs_refinement_type refine_always 499 500 test: 501 suffix: 2 502 args: -snes_monitor_short -snes_type newtontr -ksp_gmres_cgs_refinement_type refine_always 503 504 test: 505 suffix: 2a 506 filter: grep -i KSPConvergedDefault > /dev/null && echo "Found KSPConvergedDefault" 507 args: -snes_monitor_short -snes_type newtontr -ksp_gmres_cgs_refinement_type refine_always -info 508 requires: defined(PETSC_USE_INFO) 509 510 test: 511 suffix: 2b 512 filter: grep -i "User provided convergence test" > /dev/null && echo "Found User provided convergence test" 513 args: -snes_monitor_short -snes_type newtontr -ksp_gmres_cgs_refinement_type refine_always -use_convergence_test -info 514 requires: defined(PETSC_USE_INFO) 515 516 test: 517 suffix: 3 518 args: -snes_monitor_short -mat_coloring_type sl -snes_fd_coloring -mx 8 -my 11 -ksp_gmres_cgs_refinement_type refine_always 519 520 test: 521 suffix: 4 522 args: -pc -par 6.807 -snes_monitor -snes_converged_reason 523 524 TEST*/ 525