1c4762a1bSJed Brown! 242ce371bSBarry Smith! This example shows how to avoid Fortran line lengths larger than 132 characters. 3dcb3e689SBarry Smith! It avoids used of certain macros such as PetscCallA() and PetscCheckA() that 4dcb3e689SBarry Smith! generate very long lines 5dcb3e689SBarry Smith! 642ce371bSBarry Smith! We recommend starting from src/snes/tutorials/ex5f90.F90 instead of this example 7dcb3e689SBarry Smith! because that does not have the restricted formatting that makes this version 8dcb3e689SBarry Smith! more difficult to read 942ce371bSBarry Smith! 10c4762a1bSJed Brown! Description: This example solves a nonlinear system in parallel with SNES. 11c4762a1bSJed Brown! We solve the Bratu (SFI - solid fuel ignition) problem in a 2D rectangular 12c4762a1bSJed Brown! domain, using distributed arrays (DMDAs) to partition the parallel grid. 13c4762a1bSJed Brown! The command line options include: 14c4762a1bSJed Brown! -par <param>, where <param> indicates the nonlinearity of the problem 15c4762a1bSJed Brown! problem SFI: <parameter> = Bratu parameter (0 <= par <= 6.81) 16c4762a1bSJed Brown! 17c4762a1bSJed Brown! -------------------------------------------------------------------------- 18c4762a1bSJed Brown! 19c4762a1bSJed Brown! Solid Fuel Ignition (SFI) problem. This problem is modeled by 20c4762a1bSJed Brown! the partial differential equation 21c4762a1bSJed Brown! 22c4762a1bSJed Brown! -Laplacian u - lambda*exp(u) = 0, 0 < x,y < 1, 23c4762a1bSJed Brown! 24c4762a1bSJed Brown! with boundary conditions 25c4762a1bSJed Brown! 26c4762a1bSJed Brown! u = 0 for x = 0, x = 1, y = 0, y = 1. 27c4762a1bSJed Brown! 28c4762a1bSJed Brown! A finite difference approximation with the usual 5-point stencil 29c4762a1bSJed Brown! is used to discretize the boundary value problem to obtain a nonlinear 30c4762a1bSJed Brown! system of equations. 31c4762a1bSJed Brown! 32c4762a1bSJed Brown! -------------------------------------------------------------------------- 33c5e229c2SMartin Diehl#include <petsc/finclude/petscsnes.h> 34c5e229c2SMartin Diehl#include <petsc/finclude/petscdmda.h> 35*01fa2b5aSMartin Diehlmodule ex5fmodule 36dfbbaf82SBarry Smith use petscsnes 37dfbbaf82SBarry Smith use petscdmda 38e7a95102SMartin Diehl implicit none 39dfbbaf82SBarry Smith PetscInt xs, xe, xm, gxs, gxe, gxm 40dfbbaf82SBarry Smith PetscInt ys, ye, ym, gys, gye, gym 41dfbbaf82SBarry Smith PetscInt mx, my 42dfbbaf82SBarry Smith PetscMPIInt rank, size 43dfbbaf82SBarry Smith PetscReal lambda 44e7a95102SMartin Diehlcontains 45e7a95102SMartin Diehl! --------------------------------------------------------------------- 46e7a95102SMartin Diehl! 47e7a95102SMartin Diehl! FormInitialGuess - Forms initial approximation. 48e7a95102SMartin Diehl! 49e7a95102SMartin Diehl! Input Parameters: 50e7a95102SMartin Diehl! X - vector 51e7a95102SMartin Diehl! 52e7a95102SMartin Diehl! Output Parameter: 53e7a95102SMartin Diehl! X - vector 54e7a95102SMartin Diehl! 55e7a95102SMartin Diehl! Notes: 56e7a95102SMartin Diehl! This routine serves as a wrapper for the lower-level routine 57e7a95102SMartin Diehl! "ApplicationInitialGuess", where the actual computations are 58e7a95102SMartin Diehl! done using the standard Fortran style of treating the local 59e7a95102SMartin Diehl! vector data as a multidimensional array over the local mesh. 60e7a95102SMartin Diehl! This routine merely handles ghost point scatters and accesses 61e7a95102SMartin Diehl! the local vector data via VecGetArray() and VecRestoreArray(). 62e7a95102SMartin Diehl! 63e7a95102SMartin Diehl subroutine FormInitialGuess(X, ierr) 64e7a95102SMartin Diehl 65e7a95102SMartin Diehl! Input/output variables: 66e7a95102SMartin Diehl Vec X 67e7a95102SMartin Diehl PetscErrorCode ierr 68e7a95102SMartin Diehl! Declarations for use with local arrays: 69e7a95102SMartin Diehl PetscScalar, pointer :: lx_v(:) 70e7a95102SMartin Diehl 71e7a95102SMartin Diehl ierr = 0 72e7a95102SMartin Diehl 73e7a95102SMartin Diehl! Get a pointer to vector data. 74e7a95102SMartin Diehl! - For default PETSc vectors, VecGetArray() returns a pointer to 75e7a95102SMartin Diehl! the data array. Otherwise, the routine is implementation dependent. 76e7a95102SMartin Diehl! - You MUST call VecRestoreArray() when you no longer need access to 77e7a95102SMartin Diehl! the array. 78e7a95102SMartin Diehl! - Note that the Fortran interface to VecGetArray() differs from the 79e7a95102SMartin Diehl! C version. See the users manual for details. 80e7a95102SMartin Diehl 81e7a95102SMartin Diehl call VecGetArray(X, lx_v, ierr) 82e7a95102SMartin Diehl CHKERRQ(ierr) 83e7a95102SMartin Diehl 84e7a95102SMartin Diehl! Compute initial guess over the locally owned part of the grid 85e7a95102SMartin Diehl 86e7a95102SMartin Diehl call InitialGuessLocal(lx_v, ierr) 87e7a95102SMartin Diehl CHKERRQ(ierr) 88e7a95102SMartin Diehl 89e7a95102SMartin Diehl! Restore vector 90e7a95102SMartin Diehl 91e7a95102SMartin Diehl call VecRestoreArray(X, lx_v, ierr) 92e7a95102SMartin Diehl CHKERRQ(ierr) 93e7a95102SMartin Diehl 94e7a95102SMartin Diehl end 95e7a95102SMartin Diehl 96e7a95102SMartin Diehl! --------------------------------------------------------------------- 97e7a95102SMartin Diehl! 98e7a95102SMartin Diehl! InitialGuessLocal - Computes initial approximation, called by 99e7a95102SMartin Diehl! the higher level routine FormInitialGuess(). 100e7a95102SMartin Diehl! 101e7a95102SMartin Diehl! Input Parameter: 102e7a95102SMartin Diehl! x - local vector data 103e7a95102SMartin Diehl! 104e7a95102SMartin Diehl! Output Parameters: 105e7a95102SMartin Diehl! x - local vector data 106e7a95102SMartin Diehl! ierr - error code 107e7a95102SMartin Diehl! 108e7a95102SMartin Diehl! Notes: 109e7a95102SMartin Diehl! This routine uses standard Fortran-style computations over a 2-dim array. 110e7a95102SMartin Diehl! 111e7a95102SMartin Diehl subroutine InitialGuessLocal(x, ierr) 112e7a95102SMartin Diehl 113e7a95102SMartin Diehl! Input/output variables: 114e7a95102SMartin Diehl PetscScalar x(xs:xe, ys:ye) 115e7a95102SMartin Diehl PetscErrorCode ierr 116e7a95102SMartin Diehl 117e7a95102SMartin Diehl! Local variables: 118e7a95102SMartin Diehl PetscInt i, j 119e7a95102SMartin Diehl PetscReal temp1, temp, one, hx, hy 120e7a95102SMartin Diehl 121e7a95102SMartin Diehl! Set parameters 122e7a95102SMartin Diehl 123e7a95102SMartin Diehl ierr = 0 124e7a95102SMartin Diehl one = 1.0 125e7a95102SMartin Diehl hx = one/((real(mx) - 1)) 126e7a95102SMartin Diehl hy = one/((real(my) - 1)) 127e7a95102SMartin Diehl temp1 = lambda/(lambda + one) 128e7a95102SMartin Diehl 129e7a95102SMartin Diehl do j = ys, ye 130e7a95102SMartin Diehl temp = (real(min(j - 1, my - j)))*hy 131e7a95102SMartin Diehl do i = xs, xe 132e7a95102SMartin Diehl if (i == 1 .or. j == 1 .or. i == mx .or. j == my) then 133e7a95102SMartin Diehl x(i, j) = 0.0 134e7a95102SMartin Diehl else 135e7a95102SMartin Diehl x(i, j) = temp1*sqrt(min(real(min(i - 1, mx - i))*hx, (temp))) 136e7a95102SMartin Diehl end if 137e7a95102SMartin Diehl end do 138e7a95102SMartin Diehl end do 139e7a95102SMartin Diehl 140e7a95102SMartin Diehl end 141e7a95102SMartin Diehl 142e7a95102SMartin Diehl! --------------------------------------------------------------------- 143e7a95102SMartin Diehl! 144e7a95102SMartin Diehl! FormFunctionLocal - Computes nonlinear function, called by 145e7a95102SMartin Diehl! the higher level routine FormFunction(). 146e7a95102SMartin Diehl! 147e7a95102SMartin Diehl! Input Parameter: 148e7a95102SMartin Diehl! x - local vector data 149e7a95102SMartin Diehl! 150e7a95102SMartin Diehl! Output Parameters: 151e7a95102SMartin Diehl! f - local vector data, f(x) 152e7a95102SMartin Diehl! ierr - error code 153e7a95102SMartin Diehl! 154e7a95102SMartin Diehl! Notes: 155e7a95102SMartin Diehl! This routine uses standard Fortran-style computations over a 2-dim array. 156e7a95102SMartin Diehl! 157e7a95102SMartin Diehl! 158e7a95102SMartin Diehl subroutine FormFunctionLocal(info, x, f, da, ierr) 159e7a95102SMartin Diehl 160e7a95102SMartin Diehl DM da 161e7a95102SMartin Diehl 162e7a95102SMartin Diehl! Input/output variables: 163e7a95102SMartin Diehl DMDALocalInfo info 164e7a95102SMartin Diehl PetscScalar x(gxs:gxe, gys:gye) 165e7a95102SMartin Diehl PetscScalar f(xs:xe, ys:ye) 166e7a95102SMartin Diehl PetscErrorCode ierr 167e7a95102SMartin Diehl 168e7a95102SMartin Diehl! Local variables: 169e7a95102SMartin Diehl PetscScalar two, one, hx, hy 170e7a95102SMartin Diehl PetscScalar hxdhy, hydhx, sc 171e7a95102SMartin Diehl PetscScalar u, uxx, uyy 172e7a95102SMartin Diehl PetscInt i, j 173e7a95102SMartin Diehl 174e7a95102SMartin Diehl xs = info%XS + 1 175e7a95102SMartin Diehl xe = xs + info%XM - 1 176e7a95102SMartin Diehl ys = info%YS + 1 177e7a95102SMartin Diehl ye = ys + info%YM - 1 178e7a95102SMartin Diehl mx = info%MX 179e7a95102SMartin Diehl my = info%MY 180e7a95102SMartin Diehl 181e7a95102SMartin Diehl one = 1.0 182e7a95102SMartin Diehl two = 2.0 183e7a95102SMartin Diehl hx = one/(real(mx) - 1) 184e7a95102SMartin Diehl hy = one/(real(my) - 1) 185e7a95102SMartin Diehl sc = hx*hy*lambda 186e7a95102SMartin Diehl hxdhy = hx/hy 187e7a95102SMartin Diehl hydhx = hy/hx 188e7a95102SMartin Diehl 189e7a95102SMartin Diehl! Compute function over the locally owned part of the grid 190e7a95102SMartin Diehl 191e7a95102SMartin Diehl do j = ys, ye 192e7a95102SMartin Diehl do i = xs, xe 193e7a95102SMartin Diehl if (i == 1 .or. j == 1 .or. i == mx .or. j == my) then 194e7a95102SMartin Diehl f(i, j) = x(i, j) 195e7a95102SMartin Diehl else 196e7a95102SMartin Diehl u = x(i, j) 197e7a95102SMartin Diehl uxx = hydhx*(two*u - x(i - 1, j) - x(i + 1, j)) 198e7a95102SMartin Diehl uyy = hxdhy*(two*u - x(i, j - 1) - x(i, j + 1)) 199e7a95102SMartin Diehl f(i, j) = uxx + uyy - sc*exp(u) 200e7a95102SMartin Diehl end if 201e7a95102SMartin Diehl end do 202e7a95102SMartin Diehl end do 203e7a95102SMartin Diehl 204e7a95102SMartin Diehl call PetscLogFlops(11.0d0*ym*xm, ierr) 205e7a95102SMartin Diehl CHKERRQ(ierr) 206e7a95102SMartin Diehl 207e7a95102SMartin Diehl end 208e7a95102SMartin Diehl 209e7a95102SMartin Diehl! --------------------------------------------------------------------- 210e7a95102SMartin Diehl! 211e7a95102SMartin Diehl! FormJacobianLocal - Computes Jacobian matrix, called by 212e7a95102SMartin Diehl! the higher level routine FormJacobian(). 213e7a95102SMartin Diehl! 214e7a95102SMartin Diehl! Input Parameters: 215e7a95102SMartin Diehl! x - local vector data 216e7a95102SMartin Diehl! 217e7a95102SMartin Diehl! Output Parameters: 218e7a95102SMartin Diehl! jac - Jacobian matrix 219e7a95102SMartin Diehl! jac_prec - optionally different matrix used to construct the preconditioner (not used here) 220e7a95102SMartin Diehl! ierr - error code 221e7a95102SMartin Diehl! 222e7a95102SMartin Diehl! Notes: 223e7a95102SMartin Diehl! This routine uses standard Fortran-style computations over a 2-dim array. 224e7a95102SMartin Diehl! 225e7a95102SMartin Diehl! Notes: 226e7a95102SMartin Diehl! Due to grid point reordering with DMDAs, we must always work 227e7a95102SMartin Diehl! with the local grid points, and then transform them to the new 228e7a95102SMartin Diehl! global numbering with the "ltog" mapping 229e7a95102SMartin Diehl! We cannot work directly with the global numbers for the original 230e7a95102SMartin Diehl! uniprocessor grid! 231e7a95102SMartin Diehl! 232e7a95102SMartin Diehl! Two methods are available for imposing this transformation 233e7a95102SMartin Diehl! when setting matrix entries: 234e7a95102SMartin Diehl! (A) MatSetValuesLocal(), using the local ordering (including 235e7a95102SMartin Diehl! ghost points!) 236e7a95102SMartin Diehl! by calling MatSetValuesLocal() 237e7a95102SMartin Diehl! (B) MatSetValues(), using the global ordering 238e7a95102SMartin Diehl! - Use DMDAGetGlobalIndices() to extract the local-to-global map 239e7a95102SMartin Diehl! - Then apply this map explicitly yourself 240e7a95102SMartin Diehl! - Set matrix entries using the global ordering by calling 241e7a95102SMartin Diehl! MatSetValues() 242e7a95102SMartin Diehl! Option (A) seems cleaner/easier in many cases, and is the procedure 243e7a95102SMartin Diehl! used in this example. 244e7a95102SMartin Diehl! 245e7a95102SMartin Diehl subroutine FormJacobianLocal(info, x, A, jac, da, ierr) 246e7a95102SMartin Diehl 247e7a95102SMartin Diehl DM da 248e7a95102SMartin Diehl 249e7a95102SMartin Diehl! Input/output variables: 250e7a95102SMartin Diehl PetscScalar x(gxs:gxe, gys:gye) 251e7a95102SMartin Diehl Mat A, jac 252e7a95102SMartin Diehl PetscErrorCode ierr 253e7a95102SMartin Diehl DMDALocalInfo info 254e7a95102SMartin Diehl 255e7a95102SMartin Diehl! Local variables: 256e7a95102SMartin Diehl PetscInt row, col(5), i, j, i1, i5 257e7a95102SMartin Diehl PetscScalar two, one, hx, hy, v(5) 258e7a95102SMartin Diehl PetscScalar hxdhy, hydhx, sc 259e7a95102SMartin Diehl 260e7a95102SMartin Diehl! Set parameters 261e7a95102SMartin Diehl 262e7a95102SMartin Diehl i1 = 1 263e7a95102SMartin Diehl i5 = 5 264e7a95102SMartin Diehl one = 1.0 265e7a95102SMartin Diehl two = 2.0 266e7a95102SMartin Diehl hx = one/(real(mx) - 1) 267e7a95102SMartin Diehl hy = one/(real(my) - 1) 268e7a95102SMartin Diehl sc = hx*hy 269e7a95102SMartin Diehl hxdhy = hx/hy 270e7a95102SMartin Diehl hydhx = hy/hx 271e7a95102SMartin Diehl! -Wmaybe-uninitialized 272e7a95102SMartin Diehl v = 0.0 273e7a95102SMartin Diehl col = 0 274e7a95102SMartin Diehl 275e7a95102SMartin Diehl! Compute entries for the locally owned part of the Jacobian. 276e7a95102SMartin Diehl! - Currently, all PETSc parallel matrix formats are partitioned by 277e7a95102SMartin Diehl! contiguous chunks of rows across the processors. 278e7a95102SMartin Diehl! - Each processor needs to insert only elements that it owns 279e7a95102SMartin Diehl! locally (but any non-local elements will be sent to the 280e7a95102SMartin Diehl! appropriate processor during matrix assembly). 281e7a95102SMartin Diehl! - Here, we set all entries for a particular row at once. 282e7a95102SMartin Diehl! - We can set matrix entries either using either 283e7a95102SMartin Diehl! MatSetValuesLocal() or MatSetValues(), as discussed above. 284e7a95102SMartin Diehl! - Note that MatSetValues() uses 0-based row and column numbers 285e7a95102SMartin Diehl! in Fortran as well as in C. 286e7a95102SMartin Diehl 287e7a95102SMartin Diehl do j = ys, ye 288e7a95102SMartin Diehl row = (j - gys)*gxm + xs - gxs - 1 289e7a95102SMartin Diehl do i = xs, xe 290e7a95102SMartin Diehl row = row + 1 291e7a95102SMartin Diehl! boundary points 292e7a95102SMartin Diehl if (i == 1 .or. j == 1 .or. i == mx .or. j == my) then 293e7a95102SMartin Diehl! Some f90 compilers need 4th arg to be of same type in both calls 294e7a95102SMartin Diehl col(1) = row 295e7a95102SMartin Diehl v(1) = one 296e7a95102SMartin Diehl call MatSetValuesLocal(jac, i1, [row], i1, [col], [v], INSERT_VALUES, ierr) 297e7a95102SMartin Diehl CHKERRQ(ierr) 298e7a95102SMartin Diehl! interior grid points 299e7a95102SMartin Diehl else 300e7a95102SMartin Diehl v(1) = -hxdhy 301e7a95102SMartin Diehl v(2) = -hydhx 302e7a95102SMartin Diehl v(3) = two*(hydhx + hxdhy) - sc*lambda*exp(x(i, j)) 303e7a95102SMartin Diehl v(4) = -hydhx 304e7a95102SMartin Diehl v(5) = -hxdhy 305e7a95102SMartin Diehl col(1) = row - gxm 306e7a95102SMartin Diehl col(2) = row - 1 307e7a95102SMartin Diehl col(3) = row 308e7a95102SMartin Diehl col(4) = row + 1 309e7a95102SMartin Diehl col(5) = row + gxm 310e7a95102SMartin Diehl call MatSetValuesLocal(jac, i1, [row], i5, [col], [v], INSERT_VALUES, ierr) 311e7a95102SMartin Diehl CHKERRQ(ierr) 312e7a95102SMartin Diehl end if 313e7a95102SMartin Diehl end do 314e7a95102SMartin Diehl end do 315e7a95102SMartin Diehl call MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY, ierr) 316e7a95102SMartin Diehl CHKERRQ(ierr) 317e7a95102SMartin Diehl call MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY, ierr) 318e7a95102SMartin Diehl CHKERRQ(ierr) 319e7a95102SMartin Diehl if (A /= jac) then 320e7a95102SMartin Diehl call MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY, ierr) 321e7a95102SMartin Diehl CHKERRQ(ierr) 322e7a95102SMartin Diehl call MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY, ierr) 323e7a95102SMartin Diehl CHKERRQ(ierr) 324e7a95102SMartin Diehl end if 325e7a95102SMartin Diehl end 326e7a95102SMartin Diehl 327e7a95102SMartin Diehl! 328e7a95102SMartin Diehl! Simple convergence test based on the infinity norm of the residual being small 329e7a95102SMartin Diehl! 330e7a95102SMartin Diehl subroutine MySNESConverged(snes, it, xnorm, snorm, fnorm, reason, dummy, ierr) 331e7a95102SMartin Diehl 332e7a95102SMartin Diehl SNES snes 333e7a95102SMartin Diehl PetscInt it, dummy 334e7a95102SMartin Diehl PetscReal xnorm, snorm, fnorm, nrm 335e7a95102SMartin Diehl SNESConvergedReason reason 336e7a95102SMartin Diehl Vec f 337e7a95102SMartin Diehl PetscErrorCode ierr 338e7a95102SMartin Diehl 339e7a95102SMartin Diehl call SNESGetFunction(snes, f, PETSC_NULL_FUNCTION, dummy, ierr) 340e7a95102SMartin Diehl CHKERRQ(ierr) 341e7a95102SMartin Diehl call VecNorm(f, NORM_INFINITY, nrm, ierr) 342e7a95102SMartin Diehl CHKERRQ(ierr) 343e7a95102SMartin Diehl if (nrm <= 1.e-5) reason = SNES_CONVERGED_FNORM_ABS 344e7a95102SMartin Diehl 345e7a95102SMartin Diehl end 346e7a95102SMartin Diehl 347*01fa2b5aSMartin Diehlend module ex5fmodule 348c4762a1bSJed Brown 349c4762a1bSJed Brownprogram main 350*01fa2b5aSMartin Diehl use ex5fmodule 351c4762a1bSJed Brown implicit none 352c4762a1bSJed Brown 353c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 354c4762a1bSJed Brown! Variable declarations 355c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 356c4762a1bSJed Brown! 357c4762a1bSJed Brown! Variables: 358c4762a1bSJed Brown! snes - nonlinear solver 359c4762a1bSJed Brown! x, r - solution, residual vectors 360c4762a1bSJed Brown! its - iterations for convergence 361c4762a1bSJed Brown! 362c4762a1bSJed Brown! See additional variable declarations in the file ex5f.h 363c4762a1bSJed Brown! 364c4762a1bSJed Brown SNES snes 365c4762a1bSJed Brown Vec x, r 366c4762a1bSJed Brown PetscInt its, i1, i4 367c4762a1bSJed Brown PetscErrorCode ierr 368c4762a1bSJed Brown PetscReal lambda_max, lambda_min 369c4762a1bSJed Brown PetscBool flg 370c4762a1bSJed Brown DM da 371c4762a1bSJed Brown 372c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 373c4762a1bSJed Brown! Initialize program 374c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 375c4762a1bSJed Brown 37665afa91aSSatish Balay call PetscInitialize(ierr) 37765afa91aSSatish Balay CHKERRA(ierr) 37865afa91aSSatish Balay call MPI_Comm_size(PETSC_COMM_WORLD, size, ierr) 37965afa91aSSatish Balay CHKERRMPIA(ierr) 38065afa91aSSatish Balay call MPI_Comm_rank(PETSC_COMM_WORLD, rank, ierr) 38165afa91aSSatish Balay CHKERRMPIA(ierr) 382c4762a1bSJed Brown! Initialize problem parameters 383c4762a1bSJed Brown 384c4762a1bSJed Brown i1 = 1 385c4762a1bSJed Brown i4 = 4 386c4762a1bSJed Brown lambda_max = 6.81 387c4762a1bSJed Brown lambda_min = 0.0 388c4762a1bSJed Brown lambda = 6.0 38965afa91aSSatish Balay call PetscOptionsGetReal(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-par', lambda, PETSC_NULL_BOOL, ierr) 39065afa91aSSatish Balay CHKERRA(ierr) 39165afa91aSSatish Balay 392c4762a1bSJed Brown! this statement is split into multiple-lines to keep lines under 132 char limit - required by 'make check' 3934820e4eaSBarry Smith if (lambda >= lambda_max .or. lambda <= lambda_min) then 394ccfd86f1SBarry Smith ierr = PETSC_ERR_ARG_OUTOFRANGE 395dcb3e689SBarry Smith SETERRA(PETSC_COMM_WORLD, ierr, 'Lambda') 396c4762a1bSJed Brown end if 397c4762a1bSJed Brown 398c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 399c4762a1bSJed Brown! Create nonlinear solver context 400c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 401c4762a1bSJed Brown 40265afa91aSSatish Balay call SNESCreate(PETSC_COMM_WORLD, snes, ierr) 40365afa91aSSatish Balay CHKERRA(ierr) 404c4762a1bSJed Brown 405c4762a1bSJed Brown! Set convergence test routine if desired 406c4762a1bSJed Brown 40765afa91aSSatish Balay call PetscOptionsHasName(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-my_snes_convergence', flg, ierr) 40865afa91aSSatish Balay CHKERRA(ierr) 409c4762a1bSJed Brown if (flg) then 41065afa91aSSatish Balay call SNESSetConvergenceTest(snes, MySNESConverged, 0, PETSC_NULL_FUNCTION, ierr) 41165afa91aSSatish Balay CHKERRA(ierr) 412c4762a1bSJed Brown end if 413c4762a1bSJed Brown 414c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 415c4762a1bSJed Brown! Create vector data structures; set function evaluation routine 416c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 417c4762a1bSJed Brown 418c4762a1bSJed Brown! Create distributed array (DMDA) to manage parallel grid and vectors 419c4762a1bSJed Brown 42042ce371bSBarry Smith! This really needs only the star-type stencil, but we use the box stencil 42160cf0239SBarry Smith 42265afa91aSSatish Balay call DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_STAR, i4, i4, PETSC_DECIDE, PETSC_DECIDE, & 4235d83a8b1SBarry Smith i1, i1, PETSC_NULL_INTEGER_ARRAY, PETSC_NULL_INTEGER_ARRAY, da, ierr) 42465afa91aSSatish Balay CHKERRA(ierr) 42565afa91aSSatish Balay call DMSetFromOptions(da, ierr) 42665afa91aSSatish Balay CHKERRA(ierr) 42765afa91aSSatish Balay call DMSetUp(da, ierr) 42865afa91aSSatish Balay CHKERRA(ierr) 429c4762a1bSJed Brown 430c4762a1bSJed Brown! Extract global and local vectors from DMDA; then duplicate for remaining 431c4762a1bSJed Brown! vectors that are the same types 432c4762a1bSJed Brown 43365afa91aSSatish Balay call DMCreateGlobalVector(da, x, ierr) 43465afa91aSSatish Balay CHKERRA(ierr) 43565afa91aSSatish Balay call VecDuplicate(x, r, ierr) 43665afa91aSSatish Balay CHKERRA(ierr) 437c4762a1bSJed Brown 438c4762a1bSJed Brown! Get local grid boundaries (for 2-dimensional DMDA) 439c4762a1bSJed Brown 44060cf0239SBarry Smith call DMDAGetInfo(da, PETSC_NULL_INTEGER, mx, my, PETSC_NULL_INTEGER, PETSC_NULL_INTEGER, PETSC_NULL_INTEGER, & 441ce78bad3SBarry Smith PETSC_NULL_INTEGER, PETSC_NULL_INTEGER, PETSC_NULL_INTEGER, PETSC_NULL_DMBOUNDARYTYPE, PETSC_NULL_DMBOUNDARYTYPE, & 442ce78bad3SBarry Smith PETSC_NULL_DMBOUNDARYTYPE, PETSC_NULL_DMDASTENCILTYPE, ierr) 44365afa91aSSatish Balay CHKERRA(ierr) 44465afa91aSSatish Balay call DMDAGetCorners(da, xs, ys, PETSC_NULL_INTEGER, xm, ym, PETSC_NULL_INTEGER, ierr) 44565afa91aSSatish Balay CHKERRA(ierr) 44665afa91aSSatish Balay call DMDAGetGhostCorners(da, gxs, gys, PETSC_NULL_INTEGER, gxm, gym, PETSC_NULL_INTEGER, ierr) 44765afa91aSSatish Balay CHKERRA(ierr) 448c4762a1bSJed Brown 449c4762a1bSJed Brown! Here we shift the starting indices up by one so that we can easily 450c4762a1bSJed Brown! use the Fortran convention of 1-based indices (rather 0-based indices). 451c4762a1bSJed Brown 452c4762a1bSJed Brown xs = xs + 1 453c4762a1bSJed Brown ys = ys + 1 454c4762a1bSJed Brown gxs = gxs + 1 455c4762a1bSJed Brown gys = gys + 1 456c4762a1bSJed Brown 457c4762a1bSJed Brown ye = ys + ym - 1 458c4762a1bSJed Brown xe = xs + xm - 1 459c4762a1bSJed Brown gye = gys + gym - 1 460c4762a1bSJed Brown gxe = gxs + gxm - 1 461c4762a1bSJed Brown 462c4762a1bSJed Brown! Set function evaluation routine and vector 463c4762a1bSJed Brown 46465afa91aSSatish Balay call DMDASNESSetFunctionLocal(da, INSERT_VALUES, FormFunctionLocal, da, ierr) 46565afa91aSSatish Balay CHKERRA(ierr) 46665afa91aSSatish Balay call DMDASNESSetJacobianLocal(da, FormJacobianLocal, da, ierr) 46765afa91aSSatish Balay CHKERRA(ierr) 46865afa91aSSatish Balay call SNESSetDM(snes, da, ierr) 46965afa91aSSatish Balay CHKERRA(ierr) 470c4762a1bSJed Brown 471c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 472c4762a1bSJed Brown! Customize nonlinear solver; set runtime options 473c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 474c4762a1bSJed Brown 475c4762a1bSJed Brown! Set runtime options (e.g., -snes_monitor -snes_rtol <rtol> -ksp_type <type>) 476c4762a1bSJed Brown 47765afa91aSSatish Balay call SNESSetFromOptions(snes, ierr) 47865afa91aSSatish Balay CHKERRA(ierr) 479c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 480c4762a1bSJed Brown! Evaluate initial guess; then solve nonlinear system. 481c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 482c4762a1bSJed Brown 483c4762a1bSJed Brown! Note: The user should initialize the vector, x, with the initial guess 484c4762a1bSJed Brown! for the nonlinear solver prior to calling SNESSolve(). In particular, 485c4762a1bSJed Brown! to employ an initial guess of zero, the user should explicitly set 486c4762a1bSJed Brown! this vector to zero by calling VecSet(). 487c4762a1bSJed Brown 48865afa91aSSatish Balay call FormInitialGuess(x, ierr) 48965afa91aSSatish Balay CHKERRA(ierr) 49065afa91aSSatish Balay call SNESSolve(snes, PETSC_NULL_VEC, x, ierr) 49165afa91aSSatish Balay CHKERRA(ierr) 49265afa91aSSatish Balay call SNESGetIterationNumber(snes, its, ierr) 49365afa91aSSatish Balay CHKERRA(ierr) 4944820e4eaSBarry Smith if (rank == 0) then 495c4762a1bSJed Brown write (6, 100) its 496c4762a1bSJed Brown end if 497c4762a1bSJed Brown100 format('Number of SNES iterations = ', i5) 498c4762a1bSJed Brown 499c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 500c4762a1bSJed Brown! Free work space. All PETSc objects should be destroyed when they 501c4762a1bSJed Brown! are no longer needed. 502c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 503c4762a1bSJed Brown 50465afa91aSSatish Balay call VecDestroy(x, ierr) 50565afa91aSSatish Balay CHKERRA(ierr) 50665afa91aSSatish Balay call VecDestroy(r, ierr) 50765afa91aSSatish Balay CHKERRA(ierr) 50865afa91aSSatish Balay call SNESDestroy(snes, ierr) 50965afa91aSSatish Balay CHKERRA(ierr) 51065afa91aSSatish Balay call DMDestroy(da, ierr) 51165afa91aSSatish Balay CHKERRA(ierr) 51265afa91aSSatish Balay call PetscFinalize(ierr) 51365afa91aSSatish Balay CHKERRA(ierr) 514c4762a1bSJed Brownend 515c4762a1bSJed Brown!/*TEST 516c4762a1bSJed Brown! 517c4762a1bSJed Brown! build: 518c4762a1bSJed Brown! requires: !complex !single 519c4762a1bSJed Brown! 520c4762a1bSJed Brown! test: 521c4762a1bSJed Brown! nsize: 4 5228f8b3c79SStefano Zampini! args: -snes_mf -pc_type none -da_processors_x 4 -da_processors_y 1 -snes_monitor_short \ 5238f8b3c79SStefano Zampini! -ksp_gmres_cgs_refinement_type refine_always 524c4762a1bSJed Brown! 525c4762a1bSJed Brown! test: 526c4762a1bSJed Brown! suffix: 2 527c4762a1bSJed Brown! nsize: 4 528c4762a1bSJed Brown! args: -da_processors_x 2 -da_processors_y 2 -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always 529c4762a1bSJed Brown! 530c4762a1bSJed Brown! test: 531c4762a1bSJed Brown! suffix: 3 532c4762a1bSJed Brown! nsize: 3 533c4762a1bSJed Brown! args: -snes_fd -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always 534c4762a1bSJed Brown! 535c4762a1bSJed Brown! test: 536c4762a1bSJed Brown! suffix: 6 537c4762a1bSJed Brown! nsize: 1 538c4762a1bSJed Brown! args: -snes_monitor_short -my_snes_convergence 539c4762a1bSJed Brown! 540c4762a1bSJed Brown!TEST*/ 541