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