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