xref: /petsc/src/ts/tutorials/ex22f.F90 (revision d8606c274c09e255c003062beb17b1be973467bc)
1*d8606c27SBarry Smith! Time-dependent advection-reaction PDE in 1d. Demonstrates IMEX methods
2*d8606c27SBarry Smith!
3*d8606c27SBarry Smith!   u_t + a1*u_x = -k1*u + k2*v + s1
4*d8606c27SBarry Smith!   v_t + a2*v_x = k1*u - k2*v + s2
5*d8606c27SBarry Smith!   0 < x < 1;
6*d8606c27SBarry Smith!   a1 = 1, k1 = 10^6, s1 = 0,
7*d8606c27SBarry Smith!   a2 = 0, k2 = 2*k1, s2 = 1
8*d8606c27SBarry Smith!
9*d8606c27SBarry Smith!   Initial conditions:
10*d8606c27SBarry Smith!   u(x,0) = 1 + s2*x
11*d8606c27SBarry Smith!   v(x,0) = k0/k1*u(x,0) + s1/k1
12*d8606c27SBarry Smith!
13*d8606c27SBarry Smith!   Upstream boundary conditions:
14*d8606c27SBarry Smith!   u(0,t) = 1-sin(12*t)^4
15*d8606c27SBarry Smith!
16*d8606c27SBarry Smith
17*d8606c27SBarry Smith      program main
18*d8606c27SBarry Smith#include <petsc/finclude/petscts.h>
19*d8606c27SBarry Smith#include <petsc/finclude/petscdmda.h>
20*d8606c27SBarry Smith      use petscts
21*d8606c27SBarry Smith      implicit none
22*d8606c27SBarry Smith
23*d8606c27SBarry Smith!
24*d8606c27SBarry Smith!  Create an application context to contain data needed by the
25*d8606c27SBarry Smith!  application-provided call-back routines, FormJacobian() and
26*d8606c27SBarry Smith!  FormFunction(). We use a double precision array with six
27*d8606c27SBarry Smith!  entries, two for each problem parameter a, k, s.
28*d8606c27SBarry Smith!
29*d8606c27SBarry Smith      PetscReal user(6)
30*d8606c27SBarry Smith      integer user_a,user_k,user_s
31*d8606c27SBarry Smith      parameter (user_a = 0,user_k = 2,user_s = 4)
32*d8606c27SBarry Smith
33*d8606c27SBarry Smith      external FormRHSFunction,FormIFunction,FormIJacobian
34*d8606c27SBarry Smith      external FormInitialSolution
35*d8606c27SBarry Smith
36*d8606c27SBarry Smith      TS             ts
37*d8606c27SBarry Smith      SNES           snes
38*d8606c27SBarry Smith      SNESLineSearch linesearch
39*d8606c27SBarry Smith      Vec            X
40*d8606c27SBarry Smith      Mat            J
41*d8606c27SBarry Smith      PetscInt       mx
42*d8606c27SBarry Smith      PetscErrorCode ierr
43*d8606c27SBarry Smith      DM             da
44*d8606c27SBarry Smith      PetscReal      ftime,dt
45*d8606c27SBarry Smith      PetscReal      one,pone
46*d8606c27SBarry Smith      PetscInt       im11,i2
47*d8606c27SBarry Smith      PetscBool      flg
48*d8606c27SBarry Smith
49*d8606c27SBarry Smith      im11 = 11
50*d8606c27SBarry Smith      i2   = 2
51*d8606c27SBarry Smith      one = 1.0
52*d8606c27SBarry Smith      pone = one / 10
53*d8606c27SBarry Smith
54*d8606c27SBarry Smith      PetscCallA(PetscInitialize(ierr))
55*d8606c27SBarry Smith
56*d8606c27SBarry Smith! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
57*d8606c27SBarry Smith!  Create distributed array (DMDA) to manage parallel grid and vectors
58*d8606c27SBarry Smith! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
59*d8606c27SBarry Smith      PetscCallA(DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,im11,i2,i2,PETSC_NULL_INTEGER,da,ierr))
60*d8606c27SBarry Smith      PetscCallA(DMSetFromOptions(da,ierr))
61*d8606c27SBarry Smith      PetscCallA(DMSetUp(da,ierr))
62*d8606c27SBarry Smith
63*d8606c27SBarry Smith! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
64*d8606c27SBarry Smith!    Extract global vectors from DMDA;
65*d8606c27SBarry Smith! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
66*d8606c27SBarry Smith      PetscCallA(DMCreateGlobalVector(da,X,ierr))
67*d8606c27SBarry Smith
68*d8606c27SBarry Smith! Initialize user application context
69*d8606c27SBarry Smith! Use zero-based indexing for command line parameters to match ex22.c
70*d8606c27SBarry Smith      user(user_a+1) = 1.0
71*d8606c27SBarry Smith      PetscCallA(PetscOptionsGetReal(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-a0',user(user_a+1),flg,ierr))
72*d8606c27SBarry Smith      user(user_a+2) = 0.0
73*d8606c27SBarry Smith      PetscCallA(PetscOptionsGetReal(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-a1',user(user_a+2),flg,ierr))
74*d8606c27SBarry Smith      user(user_k+1) = 1000000.0
75*d8606c27SBarry Smith      PetscCallA(PetscOptionsGetReal(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-k0',user(user_k+1),flg,ierr))
76*d8606c27SBarry Smith      user(user_k+2) = 2*user(user_k+1)
77*d8606c27SBarry Smith      PetscCallA(PetscOptionsGetReal(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-k1', user(user_k+2),flg,ierr))
78*d8606c27SBarry Smith      user(user_s+1) = 0.0
79*d8606c27SBarry Smith      PetscCallA(PetscOptionsGetReal(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-s0',user(user_s+1),flg,ierr))
80*d8606c27SBarry Smith      user(user_s+2) = 1.0
81*d8606c27SBarry Smith      PetscCallA(PetscOptionsGetReal(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-s1',user(user_s+2),flg,ierr))
82*d8606c27SBarry Smith
83*d8606c27SBarry Smith! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
84*d8606c27SBarry Smith!    Create timestepping solver context
85*d8606c27SBarry Smith! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
86*d8606c27SBarry Smith      PetscCallA(TSCreate(PETSC_COMM_WORLD,ts,ierr))
87*d8606c27SBarry Smith      PetscCallA(TSSetDM(ts,da,ierr))
88*d8606c27SBarry Smith      PetscCallA(TSSetType(ts,TSARKIMEX,ierr))
89*d8606c27SBarry Smith      PetscCallA(TSSetRHSFunction(ts,PETSC_NULL_VEC,FormRHSFunction,user,ierr))
90*d8606c27SBarry Smith      PetscCallA(TSSetIFunction(ts,PETSC_NULL_VEC,FormIFunction,user,ierr))
91*d8606c27SBarry Smith      PetscCallA(DMSetMatType(da,MATAIJ,ierr))
92*d8606c27SBarry Smith      PetscCallA(DMCreateMatrix(da,J,ierr))
93*d8606c27SBarry Smith      PetscCallA(TSSetIJacobian(ts,J,J,FormIJacobian,user,ierr))
94*d8606c27SBarry Smith
95*d8606c27SBarry Smith      PetscCallA(TSGetSNES(ts,snes,ierr))
96*d8606c27SBarry Smith      PetscCallA(SNESGetLineSearch(snes,linesearch,ierr))
97*d8606c27SBarry Smith      PetscCallA(SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC,ierr))
98*d8606c27SBarry Smith
99*d8606c27SBarry Smith      ftime = 1.0
100*d8606c27SBarry Smith      PetscCallA(TSSetMaxTime(ts,ftime,ierr))
101*d8606c27SBarry Smith      PetscCallA(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER,ierr))
102*d8606c27SBarry Smith
103*d8606c27SBarry Smith! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
104*d8606c27SBarry Smith!  Set initial conditions
105*d8606c27SBarry Smith! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
106*d8606c27SBarry Smith      PetscCallA(FormInitialSolution(ts,X,user,ierr))
107*d8606c27SBarry Smith      PetscCallA(TSSetSolution(ts,X,ierr))
108*d8606c27SBarry Smith      PetscCallA(VecGetSize(X,mx,ierr))
109*d8606c27SBarry Smith!  Advective CFL, I don't know why it needs so much safety factor.
110*d8606c27SBarry Smith      dt = pone * max(user(user_a+1),user(user_a+2)) / mx;
111*d8606c27SBarry Smith      PetscCallA(TSSetTimeStep(ts,dt,ierr))
112*d8606c27SBarry Smith
113*d8606c27SBarry Smith! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
114*d8606c27SBarry Smith!   Set runtime options
115*d8606c27SBarry Smith! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
116*d8606c27SBarry Smith      PetscCallA(TSSetFromOptions(ts,ierr))
117*d8606c27SBarry Smith
118*d8606c27SBarry Smith! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
119*d8606c27SBarry Smith!  Solve nonlinear system
120*d8606c27SBarry Smith! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
121*d8606c27SBarry Smith      PetscCallA(TSSolve(ts,X,ierr))
122*d8606c27SBarry Smith
123*d8606c27SBarry Smith! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
124*d8606c27SBarry Smith!  Free work space.
125*d8606c27SBarry Smith! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
126*d8606c27SBarry Smith      PetscCallA(MatDestroy(J,ierr))
127*d8606c27SBarry Smith      PetscCallA(VecDestroy(X,ierr))
128*d8606c27SBarry Smith      PetscCallA(TSDestroy(ts,ierr))
129*d8606c27SBarry Smith      PetscCallA(DMDestroy(da,ierr))
130*d8606c27SBarry Smith      PetscCallA(PetscFinalize(ierr))
131*d8606c27SBarry Smith      end program
132*d8606c27SBarry Smith
133*d8606c27SBarry Smith! Small helper to extract the layout, result uses 1-based indexing.
134*d8606c27SBarry Smith      subroutine GetLayout(da,mx,xs,xe,gxs,gxe,ierr)
135*d8606c27SBarry Smith      use petscdmda
136*d8606c27SBarry Smith      implicit none
137*d8606c27SBarry Smith
138*d8606c27SBarry Smith      DM da
139*d8606c27SBarry Smith      PetscInt mx,xs,xe,gxs,gxe
140*d8606c27SBarry Smith      PetscErrorCode ierr
141*d8606c27SBarry Smith      PetscInt xm,gxm
142*d8606c27SBarry Smith      PetscCall(DMDAGetInfo(da,PETSC_NULL_INTEGER,mx,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,ierr))
143*d8606c27SBarry Smith      PetscCall(DMDAGetCorners(da,xs,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,xm,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,ierr))
144*d8606c27SBarry Smith      PetscCall(DMDAGetGhostCorners(da,gxs,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,gxm,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,ierr))
145*d8606c27SBarry Smith      xs = xs + 1
146*d8606c27SBarry Smith      gxs = gxs + 1
147*d8606c27SBarry Smith      xe = xs + xm - 1
148*d8606c27SBarry Smith      gxe = gxs + gxm - 1
149*d8606c27SBarry Smith      end subroutine
150*d8606c27SBarry Smith
151*d8606c27SBarry Smith      subroutine FormIFunctionLocal(mx,xs,xe,gxs,gxe,x,xdot,f,a,k,s,ierr)
152*d8606c27SBarry Smith      implicit none
153*d8606c27SBarry Smith      PetscInt mx,xs,xe,gxs,gxe
154*d8606c27SBarry Smith      PetscScalar x(2,xs:xe)
155*d8606c27SBarry Smith      PetscScalar xdot(2,xs:xe)
156*d8606c27SBarry Smith      PetscScalar f(2,xs:xe)
157*d8606c27SBarry Smith      PetscReal a(2),k(2),s(2)
158*d8606c27SBarry Smith      PetscErrorCode ierr
159*d8606c27SBarry Smith      PetscInt i
160*d8606c27SBarry Smith      do 10 i = xs,xe
161*d8606c27SBarry Smith         f(1,i) = xdot(1,i) + k(1)*x(1,i) - k(2)*x(2,i) - s(1)
162*d8606c27SBarry Smith         f(2,i) = xdot(2,i) - k(1)*x(1,i) + k(2)*x(2,i) - s(2)
163*d8606c27SBarry Smith 10   continue
164*d8606c27SBarry Smith      end subroutine
165*d8606c27SBarry Smith
166*d8606c27SBarry Smith      subroutine FormIFunction(ts,t,X,Xdot,F,user,ierr)
167*d8606c27SBarry Smith      use petscts
168*d8606c27SBarry Smith      implicit none
169*d8606c27SBarry Smith
170*d8606c27SBarry Smith      TS ts
171*d8606c27SBarry Smith      PetscReal t
172*d8606c27SBarry Smith      Vec X,Xdot,F
173*d8606c27SBarry Smith      PetscReal user(6)
174*d8606c27SBarry Smith      PetscErrorCode ierr
175*d8606c27SBarry Smith      integer user_a,user_k,user_s
176*d8606c27SBarry Smith      parameter (user_a = 1,user_k = 3,user_s = 5)
177*d8606c27SBarry Smith
178*d8606c27SBarry Smith      DM             da
179*d8606c27SBarry Smith      PetscInt       mx,xs,xe,gxs,gxe
180*d8606c27SBarry Smith      PetscOffset    ixx,ixxdot,iff
181*d8606c27SBarry Smith      PetscScalar    xx(0:1),xxdot(0:1),ff(0:1)
182*d8606c27SBarry Smith
183*d8606c27SBarry Smith      PetscCall(TSGetDM(ts,da,ierr))
184*d8606c27SBarry Smith      PetscCall(GetLayout(da,mx,xs,xe,gxs,gxe,ierr))
185*d8606c27SBarry Smith
186*d8606c27SBarry Smith! Get access to vector data
187*d8606c27SBarry Smith      PetscCall(VecGetArrayRead(X,xx,ixx,ierr))
188*d8606c27SBarry Smith      PetscCall(VecGetArrayRead(Xdot,xxdot,ixxdot,ierr))
189*d8606c27SBarry Smith      PetscCall(VecGetArray(F,ff,iff,ierr))
190*d8606c27SBarry Smith
191*d8606c27SBarry Smith      PetscCall(FormIFunctionLocal(mx,xs,xe,gxs,gxe,xx(ixx),xxdot(ixxdot),ff(iff),user(user_a),user(user_k),user(user_s),ierr))
192*d8606c27SBarry Smith
193*d8606c27SBarry Smith      PetscCall(VecRestoreArrayRead(X,xx,ixx,ierr))
194*d8606c27SBarry Smith      PetscCall(VecRestoreArrayRead(Xdot,xxdot,ixxdot,ierr))
195*d8606c27SBarry Smith      PetscCall(VecRestoreArray(F,ff,iff,ierr))
196*d8606c27SBarry Smith      end subroutine
197*d8606c27SBarry Smith
198*d8606c27SBarry Smith      subroutine FormRHSFunctionLocal(mx,xs,xe,gxs,gxe,t,x,f,a,k,s,ierr)
199*d8606c27SBarry Smith      implicit none
200*d8606c27SBarry Smith      PetscInt mx,xs,xe,gxs,gxe
201*d8606c27SBarry Smith      PetscReal t
202*d8606c27SBarry Smith      PetscScalar x(2,gxs:gxe),f(2,xs:xe)
203*d8606c27SBarry Smith      PetscReal a(2),k(2),s(2)
204*d8606c27SBarry Smith      PetscErrorCode ierr
205*d8606c27SBarry Smith      PetscInt i,j
206*d8606c27SBarry Smith      PetscReal hx,u0t(2)
207*d8606c27SBarry Smith      PetscReal one,two,three,four,six,twelve
208*d8606c27SBarry Smith      PetscReal half,third,twothird,sixth
209*d8606c27SBarry Smith      PetscReal twelfth
210*d8606c27SBarry Smith
211*d8606c27SBarry Smith      one = 1.0
212*d8606c27SBarry Smith      two = 2.0
213*d8606c27SBarry Smith      three = 3.0
214*d8606c27SBarry Smith      four = 4.0
215*d8606c27SBarry Smith      six = 6.0
216*d8606c27SBarry Smith      twelve = 12.0
217*d8606c27SBarry Smith      hx = one / mx
218*d8606c27SBarry Smith!     The Fortran standard only allows positive base for power functions; Nag compiler fails on this
219*d8606c27SBarry Smith      u0t(1) = one - abs(sin(twelve*t))**four
220*d8606c27SBarry Smith      u0t(2) = 0.0
221*d8606c27SBarry Smith      half = one/two
222*d8606c27SBarry Smith      third = one / three
223*d8606c27SBarry Smith      twothird = two / three
224*d8606c27SBarry Smith      sixth = one / six
225*d8606c27SBarry Smith      twelfth = one / twelve
226*d8606c27SBarry Smith      do 20 i = xs,xe
227*d8606c27SBarry Smith         do 10 j = 1,2
228*d8606c27SBarry Smith            if (i .eq. 1) then
229*d8606c27SBarry Smith               f(j,i) = a(j)/hx*(third*u0t(j) + half*x(j,i) - x(j,i+1)  &
230*d8606c27SBarry Smith     &              + sixth*x(j,i+2))
231*d8606c27SBarry Smith            else if (i .eq. 2) then
232*d8606c27SBarry Smith               f(j,i) = a(j)/hx*(-twelfth*u0t(j) + twothird*x(j,i-1)    &
233*d8606c27SBarry Smith     &              - twothird*x(j,i+1) + twelfth*x(j,i+2))
234*d8606c27SBarry Smith            else if (i .eq. mx-1) then
235*d8606c27SBarry Smith               f(j,i) = a(j)/hx*(-sixth*x(j,i-2) + x(j,i-1)             &
236*d8606c27SBarry Smith     &         - half*x(j,i) -third*x(j,i+1))
237*d8606c27SBarry Smith            else if (i .eq. mx) then
238*d8606c27SBarry Smith               f(j,i) = a(j)/hx*(-x(j,i) + x(j,i-1))
239*d8606c27SBarry Smith            else
240*d8606c27SBarry Smith               f(j,i) = a(j)/hx*(-twelfth*x(j,i-2)                      &
241*d8606c27SBarry Smith     &              + twothird*x(j,i-1)                                 &
242*d8606c27SBarry Smith     &              - twothird*x(j,i+1) + twelfth*x(j,i+2))
243*d8606c27SBarry Smith            end if
244*d8606c27SBarry Smith 10      continue
245*d8606c27SBarry Smith 20   continue
246*d8606c27SBarry Smith      end subroutine
247*d8606c27SBarry Smith
248*d8606c27SBarry Smith      subroutine FormRHSFunction(ts,t,X,F,user,ierr)
249*d8606c27SBarry Smith      use petscts
250*d8606c27SBarry Smith      implicit none
251*d8606c27SBarry Smith
252*d8606c27SBarry Smith      TS ts
253*d8606c27SBarry Smith      PetscReal t
254*d8606c27SBarry Smith      Vec X,F
255*d8606c27SBarry Smith      PetscReal user(6)
256*d8606c27SBarry Smith      PetscErrorCode ierr
257*d8606c27SBarry Smith      integer user_a,user_k,user_s
258*d8606c27SBarry Smith      parameter (user_a = 1,user_k = 3,user_s = 5)
259*d8606c27SBarry Smith      DM             da
260*d8606c27SBarry Smith      Vec            Xloc
261*d8606c27SBarry Smith      PetscInt       mx,xs,xe,gxs,gxe
262*d8606c27SBarry Smith      PetscOffset    ixx,iff
263*d8606c27SBarry Smith      PetscScalar    xx(0:1),ff(0:1)
264*d8606c27SBarry Smith
265*d8606c27SBarry Smith      PetscCall(TSGetDM(ts,da,ierr))
266*d8606c27SBarry Smith      PetscCall(GetLayout(da,mx,xs,xe,gxs,gxe,ierr))
267*d8606c27SBarry Smith
268*d8606c27SBarry Smith!     Scatter ghost points to local vector,using the 2-step process
269*d8606c27SBarry Smith!        DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
270*d8606c27SBarry Smith!     By placing code between these two statements, computations can be
271*d8606c27SBarry Smith!     done while messages are in transition.
272*d8606c27SBarry Smith      PetscCall(DMGetLocalVector(da,Xloc,ierr))
273*d8606c27SBarry Smith      PetscCall(DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc,ierr))
274*d8606c27SBarry Smith      PetscCall(DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc,ierr))
275*d8606c27SBarry Smith
276*d8606c27SBarry Smith! Get access to vector data
277*d8606c27SBarry Smith      PetscCall(VecGetArrayRead(Xloc,xx,ixx,ierr))
278*d8606c27SBarry Smith      PetscCall(VecGetArray(F,ff,iff,ierr))
279*d8606c27SBarry Smith
280*d8606c27SBarry Smith      PetscCall(FormRHSFunctionLocal(mx,xs,xe,gxs,gxe,t,xx(ixx),ff(iff),user(user_a),user(user_k),user(user_s),ierr))
281*d8606c27SBarry Smith
282*d8606c27SBarry Smith      PetscCall(VecRestoreArrayRead(Xloc,xx,ixx,ierr))
283*d8606c27SBarry Smith      PetscCall(VecRestoreArray(F,ff,iff,ierr))
284*d8606c27SBarry Smith      PetscCall(DMRestoreLocalVector(da,Xloc,ierr))
285*d8606c27SBarry Smith      end subroutine
286*d8606c27SBarry Smith
287*d8606c27SBarry Smith! ---------------------------------------------------------------------
288*d8606c27SBarry Smith!
289*d8606c27SBarry Smith!  IJacobian - Compute IJacobian = dF/dU + shift*dF/dUdot
290*d8606c27SBarry Smith!
291*d8606c27SBarry Smith      subroutine FormIJacobian(ts,t,X,Xdot,shift,J,Jpre,user,ierr)
292*d8606c27SBarry Smith      use petscts
293*d8606c27SBarry Smith      implicit none
294*d8606c27SBarry Smith
295*d8606c27SBarry Smith      TS ts
296*d8606c27SBarry Smith      PetscReal t,shift
297*d8606c27SBarry Smith      Vec X,Xdot
298*d8606c27SBarry Smith      Mat J,Jpre
299*d8606c27SBarry Smith      PetscReal user(6)
300*d8606c27SBarry Smith      PetscErrorCode ierr
301*d8606c27SBarry Smith      integer user_a,user_k,user_s
302*d8606c27SBarry Smith      parameter (user_a = 0,user_k = 2,user_s = 4)
303*d8606c27SBarry Smith
304*d8606c27SBarry Smith      DM             da
305*d8606c27SBarry Smith      PetscInt       mx,xs,xe,gxs,gxe
306*d8606c27SBarry Smith      PetscInt       i,i1,row,col
307*d8606c27SBarry Smith      PetscReal      k1,k2;
308*d8606c27SBarry Smith      PetscScalar    val(4)
309*d8606c27SBarry Smith
310*d8606c27SBarry Smith      PetscCall(TSGetDM(ts,da,ierr))
311*d8606c27SBarry Smith      PetscCall(GetLayout(da,mx,xs,xe,gxs,gxe,ierr))
312*d8606c27SBarry Smith
313*d8606c27SBarry Smith      i1 = 1
314*d8606c27SBarry Smith      k1 = user(user_k+1)
315*d8606c27SBarry Smith      k2 = user(user_k+2)
316*d8606c27SBarry Smith      do 10 i = xs,xe
317*d8606c27SBarry Smith         row = i-gxs
318*d8606c27SBarry Smith         col = i-gxs
319*d8606c27SBarry Smith         val(1) = shift + k1
320*d8606c27SBarry Smith         val(2) = -k2
321*d8606c27SBarry Smith         val(3) = -k1
322*d8606c27SBarry Smith         val(4) = shift + k2
323*d8606c27SBarry Smith         PetscCall(MatSetValuesBlockedLocal(Jpre,i1,row,i1,col,val,INSERT_VALUES,ierr))
324*d8606c27SBarry Smith 10   continue
325*d8606c27SBarry Smith      PetscCall(MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY,ierr))
326*d8606c27SBarry Smith      PetscCall(MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY,ierr))
327*d8606c27SBarry Smith      if (J /= Jpre) then
328*d8606c27SBarry Smith         PetscCall(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY,ierr))
329*d8606c27SBarry Smith         PetscCall(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY,ierr))
330*d8606c27SBarry Smith      end if
331*d8606c27SBarry Smith      end subroutine
332*d8606c27SBarry Smith
333*d8606c27SBarry Smith      subroutine FormInitialSolutionLocal(mx,xs,xe,gxs,gxe,x,a,k,s,ierr)
334*d8606c27SBarry Smith      implicit none
335*d8606c27SBarry Smith      PetscInt mx,xs,xe,gxs,gxe
336*d8606c27SBarry Smith      PetscScalar x(2,xs:xe)
337*d8606c27SBarry Smith      PetscReal a(2),k(2),s(2)
338*d8606c27SBarry Smith      PetscErrorCode ierr
339*d8606c27SBarry Smith
340*d8606c27SBarry Smith      PetscInt i
341*d8606c27SBarry Smith      PetscReal one,hx,r,ik
342*d8606c27SBarry Smith      one = 1.0
343*d8606c27SBarry Smith      hx = one / mx
344*d8606c27SBarry Smith      do 10 i=xs,xe
345*d8606c27SBarry Smith         r = i*hx
346*d8606c27SBarry Smith         if (k(2) .ne. 0.0) then
347*d8606c27SBarry Smith            ik = one/k(2)
348*d8606c27SBarry Smith         else
349*d8606c27SBarry Smith            ik = one
350*d8606c27SBarry Smith         end if
351*d8606c27SBarry Smith         x(1,i) = one + s(2)*r
352*d8606c27SBarry Smith         x(2,i) = k(1)*ik*x(1,i) + s(2)*ik
353*d8606c27SBarry Smith 10   continue
354*d8606c27SBarry Smith      end subroutine
355*d8606c27SBarry Smith
356*d8606c27SBarry Smith      subroutine FormInitialSolution(ts,X,user,ierr)
357*d8606c27SBarry Smith      use petscts
358*d8606c27SBarry Smith      implicit none
359*d8606c27SBarry Smith
360*d8606c27SBarry Smith      TS ts
361*d8606c27SBarry Smith      Vec X
362*d8606c27SBarry Smith      PetscReal user(6)
363*d8606c27SBarry Smith      PetscErrorCode ierr
364*d8606c27SBarry Smith      integer user_a,user_k,user_s
365*d8606c27SBarry Smith      parameter (user_a = 1,user_k = 3,user_s = 5)
366*d8606c27SBarry Smith
367*d8606c27SBarry Smith      DM             da
368*d8606c27SBarry Smith      PetscInt       mx,xs,xe,gxs,gxe
369*d8606c27SBarry Smith      PetscOffset    ixx
370*d8606c27SBarry Smith      PetscScalar    xx(0:1)
371*d8606c27SBarry Smith
372*d8606c27SBarry Smith      PetscCall(TSGetDM(ts,da,ierr))
373*d8606c27SBarry Smith      PetscCall(GetLayout(da,mx,xs,xe,gxs,gxe,ierr))
374*d8606c27SBarry Smith
375*d8606c27SBarry Smith! Get access to vector data
376*d8606c27SBarry Smith      PetscCall(VecGetArray(X,xx,ixx,ierr))
377*d8606c27SBarry Smith
378*d8606c27SBarry Smith      PetscCall(FormInitialSolutionLocal(mx,xs,xe,gxs,gxe,xx(ixx),user(user_a),user(user_k),user(user_s),ierr))
379*d8606c27SBarry Smith
380*d8606c27SBarry Smith      PetscCall(VecRestoreArray(X,xx,ixx,ierr))
381*d8606c27SBarry Smith      end subroutine
382*d8606c27SBarry Smith
383*d8606c27SBarry Smith!/*TEST
384*d8606c27SBarry Smith!
385*d8606c27SBarry Smith!    test:
386*d8606c27SBarry Smith!      args: -da_grid_x 200 -ts_arkimex_type 4
387*d8606c27SBarry Smith!      requires: !single
388*d8606c27SBarry Smith!
389*d8606c27SBarry Smith!TEST*/
390