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