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