Actual source code: ex2.c
petsc-3.4.2 2013-07-02
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*PetscExpScalar(7.0*PetscPowScalar(PetscSinScalar(.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: {
45: PetscScalar *u,*udot,*f;
48: VecGetArray(U,&u);
49: VecGetArray(Udot,&udot);
50: VecGetArray(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: VecRestoreArray(U,&u);
56: VecRestoreArray(Udot,&udot);
57: VecRestoreArray(F,&f);
58: return(0);
59: }
63: static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat *A,Mat *B,MatStructure *flag,AppCtx *ctx)
64: {
66: PetscInt rowcol[] = {0,1,2,3};
67: PetscScalar *u,*udot,J[4][4];
70: VecGetArray(U,&u);
71: VecGetArray(Udot,&udot);
72: J[0][0] = a + ctx->k2; J[0][1] = 0.0; J[0][2] = -k1(ctx,t); J[0][3] = 0.0;
73: 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];
74: 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];
75: 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];
76: MatSetValues(*B,4,rowcol,4,rowcol,&J[0][0],INSERT_VALUES);
77: VecRestoreArray(U,&u);
78: VecRestoreArray(Udot,&udot);
80: MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);
81: MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);
82: if (*A != *B) {
83: MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);
84: MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);
85: }
86: *flag = SAME_NONZERO_PATTERN;
87: return(0);
88: }
92: static PetscErrorCode Solution(TS ts,PetscReal t,Vec U,AppCtx *ctx)
93: {
97: VecCopy(ctx->initialsolution,U);
98: if (t > 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Solution not given");
99: return(0);
100: }
104: int main(int argc,char **argv)
105: {
106: TS ts; /* ODE integrator */
107: Vec U; /* solution */
108: Mat A; /* Jacobian matrix */
110: PetscMPIInt size;
111: PetscInt n = 4;
112: AppCtx ctx;
113: PetscScalar *u;
115: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
116: Initialize program
117: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
118: PetscInitialize(&argc,&argv,(char*)0,help);
119: MPI_Comm_size(PETSC_COMM_WORLD,&size);
120: if (size > 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Only for sequential runs");
122: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
123: Create necessary matrix and vectors
124: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
125: MatCreate(PETSC_COMM_WORLD,&A);
126: MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE);
127: MatSetFromOptions(A);
128: MatSetUp(A);
130: MatGetVecs(A,&U,NULL);
132: ctx.k1 = 1.0e-5;
133: ctx.k2 = 1.0e5;
134: ctx.k3 = 1.0e-16;
135: ctx.sigma2 = 1.0e6;
137: VecDuplicate(U,&ctx.initialsolution);
138: VecGetArray(ctx.initialsolution,&u);
139: u[0] = 0.0;
140: u[1] = 1.3e8;
141: u[2] = 5.0e11;
142: u[3] = 8.0e11;
143: VecRestoreArray(ctx.initialsolution,&u);
145: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
146: Create timestepping solver context
147: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
148: TSCreate(PETSC_COMM_WORLD,&ts);
149: TSSetProblemType(ts,TS_NONLINEAR);
150: TSSetType(ts,TSROSW);
151: TSSetIFunction(ts,NULL,(TSIFunction) IFunction,&ctx);
152: TSSetIJacobian(ts,A,A,(TSIJacobian)IJacobian,&ctx);
154: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
155: Set initial conditions
156: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
157: Solution(ts,0,U,&ctx);
158: TSSetInitialTimeStep(ts,4.0*3600,1.0);
159: TSSetSolution(ts,U);
161: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
162: Set solver options
163: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
164: TSSetDuration(ts,1000000,518400.0);
165: TSSetMaxStepRejections(ts,100);
166: TSSetMaxSNESFailures(ts,-1); /* unlimited */
167: TSSetFromOptions(ts);
169: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
170: Solve nonlinear system
171: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
172: TSSolve(ts,U);
174: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
175: Free work space. All PETSc objects should be destroyed when they
176: are no longer needed.
177: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
178: VecDestroy(&ctx.initialsolution);
179: MatDestroy(&A);
180: VecDestroy(&U);
181: TSDestroy(&ts);
183: PetscFinalize();
184: return(0);
185: }