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