xref: /petsc/src/ts/tutorials/advection-diffusion-reaction/ex2.c (revision 9566063d113dddea24716c546802770db7481bc0)
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 
35c4762a1bSJed Brown PetscScalar k1(AppCtx *ctx,PetscReal t)
36c4762a1bSJed Brown {
37c4762a1bSJed Brown   PetscReal th    = t/3600.0;
38c4762a1bSJed Brown   PetscReal barth = th - 24.0*PetscFloorReal(th/24.0);
39c4762a1bSJed Brown   if (((((PetscInt)th) % 24) < 4) || ((((PetscInt)th) % 24) >= 20)) return(1.0e-40);
40c4762a1bSJed Brown   else return(ctx->k1*PetscExpReal(7.0*PetscPowReal(PetscSinReal(.0625*PETSC_PI*(barth - 4.0)),.2)));
41c4762a1bSJed Brown }
42c4762a1bSJed Brown 
43c4762a1bSJed Brown static PetscErrorCode IFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,AppCtx *ctx)
44c4762a1bSJed Brown {
45c4762a1bSJed Brown   PetscScalar       *f;
46c4762a1bSJed Brown   const PetscScalar *u,*udot;
47c4762a1bSJed Brown 
48c4762a1bSJed Brown   PetscFunctionBegin;
49*9566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U,&u));
50*9566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Udot,&udot));
51*9566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(F,&f));
52c4762a1bSJed Brown   f[0] = udot[0] - k1(ctx,t)*u[2] + ctx->k2*u[0];
53c4762a1bSJed Brown   f[1] = udot[1] - k1(ctx,t)*u[2] + ctx->k3*u[1]*u[3] - ctx->sigma2;
54c4762a1bSJed Brown   f[2] = udot[2] - ctx->k3*u[1]*u[3] + k1(ctx,t)*u[2];
55c4762a1bSJed Brown   f[3] = udot[3] - ctx->k2*u[0] + ctx->k3*u[1]*u[3];
56*9566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U,&u));
57*9566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Udot,&udot));
58*9566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(F,&f));
59c4762a1bSJed Brown   PetscFunctionReturn(0);
60c4762a1bSJed Brown }
61c4762a1bSJed Brown 
62c4762a1bSJed Brown static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,AppCtx *ctx)
63c4762a1bSJed Brown {
64c4762a1bSJed Brown   PetscInt          rowcol[] = {0,1,2,3};
65c4762a1bSJed Brown   PetscScalar       J[4][4];
66c4762a1bSJed Brown   const PetscScalar *u,*udot;
67c4762a1bSJed Brown 
68c4762a1bSJed Brown   PetscFunctionBegin;
69*9566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U,&u));
70*9566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Udot,&udot));
71c4762a1bSJed Brown   J[0][0] = a + ctx->k2;   J[0][1] = 0.0;                J[0][2] = -k1(ctx,t);       J[0][3] = 0.0;
72c4762a1bSJed Brown   J[1][0] = 0.0;           J[1][1] = a + ctx->k3*u[3];   J[1][2] = -k1(ctx,t);       J[1][3] = ctx->k3*u[1];
73c4762a1bSJed Brown   J[2][0] = 0.0;           J[2][1] = -ctx->k3*u[3];      J[2][2] = a + k1(ctx,t);    J[2][3] =  -ctx->k3*u[1];
74c4762a1bSJed Brown   J[3][0] =  -ctx->k2;     J[3][1] = ctx->k3*u[3];       J[3][2] = 0.0;              J[3][3] = a + ctx->k3*u[1];
75*9566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B,4,rowcol,4,rowcol,&J[0][0],INSERT_VALUES));
76*9566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U,&u));
77*9566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Udot,&udot));
78c4762a1bSJed Brown 
79*9566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
80*9566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
81c4762a1bSJed Brown   if (A != B) {
82*9566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
83*9566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
84c4762a1bSJed Brown   }
85c4762a1bSJed Brown   PetscFunctionReturn(0);
86c4762a1bSJed Brown }
87c4762a1bSJed Brown 
88c4762a1bSJed Brown static PetscErrorCode Solution(TS ts,PetscReal t,Vec U,AppCtx *ctx)
89c4762a1bSJed Brown {
90c4762a1bSJed Brown   PetscFunctionBegin;
91*9566063dSJacob Faibussowitsch   PetscCall(VecCopy(ctx->initialsolution,U));
923c633725SBarry Smith   PetscCheck(t <= 0,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Solution not given");
93c4762a1bSJed Brown   PetscFunctionReturn(0);
94c4762a1bSJed Brown }
95c4762a1bSJed Brown 
96c4762a1bSJed Brown int main(int argc,char **argv)
97c4762a1bSJed Brown {
98c4762a1bSJed Brown   TS             ts;            /* ODE integrator */
99c4762a1bSJed Brown   Vec            U;             /* solution */
100c4762a1bSJed Brown   Mat            A;             /* Jacobian matrix */
101c4762a1bSJed Brown   PetscMPIInt    size;
102c4762a1bSJed Brown   PetscInt       n = 4;
103c4762a1bSJed Brown   AppCtx         ctx;
104c4762a1bSJed Brown   PetscScalar    *u;
105c4762a1bSJed Brown 
106c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
107c4762a1bSJed Brown      Initialize program
108c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
109*9566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
110*9566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
1113c633725SBarry Smith   PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Only for sequential runs");
112c4762a1bSJed Brown 
113c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
114c4762a1bSJed Brown     Create necessary matrix and vectors
115c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
116*9566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
117*9566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE));
118*9566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
119*9566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
120c4762a1bSJed Brown 
121*9566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A,&U,NULL));
122c4762a1bSJed Brown 
123c4762a1bSJed Brown   ctx.k1     = 1.0e-5;
124c4762a1bSJed Brown   ctx.k2     = 1.0e5;
125c4762a1bSJed Brown   ctx.k3     = 1.0e-16;
126c4762a1bSJed Brown   ctx.sigma2 = 1.0e6;
127c4762a1bSJed Brown 
128*9566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(U,&ctx.initialsolution));
129*9566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(ctx.initialsolution,&u));
130c4762a1bSJed Brown   u[0] = 0.0;
131c4762a1bSJed Brown   u[1] = 1.3e8;
132c4762a1bSJed Brown   u[2] = 5.0e11;
133c4762a1bSJed Brown   u[3] = 8.0e11;
134*9566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(ctx.initialsolution,&u));
135c4762a1bSJed Brown 
136c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
137c4762a1bSJed Brown      Create timestepping solver context
138c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
139*9566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD,&ts));
140*9566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(ts,TS_NONLINEAR));
141*9566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts,TSROSW));
142*9566063dSJacob Faibussowitsch   PetscCall(TSSetIFunction(ts,NULL,(TSIFunction) IFunction,&ctx));
143*9566063dSJacob Faibussowitsch   PetscCall(TSSetIJacobian(ts,A,A,(TSIJacobian)IJacobian,&ctx));
144c4762a1bSJed Brown 
145c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
146c4762a1bSJed Brown      Set initial conditions
147c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
148*9566063dSJacob Faibussowitsch   PetscCall(Solution(ts,0,U,&ctx));
149*9566063dSJacob Faibussowitsch   PetscCall(TSSetTime(ts,4.0*3600));
150*9566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts,1.0));
151*9566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts,U));
152c4762a1bSJed Brown 
153c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
154c4762a1bSJed Brown      Set solver options
155c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
156*9566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts,518400.0));
157*9566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER));
158*9566063dSJacob Faibussowitsch   PetscCall(TSSetMaxStepRejections(ts,100));
159*9566063dSJacob Faibussowitsch   PetscCall(TSSetMaxSNESFailures(ts,-1)); /* unlimited */
160*9566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
161c4762a1bSJed Brown 
162c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
163c4762a1bSJed Brown      Solve nonlinear system
164c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
165*9566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts,U));
166c4762a1bSJed Brown 
167c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
168c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
169c4762a1bSJed Brown      are no longer needed.
170c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
171*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ctx.initialsolution));
172*9566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
173*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&U));
174*9566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
175c4762a1bSJed Brown 
176*9566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
177b122ec5aSJacob Faibussowitsch   return 0;
178c4762a1bSJed Brown }
179c4762a1bSJed Brown 
180c4762a1bSJed Brown /*TEST
181c4762a1bSJed Brown 
182c4762a1bSJed Brown    test:
183c4762a1bSJed Brown      args: -ts_view -ts_max_time 2.e4
184c4762a1bSJed Brown      timeoutfactor: 15
185c4762a1bSJed Brown      requires: !single
186c4762a1bSJed Brown 
187c4762a1bSJed Brown TEST*/
188