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