1*c4762a1bSJed Brown! Time-dependent advection-reaction PDE in 1d. Demonstrates IMEX methods 2*c4762a1bSJed Brown! 3*c4762a1bSJed Brown! u_t + a1*u_x = -k1*u + k2*v + s1 4*c4762a1bSJed Brown! v_t + a2*v_x = k1*u - k2*v + s2 5*c4762a1bSJed Brown! 0 < x < 1; 6*c4762a1bSJed Brown! a1 = 1, k1 = 10^6, s1 = 0, 7*c4762a1bSJed Brown! a2 = 0, k2 = 2*k1, s2 = 1 8*c4762a1bSJed Brown! 9*c4762a1bSJed Brown! Initial conditions: 10*c4762a1bSJed Brown! u(x,0) = 1 + s2*x 11*c4762a1bSJed Brown! v(x,0) = k0/k1*u(x,0) + s1/k1 12*c4762a1bSJed Brown! 13*c4762a1bSJed Brown! Upstream boundary conditions: 14*c4762a1bSJed Brown! u(0,t) = 1-sin(12*t)^4 15*c4762a1bSJed Brown! 16*c4762a1bSJed Brown 17*c4762a1bSJed Brown module PETScShiftMod 18*c4762a1bSJed Brown#include <petsc/finclude/petscts.h> 19*c4762a1bSJed Brown use petscts 20*c4762a1bSJed Brown PetscScalar::PETSC_SHIFT 21*c4762a1bSJed Brown TS::tscontext 22*c4762a1bSJed Brown Mat::Jmat 23*c4762a1bSJed Brown PetscReal::MFuser(6) 24*c4762a1bSJed Brown end module PETScShiftMod 25*c4762a1bSJed Brown 26*c4762a1bSJed Brown 27*c4762a1bSJed Brownprogram main 28*c4762a1bSJed Brown use PETScShiftMod 29*c4762a1bSJed Brown use petscdmda 30*c4762a1bSJed Brown implicit none 31*c4762a1bSJed Brown 32*c4762a1bSJed Brown ! 33*c4762a1bSJed Brown ! Create an application context to contain data needed by the 34*c4762a1bSJed Brown ! application-provided call-back routines, FormJacobian() and 35*c4762a1bSJed Brown ! FormFunction(). We use a double precision array with six 36*c4762a1bSJed Brown ! entries, two for each problem parameter a, k, s. 37*c4762a1bSJed Brown ! 38*c4762a1bSJed Brown PetscReal user(6) 39*c4762a1bSJed Brown integer user_a,user_k,user_s 40*c4762a1bSJed Brown parameter (user_a = 0,user_k = 2,user_s = 4) 41*c4762a1bSJed Brown 42*c4762a1bSJed Brown external FormRHSFunction,FormIFunction 43*c4762a1bSJed Brown external FormInitialSolution 44*c4762a1bSJed Brown external FormIJacobian 45*c4762a1bSJed Brown external MyMult,FormIJacobianMF 46*c4762a1bSJed Brown 47*c4762a1bSJed Brown TS ts 48*c4762a1bSJed Brown Vec X 49*c4762a1bSJed Brown Mat J 50*c4762a1bSJed Brown PetscInt mx 51*c4762a1bSJed Brown PetscBool OptionSaveToDisk 52*c4762a1bSJed Brown PetscErrorCode ierr 53*c4762a1bSJed Brown DM da 54*c4762a1bSJed Brown PetscReal ftime,dt 55*c4762a1bSJed Brown PetscReal one,pone 56*c4762a1bSJed Brown PetscInt im11,i2 57*c4762a1bSJed Brown PetscBool flg 58*c4762a1bSJed Brown 59*c4762a1bSJed Brown 60*c4762a1bSJed Brown PetscInt xs,xe,gxs,gxe,dof,gdof 61*c4762a1bSJed Brown PetscScalar shell_shift 62*c4762a1bSJed Brown Mat A 63*c4762a1bSJed Brown 64*c4762a1bSJed Brown im11 = 11 65*c4762a1bSJed Brown i2 = 2 66*c4762a1bSJed Brown one = 1.0 67*c4762a1bSJed Brown pone = one / 10 68*c4762a1bSJed Brown 69*c4762a1bSJed Brown call PetscInitialize(PETSC_NULL_CHARACTER,ierr) 70*c4762a1bSJed Brown if (ierr .ne. 0) then 71*c4762a1bSJed Brown print*,'PetscInitialize failed' 72*c4762a1bSJed Brown stop 73*c4762a1bSJed Brown endif 74*c4762a1bSJed Brown 75*c4762a1bSJed Brown ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 76*c4762a1bSJed Brown ! Create distributed array (DMDA) to manage parallel grid and vectors 77*c4762a1bSJed Brown ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 78*c4762a1bSJed Brown call DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,im11,i2,i2,PETSC_NULL_INTEGER,da,ierr);CHKERRA(ierr) 79*c4762a1bSJed Brown call DMSetFromOptions(da,ierr);CHKERRA(ierr) 80*c4762a1bSJed Brown call DMSetUp(da,ierr);CHKERRA(ierr) 81*c4762a1bSJed Brown 82*c4762a1bSJed Brown ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 83*c4762a1bSJed Brown ! Extract global vectors from DMDA; 84*c4762a1bSJed Brown ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 85*c4762a1bSJed Brown call DMCreateGlobalVector(da,X,ierr);CHKERRA(ierr) 86*c4762a1bSJed Brown 87*c4762a1bSJed Brown ! Initialize user application context 88*c4762a1bSJed Brown ! Use zero-based indexing for command line parameters to match ex22.c 89*c4762a1bSJed Brown user(user_a+1) = 1.0 90*c4762a1bSJed Brown call PetscOptionsGetReal(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-a0',user(user_a+1),flg,ierr);CHKERRA(ierr) 91*c4762a1bSJed Brown user(user_a+2) = 0.0 92*c4762a1bSJed Brown call PetscOptionsGetReal(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-a1',user(user_a+2),flg,ierr);CHKERRA(ierr) 93*c4762a1bSJed Brown user(user_k+1) = 1000000.0 94*c4762a1bSJed Brown call PetscOptionsGetReal(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-k0', user(user_k+1),flg,ierr);CHKERRA(ierr) 95*c4762a1bSJed Brown user(user_k+2) = 2*user(user_k+1) 96*c4762a1bSJed Brown call PetscOptionsGetReal(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-k1', user(user_k+2),flg,ierr);CHKERRA(ierr) 97*c4762a1bSJed Brown user(user_s+1) = 0.0 98*c4762a1bSJed Brown call PetscOptionsGetReal(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-s0',user(user_s+1),flg,ierr);CHKERRA(ierr) 99*c4762a1bSJed Brown user(user_s+2) = 1.0 100*c4762a1bSJed Brown call PetscOptionsGetReal(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-s1',user(user_s+2),flg,ierr);CHKERRA(ierr) 101*c4762a1bSJed Brown 102*c4762a1bSJed Brown OptionSaveToDisk=.FALSE. 103*c4762a1bSJed Brown call PetscOptionsGetBool(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-sdisk',OptionSaveToDisk,flg,ierr);CHKERRA(ierr) 104*c4762a1bSJed Brown ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 105*c4762a1bSJed Brown ! Create timestepping solver context 106*c4762a1bSJed Brown ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 107*c4762a1bSJed Brown call TSCreate(PETSC_COMM_WORLD,ts,ierr);CHKERRA(ierr) 108*c4762a1bSJed Brown tscontext=ts 109*c4762a1bSJed Brown call TSSetDM(ts,da,ierr);CHKERRA(ierr) 110*c4762a1bSJed Brown call TSSetType(ts,TSARKIMEX,ierr);CHKERRA(ierr) 111*c4762a1bSJed Brown call TSSetRHSFunction(ts,PETSC_NULL_VEC,FormRHSFunction,user,ierr);CHKERRA(ierr) 112*c4762a1bSJed Brown 113*c4762a1bSJed Brown ! - - - - - - - - -- - - - - 114*c4762a1bSJed Brown ! Matrix free setup 115*c4762a1bSJed Brown call GetLayout(da,mx,xs,xe,gxs,gxe,ierr);CHKERRA(ierr) 116*c4762a1bSJed Brown dof=i2*(xe-xs+1) 117*c4762a1bSJed Brown gdof=i2*(gxe-gxs+1) 118*c4762a1bSJed Brown call MatCreateShell(PETSC_COMM_WORLD,dof,dof,gdof,gdof,shell_shift,A,ierr);CHKERRA(ierr) 119*c4762a1bSJed Brown 120*c4762a1bSJed Brown call MatShellSetOperation(A,MATOP_MULT,MyMult,ierr);CHKERRA(ierr) 121*c4762a1bSJed Brown ! - - - - - - - - - - - - 122*c4762a1bSJed Brown 123*c4762a1bSJed Brown call TSSetIFunction(ts,PETSC_NULL_VEC,FormIFunction,user,ierr);CHKERRA(ierr) 124*c4762a1bSJed Brown call DMSetMatType(da,MATAIJ,ierr);CHKERRA(ierr) 125*c4762a1bSJed Brown call DMCreateMatrix(da,J,ierr);CHKERRA(ierr) 126*c4762a1bSJed Brown 127*c4762a1bSJed Brown Jmat=J 128*c4762a1bSJed Brown 129*c4762a1bSJed Brown call TSSetIJacobian(ts,J,J,FormIJacobian,user,ierr);CHKERRA(ierr) 130*c4762a1bSJed Brown call TSSetIJacobian(ts,A,A,FormIJacobianMF,user,ierr);CHKERRA(ierr) 131*c4762a1bSJed Brown 132*c4762a1bSJed Brown 133*c4762a1bSJed Brown ftime = 1.0 134*c4762a1bSJed Brown call TSSetMaxTime(ts,ftime,ierr);CHKERRA(ierr) 135*c4762a1bSJed Brown call TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER,ierr);CHKERRA(ierr) 136*c4762a1bSJed Brown 137*c4762a1bSJed Brown ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 138*c4762a1bSJed Brown ! Set initial conditions 139*c4762a1bSJed Brown ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 140*c4762a1bSJed Brown call FormInitialSolution(ts,X,user,ierr);CHKERRA(ierr) 141*c4762a1bSJed Brown call TSSetSolution(ts,X,ierr);CHKERRA(ierr) 142*c4762a1bSJed Brown call VecGetSize(X,mx,ierr);CHKERRA(ierr) 143*c4762a1bSJed Brown ! Advective CFL, I don't know why it needs so much safety factor. 144*c4762a1bSJed Brown dt = pone * max(user(user_a+1),user(user_a+2)) / mx; 145*c4762a1bSJed Brown call TSSetTimeStep(ts,dt,ierr);CHKERRA(ierr) 146*c4762a1bSJed Brown 147*c4762a1bSJed Brown ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 148*c4762a1bSJed Brown ! Set runtime options 149*c4762a1bSJed Brown ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 150*c4762a1bSJed Brown call TSSetFromOptions(ts,ierr);CHKERRA(ierr) 151*c4762a1bSJed Brown 152*c4762a1bSJed Brown ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 153*c4762a1bSJed Brown ! Solve nonlinear system 154*c4762a1bSJed Brown ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 155*c4762a1bSJed Brown call TSSolve(ts,X,ierr);CHKERRA(ierr) 156*c4762a1bSJed Brown 157*c4762a1bSJed Brown if (OptionSaveToDisk) then 158*c4762a1bSJed Brown call GetLayout(da,mx,xs,xe,gxs,gxe,ierr);CHKERRA(ierr) 159*c4762a1bSJed Brown dof=i2*(xe-xs+1) 160*c4762a1bSJed Brown gdof=i2*(gxe-gxs+1) 161*c4762a1bSJed Brown call SaveSolutionToDisk(da,X,gdof,xs,xe);CHKERRA(ierr) 162*c4762a1bSJed Brown end if 163*c4762a1bSJed Brown 164*c4762a1bSJed Brown ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 165*c4762a1bSJed Brown ! Free work space. 166*c4762a1bSJed Brown! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 167*c4762a1bSJed Brown call MatDestroy(A,ierr);CHKERRA(ierr) 168*c4762a1bSJed Brown call MatDestroy(J,ierr);CHKERRA(ierr) 169*c4762a1bSJed Brown call VecDestroy(X,ierr);CHKERRA(ierr) 170*c4762a1bSJed Brown call TSDestroy(ts,ierr);CHKERRA(ierr) 171*c4762a1bSJed Brown call DMDestroy(da,ierr);CHKERRA(ierr) 172*c4762a1bSJed Brown call PetscFinalize(ierr) 173*c4762a1bSJed Brownend program main 174*c4762a1bSJed Brown 175*c4762a1bSJed Brown! Small helper to extract the layout, result uses 1-based indexing. 176*c4762a1bSJed Brown subroutine GetLayout(da,mx,xs,xe,gxs,gxe,ierr) 177*c4762a1bSJed Brown use petscdmda 178*c4762a1bSJed Brown implicit none 179*c4762a1bSJed Brown 180*c4762a1bSJed Brown DM da 181*c4762a1bSJed Brown PetscInt mx,xs,xe,gxs,gxe 182*c4762a1bSJed Brown PetscErrorCode ierr 183*c4762a1bSJed Brown PetscInt xm,gxm 184*c4762a1bSJed Brown call DMDAGetInfo(da,PETSC_NULL_INTEGER,mx,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER, & 185*c4762a1bSJed Brown PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,ierr);CHKERRQ(ierr) 186*c4762a1bSJed Brown call DMDAGetCorners(da,xs,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,xm,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,ierr);CHKERRQ(ierr) 187*c4762a1bSJed Brown call DMDAGetGhostCorners(da,gxs,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,gxm,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,ierr);CHKERRQ(ierr) 188*c4762a1bSJed Brown xs = xs + 1 189*c4762a1bSJed Brown gxs = gxs + 1 190*c4762a1bSJed Brown xe = xs + xm - 1 191*c4762a1bSJed Brown gxe = gxs + gxm - 1 192*c4762a1bSJed Brownend subroutine GetLayout 193*c4762a1bSJed Brown 194*c4762a1bSJed Brownsubroutine FormIFunctionLocal(mx,xs,xe,gxs,gxe,x,xdot,f,a,k,s,ierr) 195*c4762a1bSJed Brown implicit none 196*c4762a1bSJed Brown PetscInt mx,xs,xe,gxs,gxe 197*c4762a1bSJed Brown PetscScalar x(2,xs:xe) 198*c4762a1bSJed Brown PetscScalar xdot(2,xs:xe) 199*c4762a1bSJed Brown PetscScalar f(2,xs:xe) 200*c4762a1bSJed Brown PetscReal a(2),k(2),s(2) 201*c4762a1bSJed Brown PetscErrorCode ierr 202*c4762a1bSJed Brown PetscInt i 203*c4762a1bSJed Brown do i = xs,xe 204*c4762a1bSJed Brown f(1,i) = xdot(1,i) + k(1)*x(1,i) - k(2)*x(2,i) - s(1) 205*c4762a1bSJed Brown f(2,i) = xdot(2,i) - k(1)*x(1,i) + k(2)*x(2,i) - s(2) 206*c4762a1bSJed Brown end do 207*c4762a1bSJed Brownend subroutine FormIFunctionLocal 208*c4762a1bSJed Brown 209*c4762a1bSJed Brownsubroutine FormIFunction(ts,t,X,Xdot,F,user,ierr) 210*c4762a1bSJed Brown use petscdmda 211*c4762a1bSJed Brown use petscts 212*c4762a1bSJed Brown implicit none 213*c4762a1bSJed Brown 214*c4762a1bSJed Brown TS ts 215*c4762a1bSJed Brown PetscReal t 216*c4762a1bSJed Brown Vec X,Xdot,F 217*c4762a1bSJed Brown PetscReal user(6) 218*c4762a1bSJed Brown PetscErrorCode ierr 219*c4762a1bSJed Brown integer user_a,user_k,user_s 220*c4762a1bSJed Brown parameter (user_a = 1,user_k = 3,user_s = 5) 221*c4762a1bSJed Brown 222*c4762a1bSJed Brown DM da 223*c4762a1bSJed Brown PetscInt mx,xs,xe,gxs,gxe 224*c4762a1bSJed Brown PetscOffset ixx,ixxdot,iff 225*c4762a1bSJed Brown PetscScalar xx(0:1),xxdot(0:1),ff(0:1) 226*c4762a1bSJed Brown 227*c4762a1bSJed Brown call TSGetDM(ts,da,ierr);CHKERRQ(ierr) 228*c4762a1bSJed Brown call GetLayout(da,mx,xs,xe,gxs,gxe,ierr);CHKERRQ(ierr) 229*c4762a1bSJed Brown 230*c4762a1bSJed Brown ! Get access to vector data 231*c4762a1bSJed Brown call VecGetArrayRead(X,xx,ixx,ierr);CHKERRQ(ierr) 232*c4762a1bSJed Brown call VecGetArrayRead(Xdot,xxdot,ixxdot,ierr);CHKERRQ(ierr) 233*c4762a1bSJed Brown call VecGetArray(F,ff,iff,ierr);CHKERRQ(ierr) 234*c4762a1bSJed Brown 235*c4762a1bSJed Brown call FormIFunctionLocal(mx,xs,xe,gxs,gxe,xx(ixx),xxdot(ixxdot),ff(iff),user(user_a),user(user_k),user(user_s),ierr);CHKERRQ(ierr) 236*c4762a1bSJed Brown 237*c4762a1bSJed Brown call VecRestoreArrayRead(X,xx,ixx,ierr);CHKERRQ(ierr) 238*c4762a1bSJed Brown call VecRestoreArrayRead(Xdot,xxdot,ixxdot,ierr);CHKERRQ(ierr) 239*c4762a1bSJed Brown call VecRestoreArray(F,ff,iff,ierr);CHKERRQ(ierr) 240*c4762a1bSJed Brownend subroutine FormIFunction 241*c4762a1bSJed Brown 242*c4762a1bSJed Brownsubroutine FormRHSFunctionLocal(mx,xs,xe,gxs,gxe,t,x,f,a,k,s,ierr) 243*c4762a1bSJed Brown implicit none 244*c4762a1bSJed Brown PetscInt mx,xs,xe,gxs,gxe 245*c4762a1bSJed Brown PetscReal t 246*c4762a1bSJed Brown PetscScalar x(2,gxs:gxe),f(2,xs:xe) 247*c4762a1bSJed Brown PetscReal a(2),k(2),s(2) 248*c4762a1bSJed Brown PetscErrorCode ierr 249*c4762a1bSJed Brown PetscInt i,j 250*c4762a1bSJed Brown PetscReal hx,u0t(2) 251*c4762a1bSJed Brown PetscReal one,two,three,four,six,twelve 252*c4762a1bSJed Brown PetscReal half,third,twothird,sixth 253*c4762a1bSJed Brown PetscReal twelfth 254*c4762a1bSJed Brown 255*c4762a1bSJed Brown one = 1.0 256*c4762a1bSJed Brown two = 2.0 257*c4762a1bSJed Brown three = 3.0 258*c4762a1bSJed Brown four = 4.0 259*c4762a1bSJed Brown six = 6.0 260*c4762a1bSJed Brown twelve = 12.0 261*c4762a1bSJed Brown hx = one / mx 262*c4762a1bSJed Brown u0t(1) = one - sin(twelve*t)**four 263*c4762a1bSJed Brown u0t(2) = 0.0 264*c4762a1bSJed Brown half = one/two 265*c4762a1bSJed Brown third = one / three 266*c4762a1bSJed Brown twothird = two / three 267*c4762a1bSJed Brown sixth = one / six 268*c4762a1bSJed Brown twelfth = one / twelve 269*c4762a1bSJed Brown do i = xs,xe 270*c4762a1bSJed Brown do j = 1,2 271*c4762a1bSJed Brown if (i .eq. 1) then 272*c4762a1bSJed Brown f(j,i) = a(j)/hx*(third*u0t(j) + half*x(j,i) - x(j,i+1) + sixth*x(j,i+2)) 273*c4762a1bSJed Brown else if (i .eq. 2) then 274*c4762a1bSJed Brown f(j,i) = a(j)/hx*(-twelfth*u0t(j) + twothird*x(j,i-1) - twothird*x(j,i+1) + twelfth*x(j,i+2)) 275*c4762a1bSJed Brown else if (i .eq. mx-1) then 276*c4762a1bSJed Brown f(j,i) = a(j)/hx*(-sixth*x(j,i-2) + x(j,i-1) - half*x(j,i) -third*x(j,i+1)) 277*c4762a1bSJed Brown else if (i .eq. mx) then 278*c4762a1bSJed Brown f(j,i) = a(j)/hx*(-x(j,i) + x(j,i-1)) 279*c4762a1bSJed Brown else 280*c4762a1bSJed Brown f(j,i) = a(j)/hx*(-twelfth*x(j,i-2) + twothird*x(j,i-1) - twothird*x(j,i+1) + twelfth*x(j,i+2)) 281*c4762a1bSJed Brown end if 282*c4762a1bSJed Brown end do 283*c4762a1bSJed Brown end do 284*c4762a1bSJed Brown 285*c4762a1bSJed Brown#ifdef EXPLICIT_INTEGRATOR22 286*c4762a1bSJed Brown do i = xs,xe 287*c4762a1bSJed Brown f(1,i) = f(1,i) -( k(1)*x(1,i) - k(2)*x(2,i) - s(1)) 288*c4762a1bSJed Brown f(2,i) = f(2,i) -(- k(1)*x(1,i) + k(2)*x(2,i) - s(2)) 289*c4762a1bSJed Brown end do 290*c4762a1bSJed Brown#endif 291*c4762a1bSJed Brown 292*c4762a1bSJed Brownend subroutine FormRHSFunctionLocal 293*c4762a1bSJed Brown 294*c4762a1bSJed Brownsubroutine FormRHSFunction(ts,t,X,F,user,ierr) 295*c4762a1bSJed Brown use petscts 296*c4762a1bSJed Brown use petscdmda 297*c4762a1bSJed Brown implicit none 298*c4762a1bSJed Brown 299*c4762a1bSJed Brown TS ts 300*c4762a1bSJed Brown PetscReal t 301*c4762a1bSJed Brown Vec X,F 302*c4762a1bSJed Brown PetscReal user(6) 303*c4762a1bSJed Brown PetscErrorCode ierr 304*c4762a1bSJed Brown integer user_a,user_k,user_s 305*c4762a1bSJed Brown parameter (user_a = 1,user_k = 3,user_s = 5) 306*c4762a1bSJed Brown DM da 307*c4762a1bSJed Brown Vec Xloc 308*c4762a1bSJed Brown PetscInt mx,xs,xe,gxs,gxe 309*c4762a1bSJed Brown PetscOffset ixx,iff 310*c4762a1bSJed Brown PetscScalar xx(0:1),ff(0:1) 311*c4762a1bSJed Brown 312*c4762a1bSJed Brown call TSGetDM(ts,da,ierr);CHKERRQ(ierr) 313*c4762a1bSJed Brown call GetLayout(da,mx,xs,xe,gxs,gxe,ierr);CHKERRQ(ierr) 314*c4762a1bSJed Brown 315*c4762a1bSJed Brown ! Scatter ghost points to local vector,using the 2-step process 316*c4762a1bSJed Brown ! DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). 317*c4762a1bSJed Brown ! By placing code between these two statements, computations can be 318*c4762a1bSJed Brown ! done while messages are in transition. 319*c4762a1bSJed Brown call DMGetLocalVector(da,Xloc,ierr);CHKERRQ(ierr) 320*c4762a1bSJed Brown call DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc,ierr);CHKERRQ(ierr) 321*c4762a1bSJed Brown call DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc,ierr);CHKERRQ(ierr) 322*c4762a1bSJed Brown 323*c4762a1bSJed Brown ! Get access to vector data 324*c4762a1bSJed Brown call VecGetArrayRead(Xloc,xx,ixx,ierr);CHKERRQ(ierr) 325*c4762a1bSJed Brown call VecGetArray(F,ff,iff,ierr);CHKERRQ(ierr) 326*c4762a1bSJed Brown 327*c4762a1bSJed Brown call FormRHSFunctionLocal(mx,xs,xe,gxs,gxe,t,xx(ixx),ff(iff),user(user_a),user(user_k),user(user_s),ierr);CHKERRQ(ierr) 328*c4762a1bSJed Brown 329*c4762a1bSJed Brown call VecRestoreArrayRead(Xloc,xx,ixx,ierr);CHKERRQ(ierr) 330*c4762a1bSJed Brown call VecRestoreArray(F,ff,iff,ierr);CHKERRQ(ierr) 331*c4762a1bSJed Brown call DMRestoreLocalVector(da,Xloc,ierr);CHKERRQ(ierr) 332*c4762a1bSJed Brownend subroutine FormRHSFunction 333*c4762a1bSJed Brown 334*c4762a1bSJed Brown! --------------------------------------------------------------------- 335*c4762a1bSJed Brown! 336*c4762a1bSJed Brown! IJacobian - Compute IJacobian = dF/dU + shift*dF/dUdot 337*c4762a1bSJed Brown! 338*c4762a1bSJed Brownsubroutine FormIJacobian(ts,t,X,Xdot,shift,J,Jpre,user,ierr) 339*c4762a1bSJed Brown use petscts 340*c4762a1bSJed Brown use petscdmda 341*c4762a1bSJed Brown implicit none 342*c4762a1bSJed Brown 343*c4762a1bSJed Brown TS ts 344*c4762a1bSJed Brown PetscReal t,shift 345*c4762a1bSJed Brown Vec X,Xdot 346*c4762a1bSJed Brown Mat J,Jpre 347*c4762a1bSJed Brown PetscReal user(6) 348*c4762a1bSJed Brown PetscErrorCode ierr 349*c4762a1bSJed Brown integer user_a,user_k,user_s 350*c4762a1bSJed Brown parameter (user_a = 0,user_k = 2,user_s = 4) 351*c4762a1bSJed Brown 352*c4762a1bSJed Brown DM da 353*c4762a1bSJed Brown PetscInt mx,xs,xe,gxs,gxe 354*c4762a1bSJed Brown PetscInt i,i1,row,col 355*c4762a1bSJed Brown PetscReal k1,k2; 356*c4762a1bSJed Brown PetscScalar val(4) 357*c4762a1bSJed Brown 358*c4762a1bSJed Brown call TSGetDM(ts,da,ierr);CHKERRQ(ierr) 359*c4762a1bSJed Brown call GetLayout(da,mx,xs,xe,gxs,gxe,ierr);CHKERRQ(ierr) 360*c4762a1bSJed Brown 361*c4762a1bSJed Brown i1 = 1 362*c4762a1bSJed Brown k1 = user(user_k+1) 363*c4762a1bSJed Brown k2 = user(user_k+2) 364*c4762a1bSJed Brown do i = xs,xe 365*c4762a1bSJed Brown row = i-gxs 366*c4762a1bSJed Brown col = i-gxs 367*c4762a1bSJed Brown val(1) = shift + k1 368*c4762a1bSJed Brown val(2) = -k2 369*c4762a1bSJed Brown val(3) = -k1 370*c4762a1bSJed Brown val(4) = shift + k2 371*c4762a1bSJed Brown call MatSetValuesBlockedLocal(Jpre,i1,row,i1,col,val,INSERT_VALUES,ierr);CHKERRQ(ierr) 372*c4762a1bSJed Brown end do 373*c4762a1bSJed Brown call MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY,ierr);CHKERRQ(ierr) 374*c4762a1bSJed Brown call MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY,ierr);CHKERRQ(ierr) 375*c4762a1bSJed Brown if (J /= Jpre) then 376*c4762a1bSJed Brown call MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY,ierr);CHKERRQ(ierr) 377*c4762a1bSJed Brown call MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY,ierr);CHKERRQ(ierr) 378*c4762a1bSJed Brown end if 379*c4762a1bSJed Brownend subroutine FormIJacobian 380*c4762a1bSJed Brown 381*c4762a1bSJed Brownsubroutine FormInitialSolutionLocal(mx,xs,xe,gxs,gxe,x,a,k,s,ierr) 382*c4762a1bSJed Brown implicit none 383*c4762a1bSJed Brown PetscInt mx,xs,xe,gxs,gxe 384*c4762a1bSJed Brown PetscScalar x(2,xs:xe) 385*c4762a1bSJed Brown PetscReal a(2),k(2),s(2) 386*c4762a1bSJed Brown PetscErrorCode ierr 387*c4762a1bSJed Brown 388*c4762a1bSJed Brown PetscInt i 389*c4762a1bSJed Brown PetscReal one,hx,r,ik 390*c4762a1bSJed Brown one = 1.0 391*c4762a1bSJed Brown hx = one / mx 392*c4762a1bSJed Brown do i=xs,xe 393*c4762a1bSJed Brown r = i*hx 394*c4762a1bSJed Brown if (k(2) .ne. 0.0) then 395*c4762a1bSJed Brown ik = one/k(2) 396*c4762a1bSJed Brown else 397*c4762a1bSJed Brown ik = one 398*c4762a1bSJed Brown end if 399*c4762a1bSJed Brown x(1,i) = one + s(2)*r 400*c4762a1bSJed Brown x(2,i) = k(1)*ik*x(1,i) + s(2)*ik 401*c4762a1bSJed Brown end do 402*c4762a1bSJed Brownend subroutine FormInitialSolutionLocal 403*c4762a1bSJed Brown 404*c4762a1bSJed Brownsubroutine FormInitialSolution(ts,X,user,ierr) 405*c4762a1bSJed Brown use petscts 406*c4762a1bSJed Brown use petscdmda 407*c4762a1bSJed Brown implicit none 408*c4762a1bSJed Brown 409*c4762a1bSJed Brown TS ts 410*c4762a1bSJed Brown Vec X 411*c4762a1bSJed Brown PetscReal user(6) 412*c4762a1bSJed Brown PetscErrorCode ierr 413*c4762a1bSJed Brown integer user_a,user_k,user_s 414*c4762a1bSJed Brown parameter (user_a = 1,user_k = 3,user_s = 5) 415*c4762a1bSJed Brown 416*c4762a1bSJed Brown DM da 417*c4762a1bSJed Brown PetscInt mx,xs,xe,gxs,gxe 418*c4762a1bSJed Brown PetscOffset ixx 419*c4762a1bSJed Brown PetscScalar xx(0:1) 420*c4762a1bSJed Brown 421*c4762a1bSJed Brown call TSGetDM(ts,da,ierr);CHKERRQ(ierr) 422*c4762a1bSJed Brown call GetLayout(da,mx,xs,xe,gxs,gxe,ierr);CHKERRQ(ierr) 423*c4762a1bSJed Brown 424*c4762a1bSJed Brown ! Get access to vector data 425*c4762a1bSJed Brown call VecGetArray(X,xx,ixx,ierr);CHKERRQ(ierr) 426*c4762a1bSJed Brown 427*c4762a1bSJed Brown call FormInitialSolutionLocal(mx,xs,xe,gxs,gxe,xx(ixx),user(user_a),user(user_k),user(user_s),ierr);CHKERRQ(ierr) 428*c4762a1bSJed Brown 429*c4762a1bSJed Brown call VecRestoreArray(X,xx,ixx,ierr);CHKERRQ(ierr) 430*c4762a1bSJed Brownend subroutine FormInitialSolution 431*c4762a1bSJed Brown 432*c4762a1bSJed Brown 433*c4762a1bSJed Brown! --------------------------------------------------------------------- 434*c4762a1bSJed Brown! 435*c4762a1bSJed Brown! IJacobian - Compute IJacobian = dF/dU + shift*dF/dUdot 436*c4762a1bSJed Brown! 437*c4762a1bSJed Brownsubroutine FormIJacobianMF(ts,t,X,Xdot,shift,J,Jpre,user,ierr) 438*c4762a1bSJed Brown use PETScShiftMod 439*c4762a1bSJed Brown implicit none 440*c4762a1bSJed Brown TS ts 441*c4762a1bSJed Brown PetscReal t,shift 442*c4762a1bSJed Brown Vec X,Xdot 443*c4762a1bSJed Brown Mat J,Jpre 444*c4762a1bSJed Brown PetscReal user(6) 445*c4762a1bSJed Brown PetscErrorCode ierr 446*c4762a1bSJed Brown 447*c4762a1bSJed Brown ! call MatShellSetContext(J,shift,ierr) 448*c4762a1bSJed Brown PETSC_SHIFT=shift 449*c4762a1bSJed Brown MFuser=user 450*c4762a1bSJed Brown 451*c4762a1bSJed Brownend subroutine FormIJacobianMF 452*c4762a1bSJed Brown 453*c4762a1bSJed Brown! ------------------------------------------------------------------- 454*c4762a1bSJed Brown! 455*c4762a1bSJed Brown! MyMult - user provided matrix multiply 456*c4762a1bSJed Brown! 457*c4762a1bSJed Brown! Input Parameters: 458*c4762a1bSJed Brown!. X - input vector 459*c4762a1bSJed Brown! 460*c4762a1bSJed Brown! Output Parameter: 461*c4762a1bSJed Brown!. F - function vector 462*c4762a1bSJed Brown! 463*c4762a1bSJed Brownsubroutine MyMult(A,X,F,ierr) 464*c4762a1bSJed Brown use PETScShiftMod 465*c4762a1bSJed Brown implicit none 466*c4762a1bSJed Brown 467*c4762a1bSJed Brown Mat A 468*c4762a1bSJed Brown Vec X,F 469*c4762a1bSJed Brown 470*c4762a1bSJed Brown PetscErrorCode ierr 471*c4762a1bSJed Brown PetscScalar shift 472*c4762a1bSJed Brown 473*c4762a1bSJed Brown! Mat J,Jpre 474*c4762a1bSJed Brown 475*c4762a1bSJed Brown PetscReal user(6) 476*c4762a1bSJed Brown 477*c4762a1bSJed Brown integer user_a,user_k,user_s 478*c4762a1bSJed Brown parameter (user_a = 0,user_k = 2,user_s = 4) 479*c4762a1bSJed Brown 480*c4762a1bSJed Brown DM da 481*c4762a1bSJed Brown PetscInt mx,xs,xe,gxs,gxe 482*c4762a1bSJed Brown PetscInt i,i1,row,col 483*c4762a1bSJed Brown PetscReal k1,k2; 484*c4762a1bSJed Brown PetscScalar val(4) 485*c4762a1bSJed Brown 486*c4762a1bSJed Brown !call MatShellGetContext(A,shift,ierr) 487*c4762a1bSJed Brown shift=PETSC_SHIFT 488*c4762a1bSJed Brown user=MFuser 489*c4762a1bSJed Brown 490*c4762a1bSJed Brown call TSGetDM(tscontext,da,ierr) 491*c4762a1bSJed Brown call GetLayout(da,mx,xs,xe,gxs,gxe,ierr) 492*c4762a1bSJed Brown 493*c4762a1bSJed Brown i1 = 1 494*c4762a1bSJed Brown k1 = user(user_k+1) 495*c4762a1bSJed Brown k2 = user(user_k+2) 496*c4762a1bSJed Brown 497*c4762a1bSJed Brown do i = xs,xe 498*c4762a1bSJed Brown row = i-gxs 499*c4762a1bSJed Brown col = i-gxs 500*c4762a1bSJed Brown val(1) = shift + k1 501*c4762a1bSJed Brown val(2) = -k2 502*c4762a1bSJed Brown val(3) = -k1 503*c4762a1bSJed Brown val(4) = shift + k2 504*c4762a1bSJed Brown call MatSetValuesBlockedLocal(Jmat,i1,row,i1,col,val,INSERT_VALUES,ierr) 505*c4762a1bSJed Brown end do 506*c4762a1bSJed Brown 507*c4762a1bSJed Brown! call MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY,ierr) 508*c4762a1bSJed Brown! call MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY,ierr) 509*c4762a1bSJed Brown! if (J /= Jpre) then 510*c4762a1bSJed Brown call MatAssemblyBegin(Jmat,MAT_FINAL_ASSEMBLY,ierr) 511*c4762a1bSJed Brown call MatAssemblyEnd(Jmat,MAT_FINAL_ASSEMBLY,ierr) 512*c4762a1bSJed Brown! end if 513*c4762a1bSJed Brown 514*c4762a1bSJed Brown call MatMult(Jmat,X,F,ierr) 515*c4762a1bSJed Brown 516*c4762a1bSJed Brown return 517*c4762a1bSJed Brownend subroutine MyMult 518*c4762a1bSJed Brown 519*c4762a1bSJed Brown 520*c4762a1bSJed Brown! 521*c4762a1bSJed Brownsubroutine SaveSolutionToDisk(da,X,gdof,xs,xe) 522*c4762a1bSJed Brown use petscdmda 523*c4762a1bSJed Brown implicit none 524*c4762a1bSJed Brown 525*c4762a1bSJed Brown Vec X 526*c4762a1bSJed Brown DM da 527*c4762a1bSJed Brown PetscInt xs,xe,two 528*c4762a1bSJed Brown PetscInt gdof,i 529*c4762a1bSJed Brown PetscErrorCode ierr 530*c4762a1bSJed Brown PetscOffset ixx 531*c4762a1bSJed Brown PetscScalar data2(2,xs:xe),data(gdof) 532*c4762a1bSJed Brown PetscScalar xx(0:1) 533*c4762a1bSJed Brown 534*c4762a1bSJed Brown call VecGetArrayRead(X,xx,ixx,ierr) 535*c4762a1bSJed Brown 536*c4762a1bSJed Brown two = 2 537*c4762a1bSJed Brown data2=reshape(xx(ixx:ixx+gdof),(/two,xe-xs+1/)) 538*c4762a1bSJed Brown data=reshape(data2,(/gdof/)) 539*c4762a1bSJed Brown open(1020,file='solution_out_ex22f_mf.txt') 540*c4762a1bSJed Brown do i=1,gdof 541*c4762a1bSJed Brown write(1020,'(e24.16,1x)') data(i) 542*c4762a1bSJed Brown end do 543*c4762a1bSJed Brown close(1020) 544*c4762a1bSJed Brown 545*c4762a1bSJed Brown call VecRestoreArrayRead(X,xx,ixx,ierr) 546*c4762a1bSJed Brownend subroutine SaveSolutionToDisk 547*c4762a1bSJed Brown 548*c4762a1bSJed Brown!/*TEST 549*c4762a1bSJed Brown! 550*c4762a1bSJed Brown! test: 551*c4762a1bSJed Brown! args: -da_grid_x 200 -ts_arkimex_type 4 552*c4762a1bSJed Brown! requires: !single 553*c4762a1bSJed Brown! output_file: output/ex22f_mf_1.out 554*c4762a1bSJed Brown! 555*c4762a1bSJed Brown!TEST*/ 556