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