1*c4762a1bSJed Brown! 2*c4762a1bSJed Brown! 3*c4762a1bSJed Brown! Description: Uses the Newton method to solve a two-variable system. 4*c4762a1bSJed Brown! 5*c4762a1bSJed Brown!!/*T 6*c4762a1bSJed Brown! Concepts: SNES^basic uniprocessor example 7*c4762a1bSJed Brown! Processors: 1 8*c4762a1bSJed Brown!T*/ 9*c4762a1bSJed Brown 10*c4762a1bSJed Brown 11*c4762a1bSJed Brown! 12*c4762a1bSJed Brown! ----------------------------------------------------------------------- 13*c4762a1bSJed Brown 14*c4762a1bSJed Brown program main 15*c4762a1bSJed Brown#include <petsc/finclude/petsc.h> 16*c4762a1bSJed Brown use petsc 17*c4762a1bSJed Brown implicit none 18*c4762a1bSJed Brown 19*c4762a1bSJed Brown 20*c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 21*c4762a1bSJed Brown! Variable declarations 22*c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 23*c4762a1bSJed Brown! 24*c4762a1bSJed Brown! Variables: 25*c4762a1bSJed Brown! snes - nonlinear solver 26*c4762a1bSJed Brown! ksp - linear solver 27*c4762a1bSJed Brown! pc - preconditioner context 28*c4762a1bSJed Brown! ksp - Krylov subspace method context 29*c4762a1bSJed Brown! x, r - solution, residual vectors 30*c4762a1bSJed Brown! J - Jacobian matrix 31*c4762a1bSJed Brown! its - iterations for convergence 32*c4762a1bSJed Brown! 33*c4762a1bSJed Brown SNES snes 34*c4762a1bSJed Brown PC pc 35*c4762a1bSJed Brown KSP ksp 36*c4762a1bSJed Brown Vec x,r 37*c4762a1bSJed Brown Mat J 38*c4762a1bSJed Brown SNESLineSearch linesearch 39*c4762a1bSJed Brown PetscErrorCode ierr 40*c4762a1bSJed Brown PetscInt its,i2,i20 41*c4762a1bSJed Brown PetscMPIInt size,rank 42*c4762a1bSJed Brown PetscScalar pfive 43*c4762a1bSJed Brown PetscReal tol 44*c4762a1bSJed Brown PetscBool setls 45*c4762a1bSJed Brown PetscViewer viewer 46*c4762a1bSJed Brown double precision threshold,oldthreshold 47*c4762a1bSJed Brown 48*c4762a1bSJed Brown! Note: Any user-defined Fortran routines (such as FormJacobian) 49*c4762a1bSJed Brown! MUST be declared as external. 50*c4762a1bSJed Brown 51*c4762a1bSJed Brown external FormFunction, FormJacobian, MyLineSearch 52*c4762a1bSJed Brown 53*c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 54*c4762a1bSJed Brown! Macro definitions 55*c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 56*c4762a1bSJed Brown! 57*c4762a1bSJed Brown! Macros to make clearer the process of setting values in vectors and 58*c4762a1bSJed Brown! getting values from vectors. These vectors are used in the routines 59*c4762a1bSJed Brown! FormFunction() and FormJacobian(). 60*c4762a1bSJed Brown! - The element lx_a(ib) is element ib in the vector x 61*c4762a1bSJed Brown! 62*c4762a1bSJed Brown#define lx_a(ib) lx_v(lx_i + (ib)) 63*c4762a1bSJed Brown#define lf_a(ib) lf_v(lf_i + (ib)) 64*c4762a1bSJed Brown! 65*c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 66*c4762a1bSJed Brown! Beginning of program 67*c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 68*c4762a1bSJed Brown 69*c4762a1bSJed Brown call PetscInitialize(PETSC_NULL_CHARACTER,ierr) 70*c4762a1bSJed Brown if (ierr .ne. 0) then 71*c4762a1bSJed Brown print*,'Unable to initialize PETSc' 72*c4762a1bSJed Brown stop 73*c4762a1bSJed Brown endif 74*c4762a1bSJed Brown call PetscLogNestedBegin(ierr);CHKERRA(ierr); 75*c4762a1bSJed Brown threshold = 1.0 76*c4762a1bSJed Brown call PetscLogSetThreshold(threshold,oldthreshold,ierr) 77*c4762a1bSJed Brown! dummy test of logging a reduction 78*c4762a1bSJed Brown ierr = PetscAReduce() 79*c4762a1bSJed Brown call MPI_Comm_size(PETSC_COMM_WORLD,size,ierr) 80*c4762a1bSJed Brown call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr) 81*c4762a1bSJed Brown if (size .ne. 1) then; SETERRA(PETSC_COMM_SELF,PETSC_ERR_WRONG_MPI_SIZE,'Uniprocessor example'); endif 82*c4762a1bSJed Brown 83*c4762a1bSJed Brown i2 = 2 84*c4762a1bSJed Brown i20 = 20 85*c4762a1bSJed Brown! - - - - - - - - - -- - - - - - - - - - - - - - - - - - - - - - - - - - 86*c4762a1bSJed Brown! Create nonlinear solver context 87*c4762a1bSJed Brown! - - - - - - - - - -- - - - - - - - - - - - - - - - - - - - - - - - - - 88*c4762a1bSJed Brown 89*c4762a1bSJed Brown call SNESCreate(PETSC_COMM_WORLD,snes,ierr) 90*c4762a1bSJed Brown 91*c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 92*c4762a1bSJed Brown! Create matrix and vector data structures; set corresponding routines 93*c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 94*c4762a1bSJed Brown 95*c4762a1bSJed Brown! Create vectors for solution and nonlinear function 96*c4762a1bSJed Brown 97*c4762a1bSJed Brown call VecCreateSeq(PETSC_COMM_SELF,i2,x,ierr) 98*c4762a1bSJed Brown call VecDuplicate(x,r,ierr) 99*c4762a1bSJed Brown 100*c4762a1bSJed Brown! Create Jacobian matrix data structure 101*c4762a1bSJed Brown 102*c4762a1bSJed Brown call MatCreate(PETSC_COMM_SELF,J,ierr) 103*c4762a1bSJed Brown call MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,i2,i2,ierr) 104*c4762a1bSJed Brown call MatSetFromOptions(J,ierr) 105*c4762a1bSJed Brown call MatSetUp(J,ierr) 106*c4762a1bSJed Brown 107*c4762a1bSJed Brown! Set function evaluation routine and vector 108*c4762a1bSJed Brown 109*c4762a1bSJed Brown call SNESSetFunction(snes,r,FormFunction,0,ierr) 110*c4762a1bSJed Brown 111*c4762a1bSJed Brown! Set Jacobian matrix data structure and Jacobian evaluation routine 112*c4762a1bSJed Brown 113*c4762a1bSJed Brown call SNESSetJacobian(snes,J,J,FormJacobian,0,ierr) 114*c4762a1bSJed Brown 115*c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 116*c4762a1bSJed Brown! Customize nonlinear solver; set runtime options 117*c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 118*c4762a1bSJed Brown 119*c4762a1bSJed Brown! Set linear solver defaults for this problem. By extracting the 120*c4762a1bSJed Brown! KSP, KSP, and PC contexts from the SNES context, we can then 121*c4762a1bSJed Brown! directly call any KSP, KSP, and PC routines to set various options. 122*c4762a1bSJed Brown 123*c4762a1bSJed Brown call SNESGetKSP(snes,ksp,ierr) 124*c4762a1bSJed Brown call KSPGetPC(ksp,pc,ierr) 125*c4762a1bSJed Brown call PCSetType(pc,PCNONE,ierr) 126*c4762a1bSJed Brown tol = 1.e-4 127*c4762a1bSJed Brown call KSPSetTolerances(ksp,tol,PETSC_DEFAULT_REAL, & 128*c4762a1bSJed Brown & PETSC_DEFAULT_REAL,i20,ierr) 129*c4762a1bSJed Brown 130*c4762a1bSJed Brown! Set SNES/KSP/KSP/PC runtime options, e.g., 131*c4762a1bSJed Brown! -snes_view -snes_monitor -ksp_type <ksp> -pc_type <pc> 132*c4762a1bSJed Brown! These options will override those specified above as long as 133*c4762a1bSJed Brown! SNESSetFromOptions() is called _after_ any other customization 134*c4762a1bSJed Brown! routines. 135*c4762a1bSJed Brown 136*c4762a1bSJed Brown 137*c4762a1bSJed Brown call SNESSetFromOptions(snes,ierr) 138*c4762a1bSJed Brown 139*c4762a1bSJed Brown call PetscOptionsHasName(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER, & 140*c4762a1bSJed Brown & '-setls',setls,ierr) 141*c4762a1bSJed Brown 142*c4762a1bSJed Brown if (setls) then 143*c4762a1bSJed Brown call SNESGetLineSearch(snes, linesearch, ierr) 144*c4762a1bSJed Brown call SNESLineSearchSetType(linesearch, 'shell', ierr) 145*c4762a1bSJed Brown call SNESLineSearchShellSetUserFunc(linesearch, MyLineSearch, & 146*c4762a1bSJed Brown & 0, ierr) 147*c4762a1bSJed Brown endif 148*c4762a1bSJed Brown 149*c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 150*c4762a1bSJed Brown! Evaluate initial guess; then solve nonlinear system 151*c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 152*c4762a1bSJed Brown 153*c4762a1bSJed Brown! Note: The user should initialize the vector, x, with the initial guess 154*c4762a1bSJed Brown! for the nonlinear solver prior to calling SNESSolve(). In particular, 155*c4762a1bSJed Brown! to employ an initial guess of zero, the user should explicitly set 156*c4762a1bSJed Brown! this vector to zero by calling VecSet(). 157*c4762a1bSJed Brown 158*c4762a1bSJed Brown pfive = 0.5 159*c4762a1bSJed Brown call VecSet(x,pfive,ierr) 160*c4762a1bSJed Brown call SNESSolve(snes,PETSC_NULL_VEC,x,ierr) 161*c4762a1bSJed Brown call SNESGetIterationNumber(snes,its,ierr); 162*c4762a1bSJed Brown if (rank .eq. 0) then 163*c4762a1bSJed Brown write(6,100) its 164*c4762a1bSJed Brown endif 165*c4762a1bSJed Brown 100 format('Number of SNES iterations = ',i5) 166*c4762a1bSJed Brown 167*c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 168*c4762a1bSJed Brown! Free work space. All PETSc objects should be destroyed when they 169*c4762a1bSJed Brown! are no longer needed. 170*c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 171*c4762a1bSJed Brown 172*c4762a1bSJed Brown call VecDestroy(x,ierr) 173*c4762a1bSJed Brown call VecDestroy(r,ierr) 174*c4762a1bSJed Brown call MatDestroy(J,ierr) 175*c4762a1bSJed Brown call SNESDestroy(snes,ierr) 176*c4762a1bSJed Brown call PetscViewerASCIIOpen(PETSC_COMM_WORLD,'filename.xml',viewer,ierr) 177*c4762a1bSJed Brown call PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_XML,ierr) 178*c4762a1bSJed Brown call PetscLogView(viewer,ierr) 179*c4762a1bSJed Brown call PetscViewerDestroy(viewer,ierr) 180*c4762a1bSJed Brown call PetscFinalize(ierr) 181*c4762a1bSJed Brown end 182*c4762a1bSJed Brown! 183*c4762a1bSJed Brown! ------------------------------------------------------------------------ 184*c4762a1bSJed Brown! 185*c4762a1bSJed Brown! FormFunction - Evaluates nonlinear function, F(x). 186*c4762a1bSJed Brown! 187*c4762a1bSJed Brown! Input Parameters: 188*c4762a1bSJed Brown! snes - the SNES context 189*c4762a1bSJed Brown! x - input vector 190*c4762a1bSJed Brown! dummy - optional user-defined context (not used here) 191*c4762a1bSJed Brown! 192*c4762a1bSJed Brown! Output Parameter: 193*c4762a1bSJed Brown! f - function vector 194*c4762a1bSJed Brown! 195*c4762a1bSJed Brown subroutine FormFunction(snes,x,f,dummy,ierr) 196*c4762a1bSJed Brown use petscsnes 197*c4762a1bSJed Brown implicit none 198*c4762a1bSJed Brown 199*c4762a1bSJed Brown SNES snes 200*c4762a1bSJed Brown Vec x,f 201*c4762a1bSJed Brown PetscErrorCode ierr 202*c4762a1bSJed Brown integer dummy(*) 203*c4762a1bSJed Brown 204*c4762a1bSJed Brown! Declarations for use with local arrays 205*c4762a1bSJed Brown 206*c4762a1bSJed Brown PetscScalar lx_v(2),lf_v(2) 207*c4762a1bSJed Brown PetscOffset lx_i,lf_i 208*c4762a1bSJed Brown 209*c4762a1bSJed Brown! Get pointers to vector data. 210*c4762a1bSJed Brown! - For default PETSc vectors, VecGetArray() returns a pointer to 211*c4762a1bSJed Brown! the data array. Otherwise, the routine is implementation dependent. 212*c4762a1bSJed Brown! - You MUST call VecRestoreArray() when you no longer need access to 213*c4762a1bSJed Brown! the array. 214*c4762a1bSJed Brown! - Note that the Fortran interface to VecGetArray() differs from the 215*c4762a1bSJed Brown! C version. See the Fortran chapter of the users manual for details. 216*c4762a1bSJed Brown 217*c4762a1bSJed Brown call VecGetArrayRead(x,lx_v,lx_i,ierr) 218*c4762a1bSJed Brown call VecGetArray(f,lf_v,lf_i,ierr) 219*c4762a1bSJed Brown 220*c4762a1bSJed Brown! Compute function 221*c4762a1bSJed Brown 222*c4762a1bSJed Brown lf_a(1) = lx_a(1)*lx_a(1) & 223*c4762a1bSJed Brown & + lx_a(1)*lx_a(2) - 3.0 224*c4762a1bSJed Brown lf_a(2) = lx_a(1)*lx_a(2) & 225*c4762a1bSJed Brown & + lx_a(2)*lx_a(2) - 6.0 226*c4762a1bSJed Brown 227*c4762a1bSJed Brown! Restore vectors 228*c4762a1bSJed Brown 229*c4762a1bSJed Brown call VecRestoreArrayRead(x,lx_v,lx_i,ierr) 230*c4762a1bSJed Brown call VecRestoreArray(f,lf_v,lf_i,ierr) 231*c4762a1bSJed Brown 232*c4762a1bSJed Brown return 233*c4762a1bSJed Brown end 234*c4762a1bSJed Brown 235*c4762a1bSJed Brown! --------------------------------------------------------------------- 236*c4762a1bSJed Brown! 237*c4762a1bSJed Brown! FormJacobian - Evaluates Jacobian matrix. 238*c4762a1bSJed Brown! 239*c4762a1bSJed Brown! Input Parameters: 240*c4762a1bSJed Brown! snes - the SNES context 241*c4762a1bSJed Brown! x - input vector 242*c4762a1bSJed Brown! dummy - optional user-defined context (not used here) 243*c4762a1bSJed Brown! 244*c4762a1bSJed Brown! Output Parameters: 245*c4762a1bSJed Brown! A - Jacobian matrix 246*c4762a1bSJed Brown! B - optionally different preconditioning matrix 247*c4762a1bSJed Brown! flag - flag indicating matrix structure 248*c4762a1bSJed Brown! 249*c4762a1bSJed Brown subroutine FormJacobian(snes,X,jac,B,dummy,ierr) 250*c4762a1bSJed Brown use petscsnes 251*c4762a1bSJed Brown implicit none 252*c4762a1bSJed Brown 253*c4762a1bSJed Brown SNES snes 254*c4762a1bSJed Brown Vec X 255*c4762a1bSJed Brown Mat jac,B 256*c4762a1bSJed Brown PetscScalar A(4) 257*c4762a1bSJed Brown PetscErrorCode ierr 258*c4762a1bSJed Brown PetscInt idx(2),i2 259*c4762a1bSJed Brown integer dummy(*) 260*c4762a1bSJed Brown 261*c4762a1bSJed Brown! Declarations for use with local arrays 262*c4762a1bSJed Brown 263*c4762a1bSJed Brown PetscScalar lx_v(2) 264*c4762a1bSJed Brown PetscOffset lx_i 265*c4762a1bSJed Brown 266*c4762a1bSJed Brown! Get pointer to vector data 267*c4762a1bSJed Brown 268*c4762a1bSJed Brown i2 = 2 269*c4762a1bSJed Brown call VecGetArrayRead(x,lx_v,lx_i,ierr) 270*c4762a1bSJed Brown 271*c4762a1bSJed Brown! Compute Jacobian entries and insert into matrix. 272*c4762a1bSJed Brown! - Since this is such a small problem, we set all entries for 273*c4762a1bSJed Brown! the matrix at once. 274*c4762a1bSJed Brown! - Note that MatSetValues() uses 0-based row and column numbers 275*c4762a1bSJed Brown! in Fortran as well as in C (as set here in the array idx). 276*c4762a1bSJed Brown 277*c4762a1bSJed Brown idx(1) = 0 278*c4762a1bSJed Brown idx(2) = 1 279*c4762a1bSJed Brown A(1) = 2.0*lx_a(1) + lx_a(2) 280*c4762a1bSJed Brown A(2) = lx_a(1) 281*c4762a1bSJed Brown A(3) = lx_a(2) 282*c4762a1bSJed Brown A(4) = lx_a(1) + 2.0*lx_a(2) 283*c4762a1bSJed Brown call MatSetValues(B,i2,idx,i2,idx,A,INSERT_VALUES,ierr) 284*c4762a1bSJed Brown 285*c4762a1bSJed Brown! Restore vector 286*c4762a1bSJed Brown 287*c4762a1bSJed Brown call VecRestoreArrayRead(x,lx_v,lx_i,ierr) 288*c4762a1bSJed Brown 289*c4762a1bSJed Brown! Assemble matrix 290*c4762a1bSJed Brown 291*c4762a1bSJed Brown call MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY,ierr) 292*c4762a1bSJed Brown call MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY,ierr) 293*c4762a1bSJed Brown if (B .ne. jac) then 294*c4762a1bSJed Brown call MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY,ierr) 295*c4762a1bSJed Brown call MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY,ierr) 296*c4762a1bSJed Brown endif 297*c4762a1bSJed Brown 298*c4762a1bSJed Brown return 299*c4762a1bSJed Brown end 300*c4762a1bSJed Brown 301*c4762a1bSJed Brown 302*c4762a1bSJed Brown subroutine MyLineSearch(linesearch, lctx, ierr) 303*c4762a1bSJed Brown use petscsnes 304*c4762a1bSJed Brown implicit none 305*c4762a1bSJed Brown 306*c4762a1bSJed Brown SNESLineSearch linesearch 307*c4762a1bSJed Brown SNES snes 308*c4762a1bSJed Brown integer lctx 309*c4762a1bSJed Brown Vec x, f,g, y, w 310*c4762a1bSJed Brown PetscReal ynorm,gnorm,xnorm 311*c4762a1bSJed Brown PetscBool flag 312*c4762a1bSJed Brown PetscErrorCode ierr 313*c4762a1bSJed Brown 314*c4762a1bSJed Brown PetscScalar mone 315*c4762a1bSJed Brown 316*c4762a1bSJed Brown mone = -1.0 317*c4762a1bSJed Brown call SNESLineSearchGetSNES(linesearch, snes, ierr) 318*c4762a1bSJed Brown call SNESLineSearchGetVecs(linesearch, x, f, y, w, g, ierr) 319*c4762a1bSJed Brown call VecNorm(y,NORM_2,ynorm,ierr) 320*c4762a1bSJed Brown call VecAXPY(x,mone,y,ierr) 321*c4762a1bSJed Brown call SNESComputeFunction(snes,x,f,ierr) 322*c4762a1bSJed Brown call VecNorm(f,NORM_2,gnorm,ierr) 323*c4762a1bSJed Brown call VecNorm(x,NORM_2,xnorm,ierr) 324*c4762a1bSJed Brown call VecNorm(y,NORM_2,ynorm,ierr) 325*c4762a1bSJed Brown call SNESLineSearchSetNorms(linesearch, xnorm, gnorm, ynorm, & 326*c4762a1bSJed Brown & ierr) 327*c4762a1bSJed Brown flag = PETSC_FALSE 328*c4762a1bSJed Brown return 329*c4762a1bSJed Brown end 330*c4762a1bSJed Brown 331*c4762a1bSJed Brown!/*TEST 332*c4762a1bSJed Brown! 333*c4762a1bSJed Brown! test: 334*c4762a1bSJed Brown! args: -ksp_gmres_cgs_refinement_type refine_always -snes_monitor_short 335*c4762a1bSJed Brown! requires: !single 336*c4762a1bSJed Brown! 337*c4762a1bSJed Brown!TEST*/ 338