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