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