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