Actual source code: ex2.c

petsc-3.8.4 2018-03-24
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*PetscFloorReal(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: }

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

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

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

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

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

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

 92:   VecCopy(ctx->initialsolution,U);
 93:   if (t > 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Solution not given");
 94:   return(0);
 95: }

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

108:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
109:      Initialize program
110:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
111:   PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
112:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
113:   if (size > 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Only for sequential runs");

115:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
116:     Create necessary matrix and vectors
117:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
118:   MatCreate(PETSC_COMM_WORLD,&A);
119:   MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE);
120:   MatSetFromOptions(A);
121:   MatSetUp(A);

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

125:   ctx.k1     = 1.0e-5;
126:   ctx.k2     = 1.0e5;
127:   ctx.k3     = 1.0e-16;
128:   ctx.sigma2 = 1.0e6;

130:   VecDuplicate(U,&ctx.initialsolution);
131:   VecGetArray(ctx.initialsolution,&u);
132:   u[0] = 0.0;
133:   u[1] = 1.3e8;
134:   u[2] = 5.0e11;
135:   u[3] = 8.0e11;
136:   VecRestoreArray(ctx.initialsolution,&u);

138:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
139:      Create timestepping solver context
140:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
141:   TSCreate(PETSC_COMM_WORLD,&ts);
142:   TSSetProblemType(ts,TS_NONLINEAR);
143:   TSSetType(ts,TSROSW);
144:   TSSetIFunction(ts,NULL,(TSIFunction) IFunction,&ctx);
145:   TSSetIJacobian(ts,A,A,(TSIJacobian)IJacobian,&ctx);

147:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
148:      Set initial conditions
149:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
150:   Solution(ts,0,U,&ctx);
151:   TSSetTime(ts,4.0*3600);
152:   TSSetTimeStep(ts,1.0);
153:   TSSetSolution(ts,U);

155:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
156:      Set solver options
157:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
158:   TSSetMaxTime(ts,518400.0);
159:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);
160:   TSSetMaxStepRejections(ts,100);
161:   TSSetMaxSNESFailures(ts,-1); /* unlimited */
162:   TSSetFromOptions(ts);

164:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
165:      Solve nonlinear system
166:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
167:   TSSolve(ts,U);

169:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
170:      Free work space.  All PETSc objects should be destroyed when they
171:      are no longer needed.
172:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
173:   VecDestroy(&ctx.initialsolution);
174:   MatDestroy(&A);
175:   VecDestroy(&U);
176:   TSDestroy(&ts);

178:   PetscFinalize();
179:   return ierr;
180: }