13a2a065fSHong Zhang static char help[] = "A fast-slow system for testing ARKIMEX.\n";
23a2a065fSHong Zhang
33a2a065fSHong Zhang /*F
43a2a065fSHong Zhang
53a2a065fSHong Zhang \begin{eqnarray}
63a2a065fSHong Zhang ys' = -2.0*\frac{-1.0+ys^2.0-\cos(t)}{2.0*ys}+0.05*\frac{-2.0+yf^2-\cos(5.0*t)}{2.0*yf}-\frac{\sin(t)}{2.0*ys}\\
73a2a065fSHong Zhang yf' = 0.05*\frac{-1.0+ys^2-\cos(t)}{2.0*ys}-\frac{-2.0+yf^2-\cos(5.0*t)}{2.0*yf}-5.0*\frac{\sin(5.0*t)}{2.0*yf}\\
83a2a065fSHong Zhang \end{eqnarray}
93a2a065fSHong Zhang
103a2a065fSHong Zhang F*/
113a2a065fSHong Zhang
123a2a065fSHong Zhang /*
133a2a065fSHong Zhang This example demonstrates how to use ARKIMEX for solving a fast-slow system. The system is partitioned additively and component-wise at the same time.
143a2a065fSHong Zhang ys stands for the slow component and yf stands for the fast component. On the RHS for yf, only the term -\frac{-2.0+yf^2-\cos(5.0*t)}{2.0*yf} is treated implicitly while the rest is treated explicitly.
153a2a065fSHong Zhang */
163a2a065fSHong Zhang
173a2a065fSHong Zhang #include <petscts.h>
183a2a065fSHong Zhang
193a2a065fSHong Zhang typedef struct {
203a2a065fSHong Zhang PetscReal Tf, dt;
213a2a065fSHong Zhang } AppCtx;
223a2a065fSHong Zhang
RHSFunctionslow(TS ts,PetscReal t,Vec U,Vec F,AppCtx * ctx)233a2a065fSHong Zhang static PetscErrorCode RHSFunctionslow(TS ts, PetscReal t, Vec U, Vec F, AppCtx *ctx)
243a2a065fSHong Zhang {
253a2a065fSHong Zhang const PetscScalar *u;
263a2a065fSHong Zhang PetscScalar *f;
273a2a065fSHong Zhang
283a2a065fSHong Zhang PetscFunctionBegin;
293a2a065fSHong Zhang PetscCall(VecGetArrayRead(U, &u));
303a2a065fSHong Zhang PetscCall(VecGetArray(F, &f));
313a2a065fSHong Zhang f[0] = -2.0 * (-1.0 + u[0] * u[0] - PetscCosScalar(t)) / (2.0 * u[0]) + 0.05 * (-2.0 + u[1] * u[1] - PetscCosScalar(5.0 * t)) / (2.0 * u[1]) - PetscSinScalar(t) / (2.0 * u[0]);
323a2a065fSHong Zhang PetscCall(VecRestoreArrayRead(U, &u));
333a2a065fSHong Zhang PetscCall(VecRestoreArray(F, &f));
343a2a065fSHong Zhang PetscFunctionReturn(PETSC_SUCCESS);
353a2a065fSHong Zhang }
363a2a065fSHong Zhang
RHSFunctionfast(TS ts,PetscReal t,Vec U,Vec F,AppCtx * ctx)373a2a065fSHong Zhang static PetscErrorCode RHSFunctionfast(TS ts, PetscReal t, Vec U, Vec F, AppCtx *ctx)
383a2a065fSHong Zhang {
393a2a065fSHong Zhang const PetscScalar *u;
403a2a065fSHong Zhang PetscScalar *f;
413a2a065fSHong Zhang
423a2a065fSHong Zhang PetscFunctionBegin;
433a2a065fSHong Zhang PetscCall(VecGetArrayRead(U, &u));
443a2a065fSHong Zhang PetscCall(VecGetArray(F, &f));
453a2a065fSHong Zhang f[0] = 0.05 * (-1.0 + u[0] * u[0] - PetscCosScalar(t)) / (2.0 * u[0]) - 5.0 * PetscSinScalar(5.0 * t) / (2.0 * u[1]);
463a2a065fSHong Zhang PetscCall(VecRestoreArrayRead(U, &u));
473a2a065fSHong Zhang PetscCall(VecRestoreArray(F, &f));
483a2a065fSHong Zhang PetscFunctionReturn(PETSC_SUCCESS);
493a2a065fSHong Zhang }
503a2a065fSHong Zhang
IFunctionfast(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,AppCtx * ctx)513a2a065fSHong Zhang static PetscErrorCode IFunctionfast(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, AppCtx *ctx)
523a2a065fSHong Zhang {
533a2a065fSHong Zhang const PetscScalar *u, *udot;
543a2a065fSHong Zhang PetscScalar *f;
553a2a065fSHong Zhang
563a2a065fSHong Zhang PetscFunctionBegin;
573a2a065fSHong Zhang PetscCall(VecGetArrayRead(U, &u));
583a2a065fSHong Zhang PetscCall(VecGetArrayRead(Udot, &udot));
593a2a065fSHong Zhang PetscCall(VecGetArray(F, &f));
603a2a065fSHong Zhang f[0] = udot[0] + (-2.0 + u[1] * u[1] - PetscCosScalar(5.0 * t)) / (2.0 * u[1]);
613a2a065fSHong Zhang PetscCall(VecRestoreArrayRead(Udot, &udot));
623a2a065fSHong Zhang PetscCall(VecRestoreArrayRead(U, &u));
633a2a065fSHong Zhang PetscCall(VecRestoreArray(F, &f));
643a2a065fSHong Zhang PetscFunctionReturn(PETSC_SUCCESS);
653a2a065fSHong Zhang }
663a2a065fSHong Zhang
IJacobianfast(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,AppCtx * ctx)673a2a065fSHong Zhang static PetscErrorCode IJacobianfast(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal a, Mat A, Mat B, AppCtx *ctx)
683a2a065fSHong Zhang {
693a2a065fSHong Zhang PetscInt rowcol[] = {0};
703a2a065fSHong Zhang const PetscScalar *u;
713a2a065fSHong Zhang PetscScalar J[1][1];
723a2a065fSHong Zhang
733a2a065fSHong Zhang PetscFunctionBeginUser;
743a2a065fSHong Zhang PetscCall(VecGetArrayRead(U, &u));
753a2a065fSHong Zhang J[0][0] = a + (2.0 + PetscCosScalar(5.0 * t)) / (2.0 * u[1] * u[1]) + 0.5;
763a2a065fSHong Zhang PetscCall(MatSetValues(B, 1, rowcol, 1, rowcol, &J[0][0], INSERT_VALUES));
773a2a065fSHong Zhang PetscCall(VecRestoreArrayRead(U, &u));
783a2a065fSHong Zhang
793a2a065fSHong Zhang PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
803a2a065fSHong Zhang PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
813a2a065fSHong Zhang if (A != B) {
823a2a065fSHong Zhang PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
833a2a065fSHong Zhang PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
843a2a065fSHong Zhang }
853a2a065fSHong Zhang PetscFunctionReturn(PETSC_SUCCESS);
863a2a065fSHong Zhang }
873a2a065fSHong Zhang
883a2a065fSHong Zhang /*
893a2a065fSHong Zhang Define the analytic solution for check method easily
903a2a065fSHong Zhang */
sol_true(PetscReal t,Vec U)913a2a065fSHong Zhang static PetscErrorCode sol_true(PetscReal t, Vec U)
923a2a065fSHong Zhang {
933a2a065fSHong Zhang PetscScalar *u;
943a2a065fSHong Zhang
953a2a065fSHong Zhang PetscFunctionBegin;
963a2a065fSHong Zhang PetscCall(VecGetArray(U, &u));
973a2a065fSHong Zhang u[0] = PetscSqrtScalar(1.0 + PetscCosScalar(t));
983a2a065fSHong Zhang u[1] = PetscSqrtScalar(2.0 + PetscCosScalar(5.0 * t));
993a2a065fSHong Zhang PetscCall(VecRestoreArray(U, &u));
1003a2a065fSHong Zhang PetscFunctionReturn(PETSC_SUCCESS);
1013a2a065fSHong Zhang }
1023a2a065fSHong Zhang
main(int argc,char ** argv)1033a2a065fSHong Zhang int main(int argc, char **argv)
1043a2a065fSHong Zhang {
1053a2a065fSHong Zhang TS ts; /* ODE integrator */
1063a2a065fSHong Zhang Vec U; /* solution will be stored here */
1073a2a065fSHong Zhang Vec Utrue;
1083a2a065fSHong Zhang Mat A;
1093a2a065fSHong Zhang PetscMPIInt size;
1103a2a065fSHong Zhang AppCtx ctx;
1113a2a065fSHong Zhang PetscScalar *u;
1123a2a065fSHong Zhang IS iss;
1133a2a065fSHong Zhang IS isf;
1143a2a065fSHong Zhang PetscInt *indicess;
1153a2a065fSHong Zhang PetscInt *indicesf;
1163a2a065fSHong Zhang PetscInt n = 2;
1173a2a065fSHong Zhang PetscScalar error, tt;
1183a2a065fSHong Zhang
1193a2a065fSHong Zhang /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1203a2a065fSHong Zhang Initialize program
1213a2a065fSHong Zhang - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1223a2a065fSHong Zhang PetscFunctionBeginUser;
123c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help));
1243a2a065fSHong Zhang PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
1253a2a065fSHong Zhang PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only for sequential runs");
1263a2a065fSHong Zhang
1273a2a065fSHong Zhang /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1283a2a065fSHong Zhang Create index for slow part and fast part
1293a2a065fSHong Zhang - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1303a2a065fSHong Zhang PetscCall(PetscMalloc1(1, &indicess));
1313a2a065fSHong Zhang indicess[0] = 0;
1323a2a065fSHong Zhang PetscCall(PetscMalloc1(1, &indicesf));
1333a2a065fSHong Zhang indicesf[0] = 1;
1343a2a065fSHong Zhang PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 1, indicess, PETSC_COPY_VALUES, &iss));
1353a2a065fSHong Zhang PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 1, indicesf, PETSC_COPY_VALUES, &isf));
1363a2a065fSHong Zhang
1373a2a065fSHong Zhang /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1383a2a065fSHong Zhang Create necessary vector
1393a2a065fSHong Zhang - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1403a2a065fSHong Zhang PetscCall(VecCreate(PETSC_COMM_WORLD, &U));
1413a2a065fSHong Zhang PetscCall(VecSetSizes(U, n, PETSC_DETERMINE));
1423a2a065fSHong Zhang PetscCall(VecSetFromOptions(U));
1433a2a065fSHong Zhang PetscCall(VecDuplicate(U, &Utrue));
1443a2a065fSHong Zhang PetscCall(VecCopy(U, Utrue));
1453a2a065fSHong Zhang
1463a2a065fSHong Zhang /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1473a2a065fSHong Zhang Set runtime options
1483a2a065fSHong Zhang - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1493a2a065fSHong Zhang PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "ODE options", "");
1503a2a065fSHong Zhang ctx.Tf = 0.3;
1513a2a065fSHong Zhang ctx.dt = 0.01;
1523a2a065fSHong Zhang PetscCall(PetscOptionsScalar("-Tf", "", "", ctx.Tf, &ctx.Tf, NULL));
1533a2a065fSHong Zhang PetscCall(PetscOptionsScalar("-dt", "", "", ctx.dt, &ctx.dt, NULL));
1543a2a065fSHong Zhang PetscOptionsEnd();
1553a2a065fSHong Zhang
1563a2a065fSHong Zhang /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1573a2a065fSHong Zhang Initialize the solution
1583a2a065fSHong Zhang - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1593a2a065fSHong Zhang PetscCall(VecGetArray(U, &u));
1603a2a065fSHong Zhang u[0] = PetscSqrtScalar(2.0);
1613a2a065fSHong Zhang u[1] = PetscSqrtScalar(3.0);
1623a2a065fSHong Zhang PetscCall(VecRestoreArray(U, &u));
1633a2a065fSHong Zhang
1643a2a065fSHong Zhang /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1653a2a065fSHong Zhang Create timestepping solver context
1663a2a065fSHong Zhang - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1673a2a065fSHong Zhang PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
1683a2a065fSHong Zhang PetscCall(TSSetType(ts, TSARKIMEX));
1693a2a065fSHong Zhang PetscCall(TSARKIMEXSetFastSlowSplit(ts, PETSC_TRUE));
1703a2a065fSHong Zhang
1713a2a065fSHong Zhang PetscCall(TSRHSSplitSetIS(ts, "slow", iss));
1723a2a065fSHong Zhang PetscCall(TSRHSSplitSetIS(ts, "fast", isf));
1733a2a065fSHong Zhang PetscCall(TSRHSSplitSetRHSFunction(ts, "slow", NULL, (TSRHSFunctionFn *)RHSFunctionslow, &ctx));
1743a2a065fSHong Zhang PetscCall(TSRHSSplitSetRHSFunction(ts, "fast", NULL, (TSRHSFunctionFn *)RHSFunctionfast, &ctx));
1753a2a065fSHong Zhang PetscCall(TSRHSSplitSetIFunction(ts, "fast", NULL, (TSIFunctionFn *)IFunctionfast, &ctx));
1763a2a065fSHong Zhang PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
1773a2a065fSHong Zhang PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, 1, 1));
1783a2a065fSHong Zhang PetscCall(MatSetFromOptions(A));
1793a2a065fSHong Zhang PetscCall(MatSetUp(A));
1803a2a065fSHong Zhang PetscCall(TSRHSSplitSetIJacobian(ts, "fast", A, A, (TSIJacobianFn *)IJacobianfast, &ctx));
1813a2a065fSHong Zhang
1823a2a065fSHong Zhang /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1833a2a065fSHong Zhang Set initial conditions
1843a2a065fSHong Zhang - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1853a2a065fSHong Zhang PetscCall(TSSetSolution(ts, U));
1863a2a065fSHong Zhang
1873a2a065fSHong Zhang /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1883a2a065fSHong Zhang Set solver options
1893a2a065fSHong Zhang - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1903a2a065fSHong Zhang PetscCall(TSSetMaxTime(ts, ctx.Tf));
1913a2a065fSHong Zhang PetscCall(TSSetTimeStep(ts, ctx.dt));
1923a2a065fSHong Zhang PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
1933a2a065fSHong Zhang PetscCall(TSSetFromOptions(ts));
1943a2a065fSHong Zhang
1953a2a065fSHong Zhang /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1963a2a065fSHong Zhang Solve linear system
1973a2a065fSHong Zhang - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1983a2a065fSHong Zhang PetscCall(TSSolve(ts, U));
1993a2a065fSHong Zhang PetscCall(VecView(U, PETSC_VIEWER_STDOUT_WORLD));
2003a2a065fSHong Zhang
2013a2a065fSHong Zhang /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
202*f0b74427SPierre Jolivet Check the error of the PETSc solution
2033a2a065fSHong Zhang - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2043a2a065fSHong Zhang PetscCall(TSGetTime(ts, &tt));
2053a2a065fSHong Zhang PetscCall(sol_true(tt, Utrue));
2063a2a065fSHong Zhang PetscCall(VecAXPY(Utrue, -1.0, U));
2073a2a065fSHong Zhang PetscCall(VecNorm(Utrue, NORM_2, &error));
2083a2a065fSHong Zhang
2093a2a065fSHong Zhang /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
2103a2a065fSHong Zhang Print norm2 error
2113a2a065fSHong Zhang - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2123a2a065fSHong Zhang PetscCall(PetscPrintf(PETSC_COMM_WORLD, "l2 error norm: %g\n", (double)PetscAbsScalar(error)));
2133a2a065fSHong Zhang
2143a2a065fSHong Zhang /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
2153a2a065fSHong Zhang Free work space. All PETSc objects should be destroyed when they are no longer needed.
2163a2a065fSHong Zhang - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2173a2a065fSHong Zhang PetscCall(VecDestroy(&U));
2183a2a065fSHong Zhang PetscCall(TSDestroy(&ts));
2193a2a065fSHong Zhang PetscCall(VecDestroy(&Utrue));
2203a2a065fSHong Zhang PetscCall(ISDestroy(&iss));
2213a2a065fSHong Zhang PetscCall(ISDestroy(&isf));
2223a2a065fSHong Zhang PetscCall(PetscFree(indicess));
2233a2a065fSHong Zhang PetscCall(PetscFree(indicesf));
2243a2a065fSHong Zhang PetscCall(MatDestroy(&A));
2253a2a065fSHong Zhang PetscCall(PetscFinalize());
2263a2a065fSHong Zhang return 0;
2273a2a065fSHong Zhang }
2283a2a065fSHong Zhang
2293a2a065fSHong Zhang /*TEST
2303a2a065fSHong Zhang build:
2313a2a065fSHong Zhang requires: !complex
2323a2a065fSHong Zhang
2333a2a065fSHong Zhang test:
2343a2a065fSHong Zhang
2353a2a065fSHong Zhang test:
2363a2a065fSHong Zhang suffix: 2
2373a2a065fSHong Zhang args: -ts_arkimex_initial_guess_extrapolate 1
2383a2a065fSHong Zhang output_file: output/ex3fastslowsplit_1.out
2393a2a065fSHong Zhang
2403a2a065fSHong Zhang TEST*/
241