1c4762a1bSJed Brown! 2c4762a1bSJed Brown! Description: This example solves a nonlinear system in parallel with SNES. 3c4762a1bSJed Brown! We solve the Bratu (SFI - solid fuel ignition) problem in a 2D rectangular 4c4762a1bSJed Brown! domain, using distributed arrays (DMDAs) to partition the parallel grid. 5c4762a1bSJed Brown! The command line options include: 6c4762a1bSJed Brown! -par <param>, where <param> indicates the nonlinearity of the problem 7c4762a1bSJed Brown! problem SFI: <parameter> = Bratu parameter (0 <= par <= 6.81) 8c4762a1bSJed Brown! 9c4762a1bSJed Brown! 10c4762a1bSJed Brown 11c4762a1bSJed Brown! 12c4762a1bSJed Brown! -------------------------------------------------------------------------- 13c4762a1bSJed Brown! 14c4762a1bSJed Brown! Solid Fuel Ignition (SFI) problem. This problem is modeled by 15c4762a1bSJed Brown! the partial differential equation 16c4762a1bSJed Brown! 17c4762a1bSJed Brown! -Laplacian u - lambda*exp(u) = 0, 0 < x,y < 1, 18c4762a1bSJed Brown! 19c4762a1bSJed Brown! with boundary conditions 20c4762a1bSJed Brown! 21c4762a1bSJed Brown! u = 0 for x = 0, x = 1, y = 0, y = 1. 22c4762a1bSJed Brown! 23c4762a1bSJed Brown! A finite difference approximation with the usual 5-point stencil 24c4762a1bSJed Brown! is used to discretize the boundary value problem to obtain a nonlinear 25c4762a1bSJed Brown! system of equations. 26c4762a1bSJed Brown! 27c4762a1bSJed Brown! -------------------------------------------------------------------------- 28c4762a1bSJed Brown 29c4762a1bSJed Brown program main 30c4762a1bSJed Brown#include <petsc/finclude/petscsnes.h> 31c4762a1bSJed Brown use petscdmda 32c4762a1bSJed Brown use petscsnes 33c4762a1bSJed Brown implicit none 34c4762a1bSJed Brown! 35c4762a1bSJed Brown! We place common blocks, variable declarations, and other include files 36c4762a1bSJed Brown! needed for this code in the single file ex5f.h. We then need to include 37c4762a1bSJed Brown! only this file throughout the various routines in this program. See 38c4762a1bSJed Brown! additional comments in the file ex5f.h. 39c4762a1bSJed Brown! 40c4762a1bSJed Brown#include "ex5f.h" 41c4762a1bSJed Brown 42c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 43c4762a1bSJed Brown! Variable declarations 44c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 45c4762a1bSJed Brown! 46c4762a1bSJed Brown! Variables: 47c4762a1bSJed Brown! snes - nonlinear solver 48c4762a1bSJed Brown! x, r - solution, residual vectors 49c4762a1bSJed Brown! its - iterations for convergence 50c4762a1bSJed Brown! 51c4762a1bSJed Brown! See additional variable declarations in the file ex5f.h 52c4762a1bSJed Brown! 53c4762a1bSJed Brown SNES snes 54c4762a1bSJed Brown Vec x,r 55c4762a1bSJed Brown PetscInt its,i1,i4 56c4762a1bSJed Brown PetscErrorCode ierr 57c4762a1bSJed Brown PetscReal lambda_max,lambda_min 58c4762a1bSJed Brown PetscBool flg 59c4762a1bSJed Brown DM da 60c4762a1bSJed Brown 61c4762a1bSJed Brown! Note: Any user-defined Fortran routines (such as FormJacobianLocal) 62c4762a1bSJed Brown! MUST be declared as external. 63c4762a1bSJed Brown 64c4762a1bSJed Brown external FormInitialGuess 65c4762a1bSJed Brown external FormFunctionLocal,FormJacobianLocal 66c4762a1bSJed Brown external MySNESConverged 67c4762a1bSJed Brown 68c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 69c4762a1bSJed Brown! Initialize program 70c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 71c4762a1bSJed Brown 72*d8606c27SBarry Smith PetscCallA(PetscInitialize(ierr)) 73*d8606c27SBarry Smith PetscCallMPIA(MPI_Comm_size(PETSC_COMM_WORLD,size,ierr)) 74*d8606c27SBarry Smith PetscCallMPIA(MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)) 75c4762a1bSJed Brown 76c4762a1bSJed Brown! Initialize problem parameters 77c4762a1bSJed Brown 78c4762a1bSJed Brown i1 = 1 79c4762a1bSJed Brown i4 = 4 80c4762a1bSJed Brown lambda_max = 6.81 81c4762a1bSJed Brown lambda_min = 0.0 82c4762a1bSJed Brown lambda = 6.0 83*d8606c27SBarry Smith PetscCallA(PetscOptionsGetReal(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-par',lambda,PETSC_NULL_BOOL,ierr)) 84c4762a1bSJed Brown! this statement is split into multiple-lines to keep lines under 132 char limit - required by 'make check' 85c4762a1bSJed Brown if (lambda .ge. lambda_max .or. lambda .le. lambda_min) then 86c4762a1bSJed Brown ierr = PETSC_ERR_ARG_OUTOFRANGE; SETERRA(PETSC_COMM_WORLD,ierr,'Lambda') 87c4762a1bSJed Brown endif 88c4762a1bSJed Brown 89c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 90c4762a1bSJed Brown! Create nonlinear solver context 91c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 92c4762a1bSJed Brown 93*d8606c27SBarry Smith PetscCallA(SNESCreate(PETSC_COMM_WORLD,snes,ierr)) 94c4762a1bSJed Brown 95c4762a1bSJed Brown! Set convergence test routine if desired 96c4762a1bSJed Brown 97*d8606c27SBarry Smith PetscCallA(PetscOptionsHasName(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-my_snes_convergence',flg,ierr)) 98c4762a1bSJed Brown if (flg) then 99*d8606c27SBarry Smith PetscCallA(SNESSetConvergenceTest(snes,MySNESConverged,0,PETSC_NULL_FUNCTION,ierr)) 100c4762a1bSJed Brown endif 101c4762a1bSJed Brown 102c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 103c4762a1bSJed Brown! Create vector data structures; set function evaluation routine 104c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 105c4762a1bSJed Brown 106c4762a1bSJed Brown! Create distributed array (DMDA) to manage parallel grid and vectors 107c4762a1bSJed Brown 108c4762a1bSJed Brown! This really needs only the star-type stencil, but we use the box 109c4762a1bSJed Brown! stencil temporarily. 110*d8606c27SBarry Smith PetscCallA(DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,i4,i4,PETSC_DECIDE,PETSC_DECIDE,i1,i1,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,da,ierr)) 111*d8606c27SBarry Smith PetscCallA(DMSetFromOptions(da,ierr)) 112*d8606c27SBarry Smith PetscCallA(DMSetUp(da,ierr)) 113c4762a1bSJed Brown 114c4762a1bSJed Brown! Extract global and local vectors from DMDA; then duplicate for remaining 115c4762a1bSJed Brown! vectors that are the same types 116c4762a1bSJed Brown 117*d8606c27SBarry Smith PetscCallA(DMCreateGlobalVector(da,x,ierr)) 118*d8606c27SBarry Smith PetscCallA(VecDuplicate(x,r,ierr)) 119c4762a1bSJed Brown 120c4762a1bSJed Brown! Get local grid boundaries (for 2-dimensional DMDA) 121c4762a1bSJed Brown 122*d8606c27SBarry Smith PetscCallA(DMDAGetInfo(da,PETSC_NULL_INTEGER,mx,my,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,ierr)) 123*d8606c27SBarry Smith PetscCallA(DMDAGetCorners(da,xs,ys,PETSC_NULL_INTEGER,xm,ym,PETSC_NULL_INTEGER,ierr)) 124*d8606c27SBarry Smith PetscCallA(DMDAGetGhostCorners(da,gxs,gys,PETSC_NULL_INTEGER,gxm,gym,PETSC_NULL_INTEGER,ierr)) 125c4762a1bSJed Brown 126c4762a1bSJed Brown! Here we shift the starting indices up by one so that we can easily 127c4762a1bSJed Brown! use the Fortran convention of 1-based indices (rather 0-based indices). 128c4762a1bSJed Brown 129c4762a1bSJed Brown xs = xs+1 130c4762a1bSJed Brown ys = ys+1 131c4762a1bSJed Brown gxs = gxs+1 132c4762a1bSJed Brown gys = gys+1 133c4762a1bSJed Brown 134c4762a1bSJed Brown ye = ys+ym-1 135c4762a1bSJed Brown xe = xs+xm-1 136c4762a1bSJed Brown gye = gys+gym-1 137c4762a1bSJed Brown gxe = gxs+gxm-1 138c4762a1bSJed Brown 139c4762a1bSJed Brown! Set function evaluation routine and vector 140c4762a1bSJed Brown 141*d8606c27SBarry Smith PetscCallA(DMDASNESSetFunctionLocal(da,INSERT_VALUES,FormFunctionLocal,da,ierr)) 142*d8606c27SBarry Smith PetscCallA(DMDASNESSetJacobianLocal(da,FormJacobianLocal,da,ierr)) 143*d8606c27SBarry Smith PetscCallA(SNESSetDM(snes,da,ierr)) 144c4762a1bSJed Brown 145c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 146c4762a1bSJed Brown! Customize nonlinear solver; set runtime options 147c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 148c4762a1bSJed Brown 149c4762a1bSJed Brown! Set runtime options (e.g., -snes_monitor -snes_rtol <rtol> -ksp_type <type>) 150c4762a1bSJed Brown 151*d8606c27SBarry Smith PetscCallA(SNESSetFromOptions(snes,ierr)) 152c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 153c4762a1bSJed Brown! Evaluate initial guess; then solve nonlinear system. 154c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 155c4762a1bSJed Brown 156c4762a1bSJed Brown! Note: The user should initialize the vector, x, with the initial guess 157c4762a1bSJed Brown! for the nonlinear solver prior to calling SNESSolve(). In particular, 158c4762a1bSJed Brown! to employ an initial guess of zero, the user should explicitly set 159c4762a1bSJed Brown! this vector to zero by calling VecSet(). 160c4762a1bSJed Brown 161*d8606c27SBarry Smith PetscCallA(FormInitialGuess(x,ierr)) 162*d8606c27SBarry Smith PetscCallA(SNESSolve(snes,PETSC_NULL_VEC,x,ierr)) 163*d8606c27SBarry Smith PetscCallA(SNESGetIterationNumber(snes,its,ierr)) 164c4762a1bSJed Brown if (rank .eq. 0) then 165c4762a1bSJed Brown write(6,100) its 166c4762a1bSJed Brown endif 167c4762a1bSJed Brown 100 format('Number of SNES iterations = ',i5) 168c4762a1bSJed Brown 169c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 170c4762a1bSJed Brown! Free work space. All PETSc objects should be destroyed when they 171c4762a1bSJed Brown! are no longer needed. 172c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 173c4762a1bSJed Brown 174*d8606c27SBarry Smith PetscCallA(VecDestroy(x,ierr)) 175*d8606c27SBarry Smith PetscCallA(VecDestroy(r,ierr)) 176*d8606c27SBarry Smith PetscCallA(SNESDestroy(snes,ierr)) 177*d8606c27SBarry Smith PetscCallA(DMDestroy(da,ierr)) 178*d8606c27SBarry Smith PetscCallA(PetscFinalize(ierr)) 179c4762a1bSJed Brown end 180c4762a1bSJed Brown 181c4762a1bSJed Brown! --------------------------------------------------------------------- 182c4762a1bSJed Brown! 183c4762a1bSJed Brown! FormInitialGuess - Forms initial approximation. 184c4762a1bSJed Brown! 185c4762a1bSJed Brown! Input Parameters: 186c4762a1bSJed Brown! X - vector 187c4762a1bSJed Brown! 188c4762a1bSJed Brown! Output Parameter: 189c4762a1bSJed Brown! X - vector 190c4762a1bSJed Brown! 191c4762a1bSJed Brown! Notes: 192c4762a1bSJed Brown! This routine serves as a wrapper for the lower-level routine 193c4762a1bSJed Brown! "ApplicationInitialGuess", where the actual computations are 194c4762a1bSJed Brown! done using the standard Fortran style of treating the local 195c4762a1bSJed Brown! vector data as a multidimensional array over the local mesh. 196c4762a1bSJed Brown! This routine merely handles ghost point scatters and accesses 197c4762a1bSJed Brown! the local vector data via VecGetArray() and VecRestoreArray(). 198c4762a1bSJed Brown! 199c4762a1bSJed Brown subroutine FormInitialGuess(X,ierr) 200c4762a1bSJed Brown use petscsnes 201c4762a1bSJed Brown implicit none 202c4762a1bSJed Brown 203c4762a1bSJed Brown#include "ex5f.h" 204c4762a1bSJed Brown 205c4762a1bSJed Brown! Input/output variables: 206c4762a1bSJed Brown Vec X 207c4762a1bSJed Brown PetscErrorCode ierr 208c4762a1bSJed Brown 209c4762a1bSJed Brown! Declarations for use with local arrays: 210c4762a1bSJed Brown PetscScalar lx_v(0:1) 211c4762a1bSJed Brown PetscOffset lx_i 212c4762a1bSJed Brown 213c4762a1bSJed Brown ierr = 0 214c4762a1bSJed Brown 215c4762a1bSJed Brown! Get a pointer to vector data. 216c4762a1bSJed Brown! - For default PETSc vectors, VecGetArray() returns a pointer to 217c4762a1bSJed Brown! the data array. Otherwise, the routine is implementation dependent. 218c4762a1bSJed Brown! - You MUST call VecRestoreArray() when you no longer need access to 219c4762a1bSJed Brown! the array. 220c4762a1bSJed Brown! - Note that the Fortran interface to VecGetArray() differs from the 221c4762a1bSJed Brown! C version. See the users manual for details. 222c4762a1bSJed Brown 223*d8606c27SBarry Smith PetscCall(VecGetArray(X,lx_v,lx_i,ierr)) 224c4762a1bSJed Brown 225c4762a1bSJed Brown! Compute initial guess over the locally owned part of the grid 226c4762a1bSJed Brown 227*d8606c27SBarry Smith PetscCall(InitialGuessLocal(lx_v(lx_i),ierr)) 228c4762a1bSJed Brown 229c4762a1bSJed Brown! Restore vector 230c4762a1bSJed Brown 231*d8606c27SBarry Smith PetscCall(VecRestoreArray(X,lx_v,lx_i,ierr)) 232c4762a1bSJed Brown 233c4762a1bSJed Brown return 234c4762a1bSJed Brown end 235c4762a1bSJed Brown 236c4762a1bSJed Brown! --------------------------------------------------------------------- 237c4762a1bSJed Brown! 238c4762a1bSJed Brown! InitialGuessLocal - Computes initial approximation, called by 239c4762a1bSJed Brown! the higher level routine FormInitialGuess(). 240c4762a1bSJed Brown! 241c4762a1bSJed Brown! Input Parameter: 242c4762a1bSJed Brown! x - local vector data 243c4762a1bSJed Brown! 244c4762a1bSJed Brown! Output Parameters: 245c4762a1bSJed Brown! x - local vector data 246c4762a1bSJed Brown! ierr - error code 247c4762a1bSJed Brown! 248c4762a1bSJed Brown! Notes: 249c4762a1bSJed Brown! This routine uses standard Fortran-style computations over a 2-dim array. 250c4762a1bSJed Brown! 251c4762a1bSJed Brown subroutine InitialGuessLocal(x,ierr) 252c4762a1bSJed Brown use petscsnes 253c4762a1bSJed Brown implicit none 254c4762a1bSJed Brown 255c4762a1bSJed Brown#include "ex5f.h" 256c4762a1bSJed Brown 257c4762a1bSJed Brown! Input/output variables: 258c4762a1bSJed Brown PetscScalar x(xs:xe,ys:ye) 259c4762a1bSJed Brown PetscErrorCode ierr 260c4762a1bSJed Brown 261c4762a1bSJed Brown! Local variables: 262c4762a1bSJed Brown PetscInt i,j 263c4762a1bSJed Brown PetscReal temp1,temp,one,hx,hy 264c4762a1bSJed Brown 265c4762a1bSJed Brown! Set parameters 266c4762a1bSJed Brown 267c4762a1bSJed Brown ierr = 0 268c4762a1bSJed Brown one = 1.0 269c4762a1bSJed Brown hx = one/((mx-1)) 270c4762a1bSJed Brown hy = one/((my-1)) 271c4762a1bSJed Brown temp1 = lambda/(lambda + one) 272c4762a1bSJed Brown 273c4762a1bSJed Brown do 20 j=ys,ye 274c4762a1bSJed Brown temp = (min(j-1,my-j))*hy 275c4762a1bSJed Brown do 10 i=xs,xe 276c4762a1bSJed Brown if (i .eq. 1 .or. j .eq. 1 .or. i .eq. mx .or. j .eq. my) then 277c4762a1bSJed Brown x(i,j) = 0.0 278c4762a1bSJed Brown else 279c4762a1bSJed Brown x(i,j) = temp1 * sqrt(min(min(i-1,mx-i)*hx,(temp))) 280c4762a1bSJed Brown endif 281c4762a1bSJed Brown 10 continue 282c4762a1bSJed Brown 20 continue 283c4762a1bSJed Brown 284c4762a1bSJed Brown return 285c4762a1bSJed Brown end 286c4762a1bSJed Brown 287c4762a1bSJed Brown! --------------------------------------------------------------------- 288c4762a1bSJed Brown! 289c4762a1bSJed Brown! FormFunctionLocal - Computes nonlinear function, called by 290c4762a1bSJed Brown! the higher level routine FormFunction(). 291c4762a1bSJed Brown! 292c4762a1bSJed Brown! Input Parameter: 293c4762a1bSJed Brown! x - local vector data 294c4762a1bSJed Brown! 295c4762a1bSJed Brown! Output Parameters: 296c4762a1bSJed Brown! f - local vector data, f(x) 297c4762a1bSJed Brown! ierr - error code 298c4762a1bSJed Brown! 299c4762a1bSJed Brown! Notes: 300c4762a1bSJed Brown! This routine uses standard Fortran-style computations over a 2-dim array. 301c4762a1bSJed Brown! 302c4762a1bSJed Brown! 303c4762a1bSJed Brown subroutine FormFunctionLocal(info,x,f,da,ierr) 304c4762a1bSJed Brown#include <petsc/finclude/petscdmda.h> 305c4762a1bSJed Brown use petscsnes 306c4762a1bSJed Brown implicit none 307c4762a1bSJed Brown 308c4762a1bSJed Brown#include "ex5f.h" 309c4762a1bSJed Brown DM da 310c4762a1bSJed Brown 311c4762a1bSJed Brown! Input/output variables: 312c4762a1bSJed Brown DMDALocalInfo info(DMDA_LOCAL_INFO_SIZE) 313c4762a1bSJed Brown PetscScalar x(gxs:gxe,gys:gye) 314c4762a1bSJed Brown PetscScalar f(xs:xe,ys:ye) 315c4762a1bSJed Brown PetscErrorCode ierr 316c4762a1bSJed Brown 317c4762a1bSJed Brown! Local variables: 318c4762a1bSJed Brown PetscScalar two,one,hx,hy 319c4762a1bSJed Brown PetscScalar hxdhy,hydhx,sc 320c4762a1bSJed Brown PetscScalar u,uxx,uyy 321c4762a1bSJed Brown PetscInt i,j 322c4762a1bSJed Brown 323c4762a1bSJed Brown xs = info(DMDA_LOCAL_INFO_XS)+1 324c4762a1bSJed Brown xe = xs+info(DMDA_LOCAL_INFO_XM)-1 325c4762a1bSJed Brown ys = info(DMDA_LOCAL_INFO_YS)+1 326c4762a1bSJed Brown ye = ys+info(DMDA_LOCAL_INFO_YM)-1 327c4762a1bSJed Brown mx = info(DMDA_LOCAL_INFO_MX) 328c4762a1bSJed Brown my = info(DMDA_LOCAL_INFO_MY) 329c4762a1bSJed Brown 330c4762a1bSJed Brown one = 1.0 331c4762a1bSJed Brown two = 2.0 332c4762a1bSJed Brown hx = one/(mx-1) 333c4762a1bSJed Brown hy = one/(my-1) 334c4762a1bSJed Brown sc = hx*hy*lambda 335c4762a1bSJed Brown hxdhy = hx/hy 336c4762a1bSJed Brown hydhx = hy/hx 337c4762a1bSJed Brown 338c4762a1bSJed Brown! Compute function over the locally owned part of the grid 339c4762a1bSJed Brown 340c4762a1bSJed Brown do 20 j=ys,ye 341c4762a1bSJed Brown do 10 i=xs,xe 342c4762a1bSJed Brown if (i .eq. 1 .or. j .eq. 1 .or. i .eq. mx .or. j .eq. my) then 343c4762a1bSJed Brown f(i,j) = x(i,j) 344c4762a1bSJed Brown else 345c4762a1bSJed Brown u = x(i,j) 346c4762a1bSJed Brown uxx = hydhx * (two*u - x(i-1,j) - x(i+1,j)) 347c4762a1bSJed Brown uyy = hxdhy * (two*u - x(i,j-1) - x(i,j+1)) 348c4762a1bSJed Brown f(i,j) = uxx + uyy - sc*exp(u) 349c4762a1bSJed Brown endif 350c4762a1bSJed Brown 10 continue 351c4762a1bSJed Brown 20 continue 352c4762a1bSJed Brown 353*d8606c27SBarry Smith PetscCall(PetscLogFlops(11.0d0*ym*xm,ierr)) 354c4762a1bSJed Brown 355c4762a1bSJed Brown return 356c4762a1bSJed Brown end 357c4762a1bSJed Brown 358c4762a1bSJed Brown! --------------------------------------------------------------------- 359c4762a1bSJed Brown! 360c4762a1bSJed Brown! FormJacobianLocal - Computes Jacobian matrix, called by 361c4762a1bSJed Brown! the higher level routine FormJacobian(). 362c4762a1bSJed Brown! 363c4762a1bSJed Brown! Input Parameters: 364c4762a1bSJed Brown! x - local vector data 365c4762a1bSJed Brown! 366c4762a1bSJed Brown! Output Parameters: 367c4762a1bSJed Brown! jac - Jacobian matrix 368c4762a1bSJed Brown! jac_prec - optionally different preconditioning matrix (not used here) 369c4762a1bSJed Brown! ierr - error code 370c4762a1bSJed Brown! 371c4762a1bSJed Brown! Notes: 372c4762a1bSJed Brown! This routine uses standard Fortran-style computations over a 2-dim array. 373c4762a1bSJed Brown! 374c4762a1bSJed Brown! Notes: 375c4762a1bSJed Brown! Due to grid point reordering with DMDAs, we must always work 376c4762a1bSJed Brown! with the local grid points, and then transform them to the new 377c4762a1bSJed Brown! global numbering with the "ltog" mapping 378c4762a1bSJed Brown! We cannot work directly with the global numbers for the original 379c4762a1bSJed Brown! uniprocessor grid! 380c4762a1bSJed Brown! 381c4762a1bSJed Brown! Two methods are available for imposing this transformation 382c4762a1bSJed Brown! when setting matrix entries: 383c4762a1bSJed Brown! (A) MatSetValuesLocal(), using the local ordering (including 384c4762a1bSJed Brown! ghost points!) 385c4762a1bSJed Brown! by calling MatSetValuesLocal() 386c4762a1bSJed Brown! (B) MatSetValues(), using the global ordering 387c4762a1bSJed Brown! - Use DMDAGetGlobalIndices() to extract the local-to-global map 388c4762a1bSJed Brown! - Then apply this map explicitly yourself 389c4762a1bSJed Brown! - Set matrix entries using the global ordering by calling 390c4762a1bSJed Brown! MatSetValues() 391c4762a1bSJed Brown! Option (A) seems cleaner/easier in many cases, and is the procedure 392c4762a1bSJed Brown! used in this example. 393c4762a1bSJed Brown! 394c4762a1bSJed Brown subroutine FormJacobianLocal(info,x,A,jac,da,ierr) 395c4762a1bSJed Brown use petscsnes 396c4762a1bSJed Brown implicit none 397c4762a1bSJed Brown 398c4762a1bSJed Brown#include "ex5f.h" 399c4762a1bSJed Brown DM da 400c4762a1bSJed Brown 401c4762a1bSJed Brown! Input/output variables: 402c4762a1bSJed Brown PetscScalar x(gxs:gxe,gys:gye) 403c4762a1bSJed Brown Mat A,jac 404c4762a1bSJed Brown PetscErrorCode ierr 405c4762a1bSJed Brown DMDALocalInfo info(DMDA_LOCAL_INFO_SIZE) 406c4762a1bSJed Brown 407c4762a1bSJed Brown! Local variables: 408c4762a1bSJed Brown PetscInt row,col(5),i,j,i1,i5 409c4762a1bSJed Brown PetscScalar two,one,hx,hy,v(5) 410c4762a1bSJed Brown PetscScalar hxdhy,hydhx,sc 411c4762a1bSJed Brown 412c4762a1bSJed Brown! Set parameters 413c4762a1bSJed Brown 414c4762a1bSJed Brown i1 = 1 415c4762a1bSJed Brown i5 = 5 416c4762a1bSJed Brown one = 1.0 417c4762a1bSJed Brown two = 2.0 418c4762a1bSJed Brown hx = one/(mx-1) 419c4762a1bSJed Brown hy = one/(my-1) 420c4762a1bSJed Brown sc = hx*hy 421c4762a1bSJed Brown hxdhy = hx/hy 422c4762a1bSJed Brown hydhx = hy/hx 423c4762a1bSJed Brown 424c4762a1bSJed Brown! Compute entries for the locally owned part of the Jacobian. 425c4762a1bSJed Brown! - Currently, all PETSc parallel matrix formats are partitioned by 426c4762a1bSJed Brown! contiguous chunks of rows across the processors. 427c4762a1bSJed Brown! - Each processor needs to insert only elements that it owns 428c4762a1bSJed Brown! locally (but any non-local elements will be sent to the 429c4762a1bSJed Brown! appropriate processor during matrix assembly). 430c4762a1bSJed Brown! - Here, we set all entries for a particular row at once. 431c4762a1bSJed Brown! - We can set matrix entries either using either 432c4762a1bSJed Brown! MatSetValuesLocal() or MatSetValues(), as discussed above. 433c4762a1bSJed Brown! - Note that MatSetValues() uses 0-based row and column numbers 434c4762a1bSJed Brown! in Fortran as well as in C. 435c4762a1bSJed Brown 436c4762a1bSJed Brown do 20 j=ys,ye 437c4762a1bSJed Brown row = (j - gys)*gxm + xs - gxs - 1 438c4762a1bSJed Brown do 10 i=xs,xe 439c4762a1bSJed Brown row = row + 1 440c4762a1bSJed Brown! boundary points 441c4762a1bSJed Brown if (i .eq. 1 .or. j .eq. 1 .or. i .eq. mx .or. j .eq. my) then 442c4762a1bSJed Brown! Some f90 compilers need 4th arg to be of same type in both calls 443c4762a1bSJed Brown col(1) = row 444c4762a1bSJed Brown v(1) = one 445*d8606c27SBarry Smith PetscCall(MatSetValuesLocal(jac,i1,row,i1,col,v,INSERT_VALUES,ierr)) 446c4762a1bSJed Brown! interior grid points 447c4762a1bSJed Brown else 448c4762a1bSJed Brown v(1) = -hxdhy 449c4762a1bSJed Brown v(2) = -hydhx 450c4762a1bSJed Brown v(3) = two*(hydhx + hxdhy) - sc*lambda*exp(x(i,j)) 451c4762a1bSJed Brown v(4) = -hydhx 452c4762a1bSJed Brown v(5) = -hxdhy 453c4762a1bSJed Brown col(1) = row - gxm 454c4762a1bSJed Brown col(2) = row - 1 455c4762a1bSJed Brown col(3) = row 456c4762a1bSJed Brown col(4) = row + 1 457c4762a1bSJed Brown col(5) = row + gxm 458*d8606c27SBarry Smith PetscCall(MatSetValuesLocal(jac,i1,row,i5,col,v, INSERT_VALUES,ierr)) 459c4762a1bSJed Brown endif 460c4762a1bSJed Brown 10 continue 461c4762a1bSJed Brown 20 continue 462*d8606c27SBarry Smith PetscCall(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY,ierr)) 463*d8606c27SBarry Smith PetscCall(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY,ierr)) 464c4762a1bSJed Brown if (A .ne. jac) then 465*d8606c27SBarry Smith PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY,ierr)) 466*d8606c27SBarry Smith PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY,ierr)) 467c4762a1bSJed Brown endif 468c4762a1bSJed Brown return 469c4762a1bSJed Brown end 470c4762a1bSJed Brown 471c4762a1bSJed Brown! 472c4762a1bSJed Brown! Simple convergence test based on the infinity norm of the residual being small 473c4762a1bSJed Brown! 474c4762a1bSJed Brown subroutine MySNESConverged(snes,it,xnorm,snorm,fnorm,reason,dummy,ierr) 475c4762a1bSJed Brown use petscsnes 476c4762a1bSJed Brown implicit none 477c4762a1bSJed Brown 478c4762a1bSJed Brown SNES snes 479c4762a1bSJed Brown PetscInt it,dummy 480c4762a1bSJed Brown PetscReal xnorm,snorm,fnorm,nrm 481c4762a1bSJed Brown SNESConvergedReason reason 482c4762a1bSJed Brown Vec f 483c4762a1bSJed Brown PetscErrorCode ierr 484c4762a1bSJed Brown 485*d8606c27SBarry Smith PetscCall(SNESGetFunction(snes,f,PETSC_NULL_FUNCTION,dummy,ierr)) 486*d8606c27SBarry Smith PetscCall(VecNorm(f,NORM_INFINITY,nrm,ierr)) 487c4762a1bSJed Brown if (nrm .le. 1.e-5) reason = SNES_CONVERGED_FNORM_ABS 488c4762a1bSJed Brown 489c4762a1bSJed Brown end 490c4762a1bSJed Brown 491c4762a1bSJed Brown!/*TEST 492c4762a1bSJed Brown! 493c4762a1bSJed Brown! build: 494c4762a1bSJed Brown! requires: !complex !single 495c4762a1bSJed Brown! 496c4762a1bSJed Brown! test: 497c4762a1bSJed Brown! nsize: 4 4988f8b3c79SStefano Zampini! args: -snes_mf -pc_type none -da_processors_x 4 -da_processors_y 1 -snes_monitor_short \ 4998f8b3c79SStefano Zampini! -ksp_gmres_cgs_refinement_type refine_always 500c4762a1bSJed Brown! 501c4762a1bSJed Brown! test: 502c4762a1bSJed Brown! suffix: 2 503c4762a1bSJed Brown! nsize: 4 504c4762a1bSJed Brown! args: -da_processors_x 2 -da_processors_y 2 -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always 505c4762a1bSJed Brown! 506c4762a1bSJed Brown! test: 507c4762a1bSJed Brown! suffix: 3 508c4762a1bSJed Brown! nsize: 3 509c4762a1bSJed Brown! args: -snes_fd -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always 510c4762a1bSJed Brown! 511c4762a1bSJed Brown! test: 512c4762a1bSJed Brown! suffix: 6 513c4762a1bSJed Brown! nsize: 1 514c4762a1bSJed Brown! args: -snes_monitor_short -my_snes_convergence 515c4762a1bSJed Brown! 516c4762a1bSJed Brown!TEST*/ 517