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! -------------------------------------------------------------------------- 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 98dcb3e689SBarry Smith ierr = PETSC_ERR_ARG_OUTOFRANGE; 99dcb3e689SBarry 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, & 1275d83a8b1SBarry Smith i1,i1, PETSC_NULL_INTEGER_ARRAY,PETSC_NULL_INTEGER_ARRAY,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, & 1455d83a8b1SBarry Smith PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_ENUM,PETSC_NULL_ENUM, & 1465d83a8b1SBarry Smith PETSC_NULL_ENUM,PETSC_NULL_ENUM,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 end 272c4762a1bSJed Brown 273c4762a1bSJed Brown! --------------------------------------------------------------------- 274c4762a1bSJed Brown! 275c4762a1bSJed Brown! InitialGuessLocal - Computes initial approximation, called by 276c4762a1bSJed Brown! the higher level routine FormInitialGuess(). 277c4762a1bSJed Brown! 278c4762a1bSJed Brown! Input Parameter: 279c4762a1bSJed Brown! x - local vector data 280c4762a1bSJed Brown! 281c4762a1bSJed Brown! Output Parameters: 282c4762a1bSJed Brown! x - local vector data 283c4762a1bSJed Brown! ierr - error code 284c4762a1bSJed Brown! 285c4762a1bSJed Brown! Notes: 286c4762a1bSJed Brown! This routine uses standard Fortran-style computations over a 2-dim array. 287c4762a1bSJed Brown! 288c4762a1bSJed Brown subroutine InitialGuessLocal(x,ierr) 289dfbbaf82SBarry Smith use ex5fmodule 290c4762a1bSJed Brown implicit none 291c4762a1bSJed Brown 292c4762a1bSJed Brown! Input/output variables: 293c4762a1bSJed Brown PetscScalar x(xs:xe,ys:ye) 294c4762a1bSJed Brown PetscErrorCode ierr 295c4762a1bSJed Brown 296c4762a1bSJed Brown! Local variables: 297c4762a1bSJed Brown PetscInt i,j 298c4762a1bSJed Brown PetscReal temp1,temp,one,hx,hy 299c4762a1bSJed Brown 300c4762a1bSJed Brown! Set parameters 301c4762a1bSJed Brown 302c4762a1bSJed Brown ierr = 0 303c4762a1bSJed Brown one = 1.0 304bb060b03SBarry Smith hx = one/((real(mx)-1)) 305bb060b03SBarry Smith hy = one/((real(my)-1)) 306c4762a1bSJed Brown temp1 = lambda/(lambda + one) 307c4762a1bSJed Brown 308c4762a1bSJed Brown do 20 j=ys,ye 309bb060b03SBarry Smith temp = (real(min(j-1,my-j)))*hy 310c4762a1bSJed Brown do 10 i=xs,xe 311c4762a1bSJed Brown if (i .eq. 1 .or. j .eq. 1 .or. i .eq. mx .or. j .eq. my) then 312c4762a1bSJed Brown x(i,j) = 0.0 313c4762a1bSJed Brown else 314bb060b03SBarry Smith x(i,j) = temp1 * sqrt(min(real(min(i-1,mx-i))*hx,(temp))) 315c4762a1bSJed Brown endif 316c4762a1bSJed Brown 10 continue 317c4762a1bSJed Brown 20 continue 318c4762a1bSJed Brown 319c4762a1bSJed Brown end 320c4762a1bSJed Brown 321c4762a1bSJed Brown! --------------------------------------------------------------------- 322c4762a1bSJed Brown! 323c4762a1bSJed Brown! FormFunctionLocal - Computes nonlinear function, called by 324c4762a1bSJed Brown! the higher level routine FormFunction(). 325c4762a1bSJed Brown! 326c4762a1bSJed Brown! Input Parameter: 327c4762a1bSJed Brown! x - local vector data 328c4762a1bSJed Brown! 329c4762a1bSJed Brown! Output Parameters: 330c4762a1bSJed Brown! f - local vector data, f(x) 331c4762a1bSJed Brown! ierr - error code 332c4762a1bSJed Brown! 333c4762a1bSJed Brown! Notes: 334c4762a1bSJed Brown! This routine uses standard Fortran-style computations over a 2-dim array. 335c4762a1bSJed Brown! 336c4762a1bSJed Brown! 337c4762a1bSJed Brown subroutine FormFunctionLocal(info,x,f,da,ierr) 338dfbbaf82SBarry Smith use ex5fmodule 339c4762a1bSJed Brown implicit none 340c4762a1bSJed Brown 341c4762a1bSJed Brown DM da 342c4762a1bSJed Brown 343c4762a1bSJed Brown! Input/output variables: 344c4762a1bSJed Brown DMDALocalInfo info(DMDA_LOCAL_INFO_SIZE) 345c4762a1bSJed Brown PetscScalar x(gxs:gxe,gys:gye) 346c4762a1bSJed Brown PetscScalar f(xs:xe,ys:ye) 347c4762a1bSJed Brown PetscErrorCode ierr 348c4762a1bSJed Brown 349c4762a1bSJed Brown! Local variables: 350c4762a1bSJed Brown PetscScalar two,one,hx,hy 351c4762a1bSJed Brown PetscScalar hxdhy,hydhx,sc 352c4762a1bSJed Brown PetscScalar u,uxx,uyy 353c4762a1bSJed Brown PetscInt i,j 354c4762a1bSJed Brown 355c4762a1bSJed Brown xs = info(DMDA_LOCAL_INFO_XS)+1 356c4762a1bSJed Brown xe = xs+info(DMDA_LOCAL_INFO_XM)-1 357c4762a1bSJed Brown ys = info(DMDA_LOCAL_INFO_YS)+1 358c4762a1bSJed Brown ye = ys+info(DMDA_LOCAL_INFO_YM)-1 359c4762a1bSJed Brown mx = info(DMDA_LOCAL_INFO_MX) 360c4762a1bSJed Brown my = info(DMDA_LOCAL_INFO_MY) 361c4762a1bSJed Brown 362c4762a1bSJed Brown one = 1.0 363c4762a1bSJed Brown two = 2.0 364bb060b03SBarry Smith hx = one/(real(mx)-1) 365bb060b03SBarry Smith hy = one/(real(my)-1) 366c4762a1bSJed Brown sc = hx*hy*lambda 367c4762a1bSJed Brown hxdhy = hx/hy 368c4762a1bSJed Brown hydhx = hy/hx 369c4762a1bSJed Brown 370c4762a1bSJed Brown! Compute function over the locally owned part of the grid 371c4762a1bSJed Brown 372c4762a1bSJed Brown do 20 j=ys,ye 373c4762a1bSJed Brown do 10 i=xs,xe 374c4762a1bSJed Brown if (i .eq. 1 .or. j .eq. 1 .or. i .eq. mx .or. j .eq. my) then 375c4762a1bSJed Brown f(i,j) = x(i,j) 376c4762a1bSJed Brown else 377c4762a1bSJed Brown u = x(i,j) 378c4762a1bSJed Brown uxx = hydhx * (two*u - x(i-1,j) - x(i+1,j)) 379c4762a1bSJed Brown uyy = hxdhy * (two*u - x(i,j-1) - x(i,j+1)) 380c4762a1bSJed Brown f(i,j) = uxx + uyy - sc*exp(u) 381c4762a1bSJed Brown endif 382c4762a1bSJed Brown 10 continue 383c4762a1bSJed Brown 20 continue 384c4762a1bSJed Brown 38565afa91aSSatish Balay call PetscLogFlops(11.0d0*ym*xm,ierr) 38665afa91aSSatish Balay CHKERRQ(ierr) 387c4762a1bSJed Brown 388c4762a1bSJed Brown end 389c4762a1bSJed Brown 390c4762a1bSJed Brown! --------------------------------------------------------------------- 391c4762a1bSJed Brown! 392c4762a1bSJed Brown! FormJacobianLocal - Computes Jacobian matrix, called by 393c4762a1bSJed Brown! the higher level routine FormJacobian(). 394c4762a1bSJed Brown! 395c4762a1bSJed Brown! Input Parameters: 396c4762a1bSJed Brown! x - local vector data 397c4762a1bSJed Brown! 398c4762a1bSJed Brown! Output Parameters: 399c4762a1bSJed Brown! jac - Jacobian matrix 400c4762a1bSJed Brown! jac_prec - optionally different preconditioning matrix (not used here) 401c4762a1bSJed Brown! ierr - error code 402c4762a1bSJed Brown! 403c4762a1bSJed Brown! Notes: 404c4762a1bSJed Brown! This routine uses standard Fortran-style computations over a 2-dim array. 405c4762a1bSJed Brown! 406c4762a1bSJed Brown! Notes: 407c4762a1bSJed Brown! Due to grid point reordering with DMDAs, we must always work 408c4762a1bSJed Brown! with the local grid points, and then transform them to the new 409c4762a1bSJed Brown! global numbering with the "ltog" mapping 410c4762a1bSJed Brown! We cannot work directly with the global numbers for the original 411c4762a1bSJed Brown! uniprocessor grid! 412c4762a1bSJed Brown! 413c4762a1bSJed Brown! Two methods are available for imposing this transformation 414c4762a1bSJed Brown! when setting matrix entries: 415c4762a1bSJed Brown! (A) MatSetValuesLocal(), using the local ordering (including 416c4762a1bSJed Brown! ghost points!) 417c4762a1bSJed Brown! by calling MatSetValuesLocal() 418c4762a1bSJed Brown! (B) MatSetValues(), using the global ordering 419c4762a1bSJed Brown! - Use DMDAGetGlobalIndices() to extract the local-to-global map 420c4762a1bSJed Brown! - Then apply this map explicitly yourself 421c4762a1bSJed Brown! - Set matrix entries using the global ordering by calling 422c4762a1bSJed Brown! MatSetValues() 423c4762a1bSJed Brown! Option (A) seems cleaner/easier in many cases, and is the procedure 424c4762a1bSJed Brown! used in this example. 425c4762a1bSJed Brown! 426c4762a1bSJed Brown subroutine FormJacobianLocal(info,x,A,jac,da,ierr) 427dfbbaf82SBarry Smith use ex5fmodule 428c4762a1bSJed Brown implicit none 429c4762a1bSJed Brown 430c4762a1bSJed Brown DM da 431c4762a1bSJed Brown 432c4762a1bSJed Brown! Input/output variables: 433c4762a1bSJed Brown PetscScalar x(gxs:gxe,gys:gye) 434c4762a1bSJed Brown Mat A,jac 435c4762a1bSJed Brown PetscErrorCode ierr 436c4762a1bSJed Brown DMDALocalInfo info(DMDA_LOCAL_INFO_SIZE) 437c4762a1bSJed Brown 438c4762a1bSJed Brown! Local variables: 439c4762a1bSJed Brown PetscInt row,col(5),i,j,i1,i5 440c4762a1bSJed Brown PetscScalar two,one,hx,hy,v(5) 441c4762a1bSJed Brown PetscScalar hxdhy,hydhx,sc 442c4762a1bSJed Brown 443c4762a1bSJed Brown! Set parameters 444c4762a1bSJed Brown 445c4762a1bSJed Brown i1 = 1 446c4762a1bSJed Brown i5 = 5 447c4762a1bSJed Brown one = 1.0 448c4762a1bSJed Brown two = 2.0 449bb060b03SBarry Smith hx = one/(real(mx)-1) 450bb060b03SBarry Smith hy = one/(real(my)-1) 451c4762a1bSJed Brown sc = hx*hy 452c4762a1bSJed Brown hxdhy = hx/hy 453c4762a1bSJed Brown hydhx = hy/hx 454*e0b71005SStefano Zampini! -Wmaybe-uninitialized 455*e0b71005SStefano Zampini v = 0.0 456*e0b71005SStefano Zampini col = 0 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 4795d83a8b1SBarry Smith 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 4935d83a8b1SBarry Smith 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 end 509c4762a1bSJed Brown 510c4762a1bSJed Brown! 511c4762a1bSJed Brown! Simple convergence test based on the infinity norm of the residual being small 512c4762a1bSJed Brown! 513c4762a1bSJed Brown subroutine MySNESConverged(snes,it,xnorm,snorm,fnorm,reason,dummy,ierr) 514dfbbaf82SBarry Smith use ex5fmodule 515c4762a1bSJed Brown implicit none 516c4762a1bSJed Brown 517c4762a1bSJed Brown SNES snes 518c4762a1bSJed Brown PetscInt it,dummy 519c4762a1bSJed Brown PetscReal xnorm,snorm,fnorm,nrm 520c4762a1bSJed Brown SNESConvergedReason reason 521c4762a1bSJed Brown Vec f 522c4762a1bSJed Brown PetscErrorCode ierr 523c4762a1bSJed Brown 52465afa91aSSatish Balay call SNESGetFunction(snes,f,PETSC_NULL_FUNCTION,dummy,ierr) 52565afa91aSSatish Balay CHKERRQ(ierr) 52665afa91aSSatish Balay call VecNorm(f,NORM_INFINITY,nrm,ierr) 52765afa91aSSatish Balay CHKERRQ(ierr) 528c4762a1bSJed Brown if (nrm .le. 1.e-5) reason = SNES_CONVERGED_FNORM_ABS 529c4762a1bSJed Brown 530c4762a1bSJed Brown end 531c4762a1bSJed Brown 532c4762a1bSJed Brown!/*TEST 533c4762a1bSJed Brown! 534c4762a1bSJed Brown! build: 535c4762a1bSJed Brown! requires: !complex !single 536c4762a1bSJed Brown! 537c4762a1bSJed Brown! test: 538c4762a1bSJed Brown! nsize: 4 5398f8b3c79SStefano Zampini! args: -snes_mf -pc_type none -da_processors_x 4 -da_processors_y 1 -snes_monitor_short \ 5408f8b3c79SStefano Zampini! -ksp_gmres_cgs_refinement_type refine_always 541c4762a1bSJed Brown! 542c4762a1bSJed Brown! test: 543c4762a1bSJed Brown! suffix: 2 544c4762a1bSJed Brown! nsize: 4 545c4762a1bSJed Brown! args: -da_processors_x 2 -da_processors_y 2 -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always 546c4762a1bSJed Brown! 547c4762a1bSJed Brown! test: 548c4762a1bSJed Brown! suffix: 3 549c4762a1bSJed Brown! nsize: 3 550c4762a1bSJed Brown! args: -snes_fd -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always 551c4762a1bSJed Brown! 552c4762a1bSJed Brown! test: 553c4762a1bSJed Brown! suffix: 6 554c4762a1bSJed Brown! nsize: 1 555c4762a1bSJed Brown! args: -snes_monitor_short -my_snes_convergence 556c4762a1bSJed Brown! 557c4762a1bSJed Brown!TEST*/ 558