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