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