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