xref: /petsc/src/ts/tutorials/ex27.c (revision f11a936e8bf7fbe3c93dcc91a0846fab304f7fbd)
1 static char help[] = "Time-Dependent Reactive Flow example in 2D with Darcy Flow";
2 
3 /*T
4    Concepts: Solving a multicomponent time-dependent reactive flow system
5    Concepts: DMDA with timestepping
6    Processors: n
7 T*/
8 
9 /*
10 
11 This example solves the elementary chemical reaction:
12 
13 SP_A + 2SP_B = SP_C
14 
15 Subject to predetermined flow modeled as though it were in a porous media.
16 This flow is modeled as advection and diffusion of the three species as
17 
18 v = porosity*saturation*grad(q)
19 
20 and the time-dependent equation solved for each species as
21 advection + diffusion + reaction:
22 
23 v dot grad SP_i + dSP_i / dt + dispersivity*div*grad(SP_i) + R(SP_i) = 0
24 
25 The following options are available:
26 
27 -length_x - The length of the domain in the direction of the flow
28 -length_y - The length of the domain in the direction orthogonal to the flow
29 
30 -gradq_inflow - The inflow speed as if the saturation and porosity were 1.
31 -saturation - The saturation of the porous medium.
32 -porosity - The porosity of the medium.
33 -dispersivity - The dispersivity of the flow.
34 -rate_constant - The rate constant for the chemical reaction
35 -stoich - The stoichiometry matrix for the reaction
36 -sp_inflow - The species concentrations at the inflow
37 -sp_0 - The species concentrations at time 0
38 
39  */
40 
41 #include <petscdm.h>
42 #include <petscdmda.h>
43 #include <petscsnes.h>
44 #include <petscts.h>
45 
46 #define N_SPECIES 3
47 #define N_REACTIONS 1
48 #define DIM 2
49 
50 #define stoich(i, j)  ctx->stoichiometry[N_SPECIES*i + j]
51 
52 typedef struct {
53   PetscScalar sp[N_SPECIES];
54 } Field;
55 
56 typedef struct {
57 
58   Field     x_inflow;
59   Field     x_0;
60   PetscReal stoichiometry[N_SPECIES*N_REACTIONS];
61   PetscReal porosity;
62   PetscReal dispersivity;
63   PetscReal saturation;
64   PetscReal rate_constant[N_REACTIONS];
65   PetscReal gradq_inflow;
66   PetscReal length[DIM];
67 } AppCtx;
68 
69 extern PetscErrorCode FormInitialGuess(DM da,AppCtx *ctx,Vec X);
70 extern PetscErrorCode FormIFunctionLocal(DMDALocalInfo*,PetscReal,Field**,Field**,Field**,AppCtx*);
71 extern PetscErrorCode FormIFunction(TS,PetscReal,Vec,Vec,Vec,void*);
72 extern PetscErrorCode ReactingFlowPostCheck(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*);
73 
74 PetscErrorCode SetFromOptions(AppCtx * ctx)
75 {
76   PetscErrorCode ierr;
77   PetscInt       i,j;
78 
79   PetscFunctionBeginUser;
80   ctx->dispersivity     = 0.5;
81   ctx->porosity         = 0.25;
82   ctx->saturation       = 1.0;
83   ctx->gradq_inflow     = 1.0;
84   ctx->rate_constant[0] = 0.5;
85 
86   ctx->length[0] = 100.0;
87   ctx->length[1] = 100.0;
88 
89   ctx->x_0.sp[0] = 0.0;
90   ctx->x_0.sp[1] = 0.0;
91   ctx->x_0.sp[2] = 0.0;
92 
93   ctx->x_inflow.sp[0] = 0.05;
94   ctx->x_inflow.sp[1] = 0.05;
95   ctx->x_inflow.sp[2] = 0.0;
96 
97   for (i = 0; i < N_REACTIONS; i++) {
98     for (j = 0; j < N_SPECIES; j++) stoich(i, j) = 0.;
99   }
100 
101   /* set up a pretty easy example */
102   stoich(0, 0) = -1.;
103   stoich(0, 1) = -2.;
104   stoich(0, 2) = 1.;
105 
106   PetscInt as = N_SPECIES;
107   ierr = PetscOptionsGetReal(NULL,NULL,"-length_x",&ctx->length[0],NULL);CHKERRQ(ierr);
108   ierr = PetscOptionsGetReal(NULL,NULL,"-length_y",&ctx->length[1],NULL);CHKERRQ(ierr);
109   ierr = PetscOptionsGetReal(NULL,NULL,"-porosity",&ctx->porosity,NULL);CHKERRQ(ierr);
110   ierr = PetscOptionsGetReal(NULL,NULL,"-saturation",&ctx->saturation,NULL);CHKERRQ(ierr);
111   ierr = PetscOptionsGetReal(NULL,NULL,"-dispersivity",&ctx->dispersivity,NULL);CHKERRQ(ierr);
112   ierr = PetscOptionsGetReal(NULL,NULL,"-gradq_inflow",&ctx->gradq_inflow,NULL);CHKERRQ(ierr);
113   ierr = PetscOptionsGetReal(NULL,NULL,"-rate_constant",&ctx->rate_constant[0],NULL);CHKERRQ(ierr);
114   ierr = PetscOptionsGetRealArray(NULL,NULL,"-sp_inflow",ctx->x_inflow.sp,&as,NULL);CHKERRQ(ierr);
115   ierr = PetscOptionsGetRealArray(NULL,NULL,"-sp_0",ctx->x_0.sp,&as,NULL);CHKERRQ(ierr);
116   as   = N_SPECIES;
117   ierr = PetscOptionsGetRealArray(NULL,NULL,"-stoich",ctx->stoichiometry,&as,NULL);CHKERRQ(ierr);
118   PetscFunctionReturn(0);
119 }
120 
121 int main(int argc,char **argv)
122 {
123   TS             ts;
124   SNES           snes;
125   SNESLineSearch linesearch;
126   Vec            x;
127   AppCtx         ctx;
128   PetscErrorCode ierr;
129   DM             da;
130 
131   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
132   ierr = SetFromOptions(&ctx);CHKERRQ(ierr);
133   ierr = TSCreate(PETSC_COMM_WORLD, &ts);CHKERRQ(ierr);
134   ierr = TSSetType(ts,TSCN);CHKERRQ(ierr);
135   ierr = TSSetProblemType(ts,TS_NONLINEAR);CHKERRQ(ierr);
136   ierr = TSSetIFunction(ts, NULL, FormIFunction, &ctx);CHKERRQ(ierr);
137   ierr = DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,-4,-4,PETSC_DECIDE,PETSC_DECIDE,N_SPECIES,1,NULL,NULL,&da);CHKERRQ(ierr);
138   ierr = DMSetFromOptions(da);CHKERRQ(ierr);
139   ierr = DMSetUp(da);CHKERRQ(ierr);
140   ierr = DMDASetUniformCoordinates(da, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0);CHKERRQ(ierr);
141   ierr = DMDASetFieldName(da,0,"species A");CHKERRQ(ierr);
142   ierr = DMDASetFieldName(da,1,"species B");CHKERRQ(ierr);
143   ierr = DMDASetFieldName(da,2,"species C");CHKERRQ(ierr);
144   ierr = DMSetApplicationContext(da,&ctx);CHKERRQ(ierr);
145   ierr = DMCreateGlobalVector(da,&x);CHKERRQ(ierr);
146   ierr = FormInitialGuess(da, &ctx, x);CHKERRQ(ierr);
147 
148   ierr = TSSetDM(ts, da);CHKERRQ(ierr);
149   ierr = TSSetMaxTime(ts,1000.0);CHKERRQ(ierr);
150   ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr);
151   ierr = TSSetTimeStep(ts,1.0);CHKERRQ(ierr);
152   ierr = TSSetSolution(ts,x);CHKERRQ(ierr);
153   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
154 
155   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
156   ierr = SNESGetLineSearch(snes,&linesearch);CHKERRQ(ierr);
157   ierr = SNESLineSearchSetPostCheck(linesearch, ReactingFlowPostCheck, (void*)&ctx);CHKERRQ(ierr);
158   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
159   ierr = TSSolve(ts,x);CHKERRQ(ierr);
160 
161   ierr = VecDestroy(&x);CHKERRQ(ierr);
162   ierr = TSDestroy(&ts);CHKERRQ(ierr);
163   ierr = DMDestroy(&da);CHKERRQ(ierr);
164   ierr = PetscFinalize();
165   PetscFunctionReturn(0);
166 }
167 
168 /* ------------------------------------------------------------------- */
169 
170 PetscErrorCode FormInitialGuess(DM da,AppCtx *ctx,Vec X)
171 {
172   PetscInt       i,j,l,Mx,My,xs,ys,xm,ym;
173   PetscErrorCode ierr;
174   Field          **x;
175 
176   PetscFunctionBeginUser;
177   ierr = DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
178                      PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
179 
180   ierr = DMDAVecGetArray(da,X,&x);CHKERRQ(ierr);
181   ierr = DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr);
182 
183   for (j=ys; j<ys+ym; j++) {
184     for (i=xs; i<xs+xm; i++) {
185       for (l = 0; l < N_SPECIES; l++) {
186         if (i == 0) {
187           if (l == 0)      x[j][i].sp[l] = (ctx->x_inflow.sp[l]*((PetscScalar)j) / (My - 1));
188           else if (l == 1) x[j][i].sp[l] = (ctx->x_inflow.sp[l]*(1. - ((PetscScalar)j) / (My - 1)));
189           else             x[j][i].sp[l] = ctx->x_0.sp[l];
190         }
191       }
192     }
193   }
194   ierr = DMDAVecRestoreArray(da,X,&x);CHKERRQ(ierr);
195   PetscFunctionReturn(0);
196 
197 }
198 
199 PetscErrorCode FormIFunctionLocal(DMDALocalInfo *info,PetscScalar ptime,Field **x,Field **xt,Field **f,AppCtx *ctx)
200 {
201   PetscErrorCode ierr;
202   PetscInt       i,j,l,m;
203   PetscReal      hx,hy,dhx,dhy,hxdhy,hydhx,scale;
204   PetscScalar    u,uxx,uyy;
205   PetscScalar    vx, vy,sxp,syp,sxm,sym,avx,vxp,vxm,avy,vyp,vym,f_advect;
206   PetscScalar    rate;
207 
208   PetscFunctionBeginUser;
209   hx = ctx->length[0]/((PetscReal)(info->mx-1));
210   hy = ctx->length[1]/((PetscReal)(info->my-1));
211 
212   dhx   =     1. / hx;
213   dhy   =     1. / hy;
214   hxdhy =  hx/hy;
215   hydhx =  hy/hx;
216   scale =   hx*hy;
217 
218   for (j=info->ys; j<info->ys+info->ym; j++) {
219     for (i=info->xs; i<info->xs+info->xm; i++) {
220       vx  = ctx->gradq_inflow*ctx->porosity*ctx->saturation;
221       vy  = 0.0*dhy;
222       avx = PetscAbsScalar(vx);
223       vxp = .5*(vx+avx); vxm = .5*(vx-avx);
224       avy = PetscAbsScalar(vy);
225       vyp = .5*(vy+avy); vym = .5*(vy-avy);
226       /* chemical reactions */
227       for (l = 0; l < N_SPECIES; l++) {
228         /* determine the velocites as the gradients of the pressure */
229         if (i == 0) {
230           sxp = (x[j][i+1].sp[l] - x[j][i].sp[l])*dhx;
231           sxm = sxp;
232         } else if (i == info->mx - 1) {
233           sxp = (x[j][i].sp[l] - x[j][i-1].sp[l])*dhx;
234           sxm = sxp;
235         } else {
236           sxm = (x[j][i+1].sp[l] - x[j][i].sp[l])*dhx;
237           sxp = (x[j][i].sp[l] - x[j][i-1].sp[l])*dhx;
238         }
239         if (j == 0) {
240           syp = (x[j+1][i].sp[l] - x[j][i].sp[l])*dhy;
241           sym = syp;
242         } else if (j == info->my - 1) {
243           syp = (x[j][i].sp[l] - x[j-1][i].sp[l])*dhy;
244           sym = syp;
245         } else {
246           sym = (x[j+1][i].sp[l] - x[j][i].sp[l])*dhy;
247           syp = (x[j][i].sp[l] - x[j-1][i].sp[l])*dhy;
248         } /* 4 flops per species*point */
249 
250         if (i == 0) {
251           if (l == 0)      f[j][i].sp[l] = (x[j][i].sp[l] - ctx->x_inflow.sp[l]*((PetscScalar)j) / (info->my - 1));
252           else if (l == 1) f[j][i].sp[l] = (x[j][i].sp[l] - ctx->x_inflow.sp[l]*(1. - ((PetscScalar)j) / (info->my - 1)));
253           else             f[j][i].sp[l] = x[j][i].sp[l];
254 
255         } else {
256           f[j][i].sp[l] = xt[j][i].sp[l]*scale;
257           u       = x[j][i].sp[l];
258           if (j == 0) uyy = u - x[j+1][i].sp[l];
259           else if (j == info->my - 1) uyy = u - x[j-1][i].sp[l];
260           else                        uyy = (2.0*u - x[j-1][i].sp[l] - x[j+1][i].sp[l])*hxdhy;
261 
262           if (i != info->mx - 1) uxx = (2.0*u - x[j][i-1].sp[l] - x[j][i+1].sp[l])*hydhx;
263           else                   uxx = u - x[j][i-1].sp[l];
264           /* 10 flops per species*point */
265 
266           f_advect       = 0.;
267           f_advect      += scale*(vxp*sxp + vxm*sxm);
268           f_advect      += scale*(vyp*syp + vym*sym);
269           f[j][i].sp[l] += f_advect + ctx->dispersivity*(uxx + uyy);
270           /* 14 flops per species*point */
271         }
272       }
273       /* reaction */
274       if (i != 0) {
275         for (m = 0; m < N_REACTIONS; m++) {
276           rate = ctx->rate_constant[m];
277           for (l = 0; l < N_SPECIES; l++) {
278             if (stoich(m, l) < 0) {
279               /* assume an elementary reaction */
280               rate *= PetscPowScalar(x[j][i].sp[l], PetscAbsScalar(stoich(m, l)));
281               /* ~10 flops per reaction per species per point */
282             }
283           }
284           for (l = 0; l < N_SPECIES; l++) {
285               f[j][i].sp[l] += -scale*stoich(m, l)*rate;  /* Reaction term */
286               /* 3 flops per reaction per species per point */
287           }
288         }
289       }
290     }
291   }
292   ierr = PetscLogFlops((N_SPECIES*(28.0 + 13.0*N_REACTIONS))*info->ym*info->xm);CHKERRQ(ierr);
293   PetscFunctionReturn(0);
294 }
295 
296 PetscErrorCode ReactingFlowPostCheck(SNESLineSearch linesearch, Vec X, Vec Y, Vec W, PetscBool *changed_y, PetscBool *changed_w, void *vctx)
297 {
298   PetscInt       i,j,l,Mx,My,xs,ys,xm,ym;
299   PetscErrorCode ierr;
300   Field          **x;
301   SNES           snes;
302   DM             da;
303   PetscScalar    min;
304 
305   PetscFunctionBeginUser;
306    *changed_w = PETSC_FALSE;
307   ierr = VecMin(X,NULL,&min);CHKERRQ(ierr);
308   if (min >= 0.) PetscFunctionReturn(0);
309 
310   *changed_w = PETSC_TRUE;
311   ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr);
312   ierr = SNESGetDM(snes,&da);CHKERRQ(ierr);
313   ierr = DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
314                      PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
315   ierr = DMDAVecGetArray(da,W,&x);CHKERRQ(ierr);
316   ierr = DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr);
317   for (j=ys; j<ys+ym; j++) {
318     for (i=xs; i<xs+xm; i++) {
319       for (l = 0; l < N_SPECIES; l++) {
320         if (x[j][i].sp[l] < 0.) x[j][i].sp[l] = 0.;
321       }
322     }
323   }
324   ierr = DMDAVecRestoreArray(da,W,&x);CHKERRQ(ierr);
325   PetscFunctionReturn(0);
326 }
327 
328 PetscErrorCode FormIFunction(TS ts,PetscReal ptime,Vec X,Vec Xdot,Vec F,void *user)
329 {
330   DMDALocalInfo  info;
331   Field          **u,**udot,**fu;
332   PetscErrorCode ierr;
333   Vec            localX;
334   DM             da;
335 
336   PetscFunctionBeginUser;
337   ierr = TSGetDM(ts,(DM*)&da);CHKERRQ(ierr);
338   ierr = DMGetLocalVector(da,&localX);CHKERRQ(ierr);
339   ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
340   ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
341   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
342   ierr = DMDAVecGetArrayRead(da,localX,&u);CHKERRQ(ierr);
343   ierr = DMDAVecGetArrayRead(da,Xdot,&udot);CHKERRQ(ierr);
344   ierr = DMDAVecGetArray(da,F,&fu);CHKERRQ(ierr);
345   ierr = FormIFunctionLocal(&info,ptime,u,udot,fu,(AppCtx*)user);CHKERRQ(ierr);
346   ierr = DMDAVecRestoreArrayRead(da,localX,&u);CHKERRQ(ierr);
347   ierr = DMDAVecRestoreArrayRead(da,Xdot,&udot);CHKERRQ(ierr);
348   ierr = DMDAVecRestoreArray(da,F,&fu);CHKERRQ(ierr);
349   ierr = DMRestoreLocalVector(da,&localX);CHKERRQ(ierr);
350   PetscFunctionReturn(0);
351 }
352