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