1c4762a1bSJed Brown static char help[] = "Reaction Equation from Chemistry\n";
2c4762a1bSJed Brown
3c4762a1bSJed Brown /*
4c4762a1bSJed Brown
5c4762a1bSJed Brown Page 6, An example from Atomospheric Chemistry
6c4762a1bSJed Brown
7c4762a1bSJed Brown u_1_t =
8c4762a1bSJed Brown u_2_t =
9c4762a1bSJed Brown u_3_t =
10c4762a1bSJed Brown u_4_t =
11c4762a1bSJed Brown
12c4762a1bSJed Brown -ts_monitor_lg_error -ts_monitor_lg_solution -ts_view -ts_max_time 2.e4
13c4762a1bSJed Brown
14c4762a1bSJed Brown */
15c4762a1bSJed Brown
16c4762a1bSJed Brown /*
17c4762a1bSJed Brown Include "petscts.h" so that we can use TS solvers. Note that this
18c4762a1bSJed Brown file automatically includes:
19c4762a1bSJed Brown petscsys.h - base PETSc routines petscvec.h - vectors
20c4762a1bSJed Brown petscmat.h - matrices
21c4762a1bSJed Brown petscis.h - index sets petscksp.h - Krylov subspace methods
22c4762a1bSJed Brown petscviewer.h - viewers petscpc.h - preconditioners
23c4762a1bSJed Brown petscksp.h - linear solvers
24c4762a1bSJed Brown */
25c4762a1bSJed Brown
26c4762a1bSJed Brown #include <petscts.h>
27c4762a1bSJed Brown
28c4762a1bSJed Brown typedef struct {
29c4762a1bSJed Brown PetscScalar k1, k2, k3;
30c4762a1bSJed Brown PetscScalar sigma2;
31c4762a1bSJed Brown Vec initialsolution;
32c4762a1bSJed Brown } AppCtx;
33c4762a1bSJed Brown
k1(AppCtx * ctx,PetscReal t)34d71ae5a4SJacob Faibussowitsch PetscScalar k1(AppCtx *ctx, PetscReal t)
35d71ae5a4SJacob Faibussowitsch {
36c4762a1bSJed Brown PetscReal th = t / 3600.0;
37c4762a1bSJed Brown PetscReal barth = th - 24.0 * PetscFloorReal(th / 24.0);
384ad8454bSPierre Jolivet if (((((PetscInt)th) % 24) < 4) || ((((PetscInt)th) % 24) >= 20)) return 1.0e-40;
394ad8454bSPierre Jolivet else return ctx->k1 * PetscExpReal(7.0 * PetscPowReal(PetscSinReal(.0625 * PETSC_PI * (barth - 4.0)), .2));
40c4762a1bSJed Brown }
41c4762a1bSJed Brown
IFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,AppCtx * ctx)42d71ae5a4SJacob Faibussowitsch static PetscErrorCode IFunction(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, AppCtx *ctx)
43d71ae5a4SJacob Faibussowitsch {
44c4762a1bSJed Brown PetscScalar *f;
45c4762a1bSJed Brown const PetscScalar *u, *udot;
46c4762a1bSJed Brown
47c4762a1bSJed Brown PetscFunctionBegin;
489566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u));
499566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Udot, &udot));
509566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(F, &f));
51c4762a1bSJed Brown f[0] = udot[0] - k1(ctx, t) * u[2] + ctx->k2 * u[0];
52c4762a1bSJed Brown f[1] = udot[1] - k1(ctx, t) * u[2] + ctx->k3 * u[1] * u[3] - ctx->sigma2;
53c4762a1bSJed Brown f[2] = udot[2] - ctx->k3 * u[1] * u[3] + k1(ctx, t) * u[2];
54c4762a1bSJed Brown f[3] = udot[3] - ctx->k2 * u[0] + ctx->k3 * u[1] * u[3];
559566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u));
569566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Udot, &udot));
579566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(F, &f));
583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
59c4762a1bSJed Brown }
60c4762a1bSJed Brown
IJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,AppCtx * ctx)61d71ae5a4SJacob Faibussowitsch static PetscErrorCode IJacobian(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal a, Mat A, Mat B, AppCtx *ctx)
62d71ae5a4SJacob Faibussowitsch {
63c4762a1bSJed Brown PetscInt rowcol[] = {0, 1, 2, 3};
64c4762a1bSJed Brown PetscScalar J[4][4];
65c4762a1bSJed Brown const PetscScalar *u, *udot;
66c4762a1bSJed Brown
67c4762a1bSJed Brown PetscFunctionBegin;
689566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u));
699566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Udot, &udot));
709371c9d4SSatish Balay J[0][0] = a + ctx->k2;
719371c9d4SSatish Balay J[0][1] = 0.0;
729371c9d4SSatish Balay J[0][2] = -k1(ctx, t);
739371c9d4SSatish Balay J[0][3] = 0.0;
749371c9d4SSatish Balay J[1][0] = 0.0;
759371c9d4SSatish Balay J[1][1] = a + ctx->k3 * u[3];
769371c9d4SSatish Balay J[1][2] = -k1(ctx, t);
779371c9d4SSatish Balay J[1][3] = ctx->k3 * u[1];
789371c9d4SSatish Balay J[2][0] = 0.0;
799371c9d4SSatish Balay J[2][1] = -ctx->k3 * u[3];
809371c9d4SSatish Balay J[2][2] = a + k1(ctx, t);
819371c9d4SSatish Balay J[2][3] = -ctx->k3 * u[1];
829371c9d4SSatish Balay J[3][0] = -ctx->k2;
839371c9d4SSatish Balay J[3][1] = ctx->k3 * u[3];
849371c9d4SSatish Balay J[3][2] = 0.0;
859371c9d4SSatish Balay J[3][3] = a + ctx->k3 * u[1];
869566063dSJacob Faibussowitsch PetscCall(MatSetValues(B, 4, rowcol, 4, rowcol, &J[0][0], INSERT_VALUES));
879566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u));
889566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Udot, &udot));
89c4762a1bSJed Brown
909566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
919566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
92c4762a1bSJed Brown if (A != B) {
939566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
949566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
95c4762a1bSJed Brown }
963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
97c4762a1bSJed Brown }
98c4762a1bSJed Brown
Solution(TS ts,PetscReal t,Vec U,AppCtx * ctx)99d71ae5a4SJacob Faibussowitsch static PetscErrorCode Solution(TS ts, PetscReal t, Vec U, AppCtx *ctx)
100d71ae5a4SJacob Faibussowitsch {
101c4762a1bSJed Brown PetscFunctionBegin;
1029566063dSJacob Faibussowitsch PetscCall(VecCopy(ctx->initialsolution, U));
1033c633725SBarry Smith PetscCheck(t <= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Solution not given");
1043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
105c4762a1bSJed Brown }
106c4762a1bSJed Brown
main(int argc,char ** argv)107d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
108d71ae5a4SJacob Faibussowitsch {
109c4762a1bSJed Brown TS ts; /* ODE integrator */
110c4762a1bSJed Brown Vec U; /* solution */
111c4762a1bSJed Brown Mat A; /* Jacobian matrix */
112c4762a1bSJed Brown PetscMPIInt size;
113c4762a1bSJed Brown PetscInt n = 4;
114c4762a1bSJed Brown AppCtx ctx;
115c4762a1bSJed Brown PetscScalar *u;
116c4762a1bSJed Brown
117c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
118c4762a1bSJed Brown Initialize program
119c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
120327415f7SBarry Smith PetscFunctionBeginUser;
121*c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help));
1229566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
1233c633725SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only for sequential runs");
124c4762a1bSJed Brown
125c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
126c4762a1bSJed Brown Create necessary matrix and vectors
127c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1289566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
1299566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, n, n, PETSC_DETERMINE, PETSC_DETERMINE));
1309566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A));
1319566063dSJacob Faibussowitsch PetscCall(MatSetUp(A));
132c4762a1bSJed Brown
1339566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &U, NULL));
134c4762a1bSJed Brown
135c4762a1bSJed Brown ctx.k1 = 1.0e-5;
136c4762a1bSJed Brown ctx.k2 = 1.0e5;
137c4762a1bSJed Brown ctx.k3 = 1.0e-16;
138c4762a1bSJed Brown ctx.sigma2 = 1.0e6;
139c4762a1bSJed Brown
1409566063dSJacob Faibussowitsch PetscCall(VecDuplicate(U, &ctx.initialsolution));
1419566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(ctx.initialsolution, &u));
142c4762a1bSJed Brown u[0] = 0.0;
143c4762a1bSJed Brown u[1] = 1.3e8;
144c4762a1bSJed Brown u[2] = 5.0e11;
145c4762a1bSJed Brown u[3] = 8.0e11;
1469566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(ctx.initialsolution, &u));
147c4762a1bSJed Brown
148c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
149c4762a1bSJed Brown Create timestepping solver context
150c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1519566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
1529566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
1539566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSROSW));
1548434afd1SBarry Smith PetscCall(TSSetIFunction(ts, NULL, (TSIFunctionFn *)IFunction, &ctx));
1558434afd1SBarry Smith PetscCall(TSSetIJacobian(ts, A, A, (TSIJacobianFn *)IJacobian, &ctx));
156c4762a1bSJed Brown
157c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
158c4762a1bSJed Brown Set initial conditions
159c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1609566063dSJacob Faibussowitsch PetscCall(Solution(ts, 0, U, &ctx));
1619566063dSJacob Faibussowitsch PetscCall(TSSetTime(ts, 4.0 * 3600));
1629566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, 1.0));
1639566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts, U));
164c4762a1bSJed Brown
165c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
166c4762a1bSJed Brown Set solver options
167c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1689566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, 518400.0));
1699566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
1709566063dSJacob Faibussowitsch PetscCall(TSSetMaxStepRejections(ts, 100));
1719566063dSJacob Faibussowitsch PetscCall(TSSetMaxSNESFailures(ts, -1)); /* unlimited */
1729566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts));
173c4762a1bSJed Brown
174c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
175c4762a1bSJed Brown Solve nonlinear system
176c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1779566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, U));
178c4762a1bSJed Brown
179c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
180c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they
181c4762a1bSJed Brown are no longer needed.
182c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1839566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ctx.initialsolution));
1849566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A));
1859566063dSJacob Faibussowitsch PetscCall(VecDestroy(&U));
1869566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts));
187c4762a1bSJed Brown
1889566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
189b122ec5aSJacob Faibussowitsch return 0;
190c4762a1bSJed Brown }
191c4762a1bSJed Brown
192c4762a1bSJed Brown /*TEST
193c4762a1bSJed Brown
194c4762a1bSJed Brown test:
195c4762a1bSJed Brown args: -ts_view -ts_max_time 2.e4
196c4762a1bSJed Brown timeoutfactor: 15
197c4762a1bSJed Brown requires: !single
198c4762a1bSJed Brown
199c4762a1bSJed Brown TEST*/
200