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