xref: /petsc/src/ts/tutorials/ex27.c (revision 76d901e46dda72c1afe96306c7cb4731c47d4e87)
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   PetscFunctionBeginUser;
124   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
125   PetscCall(SetFromOptions(&ctx));
126   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
127   PetscCall(TSSetType(ts,TSCN));
128   PetscCall(TSSetProblemType(ts,TS_NONLINEAR));
129   PetscCall(TSSetIFunction(ts, NULL, FormIFunction, &ctx));
130   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));
131   PetscCall(DMSetFromOptions(da));
132   PetscCall(DMSetUp(da));
133   PetscCall(DMDASetUniformCoordinates(da, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0));
134   PetscCall(DMDASetFieldName(da,0,"species A"));
135   PetscCall(DMDASetFieldName(da,1,"species B"));
136   PetscCall(DMDASetFieldName(da,2,"species C"));
137   PetscCall(DMSetApplicationContext(da,&ctx));
138   PetscCall(DMCreateGlobalVector(da,&x));
139   PetscCall(FormInitialGuess(da, &ctx, x));
140 
141   PetscCall(TSSetDM(ts, da));
142   PetscCall(TSSetMaxTime(ts,1000.0));
143   PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER));
144   PetscCall(TSSetTimeStep(ts,1.0));
145   PetscCall(TSSetSolution(ts,x));
146   PetscCall(TSSetFromOptions(ts));
147 
148   PetscCall(TSGetSNES(ts,&snes));
149   PetscCall(SNESGetLineSearch(snes,&linesearch));
150   PetscCall(SNESLineSearchSetPostCheck(linesearch, ReactingFlowPostCheck, (void*)&ctx));
151   PetscCall(SNESSetFromOptions(snes));
152   PetscCall(TSSolve(ts,x));
153 
154   PetscCall(VecDestroy(&x));
155   PetscCall(TSDestroy(&ts));
156   PetscCall(DMDestroy(&da));
157   PetscCall(PetscFinalize());
158   PetscFunctionReturn(0);
159 }
160 
161 /* ------------------------------------------------------------------- */
162 
163 PetscErrorCode FormInitialGuess(DM da,AppCtx *ctx,Vec X)
164 {
165   PetscInt       i,j,l,Mx,My,xs,ys,xm,ym;
166   Field          **x;
167 
168   PetscFunctionBeginUser;
169   PetscCall(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   PetscCall(DMDAVecGetArray(da,X,&x));
172   PetscCall(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL));
173 
174   for (j=ys; j<ys+ym; j++) {
175     for (i=xs; i<xs+xm; i++) {
176       for (l = 0; l < N_SPECIES; l++) {
177         if (i == 0) {
178           if (l == 0)      x[j][i].sp[l] = (ctx->x_inflow.sp[l]*((PetscScalar)j) / (My - 1));
179           else if (l == 1) x[j][i].sp[l] = (ctx->x_inflow.sp[l]*(1. - ((PetscScalar)j) / (My - 1)));
180           else             x[j][i].sp[l] = ctx->x_0.sp[l];
181         }
182       }
183     }
184   }
185   PetscCall(DMDAVecRestoreArray(da,X,&x));
186   PetscFunctionReturn(0);
187 
188 }
189 
190 PetscErrorCode FormIFunctionLocal(DMDALocalInfo *info,PetscScalar ptime,Field **x,Field **xt,Field **f,AppCtx *ctx)
191 {
192   PetscInt       i,j,l,m;
193   PetscReal      hx,hy,dhx,dhy,hxdhy,hydhx,scale;
194   PetscScalar    u,uxx,uyy;
195   PetscScalar    vx, vy,sxp,syp,sxm,sym,avx,vxp,vxm,avy,vyp,vym,f_advect;
196   PetscScalar    rate;
197 
198   PetscFunctionBeginUser;
199   hx = ctx->length[0]/((PetscReal)(info->mx-1));
200   hy = ctx->length[1]/((PetscReal)(info->my-1));
201 
202   dhx   =     1. / hx;
203   dhy   =     1. / hy;
204   hxdhy =  hx/hy;
205   hydhx =  hy/hx;
206   scale =   hx*hy;
207 
208   for (j=info->ys; j<info->ys+info->ym; j++) {
209     for (i=info->xs; i<info->xs+info->xm; i++) {
210       vx  = ctx->gradq_inflow*ctx->porosity*ctx->saturation;
211       vy  = 0.0*dhy;
212       avx = PetscAbsScalar(vx);
213       vxp = .5*(vx+avx); vxm = .5*(vx-avx);
214       avy = PetscAbsScalar(vy);
215       vyp = .5*(vy+avy); vym = .5*(vy-avy);
216       /* chemical reactions */
217       for (l = 0; l < N_SPECIES; l++) {
218         /* determine the velocites as the gradients of the pressure */
219         if (i == 0) {
220           sxp = (x[j][i+1].sp[l] - x[j][i].sp[l])*dhx;
221           sxm = sxp;
222         } else if (i == info->mx - 1) {
223           sxp = (x[j][i].sp[l] - x[j][i-1].sp[l])*dhx;
224           sxm = sxp;
225         } else {
226           sxm = (x[j][i+1].sp[l] - x[j][i].sp[l])*dhx;
227           sxp = (x[j][i].sp[l] - x[j][i-1].sp[l])*dhx;
228         }
229         if (j == 0) {
230           syp = (x[j+1][i].sp[l] - x[j][i].sp[l])*dhy;
231           sym = syp;
232         } else if (j == info->my - 1) {
233           syp = (x[j][i].sp[l] - x[j-1][i].sp[l])*dhy;
234           sym = syp;
235         } else {
236           sym = (x[j+1][i].sp[l] - x[j][i].sp[l])*dhy;
237           syp = (x[j][i].sp[l] - x[j-1][i].sp[l])*dhy;
238         } /* 4 flops per species*point */
239 
240         if (i == 0) {
241           if (l == 0)      f[j][i].sp[l] = (x[j][i].sp[l] - ctx->x_inflow.sp[l]*((PetscScalar)j) / (info->my - 1));
242           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)));
243           else             f[j][i].sp[l] = x[j][i].sp[l];
244 
245         } else {
246           f[j][i].sp[l] = xt[j][i].sp[l]*scale;
247           u       = x[j][i].sp[l];
248           if (j == 0) uyy = u - x[j+1][i].sp[l];
249           else if (j == info->my - 1) uyy = u - x[j-1][i].sp[l];
250           else                        uyy = (2.0*u - x[j-1][i].sp[l] - x[j+1][i].sp[l])*hxdhy;
251 
252           if (i != info->mx - 1) uxx = (2.0*u - x[j][i-1].sp[l] - x[j][i+1].sp[l])*hydhx;
253           else                   uxx = u - x[j][i-1].sp[l];
254           /* 10 flops per species*point */
255 
256           f_advect       = 0.;
257           f_advect      += scale*(vxp*sxp + vxm*sxm);
258           f_advect      += scale*(vyp*syp + vym*sym);
259           f[j][i].sp[l] += f_advect + ctx->dispersivity*(uxx + uyy);
260           /* 14 flops per species*point */
261         }
262       }
263       /* reaction */
264       if (i != 0) {
265         for (m = 0; m < N_REACTIONS; m++) {
266           rate = ctx->rate_constant[m];
267           for (l = 0; l < N_SPECIES; l++) {
268             if (stoich(m, l) < 0) {
269               /* assume an elementary reaction */
270               rate *= PetscPowScalar(x[j][i].sp[l], PetscAbsScalar(stoich(m, l)));
271               /* ~10 flops per reaction per species per point */
272             }
273           }
274           for (l = 0; l < N_SPECIES; l++) {
275               f[j][i].sp[l] += -scale*stoich(m, l)*rate;  /* Reaction term */
276               /* 3 flops per reaction per species per point */
277           }
278         }
279       }
280     }
281   }
282   PetscCall(PetscLogFlops((N_SPECIES*(28.0 + 13.0*N_REACTIONS))*info->ym*info->xm));
283   PetscFunctionReturn(0);
284 }
285 
286 PetscErrorCode ReactingFlowPostCheck(SNESLineSearch linesearch, Vec X, Vec Y, Vec W, PetscBool *changed_y, PetscBool *changed_w, void *vctx)
287 {
288   PetscInt       i,j,l,Mx,My,xs,ys,xm,ym;
289   Field          **x;
290   SNES           snes;
291   DM             da;
292   PetscScalar    min;
293 
294   PetscFunctionBeginUser;
295    *changed_w = PETSC_FALSE;
296   PetscCall(VecMin(X,NULL,&min));
297   if (min >= 0.) PetscFunctionReturn(0);
298 
299   *changed_w = PETSC_TRUE;
300   PetscCall(SNESLineSearchGetSNES(linesearch, &snes));
301   PetscCall(SNESGetDM(snes,&da));
302   PetscCall(DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
303                         PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE));
304   PetscCall(DMDAVecGetArray(da,W,&x));
305   PetscCall(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL));
306   for (j=ys; j<ys+ym; j++) {
307     for (i=xs; i<xs+xm; i++) {
308       for (l = 0; l < N_SPECIES; l++) {
309         if (x[j][i].sp[l] < 0.) x[j][i].sp[l] = 0.;
310       }
311     }
312   }
313   PetscCall(DMDAVecRestoreArray(da,W,&x));
314   PetscFunctionReturn(0);
315 }
316 
317 PetscErrorCode FormIFunction(TS ts,PetscReal ptime,Vec X,Vec Xdot,Vec F,void *user)
318 {
319   DMDALocalInfo  info;
320   Field          **u,**udot,**fu;
321   Vec            localX;
322   DM             da;
323 
324   PetscFunctionBeginUser;
325   PetscCall(TSGetDM(ts,(DM*)&da));
326   PetscCall(DMGetLocalVector(da,&localX));
327   PetscCall(DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX));
328   PetscCall(DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX));
329   PetscCall(DMDAGetLocalInfo(da,&info));
330   PetscCall(DMDAVecGetArrayRead(da,localX,&u));
331   PetscCall(DMDAVecGetArrayRead(da,Xdot,&udot));
332   PetscCall(DMDAVecGetArray(da,F,&fu));
333   PetscCall(FormIFunctionLocal(&info,ptime,u,udot,fu,(AppCtx*)user));
334   PetscCall(DMDAVecRestoreArrayRead(da,localX,&u));
335   PetscCall(DMDAVecRestoreArrayRead(da,Xdot,&udot));
336   PetscCall(DMDAVecRestoreArray(da,F,&fu));
337   PetscCall(DMRestoreLocalVector(da,&localX));
338   PetscFunctionReturn(0);
339 }
340