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