1c4762a1bSJed Brown! 2c4762a1bSJed Brown! Description: This example solves a nonlinear system on 1 processor with SNES. 3c4762a1bSJed Brown! We solve the Bratu (SFI - solid fuel ignition) problem in a 2D rectangular 4c4762a1bSJed Brown! domain. The command line options include: 5c4762a1bSJed Brown! -par <parameter>, where <parameter> indicates the nonlinearity of the problem 6c4762a1bSJed Brown! problem SFI: <parameter> = Bratu parameter (0 <= par <= 6.81) 7c4762a1bSJed Brown! -mx <xg>, where <xg> = number of grid points in the x-direction 8c4762a1bSJed Brown! -my <yg>, where <yg> = number of grid points in the y-direction 9c4762a1bSJed Brown! 10c4762a1bSJed Brown 11c4762a1bSJed Brown! 12c4762a1bSJed Brown! -------------------------------------------------------------------------- 13c4762a1bSJed Brown! 14c4762a1bSJed Brown! Solid Fuel Ignition (SFI) problem. This problem is modeled by 15c4762a1bSJed Brown! the partial differential equation 16c4762a1bSJed Brown! 17c4762a1bSJed Brown! -Laplacian u - lambda*exp(u) = 0, 0 < x,y < 1, 18c4762a1bSJed Brown! 19c4762a1bSJed Brown! with boundary conditions 20c4762a1bSJed Brown! 21c4762a1bSJed Brown! u = 0 for x = 0, x = 1, y = 0, y = 1. 22c4762a1bSJed Brown! 23c4762a1bSJed Brown! A finite difference approximation with the usual 5-point stencil 24c4762a1bSJed Brown! is used to discretize the boundary value problem to obtain a nonlinear 25c4762a1bSJed Brown! system of equations. 26c4762a1bSJed Brown! 27c4762a1bSJed Brown! The parallel version of this code is snes/tutorials/ex5f.F 28c4762a1bSJed Brown! 29c4762a1bSJed Brown! -------------------------------------------------------------------------- 30c4762a1bSJed Brown subroutine postcheck(snes,x,y,w,changed_y,changed_w,ctx,ierr) 31c4762a1bSJed Brown#include <petsc/finclude/petscsnes.h> 32c4762a1bSJed Brown use petscsnes 33c4762a1bSJed Brown implicit none 34c4762a1bSJed Brown SNES snes 35c4762a1bSJed Brown PetscReal norm 36c4762a1bSJed Brown Vec tmp,x,y,w 37c4762a1bSJed Brown PetscBool changed_w,changed_y 38c4762a1bSJed Brown PetscErrorCode ierr 39c4762a1bSJed Brown PetscInt ctx 40c4762a1bSJed Brown PetscScalar mone 4124fb275aSStefano Zampini MPI_Comm comm 42c4762a1bSJed Brown 4324fb275aSStefano Zampini character(len=PETSC_MAX_PATH_LEN) :: outputString 4424fb275aSStefano Zampini 4524fb275aSStefano Zampini PetscCallA(PetscObjectGetComm(snes,comm,ierr)) 46d8606c27SBarry Smith PetscCallA(VecDuplicate(x,tmp,ierr)) 47c4762a1bSJed Brown mone = -1.0 48d8606c27SBarry Smith PetscCallA(VecWAXPY(tmp,mone,x,w,ierr)) 49d8606c27SBarry Smith PetscCallA(VecNorm(tmp,NORM_2,norm,ierr)) 50d8606c27SBarry Smith PetscCallA(VecDestroy(tmp,ierr)) 5124fb275aSStefano Zampini write(outputString,*) norm 5224fb275aSStefano Zampini PetscCallA(PetscPrintf(comm,'Norm of search step '//trim(outputString)//'\n',ierr)) 53c4762a1bSJed Brown end 54c4762a1bSJed Brown 55c4762a1bSJed Brown program main 56c4762a1bSJed Brown#include <petsc/finclude/petscdraw.h> 57c4762a1bSJed Brown use petscsnes 58c4762a1bSJed Brown implicit none 59c4762a1bSJed Brown! 60c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 61c4762a1bSJed Brown! Variable declarations 62c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 63c4762a1bSJed Brown! 64c4762a1bSJed Brown! Variables: 65c4762a1bSJed Brown! snes - nonlinear solver 66c4762a1bSJed Brown! x, r - solution, residual vectors 67c4762a1bSJed Brown! J - Jacobian matrix 68c4762a1bSJed Brown! its - iterations for convergence 69c4762a1bSJed Brown! matrix_free - flag - 1 indicates matrix-free version 70c4762a1bSJed Brown! lambda - nonlinearity parameter 71c4762a1bSJed Brown! draw - drawing context 72c4762a1bSJed Brown! 73c4762a1bSJed Brown SNES snes 74c4762a1bSJed Brown MatColoring mc 75c4762a1bSJed Brown Vec x,r 76c4762a1bSJed Brown PetscDraw draw 77c4762a1bSJed Brown Mat J 78c4762a1bSJed Brown PetscBool matrix_free,flg,fd_coloring 79c4762a1bSJed Brown PetscErrorCode ierr 80c4762a1bSJed Brown PetscInt its,N, mx,my,i5 81c4762a1bSJed Brown PetscMPIInt size,rank 82c4762a1bSJed Brown PetscReal lambda_max,lambda_min,lambda 83c4762a1bSJed Brown MatFDColoring fdcoloring 84c4762a1bSJed Brown ISColoring iscoloring 85c4762a1bSJed Brown PetscBool pc 86c4762a1bSJed Brown external postcheck 87c4762a1bSJed Brown 8824fb275aSStefano Zampini character(len=PETSC_MAX_PATH_LEN) :: outputString 8924fb275aSStefano Zampini 9042ce371bSBarry Smith PetscScalar,pointer :: lx_v(:) 91c4762a1bSJed Brown 92c4762a1bSJed Brown! Store parameters in common block 93c4762a1bSJed Brown 94c4762a1bSJed Brown common /params/ lambda,mx,my,fd_coloring 95c4762a1bSJed Brown 96c4762a1bSJed Brown! Note: Any user-defined Fortran routines (such as FormJacobian) 97c4762a1bSJed Brown! MUST be declared as external. 98c4762a1bSJed Brown 99c4762a1bSJed Brown external FormFunction,FormInitialGuess,FormJacobian 100c4762a1bSJed Brown 101c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 102c4762a1bSJed Brown! Initialize program 103c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 104c4762a1bSJed Brown 105d8606c27SBarry Smith PetscCallA(PetscInitialize(ierr)) 106d8606c27SBarry Smith PetscCallMPIA(MPI_Comm_size(PETSC_COMM_WORLD,size,ierr)) 107d8606c27SBarry Smith PetscCallMPIA(MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)) 108c4762a1bSJed Brown 109dcb3e689SBarry Smith PetscCheckA(size .eq. 1,PETSC_COMM_SELF,PETSC_ERR_WRONG_MPI_SIZE,'This is a uniprocessor example only') 110c4762a1bSJed Brown 111c4762a1bSJed Brown! Initialize problem parameters 112c4762a1bSJed Brown i5 = 5 113c4762a1bSJed Brown lambda_max = 6.81 114c4762a1bSJed Brown lambda_min = 0.0 115c4762a1bSJed Brown lambda = 6.0 116c4762a1bSJed Brown mx = 4 117c4762a1bSJed Brown my = 4 118d8606c27SBarry Smith PetscCallA(PetscOptionsGetInt(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-mx',mx,flg,ierr)) 119d8606c27SBarry Smith PetscCallA(PetscOptionsGetInt(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-my',my,flg,ierr)) 120d8606c27SBarry Smith PetscCallA(PetscOptionsGetReal(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-par',lambda,flg,ierr)) 121dcb3e689SBarry Smith PetscCheckA(lambda .lt. lambda_max .and. lambda .gt. lambda_min,PETSC_COMM_SELF,PETSC_ERR_USER,'Lambda out of range ') 122c4762a1bSJed Brown N = mx*my 123c4762a1bSJed Brown pc = PETSC_FALSE 124d8606c27SBarry Smith PetscCallA(PetscOptionsGetBool(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-pc',pc,PETSC_NULL_BOOL,ierr)) 125c4762a1bSJed Brown 126c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 127c4762a1bSJed Brown! Create nonlinear solver context 128c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 129c4762a1bSJed Brown 130d8606c27SBarry Smith PetscCallA(SNESCreate(PETSC_COMM_WORLD,snes,ierr)) 131c4762a1bSJed Brown 132c4762a1bSJed Brown if (pc .eqv. PETSC_TRUE) then 133d8606c27SBarry Smith PetscCallA(SNESSetType(snes,SNESNEWTONTR,ierr)) 134d8606c27SBarry Smith PetscCallA(SNESNewtonTRSetPostCheck(snes, postcheck,snes,ierr)) 135c4762a1bSJed Brown endif 136c4762a1bSJed Brown 137c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 138c4762a1bSJed Brown! Create vector data structures; set function evaluation routine 139c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 140c4762a1bSJed Brown 141d8606c27SBarry Smith PetscCallA(VecCreate(PETSC_COMM_WORLD,x,ierr)) 142d8606c27SBarry Smith PetscCallA(VecSetSizes(x,PETSC_DECIDE,N,ierr)) 143d8606c27SBarry Smith PetscCallA(VecSetFromOptions(x,ierr)) 144d8606c27SBarry Smith PetscCallA(VecDuplicate(x,r,ierr)) 145c4762a1bSJed Brown 146c4762a1bSJed Brown! Set function evaluation routine and vector. Whenever the nonlinear 147c4762a1bSJed Brown! solver needs to evaluate the nonlinear function, it will call this 148c4762a1bSJed Brown! routine. 149c4762a1bSJed Brown! - Note that the final routine argument is the user-defined 150c4762a1bSJed Brown! context that provides application-specific data for the 151c4762a1bSJed Brown! function evaluation routine. 152c4762a1bSJed Brown 153d8606c27SBarry Smith PetscCallA(SNESSetFunction(snes,r,FormFunction,fdcoloring,ierr)) 154c4762a1bSJed Brown 155c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 156c4762a1bSJed Brown! Create matrix data structure; set Jacobian evaluation routine 157c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 158c4762a1bSJed Brown 159c4762a1bSJed Brown! Create matrix. Here we only approximately preallocate storage space 160c4762a1bSJed Brown! for the Jacobian. See the users manual for a discussion of better 161c4762a1bSJed Brown! techniques for preallocating matrix memory. 162c4762a1bSJed Brown 163d8606c27SBarry Smith PetscCallA(PetscOptionsHasName(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-snes_mf',matrix_free,ierr)) 164c4762a1bSJed Brown if (.not. matrix_free) then 1655d83a8b1SBarry Smith PetscCallA(MatCreateSeqAIJ(PETSC_COMM_WORLD,N,N,i5,PETSC_NULL_INTEGER_ARRAY,J,ierr)) 166c4762a1bSJed Brown endif 167c4762a1bSJed Brown 168c4762a1bSJed Brown! 169c4762a1bSJed Brown! This option will cause the Jacobian to be computed via finite differences 170c4762a1bSJed Brown! efficiently using a coloring of the columns of the matrix. 171c4762a1bSJed Brown! 172c4762a1bSJed Brown fd_coloring = .false. 173d8606c27SBarry Smith PetscCallA(PetscOptionsHasName(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-snes_fd_coloring',fd_coloring,ierr)) 174c4762a1bSJed Brown if (fd_coloring) then 175c4762a1bSJed Brown 176c4762a1bSJed Brown! 177c4762a1bSJed Brown! This initializes the nonzero structure of the Jacobian. This is artificial 178c4762a1bSJed Brown! because clearly if we had a routine to compute the Jacobian we won't need 179c4762a1bSJed Brown! to use finite differences. 180c4762a1bSJed Brown! 181d8606c27SBarry Smith PetscCallA(FormJacobian(snes,x,J,J,0,ierr)) 182c4762a1bSJed Brown! 183c4762a1bSJed Brown! Color the matrix, i.e. determine groups of columns that share no common 184a5b23f4aSJose E. Roman! rows. These columns in the Jacobian can all be computed simultaneously. 185c4762a1bSJed Brown! 186d8606c27SBarry Smith PetscCallA(MatColoringCreate(J,mc,ierr)) 187d8606c27SBarry Smith PetscCallA(MatColoringSetType(mc,MATCOLORINGNATURAL,ierr)) 188d8606c27SBarry Smith PetscCallA(MatColoringSetFromOptions(mc,ierr)) 189d8606c27SBarry Smith PetscCallA(MatColoringApply(mc,iscoloring,ierr)) 190d8606c27SBarry Smith PetscCallA(MatColoringDestroy(mc,ierr)) 191c4762a1bSJed Brown! 192c4762a1bSJed Brown! Create the data structure that SNESComputeJacobianDefaultColor() uses 193c4762a1bSJed Brown! to compute the actual Jacobians via finite differences. 194c4762a1bSJed Brown! 195d8606c27SBarry Smith PetscCallA(MatFDColoringCreate(J,iscoloring,fdcoloring,ierr)) 196d8606c27SBarry Smith PetscCallA(MatFDColoringSetFunction(fdcoloring,FormFunction,fdcoloring,ierr)) 197d8606c27SBarry Smith PetscCallA(MatFDColoringSetFromOptions(fdcoloring,ierr)) 198d8606c27SBarry Smith PetscCallA(MatFDColoringSetUp(J,iscoloring,fdcoloring,ierr)) 199c4762a1bSJed Brown! 200c4762a1bSJed Brown! Tell SNES to use the routine SNESComputeJacobianDefaultColor() 201c4762a1bSJed Brown! to compute Jacobians. 202c4762a1bSJed Brown! 203d8606c27SBarry Smith PetscCallA(SNESSetJacobian(snes,J,J,SNESComputeJacobianDefaultColor,fdcoloring,ierr)) 204*f51a5268SBarry Smith PetscCallA(SNESGetJacobian(snes,J,PETSC_NULL_MAT,PETSC_NULL_FUNCTION,PETSC_NULL_INTEGER,ierr)) 205*f51a5268SBarry Smith PetscCallA(SNESGetJacobian(snes,J,PETSC_NULL_MAT,PETSC_NULL_FUNCTION,fdcoloring,ierr)) 206d8606c27SBarry Smith PetscCallA(ISColoringDestroy(iscoloring,ierr)) 207c4762a1bSJed Brown 208c4762a1bSJed Brown else if (.not. matrix_free) then 209c4762a1bSJed Brown 210c4762a1bSJed Brown! Set Jacobian matrix data structure and default Jacobian evaluation 211c4762a1bSJed Brown! routine. Whenever the nonlinear solver needs to compute the 212c4762a1bSJed Brown! Jacobian matrix, it will call this routine. 213c4762a1bSJed Brown! - Note that the final routine argument is the user-defined 214c4762a1bSJed Brown! context that provides application-specific data for the 215c4762a1bSJed Brown! Jacobian evaluation routine. 216c4762a1bSJed Brown! - The user can override with: 217c4762a1bSJed Brown! -snes_fd : default finite differencing approximation of Jacobian 218c4762a1bSJed Brown! -snes_mf : matrix-free Newton-Krylov method with no preconditioning 219c4762a1bSJed Brown! (unless user explicitly sets preconditioner) 220c4762a1bSJed Brown! -snes_mf_operator : form preconditioning matrix as set by the user, 221c4762a1bSJed Brown! but use matrix-free approx for Jacobian-vector 222c4762a1bSJed Brown! products within Newton-Krylov method 223c4762a1bSJed Brown! 224d8606c27SBarry Smith PetscCallA(SNESSetJacobian(snes,J,J,FormJacobian,0,ierr)) 225c4762a1bSJed Brown endif 226c4762a1bSJed Brown 227c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 228c4762a1bSJed Brown! Customize nonlinear solver; set runtime options 229c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 230c4762a1bSJed Brown 231c4762a1bSJed Brown! Set runtime options (e.g., -snes_monitor -snes_rtol <rtol> -ksp_type <type>) 232c4762a1bSJed Brown 233d8606c27SBarry Smith PetscCallA(SNESSetFromOptions(snes,ierr)) 234c4762a1bSJed Brown 235c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 236c4762a1bSJed Brown! Evaluate initial guess; then solve nonlinear system. 237c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 238c4762a1bSJed Brown 239c4762a1bSJed Brown! Note: The user should initialize the vector, x, with the initial guess 240c4762a1bSJed Brown! for the nonlinear solver prior to calling SNESSolve(). In particular, 241c4762a1bSJed Brown! to employ an initial guess of zero, the user should explicitly set 242c4762a1bSJed Brown! this vector to zero by calling VecSet(). 243c4762a1bSJed Brown 244d8606c27SBarry Smith PetscCallA(FormInitialGuess(x,ierr)) 245d8606c27SBarry Smith PetscCallA(SNESSolve(snes,PETSC_NULL_VEC,x,ierr)) 246d8606c27SBarry Smith PetscCallA(SNESGetIterationNumber(snes,its,ierr)) 24724fb275aSStefano Zampini write(outputString,*) its 24824fb275aSStefano Zampini PetscCallA(PetscPrintf(PETSC_COMM_WORLD,'Number of SNES iterations = '//trim(outputString)//'\n',ierr)) 249c4762a1bSJed Brown 250c4762a1bSJed Brown! PetscDraw contour plot of solution 251c4762a1bSJed Brown 252d8606c27SBarry Smith PetscCallA(PetscDrawCreate(PETSC_COMM_WORLD,PETSC_NULL_CHARACTER,'Solution',300,0,300,300,draw,ierr)) 253d8606c27SBarry Smith PetscCallA(PetscDrawSetFromOptions(draw,ierr)) 254c4762a1bSJed Brown 25542ce371bSBarry Smith PetscCallA(VecGetArrayReadF90(x,lx_v,ierr)) 25642ce371bSBarry Smith PetscCallA(PetscDrawTensorContour(draw,mx,my,PETSC_NULL_REAL,PETSC_NULL_REAL,lx_v,ierr)) 25742ce371bSBarry Smith PetscCallA(VecRestoreArrayReadF90(x,lx_v,ierr)) 258c4762a1bSJed Brown 259c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 260c4762a1bSJed Brown! Free work space. All PETSc objects should be destroyed when they 261c4762a1bSJed Brown! are no longer needed. 262c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 263c4762a1bSJed Brown 264d8606c27SBarry Smith if (.not. matrix_free) PetscCallA(MatDestroy(J,ierr)) 265d8606c27SBarry Smith if (fd_coloring) PetscCallA(MatFDColoringDestroy(fdcoloring,ierr)) 266c4762a1bSJed Brown 267d8606c27SBarry Smith PetscCallA(VecDestroy(x,ierr)) 268d8606c27SBarry Smith PetscCallA(VecDestroy(r,ierr)) 269d8606c27SBarry Smith PetscCallA(SNESDestroy(snes,ierr)) 270d8606c27SBarry Smith PetscCallA(PetscDrawDestroy(draw,ierr)) 271d8606c27SBarry Smith PetscCallA(PetscFinalize(ierr)) 272c4762a1bSJed Brown end 273c4762a1bSJed Brown 274c4762a1bSJed Brown! --------------------------------------------------------------------- 275c4762a1bSJed Brown! 276c4762a1bSJed Brown! FormInitialGuess - Forms initial approximation. 277c4762a1bSJed Brown! 278c4762a1bSJed Brown! Input Parameter: 279c4762a1bSJed Brown! X - vector 280c4762a1bSJed Brown! 281c4762a1bSJed Brown! Output Parameters: 282c4762a1bSJed Brown! X - vector 283c4762a1bSJed Brown! ierr - error code 284c4762a1bSJed Brown! 285c4762a1bSJed Brown! Notes: 286c4762a1bSJed Brown! This routine serves as a wrapper for the lower-level routine 287c4762a1bSJed Brown! "ApplicationInitialGuess", where the actual computations are 288c4762a1bSJed Brown! done using the standard Fortran style of treating the local 289c4762a1bSJed Brown! vector data as a multidimensional array over the local mesh. 290c4762a1bSJed Brown! This routine merely accesses the local vector data via 29142ce371bSBarry Smith! VecGetArrayF90() and VecRestoreArrayF90(). 292c4762a1bSJed Brown! 293c4762a1bSJed Brown subroutine FormInitialGuess(X,ierr) 294c4762a1bSJed Brown use petscsnes 295c4762a1bSJed Brown implicit none 296c4762a1bSJed Brown 297c4762a1bSJed Brown! Input/output variables: 298c4762a1bSJed Brown Vec X 299c4762a1bSJed Brown PetscErrorCode ierr 300c4762a1bSJed Brown 301c4762a1bSJed Brown! Declarations for use with local arrays: 30242ce371bSBarry Smith PetscScalar,pointer :: lx_v(:) 303c4762a1bSJed Brown 304c4762a1bSJed Brown ierr = 0 305c4762a1bSJed Brown 306c4762a1bSJed Brown! Get a pointer to vector data. 30742ce371bSBarry Smith! - VecGetArrayF90() returns a pointer to the data array. 30842ce371bSBarry Smith! - You MUST call VecRestoreArrayF90() when you no longer need access to 309c4762a1bSJed Brown! the array. 310c4762a1bSJed Brown 31142ce371bSBarry Smith PetscCallA(VecGetArrayF90(X,lx_v,ierr)) 312c4762a1bSJed Brown 313c4762a1bSJed Brown! Compute initial guess 314c4762a1bSJed Brown 31542ce371bSBarry Smith PetscCallA(ApplicationInitialGuess(lx_v,ierr)) 316c4762a1bSJed Brown 317c4762a1bSJed Brown! Restore vector 318c4762a1bSJed Brown 31942ce371bSBarry Smith PetscCallA(VecRestoreArrayF90(X,lx_v,ierr)) 320c4762a1bSJed Brown 321c4762a1bSJed Brown end 322c4762a1bSJed Brown 323c4762a1bSJed Brown! ApplicationInitialGuess - Computes initial approximation, called by 324c4762a1bSJed Brown! the higher level routine FormInitialGuess(). 325c4762a1bSJed Brown! 326c4762a1bSJed Brown! Input Parameter: 327c4762a1bSJed Brown! x - local vector data 328c4762a1bSJed Brown! 329c4762a1bSJed Brown! Output Parameters: 330c4762a1bSJed Brown! f - local vector data, f(x) 331c4762a1bSJed Brown! ierr - error code 332c4762a1bSJed Brown! 333c4762a1bSJed Brown! Notes: 334c4762a1bSJed Brown! This routine uses standard Fortran-style computations over a 2-dim array. 335c4762a1bSJed Brown! 336c4762a1bSJed Brown subroutine ApplicationInitialGuess(x,ierr) 337c4762a1bSJed Brown use petscksp 338c4762a1bSJed Brown implicit none 339c4762a1bSJed Brown 340c4762a1bSJed Brown! Common blocks: 341c4762a1bSJed Brown PetscReal lambda 342c4762a1bSJed Brown PetscInt mx,my 343c4762a1bSJed Brown PetscBool fd_coloring 344c4762a1bSJed Brown common /params/ lambda,mx,my,fd_coloring 345c4762a1bSJed Brown 346c4762a1bSJed Brown! Input/output variables: 347c4762a1bSJed Brown PetscScalar x(mx,my) 348c4762a1bSJed Brown PetscErrorCode ierr 349c4762a1bSJed Brown 350c4762a1bSJed Brown! Local variables: 351c4762a1bSJed Brown PetscInt i,j 352c4762a1bSJed Brown PetscReal temp1,temp,hx,hy,one 353c4762a1bSJed Brown 354c4762a1bSJed Brown! Set parameters 355c4762a1bSJed Brown 356c4762a1bSJed Brown ierr = 0 357c4762a1bSJed Brown one = 1.0 358c4762a1bSJed Brown hx = one/(mx-1) 359c4762a1bSJed Brown hy = one/(my-1) 360c4762a1bSJed Brown temp1 = lambda/(lambda + one) 361c4762a1bSJed Brown 362c4762a1bSJed Brown do 20 j=1,my 363c4762a1bSJed Brown temp = min(j-1,my-j)*hy 364c4762a1bSJed Brown do 10 i=1,mx 365c4762a1bSJed Brown if (i .eq. 1 .or. j .eq. 1 .or. i .eq. mx .or. j .eq. my) then 366c4762a1bSJed Brown x(i,j) = 0.0 367c4762a1bSJed Brown else 368c4762a1bSJed Brown x(i,j) = temp1 * sqrt(min(min(i-1,mx-i)*hx,temp)) 369c4762a1bSJed Brown endif 370c4762a1bSJed Brown 10 continue 371c4762a1bSJed Brown 20 continue 372c4762a1bSJed Brown 373c4762a1bSJed Brown end 374c4762a1bSJed Brown 375c4762a1bSJed Brown! --------------------------------------------------------------------- 376c4762a1bSJed Brown! 377c4762a1bSJed Brown! FormFunction - Evaluates nonlinear function, F(x). 378c4762a1bSJed Brown! 379c4762a1bSJed Brown! Input Parameters: 380c4762a1bSJed Brown! snes - the SNES context 381c4762a1bSJed Brown! X - input vector 382c4762a1bSJed Brown! dummy - optional user-defined context, as set by SNESSetFunction() 383c4762a1bSJed Brown! (not used here) 384c4762a1bSJed Brown! 385c4762a1bSJed Brown! Output Parameter: 386c4762a1bSJed Brown! F - vector with newly computed function 387c4762a1bSJed Brown! 388c4762a1bSJed Brown! Notes: 389c4762a1bSJed Brown! This routine serves as a wrapper for the lower-level routine 390c4762a1bSJed Brown! "ApplicationFunction", where the actual computations are 391c4762a1bSJed Brown! done using the standard Fortran style of treating the local 392c4762a1bSJed Brown! vector data as a multidimensional array over the local mesh. 393c4762a1bSJed Brown! This routine merely accesses the local vector data via 39442ce371bSBarry Smith! VecGetArrayF90() and VecRestoreArrayF90(). 395c4762a1bSJed Brown! 396c4762a1bSJed Brown subroutine FormFunction(snes,X,F,fdcoloring,ierr) 397c4762a1bSJed Brown use petscsnes 398c4762a1bSJed Brown implicit none 399c4762a1bSJed Brown 400c4762a1bSJed Brown! Input/output variables: 401c4762a1bSJed Brown SNES snes 402c4762a1bSJed Brown Vec X,F 403c4762a1bSJed Brown PetscErrorCode ierr 404c4762a1bSJed Brown MatFDColoring fdcoloring 405c4762a1bSJed Brown 406c4762a1bSJed Brown! Common blocks: 407c4762a1bSJed Brown PetscReal lambda 408c4762a1bSJed Brown PetscInt mx,my 409c4762a1bSJed Brown PetscBool fd_coloring 410c4762a1bSJed Brown common /params/ lambda,mx,my,fd_coloring 411c4762a1bSJed Brown 412c4762a1bSJed Brown! Declarations for use with local arrays: 41342ce371bSBarry Smith PetscScalar,pointer :: lx_v(:), lf_v(:) 414c4762a1bSJed Brown PetscInt, pointer :: indices(:) 415c4762a1bSJed Brown 416c4762a1bSJed Brown! Get pointers to vector data. 41742ce371bSBarry Smith! - VecGetArrayF90() returns a pointer to the data array. 41842ce371bSBarry Smith! - You MUST call VecRestoreArrayF90() when you no longer need access to 419c4762a1bSJed Brown! the array. 420c4762a1bSJed Brown 42142ce371bSBarry Smith PetscCallA(VecGetArrayReadF90(X,lx_v,ierr)) 42242ce371bSBarry Smith PetscCallA(VecGetArrayF90(F,lf_v,ierr)) 423c4762a1bSJed Brown 424c4762a1bSJed Brown! Compute function 425c4762a1bSJed Brown 42642ce371bSBarry Smith PetscCallA(ApplicationFunction(lx_v,lf_v,ierr)) 427c4762a1bSJed Brown 428c4762a1bSJed Brown! Restore vectors 429c4762a1bSJed Brown 43042ce371bSBarry Smith PetscCallA(VecRestoreArrayReadF90(X,lx_v,ierr)) 43142ce371bSBarry Smith PetscCallA(VecRestoreArrayF90(F,lf_v,ierr)) 432c4762a1bSJed Brown 433d8606c27SBarry Smith PetscCallA(PetscLogFlops(11.0d0*mx*my,ierr)) 434c4762a1bSJed Brown! 435c4762a1bSJed Brown! fdcoloring is in the common block and used here ONLY to test the 436c4762a1bSJed Brown! calls to MatFDColoringGetPerturbedColumnsF90() and MatFDColoringRestorePerturbedColumnsF90() 437c4762a1bSJed Brown! 438c4762a1bSJed Brown if (fd_coloring) then 439d8606c27SBarry Smith PetscCallA(MatFDColoringGetPerturbedColumnsF90(fdcoloring,indices,ierr)) 440c4762a1bSJed Brown print*,'Indices from GetPerturbedColumnsF90' 441c4762a1bSJed Brown write(*,1000) indices 442c4762a1bSJed Brown 1000 format(50i4) 443d8606c27SBarry Smith PetscCallA(MatFDColoringRestorePerturbedColumnsF90(fdcoloring,indices,ierr)) 444c4762a1bSJed Brown endif 445c4762a1bSJed Brown end 446c4762a1bSJed Brown 447c4762a1bSJed Brown! --------------------------------------------------------------------- 448c4762a1bSJed Brown! 449c4762a1bSJed Brown! ApplicationFunction - Computes nonlinear function, called by 450c4762a1bSJed Brown! the higher level routine FormFunction(). 451c4762a1bSJed Brown! 452c4762a1bSJed Brown! Input Parameter: 453c4762a1bSJed Brown! x - local vector data 454c4762a1bSJed Brown! 455c4762a1bSJed Brown! Output Parameters: 456c4762a1bSJed Brown! f - local vector data, f(x) 457c4762a1bSJed Brown! ierr - error code 458c4762a1bSJed Brown! 459c4762a1bSJed Brown! Notes: 460c4762a1bSJed Brown! This routine uses standard Fortran-style computations over a 2-dim array. 461c4762a1bSJed Brown! 462c4762a1bSJed Brown subroutine ApplicationFunction(x,f,ierr) 463c4762a1bSJed Brown use petscsnes 464c4762a1bSJed Brown implicit none 465c4762a1bSJed Brown 466c4762a1bSJed Brown! Common blocks: 467c4762a1bSJed Brown PetscReal lambda 468c4762a1bSJed Brown PetscInt mx,my 469c4762a1bSJed Brown PetscBool fd_coloring 470c4762a1bSJed Brown common /params/ lambda,mx,my,fd_coloring 471c4762a1bSJed Brown 472c4762a1bSJed Brown! Input/output variables: 473c4762a1bSJed Brown PetscScalar x(mx,my),f(mx,my) 474c4762a1bSJed Brown PetscErrorCode ierr 475c4762a1bSJed Brown 476c4762a1bSJed Brown! Local variables: 477c4762a1bSJed Brown PetscScalar two,one,hx,hy 478c4762a1bSJed Brown PetscScalar hxdhy,hydhx,sc 479c4762a1bSJed Brown PetscScalar u,uxx,uyy 480c4762a1bSJed Brown PetscInt i,j 481c4762a1bSJed Brown 482c4762a1bSJed Brown ierr = 0 483c4762a1bSJed Brown one = 1.0 484c4762a1bSJed Brown two = 2.0 485c4762a1bSJed Brown hx = one/(mx-1) 486c4762a1bSJed Brown hy = one/(my-1) 487c4762a1bSJed Brown sc = hx*hy*lambda 488c4762a1bSJed Brown hxdhy = hx/hy 489c4762a1bSJed Brown hydhx = hy/hx 490c4762a1bSJed Brown 491c4762a1bSJed Brown! Compute function 492c4762a1bSJed Brown 493c4762a1bSJed Brown do 20 j=1,my 494c4762a1bSJed Brown do 10 i=1,mx 495c4762a1bSJed Brown if (i .eq. 1 .or. j .eq. 1 .or. i .eq. mx .or. j .eq. my) then 496c4762a1bSJed Brown f(i,j) = x(i,j) 497c4762a1bSJed Brown else 498c4762a1bSJed Brown u = x(i,j) 499c4762a1bSJed Brown uxx = hydhx * (two*u - x(i-1,j) - x(i+1,j)) 500c4762a1bSJed Brown uyy = hxdhy * (two*u - x(i,j-1) - x(i,j+1)) 501c4762a1bSJed Brown f(i,j) = uxx + uyy - sc*exp(u) 502c4762a1bSJed Brown endif 503c4762a1bSJed Brown 10 continue 504c4762a1bSJed Brown 20 continue 505c4762a1bSJed Brown 506c4762a1bSJed Brown end 507c4762a1bSJed Brown 508c4762a1bSJed Brown! --------------------------------------------------------------------- 509c4762a1bSJed Brown! 510c4762a1bSJed Brown! FormJacobian - Evaluates Jacobian matrix. 511c4762a1bSJed Brown! 512c4762a1bSJed Brown! Input Parameters: 513c4762a1bSJed Brown! snes - the SNES context 514c4762a1bSJed Brown! x - input vector 515c4762a1bSJed Brown! dummy - optional user-defined context, as set by SNESSetJacobian() 516c4762a1bSJed Brown! (not used here) 517c4762a1bSJed Brown! 518c4762a1bSJed Brown! Output Parameters: 519c4762a1bSJed Brown! jac - Jacobian matrix 520c4762a1bSJed Brown! jac_prec - optionally different preconditioning matrix (not used here) 521c4762a1bSJed Brown! flag - flag indicating matrix structure 522c4762a1bSJed Brown! 523c4762a1bSJed Brown! Notes: 524c4762a1bSJed Brown! This routine serves as a wrapper for the lower-level routine 525c4762a1bSJed Brown! "ApplicationJacobian", where the actual computations are 526c4762a1bSJed Brown! done using the standard Fortran style of treating the local 527c4762a1bSJed Brown! vector data as a multidimensional array over the local mesh. 528c4762a1bSJed Brown! This routine merely accesses the local vector data via 52942ce371bSBarry Smith! VecGetArrayF90() and VecRestoreArrayF90(). 530c4762a1bSJed Brown! 531c4762a1bSJed Brown subroutine FormJacobian(snes,X,jac,jac_prec,dummy,ierr) 532c4762a1bSJed Brown use petscsnes 533c4762a1bSJed Brown implicit none 534c4762a1bSJed Brown 535c4762a1bSJed Brown! Input/output variables: 536c4762a1bSJed Brown SNES snes 537c4762a1bSJed Brown Vec X 538c4762a1bSJed Brown Mat jac,jac_prec 539c4762a1bSJed Brown PetscErrorCode ierr 540c4762a1bSJed Brown integer dummy 541c4762a1bSJed Brown 542c4762a1bSJed Brown! Common blocks: 543c4762a1bSJed Brown PetscReal lambda 544c4762a1bSJed Brown PetscInt mx,my 545c4762a1bSJed Brown PetscBool fd_coloring 546c4762a1bSJed Brown common /params/ lambda,mx,my,fd_coloring 547c4762a1bSJed Brown 548c4762a1bSJed Brown! Declarations for use with local array: 54942ce371bSBarry Smith PetscScalar,pointer :: lx_v(:) 550c4762a1bSJed Brown 551c4762a1bSJed Brown! Get a pointer to vector data 552c4762a1bSJed Brown 55342ce371bSBarry Smith PetscCallA(VecGetArrayReadF90(X,lx_v,ierr)) 554c4762a1bSJed Brown 555c4762a1bSJed Brown! Compute Jacobian entries 556c4762a1bSJed Brown 55742ce371bSBarry Smith PetscCallA(ApplicationJacobian(lx_v,jac,jac_prec,ierr)) 558c4762a1bSJed Brown 559c4762a1bSJed Brown! Restore vector 560c4762a1bSJed Brown 56142ce371bSBarry Smith PetscCallA(VecRestoreArrayReadF90(X,lx_v,ierr)) 562c4762a1bSJed Brown 563c4762a1bSJed Brown! Assemble matrix 564c4762a1bSJed Brown 565d8606c27SBarry Smith PetscCallA(MatAssemblyBegin(jac_prec,MAT_FINAL_ASSEMBLY,ierr)) 566d8606c27SBarry Smith PetscCallA(MatAssemblyEnd(jac_prec,MAT_FINAL_ASSEMBLY,ierr)) 567c4762a1bSJed Brown 568c4762a1bSJed Brown end 569c4762a1bSJed Brown 570c4762a1bSJed Brown! --------------------------------------------------------------------- 571c4762a1bSJed Brown! 572c4762a1bSJed Brown! ApplicationJacobian - Computes Jacobian matrix, called by 573c4762a1bSJed Brown! the higher level routine FormJacobian(). 574c4762a1bSJed Brown! 575c4762a1bSJed Brown! Input Parameters: 576c4762a1bSJed Brown! x - local vector data 577c4762a1bSJed Brown! 578c4762a1bSJed Brown! Output Parameters: 579c4762a1bSJed Brown! jac - Jacobian matrix 580c4762a1bSJed Brown! jac_prec - optionally different preconditioning matrix (not used here) 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 subroutine ApplicationJacobian(x,jac,jac_prec,ierr) 587c4762a1bSJed Brown use petscsnes 588c4762a1bSJed Brown implicit none 589c4762a1bSJed Brown 590c4762a1bSJed Brown! Common blocks: 591c4762a1bSJed Brown PetscReal lambda 592c4762a1bSJed Brown PetscInt mx,my 593c4762a1bSJed Brown PetscBool fd_coloring 594c4762a1bSJed Brown common /params/ lambda,mx,my,fd_coloring 595c4762a1bSJed Brown 596c4762a1bSJed Brown! Input/output variables: 597c4762a1bSJed Brown PetscScalar x(mx,my) 598c4762a1bSJed Brown Mat jac,jac_prec 599c4762a1bSJed Brown PetscErrorCode ierr 600c4762a1bSJed Brown 601c4762a1bSJed Brown! Local variables: 602c4762a1bSJed Brown PetscInt i,j,row(1),col(5),i1,i5 603c4762a1bSJed Brown PetscScalar two,one, hx,hy 604c4762a1bSJed Brown PetscScalar hxdhy,hydhx,sc,v(5) 605c4762a1bSJed Brown 606c4762a1bSJed Brown! Set parameters 607c4762a1bSJed Brown 608c4762a1bSJed Brown i1 = 1 609c4762a1bSJed Brown i5 = 5 610c4762a1bSJed Brown one = 1.0 611c4762a1bSJed Brown two = 2.0 612c4762a1bSJed Brown hx = one/(mx-1) 613c4762a1bSJed Brown hy = one/(my-1) 614c4762a1bSJed Brown sc = hx*hy 615c4762a1bSJed Brown hxdhy = hx/hy 616c4762a1bSJed Brown hydhx = hy/hx 617c4762a1bSJed Brown 618c4762a1bSJed Brown! Compute entries of the Jacobian matrix 619c4762a1bSJed Brown! - Here, we set all entries for a particular row at once. 620c4762a1bSJed Brown! - Note that MatSetValues() uses 0-based row and column numbers 621c4762a1bSJed Brown! in Fortran as well as in C. 622c4762a1bSJed Brown 623c4762a1bSJed Brown do 20 j=1,my 624c4762a1bSJed Brown row(1) = (j-1)*mx - 1 625c4762a1bSJed Brown do 10 i=1,mx 626c4762a1bSJed Brown row(1) = row(1) + 1 627c4762a1bSJed Brown! boundary points 628c4762a1bSJed Brown if (i .eq. 1 .or. j .eq. 1 .or. i .eq. mx .or. j .eq. my) then 6295d83a8b1SBarry Smith PetscCallA(MatSetValues(jac_prec,i1,row,i1,row,[one],INSERT_VALUES,ierr)) 630c4762a1bSJed Brown! interior grid points 631c4762a1bSJed Brown else 632c4762a1bSJed Brown v(1) = -hxdhy 633c4762a1bSJed Brown v(2) = -hydhx 634c4762a1bSJed Brown v(3) = two*(hydhx + hxdhy) - sc*lambda*exp(x(i,j)) 635c4762a1bSJed Brown v(4) = -hydhx 636c4762a1bSJed Brown v(5) = -hxdhy 637c4762a1bSJed Brown col(1) = row(1) - mx 638c4762a1bSJed Brown col(2) = row(1) - 1 639c4762a1bSJed Brown col(3) = row(1) 640c4762a1bSJed Brown col(4) = row(1) + 1 641c4762a1bSJed Brown col(5) = row(1) + mx 642d8606c27SBarry Smith PetscCallA(MatSetValues(jac_prec,i1,row,i5,col,v,INSERT_VALUES,ierr)) 643c4762a1bSJed Brown endif 644c4762a1bSJed Brown 10 continue 645c4762a1bSJed Brown 20 continue 646c4762a1bSJed Brown 647c4762a1bSJed Brown end 648c4762a1bSJed Brown 649c4762a1bSJed Brown! 650c4762a1bSJed Brown!/*TEST 651c4762a1bSJed Brown! 652c4762a1bSJed Brown! build: 653c4762a1bSJed Brown! requires: !single 654c4762a1bSJed Brown! 655c4762a1bSJed Brown! test: 656c4762a1bSJed Brown! args: -snes_monitor_short -nox -snes_type newtontr -ksp_gmres_cgs_refinement_type refine_always 657c4762a1bSJed Brown! 658c4762a1bSJed Brown! test: 659c4762a1bSJed Brown! suffix: 2 660c4762a1bSJed Brown! args: -snes_monitor_short -nox -snes_fd -ksp_gmres_cgs_refinement_type refine_always 661c4762a1bSJed Brown! 662c4762a1bSJed Brown! test: 663c4762a1bSJed Brown! suffix: 3 664c4762a1bSJed Brown! args: -snes_monitor_short -nox -snes_fd_coloring -mat_coloring_type sl -ksp_gmres_cgs_refinement_type refine_always 665c4762a1bSJed Brown! filter: sort -b 666c4762a1bSJed Brown! filter_output: sort -b 667c4762a1bSJed Brown! 668c4762a1bSJed Brown! test: 669c4762a1bSJed Brown! suffix: 4 670c4762a1bSJed Brown! args: -pc -par 6.807 -nox 671c4762a1bSJed Brown! 672c4762a1bSJed Brown!TEST*/ 673