xref: /petsc/src/ts/tutorials/advection-diffusion-reaction/ex2.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Reaction Equation from Chemistry\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown /*
5c4762a1bSJed Brown 
6c4762a1bSJed Brown      Page 6, An example from Atomospheric Chemistry
7c4762a1bSJed Brown 
8c4762a1bSJed Brown                  u_1_t =
9c4762a1bSJed Brown                  u_2_t =
10c4762a1bSJed Brown                  u_3_t =
11c4762a1bSJed Brown                  u_4_t =
12c4762a1bSJed Brown 
13c4762a1bSJed Brown   -ts_monitor_lg_error -ts_monitor_lg_solution  -ts_view -ts_max_time 2.e4
14c4762a1bSJed Brown 
15c4762a1bSJed Brown */
16c4762a1bSJed Brown 
17c4762a1bSJed Brown /*
18c4762a1bSJed Brown    Include "petscts.h" so that we can use TS solvers.  Note that this
19c4762a1bSJed Brown    file automatically includes:
20c4762a1bSJed Brown      petscsys.h       - base PETSc routines   petscvec.h - vectors
21c4762a1bSJed Brown      petscmat.h - matrices
22c4762a1bSJed Brown      petscis.h     - index sets            petscksp.h - Krylov subspace methods
23c4762a1bSJed Brown      petscviewer.h - viewers               petscpc.h  - preconditioners
24c4762a1bSJed Brown      petscksp.h   - linear solvers
25c4762a1bSJed Brown */
26c4762a1bSJed Brown 
27c4762a1bSJed Brown #include <petscts.h>
28c4762a1bSJed Brown 
29c4762a1bSJed Brown typedef struct {
30c4762a1bSJed Brown   PetscScalar k1, k2, k3;
31c4762a1bSJed Brown   PetscScalar sigma2;
32c4762a1bSJed Brown   Vec         initialsolution;
33c4762a1bSJed Brown } AppCtx;
34c4762a1bSJed Brown 
35*9371c9d4SSatish Balay PetscScalar k1(AppCtx *ctx, PetscReal t) {
36c4762a1bSJed Brown   PetscReal th    = t / 3600.0;
37c4762a1bSJed Brown   PetscReal barth = th - 24.0 * PetscFloorReal(th / 24.0);
38c4762a1bSJed Brown   if (((((PetscInt)th) % 24) < 4) || ((((PetscInt)th) % 24) >= 20)) return (1.0e-40);
39c4762a1bSJed Brown   else return (ctx->k1 * PetscExpReal(7.0 * PetscPowReal(PetscSinReal(.0625 * PETSC_PI * (barth - 4.0)), .2)));
40c4762a1bSJed Brown }
41c4762a1bSJed Brown 
42*9371c9d4SSatish Balay static PetscErrorCode IFunction(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, AppCtx *ctx) {
43c4762a1bSJed Brown   PetscScalar       *f;
44c4762a1bSJed Brown   const PetscScalar *u, *udot;
45c4762a1bSJed Brown 
46c4762a1bSJed Brown   PetscFunctionBegin;
479566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
489566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Udot, &udot));
499566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(F, &f));
50c4762a1bSJed Brown   f[0] = udot[0] - k1(ctx, t) * u[2] + ctx->k2 * u[0];
51c4762a1bSJed Brown   f[1] = udot[1] - k1(ctx, t) * u[2] + ctx->k3 * u[1] * u[3] - ctx->sigma2;
52c4762a1bSJed Brown   f[2] = udot[2] - ctx->k3 * u[1] * u[3] + k1(ctx, t) * u[2];
53c4762a1bSJed Brown   f[3] = udot[3] - ctx->k2 * u[0] + ctx->k3 * u[1] * u[3];
549566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
559566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Udot, &udot));
569566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(F, &f));
57c4762a1bSJed Brown   PetscFunctionReturn(0);
58c4762a1bSJed Brown }
59c4762a1bSJed Brown 
60*9371c9d4SSatish Balay static PetscErrorCode IJacobian(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal a, Mat A, Mat B, AppCtx *ctx) {
61c4762a1bSJed Brown   PetscInt           rowcol[] = {0, 1, 2, 3};
62c4762a1bSJed Brown   PetscScalar        J[4][4];
63c4762a1bSJed Brown   const PetscScalar *u, *udot;
64c4762a1bSJed Brown 
65c4762a1bSJed Brown   PetscFunctionBegin;
669566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
679566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Udot, &udot));
68*9371c9d4SSatish Balay   J[0][0] = a + ctx->k2;
69*9371c9d4SSatish Balay   J[0][1] = 0.0;
70*9371c9d4SSatish Balay   J[0][2] = -k1(ctx, t);
71*9371c9d4SSatish Balay   J[0][3] = 0.0;
72*9371c9d4SSatish Balay   J[1][0] = 0.0;
73*9371c9d4SSatish Balay   J[1][1] = a + ctx->k3 * u[3];
74*9371c9d4SSatish Balay   J[1][2] = -k1(ctx, t);
75*9371c9d4SSatish Balay   J[1][3] = ctx->k3 * u[1];
76*9371c9d4SSatish Balay   J[2][0] = 0.0;
77*9371c9d4SSatish Balay   J[2][1] = -ctx->k3 * u[3];
78*9371c9d4SSatish Balay   J[2][2] = a + k1(ctx, t);
79*9371c9d4SSatish Balay   J[2][3] = -ctx->k3 * u[1];
80*9371c9d4SSatish Balay   J[3][0] = -ctx->k2;
81*9371c9d4SSatish Balay   J[3][1] = ctx->k3 * u[3];
82*9371c9d4SSatish Balay   J[3][2] = 0.0;
83*9371c9d4SSatish Balay   J[3][3] = a + ctx->k3 * u[1];
849566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B, 4, rowcol, 4, rowcol, &J[0][0], INSERT_VALUES));
859566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
869566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Udot, &udot));
87c4762a1bSJed Brown 
889566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
899566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
90c4762a1bSJed Brown   if (A != B) {
919566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
929566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
93c4762a1bSJed Brown   }
94c4762a1bSJed Brown   PetscFunctionReturn(0);
95c4762a1bSJed Brown }
96c4762a1bSJed Brown 
97*9371c9d4SSatish Balay static PetscErrorCode Solution(TS ts, PetscReal t, Vec U, AppCtx *ctx) {
98c4762a1bSJed Brown   PetscFunctionBegin;
999566063dSJacob Faibussowitsch   PetscCall(VecCopy(ctx->initialsolution, U));
1003c633725SBarry Smith   PetscCheck(t <= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Solution not given");
101c4762a1bSJed Brown   PetscFunctionReturn(0);
102c4762a1bSJed Brown }
103c4762a1bSJed Brown 
104*9371c9d4SSatish Balay int main(int argc, char **argv) {
105c4762a1bSJed Brown   TS           ts; /* ODE integrator */
106c4762a1bSJed Brown   Vec          U;  /* solution */
107c4762a1bSJed Brown   Mat          A;  /* Jacobian matrix */
108c4762a1bSJed Brown   PetscMPIInt  size;
109c4762a1bSJed Brown   PetscInt     n = 4;
110c4762a1bSJed Brown   AppCtx       ctx;
111c4762a1bSJed Brown   PetscScalar *u;
112c4762a1bSJed Brown 
113c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
114c4762a1bSJed Brown      Initialize program
115c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
116327415f7SBarry Smith   PetscFunctionBeginUser;
1179566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
1189566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
1193c633725SBarry Smith   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only for sequential runs");
120c4762a1bSJed Brown 
121c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
122c4762a1bSJed Brown     Create necessary matrix and vectors
123c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1249566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
1259566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, n, n, PETSC_DETERMINE, PETSC_DETERMINE));
1269566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
1279566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
128c4762a1bSJed Brown 
1299566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A, &U, NULL));
130c4762a1bSJed Brown 
131c4762a1bSJed Brown   ctx.k1     = 1.0e-5;
132c4762a1bSJed Brown   ctx.k2     = 1.0e5;
133c4762a1bSJed Brown   ctx.k3     = 1.0e-16;
134c4762a1bSJed Brown   ctx.sigma2 = 1.0e6;
135c4762a1bSJed Brown 
1369566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(U, &ctx.initialsolution));
1379566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(ctx.initialsolution, &u));
138c4762a1bSJed Brown   u[0] = 0.0;
139c4762a1bSJed Brown   u[1] = 1.3e8;
140c4762a1bSJed Brown   u[2] = 5.0e11;
141c4762a1bSJed Brown   u[3] = 8.0e11;
1429566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(ctx.initialsolution, &u));
143c4762a1bSJed Brown 
144c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
145c4762a1bSJed Brown      Create timestepping solver context
146c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1479566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
1489566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
1499566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts, TSROSW));
1509566063dSJacob Faibussowitsch   PetscCall(TSSetIFunction(ts, NULL, (TSIFunction)IFunction, &ctx));
1519566063dSJacob Faibussowitsch   PetscCall(TSSetIJacobian(ts, A, A, (TSIJacobian)IJacobian, &ctx));
152c4762a1bSJed Brown 
153c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
154c4762a1bSJed Brown      Set initial conditions
155c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1569566063dSJacob Faibussowitsch   PetscCall(Solution(ts, 0, U, &ctx));
1579566063dSJacob Faibussowitsch   PetscCall(TSSetTime(ts, 4.0 * 3600));
1589566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, 1.0));
1599566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts, U));
160c4762a1bSJed Brown 
161c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
162c4762a1bSJed Brown      Set solver options
163c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1649566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, 518400.0));
1659566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
1669566063dSJacob Faibussowitsch   PetscCall(TSSetMaxStepRejections(ts, 100));
1679566063dSJacob Faibussowitsch   PetscCall(TSSetMaxSNESFailures(ts, -1)); /* unlimited */
1689566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
169c4762a1bSJed Brown 
170c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
171c4762a1bSJed Brown      Solve nonlinear system
172c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1739566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, U));
174c4762a1bSJed Brown 
175c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
176c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
177c4762a1bSJed Brown      are no longer needed.
178c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1799566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ctx.initialsolution));
1809566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
1819566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&U));
1829566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
183c4762a1bSJed Brown 
1849566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
185b122ec5aSJacob Faibussowitsch   return 0;
186c4762a1bSJed Brown }
187c4762a1bSJed Brown 
188c4762a1bSJed Brown /*TEST
189c4762a1bSJed Brown 
190c4762a1bSJed Brown    test:
191c4762a1bSJed Brown      args: -ts_view -ts_max_time 2.e4
192c4762a1bSJed Brown      timeoutfactor: 15
193c4762a1bSJed Brown      requires: !single
194c4762a1bSJed Brown 
195c4762a1bSJed Brown TEST*/
196