1c4762a1bSJed Brown! 242ce371bSBarry Smith! This example shows how to avoid Fortran line lengths larger than 132 characters. 3*dcb3e689SBarry Smith! It avoids used of certain macros such as PetscCallA() and PetscCheckA() that 4*dcb3e689SBarry Smith! generate very long lines 5*dcb3e689SBarry Smith! 642ce371bSBarry Smith! We recommend starting from src/snes/tutorials/ex5f90.F90 instead of this example 7*dcb3e689SBarry Smith! because that does not have the restricted formatting that makes this version 8*dcb3e689SBarry 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! -------------------------------------------------------------------------- 33dfbbaf82SBarry Smith module ex5fmodule 34dfbbaf82SBarry Smith use petscsnes 35dfbbaf82SBarry Smith use petscdmda 36dfbbaf82SBarry Smith#include <petsc/finclude/petscsnes.h> 37dfbbaf82SBarry Smith#include <petsc/finclude/petscdm.h> 38dfbbaf82SBarry Smith#include <petsc/finclude/petscdmda.h> 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 44dfbbaf82SBarry Smith end module ex5fmodule 45c4762a1bSJed Brown 46c4762a1bSJed Brown program main 47dfbbaf82SBarry Smith use ex5fmodule 48c4762a1bSJed Brown implicit none 49c4762a1bSJed Brown 50c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 51c4762a1bSJed Brown! Variable declarations 52c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 53c4762a1bSJed Brown! 54c4762a1bSJed Brown! Variables: 55c4762a1bSJed Brown! snes - nonlinear solver 56c4762a1bSJed Brown! x, r - solution, residual vectors 57c4762a1bSJed Brown! its - iterations for convergence 58c4762a1bSJed Brown! 59c4762a1bSJed Brown! See additional variable declarations in the file ex5f.h 60c4762a1bSJed Brown! 61c4762a1bSJed Brown SNES snes 62c4762a1bSJed Brown Vec x,r 63c4762a1bSJed Brown PetscInt its,i1,i4 64c4762a1bSJed Brown PetscErrorCode ierr 65c4762a1bSJed Brown PetscReal lambda_max,lambda_min 66c4762a1bSJed Brown PetscBool flg 67c4762a1bSJed Brown DM da 68c4762a1bSJed Brown 69c4762a1bSJed Brown! Note: Any user-defined Fortran routines (such as FormJacobianLocal) 70c4762a1bSJed Brown! MUST be declared as external. 71c4762a1bSJed Brown 72c4762a1bSJed Brown external FormInitialGuess 73c4762a1bSJed Brown external FormFunctionLocal,FormJacobianLocal 74c4762a1bSJed Brown external MySNESConverged 75c4762a1bSJed Brown 76c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 77c4762a1bSJed Brown! Initialize program 78c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 79c4762a1bSJed Brown 8065afa91aSSatish Balay call PetscInitialize(ierr) 8165afa91aSSatish Balay CHKERRA(ierr) 8265afa91aSSatish Balay call MPI_Comm_size(PETSC_COMM_WORLD,size,ierr) 8365afa91aSSatish Balay CHKERRMPIA(ierr) 8465afa91aSSatish Balay call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr) 8565afa91aSSatish Balay CHKERRMPIA(ierr) 86c4762a1bSJed Brown! Initialize problem parameters 87c4762a1bSJed Brown 88c4762a1bSJed Brown i1 = 1 89c4762a1bSJed Brown i4 = 4 90c4762a1bSJed Brown lambda_max = 6.81 91c4762a1bSJed Brown lambda_min = 0.0 92c4762a1bSJed Brown lambda = 6.0 9365afa91aSSatish Balay call PetscOptionsGetReal(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-par',lambda,PETSC_NULL_BOOL,ierr) 9465afa91aSSatish Balay CHKERRA(ierr) 9565afa91aSSatish Balay 96c4762a1bSJed Brown! this statement is split into multiple-lines to keep lines under 132 char limit - required by 'make check' 97c4762a1bSJed Brown if (lambda .ge. lambda_max .or. lambda .le. lambda_min) then 98*dcb3e689SBarry Smith ierr = PETSC_ERR_ARG_OUTOFRANGE; 99*dcb3e689SBarry Smith SETERRA(PETSC_COMM_WORLD,ierr,'Lambda') 100c4762a1bSJed Brown endif 101c4762a1bSJed Brown 102c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 103c4762a1bSJed Brown! Create nonlinear solver context 104c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 105c4762a1bSJed Brown 10665afa91aSSatish Balay call SNESCreate(PETSC_COMM_WORLD,snes,ierr) 10765afa91aSSatish Balay CHKERRA(ierr) 108c4762a1bSJed Brown 109c4762a1bSJed Brown! Set convergence test routine if desired 110c4762a1bSJed Brown 11165afa91aSSatish Balay call PetscOptionsHasName(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-my_snes_convergence',flg,ierr) 11265afa91aSSatish Balay CHKERRA(ierr) 113c4762a1bSJed Brown if (flg) then 11465afa91aSSatish Balay call SNESSetConvergenceTest(snes,MySNESConverged,0,PETSC_NULL_FUNCTION,ierr) 11565afa91aSSatish Balay CHKERRA(ierr) 116c4762a1bSJed Brown endif 117c4762a1bSJed Brown 118c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 119c4762a1bSJed Brown! Create vector data structures; set function evaluation routine 120c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 121c4762a1bSJed Brown 122c4762a1bSJed Brown! Create distributed array (DMDA) to manage parallel grid and vectors 123c4762a1bSJed Brown 12442ce371bSBarry Smith! This really needs only the star-type stencil, but we use the box stencil 12560cf0239SBarry Smith 12665afa91aSSatish Balay call DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,i4,i4,PETSC_DECIDE,PETSC_DECIDE, & 12765afa91aSSatish Balay i1,i1, PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,da,ierr) 12865afa91aSSatish Balay CHKERRA(ierr) 12965afa91aSSatish Balay call DMSetFromOptions(da,ierr) 13065afa91aSSatish Balay CHKERRA(ierr) 13165afa91aSSatish Balay call DMSetUp(da,ierr) 13265afa91aSSatish Balay CHKERRA(ierr) 133c4762a1bSJed Brown 134c4762a1bSJed Brown! Extract global and local vectors from DMDA; then duplicate for remaining 135c4762a1bSJed Brown! vectors that are the same types 136c4762a1bSJed Brown 13765afa91aSSatish Balay call DMCreateGlobalVector(da,x,ierr) 13865afa91aSSatish Balay CHKERRA(ierr) 13965afa91aSSatish Balay call VecDuplicate(x,r,ierr) 14065afa91aSSatish Balay CHKERRA(ierr) 141c4762a1bSJed Brown 142c4762a1bSJed Brown! Get local grid boundaries (for 2-dimensional DMDA) 143c4762a1bSJed Brown 14460cf0239SBarry Smith call DMDAGetInfo(da,PETSC_NULL_INTEGER,mx,my,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER, & 14560cf0239SBarry Smith PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER, & 14660cf0239SBarry Smith PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,ierr) 14765afa91aSSatish Balay CHKERRA(ierr) 14865afa91aSSatish Balay call DMDAGetCorners(da,xs,ys,PETSC_NULL_INTEGER,xm,ym,PETSC_NULL_INTEGER,ierr) 14965afa91aSSatish Balay CHKERRA(ierr) 15065afa91aSSatish Balay call DMDAGetGhostCorners(da,gxs,gys,PETSC_NULL_INTEGER,gxm,gym,PETSC_NULL_INTEGER,ierr) 15165afa91aSSatish Balay CHKERRA(ierr) 152c4762a1bSJed Brown 153c4762a1bSJed Brown! Here we shift the starting indices up by one so that we can easily 154c4762a1bSJed Brown! use the Fortran convention of 1-based indices (rather 0-based indices). 155c4762a1bSJed Brown 156c4762a1bSJed Brown xs = xs+1 157c4762a1bSJed Brown ys = ys+1 158c4762a1bSJed Brown gxs = gxs+1 159c4762a1bSJed Brown gys = gys+1 160c4762a1bSJed Brown 161c4762a1bSJed Brown ye = ys+ym-1 162c4762a1bSJed Brown xe = xs+xm-1 163c4762a1bSJed Brown gye = gys+gym-1 164c4762a1bSJed Brown gxe = gxs+gxm-1 165c4762a1bSJed Brown 166c4762a1bSJed Brown! Set function evaluation routine and vector 167c4762a1bSJed Brown 16865afa91aSSatish Balay call DMDASNESSetFunctionLocal(da,INSERT_VALUES,FormFunctionLocal,da,ierr) 16965afa91aSSatish Balay CHKERRA(ierr) 17065afa91aSSatish Balay call DMDASNESSetJacobianLocal(da,FormJacobianLocal,da,ierr) 17165afa91aSSatish Balay CHKERRA(ierr) 17265afa91aSSatish Balay call SNESSetDM(snes,da,ierr) 17365afa91aSSatish Balay CHKERRA(ierr) 174c4762a1bSJed Brown 175c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 176c4762a1bSJed Brown! Customize nonlinear solver; set runtime options 177c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 178c4762a1bSJed Brown 179c4762a1bSJed Brown! Set runtime options (e.g., -snes_monitor -snes_rtol <rtol> -ksp_type <type>) 180c4762a1bSJed Brown 18165afa91aSSatish Balay call SNESSetFromOptions(snes,ierr) 18265afa91aSSatish Balay CHKERRA(ierr) 183c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 184c4762a1bSJed Brown! Evaluate initial guess; then solve nonlinear system. 185c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 186c4762a1bSJed Brown 187c4762a1bSJed Brown! Note: The user should initialize the vector, x, with the initial guess 188c4762a1bSJed Brown! for the nonlinear solver prior to calling SNESSolve(). In particular, 189c4762a1bSJed Brown! to employ an initial guess of zero, the user should explicitly set 190c4762a1bSJed Brown! this vector to zero by calling VecSet(). 191c4762a1bSJed Brown 19265afa91aSSatish Balay call FormInitialGuess(x,ierr) 19365afa91aSSatish Balay CHKERRA(ierr) 19465afa91aSSatish Balay call SNESSolve(snes,PETSC_NULL_VEC,x,ierr) 19565afa91aSSatish Balay CHKERRA(ierr) 19665afa91aSSatish Balay call SNESGetIterationNumber(snes,its,ierr) 19765afa91aSSatish Balay CHKERRA(ierr) 198c4762a1bSJed Brown if (rank .eq. 0) then 199c4762a1bSJed Brown write(6,100) its 200c4762a1bSJed Brown endif 201c4762a1bSJed Brown 100 format('Number of SNES iterations = ',i5) 202c4762a1bSJed Brown 203c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 204c4762a1bSJed Brown! Free work space. All PETSc objects should be destroyed when they 205c4762a1bSJed Brown! are no longer needed. 206c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 207c4762a1bSJed Brown 20865afa91aSSatish Balay call VecDestroy(x,ierr) 20965afa91aSSatish Balay CHKERRA(ierr) 21065afa91aSSatish Balay call VecDestroy(r,ierr) 21165afa91aSSatish Balay CHKERRA(ierr) 21265afa91aSSatish Balay call SNESDestroy(snes,ierr) 21365afa91aSSatish Balay CHKERRA(ierr) 21465afa91aSSatish Balay call DMDestroy(da,ierr) 21565afa91aSSatish Balay CHKERRA(ierr) 21665afa91aSSatish Balay call PetscFinalize(ierr) 21765afa91aSSatish Balay CHKERRA(ierr) 218c4762a1bSJed Brown end 219c4762a1bSJed Brown 220c4762a1bSJed Brown! --------------------------------------------------------------------- 221c4762a1bSJed Brown! 222c4762a1bSJed Brown! FormInitialGuess - Forms initial approximation. 223c4762a1bSJed Brown! 224c4762a1bSJed Brown! Input Parameters: 225c4762a1bSJed Brown! X - vector 226c4762a1bSJed Brown! 227c4762a1bSJed Brown! Output Parameter: 228c4762a1bSJed Brown! X - vector 229c4762a1bSJed Brown! 230c4762a1bSJed Brown! Notes: 231c4762a1bSJed Brown! This routine serves as a wrapper for the lower-level routine 232c4762a1bSJed Brown! "ApplicationInitialGuess", where the actual computations are 233c4762a1bSJed Brown! done using the standard Fortran style of treating the local 234c4762a1bSJed Brown! vector data as a multidimensional array over the local mesh. 235c4762a1bSJed Brown! This routine merely handles ghost point scatters and accesses 23642ce371bSBarry Smith! the local vector data via VecGetArrayF90() and VecRestoreArrayF90(). 237c4762a1bSJed Brown! 238c4762a1bSJed Brown subroutine FormInitialGuess(X,ierr) 239dfbbaf82SBarry Smith use ex5fmodule 240c4762a1bSJed Brown implicit none 241c4762a1bSJed Brown 242c4762a1bSJed Brown! Input/output variables: 243c4762a1bSJed Brown Vec X 244c4762a1bSJed Brown PetscErrorCode ierr 245c4762a1bSJed Brown! Declarations for use with local arrays: 24642ce371bSBarry Smith PetscScalar, pointer :: lx_v(:) 247c4762a1bSJed Brown 248c4762a1bSJed Brown ierr = 0 249c4762a1bSJed Brown 250c4762a1bSJed Brown! Get a pointer to vector data. 251c4762a1bSJed Brown! - For default PETSc vectors, VecGetArray() returns a pointer to 252c4762a1bSJed Brown! the data array. Otherwise, the routine is implementation dependent. 253c4762a1bSJed Brown! - You MUST call VecRestoreArray() when you no longer need access to 254c4762a1bSJed Brown! the array. 255c4762a1bSJed Brown! - Note that the Fortran interface to VecGetArray() differs from the 256c4762a1bSJed Brown! C version. See the users manual for details. 257c4762a1bSJed Brown 25842ce371bSBarry Smith call VecGetArrayF90(X,lx_v,ierr) 25965afa91aSSatish Balay CHKERRQ(ierr) 260c4762a1bSJed Brown 261c4762a1bSJed Brown! Compute initial guess over the locally owned part of the grid 262c4762a1bSJed Brown 26342ce371bSBarry Smith call InitialGuessLocal(lx_v,ierr) 26465afa91aSSatish Balay CHKERRQ(ierr) 265c4762a1bSJed Brown 266c4762a1bSJed Brown! Restore vector 267c4762a1bSJed Brown 26842ce371bSBarry Smith call VecRestoreArrayF90(X,lx_v,ierr) 26965afa91aSSatish Balay CHKERRQ(ierr) 270c4762a1bSJed Brown 271c4762a1bSJed Brown return 272c4762a1bSJed Brown end 273c4762a1bSJed Brown 274c4762a1bSJed Brown! --------------------------------------------------------------------- 275c4762a1bSJed Brown! 276c4762a1bSJed Brown! InitialGuessLocal - Computes initial approximation, called by 277c4762a1bSJed Brown! the higher level routine FormInitialGuess(). 278c4762a1bSJed Brown! 279c4762a1bSJed Brown! Input Parameter: 280c4762a1bSJed Brown! x - local vector data 281c4762a1bSJed Brown! 282c4762a1bSJed Brown! Output Parameters: 283c4762a1bSJed Brown! x - local vector data 284c4762a1bSJed Brown! ierr - error code 285c4762a1bSJed Brown! 286c4762a1bSJed Brown! Notes: 287c4762a1bSJed Brown! This routine uses standard Fortran-style computations over a 2-dim array. 288c4762a1bSJed Brown! 289c4762a1bSJed Brown subroutine InitialGuessLocal(x,ierr) 290dfbbaf82SBarry Smith use ex5fmodule 291c4762a1bSJed Brown implicit none 292c4762a1bSJed Brown 293c4762a1bSJed Brown! Input/output variables: 294c4762a1bSJed Brown PetscScalar x(xs:xe,ys:ye) 295c4762a1bSJed Brown PetscErrorCode ierr 296c4762a1bSJed Brown 297c4762a1bSJed Brown! Local variables: 298c4762a1bSJed Brown PetscInt i,j 299c4762a1bSJed Brown PetscReal temp1,temp,one,hx,hy 300c4762a1bSJed Brown 301c4762a1bSJed Brown! Set parameters 302c4762a1bSJed Brown 303c4762a1bSJed Brown ierr = 0 304c4762a1bSJed Brown one = 1.0 305c4762a1bSJed Brown hx = one/((mx-1)) 306c4762a1bSJed Brown hy = one/((my-1)) 307c4762a1bSJed Brown temp1 = lambda/(lambda + one) 308c4762a1bSJed Brown 309c4762a1bSJed Brown do 20 j=ys,ye 310c4762a1bSJed Brown temp = (min(j-1,my-j))*hy 311c4762a1bSJed Brown do 10 i=xs,xe 312c4762a1bSJed Brown if (i .eq. 1 .or. j .eq. 1 .or. i .eq. mx .or. j .eq. my) then 313c4762a1bSJed Brown x(i,j) = 0.0 314c4762a1bSJed Brown else 315c4762a1bSJed Brown x(i,j) = temp1 * sqrt(min(min(i-1,mx-i)*hx,(temp))) 316c4762a1bSJed Brown endif 317c4762a1bSJed Brown 10 continue 318c4762a1bSJed Brown 20 continue 319c4762a1bSJed Brown 320c4762a1bSJed Brown return 321c4762a1bSJed Brown end 322c4762a1bSJed Brown 323c4762a1bSJed Brown! --------------------------------------------------------------------- 324c4762a1bSJed Brown! 325c4762a1bSJed Brown! FormFunctionLocal - Computes nonlinear function, called by 326c4762a1bSJed Brown! the higher level routine FormFunction(). 327c4762a1bSJed Brown! 328c4762a1bSJed Brown! Input Parameter: 329c4762a1bSJed Brown! x - local vector data 330c4762a1bSJed Brown! 331c4762a1bSJed Brown! Output Parameters: 332c4762a1bSJed Brown! f - local vector data, f(x) 333c4762a1bSJed Brown! ierr - error code 334c4762a1bSJed Brown! 335c4762a1bSJed Brown! Notes: 336c4762a1bSJed Brown! This routine uses standard Fortran-style computations over a 2-dim array. 337c4762a1bSJed Brown! 338c4762a1bSJed Brown! 339c4762a1bSJed Brown subroutine FormFunctionLocal(info,x,f,da,ierr) 340dfbbaf82SBarry Smith use ex5fmodule 341c4762a1bSJed Brown implicit none 342c4762a1bSJed Brown 343c4762a1bSJed Brown DM da 344c4762a1bSJed Brown 345c4762a1bSJed Brown! Input/output variables: 346c4762a1bSJed Brown DMDALocalInfo info(DMDA_LOCAL_INFO_SIZE) 347c4762a1bSJed Brown PetscScalar x(gxs:gxe,gys:gye) 348c4762a1bSJed Brown PetscScalar f(xs:xe,ys:ye) 349c4762a1bSJed Brown PetscErrorCode ierr 350c4762a1bSJed Brown 351c4762a1bSJed Brown! Local variables: 352c4762a1bSJed Brown PetscScalar two,one,hx,hy 353c4762a1bSJed Brown PetscScalar hxdhy,hydhx,sc 354c4762a1bSJed Brown PetscScalar u,uxx,uyy 355c4762a1bSJed Brown PetscInt i,j 356c4762a1bSJed Brown 357c4762a1bSJed Brown xs = info(DMDA_LOCAL_INFO_XS)+1 358c4762a1bSJed Brown xe = xs+info(DMDA_LOCAL_INFO_XM)-1 359c4762a1bSJed Brown ys = info(DMDA_LOCAL_INFO_YS)+1 360c4762a1bSJed Brown ye = ys+info(DMDA_LOCAL_INFO_YM)-1 361c4762a1bSJed Brown mx = info(DMDA_LOCAL_INFO_MX) 362c4762a1bSJed Brown my = info(DMDA_LOCAL_INFO_MY) 363c4762a1bSJed Brown 364c4762a1bSJed Brown one = 1.0 365c4762a1bSJed Brown two = 2.0 366c4762a1bSJed Brown hx = one/(mx-1) 367c4762a1bSJed Brown hy = one/(my-1) 368c4762a1bSJed Brown sc = hx*hy*lambda 369c4762a1bSJed Brown hxdhy = hx/hy 370c4762a1bSJed Brown hydhx = hy/hx 371c4762a1bSJed Brown 372c4762a1bSJed Brown! Compute function over the locally owned part of the grid 373c4762a1bSJed Brown 374c4762a1bSJed Brown do 20 j=ys,ye 375c4762a1bSJed Brown do 10 i=xs,xe 376c4762a1bSJed Brown if (i .eq. 1 .or. j .eq. 1 .or. i .eq. mx .or. j .eq. my) then 377c4762a1bSJed Brown f(i,j) = x(i,j) 378c4762a1bSJed Brown else 379c4762a1bSJed Brown u = x(i,j) 380c4762a1bSJed Brown uxx = hydhx * (two*u - x(i-1,j) - x(i+1,j)) 381c4762a1bSJed Brown uyy = hxdhy * (two*u - x(i,j-1) - x(i,j+1)) 382c4762a1bSJed Brown f(i,j) = uxx + uyy - sc*exp(u) 383c4762a1bSJed Brown endif 384c4762a1bSJed Brown 10 continue 385c4762a1bSJed Brown 20 continue 386c4762a1bSJed Brown 38765afa91aSSatish Balay call PetscLogFlops(11.0d0*ym*xm,ierr) 38865afa91aSSatish Balay CHKERRQ(ierr) 389c4762a1bSJed Brown 390c4762a1bSJed Brown return 391c4762a1bSJed Brown end 392c4762a1bSJed Brown 393c4762a1bSJed Brown! --------------------------------------------------------------------- 394c4762a1bSJed Brown! 395c4762a1bSJed Brown! FormJacobianLocal - Computes Jacobian matrix, called by 396c4762a1bSJed Brown! the higher level routine FormJacobian(). 397c4762a1bSJed Brown! 398c4762a1bSJed Brown! Input Parameters: 399c4762a1bSJed Brown! x - local vector data 400c4762a1bSJed Brown! 401c4762a1bSJed Brown! Output Parameters: 402c4762a1bSJed Brown! jac - Jacobian matrix 403c4762a1bSJed Brown! jac_prec - optionally different preconditioning matrix (not used here) 404c4762a1bSJed Brown! ierr - error code 405c4762a1bSJed Brown! 406c4762a1bSJed Brown! Notes: 407c4762a1bSJed Brown! This routine uses standard Fortran-style computations over a 2-dim array. 408c4762a1bSJed Brown! 409c4762a1bSJed Brown! Notes: 410c4762a1bSJed Brown! Due to grid point reordering with DMDAs, we must always work 411c4762a1bSJed Brown! with the local grid points, and then transform them to the new 412c4762a1bSJed Brown! global numbering with the "ltog" mapping 413c4762a1bSJed Brown! We cannot work directly with the global numbers for the original 414c4762a1bSJed Brown! uniprocessor grid! 415c4762a1bSJed Brown! 416c4762a1bSJed Brown! Two methods are available for imposing this transformation 417c4762a1bSJed Brown! when setting matrix entries: 418c4762a1bSJed Brown! (A) MatSetValuesLocal(), using the local ordering (including 419c4762a1bSJed Brown! ghost points!) 420c4762a1bSJed Brown! by calling MatSetValuesLocal() 421c4762a1bSJed Brown! (B) MatSetValues(), using the global ordering 422c4762a1bSJed Brown! - Use DMDAGetGlobalIndices() to extract the local-to-global map 423c4762a1bSJed Brown! - Then apply this map explicitly yourself 424c4762a1bSJed Brown! - Set matrix entries using the global ordering by calling 425c4762a1bSJed Brown! MatSetValues() 426c4762a1bSJed Brown! Option (A) seems cleaner/easier in many cases, and is the procedure 427c4762a1bSJed Brown! used in this example. 428c4762a1bSJed Brown! 429c4762a1bSJed Brown subroutine FormJacobianLocal(info,x,A,jac,da,ierr) 430dfbbaf82SBarry Smith use ex5fmodule 431c4762a1bSJed Brown implicit none 432c4762a1bSJed Brown 433c4762a1bSJed Brown DM da 434c4762a1bSJed Brown 435c4762a1bSJed Brown! Input/output variables: 436c4762a1bSJed Brown PetscScalar x(gxs:gxe,gys:gye) 437c4762a1bSJed Brown Mat A,jac 438c4762a1bSJed Brown PetscErrorCode ierr 439c4762a1bSJed Brown DMDALocalInfo info(DMDA_LOCAL_INFO_SIZE) 440c4762a1bSJed Brown 441c4762a1bSJed Brown! Local variables: 442c4762a1bSJed Brown PetscInt row,col(5),i,j,i1,i5 443c4762a1bSJed Brown PetscScalar two,one,hx,hy,v(5) 444c4762a1bSJed Brown PetscScalar hxdhy,hydhx,sc 445c4762a1bSJed Brown 446c4762a1bSJed Brown! Set parameters 447c4762a1bSJed Brown 448c4762a1bSJed Brown i1 = 1 449c4762a1bSJed Brown i5 = 5 450c4762a1bSJed Brown one = 1.0 451c4762a1bSJed Brown two = 2.0 452c4762a1bSJed Brown hx = one/(mx-1) 453c4762a1bSJed Brown hy = one/(my-1) 454c4762a1bSJed Brown sc = hx*hy 455c4762a1bSJed Brown hxdhy = hx/hy 456c4762a1bSJed Brown hydhx = hy/hx 457c4762a1bSJed Brown 458c4762a1bSJed Brown! Compute entries for the locally owned part of the Jacobian. 459c4762a1bSJed Brown! - Currently, all PETSc parallel matrix formats are partitioned by 460c4762a1bSJed Brown! contiguous chunks of rows across the processors. 461c4762a1bSJed Brown! - Each processor needs to insert only elements that it owns 462c4762a1bSJed Brown! locally (but any non-local elements will be sent to the 463c4762a1bSJed Brown! appropriate processor during matrix assembly). 464c4762a1bSJed Brown! - Here, we set all entries for a particular row at once. 465c4762a1bSJed Brown! - We can set matrix entries either using either 466c4762a1bSJed Brown! MatSetValuesLocal() or MatSetValues(), as discussed above. 467c4762a1bSJed Brown! - Note that MatSetValues() uses 0-based row and column numbers 468c4762a1bSJed Brown! in Fortran as well as in C. 469c4762a1bSJed Brown 470c4762a1bSJed Brown do 20 j=ys,ye 471c4762a1bSJed Brown row = (j - gys)*gxm + xs - gxs - 1 472c4762a1bSJed Brown do 10 i=xs,xe 473c4762a1bSJed Brown row = row + 1 474c4762a1bSJed Brown! boundary points 475c4762a1bSJed Brown if (i .eq. 1 .or. j .eq. 1 .or. i .eq. mx .or. j .eq. my) then 476c4762a1bSJed Brown! Some f90 compilers need 4th arg to be of same type in both calls 477c4762a1bSJed Brown col(1) = row 478c4762a1bSJed Brown v(1) = one 47965afa91aSSatish Balay call MatSetValuesLocal(jac,i1,row,i1,col,v,INSERT_VALUES,ierr) 48065afa91aSSatish Balay CHKERRQ(ierr) 481c4762a1bSJed Brown! interior grid points 482c4762a1bSJed Brown else 483c4762a1bSJed Brown v(1) = -hxdhy 484c4762a1bSJed Brown v(2) = -hydhx 485c4762a1bSJed Brown v(3) = two*(hydhx + hxdhy) - sc*lambda*exp(x(i,j)) 486c4762a1bSJed Brown v(4) = -hydhx 487c4762a1bSJed Brown v(5) = -hxdhy 488c4762a1bSJed Brown col(1) = row - gxm 489c4762a1bSJed Brown col(2) = row - 1 490c4762a1bSJed Brown col(3) = row 491c4762a1bSJed Brown col(4) = row + 1 492c4762a1bSJed Brown col(5) = row + gxm 49365afa91aSSatish Balay call MatSetValuesLocal(jac,i1,row,i5,col,v, INSERT_VALUES,ierr) 49465afa91aSSatish Balay CHKERRQ(ierr) 495c4762a1bSJed Brown endif 496c4762a1bSJed Brown 10 continue 497c4762a1bSJed Brown 20 continue 49865afa91aSSatish Balay call MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY,ierr) 49965afa91aSSatish Balay CHKERRQ(ierr) 50065afa91aSSatish Balay call MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY,ierr) 50165afa91aSSatish Balay CHKERRQ(ierr) 502c4762a1bSJed Brown if (A .ne. jac) then 50365afa91aSSatish Balay call MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY,ierr) 50465afa91aSSatish Balay CHKERRQ(ierr) 50565afa91aSSatish Balay call MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY,ierr) 50665afa91aSSatish Balay CHKERRQ(ierr) 507c4762a1bSJed Brown endif 508c4762a1bSJed Brown return 509c4762a1bSJed Brown end 510c4762a1bSJed Brown 511c4762a1bSJed Brown! 512c4762a1bSJed Brown! Simple convergence test based on the infinity norm of the residual being small 513c4762a1bSJed Brown! 514c4762a1bSJed Brown subroutine MySNESConverged(snes,it,xnorm,snorm,fnorm,reason,dummy,ierr) 515dfbbaf82SBarry Smith use ex5fmodule 516c4762a1bSJed Brown implicit none 517c4762a1bSJed Brown 518c4762a1bSJed Brown SNES snes 519c4762a1bSJed Brown PetscInt it,dummy 520c4762a1bSJed Brown PetscReal xnorm,snorm,fnorm,nrm 521c4762a1bSJed Brown SNESConvergedReason reason 522c4762a1bSJed Brown Vec f 523c4762a1bSJed Brown PetscErrorCode ierr 524c4762a1bSJed Brown 52565afa91aSSatish Balay call SNESGetFunction(snes,f,PETSC_NULL_FUNCTION,dummy,ierr) 52665afa91aSSatish Balay CHKERRQ(ierr) 52765afa91aSSatish Balay call VecNorm(f,NORM_INFINITY,nrm,ierr) 52865afa91aSSatish Balay CHKERRQ(ierr) 529c4762a1bSJed Brown if (nrm .le. 1.e-5) reason = SNES_CONVERGED_FNORM_ABS 530c4762a1bSJed Brown 531c4762a1bSJed Brown end 532c4762a1bSJed Brown 533c4762a1bSJed Brown!/*TEST 534c4762a1bSJed Brown! 535c4762a1bSJed Brown! build: 536c4762a1bSJed Brown! requires: !complex !single 537c4762a1bSJed Brown! 538c4762a1bSJed Brown! test: 539c4762a1bSJed Brown! nsize: 4 5408f8b3c79SStefano Zampini! args: -snes_mf -pc_type none -da_processors_x 4 -da_processors_y 1 -snes_monitor_short \ 5418f8b3c79SStefano Zampini! -ksp_gmres_cgs_refinement_type refine_always 542c4762a1bSJed Brown! 543c4762a1bSJed Brown! test: 544c4762a1bSJed Brown! suffix: 2 545c4762a1bSJed Brown! nsize: 4 546c4762a1bSJed Brown! args: -da_processors_x 2 -da_processors_y 2 -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always 547c4762a1bSJed Brown! 548c4762a1bSJed Brown! test: 549c4762a1bSJed Brown! suffix: 3 550c4762a1bSJed Brown! nsize: 3 551c4762a1bSJed Brown! args: -snes_fd -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always 552c4762a1bSJed Brown! 553c4762a1bSJed Brown! test: 554c4762a1bSJed Brown! suffix: 6 555c4762a1bSJed Brown! nsize: 1 556c4762a1bSJed Brown! args: -snes_monitor_short -my_snes_convergence 557c4762a1bSJed Brown! 558c4762a1bSJed Brown!TEST*/ 559