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