1c4762a1bSJed Brown! 2c4762a1bSJed Brown! Description: 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 <parameter>, where <parameter> indicates the nonlinearity of the problem 7c4762a1bSJed Brown! problem SFI: <parameter> = Bratu parameter (0 <= par <= 6.81) 8c4762a1bSJed Brown! 9c4762a1bSJed Brown 10c4762a1bSJed Brown! 11c4762a1bSJed Brown! -------------------------------------------------------------------------- 12c4762a1bSJed Brown! 13c4762a1bSJed Brown! Solid Fuel Ignition (SFI) problem. This problem is modeled by 14c4762a1bSJed Brown! the partial differential equation 15c4762a1bSJed Brown! 16c4762a1bSJed Brown! -Laplacian u - lambda*exp(u) = 0, 0 < x,y < 1, 17c4762a1bSJed Brown! 18c4762a1bSJed Brown! with boundary conditions 19c4762a1bSJed Brown! 20c4762a1bSJed Brown! u = 0 for x = 0, x = 1, y = 0, y = 1. 21c4762a1bSJed Brown! 22c4762a1bSJed Brown! A finite difference approximation with the usual 5-point stencil 23c4762a1bSJed Brown! is used to discretize the boundary value problem to obtain a nonlinear 24c4762a1bSJed Brown! system of equations. 25c4762a1bSJed Brown! 26c4762a1bSJed Brown! The uniprocessor version of this code is snes/tutorials/ex4f.F 27c4762a1bSJed Brown! 28c4762a1bSJed Brown! -------------------------------------------------------------------------- 29c4762a1bSJed Brown! The following define must be used before including any PETSc include files 30c4762a1bSJed Brown! into a module or interface. This is because they can't handle declarations 31c4762a1bSJed Brown! in them 32c4762a1bSJed Brown! 33c4762a1bSJed Brown 34c4762a1bSJed Brown module f90module 35c4762a1bSJed Brown use petscsys 36c4762a1bSJed Brown use petscis 37c4762a1bSJed Brown use petscvec 38c4762a1bSJed Brown use petscdm 39c4762a1bSJed Brown use petscdmda 40c4762a1bSJed Brown use petscmat 41c4762a1bSJed Brown use petscpc 42c4762a1bSJed Brown use petscksp 43c4762a1bSJed Brown use petscsnes 44c4762a1bSJed Brown#include <petsc/finclude/petscsnes.h> 45c4762a1bSJed Brown type userctx 46c4762a1bSJed Brown PetscInt xs,xe,xm,gxs,gxe,gxm 47c4762a1bSJed Brown PetscInt ys,ye,ym,gys,gye,gym 48c4762a1bSJed Brown PetscInt mx,my 49c4762a1bSJed Brown PetscMPIInt rank 50c4762a1bSJed Brown PetscReal lambda 51c4762a1bSJed Brown end type userctx 52c4762a1bSJed Brown 53c4762a1bSJed Brown contains 54c4762a1bSJed Brown! --------------------------------------------------------------------- 55c4762a1bSJed Brown! 56c4762a1bSJed Brown! FormFunction - Evaluates nonlinear function, F(x). 57c4762a1bSJed Brown! 58c4762a1bSJed Brown! Input Parameters: 59c4762a1bSJed Brown! snes - the SNES context 60c4762a1bSJed Brown! X - input vector 61c4762a1bSJed Brown! dummy - optional user-defined context, as set by SNESSetFunction() 62c4762a1bSJed Brown! (not used here) 63c4762a1bSJed Brown! 64c4762a1bSJed Brown! Output Parameter: 65c4762a1bSJed Brown! F - function vector 66c4762a1bSJed Brown! 67c4762a1bSJed Brown! Notes: 68c4762a1bSJed Brown! This routine serves as a wrapper for the lower-level routine 69c4762a1bSJed Brown! "FormFunctionLocal", where the actual computations are 70c4762a1bSJed Brown! done using the standard Fortran style of treating the local 71c4762a1bSJed Brown! vector data as a multidimensional array over the local mesh. 72c4762a1bSJed Brown! This routine merely handles ghost point scatters and accesses 73c4762a1bSJed Brown! the local vector data via VecGetArrayF90() and VecRestoreArrayF90(). 74c4762a1bSJed Brown! 75c4762a1bSJed Brown subroutine FormFunction(snes,X,F,user,ierr) 76c4762a1bSJed Brown implicit none 77c4762a1bSJed Brown 78c4762a1bSJed Brown! Input/output variables: 79c4762a1bSJed Brown SNES snes 80c4762a1bSJed Brown Vec X,F 81c4762a1bSJed Brown PetscErrorCode ierr 82c4762a1bSJed Brown type (userctx) user 83c4762a1bSJed Brown DM da 84c4762a1bSJed Brown 85c4762a1bSJed Brown! Declarations for use with local arrays: 86c4762a1bSJed Brown PetscScalar,pointer :: lx_v(:),lf_v(:) 87c4762a1bSJed Brown Vec localX 88c4762a1bSJed Brown 89c4762a1bSJed Brown! Scatter ghost points to local vector, using the 2-step process 90c4762a1bSJed Brown! DMGlobalToLocalBegin(), DMGlobalToLocalEnd(). 91c4762a1bSJed Brown! By placing code between these two statements, computations can 92c4762a1bSJed Brown! be done while messages are in transition. 93*d8606c27SBarry Smith PetscCall(SNESGetDM(snes,da,ierr)) 94*d8606c27SBarry Smith PetscCall(DMGetLocalVector(da,localX,ierr)) 95*d8606c27SBarry Smith PetscCall(DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX,ierr)) 96*d8606c27SBarry Smith PetscCall(DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX,ierr)) 97c4762a1bSJed Brown 98c4762a1bSJed Brown! Get a pointer to vector data. 99c4762a1bSJed Brown! - For default PETSc vectors, VecGetArray90() returns a pointer to 100c4762a1bSJed Brown! the data array. Otherwise, the routine is implementation dependent. 101c4762a1bSJed Brown! - You MUST call VecRestoreArrayF90() when you no longer need access to 102c4762a1bSJed Brown! the array. 103c4762a1bSJed Brown! - Note that the interface to VecGetArrayF90() differs from VecGetArray(), 104c4762a1bSJed Brown! and is useable from Fortran-90 Only. 105c4762a1bSJed Brown 106*d8606c27SBarry Smith PetscCall(VecGetArrayF90(localX,lx_v,ierr)) 107*d8606c27SBarry Smith PetscCall(VecGetArrayF90(F,lf_v,ierr)) 108c4762a1bSJed Brown 109c4762a1bSJed Brown! Compute function over the locally owned part of the grid 110*d8606c27SBarry Smith PetscCall(FormFunctionLocal(lx_v,lf_v,user,ierr)) 111c4762a1bSJed Brown 112c4762a1bSJed Brown! Restore vectors 113*d8606c27SBarry Smith PetscCall(VecRestoreArrayF90(localX,lx_v,ierr)) 114*d8606c27SBarry Smith PetscCall(VecRestoreArrayF90(F,lf_v,ierr)) 115c4762a1bSJed Brown 116c4762a1bSJed Brown! Insert values into global vector 117c4762a1bSJed Brown 118*d8606c27SBarry Smith PetscCall(DMRestoreLocalVector(da,localX,ierr)) 119*d8606c27SBarry Smith PetscCall(PetscLogFlops(11.0d0*user%ym*user%xm,ierr)) 120c4762a1bSJed Brown 121*d8606c27SBarry Smith! PetscCallA(VecView(X,PETSC_VIEWER_STDOUT_WORLD,ierr)) 122*d8606c27SBarry Smith! PetscCallA(VecView(F,PETSC_VIEWER_STDOUT_WORLD,ierr)) 123c4762a1bSJed Brown return 124c4762a1bSJed Brown end subroutine formfunction 125c4762a1bSJed Brown end module f90module 126c4762a1bSJed Brown 127c4762a1bSJed Brown module f90moduleinterfaces 128c4762a1bSJed Brown use f90module 129c4762a1bSJed Brown 130c4762a1bSJed Brown Interface SNESSetApplicationContext 131c4762a1bSJed Brown Subroutine SNESSetApplicationContext(snes,ctx,ierr) 132c4762a1bSJed Brown use f90module 133c4762a1bSJed Brown SNES snes 134c4762a1bSJed Brown type(userctx) ctx 135c4762a1bSJed Brown PetscErrorCode ierr 136c4762a1bSJed Brown End Subroutine 137c4762a1bSJed Brown End Interface SNESSetApplicationContext 138c4762a1bSJed Brown 139c4762a1bSJed Brown Interface SNESGetApplicationContext 140c4762a1bSJed Brown Subroutine SNESGetApplicationContext(snes,ctx,ierr) 141c4762a1bSJed Brown use f90module 142c4762a1bSJed Brown SNES snes 143c4762a1bSJed Brown type(userctx), pointer :: ctx 144c4762a1bSJed Brown PetscErrorCode ierr 145c4762a1bSJed Brown End Subroutine 146c4762a1bSJed Brown End Interface SNESGetApplicationContext 147c4762a1bSJed Brown end module f90moduleinterfaces 148c4762a1bSJed Brown 149c4762a1bSJed Brown program main 150c4762a1bSJed Brown use f90module 151c4762a1bSJed Brown use f90moduleinterfaces 152c4762a1bSJed Brown implicit none 153c4762a1bSJed Brown! 154c4762a1bSJed Brown 155c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 156c4762a1bSJed Brown! Variable declarations 157c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 158c4762a1bSJed Brown! 159c4762a1bSJed Brown! Variables: 160c4762a1bSJed Brown! snes - nonlinear solver 161c4762a1bSJed Brown! x, r - solution, residual vectors 162c4762a1bSJed Brown! J - Jacobian matrix 163c4762a1bSJed Brown! its - iterations for convergence 164c4762a1bSJed Brown! Nx, Ny - number of preocessors in x- and y- directions 165c4762a1bSJed Brown! matrix_free - flag - 1 indicates matrix-free version 166c4762a1bSJed Brown! 167c4762a1bSJed Brown SNES snes 168c4762a1bSJed Brown Vec x,r 169c4762a1bSJed Brown Mat J 170c4762a1bSJed Brown PetscErrorCode ierr 171c4762a1bSJed Brown PetscInt its 172c4762a1bSJed Brown PetscBool flg,matrix_free 173c4762a1bSJed Brown PetscInt ione,nfour 174c4762a1bSJed Brown PetscReal lambda_max,lambda_min 175c4762a1bSJed Brown type (userctx) user 176c4762a1bSJed Brown DM da 177c4762a1bSJed Brown 178c4762a1bSJed Brown! Note: Any user-defined Fortran routines (such as FormJacobian) 179c4762a1bSJed Brown! MUST be declared as external. 180c4762a1bSJed Brown external FormInitialGuess,FormJacobian 181c4762a1bSJed Brown 182c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 183c4762a1bSJed Brown! Initialize program 184c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 185*d8606c27SBarry Smith PetscCallA(PetscInitialize(ierr)) 186*d8606c27SBarry Smith PetscCallMPIA(MPI_Comm_rank(PETSC_COMM_WORLD,user%rank,ierr)) 187c4762a1bSJed Brown 188c4762a1bSJed Brown! Initialize problem parameters 189c4762a1bSJed Brown lambda_max = 6.81 190c4762a1bSJed Brown lambda_min = 0.0 191c4762a1bSJed Brown user%lambda = 6.0 192c4762a1bSJed Brown ione = 1 193c4762a1bSJed Brown nfour = 4 194*d8606c27SBarry Smith PetscCallA(PetscOptionsGetReal(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-par',user%lambda,flg,ierr)) 195*d8606c27SBarry Smith if (user%lambda .ge. lambda_max .or. user%lambda .le. lambda_min) then 196*d8606c27SBarry Smith SETERRA(PETSC_COMM_SELF,PETSC_ERR_USER,'Lambda provided with -par is out of range ') 197*d8606c27SBarry Smith endif 198c4762a1bSJed Brown 199c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 200c4762a1bSJed Brown! Create nonlinear solver context 201c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 202*d8606c27SBarry Smith PetscCallA(SNESCreate(PETSC_COMM_WORLD,snes,ierr)) 203c4762a1bSJed Brown 204c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 205c4762a1bSJed Brown! Create vector data structures; set function evaluation routine 206c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 207c4762a1bSJed Brown 208c4762a1bSJed Brown! Create distributed array (DMDA) to manage parallel grid and vectors 209c4762a1bSJed Brown 210c4762a1bSJed Brown! This really needs only the star-type stencil, but we use the box 211c4762a1bSJed Brown! stencil temporarily. 212*d8606c27SBarry Smith PetscCallA(DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,nfour,nfour,PETSC_DECIDE,PETSC_DECIDE,ione,ione,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,da,ierr)) 213*d8606c27SBarry Smith PetscCallA(DMSetFromOptions(da,ierr)) 214*d8606c27SBarry Smith PetscCallA(DMSetUp(da,ierr)) 215c4762a1bSJed Brown 216*d8606c27SBarry Smith PetscCallA(DMDAGetInfo(da,PETSC_NULL_INTEGER,user%mx,user%my,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,ierr)) 217c4762a1bSJed Brown 218c4762a1bSJed Brown! 219c4762a1bSJed Brown! Visualize the distribution of the array across the processors 220c4762a1bSJed Brown! 221*d8606c27SBarry Smith! PetscCallA(DMView(da,PETSC_VIEWER_DRAW_WORLD,ierr)) 222c4762a1bSJed Brown 223c4762a1bSJed Brown! Extract global and local vectors from DMDA; then duplicate for remaining 224c4762a1bSJed Brown! vectors that are the same types 225*d8606c27SBarry Smith PetscCallA(DMCreateGlobalVector(da,x,ierr)) 226*d8606c27SBarry Smith PetscCallA(VecDuplicate(x,r,ierr)) 227c4762a1bSJed Brown 228c4762a1bSJed Brown! Get local grid boundaries (for 2-dimensional DMDA) 229*d8606c27SBarry Smith PetscCallA(DMDAGetCorners(da,user%xs,user%ys,PETSC_NULL_INTEGER,user%xm,user%ym,PETSC_NULL_INTEGER,ierr)) 230*d8606c27SBarry Smith PetscCallA(DMDAGetGhostCorners(da,user%gxs,user%gys,PETSC_NULL_INTEGER,user%gxm,user%gym,PETSC_NULL_INTEGER,ierr)) 231c4762a1bSJed Brown 232c4762a1bSJed Brown! Here we shift the starting indices up by one so that we can easily 233c4762a1bSJed Brown! use the Fortran convention of 1-based indices (rather 0-based indices). 234c4762a1bSJed Brown user%xs = user%xs+1 235c4762a1bSJed Brown user%ys = user%ys+1 236c4762a1bSJed Brown user%gxs = user%gxs+1 237c4762a1bSJed Brown user%gys = user%gys+1 238c4762a1bSJed Brown 239c4762a1bSJed Brown user%ye = user%ys+user%ym-1 240c4762a1bSJed Brown user%xe = user%xs+user%xm-1 241c4762a1bSJed Brown user%gye = user%gys+user%gym-1 242c4762a1bSJed Brown user%gxe = user%gxs+user%gxm-1 243c4762a1bSJed Brown 244*d8606c27SBarry Smith PetscCallA(SNESSetApplicationContext(snes,user,ierr)) 245c4762a1bSJed Brown 246c4762a1bSJed Brown! Set function evaluation routine and vector 247*d8606c27SBarry Smith PetscCallA(SNESSetFunction(snes,r,FormFunction,user,ierr)) 248c4762a1bSJed Brown 249c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 250c4762a1bSJed Brown! Create matrix data structure; set Jacobian evaluation routine 251c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 252c4762a1bSJed Brown 253c4762a1bSJed Brown! Set Jacobian matrix data structure and default Jacobian evaluation 254c4762a1bSJed Brown! routine. User can override with: 255c4762a1bSJed Brown! -snes_fd : default finite differencing approximation of Jacobian 256c4762a1bSJed Brown! -snes_mf : matrix-free Newton-Krylov method with no preconditioning 257c4762a1bSJed Brown! (unless user explicitly sets preconditioner) 258c4762a1bSJed Brown! -snes_mf_operator : form preconditioning matrix as set by the user, 259c4762a1bSJed Brown! but use matrix-free approx for Jacobian-vector 260c4762a1bSJed Brown! products within Newton-Krylov method 261c4762a1bSJed Brown! 262c4762a1bSJed Brown! Note: For the parallel case, vectors and matrices MUST be partitioned 263c4762a1bSJed Brown! accordingly. When using distributed arrays (DMDAs) to create vectors, 264c4762a1bSJed Brown! the DMDAs determine the problem partitioning. We must explicitly 265c4762a1bSJed Brown! specify the local matrix dimensions upon its creation for compatibility 266c4762a1bSJed Brown! with the vector distribution. Thus, the generic MatCreate() routine 267c4762a1bSJed Brown! is NOT sufficient when working with distributed arrays. 268c4762a1bSJed Brown! 269c4762a1bSJed Brown! Note: Here we only approximately preallocate storage space for the 270c4762a1bSJed Brown! Jacobian. See the users manual for a discussion of better techniques 271c4762a1bSJed Brown! for preallocating matrix memory. 272c4762a1bSJed Brown 273*d8606c27SBarry Smith PetscCallA(PetscOptionsHasName(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-snes_mf',matrix_free,ierr)) 274c4762a1bSJed Brown if (.not. matrix_free) then 275*d8606c27SBarry Smith PetscCallA(DMSetMatType(da,MATAIJ,ierr)) 276*d8606c27SBarry Smith PetscCallA(DMCreateMatrix(da,J,ierr)) 277*d8606c27SBarry Smith PetscCallA(SNESSetJacobian(snes,J,J,FormJacobian,user,ierr)) 278c4762a1bSJed Brown endif 279c4762a1bSJed Brown 280c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 281c4762a1bSJed Brown! Customize nonlinear solver; set runtime options 282c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 283c4762a1bSJed Brown! Set runtime options (e.g., -snes_monitor -snes_rtol <rtol> -ksp_type <type>) 284*d8606c27SBarry Smith PetscCallA(SNESSetDM(snes,da,ierr)) 285*d8606c27SBarry Smith PetscCallA(SNESSetFromOptions(snes,ierr)) 286c4762a1bSJed Brown 287c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 288c4762a1bSJed Brown! Evaluate initial guess; then solve nonlinear system. 289c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 290c4762a1bSJed Brown! Note: The user should initialize the vector, x, with the initial guess 291c4762a1bSJed Brown! for the nonlinear solver prior to calling SNESSolve(). In particular, 292c4762a1bSJed Brown! to employ an initial guess of zero, the user should explicitly set 293c4762a1bSJed Brown! this vector to zero by calling VecSet(). 294c4762a1bSJed Brown 295*d8606c27SBarry Smith PetscCallA(FormInitialGuess(snes,x,ierr)) 296*d8606c27SBarry Smith PetscCallA(SNESSolve(snes,PETSC_NULL_VEC,x,ierr)) 297*d8606c27SBarry Smith PetscCallA(SNESGetIterationNumber(snes,its,ierr)) 298c4762a1bSJed Brown if (user%rank .eq. 0) then 299c4762a1bSJed Brown write(6,100) its 300c4762a1bSJed Brown endif 301c4762a1bSJed Brown 100 format('Number of SNES iterations = ',i5) 302c4762a1bSJed Brown 303c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 304c4762a1bSJed Brown! Free work space. All PETSc objects should be destroyed when they 305c4762a1bSJed Brown! are no longer needed. 306c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 307*d8606c27SBarry Smith if (.not. matrix_free) PetscCallA(MatDestroy(J,ierr)) 308*d8606c27SBarry Smith PetscCallA(VecDestroy(x,ierr)) 309*d8606c27SBarry Smith PetscCallA(VecDestroy(r,ierr)) 310*d8606c27SBarry Smith PetscCallA(SNESDestroy(snes,ierr)) 311*d8606c27SBarry Smith PetscCallA(DMDestroy(da,ierr)) 312c4762a1bSJed Brown 313*d8606c27SBarry Smith PetscCallA(PetscFinalize(ierr)) 314c4762a1bSJed Brown end 315c4762a1bSJed Brown 316c4762a1bSJed Brown! --------------------------------------------------------------------- 317c4762a1bSJed Brown! 318c4762a1bSJed Brown! FormInitialGuess - Forms initial approximation. 319c4762a1bSJed Brown! 320c4762a1bSJed Brown! Input Parameters: 321c4762a1bSJed Brown! X - vector 322c4762a1bSJed Brown! 323c4762a1bSJed Brown! Output Parameter: 324c4762a1bSJed Brown! X - vector 325c4762a1bSJed Brown! 326c4762a1bSJed Brown! Notes: 327c4762a1bSJed Brown! This routine serves as a wrapper for the lower-level routine 328c4762a1bSJed Brown! "InitialGuessLocal", where the actual computations are 329c4762a1bSJed Brown! done using the standard Fortran style of treating the local 330c4762a1bSJed Brown! vector data as a multidimensional array over the local mesh. 331c4762a1bSJed Brown! This routine merely handles ghost point scatters and accesses 332c4762a1bSJed Brown! the local vector data via VecGetArrayF90() and VecRestoreArrayF90(). 333c4762a1bSJed Brown! 334c4762a1bSJed Brown subroutine FormInitialGuess(snes,X,ierr) 335c4762a1bSJed Brown use f90module 336c4762a1bSJed Brown use f90moduleinterfaces 337c4762a1bSJed Brown implicit none 338c4762a1bSJed Brown 339c4762a1bSJed Brown! Input/output variables: 340c4762a1bSJed Brown SNES snes 341c4762a1bSJed Brown type(userctx), pointer:: puser 342c4762a1bSJed Brown Vec X 343c4762a1bSJed Brown PetscErrorCode ierr 344c4762a1bSJed Brown DM da 345c4762a1bSJed Brown 346c4762a1bSJed Brown! Declarations for use with local arrays: 347c4762a1bSJed Brown PetscScalar,pointer :: lx_v(:) 348c4762a1bSJed Brown 349c4762a1bSJed Brown ierr = 0 350*d8606c27SBarry Smith PetscCallA(SNESGetDM(snes,da,ierr)) 351*d8606c27SBarry Smith PetscCallA(SNESGetApplicationContext(snes,puser,ierr)) 352c4762a1bSJed Brown! Get a pointer to vector data. 353c4762a1bSJed Brown! - For default PETSc vectors, VecGetArray90() returns a pointer to 354c4762a1bSJed Brown! the data array. Otherwise, the routine is implementation dependent. 355c4762a1bSJed Brown! - You MUST call VecRestoreArrayF90() when you no longer need access to 356c4762a1bSJed Brown! the array. 357c4762a1bSJed Brown! - Note that the interface to VecGetArrayF90() differs from VecGetArray(), 358c4762a1bSJed Brown! and is useable from Fortran-90 Only. 359c4762a1bSJed Brown 360*d8606c27SBarry Smith PetscCallA(VecGetArrayF90(X,lx_v,ierr)) 361c4762a1bSJed Brown 362c4762a1bSJed Brown! Compute initial guess over the locally owned part of the grid 363*d8606c27SBarry Smith PetscCallA(InitialGuessLocal(puser,lx_v,ierr)) 364c4762a1bSJed Brown 365c4762a1bSJed Brown! Restore vector 366*d8606c27SBarry Smith PetscCallA(VecRestoreArrayF90(X,lx_v,ierr)) 367c4762a1bSJed Brown 368c4762a1bSJed Brown! Insert values into global vector 369c4762a1bSJed Brown 370c4762a1bSJed Brown return 371c4762a1bSJed Brown end 372c4762a1bSJed Brown 373c4762a1bSJed Brown! --------------------------------------------------------------------- 374c4762a1bSJed Brown! 375c4762a1bSJed Brown! InitialGuessLocal - Computes initial approximation, called by 376c4762a1bSJed Brown! the higher level routine FormInitialGuess(). 377c4762a1bSJed Brown! 378c4762a1bSJed Brown! Input Parameter: 379c4762a1bSJed Brown! x - local vector data 380c4762a1bSJed Brown! 381c4762a1bSJed Brown! Output Parameters: 382c4762a1bSJed Brown! x - local vector data 383c4762a1bSJed Brown! ierr - error code 384c4762a1bSJed Brown! 385c4762a1bSJed Brown! Notes: 386c4762a1bSJed Brown! This routine uses standard Fortran-style computations over a 2-dim array. 387c4762a1bSJed Brown! 388c4762a1bSJed Brown subroutine InitialGuessLocal(user,x,ierr) 389c4762a1bSJed Brown use f90module 390c4762a1bSJed Brown implicit none 391c4762a1bSJed Brown 392c4762a1bSJed Brown! Input/output variables: 393c4762a1bSJed Brown type (userctx) user 394c4762a1bSJed Brown PetscScalar x(user%xs:user%xe,user%ys:user%ye) 395c4762a1bSJed Brown PetscErrorCode ierr 396c4762a1bSJed Brown 397c4762a1bSJed Brown! Local variables: 398c4762a1bSJed Brown PetscInt i,j 399c4762a1bSJed Brown PetscReal temp1,temp,hx,hy 400c4762a1bSJed Brown PetscReal one 401c4762a1bSJed Brown 402c4762a1bSJed Brown! Set parameters 403c4762a1bSJed Brown 404c4762a1bSJed Brown ierr = 0 405c4762a1bSJed Brown one = 1.0 406c4762a1bSJed Brown hx = one/(user%mx-1) 407c4762a1bSJed Brown hy = one/(user%my-1) 408c4762a1bSJed Brown temp1 = user%lambda/(user%lambda + one) 409c4762a1bSJed Brown 410c4762a1bSJed Brown do 20 j=user%ys,user%ye 411c4762a1bSJed Brown temp = min(j-1,user%my-j)*hy 412c4762a1bSJed Brown do 10 i=user%xs,user%xe 413c4762a1bSJed Brown if (i .eq. 1 .or. j .eq. 1 .or. i .eq. user%mx .or. j .eq. user%my) then 414c4762a1bSJed Brown x(i,j) = 0.0 415c4762a1bSJed Brown else 416c4762a1bSJed Brown x(i,j) = temp1 * sqrt(min(hx*min(i-1,user%mx-i),temp)) 417c4762a1bSJed Brown endif 418c4762a1bSJed Brown 10 continue 419c4762a1bSJed Brown 20 continue 420c4762a1bSJed Brown 421c4762a1bSJed Brown return 422c4762a1bSJed Brown end 423c4762a1bSJed Brown 424c4762a1bSJed Brown! --------------------------------------------------------------------- 425c4762a1bSJed Brown! 426c4762a1bSJed Brown! FormFunctionLocal - Computes nonlinear function, called by 427c4762a1bSJed Brown! the higher level routine FormFunction(). 428c4762a1bSJed Brown! 429c4762a1bSJed Brown! Input Parameter: 430c4762a1bSJed Brown! x - local vector data 431c4762a1bSJed Brown! 432c4762a1bSJed Brown! Output Parameters: 433c4762a1bSJed Brown! f - local vector data, f(x) 434c4762a1bSJed Brown! ierr - error code 435c4762a1bSJed Brown! 436c4762a1bSJed Brown! Notes: 437c4762a1bSJed Brown! This routine uses standard Fortran-style computations over a 2-dim array. 438c4762a1bSJed Brown! 439c4762a1bSJed Brown subroutine FormFunctionLocal(x,f,user,ierr) 440c4762a1bSJed Brown use f90module 441c4762a1bSJed Brown 442c4762a1bSJed Brown implicit none 443c4762a1bSJed Brown 444c4762a1bSJed Brown! Input/output variables: 445c4762a1bSJed Brown type (userctx) user 446c4762a1bSJed Brown PetscScalar x(user%gxs:user%gxe,user%gys:user%gye) 447c4762a1bSJed Brown PetscScalar f(user%xs:user%xe,user%ys:user%ye) 448c4762a1bSJed Brown PetscErrorCode ierr 449c4762a1bSJed Brown 450c4762a1bSJed Brown! Local variables: 451c4762a1bSJed Brown PetscScalar two,one,hx,hy,hxdhy,hydhx,sc 452c4762a1bSJed Brown PetscScalar u,uxx,uyy 453c4762a1bSJed Brown PetscInt i,j 454c4762a1bSJed Brown 455c4762a1bSJed Brown one = 1.0 456c4762a1bSJed Brown two = 2.0 457c4762a1bSJed Brown hx = one/(user%mx-1) 458c4762a1bSJed Brown hy = one/(user%my-1) 459c4762a1bSJed Brown sc = hx*hy*user%lambda 460c4762a1bSJed Brown hxdhy = hx/hy 461c4762a1bSJed Brown hydhx = hy/hx 462c4762a1bSJed Brown 463c4762a1bSJed Brown! Compute function over the locally owned part of the grid 464c4762a1bSJed Brown 465c4762a1bSJed Brown do 20 j=user%ys,user%ye 466c4762a1bSJed Brown do 10 i=user%xs,user%xe 467c4762a1bSJed Brown if (i .eq. 1 .or. j .eq. 1 .or. i .eq. user%mx .or. j .eq. user%my) then 468c4762a1bSJed Brown f(i,j) = x(i,j) 469c4762a1bSJed Brown else 470c4762a1bSJed Brown u = x(i,j) 471c4762a1bSJed Brown uxx = hydhx * (two*u - x(i-1,j) - x(i+1,j)) 472c4762a1bSJed Brown uyy = hxdhy * (two*u - x(i,j-1) - x(i,j+1)) 473c4762a1bSJed Brown f(i,j) = uxx + uyy - sc*exp(u) 474c4762a1bSJed Brown endif 475c4762a1bSJed Brown 10 continue 476c4762a1bSJed Brown 20 continue 477c4762a1bSJed Brown 478c4762a1bSJed Brown return 479c4762a1bSJed Brown end 480c4762a1bSJed Brown 481c4762a1bSJed Brown! --------------------------------------------------------------------- 482c4762a1bSJed Brown! 483c4762a1bSJed Brown! FormJacobian - Evaluates Jacobian matrix. 484c4762a1bSJed Brown! 485c4762a1bSJed Brown! Input Parameters: 486c4762a1bSJed Brown! snes - the SNES context 487c4762a1bSJed Brown! x - input vector 488c4762a1bSJed Brown! dummy - optional user-defined context, as set by SNESSetJacobian() 489c4762a1bSJed Brown! (not used here) 490c4762a1bSJed Brown! 491c4762a1bSJed Brown! Output Parameters: 492c4762a1bSJed Brown! jac - Jacobian matrix 493c4762a1bSJed Brown! jac_prec - optionally different preconditioning matrix (not used here) 494c4762a1bSJed Brown! flag - flag indicating matrix structure 495c4762a1bSJed Brown! 496c4762a1bSJed Brown! Notes: 497c4762a1bSJed Brown! This routine serves as a wrapper for the lower-level routine 498c4762a1bSJed Brown! "FormJacobianLocal", where the actual computations are 499c4762a1bSJed Brown! done using the standard Fortran style of treating the local 500c4762a1bSJed Brown! vector data as a multidimensional array over the local mesh. 501c4762a1bSJed Brown! This routine merely accesses the local vector data via 502c4762a1bSJed Brown! VecGetArrayF90() and VecRestoreArrayF90(). 503c4762a1bSJed Brown! 504c4762a1bSJed Brown! Notes: 505c4762a1bSJed Brown! Due to grid point reordering with DMDAs, we must always work 506c4762a1bSJed Brown! with the local grid points, and then transform them to the new 507c4762a1bSJed Brown! global numbering with the "ltog" mapping 508c4762a1bSJed Brown! We cannot work directly with the global numbers for the original 509c4762a1bSJed Brown! uniprocessor grid! 510c4762a1bSJed Brown! 511c4762a1bSJed Brown! Two methods are available for imposing this transformation 512c4762a1bSJed Brown! when setting matrix entries: 513c4762a1bSJed Brown! (A) MatSetValuesLocal(), using the local ordering (including 514c4762a1bSJed Brown! ghost points!) 515c4762a1bSJed Brown! - Set matrix entries using the local ordering 516c4762a1bSJed Brown! by calling MatSetValuesLocal() 517c4762a1bSJed Brown! (B) MatSetValues(), using the global ordering 518c4762a1bSJed Brown 519c4762a1bSJed Brown! - Set matrix entries using the global ordering by calling 520c4762a1bSJed Brown! MatSetValues() 521c4762a1bSJed Brown! Option (A) seems cleaner/easier in many cases, and is the procedure 522c4762a1bSJed Brown! used in this example. 523c4762a1bSJed Brown! 524c4762a1bSJed Brown subroutine FormJacobian(snes,X,jac,jac_prec,user,ierr) 525c4762a1bSJed Brown use f90module 526c4762a1bSJed Brown implicit none 527c4762a1bSJed Brown 528c4762a1bSJed Brown! Input/output variables: 529c4762a1bSJed Brown SNES snes 530c4762a1bSJed Brown Vec X 531c4762a1bSJed Brown Mat jac,jac_prec 532c4762a1bSJed Brown type(userctx) user 533c4762a1bSJed Brown PetscErrorCode ierr 534c4762a1bSJed Brown DM da 535c4762a1bSJed Brown 536c4762a1bSJed Brown! Declarations for use with local arrays: 537c4762a1bSJed Brown PetscScalar,pointer :: lx_v(:) 538c4762a1bSJed Brown Vec localX 539c4762a1bSJed Brown 540c4762a1bSJed Brown! Scatter ghost points to local vector, using the 2-step process 541c4762a1bSJed Brown! DMGlobalToLocalBegin(), DMGlobalToLocalEnd() 542c4762a1bSJed Brown! Computations can be done while messages are in transition, 543c4762a1bSJed Brown! by placing code between these two statements. 544c4762a1bSJed Brown 545*d8606c27SBarry Smith PetscCallA(SNESGetDM(snes,da,ierr)) 546*d8606c27SBarry Smith PetscCallA(DMGetLocalVector(da,localX,ierr)) 547*d8606c27SBarry Smith PetscCallA(DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX,ierr)) 548*d8606c27SBarry Smith PetscCallA(DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX,ierr)) 549c4762a1bSJed Brown 550c4762a1bSJed Brown! Get a pointer to vector data 551*d8606c27SBarry Smith PetscCallA(VecGetArrayF90(localX,lx_v,ierr)) 552c4762a1bSJed Brown 553c4762a1bSJed Brown! Compute entries for the locally owned part of the Jacobian preconditioner. 554*d8606c27SBarry Smith PetscCallA(FormJacobianLocal(lx_v,jac_prec,user,ierr)) 555c4762a1bSJed Brown 556c4762a1bSJed Brown! Assemble matrix, using the 2-step process: 557c4762a1bSJed Brown! MatAssemblyBegin(), MatAssemblyEnd() 558c4762a1bSJed Brown! Computations can be done while messages are in transition, 559c4762a1bSJed Brown! by placing code between these two statements. 560c4762a1bSJed Brown 561*d8606c27SBarry Smith PetscCallA(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY,ierr)) 562c4762a1bSJed Brown if (jac .ne. jac_prec) then 563*d8606c27SBarry Smith PetscCallA(MatAssemblyBegin(jac_prec,MAT_FINAL_ASSEMBLY,ierr)) 564c4762a1bSJed Brown endif 565*d8606c27SBarry Smith PetscCallA(VecRestoreArrayF90(localX,lx_v,ierr)) 566*d8606c27SBarry Smith PetscCallA(DMRestoreLocalVector(da,localX,ierr)) 567*d8606c27SBarry Smith PetscCallA(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY,ierr)) 568c4762a1bSJed Brown if (jac .ne. jac_prec) then 569*d8606c27SBarry Smith PetscCallA(MatAssemblyEnd(jac_prec,MAT_FINAL_ASSEMBLY,ierr)) 570c4762a1bSJed Brown endif 571c4762a1bSJed Brown 572c4762a1bSJed Brown! Tell the matrix we will never add a new nonzero location to the 573c4762a1bSJed Brown! matrix. If we do it will generate an error. 574c4762a1bSJed Brown 575*d8606c27SBarry Smith PetscCallA(MatSetOption(jac,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE,ierr)) 576c4762a1bSJed Brown 577c4762a1bSJed Brown return 578c4762a1bSJed Brown end 579c4762a1bSJed Brown 580c4762a1bSJed Brown! --------------------------------------------------------------------- 581c4762a1bSJed Brown! 582c4762a1bSJed Brown! FormJacobianLocal - Computes Jacobian preconditioner matrix, 583c4762a1bSJed Brown! called by the higher level routine FormJacobian(). 584c4762a1bSJed Brown! 585c4762a1bSJed Brown! Input Parameters: 586c4762a1bSJed Brown! x - local vector data 587c4762a1bSJed Brown! 588c4762a1bSJed Brown! Output Parameters: 589c4762a1bSJed Brown! jac_prec - Jacobian preconditioner matrix 590c4762a1bSJed Brown! ierr - error code 591c4762a1bSJed Brown! 592c4762a1bSJed Brown! Notes: 593c4762a1bSJed Brown! This routine uses standard Fortran-style computations over a 2-dim array. 594c4762a1bSJed Brown! 595c4762a1bSJed Brown! Notes: 596c4762a1bSJed Brown! Due to grid point reordering with DMDAs, we must always work 597c4762a1bSJed Brown! with the local grid points, and then transform them to the new 598c4762a1bSJed Brown! global numbering with the "ltog" mapping 599c4762a1bSJed Brown! We cannot work directly with the global numbers for the original 600c4762a1bSJed Brown! uniprocessor grid! 601c4762a1bSJed Brown! 602c4762a1bSJed Brown! Two methods are available for imposing this transformation 603c4762a1bSJed Brown! when setting matrix entries: 604c4762a1bSJed Brown! (A) MatSetValuesLocal(), using the local ordering (including 605c4762a1bSJed Brown! ghost points!) 606c4762a1bSJed Brown! - Set matrix entries using the local ordering 607c4762a1bSJed Brown! by calling MatSetValuesLocal() 608c4762a1bSJed Brown! (B) MatSetValues(), using the global ordering 609c4762a1bSJed Brown! - Then apply this map explicitly yourself 610c4762a1bSJed Brown! - Set matrix entries using the global ordering by calling 611c4762a1bSJed Brown! MatSetValues() 612c4762a1bSJed Brown! Option (A) seems cleaner/easier in many cases, and is the procedure 613c4762a1bSJed Brown! used in this example. 614c4762a1bSJed Brown! 615c4762a1bSJed Brown subroutine FormJacobianLocal(x,jac_prec,user,ierr) 616c4762a1bSJed Brown use f90module 617c4762a1bSJed Brown implicit none 618c4762a1bSJed Brown 619c4762a1bSJed Brown! Input/output variables: 620c4762a1bSJed Brown type (userctx) user 621c4762a1bSJed Brown PetscScalar x(user%gxs:user%gxe,user%gys:user%gye) 622c4762a1bSJed Brown Mat jac_prec 623c4762a1bSJed Brown PetscErrorCode ierr 624c4762a1bSJed Brown 625c4762a1bSJed Brown! Local variables: 626c4762a1bSJed Brown PetscInt row,col(5),i,j 627c4762a1bSJed Brown PetscInt ione,ifive 628c4762a1bSJed Brown PetscScalar two,one,hx,hy,hxdhy 629c4762a1bSJed Brown PetscScalar hydhx,sc,v(5) 630c4762a1bSJed Brown 631c4762a1bSJed Brown! Set parameters 632c4762a1bSJed Brown ione = 1 633c4762a1bSJed Brown ifive = 5 634c4762a1bSJed Brown one = 1.0 635c4762a1bSJed Brown two = 2.0 636c4762a1bSJed Brown hx = one/(user%mx-1) 637c4762a1bSJed Brown hy = one/(user%my-1) 638c4762a1bSJed Brown sc = hx*hy 639c4762a1bSJed Brown hxdhy = hx/hy 640c4762a1bSJed Brown hydhx = hy/hx 641c4762a1bSJed Brown 642c4762a1bSJed Brown! Compute entries for the locally owned part of the Jacobian. 643c4762a1bSJed Brown! - Currently, all PETSc parallel matrix formats are partitioned by 644c4762a1bSJed Brown! contiguous chunks of rows across the processors. 645c4762a1bSJed Brown! - Each processor needs to insert only elements that it owns 646c4762a1bSJed Brown! locally (but any non-local elements will be sent to the 647c4762a1bSJed Brown! appropriate processor during matrix assembly). 648c4762a1bSJed Brown! - Here, we set all entries for a particular row at once. 649c4762a1bSJed Brown! - We can set matrix entries either using either 650c4762a1bSJed Brown! MatSetValuesLocal() or MatSetValues(), as discussed above. 651c4762a1bSJed Brown! - Note that MatSetValues() uses 0-based row and column numbers 652c4762a1bSJed Brown! in Fortran as well as in C. 653c4762a1bSJed Brown 654c4762a1bSJed Brown do 20 j=user%ys,user%ye 655c4762a1bSJed Brown row = (j - user%gys)*user%gxm + user%xs - user%gxs - 1 656c4762a1bSJed Brown do 10 i=user%xs,user%xe 657c4762a1bSJed Brown row = row + 1 658c4762a1bSJed Brown! boundary points 659c4762a1bSJed Brown if (i .eq. 1 .or. j .eq. 1 .or. i .eq. user%mx .or. j .eq. user%my) then 660c4762a1bSJed Brown col(1) = row 661c4762a1bSJed Brown v(1) = one 662*d8606c27SBarry Smith PetscCallA(MatSetValuesLocal(jac_prec,ione,row,ione,col,v,INSERT_VALUES,ierr)) 663c4762a1bSJed Brown! interior grid points 664c4762a1bSJed Brown else 665c4762a1bSJed Brown v(1) = -hxdhy 666c4762a1bSJed Brown v(2) = -hydhx 667c4762a1bSJed Brown v(3) = two*(hydhx + hxdhy) - sc*user%lambda*exp(x(i,j)) 668c4762a1bSJed Brown v(4) = -hydhx 669c4762a1bSJed Brown v(5) = -hxdhy 670c4762a1bSJed Brown col(1) = row - user%gxm 671c4762a1bSJed Brown col(2) = row - 1 672c4762a1bSJed Brown col(3) = row 673c4762a1bSJed Brown col(4) = row + 1 674c4762a1bSJed Brown col(5) = row + user%gxm 675*d8606c27SBarry Smith PetscCallA(MatSetValuesLocal(jac_prec,ione,row,ifive,col,v,INSERT_VALUES,ierr)) 676c4762a1bSJed Brown endif 677c4762a1bSJed Brown 10 continue 678c4762a1bSJed Brown 20 continue 679c4762a1bSJed Brown 680c4762a1bSJed Brown return 681c4762a1bSJed Brown end 682c4762a1bSJed Brown 683c4762a1bSJed Brown! 684c4762a1bSJed Brown!/*TEST 685c4762a1bSJed Brown! 686c4762a1bSJed Brown! test: 687c4762a1bSJed Brown! nsize: 4 688c4762a1bSJed 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 689c4762a1bSJed Brown! requires: !single 690c4762a1bSJed Brown! 691c4762a1bSJed Brown! test: 692c4762a1bSJed Brown! suffix: 2 693c4762a1bSJed Brown! nsize: 4 694c4762a1bSJed Brown! args: -da_processors_x 2 -da_processors_y 2 -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always 695c4762a1bSJed Brown! requires: !single 696c4762a1bSJed Brown! 697c4762a1bSJed Brown! test: 698c4762a1bSJed Brown! suffix: 3 699c4762a1bSJed Brown! nsize: 3 700c4762a1bSJed Brown! args: -snes_fd -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always 701c4762a1bSJed Brown! requires: !single 702c4762a1bSJed Brown! 703c4762a1bSJed Brown! test: 704c4762a1bSJed Brown! suffix: 4 705c4762a1bSJed Brown! nsize: 3 706c4762a1bSJed Brown! args: -snes_mf_operator -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always 707c4762a1bSJed Brown! requires: !single 708c4762a1bSJed Brown! 709c4762a1bSJed Brown! test: 710c4762a1bSJed Brown! suffix: 5 711c4762a1bSJed Brown! requires: !single 712c4762a1bSJed Brown! 713c4762a1bSJed Brown!TEST*/ 714