1*d8606c27SBarry Smith! 2*d8606c27SBarry Smith! 3*d8606c27SBarry Smith! This example demonstrates basic use of the SNES Fortran interface. 4*d8606c27SBarry Smith! 5*d8606c27SBarry Smith! 6*d8606c27SBarry Smith module UserModule 7*d8606c27SBarry Smith#include <petsc/finclude/petscsnes.h> 8*d8606c27SBarry Smith use petscsnes 9*d8606c27SBarry Smith type User 10*d8606c27SBarry Smith DM da 11*d8606c27SBarry Smith Vec F 12*d8606c27SBarry Smith Vec xl 13*d8606c27SBarry Smith MPI_Comm comm 14*d8606c27SBarry Smith PetscInt N 15*d8606c27SBarry Smith end type User 16*d8606c27SBarry Smith save 17*d8606c27SBarry Smith type monctx 18*d8606c27SBarry Smith PetscInt :: its,lag 19*d8606c27SBarry Smith end type monctx 20*d8606c27SBarry Smith end module 21*d8606c27SBarry Smith 22*d8606c27SBarry Smith! --------------------------------------------------------------------- 23*d8606c27SBarry Smith! --------------------------------------------------------------------- 24*d8606c27SBarry Smith! Subroutine FormMonitor 25*d8606c27SBarry Smith! This function lets up keep track of the SNES progress at each step 26*d8606c27SBarry Smith! In this routine, we determine when the Jacobian is rebuilt with the parameter 'jag' 27*d8606c27SBarry Smith! 28*d8606c27SBarry Smith! Input Parameters: 29*d8606c27SBarry Smith! snes - SNES nonlinear solver context 30*d8606c27SBarry Smith! its - current nonlinear iteration, starting from a call of SNESSolve() 31*d8606c27SBarry Smith! norm - 2-norm of current residual (may be approximate) 32*d8606c27SBarry Smith! snesm - monctx designed module (included in Snesmmod) 33*d8606c27SBarry Smith! --------------------------------------------------------------------- 34*d8606c27SBarry Smith subroutine FormMonitor(snes,its,norm,snesm,ierr) 35*d8606c27SBarry Smith use UserModule 36*d8606c27SBarry Smith implicit none 37*d8606c27SBarry Smith 38*d8606c27SBarry Smith SNES :: snes 39*d8606c27SBarry Smith PetscInt :: its,one,mone 40*d8606c27SBarry Smith PetscScalar :: norm 41*d8606c27SBarry Smith type(monctx) :: snesm 42*d8606c27SBarry Smith PetscErrorCode :: ierr 43*d8606c27SBarry Smith 44*d8606c27SBarry Smith! write(6,*) ' ' 45*d8606c27SBarry Smith! write(6,*) ' its ',its,snesm%its,'lag', 46*d8606c27SBarry Smith! & snesm%lag 47*d8606c27SBarry Smith! call flush(6) 48*d8606c27SBarry Smith if (mod(snesm%its,snesm%lag).eq.0) then 49*d8606c27SBarry Smith one = 1 50*d8606c27SBarry Smith PetscCall(SNESSetLagJacobian(snes,one,ierr)) ! build jacobian 51*d8606c27SBarry Smith else 52*d8606c27SBarry Smith mone = -1 53*d8606c27SBarry Smith PetscCall(SNESSetLagJacobian(snes,mone,ierr)) ! do NOT build jacobian 54*d8606c27SBarry Smith endif 55*d8606c27SBarry Smith snesm%its = snesm%its + 1 56*d8606c27SBarry Smith end subroutine FormMonitor 57*d8606c27SBarry Smith 58*d8606c27SBarry Smith! Note: Any user-defined Fortran routines (such as FormJacobian) 59*d8606c27SBarry Smith! MUST be declared as external. 60*d8606c27SBarry Smith! 61*d8606c27SBarry Smith! 62*d8606c27SBarry Smith! Macros to make setting/getting values into vector clearer. 63*d8606c27SBarry Smith! The element xx(ib) is the ibth element in the vector indicated by ctx%F 64*d8606c27SBarry Smith#define xx(ib) vxx(ixx + (ib)) 65*d8606c27SBarry Smith#define ff(ib) vff(iff + (ib)) 66*d8606c27SBarry Smith#define F2(ib) vF2(iF2 + (ib)) 67*d8606c27SBarry Smith program main 68*d8606c27SBarry Smith use UserModule 69*d8606c27SBarry Smith implicit none 70*d8606c27SBarry Smith type(User) ctx 71*d8606c27SBarry Smith PetscMPIInt rank,size 72*d8606c27SBarry Smith PetscErrorCode ierr 73*d8606c27SBarry Smith PetscInt N,start,end,nn,i 74*d8606c27SBarry Smith PetscInt ii,its,i1,i0,i3 75*d8606c27SBarry Smith PetscBool flg 76*d8606c27SBarry Smith SNES snes 77*d8606c27SBarry Smith Mat J 78*d8606c27SBarry Smith Vec x,r,u 79*d8606c27SBarry Smith PetscScalar xp,FF,UU,h 80*d8606c27SBarry Smith character*(10) matrixname 81*d8606c27SBarry Smith external FormJacobian,FormFunction 82*d8606c27SBarry Smith external formmonitor 83*d8606c27SBarry Smith type(monctx) :: snesm 84*d8606c27SBarry Smith 85*d8606c27SBarry Smith PetscCallA(PetscInitialize(ierr)) 86*d8606c27SBarry Smith i1 = 1 87*d8606c27SBarry Smith i0 = 0 88*d8606c27SBarry Smith i3 = 3 89*d8606c27SBarry Smith N = 10 90*d8606c27SBarry Smith PetscCallA(PetscOptionsGetInt(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-n',N,flg,ierr)) 91*d8606c27SBarry Smith h = 1.0/real(N-1) 92*d8606c27SBarry Smith ctx%N = N 93*d8606c27SBarry Smith ctx%comm = PETSC_COMM_WORLD 94*d8606c27SBarry Smith 95*d8606c27SBarry Smith PetscCallMPIA(MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)) 96*d8606c27SBarry Smith PetscCallMPIA(MPI_Comm_size(PETSC_COMM_WORLD,size,ierr)) 97*d8606c27SBarry Smith 98*d8606c27SBarry Smith! Set up data structures 99*d8606c27SBarry Smith PetscCallA(DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,N,i1,i1,PETSC_NULL_INTEGER,ctx%da,ierr)) 100*d8606c27SBarry Smith PetscCallA(DMSetFromOptions(ctx%da,ierr)) 101*d8606c27SBarry Smith PetscCallA(DMSetUp(ctx%da,ierr)) 102*d8606c27SBarry Smith PetscCallA(DMCreateGlobalVector(ctx%da,x,ierr)) 103*d8606c27SBarry Smith PetscCallA(DMCreateLocalVector(ctx%da,ctx%xl,ierr)) 104*d8606c27SBarry Smith 105*d8606c27SBarry Smith PetscCallA(PetscObjectSetName(x,'Approximate Solution',ierr)) 106*d8606c27SBarry Smith PetscCallA(VecDuplicate(x,r,ierr)) 107*d8606c27SBarry Smith PetscCallA(VecDuplicate(x,ctx%F,ierr)) 108*d8606c27SBarry Smith PetscCallA(VecDuplicate(x,U,ierr)) 109*d8606c27SBarry Smith PetscCallA(PetscObjectSetName(U,'Exact Solution',ierr)) 110*d8606c27SBarry Smith 111*d8606c27SBarry Smith PetscCallA(MatCreateAIJ(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,N,N,i3,PETSC_NULL_INTEGER,i0,PETSC_NULL_INTEGER,J,ierr)) 112*d8606c27SBarry Smith PetscCallA(MatSetOption(J,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE,ierr)) 113*d8606c27SBarry Smith PetscCallA(MatGetType(J,matrixname,ierr)) 114*d8606c27SBarry Smith 115*d8606c27SBarry Smith! Store right-hand-side of PDE and exact solution 116*d8606c27SBarry Smith PetscCallA(VecGetOwnershipRange(x,start,end,ierr)) 117*d8606c27SBarry Smith xp = h*start 118*d8606c27SBarry Smith nn = end - start 119*d8606c27SBarry Smith ii = start 120*d8606c27SBarry Smith do 10, i=0,nn-1 121*d8606c27SBarry Smith FF = 6.0*xp + (xp+1.e-12)**6.e0 122*d8606c27SBarry Smith UU = xp*xp*xp 123*d8606c27SBarry Smith PetscCallA(VecSetValues(ctx%F,i1,ii,FF,INSERT_VALUES,ierr)) 124*d8606c27SBarry Smith PetscCallA(VecSetValues(U,i1,ii,UU,INSERT_VALUES,ierr)) 125*d8606c27SBarry Smith xp = xp + h 126*d8606c27SBarry Smith ii = ii + 1 127*d8606c27SBarry Smith 10 continue 128*d8606c27SBarry Smith PetscCallA(VecAssemblyBegin(ctx%F,ierr)) 129*d8606c27SBarry Smith PetscCallA(VecAssemblyEnd(ctx%F,ierr)) 130*d8606c27SBarry Smith PetscCallA(VecAssemblyBegin(U,ierr)) 131*d8606c27SBarry Smith PetscCallA(VecAssemblyEnd(U,ierr)) 132*d8606c27SBarry Smith 133*d8606c27SBarry Smith! Create nonlinear solver 134*d8606c27SBarry Smith PetscCallA(SNESCreate(PETSC_COMM_WORLD,snes,ierr)) 135*d8606c27SBarry Smith 136*d8606c27SBarry Smith! Set various routines and options 137*d8606c27SBarry Smith PetscCallA(SNESSetFunction(snes,r,FormFunction,ctx,ierr)) 138*d8606c27SBarry Smith PetscCallA(SNESSetJacobian(snes,J,J,FormJacobian,ctx,ierr)) 139*d8606c27SBarry Smith 140*d8606c27SBarry Smith snesm%its = 0 141*d8606c27SBarry Smith PetscCallA(SNESGetLagJacobian(snes,snesm%lag,ierr)) 142*d8606c27SBarry Smith PetscCallA(SNESMonitorSet(snes,FormMonitor,snesm,PETSC_NULL_FUNCTION,ierr)) 143*d8606c27SBarry Smith PetscCallA(SNESSetFromOptions(snes,ierr)) 144*d8606c27SBarry Smith 145*d8606c27SBarry Smith! Solve nonlinear system 146*d8606c27SBarry Smith PetscCallA(FormInitialGuess(snes,x,ierr)) 147*d8606c27SBarry Smith PetscCallA(SNESSolve(snes,PETSC_NULL_VEC,x,ierr)) 148*d8606c27SBarry Smith PetscCallA(SNESGetIterationNumber(snes,its,ierr)) 149*d8606c27SBarry Smith 150*d8606c27SBarry Smith! Free work space. All PETSc objects should be destroyed when they 151*d8606c27SBarry Smith! are no longer needed. 152*d8606c27SBarry Smith PetscCallA(VecDestroy(x,ierr)) 153*d8606c27SBarry Smith PetscCallA(VecDestroy(ctx%xl,ierr)) 154*d8606c27SBarry Smith PetscCallA(VecDestroy(r,ierr)) 155*d8606c27SBarry Smith PetscCallA(VecDestroy(U,ierr)) 156*d8606c27SBarry Smith PetscCallA(VecDestroy(ctx%F,ierr)) 157*d8606c27SBarry Smith PetscCallA(MatDestroy(J,ierr)) 158*d8606c27SBarry Smith PetscCallA(SNESDestroy(snes,ierr)) 159*d8606c27SBarry Smith PetscCallA(DMDestroy(ctx%da,ierr)) 160*d8606c27SBarry Smith PetscCallA(PetscFinalize(ierr)) 161*d8606c27SBarry Smith end 162*d8606c27SBarry Smith 163*d8606c27SBarry Smith! -------------------- Evaluate Function F(x) --------------------- 164*d8606c27SBarry Smith 165*d8606c27SBarry Smith subroutine FormFunction(snes,x,f,ctx,ierr) 166*d8606c27SBarry Smith use UserModule 167*d8606c27SBarry Smith implicit none 168*d8606c27SBarry Smith SNES snes 169*d8606c27SBarry Smith Vec x,f 170*d8606c27SBarry Smith type(User) ctx 171*d8606c27SBarry Smith PetscMPIInt rank,size,zero 172*d8606c27SBarry Smith PetscInt i,s,n 173*d8606c27SBarry Smith PetscErrorCode ierr 174*d8606c27SBarry Smith PetscOffset ixx,iff,iF2 175*d8606c27SBarry Smith PetscScalar h,d,vf2(2),vxx(2),vff(2) 176*d8606c27SBarry Smith 177*d8606c27SBarry Smith zero = 0 178*d8606c27SBarry Smith PetscCallMPI(MPI_Comm_rank(ctx%comm,rank,ierr)) 179*d8606c27SBarry Smith PetscCallMPI(MPI_Comm_size(ctx%comm,size,ierr)) 180*d8606c27SBarry Smith h = 1.0/(real(ctx%N) - 1.0) 181*d8606c27SBarry Smith PetscCall(DMGlobalToLocalBegin(ctx%da,x,INSERT_VALUES,ctx%xl,ierr)) 182*d8606c27SBarry Smith PetscCall(DMGlobalToLocalEnd(ctx%da,x,INSERT_VALUES,ctx%xl,ierr)) 183*d8606c27SBarry Smith 184*d8606c27SBarry Smith PetscCall(VecGetLocalSize(ctx%xl,n,ierr)) 185*d8606c27SBarry Smith if (n .gt. 1000) then 186*d8606c27SBarry Smith print*, 'Local work array not big enough' 187*d8606c27SBarry Smith call MPI_Abort(PETSC_COMM_WORLD,zero,ierr) 188*d8606c27SBarry Smith endif 189*d8606c27SBarry Smith 190*d8606c27SBarry Smith! 191*d8606c27SBarry Smith! This sets the index ixx so that vxx(ixx+1) is the first local 192*d8606c27SBarry Smith! element in the vector indicated by ctx%xl. 193*d8606c27SBarry Smith! 194*d8606c27SBarry Smith PetscCall(VecGetArrayRead(ctx%xl,vxx,ixx,ierr)) 195*d8606c27SBarry Smith PetscCall(VecGetArray(f,vff,iff,ierr)) 196*d8606c27SBarry Smith PetscCall(VecGetArray(ctx%F,vF2,iF2,ierr)) 197*d8606c27SBarry Smith 198*d8606c27SBarry Smith d = h*h 199*d8606c27SBarry Smith 200*d8606c27SBarry Smith! 201*d8606c27SBarry Smith! Note that the array vxx() was obtained from a ghosted local vector 202*d8606c27SBarry Smith! ctx%xl while the array vff() was obtained from the non-ghosted parallel 203*d8606c27SBarry Smith! vector F. This is why there is a need for shift variable s. Since vff() 204*d8606c27SBarry Smith! does not have locations for the ghost variables we need to index in it 205*d8606c27SBarry Smith! slightly different then indexing into vxx(). For example on processor 206*d8606c27SBarry Smith! 1 (the second processor) 207*d8606c27SBarry Smith! 208*d8606c27SBarry Smith! xx(1) xx(2) xx(3) ..... 209*d8606c27SBarry Smith! ^^^^^^^ ^^^^^ ^^^^^ 210*d8606c27SBarry Smith! ghost value 1st local value 2nd local value 211*d8606c27SBarry Smith! 212*d8606c27SBarry Smith! ff(1) ff(2) 213*d8606c27SBarry Smith! ^^^^^^^ ^^^^^^^ 214*d8606c27SBarry Smith! 1st local value 2nd local value 215*d8606c27SBarry Smith! 216*d8606c27SBarry Smith if (rank .eq. 0) then 217*d8606c27SBarry Smith s = 0 218*d8606c27SBarry Smith ff(1) = xx(1) 219*d8606c27SBarry Smith else 220*d8606c27SBarry Smith s = 1 221*d8606c27SBarry Smith endif 222*d8606c27SBarry Smith 223*d8606c27SBarry Smith do 10 i=1,n-2 224*d8606c27SBarry Smith ff(i-s+1) = d*(xx(i) - 2.0*xx(i+1) + xx(i+2)) + xx(i+1)*xx(i+1) - F2(i-s+1) 225*d8606c27SBarry Smith 10 continue 226*d8606c27SBarry Smith 227*d8606c27SBarry Smith if (rank .eq. size-1) then 228*d8606c27SBarry Smith ff(n-s) = xx(n) - 1.0 229*d8606c27SBarry Smith endif 230*d8606c27SBarry Smith 231*d8606c27SBarry Smith PetscCall(VecRestoreArray(f,vff,iff,ierr)) 232*d8606c27SBarry Smith PetscCall(VecRestoreArrayRead(ctx%xl,vxx,ixx,ierr)) 233*d8606c27SBarry Smith PetscCall(VecRestoreArray(ctx%F,vF2,iF2,ierr)) 234*d8606c27SBarry Smith return 235*d8606c27SBarry Smith end 236*d8606c27SBarry Smith 237*d8606c27SBarry Smith! -------------------- Form initial approximation ----------------- 238*d8606c27SBarry Smith 239*d8606c27SBarry Smith subroutine FormInitialGuess(snes,x,ierr) 240*d8606c27SBarry Smith use UserModule 241*d8606c27SBarry Smith implicit none 242*d8606c27SBarry Smith 243*d8606c27SBarry Smith PetscErrorCode ierr 244*d8606c27SBarry Smith Vec x 245*d8606c27SBarry Smith SNES snes 246*d8606c27SBarry Smith PetscScalar five 247*d8606c27SBarry Smith 248*d8606c27SBarry Smith five = .5 249*d8606c27SBarry Smith PetscCall(VecSet(x,five,ierr)) 250*d8606c27SBarry Smith return 251*d8606c27SBarry Smith end 252*d8606c27SBarry Smith 253*d8606c27SBarry Smith! -------------------- Evaluate Jacobian -------------------- 254*d8606c27SBarry Smith 255*d8606c27SBarry Smith subroutine FormJacobian(snes,x,jac,B,ctx,ierr) 256*d8606c27SBarry Smith use UserModule 257*d8606c27SBarry Smith implicit none 258*d8606c27SBarry Smith 259*d8606c27SBarry Smith SNES snes 260*d8606c27SBarry Smith Vec x 261*d8606c27SBarry Smith Mat jac,B 262*d8606c27SBarry Smith type(User) ctx 263*d8606c27SBarry Smith PetscInt ii,istart,iend 264*d8606c27SBarry Smith PetscInt i,j,n,end,start,i1 265*d8606c27SBarry Smith PetscErrorCode ierr 266*d8606c27SBarry Smith PetscMPIInt rank,size 267*d8606c27SBarry Smith PetscOffset ixx 268*d8606c27SBarry Smith PetscScalar d,A,h,vxx(2) 269*d8606c27SBarry Smith 270*d8606c27SBarry Smith i1 = 1 271*d8606c27SBarry Smith h = 1.0/(real(ctx%N) - 1.0) 272*d8606c27SBarry Smith d = h*h 273*d8606c27SBarry Smith PetscCallMPI(MPI_Comm_rank(ctx%comm,rank,ierr)) 274*d8606c27SBarry Smith PetscCallMPI(MPI_Comm_size(ctx%comm,size,ierr)) 275*d8606c27SBarry Smith 276*d8606c27SBarry Smith PetscCall(VecGetArrayRead(x,vxx,ixx,ierr)) 277*d8606c27SBarry Smith PetscCall(VecGetOwnershipRange(x,start,end,ierr)) 278*d8606c27SBarry Smith n = end - start 279*d8606c27SBarry Smith 280*d8606c27SBarry Smith if (rank .eq. 0) then 281*d8606c27SBarry Smith A = 1.0 282*d8606c27SBarry Smith PetscCall(MatSetValues(jac,i1,start,i1,start,A,INSERT_VALUES,ierr)) 283*d8606c27SBarry Smith istart = 1 284*d8606c27SBarry Smith else 285*d8606c27SBarry Smith istart = 0 286*d8606c27SBarry Smith endif 287*d8606c27SBarry Smith if (rank .eq. size-1) then 288*d8606c27SBarry Smith i = INT(ctx%N-1) 289*d8606c27SBarry Smith A = 1.0 290*d8606c27SBarry Smith PetscCall(MatSetValues(jac,i1,i,i1,i,A,INSERT_VALUES,ierr)) 291*d8606c27SBarry Smith iend = n-1 292*d8606c27SBarry Smith else 293*d8606c27SBarry Smith iend = n 294*d8606c27SBarry Smith endif 295*d8606c27SBarry Smith do 10 i=istart,iend-1 296*d8606c27SBarry Smith ii = i + start 297*d8606c27SBarry Smith j = start + i - 1 298*d8606c27SBarry Smith PetscCall(MatSetValues(jac,i1,ii,i1,j,d,INSERT_VALUES,ierr)) 299*d8606c27SBarry Smith j = start + i + 1 300*d8606c27SBarry Smith PetscCall(MatSetValues(jac,i1,ii,i1,j,d,INSERT_VALUES,ierr)) 301*d8606c27SBarry Smith A = -2.0*d + 2.0*xx(i+1) 302*d8606c27SBarry Smith PetscCall(MatSetValues(jac,i1,ii,i1,ii,A,INSERT_VALUES,ierr)) 303*d8606c27SBarry Smith 10 continue 304*d8606c27SBarry Smith PetscCall(VecRestoreArrayRead(x,vxx,ixx,ierr)) 305*d8606c27SBarry Smith PetscCall(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY,ierr)) 306*d8606c27SBarry Smith PetscCall(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY,ierr)) 307*d8606c27SBarry Smith return 308*d8606c27SBarry Smith end 309*d8606c27SBarry Smith 310*d8606c27SBarry Smith!/*TEST 311*d8606c27SBarry Smith! 312*d8606c27SBarry Smith! test: 313*d8606c27SBarry Smith! nsize: 2 314*d8606c27SBarry Smith! args: -ksp_gmres_cgs_refinement_type refine_always -n 10 -snes_monitor_short 315*d8606c27SBarry Smith! output_file: output/ex12_1.out 316*d8606c27SBarry Smith! 317*d8606c27SBarry Smith!TEST*/ 318