1e77315bcSPatrick Sanan 2e77315bcSPatrick Sanan static char help[] = "Nonlinear Radiative Transport PDE with multigrid in 3d.\n\ 3e77315bcSPatrick Sanan Uses 3-dimensional distributed arrays.\n\ 4e77315bcSPatrick Sanan A 3-dim simplified Radiative Transport test problem is used, with analytic Jacobian. \n\ 5e77315bcSPatrick Sanan \n\ 6e77315bcSPatrick Sanan Solves the linear systems via multilevel methods \n\ 7e77315bcSPatrick Sanan \n\ 8e77315bcSPatrick Sanan The command line\n\ 9e77315bcSPatrick Sanan options are:\n\ 10e77315bcSPatrick Sanan -tleft <tl>, where <tl> indicates the left Diriclet BC \n\ 11e77315bcSPatrick Sanan -tright <tr>, where <tr> indicates the right Diriclet BC \n\ 12e77315bcSPatrick Sanan -beta <beta>, where <beta> indicates the exponent in T \n\n"; 13e77315bcSPatrick Sanan 14e77315bcSPatrick Sanan /* 15e77315bcSPatrick Sanan 16e77315bcSPatrick Sanan This example models the partial differential equation 17e77315bcSPatrick Sanan 18e77315bcSPatrick Sanan - Div(alpha* T^beta (GRAD T)) = 0. 19e77315bcSPatrick Sanan 20e77315bcSPatrick Sanan where beta = 2.5 and alpha = 1.0 21e77315bcSPatrick Sanan 22e77315bcSPatrick Sanan BC: T_left = 1.0, T_right = 0.1, dT/dn_top = dTdn_bottom = dT/dn_up = dT/dn_down = 0. 23e77315bcSPatrick Sanan 24e77315bcSPatrick Sanan in the unit square, which is uniformly discretized in each of x and 25e77315bcSPatrick Sanan y in this simple encoding. The degrees of freedom are cell centered. 26e77315bcSPatrick Sanan 27e77315bcSPatrick Sanan A finite volume approximation with the usual 7-point stencil 28e77315bcSPatrick Sanan is used to discretize the boundary value problem to obtain a 29e77315bcSPatrick Sanan nonlinear system of equations. 30e77315bcSPatrick Sanan 31e77315bcSPatrick Sanan This code was contributed by Nickolas Jovanovic based on ex18.c 32e77315bcSPatrick Sanan 33e77315bcSPatrick Sanan */ 34e77315bcSPatrick Sanan 35e77315bcSPatrick Sanan #include <petscsnes.h> 36e77315bcSPatrick Sanan #include <petscdm.h> 37e77315bcSPatrick Sanan #include <petscdmda.h> 38e77315bcSPatrick Sanan 39e77315bcSPatrick Sanan /* User-defined application context */ 40e77315bcSPatrick Sanan 41e77315bcSPatrick Sanan typedef struct { 42e77315bcSPatrick Sanan PetscReal tleft, tright; /* Dirichlet boundary conditions */ 43e77315bcSPatrick Sanan PetscReal beta, bm1, coef; /* nonlinear diffusivity parameterizations */ 44e77315bcSPatrick Sanan } AppCtx; 45e77315bcSPatrick Sanan 46e77315bcSPatrick Sanan #define POWFLOP 5 /* assume a pow() takes five flops */ 47e77315bcSPatrick Sanan 48e77315bcSPatrick Sanan extern PetscErrorCode FormInitialGuess(SNES, Vec, void *); 49e77315bcSPatrick Sanan extern PetscErrorCode FormFunction(SNES, Vec, Vec, void *); 50e77315bcSPatrick Sanan extern PetscErrorCode FormJacobian(SNES, Vec, Mat, Mat, void *); 51e77315bcSPatrick Sanan 529371c9d4SSatish Balay int main(int argc, char **argv) { 53e77315bcSPatrick Sanan SNES snes; 54e77315bcSPatrick Sanan AppCtx user; 55e77315bcSPatrick Sanan PetscInt its, lits; 56e77315bcSPatrick Sanan PetscReal litspit; 57e77315bcSPatrick Sanan DM da; 58e77315bcSPatrick Sanan 59327415f7SBarry Smith PetscFunctionBeginUser; 609566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 61e77315bcSPatrick Sanan /* set problem parameters */ 62e77315bcSPatrick Sanan user.tleft = 1.0; 63e77315bcSPatrick Sanan user.tright = 0.1; 64e77315bcSPatrick Sanan user.beta = 2.5; 659566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-tleft", &user.tleft, NULL)); 669566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-tright", &user.tright, NULL)); 679566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-beta", &user.beta, NULL)); 68e77315bcSPatrick Sanan user.bm1 = user.beta - 1.0; 69e77315bcSPatrick Sanan user.coef = user.beta / 2.0; 70e77315bcSPatrick Sanan 71e77315bcSPatrick Sanan /* 72e77315bcSPatrick Sanan Set the DMDA (grid structure) for the grids. 73e77315bcSPatrick Sanan */ 749566063dSJacob Faibussowitsch PetscCall(DMDACreate3d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_STAR, 5, 5, 5, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE, 1, 1, 0, 0, 0, &da)); 759566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da)); 769566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 779566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(da, &user)); 78e77315bcSPatrick Sanan 79e77315bcSPatrick Sanan /* 80e77315bcSPatrick Sanan Create the nonlinear solver 81e77315bcSPatrick Sanan */ 829566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 839566063dSJacob Faibussowitsch PetscCall(SNESSetDM(snes, da)); 849566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes, NULL, FormFunction, &user)); 859566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, NULL, NULL, FormJacobian, &user)); 869566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 879566063dSJacob Faibussowitsch PetscCall(SNESSetComputeInitialGuess(snes, FormInitialGuess, NULL)); 88e77315bcSPatrick Sanan 899566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, NULL)); 909566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(snes, &its)); 919566063dSJacob Faibussowitsch PetscCall(SNESGetLinearSolveIterations(snes, &lits)); 92e77315bcSPatrick Sanan litspit = ((PetscReal)lits) / ((PetscReal)its); 9363a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Number of SNES iterations = %" PetscInt_FMT "\n", its)); 9463a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Number of Linear iterations = %" PetscInt_FMT "\n", lits)); 959566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Average Linear its / SNES = %e\n", (double)litspit)); 96e77315bcSPatrick Sanan 979566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 989566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da)); 999566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 100b122ec5aSJacob Faibussowitsch return 0; 101e77315bcSPatrick Sanan } 102e77315bcSPatrick Sanan /* -------------------- Form initial approximation ----------------- */ 1039371c9d4SSatish Balay PetscErrorCode FormInitialGuess(SNES snes, Vec X, void *ctx) { 104e77315bcSPatrick Sanan AppCtx *user; 105e77315bcSPatrick Sanan PetscInt i, j, k, xs, ys, xm, ym, zs, zm; 106e77315bcSPatrick Sanan PetscScalar ***x; 107e77315bcSPatrick Sanan DM da; 108e77315bcSPatrick Sanan 109e77315bcSPatrick Sanan PetscFunctionBeginUser; 1109566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &da)); 1119566063dSJacob Faibussowitsch PetscCall(DMGetApplicationContext(da, &user)); 1129566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm)); 1139566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, X, &x)); 114e77315bcSPatrick Sanan 115e77315bcSPatrick Sanan /* Compute initial guess */ 116e77315bcSPatrick Sanan for (k = zs; k < zs + zm; k++) { 117e77315bcSPatrick Sanan for (j = ys; j < ys + ym; j++) { 118*ad540459SPierre Jolivet for (i = xs; i < xs + xm; i++) x[k][j][i] = user->tleft; 119e77315bcSPatrick Sanan } 120e77315bcSPatrick Sanan } 1219566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, X, &x)); 122e77315bcSPatrick Sanan PetscFunctionReturn(0); 123e77315bcSPatrick Sanan } 124e77315bcSPatrick Sanan /* -------------------- Evaluate Function F(x) --------------------- */ 1259371c9d4SSatish Balay PetscErrorCode FormFunction(SNES snes, Vec X, Vec F, void *ptr) { 126e77315bcSPatrick Sanan AppCtx *user = (AppCtx *)ptr; 127e77315bcSPatrick Sanan PetscInt i, j, k, mx, my, mz, xs, ys, zs, xm, ym, zm; 128e77315bcSPatrick Sanan PetscScalar zero = 0.0, one = 1.0; 129e77315bcSPatrick Sanan PetscScalar hx, hy, hz, hxhydhz, hyhzdhx, hzhxdhy; 130e77315bcSPatrick Sanan PetscScalar t0, tn, ts, te, tw, an, as, ae, aw, dn, ds, de, dw, fn = 0.0, fs = 0.0, fe = 0.0, fw = 0.0; 131e77315bcSPatrick Sanan PetscScalar tleft, tright, beta, td, ad, dd, fd = 0.0, tu, au, du = 0.0, fu = 0.0; 132e77315bcSPatrick Sanan PetscScalar ***x, ***f; 133e77315bcSPatrick Sanan Vec localX; 134e77315bcSPatrick Sanan DM da; 135e77315bcSPatrick Sanan 136e77315bcSPatrick Sanan PetscFunctionBeginUser; 1379566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &da)); 1389566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &localX)); 1399566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, NULL, &mx, &my, &mz, 0, 0, 0, 0, 0, 0, 0, 0, 0)); 1409371c9d4SSatish Balay hx = one / (PetscReal)(mx - 1); 1419371c9d4SSatish Balay hy = one / (PetscReal)(my - 1); 1429371c9d4SSatish Balay hz = one / (PetscReal)(mz - 1); 1439371c9d4SSatish Balay hxhydhz = hx * hy / hz; 1449371c9d4SSatish Balay hyhzdhx = hy * hz / hx; 1459371c9d4SSatish Balay hzhxdhy = hz * hx / hy; 1469371c9d4SSatish Balay tleft = user->tleft; 1479371c9d4SSatish Balay tright = user->tright; 148e77315bcSPatrick Sanan beta = user->beta; 149e77315bcSPatrick Sanan 150e77315bcSPatrick Sanan /* Get ghost points */ 1519566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, localX)); 1529566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, localX)); 1539566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm)); 1549566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, localX, &x)); 1559566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, F, &f)); 156e77315bcSPatrick Sanan 157e77315bcSPatrick Sanan /* Evaluate function */ 158e77315bcSPatrick Sanan for (k = zs; k < zs + zm; k++) { 159e77315bcSPatrick Sanan for (j = ys; j < ys + ym; j++) { 160e77315bcSPatrick Sanan for (i = xs; i < xs + xm; i++) { 161e77315bcSPatrick Sanan t0 = x[k][j][i]; 162e77315bcSPatrick Sanan 163e77315bcSPatrick Sanan if (i > 0 && i < mx - 1 && j > 0 && j < my - 1 && k > 0 && k < mz - 1) { 164e77315bcSPatrick Sanan /* general interior volume */ 165e77315bcSPatrick Sanan 166e77315bcSPatrick Sanan tw = x[k][j][i - 1]; 167e77315bcSPatrick Sanan aw = 0.5 * (t0 + tw); 168e77315bcSPatrick Sanan dw = PetscPowScalar(aw, beta); 169e77315bcSPatrick Sanan fw = dw * (t0 - tw); 170e77315bcSPatrick Sanan 171e77315bcSPatrick Sanan te = x[k][j][i + 1]; 172e77315bcSPatrick Sanan ae = 0.5 * (t0 + te); 173e77315bcSPatrick Sanan de = PetscPowScalar(ae, beta); 174e77315bcSPatrick Sanan fe = de * (te - t0); 175e77315bcSPatrick Sanan 176e77315bcSPatrick Sanan ts = x[k][j - 1][i]; 177e77315bcSPatrick Sanan as = 0.5 * (t0 + ts); 178e77315bcSPatrick Sanan ds = PetscPowScalar(as, beta); 179e77315bcSPatrick Sanan fs = ds * (t0 - ts); 180e77315bcSPatrick Sanan 181e77315bcSPatrick Sanan tn = x[k][j + 1][i]; 182e77315bcSPatrick Sanan an = 0.5 * (t0 + tn); 183e77315bcSPatrick Sanan dn = PetscPowScalar(an, beta); 184e77315bcSPatrick Sanan fn = dn * (tn - t0); 185e77315bcSPatrick Sanan 186e77315bcSPatrick Sanan td = x[k - 1][j][i]; 187e77315bcSPatrick Sanan ad = 0.5 * (t0 + td); 188e77315bcSPatrick Sanan dd = PetscPowScalar(ad, beta); 189e77315bcSPatrick Sanan fd = dd * (t0 - td); 190e77315bcSPatrick Sanan 191e77315bcSPatrick Sanan tu = x[k + 1][j][i]; 192e77315bcSPatrick Sanan au = 0.5 * (t0 + tu); 193e77315bcSPatrick Sanan du = PetscPowScalar(au, beta); 194e77315bcSPatrick Sanan fu = du * (tu - t0); 195e77315bcSPatrick Sanan 196e77315bcSPatrick Sanan } else if (i == 0) { 197e77315bcSPatrick Sanan /* left-hand (west) boundary */ 198e77315bcSPatrick Sanan tw = tleft; 199e77315bcSPatrick Sanan aw = 0.5 * (t0 + tw); 200e77315bcSPatrick Sanan dw = PetscPowScalar(aw, beta); 201e77315bcSPatrick Sanan fw = dw * (t0 - tw); 202e77315bcSPatrick Sanan 203e77315bcSPatrick Sanan te = x[k][j][i + 1]; 204e77315bcSPatrick Sanan ae = 0.5 * (t0 + te); 205e77315bcSPatrick Sanan de = PetscPowScalar(ae, beta); 206e77315bcSPatrick Sanan fe = de * (te - t0); 207e77315bcSPatrick Sanan 208e77315bcSPatrick Sanan if (j > 0) { 209e77315bcSPatrick Sanan ts = x[k][j - 1][i]; 210e77315bcSPatrick Sanan as = 0.5 * (t0 + ts); 211e77315bcSPatrick Sanan ds = PetscPowScalar(as, beta); 212e77315bcSPatrick Sanan fs = ds * (t0 - ts); 213e77315bcSPatrick Sanan } else { 214e77315bcSPatrick Sanan fs = zero; 215e77315bcSPatrick Sanan } 216e77315bcSPatrick Sanan 217e77315bcSPatrick Sanan if (j < my - 1) { 218e77315bcSPatrick Sanan tn = x[k][j + 1][i]; 219e77315bcSPatrick Sanan an = 0.5 * (t0 + tn); 220e77315bcSPatrick Sanan dn = PetscPowScalar(an, beta); 221e77315bcSPatrick Sanan fn = dn * (tn - t0); 222e77315bcSPatrick Sanan } else { 223e77315bcSPatrick Sanan fn = zero; 224e77315bcSPatrick Sanan } 225e77315bcSPatrick Sanan 226e77315bcSPatrick Sanan if (k > 0) { 227e77315bcSPatrick Sanan td = x[k - 1][j][i]; 228e77315bcSPatrick Sanan ad = 0.5 * (t0 + td); 229e77315bcSPatrick Sanan dd = PetscPowScalar(ad, beta); 230e77315bcSPatrick Sanan fd = dd * (t0 - td); 231e77315bcSPatrick Sanan } else { 232e77315bcSPatrick Sanan fd = zero; 233e77315bcSPatrick Sanan } 234e77315bcSPatrick Sanan 235e77315bcSPatrick Sanan if (k < mz - 1) { 236e77315bcSPatrick Sanan tu = x[k + 1][j][i]; 237e77315bcSPatrick Sanan au = 0.5 * (t0 + tu); 238e77315bcSPatrick Sanan du = PetscPowScalar(au, beta); 239e77315bcSPatrick Sanan fu = du * (tu - t0); 240e77315bcSPatrick Sanan } else { 241e77315bcSPatrick Sanan fu = zero; 242e77315bcSPatrick Sanan } 243e77315bcSPatrick Sanan 244e77315bcSPatrick Sanan } else if (i == mx - 1) { 245e77315bcSPatrick Sanan /* right-hand (east) boundary */ 246e77315bcSPatrick Sanan tw = x[k][j][i - 1]; 247e77315bcSPatrick Sanan aw = 0.5 * (t0 + tw); 248e77315bcSPatrick Sanan dw = PetscPowScalar(aw, beta); 249e77315bcSPatrick Sanan fw = dw * (t0 - tw); 250e77315bcSPatrick Sanan 251e77315bcSPatrick Sanan te = tright; 252e77315bcSPatrick Sanan ae = 0.5 * (t0 + te); 253e77315bcSPatrick Sanan de = PetscPowScalar(ae, beta); 254e77315bcSPatrick Sanan fe = de * (te - t0); 255e77315bcSPatrick Sanan 256e77315bcSPatrick Sanan if (j > 0) { 257e77315bcSPatrick Sanan ts = x[k][j - 1][i]; 258e77315bcSPatrick Sanan as = 0.5 * (t0 + ts); 259e77315bcSPatrick Sanan ds = PetscPowScalar(as, beta); 260e77315bcSPatrick Sanan fs = ds * (t0 - ts); 261e77315bcSPatrick Sanan } else { 262e77315bcSPatrick Sanan fs = zero; 263e77315bcSPatrick Sanan } 264e77315bcSPatrick Sanan 265e77315bcSPatrick Sanan if (j < my - 1) { 266e77315bcSPatrick Sanan tn = x[k][j + 1][i]; 267e77315bcSPatrick Sanan an = 0.5 * (t0 + tn); 268e77315bcSPatrick Sanan dn = PetscPowScalar(an, beta); 269e77315bcSPatrick Sanan fn = dn * (tn - t0); 270e77315bcSPatrick Sanan } else { 271e77315bcSPatrick Sanan fn = zero; 272e77315bcSPatrick Sanan } 273e77315bcSPatrick Sanan 274e77315bcSPatrick Sanan if (k > 0) { 275e77315bcSPatrick Sanan td = x[k - 1][j][i]; 276e77315bcSPatrick Sanan ad = 0.5 * (t0 + td); 277e77315bcSPatrick Sanan dd = PetscPowScalar(ad, beta); 278e77315bcSPatrick Sanan fd = dd * (t0 - td); 279e77315bcSPatrick Sanan } else { 280e77315bcSPatrick Sanan fd = zero; 281e77315bcSPatrick Sanan } 282e77315bcSPatrick Sanan 283e77315bcSPatrick Sanan if (k < mz - 1) { 284e77315bcSPatrick Sanan tu = x[k + 1][j][i]; 285e77315bcSPatrick Sanan au = 0.5 * (t0 + tu); 286e77315bcSPatrick Sanan du = PetscPowScalar(au, beta); 287e77315bcSPatrick Sanan fu = du * (tu - t0); 288e77315bcSPatrick Sanan } else { 289e77315bcSPatrick Sanan fu = zero; 290e77315bcSPatrick Sanan } 291e77315bcSPatrick Sanan 292e77315bcSPatrick Sanan } else if (j == 0) { 293e77315bcSPatrick Sanan /* bottom (south) boundary, and i <> 0 or mx-1 */ 294e77315bcSPatrick Sanan tw = x[k][j][i - 1]; 295e77315bcSPatrick Sanan aw = 0.5 * (t0 + tw); 296e77315bcSPatrick Sanan dw = PetscPowScalar(aw, beta); 297e77315bcSPatrick Sanan fw = dw * (t0 - tw); 298e77315bcSPatrick Sanan 299e77315bcSPatrick Sanan te = x[k][j][i + 1]; 300e77315bcSPatrick Sanan ae = 0.5 * (t0 + te); 301e77315bcSPatrick Sanan de = PetscPowScalar(ae, beta); 302e77315bcSPatrick Sanan fe = de * (te - t0); 303e77315bcSPatrick Sanan 304e77315bcSPatrick Sanan fs = zero; 305e77315bcSPatrick Sanan 306e77315bcSPatrick Sanan tn = x[k][j + 1][i]; 307e77315bcSPatrick Sanan an = 0.5 * (t0 + tn); 308e77315bcSPatrick Sanan dn = PetscPowScalar(an, beta); 309e77315bcSPatrick Sanan fn = dn * (tn - t0); 310e77315bcSPatrick Sanan 311e77315bcSPatrick Sanan if (k > 0) { 312e77315bcSPatrick Sanan td = x[k - 1][j][i]; 313e77315bcSPatrick Sanan ad = 0.5 * (t0 + td); 314e77315bcSPatrick Sanan dd = PetscPowScalar(ad, beta); 315e77315bcSPatrick Sanan fd = dd * (t0 - td); 316e77315bcSPatrick Sanan } else { 317e77315bcSPatrick Sanan fd = zero; 318e77315bcSPatrick Sanan } 319e77315bcSPatrick Sanan 320e77315bcSPatrick Sanan if (k < mz - 1) { 321e77315bcSPatrick Sanan tu = x[k + 1][j][i]; 322e77315bcSPatrick Sanan au = 0.5 * (t0 + tu); 323e77315bcSPatrick Sanan du = PetscPowScalar(au, beta); 324e77315bcSPatrick Sanan fu = du * (tu - t0); 325e77315bcSPatrick Sanan } else { 326e77315bcSPatrick Sanan fu = zero; 327e77315bcSPatrick Sanan } 328e77315bcSPatrick Sanan 329e77315bcSPatrick Sanan } else if (j == my - 1) { 330e77315bcSPatrick Sanan /* top (north) boundary, and i <> 0 or mx-1 */ 331e77315bcSPatrick Sanan tw = x[k][j][i - 1]; 332e77315bcSPatrick Sanan aw = 0.5 * (t0 + tw); 333e77315bcSPatrick Sanan dw = PetscPowScalar(aw, beta); 334e77315bcSPatrick Sanan fw = dw * (t0 - tw); 335e77315bcSPatrick Sanan 336e77315bcSPatrick Sanan te = x[k][j][i + 1]; 337e77315bcSPatrick Sanan ae = 0.5 * (t0 + te); 338e77315bcSPatrick Sanan de = PetscPowScalar(ae, beta); 339e77315bcSPatrick Sanan fe = de * (te - t0); 340e77315bcSPatrick Sanan 341e77315bcSPatrick Sanan ts = x[k][j - 1][i]; 342e77315bcSPatrick Sanan as = 0.5 * (t0 + ts); 343e77315bcSPatrick Sanan ds = PetscPowScalar(as, beta); 344e77315bcSPatrick Sanan fs = ds * (t0 - ts); 345e77315bcSPatrick Sanan 346e77315bcSPatrick Sanan fn = zero; 347e77315bcSPatrick Sanan 348e77315bcSPatrick Sanan if (k > 0) { 349e77315bcSPatrick Sanan td = x[k - 1][j][i]; 350e77315bcSPatrick Sanan ad = 0.5 * (t0 + td); 351e77315bcSPatrick Sanan dd = PetscPowScalar(ad, beta); 352e77315bcSPatrick Sanan fd = dd * (t0 - td); 353e77315bcSPatrick Sanan } else { 354e77315bcSPatrick Sanan fd = zero; 355e77315bcSPatrick Sanan } 356e77315bcSPatrick Sanan 357e77315bcSPatrick Sanan if (k < mz - 1) { 358e77315bcSPatrick Sanan tu = x[k + 1][j][i]; 359e77315bcSPatrick Sanan au = 0.5 * (t0 + tu); 360e77315bcSPatrick Sanan du = PetscPowScalar(au, beta); 361e77315bcSPatrick Sanan fu = du * (tu - t0); 362e77315bcSPatrick Sanan } else { 363e77315bcSPatrick Sanan fu = zero; 364e77315bcSPatrick Sanan } 365e77315bcSPatrick Sanan 366e77315bcSPatrick Sanan } else if (k == 0) { 367e77315bcSPatrick Sanan /* down boundary (interior only) */ 368e77315bcSPatrick Sanan tw = x[k][j][i - 1]; 369e77315bcSPatrick Sanan aw = 0.5 * (t0 + tw); 370e77315bcSPatrick Sanan dw = PetscPowScalar(aw, beta); 371e77315bcSPatrick Sanan fw = dw * (t0 - tw); 372e77315bcSPatrick Sanan 373e77315bcSPatrick Sanan te = x[k][j][i + 1]; 374e77315bcSPatrick Sanan ae = 0.5 * (t0 + te); 375e77315bcSPatrick Sanan de = PetscPowScalar(ae, beta); 376e77315bcSPatrick Sanan fe = de * (te - t0); 377e77315bcSPatrick Sanan 378e77315bcSPatrick Sanan ts = x[k][j - 1][i]; 379e77315bcSPatrick Sanan as = 0.5 * (t0 + ts); 380e77315bcSPatrick Sanan ds = PetscPowScalar(as, beta); 381e77315bcSPatrick Sanan fs = ds * (t0 - ts); 382e77315bcSPatrick Sanan 383e77315bcSPatrick Sanan tn = x[k][j + 1][i]; 384e77315bcSPatrick Sanan an = 0.5 * (t0 + tn); 385e77315bcSPatrick Sanan dn = PetscPowScalar(an, beta); 386e77315bcSPatrick Sanan fn = dn * (tn - t0); 387e77315bcSPatrick Sanan 388e77315bcSPatrick Sanan fd = zero; 389e77315bcSPatrick Sanan 390e77315bcSPatrick Sanan tu = x[k + 1][j][i]; 391e77315bcSPatrick Sanan au = 0.5 * (t0 + tu); 392e77315bcSPatrick Sanan du = PetscPowScalar(au, beta); 393e77315bcSPatrick Sanan fu = du * (tu - t0); 394e77315bcSPatrick Sanan 395e77315bcSPatrick Sanan } else if (k == mz - 1) { 396e77315bcSPatrick Sanan /* up boundary (interior only) */ 397e77315bcSPatrick Sanan tw = x[k][j][i - 1]; 398e77315bcSPatrick Sanan aw = 0.5 * (t0 + tw); 399e77315bcSPatrick Sanan dw = PetscPowScalar(aw, beta); 400e77315bcSPatrick Sanan fw = dw * (t0 - tw); 401e77315bcSPatrick Sanan 402e77315bcSPatrick Sanan te = x[k][j][i + 1]; 403e77315bcSPatrick Sanan ae = 0.5 * (t0 + te); 404e77315bcSPatrick Sanan de = PetscPowScalar(ae, beta); 405e77315bcSPatrick Sanan fe = de * (te - t0); 406e77315bcSPatrick Sanan 407e77315bcSPatrick Sanan ts = x[k][j - 1][i]; 408e77315bcSPatrick Sanan as = 0.5 * (t0 + ts); 409e77315bcSPatrick Sanan ds = PetscPowScalar(as, beta); 410e77315bcSPatrick Sanan fs = ds * (t0 - ts); 411e77315bcSPatrick Sanan 412e77315bcSPatrick Sanan tn = x[k][j + 1][i]; 413e77315bcSPatrick Sanan an = 0.5 * (t0 + tn); 414e77315bcSPatrick Sanan dn = PetscPowScalar(an, beta); 415e77315bcSPatrick Sanan fn = dn * (tn - t0); 416e77315bcSPatrick Sanan 417e77315bcSPatrick Sanan td = x[k - 1][j][i]; 418e77315bcSPatrick Sanan ad = 0.5 * (t0 + td); 419e77315bcSPatrick Sanan dd = PetscPowScalar(ad, beta); 420e77315bcSPatrick Sanan fd = dd * (t0 - td); 421e77315bcSPatrick Sanan 422e77315bcSPatrick Sanan fu = zero; 423e77315bcSPatrick Sanan } 424e77315bcSPatrick Sanan 425e77315bcSPatrick Sanan f[k][j][i] = -hyhzdhx * (fe - fw) - hzhxdhy * (fn - fs) - hxhydhz * (fu - fd); 426e77315bcSPatrick Sanan } 427e77315bcSPatrick Sanan } 428e77315bcSPatrick Sanan } 4299566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, localX, &x)); 4309566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, F, &f)); 4319566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &localX)); 4329566063dSJacob Faibussowitsch PetscCall(PetscLogFlops((22.0 + 4.0 * POWFLOP) * ym * xm)); 433e77315bcSPatrick Sanan PetscFunctionReturn(0); 434e77315bcSPatrick Sanan } 435e77315bcSPatrick Sanan /* -------------------- Evaluate Jacobian F(x) --------------------- */ 4369371c9d4SSatish Balay PetscErrorCode FormJacobian(SNES snes, Vec X, Mat J, Mat jac, void *ptr) { 437e77315bcSPatrick Sanan AppCtx *user = (AppCtx *)ptr; 438e77315bcSPatrick Sanan PetscInt i, j, k, mx, my, mz, xs, ys, zs, xm, ym, zm; 439e77315bcSPatrick Sanan PetscScalar one = 1.0; 440e77315bcSPatrick Sanan PetscScalar hx, hy, hz, hxhydhz, hyhzdhx, hzhxdhy; 441e77315bcSPatrick Sanan PetscScalar t0, tn, ts, te, tw, an, as, ae, aw, dn, ds, de, dw; 442e77315bcSPatrick Sanan PetscScalar tleft, tright, beta, td, ad, dd, tu, au, du, v[7], bm1, coef; 443e77315bcSPatrick Sanan PetscScalar ***x, bn, bs, be, bw, bu, bd, gn, gs, ge, gw, gu, gd; 444e77315bcSPatrick Sanan Vec localX; 445e77315bcSPatrick Sanan MatStencil c[7], row; 446e77315bcSPatrick Sanan DM da; 447e77315bcSPatrick Sanan 448e77315bcSPatrick Sanan PetscFunctionBeginUser; 4499566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &da)); 4509566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &localX)); 4519566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, NULL, &mx, &my, &mz, 0, 0, 0, 0, 0, 0, 0, 0, 0)); 4529371c9d4SSatish Balay hx = one / (PetscReal)(mx - 1); 4539371c9d4SSatish Balay hy = one / (PetscReal)(my - 1); 4549371c9d4SSatish Balay hz = one / (PetscReal)(mz - 1); 4559371c9d4SSatish Balay hxhydhz = hx * hy / hz; 4569371c9d4SSatish Balay hyhzdhx = hy * hz / hx; 4579371c9d4SSatish Balay hzhxdhy = hz * hx / hy; 4589371c9d4SSatish Balay tleft = user->tleft; 4599371c9d4SSatish Balay tright = user->tright; 4609371c9d4SSatish Balay beta = user->beta; 4619371c9d4SSatish Balay bm1 = user->bm1; 4629371c9d4SSatish Balay coef = user->coef; 463e77315bcSPatrick Sanan 464e77315bcSPatrick Sanan /* Get ghost points */ 4659566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, localX)); 4669566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, localX)); 4679566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm)); 4689566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, localX, &x)); 469e77315bcSPatrick Sanan 470e77315bcSPatrick Sanan /* Evaluate Jacobian of function */ 471e77315bcSPatrick Sanan for (k = zs; k < zs + zm; k++) { 472e77315bcSPatrick Sanan for (j = ys; j < ys + ym; j++) { 473e77315bcSPatrick Sanan for (i = xs; i < xs + xm; i++) { 474e77315bcSPatrick Sanan t0 = x[k][j][i]; 4759371c9d4SSatish Balay row.k = k; 4769371c9d4SSatish Balay row.j = j; 4779371c9d4SSatish Balay row.i = i; 478e77315bcSPatrick Sanan if (i > 0 && i < mx - 1 && j > 0 && j < my - 1 && k > 0 && k < mz - 1) { 479e77315bcSPatrick Sanan /* general interior volume */ 480e77315bcSPatrick Sanan 481e77315bcSPatrick Sanan tw = x[k][j][i - 1]; 482e77315bcSPatrick Sanan aw = 0.5 * (t0 + tw); 483e77315bcSPatrick Sanan bw = PetscPowScalar(aw, bm1); 484e77315bcSPatrick Sanan /* dw = bw * aw */ 485e77315bcSPatrick Sanan dw = PetscPowScalar(aw, beta); 486e77315bcSPatrick Sanan gw = coef * bw * (t0 - tw); 487e77315bcSPatrick Sanan 488e77315bcSPatrick Sanan te = x[k][j][i + 1]; 489e77315bcSPatrick Sanan ae = 0.5 * (t0 + te); 490e77315bcSPatrick Sanan be = PetscPowScalar(ae, bm1); 491e77315bcSPatrick Sanan /* de = be * ae; */ 492e77315bcSPatrick Sanan de = PetscPowScalar(ae, beta); 493e77315bcSPatrick Sanan ge = coef * be * (te - t0); 494e77315bcSPatrick Sanan 495e77315bcSPatrick Sanan ts = x[k][j - 1][i]; 496e77315bcSPatrick Sanan as = 0.5 * (t0 + ts); 497e77315bcSPatrick Sanan bs = PetscPowScalar(as, bm1); 498e77315bcSPatrick Sanan /* ds = bs * as; */ 499e77315bcSPatrick Sanan ds = PetscPowScalar(as, beta); 500e77315bcSPatrick Sanan gs = coef * bs * (t0 - ts); 501e77315bcSPatrick Sanan 502e77315bcSPatrick Sanan tn = x[k][j + 1][i]; 503e77315bcSPatrick Sanan an = 0.5 * (t0 + tn); 504e77315bcSPatrick Sanan bn = PetscPowScalar(an, bm1); 505e77315bcSPatrick Sanan /* dn = bn * an; */ 506e77315bcSPatrick Sanan dn = PetscPowScalar(an, beta); 507e77315bcSPatrick Sanan gn = coef * bn * (tn - t0); 508e77315bcSPatrick Sanan 509e77315bcSPatrick Sanan td = x[k - 1][j][i]; 510e77315bcSPatrick Sanan ad = 0.5 * (t0 + td); 511e77315bcSPatrick Sanan bd = PetscPowScalar(ad, bm1); 512e77315bcSPatrick Sanan /* dd = bd * ad; */ 513e77315bcSPatrick Sanan dd = PetscPowScalar(ad, beta); 514e77315bcSPatrick Sanan gd = coef * bd * (t0 - td); 515e77315bcSPatrick Sanan 516e77315bcSPatrick Sanan tu = x[k + 1][j][i]; 517e77315bcSPatrick Sanan au = 0.5 * (t0 + tu); 518e77315bcSPatrick Sanan bu = PetscPowScalar(au, bm1); 519e77315bcSPatrick Sanan /* du = bu * au; */ 520e77315bcSPatrick Sanan du = PetscPowScalar(au, beta); 521e77315bcSPatrick Sanan gu = coef * bu * (tu - t0); 522e77315bcSPatrick Sanan 5239371c9d4SSatish Balay c[0].k = k - 1; 5249371c9d4SSatish Balay c[0].j = j; 5259371c9d4SSatish Balay c[0].i = i; 5269371c9d4SSatish Balay v[0] = -hxhydhz * (dd - gd); 5279371c9d4SSatish Balay c[1].k = k; 5289371c9d4SSatish Balay c[1].j = j - 1; 5299371c9d4SSatish Balay c[1].i = i; 530e77315bcSPatrick Sanan v[1] = -hzhxdhy * (ds - gs); 5319371c9d4SSatish Balay c[2].k = k; 5329371c9d4SSatish Balay c[2].j = j; 5339371c9d4SSatish Balay c[2].i = i - 1; 534e77315bcSPatrick Sanan v[2] = -hyhzdhx * (dw - gw); 5359371c9d4SSatish Balay c[3].k = k; 5369371c9d4SSatish Balay c[3].j = j; 5379371c9d4SSatish Balay c[3].i = i; 538e77315bcSPatrick Sanan v[3] = hzhxdhy * (ds + dn + gs - gn) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (dd + du + gd - gu); 5399371c9d4SSatish Balay c[4].k = k; 5409371c9d4SSatish Balay c[4].j = j; 5419371c9d4SSatish Balay c[4].i = i + 1; 542e77315bcSPatrick Sanan v[4] = -hyhzdhx * (de + ge); 5439371c9d4SSatish Balay c[5].k = k; 5449371c9d4SSatish Balay c[5].j = j + 1; 5459371c9d4SSatish Balay c[5].i = i; 546e77315bcSPatrick Sanan v[5] = -hzhxdhy * (dn + gn); 5479371c9d4SSatish Balay c[6].k = k + 1; 5489371c9d4SSatish Balay c[6].j = j; 5499371c9d4SSatish Balay c[6].i = i; 550e77315bcSPatrick Sanan v[6] = -hxhydhz * (du + gu); 5519566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac, 1, &row, 7, c, v, INSERT_VALUES)); 552e77315bcSPatrick Sanan 553e77315bcSPatrick Sanan } else if (i == 0) { 554e77315bcSPatrick Sanan /* left-hand plane boundary */ 555e77315bcSPatrick Sanan tw = tleft; 556e77315bcSPatrick Sanan aw = 0.5 * (t0 + tw); 557e77315bcSPatrick Sanan bw = PetscPowScalar(aw, bm1); 558e77315bcSPatrick Sanan /* dw = bw * aw */ 559e77315bcSPatrick Sanan dw = PetscPowScalar(aw, beta); 560e77315bcSPatrick Sanan gw = coef * bw * (t0 - tw); 561e77315bcSPatrick Sanan 562e77315bcSPatrick Sanan te = x[k][j][i + 1]; 563e77315bcSPatrick Sanan ae = 0.5 * (t0 + te); 564e77315bcSPatrick Sanan be = PetscPowScalar(ae, bm1); 565e77315bcSPatrick Sanan /* de = be * ae; */ 566e77315bcSPatrick Sanan de = PetscPowScalar(ae, beta); 567e77315bcSPatrick Sanan ge = coef * be * (te - t0); 568e77315bcSPatrick Sanan 569e77315bcSPatrick Sanan /* left-hand bottom edge */ 570e77315bcSPatrick Sanan if (j == 0) { 571e77315bcSPatrick Sanan tn = x[k][j + 1][i]; 572e77315bcSPatrick Sanan an = 0.5 * (t0 + tn); 573e77315bcSPatrick Sanan bn = PetscPowScalar(an, bm1); 574e77315bcSPatrick Sanan /* dn = bn * an; */ 575e77315bcSPatrick Sanan dn = PetscPowScalar(an, beta); 576e77315bcSPatrick Sanan gn = coef * bn * (tn - t0); 577e77315bcSPatrick Sanan 578e77315bcSPatrick Sanan /* left-hand bottom down corner */ 579e77315bcSPatrick Sanan if (k == 0) { 580e77315bcSPatrick Sanan tu = x[k + 1][j][i]; 581e77315bcSPatrick Sanan au = 0.5 * (t0 + tu); 582e77315bcSPatrick Sanan bu = PetscPowScalar(au, bm1); 583e77315bcSPatrick Sanan /* du = bu * au; */ 584e77315bcSPatrick Sanan du = PetscPowScalar(au, beta); 585e77315bcSPatrick Sanan gu = coef * bu * (tu - t0); 586e77315bcSPatrick Sanan 5879371c9d4SSatish Balay c[0].k = k; 5889371c9d4SSatish Balay c[0].j = j; 5899371c9d4SSatish Balay c[0].i = i; 590e77315bcSPatrick Sanan v[0] = hzhxdhy * (dn - gn) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (du - gu); 5919371c9d4SSatish Balay c[1].k = k; 5929371c9d4SSatish Balay c[1].j = j; 5939371c9d4SSatish Balay c[1].i = i + 1; 594e77315bcSPatrick Sanan v[1] = -hyhzdhx * (de + ge); 5959371c9d4SSatish Balay c[2].k = k; 5969371c9d4SSatish Balay c[2].j = j + 1; 5979371c9d4SSatish Balay c[2].i = i; 598e77315bcSPatrick Sanan v[2] = -hzhxdhy * (dn + gn); 5999371c9d4SSatish Balay c[3].k = k + 1; 6009371c9d4SSatish Balay c[3].j = j; 6019371c9d4SSatish Balay c[3].i = i; 602e77315bcSPatrick Sanan v[3] = -hxhydhz * (du + gu); 6039566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac, 1, &row, 4, c, v, INSERT_VALUES)); 604e77315bcSPatrick Sanan 605e77315bcSPatrick Sanan /* left-hand bottom interior edge */ 606e77315bcSPatrick Sanan } else if (k < mz - 1) { 607e77315bcSPatrick Sanan tu = x[k + 1][j][i]; 608e77315bcSPatrick Sanan au = 0.5 * (t0 + tu); 609e77315bcSPatrick Sanan bu = PetscPowScalar(au, bm1); 610e77315bcSPatrick Sanan /* du = bu * au; */ 611e77315bcSPatrick Sanan du = PetscPowScalar(au, beta); 612e77315bcSPatrick Sanan gu = coef * bu * (tu - t0); 613e77315bcSPatrick Sanan 614e77315bcSPatrick Sanan td = x[k - 1][j][i]; 615e77315bcSPatrick Sanan ad = 0.5 * (t0 + td); 616e77315bcSPatrick Sanan bd = PetscPowScalar(ad, bm1); 617e77315bcSPatrick Sanan /* dd = bd * ad; */ 618e77315bcSPatrick Sanan dd = PetscPowScalar(ad, beta); 619e77315bcSPatrick Sanan gd = coef * bd * (td - t0); 620e77315bcSPatrick Sanan 6219371c9d4SSatish Balay c[0].k = k - 1; 6229371c9d4SSatish Balay c[0].j = j; 6239371c9d4SSatish Balay c[0].i = i; 624e77315bcSPatrick Sanan v[0] = -hxhydhz * (dd - gd); 6259371c9d4SSatish Balay c[1].k = k; 6269371c9d4SSatish Balay c[1].j = j; 6279371c9d4SSatish Balay c[1].i = i; 628e77315bcSPatrick Sanan v[1] = hzhxdhy * (dn - gn) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (dd + du + gd - gu); 6299371c9d4SSatish Balay c[2].k = k; 6309371c9d4SSatish Balay c[2].j = j; 6319371c9d4SSatish Balay c[2].i = i + 1; 632e77315bcSPatrick Sanan v[2] = -hyhzdhx * (de + ge); 6339371c9d4SSatish Balay c[3].k = k; 6349371c9d4SSatish Balay c[3].j = j + 1; 6359371c9d4SSatish Balay c[3].i = i; 636e77315bcSPatrick Sanan v[3] = -hzhxdhy * (dn + gn); 6379371c9d4SSatish Balay c[4].k = k + 1; 6389371c9d4SSatish Balay c[4].j = j; 6399371c9d4SSatish Balay c[4].i = i; 640e77315bcSPatrick Sanan v[4] = -hxhydhz * (du + gu); 6419566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac, 1, &row, 5, c, v, INSERT_VALUES)); 642e77315bcSPatrick Sanan 643e77315bcSPatrick Sanan /* left-hand bottom up corner */ 644e77315bcSPatrick Sanan } else { 645e77315bcSPatrick Sanan td = x[k - 1][j][i]; 646e77315bcSPatrick Sanan ad = 0.5 * (t0 + td); 647e77315bcSPatrick Sanan bd = PetscPowScalar(ad, bm1); 648e77315bcSPatrick Sanan /* dd = bd * ad; */ 649e77315bcSPatrick Sanan dd = PetscPowScalar(ad, beta); 650e77315bcSPatrick Sanan gd = coef * bd * (td - t0); 651e77315bcSPatrick Sanan 6529371c9d4SSatish Balay c[0].k = k - 1; 6539371c9d4SSatish Balay c[0].j = j; 6549371c9d4SSatish Balay c[0].i = i; 655e77315bcSPatrick Sanan v[0] = -hxhydhz * (dd - gd); 6569371c9d4SSatish Balay c[1].k = k; 6579371c9d4SSatish Balay c[1].j = j; 6589371c9d4SSatish Balay c[1].i = i; 659e77315bcSPatrick Sanan v[1] = hzhxdhy * (dn - gn) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (dd + gd); 6609371c9d4SSatish Balay c[2].k = k; 6619371c9d4SSatish Balay c[2].j = j; 6629371c9d4SSatish Balay c[2].i = i + 1; 663e77315bcSPatrick Sanan v[2] = -hyhzdhx * (de + ge); 6649371c9d4SSatish Balay c[3].k = k; 6659371c9d4SSatish Balay c[3].j = j + 1; 6669371c9d4SSatish Balay c[3].i = i; 667e77315bcSPatrick Sanan v[3] = -hzhxdhy * (dn + gn); 6689566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac, 1, &row, 4, c, v, INSERT_VALUES)); 669e77315bcSPatrick Sanan } 670e77315bcSPatrick Sanan 671e77315bcSPatrick Sanan /* left-hand top edge */ 672e77315bcSPatrick Sanan } else if (j == my - 1) { 673e77315bcSPatrick Sanan ts = x[k][j - 1][i]; 674e77315bcSPatrick Sanan as = 0.5 * (t0 + ts); 675e77315bcSPatrick Sanan bs = PetscPowScalar(as, bm1); 676e77315bcSPatrick Sanan /* ds = bs * as; */ 677e77315bcSPatrick Sanan ds = PetscPowScalar(as, beta); 678e77315bcSPatrick Sanan gs = coef * bs * (ts - t0); 679e77315bcSPatrick Sanan 680e77315bcSPatrick Sanan /* left-hand top down corner */ 681e77315bcSPatrick Sanan if (k == 0) { 682e77315bcSPatrick Sanan tu = x[k + 1][j][i]; 683e77315bcSPatrick Sanan au = 0.5 * (t0 + tu); 684e77315bcSPatrick Sanan bu = PetscPowScalar(au, bm1); 685e77315bcSPatrick Sanan /* du = bu * au; */ 686e77315bcSPatrick Sanan du = PetscPowScalar(au, beta); 687e77315bcSPatrick Sanan gu = coef * bu * (tu - t0); 688e77315bcSPatrick Sanan 6899371c9d4SSatish Balay c[0].k = k; 6909371c9d4SSatish Balay c[0].j = j - 1; 6919371c9d4SSatish Balay c[0].i = i; 692e77315bcSPatrick Sanan v[0] = -hzhxdhy * (ds - gs); 6939371c9d4SSatish Balay c[1].k = k; 6949371c9d4SSatish Balay c[1].j = j; 6959371c9d4SSatish Balay c[1].i = i; 696e77315bcSPatrick Sanan v[1] = hzhxdhy * (ds + gs) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (du - gu); 6979371c9d4SSatish Balay c[2].k = k; 6989371c9d4SSatish Balay c[2].j = j; 6999371c9d4SSatish Balay c[2].i = i + 1; 700e77315bcSPatrick Sanan v[2] = -hyhzdhx * (de + ge); 7019371c9d4SSatish Balay c[3].k = k + 1; 7029371c9d4SSatish Balay c[3].j = j; 7039371c9d4SSatish Balay c[3].i = i; 704e77315bcSPatrick Sanan v[3] = -hxhydhz * (du + gu); 7059566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac, 1, &row, 4, c, v, INSERT_VALUES)); 706e77315bcSPatrick Sanan 707e77315bcSPatrick Sanan /* left-hand top interior edge */ 708e77315bcSPatrick Sanan } else if (k < mz - 1) { 709e77315bcSPatrick Sanan tu = x[k + 1][j][i]; 710e77315bcSPatrick Sanan au = 0.5 * (t0 + tu); 711e77315bcSPatrick Sanan bu = PetscPowScalar(au, bm1); 712e77315bcSPatrick Sanan /* du = bu * au; */ 713e77315bcSPatrick Sanan du = PetscPowScalar(au, beta); 714e77315bcSPatrick Sanan gu = coef * bu * (tu - t0); 715e77315bcSPatrick Sanan 716e77315bcSPatrick Sanan td = x[k - 1][j][i]; 717e77315bcSPatrick Sanan ad = 0.5 * (t0 + td); 718e77315bcSPatrick Sanan bd = PetscPowScalar(ad, bm1); 719e77315bcSPatrick Sanan /* dd = bd * ad; */ 720e77315bcSPatrick Sanan dd = PetscPowScalar(ad, beta); 721e77315bcSPatrick Sanan gd = coef * bd * (td - t0); 722e77315bcSPatrick Sanan 7239371c9d4SSatish Balay c[0].k = k - 1; 7249371c9d4SSatish Balay c[0].j = j; 7259371c9d4SSatish Balay c[0].i = i; 726e77315bcSPatrick Sanan v[0] = -hxhydhz * (dd - gd); 7279371c9d4SSatish Balay c[1].k = k; 7289371c9d4SSatish Balay c[1].j = j - 1; 7299371c9d4SSatish Balay c[1].i = i; 730e77315bcSPatrick Sanan v[1] = -hzhxdhy * (ds - gs); 7319371c9d4SSatish Balay c[2].k = k; 7329371c9d4SSatish Balay c[2].j = j; 7339371c9d4SSatish Balay c[2].i = i; 734e77315bcSPatrick Sanan v[2] = hzhxdhy * (ds + gs) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (dd + du + gd - gu); 7359371c9d4SSatish Balay c[3].k = k; 7369371c9d4SSatish Balay c[3].j = j; 7379371c9d4SSatish Balay c[3].i = i + 1; 738e77315bcSPatrick Sanan v[3] = -hyhzdhx * (de + ge); 7399371c9d4SSatish Balay c[4].k = k + 1; 7409371c9d4SSatish Balay c[4].j = j; 7419371c9d4SSatish Balay c[4].i = i; 742e77315bcSPatrick Sanan v[4] = -hxhydhz * (du + gu); 7439566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac, 1, &row, 5, c, v, INSERT_VALUES)); 744e77315bcSPatrick Sanan 745e77315bcSPatrick Sanan /* left-hand top up corner */ 746e77315bcSPatrick Sanan } else { 747e77315bcSPatrick Sanan td = x[k - 1][j][i]; 748e77315bcSPatrick Sanan ad = 0.5 * (t0 + td); 749e77315bcSPatrick Sanan bd = PetscPowScalar(ad, bm1); 750e77315bcSPatrick Sanan /* dd = bd * ad; */ 751e77315bcSPatrick Sanan dd = PetscPowScalar(ad, beta); 752e77315bcSPatrick Sanan gd = coef * bd * (td - t0); 753e77315bcSPatrick Sanan 7549371c9d4SSatish Balay c[0].k = k - 1; 7559371c9d4SSatish Balay c[0].j = j; 7569371c9d4SSatish Balay c[0].i = i; 757e77315bcSPatrick Sanan v[0] = -hxhydhz * (dd - gd); 7589371c9d4SSatish Balay c[1].k = k; 7599371c9d4SSatish Balay c[1].j = j - 1; 7609371c9d4SSatish Balay c[1].i = i; 761e77315bcSPatrick Sanan v[1] = -hzhxdhy * (ds - gs); 7629371c9d4SSatish Balay c[2].k = k; 7639371c9d4SSatish Balay c[2].j = j; 7649371c9d4SSatish Balay c[2].i = i; 765e77315bcSPatrick Sanan v[2] = hzhxdhy * (ds + gs) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (dd + gd); 7669371c9d4SSatish Balay c[3].k = k; 7679371c9d4SSatish Balay c[3].j = j; 7689371c9d4SSatish Balay c[3].i = i + 1; 769e77315bcSPatrick Sanan v[3] = -hyhzdhx * (de + ge); 7709566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac, 1, &row, 4, c, v, INSERT_VALUES)); 771e77315bcSPatrick Sanan } 772e77315bcSPatrick Sanan 773e77315bcSPatrick Sanan } else { 774e77315bcSPatrick Sanan ts = x[k][j - 1][i]; 775e77315bcSPatrick Sanan as = 0.5 * (t0 + ts); 776e77315bcSPatrick Sanan bs = PetscPowScalar(as, bm1); 777e77315bcSPatrick Sanan /* ds = bs * as; */ 778e77315bcSPatrick Sanan ds = PetscPowScalar(as, beta); 779e77315bcSPatrick Sanan gs = coef * bs * (t0 - ts); 780e77315bcSPatrick Sanan 781e77315bcSPatrick Sanan tn = x[k][j + 1][i]; 782e77315bcSPatrick Sanan an = 0.5 * (t0 + tn); 783e77315bcSPatrick Sanan bn = PetscPowScalar(an, bm1); 784e77315bcSPatrick Sanan /* dn = bn * an; */ 785e77315bcSPatrick Sanan dn = PetscPowScalar(an, beta); 786e77315bcSPatrick Sanan gn = coef * bn * (tn - t0); 787e77315bcSPatrick Sanan 788e77315bcSPatrick Sanan /* left-hand down interior edge */ 789e77315bcSPatrick Sanan if (k == 0) { 790e77315bcSPatrick Sanan tu = x[k + 1][j][i]; 791e77315bcSPatrick Sanan au = 0.5 * (t0 + tu); 792e77315bcSPatrick Sanan bu = PetscPowScalar(au, bm1); 793e77315bcSPatrick Sanan /* du = bu * au; */ 794e77315bcSPatrick Sanan du = PetscPowScalar(au, beta); 795e77315bcSPatrick Sanan gu = coef * bu * (tu - t0); 796e77315bcSPatrick Sanan 7979371c9d4SSatish Balay c[0].k = k; 7989371c9d4SSatish Balay c[0].j = j - 1; 7999371c9d4SSatish Balay c[0].i = i; 800e77315bcSPatrick Sanan v[0] = -hzhxdhy * (ds - gs); 8019371c9d4SSatish Balay c[1].k = k; 8029371c9d4SSatish Balay c[1].j = j; 8039371c9d4SSatish Balay c[1].i = i; 804e77315bcSPatrick Sanan v[1] = hzhxdhy * (ds + dn + gs - gn) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (du - gu); 8059371c9d4SSatish Balay c[2].k = k; 8069371c9d4SSatish Balay c[2].j = j; 8079371c9d4SSatish Balay c[2].i = i + 1; 808e77315bcSPatrick Sanan v[2] = -hyhzdhx * (de + ge); 8099371c9d4SSatish Balay c[3].k = k; 8109371c9d4SSatish Balay c[3].j = j + 1; 8119371c9d4SSatish Balay c[3].i = i; 812e77315bcSPatrick Sanan v[3] = -hzhxdhy * (dn + gn); 8139371c9d4SSatish Balay c[4].k = k + 1; 8149371c9d4SSatish Balay c[4].j = j; 8159371c9d4SSatish Balay c[4].i = i; 816e77315bcSPatrick Sanan v[4] = -hxhydhz * (du + gu); 8179566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac, 1, &row, 5, c, v, INSERT_VALUES)); 818e77315bcSPatrick Sanan 819e77315bcSPatrick Sanan } else if (k == mz - 1) { /* left-hand up interior edge */ 820e77315bcSPatrick Sanan 821e77315bcSPatrick Sanan td = x[k - 1][j][i]; 822e77315bcSPatrick Sanan ad = 0.5 * (t0 + td); 823e77315bcSPatrick Sanan bd = PetscPowScalar(ad, bm1); 824e77315bcSPatrick Sanan /* dd = bd * ad; */ 825e77315bcSPatrick Sanan dd = PetscPowScalar(ad, beta); 826e77315bcSPatrick Sanan gd = coef * bd * (t0 - td); 827e77315bcSPatrick Sanan 8289371c9d4SSatish Balay c[0].k = k - 1; 8299371c9d4SSatish Balay c[0].j = j; 8309371c9d4SSatish Balay c[0].i = i; 831e77315bcSPatrick Sanan v[0] = -hxhydhz * (dd - gd); 8329371c9d4SSatish Balay c[1].k = k; 8339371c9d4SSatish Balay c[1].j = j - 1; 8349371c9d4SSatish Balay c[1].i = i; 835e77315bcSPatrick Sanan v[1] = -hzhxdhy * (ds - gs); 8369371c9d4SSatish Balay c[2].k = k; 8379371c9d4SSatish Balay c[2].j = j; 8389371c9d4SSatish Balay c[2].i = i; 839e77315bcSPatrick Sanan v[2] = hzhxdhy * (ds + dn + gs - gn) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (dd + gd); 8409371c9d4SSatish Balay c[3].k = k; 8419371c9d4SSatish Balay c[3].j = j; 8429371c9d4SSatish Balay c[3].i = i + 1; 843e77315bcSPatrick Sanan v[3] = -hyhzdhx * (de + ge); 8449371c9d4SSatish Balay c[4].k = k; 8459371c9d4SSatish Balay c[4].j = j + 1; 8469371c9d4SSatish Balay c[4].i = i; 847e77315bcSPatrick Sanan v[4] = -hzhxdhy * (dn + gn); 8489566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac, 1, &row, 5, c, v, INSERT_VALUES)); 849e77315bcSPatrick Sanan } else { /* left-hand interior plane */ 850e77315bcSPatrick Sanan 851e77315bcSPatrick Sanan td = x[k - 1][j][i]; 852e77315bcSPatrick Sanan ad = 0.5 * (t0 + td); 853e77315bcSPatrick Sanan bd = PetscPowScalar(ad, bm1); 854e77315bcSPatrick Sanan /* dd = bd * ad; */ 855e77315bcSPatrick Sanan dd = PetscPowScalar(ad, beta); 856e77315bcSPatrick Sanan gd = coef * bd * (t0 - td); 857e77315bcSPatrick Sanan 858e77315bcSPatrick Sanan tu = x[k + 1][j][i]; 859e77315bcSPatrick Sanan au = 0.5 * (t0 + tu); 860e77315bcSPatrick Sanan bu = PetscPowScalar(au, bm1); 861e77315bcSPatrick Sanan /* du = bu * au; */ 862e77315bcSPatrick Sanan du = PetscPowScalar(au, beta); 863e77315bcSPatrick Sanan gu = coef * bu * (tu - t0); 864e77315bcSPatrick Sanan 8659371c9d4SSatish Balay c[0].k = k - 1; 8669371c9d4SSatish Balay c[0].j = j; 8679371c9d4SSatish Balay c[0].i = i; 868e77315bcSPatrick Sanan v[0] = -hxhydhz * (dd - gd); 8699371c9d4SSatish Balay c[1].k = k; 8709371c9d4SSatish Balay c[1].j = j - 1; 8719371c9d4SSatish Balay c[1].i = i; 872e77315bcSPatrick Sanan v[1] = -hzhxdhy * (ds - gs); 8739371c9d4SSatish Balay c[2].k = k; 8749371c9d4SSatish Balay c[2].j = j; 8759371c9d4SSatish Balay c[2].i = i; 876e77315bcSPatrick Sanan v[2] = hzhxdhy * (ds + dn + gs - gn) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (dd + du + gd - gu); 8779371c9d4SSatish Balay c[3].k = k; 8789371c9d4SSatish Balay c[3].j = j; 8799371c9d4SSatish Balay c[3].i = i + 1; 880e77315bcSPatrick Sanan v[3] = -hyhzdhx * (de + ge); 8819371c9d4SSatish Balay c[4].k = k; 8829371c9d4SSatish Balay c[4].j = j + 1; 8839371c9d4SSatish Balay c[4].i = i; 884e77315bcSPatrick Sanan v[4] = -hzhxdhy * (dn + gn); 8859371c9d4SSatish Balay c[5].k = k + 1; 8869371c9d4SSatish Balay c[5].j = j; 8879371c9d4SSatish Balay c[5].i = i; 888e77315bcSPatrick Sanan v[5] = -hxhydhz * (du + gu); 8899566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac, 1, &row, 6, c, v, INSERT_VALUES)); 890e77315bcSPatrick Sanan } 891e77315bcSPatrick Sanan } 892e77315bcSPatrick Sanan 893e77315bcSPatrick Sanan } else if (i == mx - 1) { 894e77315bcSPatrick Sanan /* right-hand plane boundary */ 895e77315bcSPatrick Sanan tw = x[k][j][i - 1]; 896e77315bcSPatrick Sanan aw = 0.5 * (t0 + tw); 897e77315bcSPatrick Sanan bw = PetscPowScalar(aw, bm1); 898e77315bcSPatrick Sanan /* dw = bw * aw */ 899e77315bcSPatrick Sanan dw = PetscPowScalar(aw, beta); 900e77315bcSPatrick Sanan gw = coef * bw * (t0 - tw); 901e77315bcSPatrick Sanan 902e77315bcSPatrick Sanan te = tright; 903e77315bcSPatrick Sanan ae = 0.5 * (t0 + te); 904e77315bcSPatrick Sanan be = PetscPowScalar(ae, bm1); 905e77315bcSPatrick Sanan /* de = be * ae; */ 906e77315bcSPatrick Sanan de = PetscPowScalar(ae, beta); 907e77315bcSPatrick Sanan ge = coef * be * (te - t0); 908e77315bcSPatrick Sanan 909e77315bcSPatrick Sanan /* right-hand bottom edge */ 910e77315bcSPatrick Sanan if (j == 0) { 911e77315bcSPatrick Sanan tn = x[k][j + 1][i]; 912e77315bcSPatrick Sanan an = 0.5 * (t0 + tn); 913e77315bcSPatrick Sanan bn = PetscPowScalar(an, bm1); 914e77315bcSPatrick Sanan /* dn = bn * an; */ 915e77315bcSPatrick Sanan dn = PetscPowScalar(an, beta); 916e77315bcSPatrick Sanan gn = coef * bn * (tn - t0); 917e77315bcSPatrick Sanan 918e77315bcSPatrick Sanan /* right-hand bottom down corner */ 919e77315bcSPatrick Sanan if (k == 0) { 920e77315bcSPatrick Sanan tu = x[k + 1][j][i]; 921e77315bcSPatrick Sanan au = 0.5 * (t0 + tu); 922e77315bcSPatrick Sanan bu = PetscPowScalar(au, bm1); 923e77315bcSPatrick Sanan /* du = bu * au; */ 924e77315bcSPatrick Sanan du = PetscPowScalar(au, beta); 925e77315bcSPatrick Sanan gu = coef * bu * (tu - t0); 926e77315bcSPatrick Sanan 9279371c9d4SSatish Balay c[0].k = k; 9289371c9d4SSatish Balay c[0].j = j; 9299371c9d4SSatish Balay c[0].i = i - 1; 930e77315bcSPatrick Sanan v[0] = -hyhzdhx * (dw - gw); 9319371c9d4SSatish Balay c[1].k = k; 9329371c9d4SSatish Balay c[1].j = j; 9339371c9d4SSatish Balay c[1].i = i; 934e77315bcSPatrick Sanan v[1] = hzhxdhy * (dn - gn) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (du - gu); 9359371c9d4SSatish Balay c[2].k = k; 9369371c9d4SSatish Balay c[2].j = j + 1; 9379371c9d4SSatish Balay c[2].i = i; 938e77315bcSPatrick Sanan v[2] = -hzhxdhy * (dn + gn); 9399371c9d4SSatish Balay c[3].k = k + 1; 9409371c9d4SSatish Balay c[3].j = j; 9419371c9d4SSatish Balay c[3].i = i; 942e77315bcSPatrick Sanan v[3] = -hxhydhz * (du + gu); 9439566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac, 1, &row, 4, c, v, INSERT_VALUES)); 944e77315bcSPatrick Sanan 945e77315bcSPatrick Sanan /* right-hand bottom interior edge */ 946e77315bcSPatrick Sanan } else if (k < mz - 1) { 947e77315bcSPatrick Sanan tu = x[k + 1][j][i]; 948e77315bcSPatrick Sanan au = 0.5 * (t0 + tu); 949e77315bcSPatrick Sanan bu = PetscPowScalar(au, bm1); 950e77315bcSPatrick Sanan /* du = bu * au; */ 951e77315bcSPatrick Sanan du = PetscPowScalar(au, beta); 952e77315bcSPatrick Sanan gu = coef * bu * (tu - t0); 953e77315bcSPatrick Sanan 954e77315bcSPatrick Sanan td = x[k - 1][j][i]; 955e77315bcSPatrick Sanan ad = 0.5 * (t0 + td); 956e77315bcSPatrick Sanan bd = PetscPowScalar(ad, bm1); 957e77315bcSPatrick Sanan /* dd = bd * ad; */ 958e77315bcSPatrick Sanan dd = PetscPowScalar(ad, beta); 959e77315bcSPatrick Sanan gd = coef * bd * (td - t0); 960e77315bcSPatrick Sanan 9619371c9d4SSatish Balay c[0].k = k - 1; 9629371c9d4SSatish Balay c[0].j = j; 9639371c9d4SSatish Balay c[0].i = i; 964e77315bcSPatrick Sanan v[0] = -hxhydhz * (dd - gd); 9659371c9d4SSatish Balay c[1].k = k; 9669371c9d4SSatish Balay c[1].j = j; 9679371c9d4SSatish Balay c[1].i = i - 1; 968e77315bcSPatrick Sanan v[1] = -hyhzdhx * (dw - gw); 9699371c9d4SSatish Balay c[2].k = k; 9709371c9d4SSatish Balay c[2].j = j; 9719371c9d4SSatish Balay c[2].i = i; 972e77315bcSPatrick Sanan v[2] = hzhxdhy * (dn - gn) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (dd + du + gd - gu); 9739371c9d4SSatish Balay c[3].k = k; 9749371c9d4SSatish Balay c[3].j = j + 1; 9759371c9d4SSatish Balay c[3].i = i; 976e77315bcSPatrick Sanan v[3] = -hzhxdhy * (dn + gn); 9779371c9d4SSatish Balay c[4].k = k + 1; 9789371c9d4SSatish Balay c[4].j = j; 9799371c9d4SSatish Balay c[4].i = i; 980e77315bcSPatrick Sanan v[4] = -hxhydhz * (du + gu); 9819566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac, 1, &row, 5, c, v, INSERT_VALUES)); 982e77315bcSPatrick Sanan 983e77315bcSPatrick Sanan /* right-hand bottom up corner */ 984e77315bcSPatrick Sanan } else { 985e77315bcSPatrick Sanan td = x[k - 1][j][i]; 986e77315bcSPatrick Sanan ad = 0.5 * (t0 + td); 987e77315bcSPatrick Sanan bd = PetscPowScalar(ad, bm1); 988e77315bcSPatrick Sanan /* dd = bd * ad; */ 989e77315bcSPatrick Sanan dd = PetscPowScalar(ad, beta); 990e77315bcSPatrick Sanan gd = coef * bd * (td - t0); 991e77315bcSPatrick Sanan 9929371c9d4SSatish Balay c[0].k = k - 1; 9939371c9d4SSatish Balay c[0].j = j; 9949371c9d4SSatish Balay c[0].i = i; 995e77315bcSPatrick Sanan v[0] = -hxhydhz * (dd - gd); 9969371c9d4SSatish Balay c[1].k = k; 9979371c9d4SSatish Balay c[1].j = j; 9989371c9d4SSatish Balay c[1].i = i - 1; 999e77315bcSPatrick Sanan v[1] = -hyhzdhx * (dw - gw); 10009371c9d4SSatish Balay c[2].k = k; 10019371c9d4SSatish Balay c[2].j = j; 10029371c9d4SSatish Balay c[2].i = i; 1003e77315bcSPatrick Sanan v[2] = hzhxdhy * (dn - gn) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (dd + gd); 10049371c9d4SSatish Balay c[3].k = k; 10059371c9d4SSatish Balay c[3].j = j + 1; 10069371c9d4SSatish Balay c[3].i = i; 1007e77315bcSPatrick Sanan v[3] = -hzhxdhy * (dn + gn); 10089566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac, 1, &row, 4, c, v, INSERT_VALUES)); 1009e77315bcSPatrick Sanan } 1010e77315bcSPatrick Sanan 1011e77315bcSPatrick Sanan /* right-hand top edge */ 1012e77315bcSPatrick Sanan } else if (j == my - 1) { 1013e77315bcSPatrick Sanan ts = x[k][j - 1][i]; 1014e77315bcSPatrick Sanan as = 0.5 * (t0 + ts); 1015e77315bcSPatrick Sanan bs = PetscPowScalar(as, bm1); 1016e77315bcSPatrick Sanan /* ds = bs * as; */ 1017e77315bcSPatrick Sanan ds = PetscPowScalar(as, beta); 1018e77315bcSPatrick Sanan gs = coef * bs * (ts - t0); 1019e77315bcSPatrick Sanan 1020e77315bcSPatrick Sanan /* right-hand top down corner */ 1021e77315bcSPatrick Sanan if (k == 0) { 1022e77315bcSPatrick Sanan tu = x[k + 1][j][i]; 1023e77315bcSPatrick Sanan au = 0.5 * (t0 + tu); 1024e77315bcSPatrick Sanan bu = PetscPowScalar(au, bm1); 1025e77315bcSPatrick Sanan /* du = bu * au; */ 1026e77315bcSPatrick Sanan du = PetscPowScalar(au, beta); 1027e77315bcSPatrick Sanan gu = coef * bu * (tu - t0); 1028e77315bcSPatrick Sanan 10299371c9d4SSatish Balay c[0].k = k; 10309371c9d4SSatish Balay c[0].j = j - 1; 10319371c9d4SSatish Balay c[0].i = i; 1032e77315bcSPatrick Sanan v[0] = -hzhxdhy * (ds - gs); 10339371c9d4SSatish Balay c[1].k = k; 10349371c9d4SSatish Balay c[1].j = j; 10359371c9d4SSatish Balay c[1].i = i - 1; 1036e77315bcSPatrick Sanan v[1] = -hyhzdhx * (dw - gw); 10379371c9d4SSatish Balay c[2].k = k; 10389371c9d4SSatish Balay c[2].j = j; 10399371c9d4SSatish Balay c[2].i = i; 1040e77315bcSPatrick Sanan v[2] = hzhxdhy * (ds + gs) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (du - gu); 10419371c9d4SSatish Balay c[3].k = k + 1; 10429371c9d4SSatish Balay c[3].j = j; 10439371c9d4SSatish Balay c[3].i = i; 1044e77315bcSPatrick Sanan v[3] = -hxhydhz * (du + gu); 10459566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac, 1, &row, 4, c, v, INSERT_VALUES)); 1046e77315bcSPatrick Sanan 1047e77315bcSPatrick Sanan /* right-hand top interior edge */ 1048e77315bcSPatrick Sanan } else if (k < mz - 1) { 1049e77315bcSPatrick Sanan tu = x[k + 1][j][i]; 1050e77315bcSPatrick Sanan au = 0.5 * (t0 + tu); 1051e77315bcSPatrick Sanan bu = PetscPowScalar(au, bm1); 1052e77315bcSPatrick Sanan /* du = bu * au; */ 1053e77315bcSPatrick Sanan du = PetscPowScalar(au, beta); 1054e77315bcSPatrick Sanan gu = coef * bu * (tu - t0); 1055e77315bcSPatrick Sanan 1056e77315bcSPatrick Sanan td = x[k - 1][j][i]; 1057e77315bcSPatrick Sanan ad = 0.5 * (t0 + td); 1058e77315bcSPatrick Sanan bd = PetscPowScalar(ad, bm1); 1059e77315bcSPatrick Sanan /* dd = bd * ad; */ 1060e77315bcSPatrick Sanan dd = PetscPowScalar(ad, beta); 1061e77315bcSPatrick Sanan gd = coef * bd * (td - t0); 1062e77315bcSPatrick Sanan 10639371c9d4SSatish Balay c[0].k = k - 1; 10649371c9d4SSatish Balay c[0].j = j; 10659371c9d4SSatish Balay c[0].i = i; 1066e77315bcSPatrick Sanan v[0] = -hxhydhz * (dd - gd); 10679371c9d4SSatish Balay c[1].k = k; 10689371c9d4SSatish Balay c[1].j = j - 1; 10699371c9d4SSatish Balay c[1].i = i; 1070e77315bcSPatrick Sanan v[1] = -hzhxdhy * (ds - gs); 10719371c9d4SSatish Balay c[2].k = k; 10729371c9d4SSatish Balay c[2].j = j; 10739371c9d4SSatish Balay c[2].i = i - 1; 1074e77315bcSPatrick Sanan v[2] = -hyhzdhx * (dw - gw); 10759371c9d4SSatish Balay c[3].k = k; 10769371c9d4SSatish Balay c[3].j = j; 10779371c9d4SSatish Balay c[3].i = i; 1078e77315bcSPatrick Sanan v[3] = hzhxdhy * (ds + gs) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (dd + du + gd - gu); 10799371c9d4SSatish Balay c[4].k = k + 1; 10809371c9d4SSatish Balay c[4].j = j; 10819371c9d4SSatish Balay c[4].i = i; 1082e77315bcSPatrick Sanan v[4] = -hxhydhz * (du + gu); 10839566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac, 1, &row, 5, c, v, INSERT_VALUES)); 1084e77315bcSPatrick Sanan 1085e77315bcSPatrick Sanan /* right-hand top up corner */ 1086e77315bcSPatrick Sanan } else { 1087e77315bcSPatrick Sanan td = x[k - 1][j][i]; 1088e77315bcSPatrick Sanan ad = 0.5 * (t0 + td); 1089e77315bcSPatrick Sanan bd = PetscPowScalar(ad, bm1); 1090e77315bcSPatrick Sanan /* dd = bd * ad; */ 1091e77315bcSPatrick Sanan dd = PetscPowScalar(ad, beta); 1092e77315bcSPatrick Sanan gd = coef * bd * (td - t0); 1093e77315bcSPatrick Sanan 10949371c9d4SSatish Balay c[0].k = k - 1; 10959371c9d4SSatish Balay c[0].j = j; 10969371c9d4SSatish Balay c[0].i = i; 1097e77315bcSPatrick Sanan v[0] = -hxhydhz * (dd - gd); 10989371c9d4SSatish Balay c[1].k = k; 10999371c9d4SSatish Balay c[1].j = j - 1; 11009371c9d4SSatish Balay c[1].i = i; 1101e77315bcSPatrick Sanan v[1] = -hzhxdhy * (ds - gs); 11029371c9d4SSatish Balay c[2].k = k; 11039371c9d4SSatish Balay c[2].j = j; 11049371c9d4SSatish Balay c[2].i = i - 1; 1105e77315bcSPatrick Sanan v[2] = -hyhzdhx * (dw - gw); 11069371c9d4SSatish Balay c[3].k = k; 11079371c9d4SSatish Balay c[3].j = j; 11089371c9d4SSatish Balay c[3].i = i; 1109e77315bcSPatrick Sanan v[3] = hzhxdhy * (ds + gs) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (dd + gd); 11109566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac, 1, &row, 4, c, v, INSERT_VALUES)); 1111e77315bcSPatrick Sanan } 1112e77315bcSPatrick Sanan 1113e77315bcSPatrick Sanan } else { 1114e77315bcSPatrick Sanan ts = x[k][j - 1][i]; 1115e77315bcSPatrick Sanan as = 0.5 * (t0 + ts); 1116e77315bcSPatrick Sanan bs = PetscPowScalar(as, bm1); 1117e77315bcSPatrick Sanan /* ds = bs * as; */ 1118e77315bcSPatrick Sanan ds = PetscPowScalar(as, beta); 1119e77315bcSPatrick Sanan gs = coef * bs * (t0 - ts); 1120e77315bcSPatrick Sanan 1121e77315bcSPatrick Sanan tn = x[k][j + 1][i]; 1122e77315bcSPatrick Sanan an = 0.5 * (t0 + tn); 1123e77315bcSPatrick Sanan bn = PetscPowScalar(an, bm1); 1124e77315bcSPatrick Sanan /* dn = bn * an; */ 1125e77315bcSPatrick Sanan dn = PetscPowScalar(an, beta); 1126e77315bcSPatrick Sanan gn = coef * bn * (tn - t0); 1127e77315bcSPatrick Sanan 1128e77315bcSPatrick Sanan /* right-hand down interior edge */ 1129e77315bcSPatrick Sanan if (k == 0) { 1130e77315bcSPatrick Sanan tu = x[k + 1][j][i]; 1131e77315bcSPatrick Sanan au = 0.5 * (t0 + tu); 1132e77315bcSPatrick Sanan bu = PetscPowScalar(au, bm1); 1133e77315bcSPatrick Sanan /* du = bu * au; */ 1134e77315bcSPatrick Sanan du = PetscPowScalar(au, beta); 1135e77315bcSPatrick Sanan gu = coef * bu * (tu - t0); 1136e77315bcSPatrick Sanan 11379371c9d4SSatish Balay c[0].k = k; 11389371c9d4SSatish Balay c[0].j = j - 1; 11399371c9d4SSatish Balay c[0].i = i; 1140e77315bcSPatrick Sanan v[0] = -hzhxdhy * (ds - gs); 11419371c9d4SSatish Balay c[1].k = k; 11429371c9d4SSatish Balay c[1].j = j; 11439371c9d4SSatish Balay c[1].i = i - 1; 1144e77315bcSPatrick Sanan v[1] = -hyhzdhx * (dw - gw); 11459371c9d4SSatish Balay c[2].k = k; 11469371c9d4SSatish Balay c[2].j = j; 11479371c9d4SSatish Balay c[2].i = i; 1148e77315bcSPatrick Sanan v[2] = hzhxdhy * (ds + dn + gs - gn) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (du - gu); 11499371c9d4SSatish Balay c[3].k = k; 11509371c9d4SSatish Balay c[3].j = j + 1; 11519371c9d4SSatish Balay c[3].i = i; 1152e77315bcSPatrick Sanan v[3] = -hzhxdhy * (dn + gn); 11539371c9d4SSatish Balay c[4].k = k + 1; 11549371c9d4SSatish Balay c[4].j = j; 11559371c9d4SSatish Balay c[4].i = i; 1156e77315bcSPatrick Sanan v[4] = -hxhydhz * (du + gu); 11579566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac, 1, &row, 5, c, v, INSERT_VALUES)); 1158e77315bcSPatrick Sanan 1159e77315bcSPatrick Sanan } else if (k == mz - 1) { /* right-hand up interior edge */ 1160e77315bcSPatrick Sanan 1161e77315bcSPatrick Sanan td = x[k - 1][j][i]; 1162e77315bcSPatrick Sanan ad = 0.5 * (t0 + td); 1163e77315bcSPatrick Sanan bd = PetscPowScalar(ad, bm1); 1164e77315bcSPatrick Sanan /* dd = bd * ad; */ 1165e77315bcSPatrick Sanan dd = PetscPowScalar(ad, beta); 1166e77315bcSPatrick Sanan gd = coef * bd * (t0 - td); 1167e77315bcSPatrick Sanan 11689371c9d4SSatish Balay c[0].k = k - 1; 11699371c9d4SSatish Balay c[0].j = j; 11709371c9d4SSatish Balay c[0].i = i; 1171e77315bcSPatrick Sanan v[0] = -hxhydhz * (dd - gd); 11729371c9d4SSatish Balay c[1].k = k; 11739371c9d4SSatish Balay c[1].j = j - 1; 11749371c9d4SSatish Balay c[1].i = i; 1175e77315bcSPatrick Sanan v[1] = -hzhxdhy * (ds - gs); 11769371c9d4SSatish Balay c[2].k = k; 11779371c9d4SSatish Balay c[2].j = j; 11789371c9d4SSatish Balay c[2].i = i - 1; 1179e77315bcSPatrick Sanan v[2] = -hyhzdhx * (dw - gw); 11809371c9d4SSatish Balay c[3].k = k; 11819371c9d4SSatish Balay c[3].j = j; 11829371c9d4SSatish Balay c[3].i = i; 1183e77315bcSPatrick Sanan v[3] = hzhxdhy * (ds + dn + gs - gn) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (dd + gd); 11849371c9d4SSatish Balay c[4].k = k; 11859371c9d4SSatish Balay c[4].j = j + 1; 11869371c9d4SSatish Balay c[4].i = i; 1187e77315bcSPatrick Sanan v[4] = -hzhxdhy * (dn + gn); 11889566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac, 1, &row, 5, c, v, INSERT_VALUES)); 1189e77315bcSPatrick Sanan 1190e77315bcSPatrick Sanan } else { /* right-hand interior plane */ 1191e77315bcSPatrick Sanan 1192e77315bcSPatrick Sanan td = x[k - 1][j][i]; 1193e77315bcSPatrick Sanan ad = 0.5 * (t0 + td); 1194e77315bcSPatrick Sanan bd = PetscPowScalar(ad, bm1); 1195e77315bcSPatrick Sanan /* dd = bd * ad; */ 1196e77315bcSPatrick Sanan dd = PetscPowScalar(ad, beta); 1197e77315bcSPatrick Sanan gd = coef * bd * (t0 - td); 1198e77315bcSPatrick Sanan 1199e77315bcSPatrick Sanan tu = x[k + 1][j][i]; 1200e77315bcSPatrick Sanan au = 0.5 * (t0 + tu); 1201e77315bcSPatrick Sanan bu = PetscPowScalar(au, bm1); 1202e77315bcSPatrick Sanan /* du = bu * au; */ 1203e77315bcSPatrick Sanan du = PetscPowScalar(au, beta); 1204e77315bcSPatrick Sanan gu = coef * bu * (tu - t0); 1205e77315bcSPatrick Sanan 12069371c9d4SSatish Balay c[0].k = k - 1; 12079371c9d4SSatish Balay c[0].j = j; 12089371c9d4SSatish Balay c[0].i = i; 1209e77315bcSPatrick Sanan v[0] = -hxhydhz * (dd - gd); 12109371c9d4SSatish Balay c[1].k = k; 12119371c9d4SSatish Balay c[1].j = j - 1; 12129371c9d4SSatish Balay c[1].i = i; 1213e77315bcSPatrick Sanan v[1] = -hzhxdhy * (ds - gs); 12149371c9d4SSatish Balay c[2].k = k; 12159371c9d4SSatish Balay c[2].j = j; 12169371c9d4SSatish Balay c[2].i = i - 1; 1217e77315bcSPatrick Sanan v[2] = -hyhzdhx * (dw - gw); 12189371c9d4SSatish Balay c[3].k = k; 12199371c9d4SSatish Balay c[3].j = j; 12209371c9d4SSatish Balay c[3].i = i; 1221e77315bcSPatrick Sanan v[3] = hzhxdhy * (ds + dn + gs - gn) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (dd + du + gd - gu); 12229371c9d4SSatish Balay c[4].k = k; 12239371c9d4SSatish Balay c[4].j = j + 1; 12249371c9d4SSatish Balay c[4].i = i; 1225e77315bcSPatrick Sanan v[4] = -hzhxdhy * (dn + gn); 12269371c9d4SSatish Balay c[5].k = k + 1; 12279371c9d4SSatish Balay c[5].j = j; 12289371c9d4SSatish Balay c[5].i = i; 1229e77315bcSPatrick Sanan v[5] = -hxhydhz * (du + gu); 12309566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac, 1, &row, 6, c, v, INSERT_VALUES)); 1231e77315bcSPatrick Sanan } 1232e77315bcSPatrick Sanan } 1233e77315bcSPatrick Sanan 1234e77315bcSPatrick Sanan } else if (j == 0) { 1235e77315bcSPatrick Sanan tw = x[k][j][i - 1]; 1236e77315bcSPatrick Sanan aw = 0.5 * (t0 + tw); 1237e77315bcSPatrick Sanan bw = PetscPowScalar(aw, bm1); 1238e77315bcSPatrick Sanan /* dw = bw * aw */ 1239e77315bcSPatrick Sanan dw = PetscPowScalar(aw, beta); 1240e77315bcSPatrick Sanan gw = coef * bw * (t0 - tw); 1241e77315bcSPatrick Sanan 1242e77315bcSPatrick Sanan te = x[k][j][i + 1]; 1243e77315bcSPatrick Sanan ae = 0.5 * (t0 + te); 1244e77315bcSPatrick Sanan be = PetscPowScalar(ae, bm1); 1245e77315bcSPatrick Sanan /* de = be * ae; */ 1246e77315bcSPatrick Sanan de = PetscPowScalar(ae, beta); 1247e77315bcSPatrick Sanan ge = coef * be * (te - t0); 1248e77315bcSPatrick Sanan 1249e77315bcSPatrick Sanan tn = x[k][j + 1][i]; 1250e77315bcSPatrick Sanan an = 0.5 * (t0 + tn); 1251e77315bcSPatrick Sanan bn = PetscPowScalar(an, bm1); 1252e77315bcSPatrick Sanan /* dn = bn * an; */ 1253e77315bcSPatrick Sanan dn = PetscPowScalar(an, beta); 1254e77315bcSPatrick Sanan gn = coef * bn * (tn - t0); 1255e77315bcSPatrick Sanan 1256e77315bcSPatrick Sanan /* bottom down interior edge */ 1257e77315bcSPatrick Sanan if (k == 0) { 1258e77315bcSPatrick Sanan tu = x[k + 1][j][i]; 1259e77315bcSPatrick Sanan au = 0.5 * (t0 + tu); 1260e77315bcSPatrick Sanan bu = PetscPowScalar(au, bm1); 1261e77315bcSPatrick Sanan /* du = bu * au; */ 1262e77315bcSPatrick Sanan du = PetscPowScalar(au, beta); 1263e77315bcSPatrick Sanan gu = coef * bu * (tu - t0); 1264e77315bcSPatrick Sanan 12659371c9d4SSatish Balay c[0].k = k; 12669371c9d4SSatish Balay c[0].j = j; 12679371c9d4SSatish Balay c[0].i = i - 1; 1268e77315bcSPatrick Sanan v[0] = -hyhzdhx * (dw - gw); 12699371c9d4SSatish Balay c[1].k = k; 12709371c9d4SSatish Balay c[1].j = j; 12719371c9d4SSatish Balay c[1].i = i; 1272e77315bcSPatrick Sanan v[1] = hzhxdhy * (dn - gn) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (du - gu); 12739371c9d4SSatish Balay c[2].k = k; 12749371c9d4SSatish Balay c[2].j = j; 12759371c9d4SSatish Balay c[2].i = i + 1; 1276e77315bcSPatrick Sanan v[2] = -hyhzdhx * (de + ge); 12779371c9d4SSatish Balay c[3].k = k; 12789371c9d4SSatish Balay c[3].j = j + 1; 12799371c9d4SSatish Balay c[3].i = i; 1280e77315bcSPatrick Sanan v[3] = -hzhxdhy * (dn + gn); 12819371c9d4SSatish Balay c[4].k = k + 1; 12829371c9d4SSatish Balay c[4].j = j; 12839371c9d4SSatish Balay c[4].i = i; 1284e77315bcSPatrick Sanan v[4] = -hxhydhz * (du + gu); 12859566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac, 1, &row, 5, c, v, INSERT_VALUES)); 1286e77315bcSPatrick Sanan 1287e77315bcSPatrick Sanan } else if (k == mz - 1) { /* bottom up interior edge */ 1288e77315bcSPatrick Sanan 1289e77315bcSPatrick Sanan td = x[k - 1][j][i]; 1290e77315bcSPatrick Sanan ad = 0.5 * (t0 + td); 1291e77315bcSPatrick Sanan bd = PetscPowScalar(ad, bm1); 1292e77315bcSPatrick Sanan /* dd = bd * ad; */ 1293e77315bcSPatrick Sanan dd = PetscPowScalar(ad, beta); 1294e77315bcSPatrick Sanan gd = coef * bd * (td - t0); 1295e77315bcSPatrick Sanan 12969371c9d4SSatish Balay c[0].k = k - 1; 12979371c9d4SSatish Balay c[0].j = j; 12989371c9d4SSatish Balay c[0].i = i; 1299e77315bcSPatrick Sanan v[0] = -hxhydhz * (dd - gd); 13009371c9d4SSatish Balay c[1].k = k; 13019371c9d4SSatish Balay c[1].j = j; 13029371c9d4SSatish Balay c[1].i = i - 1; 1303e77315bcSPatrick Sanan v[1] = -hyhzdhx * (dw - gw); 13049371c9d4SSatish Balay c[2].k = k; 13059371c9d4SSatish Balay c[2].j = j; 13069371c9d4SSatish Balay c[2].i = i; 1307e77315bcSPatrick Sanan v[2] = hzhxdhy * (dn - gn) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (dd + gd); 13089371c9d4SSatish Balay c[3].k = k; 13099371c9d4SSatish Balay c[3].j = j; 13109371c9d4SSatish Balay c[3].i = i + 1; 1311e77315bcSPatrick Sanan v[3] = -hyhzdhx * (de + ge); 13129371c9d4SSatish Balay c[4].k = k; 13139371c9d4SSatish Balay c[4].j = j + 1; 13149371c9d4SSatish Balay c[4].i = i; 1315e77315bcSPatrick Sanan v[4] = -hzhxdhy * (dn + gn); 13169566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac, 1, &row, 5, c, v, INSERT_VALUES)); 1317e77315bcSPatrick Sanan 1318e77315bcSPatrick Sanan } else { /* bottom interior plane */ 1319e77315bcSPatrick Sanan 1320e77315bcSPatrick Sanan tu = x[k + 1][j][i]; 1321e77315bcSPatrick Sanan au = 0.5 * (t0 + tu); 1322e77315bcSPatrick Sanan bu = PetscPowScalar(au, bm1); 1323e77315bcSPatrick Sanan /* du = bu * au; */ 1324e77315bcSPatrick Sanan du = PetscPowScalar(au, beta); 1325e77315bcSPatrick Sanan gu = coef * bu * (tu - t0); 1326e77315bcSPatrick Sanan 1327e77315bcSPatrick Sanan td = x[k - 1][j][i]; 1328e77315bcSPatrick Sanan ad = 0.5 * (t0 + td); 1329e77315bcSPatrick Sanan bd = PetscPowScalar(ad, bm1); 1330e77315bcSPatrick Sanan /* dd = bd * ad; */ 1331e77315bcSPatrick Sanan dd = PetscPowScalar(ad, beta); 1332e77315bcSPatrick Sanan gd = coef * bd * (td - t0); 1333e77315bcSPatrick Sanan 13349371c9d4SSatish Balay c[0].k = k - 1; 13359371c9d4SSatish Balay c[0].j = j; 13369371c9d4SSatish Balay c[0].i = i; 1337e77315bcSPatrick Sanan v[0] = -hxhydhz * (dd - gd); 13389371c9d4SSatish Balay c[1].k = k; 13399371c9d4SSatish Balay c[1].j = j; 13409371c9d4SSatish Balay c[1].i = i - 1; 1341e77315bcSPatrick Sanan v[1] = -hyhzdhx * (dw - gw); 13429371c9d4SSatish Balay c[2].k = k; 13439371c9d4SSatish Balay c[2].j = j; 13449371c9d4SSatish Balay c[2].i = i; 1345e77315bcSPatrick Sanan v[2] = hzhxdhy * (dn - gn) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (dd + du + gd - gu); 13469371c9d4SSatish Balay c[3].k = k; 13479371c9d4SSatish Balay c[3].j = j; 13489371c9d4SSatish Balay c[3].i = i + 1; 1349e77315bcSPatrick Sanan v[3] = -hyhzdhx * (de + ge); 13509371c9d4SSatish Balay c[4].k = k; 13519371c9d4SSatish Balay c[4].j = j + 1; 13529371c9d4SSatish Balay c[4].i = i; 1353e77315bcSPatrick Sanan v[4] = -hzhxdhy * (dn + gn); 13549371c9d4SSatish Balay c[5].k = k + 1; 13559371c9d4SSatish Balay c[5].j = j; 13569371c9d4SSatish Balay c[5].i = i; 1357e77315bcSPatrick Sanan v[5] = -hxhydhz * (du + gu); 13589566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac, 1, &row, 6, c, v, INSERT_VALUES)); 1359e77315bcSPatrick Sanan } 1360e77315bcSPatrick Sanan 1361e77315bcSPatrick Sanan } else if (j == my - 1) { 1362e77315bcSPatrick Sanan tw = x[k][j][i - 1]; 1363e77315bcSPatrick Sanan aw = 0.5 * (t0 + tw); 1364e77315bcSPatrick Sanan bw = PetscPowScalar(aw, bm1); 1365e77315bcSPatrick Sanan /* dw = bw * aw */ 1366e77315bcSPatrick Sanan dw = PetscPowScalar(aw, beta); 1367e77315bcSPatrick Sanan gw = coef * bw * (t0 - tw); 1368e77315bcSPatrick Sanan 1369e77315bcSPatrick Sanan te = x[k][j][i + 1]; 1370e77315bcSPatrick Sanan ae = 0.5 * (t0 + te); 1371e77315bcSPatrick Sanan be = PetscPowScalar(ae, bm1); 1372e77315bcSPatrick Sanan /* de = be * ae; */ 1373e77315bcSPatrick Sanan de = PetscPowScalar(ae, beta); 1374e77315bcSPatrick Sanan ge = coef * be * (te - t0); 1375e77315bcSPatrick Sanan 1376e77315bcSPatrick Sanan ts = x[k][j - 1][i]; 1377e77315bcSPatrick Sanan as = 0.5 * (t0 + ts); 1378e77315bcSPatrick Sanan bs = PetscPowScalar(as, bm1); 1379e77315bcSPatrick Sanan /* ds = bs * as; */ 1380e77315bcSPatrick Sanan ds = PetscPowScalar(as, beta); 1381e77315bcSPatrick Sanan gs = coef * bs * (t0 - ts); 1382e77315bcSPatrick Sanan 1383e77315bcSPatrick Sanan /* top down interior edge */ 1384e77315bcSPatrick Sanan if (k == 0) { 1385e77315bcSPatrick Sanan tu = x[k + 1][j][i]; 1386e77315bcSPatrick Sanan au = 0.5 * (t0 + tu); 1387e77315bcSPatrick Sanan bu = PetscPowScalar(au, bm1); 1388e77315bcSPatrick Sanan /* du = bu * au; */ 1389e77315bcSPatrick Sanan du = PetscPowScalar(au, beta); 1390e77315bcSPatrick Sanan gu = coef * bu * (tu - t0); 1391e77315bcSPatrick Sanan 13929371c9d4SSatish Balay c[0].k = k; 13939371c9d4SSatish Balay c[0].j = j - 1; 13949371c9d4SSatish Balay c[0].i = i; 1395e77315bcSPatrick Sanan v[0] = -hzhxdhy * (ds - gs); 13969371c9d4SSatish Balay c[1].k = k; 13979371c9d4SSatish Balay c[1].j = j; 13989371c9d4SSatish Balay c[1].i = i - 1; 1399e77315bcSPatrick Sanan v[1] = -hyhzdhx * (dw - gw); 14009371c9d4SSatish Balay c[2].k = k; 14019371c9d4SSatish Balay c[2].j = j; 14029371c9d4SSatish Balay c[2].i = i; 1403e77315bcSPatrick Sanan v[2] = hzhxdhy * (ds + gs) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (du - gu); 14049371c9d4SSatish Balay c[3].k = k; 14059371c9d4SSatish Balay c[3].j = j; 14069371c9d4SSatish Balay c[3].i = i + 1; 1407e77315bcSPatrick Sanan v[3] = -hyhzdhx * (de + ge); 14089371c9d4SSatish Balay c[4].k = k + 1; 14099371c9d4SSatish Balay c[4].j = j; 14109371c9d4SSatish Balay c[4].i = i; 1411e77315bcSPatrick Sanan v[4] = -hxhydhz * (du + gu); 14129566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac, 1, &row, 5, c, v, INSERT_VALUES)); 1413e77315bcSPatrick Sanan 1414e77315bcSPatrick Sanan } else if (k == mz - 1) { /* top up interior edge */ 1415e77315bcSPatrick Sanan 1416e77315bcSPatrick Sanan td = x[k - 1][j][i]; 1417e77315bcSPatrick Sanan ad = 0.5 * (t0 + td); 1418e77315bcSPatrick Sanan bd = PetscPowScalar(ad, bm1); 1419e77315bcSPatrick Sanan /* dd = bd * ad; */ 1420e77315bcSPatrick Sanan dd = PetscPowScalar(ad, beta); 1421e77315bcSPatrick Sanan gd = coef * bd * (td - t0); 1422e77315bcSPatrick Sanan 14239371c9d4SSatish Balay c[0].k = k - 1; 14249371c9d4SSatish Balay c[0].j = j; 14259371c9d4SSatish Balay c[0].i = i; 1426e77315bcSPatrick Sanan v[0] = -hxhydhz * (dd - gd); 14279371c9d4SSatish Balay c[1].k = k; 14289371c9d4SSatish Balay c[1].j = j - 1; 14299371c9d4SSatish Balay c[1].i = i; 1430e77315bcSPatrick Sanan v[1] = -hzhxdhy * (ds - gs); 14319371c9d4SSatish Balay c[2].k = k; 14329371c9d4SSatish Balay c[2].j = j; 14339371c9d4SSatish Balay c[2].i = i - 1; 1434e77315bcSPatrick Sanan v[2] = -hyhzdhx * (dw - gw); 14359371c9d4SSatish Balay c[3].k = k; 14369371c9d4SSatish Balay c[3].j = j; 14379371c9d4SSatish Balay c[3].i = i; 1438e77315bcSPatrick Sanan v[3] = hzhxdhy * (ds + gs) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (dd + gd); 14399371c9d4SSatish Balay c[4].k = k; 14409371c9d4SSatish Balay c[4].j = j; 14419371c9d4SSatish Balay c[4].i = i + 1; 1442e77315bcSPatrick Sanan v[4] = -hyhzdhx * (de + ge); 14439566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac, 1, &row, 5, c, v, INSERT_VALUES)); 1444e77315bcSPatrick Sanan 1445e77315bcSPatrick Sanan } else { /* top interior plane */ 1446e77315bcSPatrick Sanan 1447e77315bcSPatrick Sanan tu = x[k + 1][j][i]; 1448e77315bcSPatrick Sanan au = 0.5 * (t0 + tu); 1449e77315bcSPatrick Sanan bu = PetscPowScalar(au, bm1); 1450e77315bcSPatrick Sanan /* du = bu * au; */ 1451e77315bcSPatrick Sanan du = PetscPowScalar(au, beta); 1452e77315bcSPatrick Sanan gu = coef * bu * (tu - t0); 1453e77315bcSPatrick Sanan 1454e77315bcSPatrick Sanan td = x[k - 1][j][i]; 1455e77315bcSPatrick Sanan ad = 0.5 * (t0 + td); 1456e77315bcSPatrick Sanan bd = PetscPowScalar(ad, bm1); 1457e77315bcSPatrick Sanan /* dd = bd * ad; */ 1458e77315bcSPatrick Sanan dd = PetscPowScalar(ad, beta); 1459e77315bcSPatrick Sanan gd = coef * bd * (td - t0); 1460e77315bcSPatrick Sanan 14619371c9d4SSatish Balay c[0].k = k - 1; 14629371c9d4SSatish Balay c[0].j = j; 14639371c9d4SSatish Balay c[0].i = i; 1464e77315bcSPatrick Sanan v[0] = -hxhydhz * (dd - gd); 14659371c9d4SSatish Balay c[1].k = k; 14669371c9d4SSatish Balay c[1].j = j - 1; 14679371c9d4SSatish Balay c[1].i = i; 1468e77315bcSPatrick Sanan v[1] = -hzhxdhy * (ds - gs); 14699371c9d4SSatish Balay c[2].k = k; 14709371c9d4SSatish Balay c[2].j = j; 14719371c9d4SSatish Balay c[2].i = i - 1; 1472e77315bcSPatrick Sanan v[2] = -hyhzdhx * (dw - gw); 14739371c9d4SSatish Balay c[3].k = k; 14749371c9d4SSatish Balay c[3].j = j; 14759371c9d4SSatish Balay c[3].i = i; 1476e77315bcSPatrick Sanan v[3] = hzhxdhy * (ds + gs) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (dd + du + gd - gu); 14779371c9d4SSatish Balay c[4].k = k; 14789371c9d4SSatish Balay c[4].j = j; 14799371c9d4SSatish Balay c[4].i = i + 1; 1480e77315bcSPatrick Sanan v[4] = -hyhzdhx * (de + ge); 14819371c9d4SSatish Balay c[5].k = k + 1; 14829371c9d4SSatish Balay c[5].j = j; 14839371c9d4SSatish Balay c[5].i = i; 1484e77315bcSPatrick Sanan v[5] = -hxhydhz * (du + gu); 14859566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac, 1, &row, 6, c, v, INSERT_VALUES)); 1486e77315bcSPatrick Sanan } 1487e77315bcSPatrick Sanan 1488e77315bcSPatrick Sanan } else if (k == 0) { 1489e77315bcSPatrick Sanan /* down interior plane */ 1490e77315bcSPatrick Sanan 1491e77315bcSPatrick Sanan tw = x[k][j][i - 1]; 1492e77315bcSPatrick Sanan aw = 0.5 * (t0 + tw); 1493e77315bcSPatrick Sanan bw = PetscPowScalar(aw, bm1); 1494e77315bcSPatrick Sanan /* dw = bw * aw */ 1495e77315bcSPatrick Sanan dw = PetscPowScalar(aw, beta); 1496e77315bcSPatrick Sanan gw = coef * bw * (t0 - tw); 1497e77315bcSPatrick Sanan 1498e77315bcSPatrick Sanan te = x[k][j][i + 1]; 1499e77315bcSPatrick Sanan ae = 0.5 * (t0 + te); 1500e77315bcSPatrick Sanan be = PetscPowScalar(ae, bm1); 1501e77315bcSPatrick Sanan /* de = be * ae; */ 1502e77315bcSPatrick Sanan de = PetscPowScalar(ae, beta); 1503e77315bcSPatrick Sanan ge = coef * be * (te - t0); 1504e77315bcSPatrick Sanan 1505e77315bcSPatrick Sanan ts = x[k][j - 1][i]; 1506e77315bcSPatrick Sanan as = 0.5 * (t0 + ts); 1507e77315bcSPatrick Sanan bs = PetscPowScalar(as, bm1); 1508e77315bcSPatrick Sanan /* ds = bs * as; */ 1509e77315bcSPatrick Sanan ds = PetscPowScalar(as, beta); 1510e77315bcSPatrick Sanan gs = coef * bs * (t0 - ts); 1511e77315bcSPatrick Sanan 1512e77315bcSPatrick Sanan tn = x[k][j + 1][i]; 1513e77315bcSPatrick Sanan an = 0.5 * (t0 + tn); 1514e77315bcSPatrick Sanan bn = PetscPowScalar(an, bm1); 1515e77315bcSPatrick Sanan /* dn = bn * an; */ 1516e77315bcSPatrick Sanan dn = PetscPowScalar(an, beta); 1517e77315bcSPatrick Sanan gn = coef * bn * (tn - t0); 1518e77315bcSPatrick Sanan 1519e77315bcSPatrick Sanan tu = x[k + 1][j][i]; 1520e77315bcSPatrick Sanan au = 0.5 * (t0 + tu); 1521e77315bcSPatrick Sanan bu = PetscPowScalar(au, bm1); 1522e77315bcSPatrick Sanan /* du = bu * au; */ 1523e77315bcSPatrick Sanan du = PetscPowScalar(au, beta); 1524e77315bcSPatrick Sanan gu = coef * bu * (tu - t0); 1525e77315bcSPatrick Sanan 15269371c9d4SSatish Balay c[0].k = k; 15279371c9d4SSatish Balay c[0].j = j - 1; 15289371c9d4SSatish Balay c[0].i = i; 1529e77315bcSPatrick Sanan v[0] = -hzhxdhy * (ds - gs); 15309371c9d4SSatish Balay c[1].k = k; 15319371c9d4SSatish Balay c[1].j = j; 15329371c9d4SSatish Balay c[1].i = i - 1; 1533e77315bcSPatrick Sanan v[1] = -hyhzdhx * (dw - gw); 15349371c9d4SSatish Balay c[2].k = k; 15359371c9d4SSatish Balay c[2].j = j; 15369371c9d4SSatish Balay c[2].i = i; 1537e77315bcSPatrick Sanan v[2] = hzhxdhy * (ds + dn + gs - gn) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (du - gu); 15389371c9d4SSatish Balay c[3].k = k; 15399371c9d4SSatish Balay c[3].j = j; 15409371c9d4SSatish Balay c[3].i = i + 1; 1541e77315bcSPatrick Sanan v[3] = -hyhzdhx * (de + ge); 15429371c9d4SSatish Balay c[4].k = k; 15439371c9d4SSatish Balay c[4].j = j + 1; 15449371c9d4SSatish Balay c[4].i = i; 1545e77315bcSPatrick Sanan v[4] = -hzhxdhy * (dn + gn); 15469371c9d4SSatish Balay c[5].k = k + 1; 15479371c9d4SSatish Balay c[5].j = j; 15489371c9d4SSatish Balay c[5].i = i; 1549e77315bcSPatrick Sanan v[5] = -hxhydhz * (du + gu); 15509566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac, 1, &row, 6, c, v, INSERT_VALUES)); 1551e77315bcSPatrick Sanan 1552e77315bcSPatrick Sanan } else if (k == mz - 1) { 1553e77315bcSPatrick Sanan /* up interior plane */ 1554e77315bcSPatrick Sanan 1555e77315bcSPatrick Sanan tw = x[k][j][i - 1]; 1556e77315bcSPatrick Sanan aw = 0.5 * (t0 + tw); 1557e77315bcSPatrick Sanan bw = PetscPowScalar(aw, bm1); 1558e77315bcSPatrick Sanan /* dw = bw * aw */ 1559e77315bcSPatrick Sanan dw = PetscPowScalar(aw, beta); 1560e77315bcSPatrick Sanan gw = coef * bw * (t0 - tw); 1561e77315bcSPatrick Sanan 1562e77315bcSPatrick Sanan te = x[k][j][i + 1]; 1563e77315bcSPatrick Sanan ae = 0.5 * (t0 + te); 1564e77315bcSPatrick Sanan be = PetscPowScalar(ae, bm1); 1565e77315bcSPatrick Sanan /* de = be * ae; */ 1566e77315bcSPatrick Sanan de = PetscPowScalar(ae, beta); 1567e77315bcSPatrick Sanan ge = coef * be * (te - t0); 1568e77315bcSPatrick Sanan 1569e77315bcSPatrick Sanan ts = x[k][j - 1][i]; 1570e77315bcSPatrick Sanan as = 0.5 * (t0 + ts); 1571e77315bcSPatrick Sanan bs = PetscPowScalar(as, bm1); 1572e77315bcSPatrick Sanan /* ds = bs * as; */ 1573e77315bcSPatrick Sanan ds = PetscPowScalar(as, beta); 1574e77315bcSPatrick Sanan gs = coef * bs * (t0 - ts); 1575e77315bcSPatrick Sanan 1576e77315bcSPatrick Sanan tn = x[k][j + 1][i]; 1577e77315bcSPatrick Sanan an = 0.5 * (t0 + tn); 1578e77315bcSPatrick Sanan bn = PetscPowScalar(an, bm1); 1579e77315bcSPatrick Sanan /* dn = bn * an; */ 1580e77315bcSPatrick Sanan dn = PetscPowScalar(an, beta); 1581e77315bcSPatrick Sanan gn = coef * bn * (tn - t0); 1582e77315bcSPatrick Sanan 1583e77315bcSPatrick Sanan td = x[k - 1][j][i]; 1584e77315bcSPatrick Sanan ad = 0.5 * (t0 + td); 1585e77315bcSPatrick Sanan bd = PetscPowScalar(ad, bm1); 1586e77315bcSPatrick Sanan /* dd = bd * ad; */ 1587e77315bcSPatrick Sanan dd = PetscPowScalar(ad, beta); 1588e77315bcSPatrick Sanan gd = coef * bd * (t0 - td); 1589e77315bcSPatrick Sanan 15909371c9d4SSatish Balay c[0].k = k - 1; 15919371c9d4SSatish Balay c[0].j = j; 15929371c9d4SSatish Balay c[0].i = i; 1593e77315bcSPatrick Sanan v[0] = -hxhydhz * (dd - gd); 15949371c9d4SSatish Balay c[1].k = k; 15959371c9d4SSatish Balay c[1].j = j - 1; 15969371c9d4SSatish Balay c[1].i = i; 1597e77315bcSPatrick Sanan v[1] = -hzhxdhy * (ds - gs); 15989371c9d4SSatish Balay c[2].k = k; 15999371c9d4SSatish Balay c[2].j = j; 16009371c9d4SSatish Balay c[2].i = i - 1; 1601e77315bcSPatrick Sanan v[2] = -hyhzdhx * (dw - gw); 16029371c9d4SSatish Balay c[3].k = k; 16039371c9d4SSatish Balay c[3].j = j; 16049371c9d4SSatish Balay c[3].i = i; 1605e77315bcSPatrick Sanan v[3] = hzhxdhy * (ds + dn + gs - gn) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (dd + gd); 16069371c9d4SSatish Balay c[4].k = k; 16079371c9d4SSatish Balay c[4].j = j; 16089371c9d4SSatish Balay c[4].i = i + 1; 1609e77315bcSPatrick Sanan v[4] = -hyhzdhx * (de + ge); 16109371c9d4SSatish Balay c[5].k = k; 16119371c9d4SSatish Balay c[5].j = j + 1; 16129371c9d4SSatish Balay c[5].i = i; 1613e77315bcSPatrick Sanan v[5] = -hzhxdhy * (dn + gn); 16149566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac, 1, &row, 6, c, v, INSERT_VALUES)); 1615e77315bcSPatrick Sanan } 1616e77315bcSPatrick Sanan } 1617e77315bcSPatrick Sanan } 1618e77315bcSPatrick Sanan } 16199566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY)); 16209566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY)); 1621e77315bcSPatrick Sanan if (jac != J) { 16229566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 16239566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 1624e77315bcSPatrick Sanan } 16259566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, localX, &x)); 16269566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &localX)); 1627e77315bcSPatrick Sanan 16289566063dSJacob Faibussowitsch PetscCall(PetscLogFlops((41.0 + 8.0 * POWFLOP) * xm * ym)); 1629e77315bcSPatrick Sanan PetscFunctionReturn(0); 1630e77315bcSPatrick Sanan } 1631e77315bcSPatrick Sanan 1632e77315bcSPatrick Sanan /*TEST 1633e77315bcSPatrick Sanan 1634e77315bcSPatrick Sanan test: 1635e77315bcSPatrick Sanan nsize: 4 1636e77315bcSPatrick Sanan args: -snes_monitor_short -pc_mg_type full -ksp_type fgmres -pc_type mg -snes_view -pc_mg_levels 2 -pc_mg_galerkin pmat 1637e77315bcSPatrick Sanan requires: !single 1638e77315bcSPatrick Sanan 1639e77315bcSPatrick Sanan TEST*/ 1640