Actual source code: ex2.c


  2: static char help[] = "Reaction Equation from Chemistry\n";

  4: /*

  6:      Page 6, An example from Atomospheric Chemistry

  8:                  u_1_t =
  9:                  u_2_t =
 10:                  u_3_t =
 11:                  u_4_t =

 13:   -ts_monitor_lg_error -ts_monitor_lg_solution  -ts_view -ts_max_time 2.e4

 15: */

 17: /*
 18:    Include "petscts.h" so that we can use TS solvers.  Note that this
 19:    file automatically includes:
 20:      petscsys.h       - base PETSc routines   petscvec.h - vectors
 21:      petscmat.h - matrices
 22:      petscis.h     - index sets            petscksp.h - Krylov subspace methods
 23:      petscviewer.h - viewers               petscpc.h  - preconditioners
 24:      petscksp.h   - linear solvers
 25: */

 27: #include <petscts.h>

 29: typedef struct {
 30:   PetscScalar k1,k2,k3;
 31:   PetscScalar sigma2;
 32:   Vec         initialsolution;
 33: } AppCtx;

 35: PetscScalar k1(AppCtx *ctx,PetscReal t)
 36: {
 37:   PetscReal th    = t/3600.0;
 38:   PetscReal barth = th - 24.0*PetscFloorReal(th/24.0);
 39:   if (((((PetscInt)th) % 24) < 4) || ((((PetscInt)th) % 24) >= 20)) return(1.0e-40);
 40:   else return(ctx->k1*PetscExpReal(7.0*PetscPowReal(PetscSinReal(.0625*PETSC_PI*(barth - 4.0)),.2)));
 41: }

 43: static PetscErrorCode IFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,AppCtx *ctx)
 44: {
 45:   PetscScalar       *f;
 46:   const PetscScalar *u,*udot;

 48:   VecGetArrayRead(U,&u);
 49:   VecGetArrayRead(Udot,&udot);
 50:   VecGetArrayWrite(F,&f);
 51:   f[0] = udot[0] - k1(ctx,t)*u[2] + ctx->k2*u[0];
 52:   f[1] = udot[1] - k1(ctx,t)*u[2] + ctx->k3*u[1]*u[3] - ctx->sigma2;
 53:   f[2] = udot[2] - ctx->k3*u[1]*u[3] + k1(ctx,t)*u[2];
 54:   f[3] = udot[3] - ctx->k2*u[0] + ctx->k3*u[1]*u[3];
 55:   VecRestoreArrayRead(U,&u);
 56:   VecRestoreArrayRead(Udot,&udot);
 57:   VecRestoreArrayWrite(F,&f);
 58:   return 0;
 59: }

 61: static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,AppCtx *ctx)
 62: {
 63:   PetscInt          rowcol[] = {0,1,2,3};
 64:   PetscScalar       J[4][4];
 65:   const PetscScalar *u,*udot;

 67:   VecGetArrayRead(U,&u);
 68:   VecGetArrayRead(Udot,&udot);
 69:   J[0][0] = a + ctx->k2;   J[0][1] = 0.0;                J[0][2] = -k1(ctx,t);       J[0][3] = 0.0;
 70:   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];
 71:   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];
 72:   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];
 73:   MatSetValues(B,4,rowcol,4,rowcol,&J[0][0],INSERT_VALUES);
 74:   VecRestoreArrayRead(U,&u);
 75:   VecRestoreArrayRead(Udot,&udot);

 77:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 78:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
 79:   if (A != B) {
 80:     MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
 81:     MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
 82:   }
 83:   return 0;
 84: }

 86: static PetscErrorCode Solution(TS ts,PetscReal t,Vec U,AppCtx *ctx)
 87: {
 88:   VecCopy(ctx->initialsolution,U);
 90:   return 0;
 91: }

 93: int main(int argc,char **argv)
 94: {
 95:   TS             ts;            /* ODE integrator */
 96:   Vec            U;             /* solution */
 97:   Mat            A;             /* Jacobian matrix */
 98:   PetscMPIInt    size;
 99:   PetscInt       n = 4;
100:   AppCtx         ctx;
101:   PetscScalar    *u;

103:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
104:      Initialize program
105:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
106:   PetscInitialize(&argc,&argv,(char*)0,help);
107:   MPI_Comm_size(PETSC_COMM_WORLD,&size);

110:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
111:     Create necessary matrix and vectors
112:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
113:   MatCreate(PETSC_COMM_WORLD,&A);
114:   MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE);
115:   MatSetFromOptions(A);
116:   MatSetUp(A);

118:   MatCreateVecs(A,&U,NULL);

120:   ctx.k1     = 1.0e-5;
121:   ctx.k2     = 1.0e5;
122:   ctx.k3     = 1.0e-16;
123:   ctx.sigma2 = 1.0e6;

125:   VecDuplicate(U,&ctx.initialsolution);
126:   VecGetArrayWrite(ctx.initialsolution,&u);
127:   u[0] = 0.0;
128:   u[1] = 1.3e8;
129:   u[2] = 5.0e11;
130:   u[3] = 8.0e11;
131:   VecRestoreArrayWrite(ctx.initialsolution,&u);

133:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
134:      Create timestepping solver context
135:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
136:   TSCreate(PETSC_COMM_WORLD,&ts);
137:   TSSetProblemType(ts,TS_NONLINEAR);
138:   TSSetType(ts,TSROSW);
139:   TSSetIFunction(ts,NULL,(TSIFunction) IFunction,&ctx);
140:   TSSetIJacobian(ts,A,A,(TSIJacobian)IJacobian,&ctx);

142:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
143:      Set initial conditions
144:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
145:   Solution(ts,0,U,&ctx);
146:   TSSetTime(ts,4.0*3600);
147:   TSSetTimeStep(ts,1.0);
148:   TSSetSolution(ts,U);

150:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
151:      Set solver options
152:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
153:   TSSetMaxTime(ts,518400.0);
154:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);
155:   TSSetMaxStepRejections(ts,100);
156:   TSSetMaxSNESFailures(ts,-1); /* unlimited */
157:   TSSetFromOptions(ts);

159:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
160:      Solve nonlinear system
161:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
162:   TSSolve(ts,U);

164:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
165:      Free work space.  All PETSc objects should be destroyed when they
166:      are no longer needed.
167:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
168:   VecDestroy(&ctx.initialsolution);
169:   MatDestroy(&A);
170:   VecDestroy(&U);
171:   TSDestroy(&ts);

173:   PetscFinalize();
174:   return 0;
175: }

177: /*TEST

179:    test:
180:      args: -ts_view -ts_max_time 2.e4
181:      timeoutfactor: 15
182:      requires: !single

184: TEST*/