xref: /petsc/src/ts/tutorials/advection-diffusion-reaction/ex2.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(U,&u));
50*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(Udot,&udot));
51*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(U,&u));
57*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(Udot,&udot));
58*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(U,&u));
70*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(B,4,rowcol,4,rowcol,&J[0][0],INSERT_VALUES));
76*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(U,&u));
77*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(Udot,&udot));
78c4762a1bSJed Brown 
79*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
80*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
81c4762a1bSJed Brown   if (A != B) {
82*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
83*5f80ce2aSJacob Faibussowitsch     CHKERRQ(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*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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   PetscErrorCode ierr;
102c4762a1bSJed Brown   PetscMPIInt    size;
103c4762a1bSJed Brown   PetscInt       n = 4;
104c4762a1bSJed Brown   AppCtx         ctx;
105c4762a1bSJed Brown   PetscScalar    *u;
106c4762a1bSJed Brown 
107c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
108c4762a1bSJed Brown      Initialize program
109c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
110c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
111*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
1123c633725SBarry Smith   PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Only for sequential runs");
113c4762a1bSJed Brown 
114c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
115c4762a1bSJed Brown     Create necessary matrix and vectors
116c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
117*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A));
118*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE));
119*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(A));
120*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(A));
121c4762a1bSJed Brown 
122*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(A,&U,NULL));
123c4762a1bSJed Brown 
124c4762a1bSJed Brown   ctx.k1     = 1.0e-5;
125c4762a1bSJed Brown   ctx.k2     = 1.0e5;
126c4762a1bSJed Brown   ctx.k3     = 1.0e-16;
127c4762a1bSJed Brown   ctx.sigma2 = 1.0e6;
128c4762a1bSJed Brown 
129*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(U,&ctx.initialsolution));
130*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayWrite(ctx.initialsolution,&u));
131c4762a1bSJed Brown   u[0] = 0.0;
132c4762a1bSJed Brown   u[1] = 1.3e8;
133c4762a1bSJed Brown   u[2] = 5.0e11;
134c4762a1bSJed Brown   u[3] = 8.0e11;
135*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayWrite(ctx.initialsolution,&u));
136c4762a1bSJed Brown 
137c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
138c4762a1bSJed Brown      Create timestepping solver context
139c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
140*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts));
141*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetProblemType(ts,TS_NONLINEAR));
142*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetType(ts,TSROSW));
143*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetIFunction(ts,NULL,(TSIFunction) IFunction,&ctx));
144*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetIJacobian(ts,A,A,(TSIJacobian)IJacobian,&ctx));
145c4762a1bSJed Brown 
146c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
147c4762a1bSJed Brown      Set initial conditions
148c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
149*5f80ce2aSJacob Faibussowitsch   CHKERRQ(Solution(ts,0,U,&ctx));
150*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetTime(ts,4.0*3600));
151*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetTimeStep(ts,1.0));
152*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetSolution(ts,U));
153c4762a1bSJed Brown 
154c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
155c4762a1bSJed Brown      Set solver options
156c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
157*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetMaxTime(ts,518400.0));
158*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER));
159*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetMaxStepRejections(ts,100));
160*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetMaxSNESFailures(ts,-1)); /* unlimited */
161*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetFromOptions(ts));
162c4762a1bSJed Brown 
163c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
164c4762a1bSJed Brown      Solve nonlinear system
165c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
166*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSolve(ts,U));
167c4762a1bSJed Brown 
168c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
169c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
170c4762a1bSJed Brown      are no longer needed.
171c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
172*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&ctx.initialsolution));
173*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
174*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&U));
175*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSDestroy(&ts));
176c4762a1bSJed Brown 
177c4762a1bSJed Brown   ierr = PetscFinalize();
178c4762a1bSJed Brown   return ierr;
179c4762a1bSJed Brown }
180c4762a1bSJed Brown 
181c4762a1bSJed Brown /*TEST
182c4762a1bSJed Brown 
183c4762a1bSJed Brown    test:
184c4762a1bSJed Brown      args: -ts_view -ts_max_time 2.e4
185c4762a1bSJed Brown      timeoutfactor: 15
186c4762a1bSJed Brown      requires: !single
187c4762a1bSJed Brown 
188c4762a1bSJed Brown TEST*/
189