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 67*9371c9d4SSatish Balay PetscErrorCode SetFromOptions(AppCtx *ctx) { 68c4762a1bSJed Brown PetscInt i, j; 69c4762a1bSJed Brown 70c4762a1bSJed Brown PetscFunctionBeginUser; 71c4762a1bSJed Brown ctx->dispersivity = 0.5; 72c4762a1bSJed Brown ctx->porosity = 0.25; 73c4762a1bSJed Brown ctx->saturation = 1.0; 74c4762a1bSJed Brown ctx->gradq_inflow = 1.0; 75c4762a1bSJed Brown ctx->rate_constant[0] = 0.5; 76c4762a1bSJed Brown 77c4762a1bSJed Brown ctx->length[0] = 100.0; 78c4762a1bSJed Brown ctx->length[1] = 100.0; 79c4762a1bSJed Brown 80c4762a1bSJed Brown ctx->x_0.sp[0] = 0.0; 81c4762a1bSJed Brown ctx->x_0.sp[1] = 0.0; 82c4762a1bSJed Brown ctx->x_0.sp[2] = 0.0; 83c4762a1bSJed Brown 84c4762a1bSJed Brown ctx->x_inflow.sp[0] = 0.05; 85c4762a1bSJed Brown ctx->x_inflow.sp[1] = 0.05; 86c4762a1bSJed Brown ctx->x_inflow.sp[2] = 0.0; 87c4762a1bSJed Brown 88c4762a1bSJed Brown for (i = 0; i < N_REACTIONS; i++) { 89c4762a1bSJed Brown for (j = 0; j < N_SPECIES; j++) stoich(i, j) = 0.; 90c4762a1bSJed Brown } 91c4762a1bSJed Brown 92c4762a1bSJed Brown /* set up a pretty easy example */ 93c4762a1bSJed Brown stoich(0, 0) = -1.; 94c4762a1bSJed Brown stoich(0, 1) = -2.; 95c4762a1bSJed Brown stoich(0, 2) = 1.; 96c4762a1bSJed Brown 97c4762a1bSJed Brown PetscInt as = N_SPECIES; 989566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-length_x", &ctx->length[0], NULL)); 999566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-length_y", &ctx->length[1], NULL)); 1009566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-porosity", &ctx->porosity, NULL)); 1019566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-saturation", &ctx->saturation, NULL)); 1029566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-dispersivity", &ctx->dispersivity, NULL)); 1039566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-gradq_inflow", &ctx->gradq_inflow, NULL)); 1049566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-rate_constant", &ctx->rate_constant[0], NULL)); 1059566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetRealArray(NULL, NULL, "-sp_inflow", ctx->x_inflow.sp, &as, NULL)); 1069566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetRealArray(NULL, NULL, "-sp_0", ctx->x_0.sp, &as, NULL)); 107c4762a1bSJed Brown as = N_SPECIES; 1089566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetRealArray(NULL, NULL, "-stoich", ctx->stoichiometry, &as, NULL)); 109c4762a1bSJed Brown PetscFunctionReturn(0); 110c4762a1bSJed Brown } 111c4762a1bSJed Brown 112*9371c9d4SSatish Balay int main(int argc, char **argv) { 113c4762a1bSJed Brown TS ts; 114c4762a1bSJed Brown SNES snes; 115c4762a1bSJed Brown SNESLineSearch linesearch; 116c4762a1bSJed Brown Vec x; 117c4762a1bSJed Brown AppCtx ctx; 118c4762a1bSJed Brown DM da; 119c4762a1bSJed Brown 120327415f7SBarry Smith PetscFunctionBeginUser; 1219566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 1229566063dSJacob Faibussowitsch PetscCall(SetFromOptions(&ctx)); 1239566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 1249566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSCN)); 1259566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ts, TS_NONLINEAR)); 1269566063dSJacob Faibussowitsch PetscCall(TSSetIFunction(ts, NULL, FormIFunction, &ctx)); 1279566063dSJacob 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)); 1289566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da)); 1299566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 1309566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(da, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0)); 1319566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(da, 0, "species A")); 1329566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(da, 1, "species B")); 1339566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(da, 2, "species C")); 1349566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(da, &ctx)); 1359566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da, &x)); 1369566063dSJacob Faibussowitsch PetscCall(FormInitialGuess(da, &ctx, x)); 137c4762a1bSJed Brown 1389566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts, da)); 1399566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, 1000.0)); 1409566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER)); 1419566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, 1.0)); 1429566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts, x)); 1439566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 144c4762a1bSJed Brown 1459566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &snes)); 1469566063dSJacob Faibussowitsch PetscCall(SNESGetLineSearch(snes, &linesearch)); 1479566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetPostCheck(linesearch, ReactingFlowPostCheck, (void *)&ctx)); 1489566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 1499566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, x)); 150c4762a1bSJed Brown 1519566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 1529566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 1539566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da)); 1549566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 155c4762a1bSJed Brown PetscFunctionReturn(0); 156c4762a1bSJed Brown } 157c4762a1bSJed Brown 158c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 159c4762a1bSJed Brown 160*9371c9d4SSatish Balay PetscErrorCode FormInitialGuess(DM da, AppCtx *ctx, Vec X) { 161c4762a1bSJed Brown PetscInt i, j, l, Mx, My, xs, ys, xm, ym; 162c4762a1bSJed Brown Field **x; 163c4762a1bSJed Brown 164c4762a1bSJed Brown PetscFunctionBeginUser; 165*9371c9d4SSatish 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)); 1669566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, X, &x)); 1679566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL)); 168c4762a1bSJed Brown 169c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) { 170c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 171c4762a1bSJed Brown for (l = 0; l < N_SPECIES; l++) { 172c4762a1bSJed Brown if (i == 0) { 173c4762a1bSJed Brown if (l == 0) x[j][i].sp[l] = (ctx->x_inflow.sp[l] * ((PetscScalar)j) / (My - 1)); 174c4762a1bSJed Brown else if (l == 1) x[j][i].sp[l] = (ctx->x_inflow.sp[l] * (1. - ((PetscScalar)j) / (My - 1))); 175c4762a1bSJed Brown else x[j][i].sp[l] = ctx->x_0.sp[l]; 176c4762a1bSJed Brown } 177c4762a1bSJed Brown } 178c4762a1bSJed Brown } 179c4762a1bSJed Brown } 1809566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, X, &x)); 181c4762a1bSJed Brown PetscFunctionReturn(0); 182c4762a1bSJed Brown } 183c4762a1bSJed Brown 184*9371c9d4SSatish Balay PetscErrorCode FormIFunctionLocal(DMDALocalInfo *info, PetscScalar ptime, Field **x, Field **xt, Field **f, AppCtx *ctx) { 185c4762a1bSJed Brown PetscInt i, j, l, m; 186c4762a1bSJed Brown PetscReal hx, hy, dhx, dhy, hxdhy, hydhx, scale; 187c4762a1bSJed Brown PetscScalar u, uxx, uyy; 188c4762a1bSJed Brown PetscScalar vx, vy, sxp, syp, sxm, sym, avx, vxp, vxm, avy, vyp, vym, f_advect; 189c4762a1bSJed Brown PetscScalar rate; 190c4762a1bSJed Brown 191c4762a1bSJed Brown PetscFunctionBeginUser; 192c4762a1bSJed Brown hx = ctx->length[0] / ((PetscReal)(info->mx - 1)); 193c4762a1bSJed Brown hy = ctx->length[1] / ((PetscReal)(info->my - 1)); 194c4762a1bSJed Brown 195c4762a1bSJed Brown dhx = 1. / hx; 196c4762a1bSJed Brown dhy = 1. / hy; 197c4762a1bSJed Brown hxdhy = hx / hy; 198c4762a1bSJed Brown hydhx = hy / hx; 199c4762a1bSJed Brown scale = hx * hy; 200c4762a1bSJed Brown 201c4762a1bSJed Brown for (j = info->ys; j < info->ys + info->ym; j++) { 202c4762a1bSJed Brown for (i = info->xs; i < info->xs + info->xm; i++) { 203c4762a1bSJed Brown vx = ctx->gradq_inflow * ctx->porosity * ctx->saturation; 204c4762a1bSJed Brown vy = 0.0 * dhy; 205c4762a1bSJed Brown avx = PetscAbsScalar(vx); 206*9371c9d4SSatish Balay vxp = .5 * (vx + avx); 207*9371c9d4SSatish Balay vxm = .5 * (vx - avx); 208c4762a1bSJed Brown avy = PetscAbsScalar(vy); 209*9371c9d4SSatish Balay vyp = .5 * (vy + avy); 210*9371c9d4SSatish Balay vym = .5 * (vy - avy); 211c4762a1bSJed Brown /* chemical reactions */ 212c4762a1bSJed Brown for (l = 0; l < N_SPECIES; l++) { 213c4762a1bSJed Brown /* determine the velocites as the gradients of the pressure */ 214c4762a1bSJed Brown if (i == 0) { 215c4762a1bSJed Brown sxp = (x[j][i + 1].sp[l] - x[j][i].sp[l]) * dhx; 216c4762a1bSJed Brown sxm = sxp; 217c4762a1bSJed Brown } else if (i == info->mx - 1) { 218c4762a1bSJed Brown sxp = (x[j][i].sp[l] - x[j][i - 1].sp[l]) * dhx; 219c4762a1bSJed Brown sxm = sxp; 220c4762a1bSJed Brown } else { 221c4762a1bSJed Brown sxm = (x[j][i + 1].sp[l] - x[j][i].sp[l]) * dhx; 222c4762a1bSJed Brown sxp = (x[j][i].sp[l] - x[j][i - 1].sp[l]) * dhx; 223c4762a1bSJed Brown } 224c4762a1bSJed Brown if (j == 0) { 225c4762a1bSJed Brown syp = (x[j + 1][i].sp[l] - x[j][i].sp[l]) * dhy; 226c4762a1bSJed Brown sym = syp; 227c4762a1bSJed Brown } else if (j == info->my - 1) { 228c4762a1bSJed Brown syp = (x[j][i].sp[l] - x[j - 1][i].sp[l]) * dhy; 229c4762a1bSJed Brown sym = syp; 230c4762a1bSJed Brown } else { 231c4762a1bSJed Brown sym = (x[j + 1][i].sp[l] - x[j][i].sp[l]) * dhy; 232c4762a1bSJed Brown syp = (x[j][i].sp[l] - x[j - 1][i].sp[l]) * dhy; 233c4762a1bSJed Brown } /* 4 flops per species*point */ 234c4762a1bSJed Brown 235c4762a1bSJed Brown if (i == 0) { 236c4762a1bSJed Brown if (l == 0) f[j][i].sp[l] = (x[j][i].sp[l] - ctx->x_inflow.sp[l] * ((PetscScalar)j) / (info->my - 1)); 237c4762a1bSJed 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))); 238c4762a1bSJed Brown else f[j][i].sp[l] = x[j][i].sp[l]; 239c4762a1bSJed Brown 240c4762a1bSJed Brown } else { 241c4762a1bSJed Brown f[j][i].sp[l] = xt[j][i].sp[l] * scale; 242c4762a1bSJed Brown u = x[j][i].sp[l]; 243c4762a1bSJed Brown if (j == 0) uyy = u - x[j + 1][i].sp[l]; 244c4762a1bSJed Brown else if (j == info->my - 1) uyy = u - x[j - 1][i].sp[l]; 245c4762a1bSJed Brown else uyy = (2.0 * u - x[j - 1][i].sp[l] - x[j + 1][i].sp[l]) * hxdhy; 246c4762a1bSJed Brown 247c4762a1bSJed Brown if (i != info->mx - 1) uxx = (2.0 * u - x[j][i - 1].sp[l] - x[j][i + 1].sp[l]) * hydhx; 248c4762a1bSJed Brown else uxx = u - x[j][i - 1].sp[l]; 249c4762a1bSJed Brown /* 10 flops per species*point */ 250c4762a1bSJed Brown 251c4762a1bSJed Brown f_advect = 0.; 252c4762a1bSJed Brown f_advect += scale * (vxp * sxp + vxm * sxm); 253c4762a1bSJed Brown f_advect += scale * (vyp * syp + vym * sym); 254c4762a1bSJed Brown f[j][i].sp[l] += f_advect + ctx->dispersivity * (uxx + uyy); 255c4762a1bSJed Brown /* 14 flops per species*point */ 256c4762a1bSJed Brown } 257c4762a1bSJed Brown } 258c4762a1bSJed Brown /* reaction */ 259c4762a1bSJed Brown if (i != 0) { 260c4762a1bSJed Brown for (m = 0; m < N_REACTIONS; m++) { 261c4762a1bSJed Brown rate = ctx->rate_constant[m]; 262c4762a1bSJed Brown for (l = 0; l < N_SPECIES; l++) { 263c4762a1bSJed Brown if (stoich(m, l) < 0) { 264c4762a1bSJed Brown /* assume an elementary reaction */ 265c4762a1bSJed Brown rate *= PetscPowScalar(x[j][i].sp[l], PetscAbsScalar(stoich(m, l))); 266c4762a1bSJed Brown /* ~10 flops per reaction per species per point */ 267c4762a1bSJed Brown } 268c4762a1bSJed Brown } 269c4762a1bSJed Brown for (l = 0; l < N_SPECIES; l++) { 270c4762a1bSJed Brown f[j][i].sp[l] += -scale * stoich(m, l) * rate; /* Reaction term */ 271c4762a1bSJed Brown /* 3 flops per reaction per species per point */ 272c4762a1bSJed Brown } 273c4762a1bSJed Brown } 274c4762a1bSJed Brown } 275c4762a1bSJed Brown } 276c4762a1bSJed Brown } 2779566063dSJacob Faibussowitsch PetscCall(PetscLogFlops((N_SPECIES * (28.0 + 13.0 * N_REACTIONS)) * info->ym * info->xm)); 278c4762a1bSJed Brown PetscFunctionReturn(0); 279c4762a1bSJed Brown } 280c4762a1bSJed Brown 281*9371c9d4SSatish Balay PetscErrorCode ReactingFlowPostCheck(SNESLineSearch linesearch, Vec X, Vec Y, Vec W, PetscBool *changed_y, PetscBool *changed_w, void *vctx) { 282c4762a1bSJed Brown PetscInt i, j, l, Mx, My, xs, ys, xm, ym; 283c4762a1bSJed Brown Field **x; 284c4762a1bSJed Brown SNES snes; 285c4762a1bSJed Brown DM da; 286c4762a1bSJed Brown PetscScalar min; 287c4762a1bSJed Brown 288c4762a1bSJed Brown PetscFunctionBeginUser; 289c4762a1bSJed Brown *changed_w = PETSC_FALSE; 2909566063dSJacob Faibussowitsch PetscCall(VecMin(X, NULL, &min)); 291c4762a1bSJed Brown if (min >= 0.) PetscFunctionReturn(0); 292c4762a1bSJed Brown 293c4762a1bSJed Brown *changed_w = PETSC_TRUE; 2949566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetSNES(linesearch, &snes)); 2959566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &da)); 296*9371c9d4SSatish 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)); 2979566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, W, &x)); 2989566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL)); 299c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) { 300c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 301c4762a1bSJed Brown for (l = 0; l < N_SPECIES; l++) { 302c4762a1bSJed Brown if (x[j][i].sp[l] < 0.) x[j][i].sp[l] = 0.; 303c4762a1bSJed Brown } 304c4762a1bSJed Brown } 305c4762a1bSJed Brown } 3069566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, W, &x)); 307c4762a1bSJed Brown PetscFunctionReturn(0); 308c4762a1bSJed Brown } 309c4762a1bSJed Brown 310*9371c9d4SSatish Balay PetscErrorCode FormIFunction(TS ts, PetscReal ptime, Vec X, Vec Xdot, Vec F, void *user) { 311c4762a1bSJed Brown DMDALocalInfo info; 312c4762a1bSJed Brown Field **u, **udot, **fu; 313c4762a1bSJed Brown Vec localX; 314c4762a1bSJed Brown DM da; 315c4762a1bSJed Brown 316c4762a1bSJed Brown PetscFunctionBeginUser; 3179566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, (DM *)&da)); 3189566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &localX)); 3199566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, localX)); 3209566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, localX)); 3219566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(da, &info)); 3229566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, localX, &u)); 3239566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, Xdot, &udot)); 3249566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, F, &fu)); 3259566063dSJacob Faibussowitsch PetscCall(FormIFunctionLocal(&info, ptime, u, udot, fu, (AppCtx *)user)); 3269566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, localX, &u)); 3279566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, Xdot, &udot)); 3289566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, F, &fu)); 3299566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &localX)); 330c4762a1bSJed Brown PetscFunctionReturn(0); 331c4762a1bSJed Brown } 332