1c4762a1bSJed Brown static char help[] = "Time-Dependent Reactive Flow example in 2D with Darcy Flow";
2c4762a1bSJed Brown
3c4762a1bSJed Brown /*
4c4762a1bSJed Brown
5c4762a1bSJed Brown This example solves the elementary chemical reaction:
6c4762a1bSJed Brown
7c4762a1bSJed Brown SP_A + 2SP_B = SP_C
8c4762a1bSJed Brown
9c4762a1bSJed Brown Subject to predetermined flow modeled as though it were in a porous media.
10c4762a1bSJed Brown This flow is modeled as advection and diffusion of the three species as
11c4762a1bSJed Brown
12c4762a1bSJed Brown v = porosity*saturation*grad(q)
13c4762a1bSJed Brown
14c4762a1bSJed Brown and the time-dependent equation solved for each species as
15c4762a1bSJed Brown advection + diffusion + reaction:
16c4762a1bSJed Brown
17c4762a1bSJed Brown v dot grad SP_i + dSP_i / dt + dispersivity*div*grad(SP_i) + R(SP_i) = 0
18c4762a1bSJed Brown
19c4762a1bSJed Brown The following options are available:
20c4762a1bSJed Brown
21c4762a1bSJed Brown -length_x - The length of the domain in the direction of the flow
22c4762a1bSJed Brown -length_y - The length of the domain in the direction orthogonal to the flow
23c4762a1bSJed Brown
24c4762a1bSJed Brown -gradq_inflow - The inflow speed as if the saturation and porosity were 1.
25c4762a1bSJed Brown -saturation - The saturation of the porous medium.
26c4762a1bSJed Brown -porosity - The porosity of the medium.
27c4762a1bSJed Brown -dispersivity - The dispersivity of the flow.
28c4762a1bSJed Brown -rate_constant - The rate constant for the chemical reaction
29c4762a1bSJed Brown -stoich - The stoichiometry matrix for the reaction
30c4762a1bSJed Brown -sp_inflow - The species concentrations at the inflow
31c4762a1bSJed Brown -sp_0 - The species concentrations at time 0
32c4762a1bSJed Brown
33c4762a1bSJed Brown */
34c4762a1bSJed Brown
35c4762a1bSJed Brown #include <petscdm.h>
36c4762a1bSJed Brown #include <petscdmda.h>
37c4762a1bSJed Brown #include <petscsnes.h>
38c4762a1bSJed Brown #include <petscts.h>
39c4762a1bSJed Brown
40c4762a1bSJed Brown #define N_SPECIES 3
41c4762a1bSJed Brown #define N_REACTIONS 1
42c4762a1bSJed Brown #define DIM 2
43c4762a1bSJed Brown
44c4762a1bSJed Brown #define stoich(i, j) ctx->stoichiometry[N_SPECIES * i + j]
45c4762a1bSJed Brown
46c4762a1bSJed Brown typedef struct {
47c4762a1bSJed Brown PetscScalar sp[N_SPECIES];
48c4762a1bSJed Brown } Field;
49c4762a1bSJed Brown
50c4762a1bSJed Brown typedef struct {
51c4762a1bSJed Brown Field x_inflow;
52c4762a1bSJed Brown Field x_0;
53c4762a1bSJed Brown PetscReal stoichiometry[N_SPECIES * N_REACTIONS];
54c4762a1bSJed Brown PetscReal porosity;
55c4762a1bSJed Brown PetscReal dispersivity;
56c4762a1bSJed Brown PetscReal saturation;
57c4762a1bSJed Brown PetscReal rate_constant[N_REACTIONS];
58c4762a1bSJed Brown PetscReal gradq_inflow;
59c4762a1bSJed Brown PetscReal length[DIM];
60c4762a1bSJed Brown } AppCtx;
61c4762a1bSJed Brown
62c4762a1bSJed Brown extern PetscErrorCode FormInitialGuess(DM da, AppCtx *ctx, Vec X);
63c4762a1bSJed Brown extern PetscErrorCode FormIFunctionLocal(DMDALocalInfo *, PetscReal, Field **, Field **, Field **, AppCtx *);
64c4762a1bSJed Brown extern PetscErrorCode FormIFunction(TS, PetscReal, Vec, Vec, Vec, void *);
65c4762a1bSJed Brown extern PetscErrorCode ReactingFlowPostCheck(SNESLineSearch, Vec, Vec, Vec, PetscBool *, PetscBool *, void *);
66c4762a1bSJed Brown
SetFromOptions(AppCtx * ctx)67d71ae5a4SJacob Faibussowitsch PetscErrorCode SetFromOptions(AppCtx *ctx)
68d71ae5a4SJacob Faibussowitsch {
69c4762a1bSJed Brown PetscInt i, j;
70c4762a1bSJed Brown
71c4762a1bSJed Brown PetscFunctionBeginUser;
72c4762a1bSJed Brown ctx->dispersivity = 0.5;
73c4762a1bSJed Brown ctx->porosity = 0.25;
74c4762a1bSJed Brown ctx->saturation = 1.0;
75c4762a1bSJed Brown ctx->gradq_inflow = 1.0;
76c4762a1bSJed Brown ctx->rate_constant[0] = 0.5;
77c4762a1bSJed Brown
78c4762a1bSJed Brown ctx->length[0] = 100.0;
79c4762a1bSJed Brown ctx->length[1] = 100.0;
80c4762a1bSJed Brown
81c4762a1bSJed Brown ctx->x_0.sp[0] = 0.0;
82c4762a1bSJed Brown ctx->x_0.sp[1] = 0.0;
83c4762a1bSJed Brown ctx->x_0.sp[2] = 0.0;
84c4762a1bSJed Brown
85c4762a1bSJed Brown ctx->x_inflow.sp[0] = 0.05;
86c4762a1bSJed Brown ctx->x_inflow.sp[1] = 0.05;
87c4762a1bSJed Brown ctx->x_inflow.sp[2] = 0.0;
88c4762a1bSJed Brown
89c4762a1bSJed Brown for (i = 0; i < N_REACTIONS; i++) {
90c4762a1bSJed Brown for (j = 0; j < N_SPECIES; j++) stoich(i, j) = 0.;
91c4762a1bSJed Brown }
92c4762a1bSJed Brown
93c4762a1bSJed Brown /* set up a pretty easy example */
94c4762a1bSJed Brown stoich(0, 0) = -1.;
95c4762a1bSJed Brown stoich(0, 1) = -2.;
96c4762a1bSJed Brown stoich(0, 2) = 1.;
97c4762a1bSJed Brown
98c4762a1bSJed Brown PetscInt as = N_SPECIES;
999566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-length_x", &ctx->length[0], NULL));
1009566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-length_y", &ctx->length[1], NULL));
1019566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-porosity", &ctx->porosity, NULL));
1029566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-saturation", &ctx->saturation, NULL));
1039566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-dispersivity", &ctx->dispersivity, NULL));
1049566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-gradq_inflow", &ctx->gradq_inflow, NULL));
1059566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-rate_constant", &ctx->rate_constant[0], NULL));
1069566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetRealArray(NULL, NULL, "-sp_inflow", ctx->x_inflow.sp, &as, NULL));
1079566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetRealArray(NULL, NULL, "-sp_0", ctx->x_0.sp, &as, NULL));
108c4762a1bSJed Brown as = N_SPECIES;
1099566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetRealArray(NULL, NULL, "-stoich", ctx->stoichiometry, &as, NULL));
1103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
111c4762a1bSJed Brown }
112c4762a1bSJed Brown
main(int argc,char ** argv)113d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
114d71ae5a4SJacob Faibussowitsch {
115c4762a1bSJed Brown TS ts;
116c4762a1bSJed Brown SNES snes;
117c4762a1bSJed Brown SNESLineSearch linesearch;
118c4762a1bSJed Brown Vec x;
119c4762a1bSJed Brown AppCtx ctx;
120c4762a1bSJed Brown DM da;
121c4762a1bSJed Brown
122327415f7SBarry Smith PetscFunctionBeginUser;
123*c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help));
1249566063dSJacob Faibussowitsch PetscCall(SetFromOptions(&ctx));
1259566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
1269566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSCN));
1279566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
1289566063dSJacob Faibussowitsch PetscCall(TSSetIFunction(ts, NULL, FormIFunction, &ctx));
1299566063dSJacob Faibussowitsch 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));
1309566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da));
1319566063dSJacob Faibussowitsch PetscCall(DMSetUp(da));
1329566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(da, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0));
1339566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(da, 0, "species A"));
1349566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(da, 1, "species B"));
1359566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(da, 2, "species C"));
1369566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(da, &ctx));
1379566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da, &x));
1389566063dSJacob Faibussowitsch PetscCall(FormInitialGuess(da, &ctx, x));
139c4762a1bSJed Brown
1409566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts, da));
1419566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, 1000.0));
1429566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
1439566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, 1.0));
1449566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts, x));
1459566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts));
146c4762a1bSJed Brown
1479566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &snes));
1489566063dSJacob Faibussowitsch PetscCall(SNESGetLineSearch(snes, &linesearch));
1499566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetPostCheck(linesearch, ReactingFlowPostCheck, (void *)&ctx));
1509566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes));
1519566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, x));
152c4762a1bSJed Brown
1539566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x));
1549566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts));
1559566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da));
1569566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
1573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
158c4762a1bSJed Brown }
159c4762a1bSJed Brown
160c4762a1bSJed Brown /* ------------------------------------------------------------------- */
161c4762a1bSJed Brown
FormInitialGuess(DM da,AppCtx * ctx,Vec X)162d71ae5a4SJacob Faibussowitsch PetscErrorCode FormInitialGuess(DM da, AppCtx *ctx, Vec X)
163d71ae5a4SJacob Faibussowitsch {
164c4762a1bSJed Brown PetscInt i, j, l, Mx, My, xs, ys, xm, ym;
165c4762a1bSJed Brown Field **x;
166c4762a1bSJed Brown
167c4762a1bSJed Brown PetscFunctionBeginUser;
1689371c9d4SSatish Balay 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));
1699566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, X, &x));
1709566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL));
171c4762a1bSJed Brown
172c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) {
173c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) {
174c4762a1bSJed Brown for (l = 0; l < N_SPECIES; l++) {
175c4762a1bSJed Brown if (i == 0) {
176c4762a1bSJed Brown if (l == 0) x[j][i].sp[l] = (ctx->x_inflow.sp[l] * ((PetscScalar)j) / (My - 1));
177c4762a1bSJed Brown else if (l == 1) x[j][i].sp[l] = (ctx->x_inflow.sp[l] * (1. - ((PetscScalar)j) / (My - 1)));
178c4762a1bSJed Brown else x[j][i].sp[l] = ctx->x_0.sp[l];
179c4762a1bSJed Brown }
180c4762a1bSJed Brown }
181c4762a1bSJed Brown }
182c4762a1bSJed Brown }
1839566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, X, &x));
1843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
185c4762a1bSJed Brown }
186c4762a1bSJed Brown
FormIFunctionLocal(DMDALocalInfo * info,PetscScalar ptime,Field ** x,Field ** xt,Field ** f,AppCtx * ctx)187d71ae5a4SJacob Faibussowitsch PetscErrorCode FormIFunctionLocal(DMDALocalInfo *info, PetscScalar ptime, Field **x, Field **xt, Field **f, AppCtx *ctx)
188d71ae5a4SJacob Faibussowitsch {
189c4762a1bSJed Brown PetscInt i, j, l, m;
190c4762a1bSJed Brown PetscReal hx, hy, dhx, dhy, hxdhy, hydhx, scale;
191c4762a1bSJed Brown PetscScalar u, uxx, uyy;
192c4762a1bSJed Brown PetscScalar vx, vy, sxp, syp, sxm, sym, avx, vxp, vxm, avy, vyp, vym, f_advect;
193c4762a1bSJed Brown PetscScalar rate;
194c4762a1bSJed Brown
195c4762a1bSJed Brown PetscFunctionBeginUser;
196c4762a1bSJed Brown hx = ctx->length[0] / ((PetscReal)(info->mx - 1));
197c4762a1bSJed Brown hy = ctx->length[1] / ((PetscReal)(info->my - 1));
198c4762a1bSJed Brown
199c4762a1bSJed Brown dhx = 1. / hx;
200c4762a1bSJed Brown dhy = 1. / hy;
201c4762a1bSJed Brown hxdhy = hx / hy;
202c4762a1bSJed Brown hydhx = hy / hx;
203c4762a1bSJed Brown scale = hx * hy;
204c4762a1bSJed Brown
205c4762a1bSJed Brown for (j = info->ys; j < info->ys + info->ym; j++) {
206c4762a1bSJed Brown for (i = info->xs; i < info->xs + info->xm; i++) {
207c4762a1bSJed Brown vx = ctx->gradq_inflow * ctx->porosity * ctx->saturation;
208c4762a1bSJed Brown vy = 0.0 * dhy;
209c4762a1bSJed Brown avx = PetscAbsScalar(vx);
2109371c9d4SSatish Balay vxp = .5 * (vx + avx);
2119371c9d4SSatish Balay vxm = .5 * (vx - avx);
212c4762a1bSJed Brown avy = PetscAbsScalar(vy);
2139371c9d4SSatish Balay vyp = .5 * (vy + avy);
2149371c9d4SSatish Balay vym = .5 * (vy - avy);
215c4762a1bSJed Brown /* chemical reactions */
216c4762a1bSJed Brown for (l = 0; l < N_SPECIES; l++) {
217c4762a1bSJed Brown /* determine the velocites as the gradients of the pressure */
218c4762a1bSJed Brown if (i == 0) {
219c4762a1bSJed Brown sxp = (x[j][i + 1].sp[l] - x[j][i].sp[l]) * dhx;
220c4762a1bSJed Brown sxm = sxp;
221c4762a1bSJed Brown } else if (i == info->mx - 1) {
222c4762a1bSJed Brown sxp = (x[j][i].sp[l] - x[j][i - 1].sp[l]) * dhx;
223c4762a1bSJed Brown sxm = sxp;
224c4762a1bSJed Brown } else {
225c4762a1bSJed Brown sxm = (x[j][i + 1].sp[l] - x[j][i].sp[l]) * dhx;
226c4762a1bSJed Brown sxp = (x[j][i].sp[l] - x[j][i - 1].sp[l]) * dhx;
227c4762a1bSJed Brown }
228c4762a1bSJed Brown if (j == 0) {
229c4762a1bSJed Brown syp = (x[j + 1][i].sp[l] - x[j][i].sp[l]) * dhy;
230c4762a1bSJed Brown sym = syp;
231c4762a1bSJed Brown } else if (j == info->my - 1) {
232c4762a1bSJed Brown syp = (x[j][i].sp[l] - x[j - 1][i].sp[l]) * dhy;
233c4762a1bSJed Brown sym = syp;
234c4762a1bSJed Brown } else {
235c4762a1bSJed Brown sym = (x[j + 1][i].sp[l] - x[j][i].sp[l]) * dhy;
236c4762a1bSJed Brown syp = (x[j][i].sp[l] - x[j - 1][i].sp[l]) * dhy;
237c4762a1bSJed Brown } /* 4 flops per species*point */
238c4762a1bSJed Brown
239c4762a1bSJed Brown if (i == 0) {
240c4762a1bSJed Brown if (l == 0) f[j][i].sp[l] = (x[j][i].sp[l] - ctx->x_inflow.sp[l] * ((PetscScalar)j) / (info->my - 1));
241c4762a1bSJed Brown 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)));
242c4762a1bSJed Brown else f[j][i].sp[l] = x[j][i].sp[l];
243c4762a1bSJed Brown
244c4762a1bSJed Brown } else {
245c4762a1bSJed Brown f[j][i].sp[l] = xt[j][i].sp[l] * scale;
246c4762a1bSJed Brown u = x[j][i].sp[l];
247c4762a1bSJed Brown if (j == 0) uyy = u - x[j + 1][i].sp[l];
248c4762a1bSJed Brown else if (j == info->my - 1) uyy = u - x[j - 1][i].sp[l];
249c4762a1bSJed Brown else uyy = (2.0 * u - x[j - 1][i].sp[l] - x[j + 1][i].sp[l]) * hxdhy;
250c4762a1bSJed Brown
251c4762a1bSJed Brown if (i != info->mx - 1) uxx = (2.0 * u - x[j][i - 1].sp[l] - x[j][i + 1].sp[l]) * hydhx;
252c4762a1bSJed Brown else uxx = u - x[j][i - 1].sp[l];
253c4762a1bSJed Brown /* 10 flops per species*point */
254c4762a1bSJed Brown
255c4762a1bSJed Brown f_advect = 0.;
256c4762a1bSJed Brown f_advect += scale * (vxp * sxp + vxm * sxm);
257c4762a1bSJed Brown f_advect += scale * (vyp * syp + vym * sym);
258c4762a1bSJed Brown f[j][i].sp[l] += f_advect + ctx->dispersivity * (uxx + uyy);
259c4762a1bSJed Brown /* 14 flops per species*point */
260c4762a1bSJed Brown }
261c4762a1bSJed Brown }
262c4762a1bSJed Brown /* reaction */
263c4762a1bSJed Brown if (i != 0) {
264c4762a1bSJed Brown for (m = 0; m < N_REACTIONS; m++) {
265c4762a1bSJed Brown rate = ctx->rate_constant[m];
266c4762a1bSJed Brown for (l = 0; l < N_SPECIES; l++) {
267c4762a1bSJed Brown if (stoich(m, l) < 0) {
268c4762a1bSJed Brown /* assume an elementary reaction */
269c4762a1bSJed Brown rate *= PetscPowScalar(x[j][i].sp[l], PetscAbsScalar(stoich(m, l)));
270c4762a1bSJed Brown /* ~10 flops per reaction per species per point */
271c4762a1bSJed Brown }
272c4762a1bSJed Brown }
273c4762a1bSJed Brown for (l = 0; l < N_SPECIES; l++) {
274c4762a1bSJed Brown f[j][i].sp[l] += -scale * stoich(m, l) * rate; /* Reaction term */
275c4762a1bSJed Brown /* 3 flops per reaction per species per point */
276c4762a1bSJed Brown }
277c4762a1bSJed Brown }
278c4762a1bSJed Brown }
279c4762a1bSJed Brown }
280c4762a1bSJed Brown }
2819566063dSJacob Faibussowitsch PetscCall(PetscLogFlops((N_SPECIES * (28.0 + 13.0 * N_REACTIONS)) * info->ym * info->xm));
2823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
283c4762a1bSJed Brown }
284c4762a1bSJed Brown
ReactingFlowPostCheck(SNESLineSearch linesearch,Vec X,Vec Y,Vec W,PetscBool * changed_y,PetscBool * changed_w,void * vctx)285d71ae5a4SJacob Faibussowitsch PetscErrorCode ReactingFlowPostCheck(SNESLineSearch linesearch, Vec X, Vec Y, Vec W, PetscBool *changed_y, PetscBool *changed_w, void *vctx)
286d71ae5a4SJacob Faibussowitsch {
287c4762a1bSJed Brown PetscInt i, j, l, Mx, My, xs, ys, xm, ym;
288c4762a1bSJed Brown Field **x;
289c4762a1bSJed Brown SNES snes;
290c4762a1bSJed Brown DM da;
291c4762a1bSJed Brown PetscScalar min;
292c4762a1bSJed Brown
293c4762a1bSJed Brown PetscFunctionBeginUser;
294c4762a1bSJed Brown *changed_w = PETSC_FALSE;
2959566063dSJacob Faibussowitsch PetscCall(VecMin(X, NULL, &min));
2963ba16761SJacob Faibussowitsch if (min >= 0.) PetscFunctionReturn(PETSC_SUCCESS);
297c4762a1bSJed Brown
298c4762a1bSJed Brown *changed_w = PETSC_TRUE;
2999566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetSNES(linesearch, &snes));
3009566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &da));
3019371c9d4SSatish Balay 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));
3029566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, W, &x));
3039566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL));
304c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) {
305c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) {
306c4762a1bSJed Brown for (l = 0; l < N_SPECIES; l++) {
307c4762a1bSJed Brown if (x[j][i].sp[l] < 0.) x[j][i].sp[l] = 0.;
308c4762a1bSJed Brown }
309c4762a1bSJed Brown }
310c4762a1bSJed Brown }
3119566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, W, &x));
3123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
313c4762a1bSJed Brown }
314c4762a1bSJed Brown
FormIFunction(TS ts,PetscReal ptime,Vec X,Vec Xdot,Vec F,void * user)315d71ae5a4SJacob Faibussowitsch PetscErrorCode FormIFunction(TS ts, PetscReal ptime, Vec X, Vec Xdot, Vec F, void *user)
316d71ae5a4SJacob Faibussowitsch {
317c4762a1bSJed Brown DMDALocalInfo info;
318c4762a1bSJed Brown Field **u, **udot, **fu;
319c4762a1bSJed Brown Vec localX;
320c4762a1bSJed Brown DM da;
321c4762a1bSJed Brown
322c4762a1bSJed Brown PetscFunctionBeginUser;
3239566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, (DM *)&da));
3249566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &localX));
3259566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, localX));
3269566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, localX));
3279566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(da, &info));
3289566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, localX, &u));
3299566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, Xdot, &udot));
3309566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, F, &fu));
3319566063dSJacob Faibussowitsch PetscCall(FormIFunctionLocal(&info, ptime, u, udot, fu, (AppCtx *)user));
3329566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, localX, &u));
3339566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, Xdot, &udot));
3349566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, F, &fu));
3359566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &localX));
3363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
337c4762a1bSJed Brown }
338