Actual source code: ex2.c

petsc-3.7.2 2016-06-05
Report Typos and Errors
  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 =
 12: */


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

 26: typedef struct {
 27:   PetscScalar k1,k2,k3;
 28:   PetscScalar sigma2;
 29:   Vec         initialsolution;
 30: } AppCtx;

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

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

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

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

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

 82:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 83:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
 84:   if (A != B) {
 85:     MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
 86:     MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
 87:   }
 88:   return(0);
 89: }

 93: static PetscErrorCode Solution(TS ts,PetscReal t,Vec U,AppCtx *ctx)
 94: {

 98:   VecCopy(ctx->initialsolution,U);
 99:   if (t > 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Solution not given");
100:   return(0);
101: }

105: int main(int argc,char **argv)
106: {
107:   TS             ts;            /* ODE integrator */
108:   Vec            U;             /* solution */
109:   Mat            A;             /* Jacobian matrix */
111:   PetscMPIInt    size;
112:   PetscInt       n = 4;
113:   AppCtx         ctx;
114:   PetscScalar    *u;

116:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
117:      Initialize program
118:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
119:   PetscInitialize(&argc,&argv,(char*)0,help);
120:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
121:   if (size > 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Only for sequential runs");

123:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
124:     Create necessary matrix and vectors
125:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
126:   MatCreate(PETSC_COMM_WORLD,&A);
127:   MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE);
128:   MatSetFromOptions(A);
129:   MatSetUp(A);

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

133:   ctx.k1     = 1.0e-5;
134:   ctx.k2     = 1.0e5;
135:   ctx.k3     = 1.0e-16;
136:   ctx.sigma2 = 1.0e6;

138:   VecDuplicate(U,&ctx.initialsolution);
139:   VecGetArray(ctx.initialsolution,&u);
140:   u[0] = 0.0;
141:   u[1] = 1.3e8;
142:   u[2] = 5.0e11;
143:   u[3] = 8.0e11;
144:   VecRestoreArray(ctx.initialsolution,&u);

146:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
147:      Create timestepping solver context
148:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
149:   TSCreate(PETSC_COMM_WORLD,&ts);
150:   TSSetProblemType(ts,TS_NONLINEAR);
151:   TSSetType(ts,TSROSW);
152:   TSSetIFunction(ts,NULL,(TSIFunction) IFunction,&ctx);
153:   TSSetIJacobian(ts,A,A,(TSIJacobian)IJacobian,&ctx);

155:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
156:      Set initial conditions
157:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
158:   Solution(ts,0,U,&ctx);
159:   TSSetInitialTimeStep(ts,4.0*3600,1.0);
160:   TSSetSolution(ts,U);

162:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
163:      Set solver options
164:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
165:   TSSetDuration(ts,1000000,518400.0);
166:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);
167:   TSSetMaxStepRejections(ts,100);
168:   TSSetMaxSNESFailures(ts,-1); /* unlimited */
169:   TSSetFromOptions(ts);

171:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
172:      Solve nonlinear system
173:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
174:   TSSolve(ts,U);

176:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
177:      Free work space.  All PETSc objects should be destroyed when they
178:      are no longer needed.
179:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
180:   VecDestroy(&ctx.initialsolution);
181:   MatDestroy(&A);
182:   VecDestroy(&U);
183:   TSDestroy(&ts);

185:   PetscFinalize();
186:   return(0);
187: }