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