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