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 52d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 53d71ae5a4SJacob Faibussowitsch { 54e77315bcSPatrick Sanan SNES snes; 55e77315bcSPatrick Sanan AppCtx user; 56e77315bcSPatrick Sanan PetscInt its, lits; 57e77315bcSPatrick Sanan PetscReal litspit; 58e77315bcSPatrick Sanan DM da; 59e77315bcSPatrick Sanan 60327415f7SBarry Smith PetscFunctionBeginUser; 619566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 62e77315bcSPatrick Sanan /* set problem parameters */ 63e77315bcSPatrick Sanan user.tleft = 1.0; 64e77315bcSPatrick Sanan user.tright = 0.1; 65e77315bcSPatrick Sanan user.beta = 2.5; 669566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-tleft", &user.tleft, NULL)); 679566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-tright", &user.tright, NULL)); 689566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-beta", &user.beta, NULL)); 69e77315bcSPatrick Sanan user.bm1 = user.beta - 1.0; 70e77315bcSPatrick Sanan user.coef = user.beta / 2.0; 71e77315bcSPatrick Sanan 72e77315bcSPatrick Sanan /* 73e77315bcSPatrick Sanan Set the DMDA (grid structure) for the grids. 74e77315bcSPatrick Sanan */ 759566063dSJacob 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)); 769566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da)); 779566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 789566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(da, &user)); 79e77315bcSPatrick Sanan 80e77315bcSPatrick Sanan /* 81e77315bcSPatrick Sanan Create the nonlinear solver 82e77315bcSPatrick Sanan */ 839566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 849566063dSJacob Faibussowitsch PetscCall(SNESSetDM(snes, da)); 859566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes, NULL, FormFunction, &user)); 869566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, NULL, NULL, FormJacobian, &user)); 879566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 889566063dSJacob Faibussowitsch PetscCall(SNESSetComputeInitialGuess(snes, FormInitialGuess, NULL)); 89e77315bcSPatrick Sanan 909566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, NULL)); 919566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(snes, &its)); 929566063dSJacob Faibussowitsch PetscCall(SNESGetLinearSolveIterations(snes, &lits)); 93e77315bcSPatrick Sanan litspit = ((PetscReal)lits) / ((PetscReal)its); 9463a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Number of SNES iterations = %" PetscInt_FMT "\n", its)); 9563a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Number of Linear iterations = %" PetscInt_FMT "\n", lits)); 969566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Average Linear its / SNES = %e\n", (double)litspit)); 97e77315bcSPatrick Sanan 989566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 999566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da)); 1009566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 101b122ec5aSJacob Faibussowitsch return 0; 102e77315bcSPatrick Sanan } 103e77315bcSPatrick Sanan /* -------------------- Form initial approximation ----------------- */ 104d71ae5a4SJacob Faibussowitsch PetscErrorCode FormInitialGuess(SNES snes, Vec X, void *ctx) 105d71ae5a4SJacob Faibussowitsch { 106e77315bcSPatrick Sanan AppCtx *user; 107e77315bcSPatrick Sanan PetscInt i, j, k, xs, ys, xm, ym, zs, zm; 108e77315bcSPatrick Sanan PetscScalar ***x; 109e77315bcSPatrick Sanan DM da; 110e77315bcSPatrick Sanan 111e77315bcSPatrick Sanan PetscFunctionBeginUser; 1129566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &da)); 1139566063dSJacob Faibussowitsch PetscCall(DMGetApplicationContext(da, &user)); 1149566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm)); 1159566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, X, &x)); 116e77315bcSPatrick Sanan 117e77315bcSPatrick Sanan /* Compute initial guess */ 118e77315bcSPatrick Sanan for (k = zs; k < zs + zm; k++) { 119e77315bcSPatrick Sanan for (j = ys; j < ys + ym; j++) { 120ad540459SPierre Jolivet for (i = xs; i < xs + xm; i++) x[k][j][i] = user->tleft; 121e77315bcSPatrick Sanan } 122e77315bcSPatrick Sanan } 1239566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, X, &x)); 124*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 125e77315bcSPatrick Sanan } 126e77315bcSPatrick Sanan /* -------------------- Evaluate Function F(x) --------------------- */ 127d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunction(SNES snes, Vec X, Vec F, void *ptr) 128d71ae5a4SJacob Faibussowitsch { 129e77315bcSPatrick Sanan AppCtx *user = (AppCtx *)ptr; 130e77315bcSPatrick Sanan PetscInt i, j, k, mx, my, mz, xs, ys, zs, xm, ym, zm; 131e77315bcSPatrick Sanan PetscScalar zero = 0.0, one = 1.0; 132e77315bcSPatrick Sanan PetscScalar hx, hy, hz, hxhydhz, hyhzdhx, hzhxdhy; 133e77315bcSPatrick 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; 134e77315bcSPatrick Sanan PetscScalar tleft, tright, beta, td, ad, dd, fd = 0.0, tu, au, du = 0.0, fu = 0.0; 135e77315bcSPatrick Sanan PetscScalar ***x, ***f; 136e77315bcSPatrick Sanan Vec localX; 137e77315bcSPatrick Sanan DM da; 138e77315bcSPatrick Sanan 139e77315bcSPatrick Sanan PetscFunctionBeginUser; 1409566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &da)); 1419566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &localX)); 1429566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, NULL, &mx, &my, &mz, 0, 0, 0, 0, 0, 0, 0, 0, 0)); 1439371c9d4SSatish Balay hx = one / (PetscReal)(mx - 1); 1449371c9d4SSatish Balay hy = one / (PetscReal)(my - 1); 1459371c9d4SSatish Balay hz = one / (PetscReal)(mz - 1); 1469371c9d4SSatish Balay hxhydhz = hx * hy / hz; 1479371c9d4SSatish Balay hyhzdhx = hy * hz / hx; 1489371c9d4SSatish Balay hzhxdhy = hz * hx / hy; 1499371c9d4SSatish Balay tleft = user->tleft; 1509371c9d4SSatish Balay tright = user->tright; 151e77315bcSPatrick Sanan beta = user->beta; 152e77315bcSPatrick Sanan 153e77315bcSPatrick Sanan /* Get ghost points */ 1549566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, localX)); 1559566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, localX)); 1569566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm)); 1579566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, localX, &x)); 1589566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, F, &f)); 159e77315bcSPatrick Sanan 160e77315bcSPatrick Sanan /* Evaluate function */ 161e77315bcSPatrick Sanan for (k = zs; k < zs + zm; k++) { 162e77315bcSPatrick Sanan for (j = ys; j < ys + ym; j++) { 163e77315bcSPatrick Sanan for (i = xs; i < xs + xm; i++) { 164e77315bcSPatrick Sanan t0 = x[k][j][i]; 165e77315bcSPatrick Sanan 166e77315bcSPatrick Sanan if (i > 0 && i < mx - 1 && j > 0 && j < my - 1 && k > 0 && k < mz - 1) { 167e77315bcSPatrick Sanan /* general interior volume */ 168e77315bcSPatrick Sanan 169e77315bcSPatrick Sanan tw = x[k][j][i - 1]; 170e77315bcSPatrick Sanan aw = 0.5 * (t0 + tw); 171e77315bcSPatrick Sanan dw = PetscPowScalar(aw, beta); 172e77315bcSPatrick Sanan fw = dw * (t0 - tw); 173e77315bcSPatrick Sanan 174e77315bcSPatrick Sanan te = x[k][j][i + 1]; 175e77315bcSPatrick Sanan ae = 0.5 * (t0 + te); 176e77315bcSPatrick Sanan de = PetscPowScalar(ae, beta); 177e77315bcSPatrick Sanan fe = de * (te - t0); 178e77315bcSPatrick Sanan 179e77315bcSPatrick Sanan ts = x[k][j - 1][i]; 180e77315bcSPatrick Sanan as = 0.5 * (t0 + ts); 181e77315bcSPatrick Sanan ds = PetscPowScalar(as, beta); 182e77315bcSPatrick Sanan fs = ds * (t0 - ts); 183e77315bcSPatrick Sanan 184e77315bcSPatrick Sanan tn = x[k][j + 1][i]; 185e77315bcSPatrick Sanan an = 0.5 * (t0 + tn); 186e77315bcSPatrick Sanan dn = PetscPowScalar(an, beta); 187e77315bcSPatrick Sanan fn = dn * (tn - t0); 188e77315bcSPatrick Sanan 189e77315bcSPatrick Sanan td = x[k - 1][j][i]; 190e77315bcSPatrick Sanan ad = 0.5 * (t0 + td); 191e77315bcSPatrick Sanan dd = PetscPowScalar(ad, beta); 192e77315bcSPatrick Sanan fd = dd * (t0 - td); 193e77315bcSPatrick Sanan 194e77315bcSPatrick Sanan tu = x[k + 1][j][i]; 195e77315bcSPatrick Sanan au = 0.5 * (t0 + tu); 196e77315bcSPatrick Sanan du = PetscPowScalar(au, beta); 197e77315bcSPatrick Sanan fu = du * (tu - t0); 198e77315bcSPatrick Sanan 199e77315bcSPatrick Sanan } else if (i == 0) { 200e77315bcSPatrick Sanan /* left-hand (west) boundary */ 201e77315bcSPatrick Sanan tw = tleft; 202e77315bcSPatrick Sanan aw = 0.5 * (t0 + tw); 203e77315bcSPatrick Sanan dw = PetscPowScalar(aw, beta); 204e77315bcSPatrick Sanan fw = dw * (t0 - tw); 205e77315bcSPatrick Sanan 206e77315bcSPatrick Sanan te = x[k][j][i + 1]; 207e77315bcSPatrick Sanan ae = 0.5 * (t0 + te); 208e77315bcSPatrick Sanan de = PetscPowScalar(ae, beta); 209e77315bcSPatrick Sanan fe = de * (te - t0); 210e77315bcSPatrick Sanan 211e77315bcSPatrick Sanan if (j > 0) { 212e77315bcSPatrick Sanan ts = x[k][j - 1][i]; 213e77315bcSPatrick Sanan as = 0.5 * (t0 + ts); 214e77315bcSPatrick Sanan ds = PetscPowScalar(as, beta); 215e77315bcSPatrick Sanan fs = ds * (t0 - ts); 216e77315bcSPatrick Sanan } else { 217e77315bcSPatrick Sanan fs = zero; 218e77315bcSPatrick Sanan } 219e77315bcSPatrick Sanan 220e77315bcSPatrick Sanan if (j < my - 1) { 221e77315bcSPatrick Sanan tn = x[k][j + 1][i]; 222e77315bcSPatrick Sanan an = 0.5 * (t0 + tn); 223e77315bcSPatrick Sanan dn = PetscPowScalar(an, beta); 224e77315bcSPatrick Sanan fn = dn * (tn - t0); 225e77315bcSPatrick Sanan } else { 226e77315bcSPatrick Sanan fn = zero; 227e77315bcSPatrick Sanan } 228e77315bcSPatrick Sanan 229e77315bcSPatrick Sanan if (k > 0) { 230e77315bcSPatrick Sanan td = x[k - 1][j][i]; 231e77315bcSPatrick Sanan ad = 0.5 * (t0 + td); 232e77315bcSPatrick Sanan dd = PetscPowScalar(ad, beta); 233e77315bcSPatrick Sanan fd = dd * (t0 - td); 234e77315bcSPatrick Sanan } else { 235e77315bcSPatrick Sanan fd = zero; 236e77315bcSPatrick Sanan } 237e77315bcSPatrick Sanan 238e77315bcSPatrick Sanan if (k < mz - 1) { 239e77315bcSPatrick Sanan tu = x[k + 1][j][i]; 240e77315bcSPatrick Sanan au = 0.5 * (t0 + tu); 241e77315bcSPatrick Sanan du = PetscPowScalar(au, beta); 242e77315bcSPatrick Sanan fu = du * (tu - t0); 243e77315bcSPatrick Sanan } else { 244e77315bcSPatrick Sanan fu = zero; 245e77315bcSPatrick Sanan } 246e77315bcSPatrick Sanan 247e77315bcSPatrick Sanan } else if (i == mx - 1) { 248e77315bcSPatrick Sanan /* right-hand (east) boundary */ 249e77315bcSPatrick Sanan tw = x[k][j][i - 1]; 250e77315bcSPatrick Sanan aw = 0.5 * (t0 + tw); 251e77315bcSPatrick Sanan dw = PetscPowScalar(aw, beta); 252e77315bcSPatrick Sanan fw = dw * (t0 - tw); 253e77315bcSPatrick Sanan 254e77315bcSPatrick Sanan te = tright; 255e77315bcSPatrick Sanan ae = 0.5 * (t0 + te); 256e77315bcSPatrick Sanan de = PetscPowScalar(ae, beta); 257e77315bcSPatrick Sanan fe = de * (te - t0); 258e77315bcSPatrick Sanan 259e77315bcSPatrick Sanan if (j > 0) { 260e77315bcSPatrick Sanan ts = x[k][j - 1][i]; 261e77315bcSPatrick Sanan as = 0.5 * (t0 + ts); 262e77315bcSPatrick Sanan ds = PetscPowScalar(as, beta); 263e77315bcSPatrick Sanan fs = ds * (t0 - ts); 264e77315bcSPatrick Sanan } else { 265e77315bcSPatrick Sanan fs = zero; 266e77315bcSPatrick Sanan } 267e77315bcSPatrick Sanan 268e77315bcSPatrick Sanan if (j < my - 1) { 269e77315bcSPatrick Sanan tn = x[k][j + 1][i]; 270e77315bcSPatrick Sanan an = 0.5 * (t0 + tn); 271e77315bcSPatrick Sanan dn = PetscPowScalar(an, beta); 272e77315bcSPatrick Sanan fn = dn * (tn - t0); 273e77315bcSPatrick Sanan } else { 274e77315bcSPatrick Sanan fn = zero; 275e77315bcSPatrick Sanan } 276e77315bcSPatrick Sanan 277e77315bcSPatrick Sanan if (k > 0) { 278e77315bcSPatrick Sanan td = x[k - 1][j][i]; 279e77315bcSPatrick Sanan ad = 0.5 * (t0 + td); 280e77315bcSPatrick Sanan dd = PetscPowScalar(ad, beta); 281e77315bcSPatrick Sanan fd = dd * (t0 - td); 282e77315bcSPatrick Sanan } else { 283e77315bcSPatrick Sanan fd = zero; 284e77315bcSPatrick Sanan } 285e77315bcSPatrick Sanan 286e77315bcSPatrick Sanan if (k < mz - 1) { 287e77315bcSPatrick Sanan tu = x[k + 1][j][i]; 288e77315bcSPatrick Sanan au = 0.5 * (t0 + tu); 289e77315bcSPatrick Sanan du = PetscPowScalar(au, beta); 290e77315bcSPatrick Sanan fu = du * (tu - t0); 291e77315bcSPatrick Sanan } else { 292e77315bcSPatrick Sanan fu = zero; 293e77315bcSPatrick Sanan } 294e77315bcSPatrick Sanan 295e77315bcSPatrick Sanan } else if (j == 0) { 296e77315bcSPatrick Sanan /* bottom (south) boundary, and i <> 0 or mx-1 */ 297e77315bcSPatrick Sanan tw = x[k][j][i - 1]; 298e77315bcSPatrick Sanan aw = 0.5 * (t0 + tw); 299e77315bcSPatrick Sanan dw = PetscPowScalar(aw, beta); 300e77315bcSPatrick Sanan fw = dw * (t0 - tw); 301e77315bcSPatrick Sanan 302e77315bcSPatrick Sanan te = x[k][j][i + 1]; 303e77315bcSPatrick Sanan ae = 0.5 * (t0 + te); 304e77315bcSPatrick Sanan de = PetscPowScalar(ae, beta); 305e77315bcSPatrick Sanan fe = de * (te - t0); 306e77315bcSPatrick Sanan 307e77315bcSPatrick Sanan fs = zero; 308e77315bcSPatrick Sanan 309e77315bcSPatrick Sanan tn = x[k][j + 1][i]; 310e77315bcSPatrick Sanan an = 0.5 * (t0 + tn); 311e77315bcSPatrick Sanan dn = PetscPowScalar(an, beta); 312e77315bcSPatrick Sanan fn = dn * (tn - t0); 313e77315bcSPatrick Sanan 314e77315bcSPatrick Sanan if (k > 0) { 315e77315bcSPatrick Sanan td = x[k - 1][j][i]; 316e77315bcSPatrick Sanan ad = 0.5 * (t0 + td); 317e77315bcSPatrick Sanan dd = PetscPowScalar(ad, beta); 318e77315bcSPatrick Sanan fd = dd * (t0 - td); 319e77315bcSPatrick Sanan } else { 320e77315bcSPatrick Sanan fd = zero; 321e77315bcSPatrick Sanan } 322e77315bcSPatrick Sanan 323e77315bcSPatrick Sanan if (k < mz - 1) { 324e77315bcSPatrick Sanan tu = x[k + 1][j][i]; 325e77315bcSPatrick Sanan au = 0.5 * (t0 + tu); 326e77315bcSPatrick Sanan du = PetscPowScalar(au, beta); 327e77315bcSPatrick Sanan fu = du * (tu - t0); 328e77315bcSPatrick Sanan } else { 329e77315bcSPatrick Sanan fu = zero; 330e77315bcSPatrick Sanan } 331e77315bcSPatrick Sanan 332e77315bcSPatrick Sanan } else if (j == my - 1) { 333e77315bcSPatrick Sanan /* top (north) boundary, and i <> 0 or mx-1 */ 334e77315bcSPatrick Sanan tw = x[k][j][i - 1]; 335e77315bcSPatrick Sanan aw = 0.5 * (t0 + tw); 336e77315bcSPatrick Sanan dw = PetscPowScalar(aw, beta); 337e77315bcSPatrick Sanan fw = dw * (t0 - tw); 338e77315bcSPatrick Sanan 339e77315bcSPatrick Sanan te = x[k][j][i + 1]; 340e77315bcSPatrick Sanan ae = 0.5 * (t0 + te); 341e77315bcSPatrick Sanan de = PetscPowScalar(ae, beta); 342e77315bcSPatrick Sanan fe = de * (te - t0); 343e77315bcSPatrick Sanan 344e77315bcSPatrick Sanan ts = x[k][j - 1][i]; 345e77315bcSPatrick Sanan as = 0.5 * (t0 + ts); 346e77315bcSPatrick Sanan ds = PetscPowScalar(as, beta); 347e77315bcSPatrick Sanan fs = ds * (t0 - ts); 348e77315bcSPatrick Sanan 349e77315bcSPatrick Sanan fn = zero; 350e77315bcSPatrick Sanan 351e77315bcSPatrick Sanan if (k > 0) { 352e77315bcSPatrick Sanan td = x[k - 1][j][i]; 353e77315bcSPatrick Sanan ad = 0.5 * (t0 + td); 354e77315bcSPatrick Sanan dd = PetscPowScalar(ad, beta); 355e77315bcSPatrick Sanan fd = dd * (t0 - td); 356e77315bcSPatrick Sanan } else { 357e77315bcSPatrick Sanan fd = zero; 358e77315bcSPatrick Sanan } 359e77315bcSPatrick Sanan 360e77315bcSPatrick Sanan if (k < mz - 1) { 361e77315bcSPatrick Sanan tu = x[k + 1][j][i]; 362e77315bcSPatrick Sanan au = 0.5 * (t0 + tu); 363e77315bcSPatrick Sanan du = PetscPowScalar(au, beta); 364e77315bcSPatrick Sanan fu = du * (tu - t0); 365e77315bcSPatrick Sanan } else { 366e77315bcSPatrick Sanan fu = zero; 367e77315bcSPatrick Sanan } 368e77315bcSPatrick Sanan 369e77315bcSPatrick Sanan } else if (k == 0) { 370e77315bcSPatrick Sanan /* down boundary (interior only) */ 371e77315bcSPatrick Sanan tw = x[k][j][i - 1]; 372e77315bcSPatrick Sanan aw = 0.5 * (t0 + tw); 373e77315bcSPatrick Sanan dw = PetscPowScalar(aw, beta); 374e77315bcSPatrick Sanan fw = dw * (t0 - tw); 375e77315bcSPatrick Sanan 376e77315bcSPatrick Sanan te = x[k][j][i + 1]; 377e77315bcSPatrick Sanan ae = 0.5 * (t0 + te); 378e77315bcSPatrick Sanan de = PetscPowScalar(ae, beta); 379e77315bcSPatrick Sanan fe = de * (te - t0); 380e77315bcSPatrick Sanan 381e77315bcSPatrick Sanan ts = x[k][j - 1][i]; 382e77315bcSPatrick Sanan as = 0.5 * (t0 + ts); 383e77315bcSPatrick Sanan ds = PetscPowScalar(as, beta); 384e77315bcSPatrick Sanan fs = ds * (t0 - ts); 385e77315bcSPatrick Sanan 386e77315bcSPatrick Sanan tn = x[k][j + 1][i]; 387e77315bcSPatrick Sanan an = 0.5 * (t0 + tn); 388e77315bcSPatrick Sanan dn = PetscPowScalar(an, beta); 389e77315bcSPatrick Sanan fn = dn * (tn - t0); 390e77315bcSPatrick Sanan 391e77315bcSPatrick Sanan fd = zero; 392e77315bcSPatrick Sanan 393e77315bcSPatrick Sanan tu = x[k + 1][j][i]; 394e77315bcSPatrick Sanan au = 0.5 * (t0 + tu); 395e77315bcSPatrick Sanan du = PetscPowScalar(au, beta); 396e77315bcSPatrick Sanan fu = du * (tu - t0); 397e77315bcSPatrick Sanan 398e77315bcSPatrick Sanan } else if (k == mz - 1) { 399e77315bcSPatrick Sanan /* up boundary (interior only) */ 400e77315bcSPatrick Sanan tw = x[k][j][i - 1]; 401e77315bcSPatrick Sanan aw = 0.5 * (t0 + tw); 402e77315bcSPatrick Sanan dw = PetscPowScalar(aw, beta); 403e77315bcSPatrick Sanan fw = dw * (t0 - tw); 404e77315bcSPatrick Sanan 405e77315bcSPatrick Sanan te = x[k][j][i + 1]; 406e77315bcSPatrick Sanan ae = 0.5 * (t0 + te); 407e77315bcSPatrick Sanan de = PetscPowScalar(ae, beta); 408e77315bcSPatrick Sanan fe = de * (te - t0); 409e77315bcSPatrick Sanan 410e77315bcSPatrick Sanan ts = x[k][j - 1][i]; 411e77315bcSPatrick Sanan as = 0.5 * (t0 + ts); 412e77315bcSPatrick Sanan ds = PetscPowScalar(as, beta); 413e77315bcSPatrick Sanan fs = ds * (t0 - ts); 414e77315bcSPatrick Sanan 415e77315bcSPatrick Sanan tn = x[k][j + 1][i]; 416e77315bcSPatrick Sanan an = 0.5 * (t0 + tn); 417e77315bcSPatrick Sanan dn = PetscPowScalar(an, beta); 418e77315bcSPatrick Sanan fn = dn * (tn - t0); 419e77315bcSPatrick Sanan 420e77315bcSPatrick Sanan td = x[k - 1][j][i]; 421e77315bcSPatrick Sanan ad = 0.5 * (t0 + td); 422e77315bcSPatrick Sanan dd = PetscPowScalar(ad, beta); 423e77315bcSPatrick Sanan fd = dd * (t0 - td); 424e77315bcSPatrick Sanan 425e77315bcSPatrick Sanan fu = zero; 426e77315bcSPatrick Sanan } 427e77315bcSPatrick Sanan 428e77315bcSPatrick Sanan f[k][j][i] = -hyhzdhx * (fe - fw) - hzhxdhy * (fn - fs) - hxhydhz * (fu - fd); 429e77315bcSPatrick Sanan } 430e77315bcSPatrick Sanan } 431e77315bcSPatrick Sanan } 4329566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, localX, &x)); 4339566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, F, &f)); 4349566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &localX)); 4359566063dSJacob Faibussowitsch PetscCall(PetscLogFlops((22.0 + 4.0 * POWFLOP) * ym * xm)); 436*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 437e77315bcSPatrick Sanan } 438e77315bcSPatrick Sanan /* -------------------- Evaluate Jacobian F(x) --------------------- */ 439d71ae5a4SJacob Faibussowitsch PetscErrorCode FormJacobian(SNES snes, Vec X, Mat J, Mat jac, void *ptr) 440d71ae5a4SJacob Faibussowitsch { 441e77315bcSPatrick Sanan AppCtx *user = (AppCtx *)ptr; 442e77315bcSPatrick Sanan PetscInt i, j, k, mx, my, mz, xs, ys, zs, xm, ym, zm; 443e77315bcSPatrick Sanan PetscScalar one = 1.0; 444e77315bcSPatrick Sanan PetscScalar hx, hy, hz, hxhydhz, hyhzdhx, hzhxdhy; 445e77315bcSPatrick Sanan PetscScalar t0, tn, ts, te, tw, an, as, ae, aw, dn, ds, de, dw; 446e77315bcSPatrick Sanan PetscScalar tleft, tright, beta, td, ad, dd, tu, au, du, v[7], bm1, coef; 447e77315bcSPatrick Sanan PetscScalar ***x, bn, bs, be, bw, bu, bd, gn, gs, ge, gw, gu, gd; 448e77315bcSPatrick Sanan Vec localX; 449e77315bcSPatrick Sanan MatStencil c[7], row; 450e77315bcSPatrick Sanan DM da; 451e77315bcSPatrick Sanan 452e77315bcSPatrick Sanan PetscFunctionBeginUser; 4539566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &da)); 4549566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &localX)); 4559566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, NULL, &mx, &my, &mz, 0, 0, 0, 0, 0, 0, 0, 0, 0)); 4569371c9d4SSatish Balay hx = one / (PetscReal)(mx - 1); 4579371c9d4SSatish Balay hy = one / (PetscReal)(my - 1); 4589371c9d4SSatish Balay hz = one / (PetscReal)(mz - 1); 4599371c9d4SSatish Balay hxhydhz = hx * hy / hz; 4609371c9d4SSatish Balay hyhzdhx = hy * hz / hx; 4619371c9d4SSatish Balay hzhxdhy = hz * hx / hy; 4629371c9d4SSatish Balay tleft = user->tleft; 4639371c9d4SSatish Balay tright = user->tright; 4649371c9d4SSatish Balay beta = user->beta; 4659371c9d4SSatish Balay bm1 = user->bm1; 4669371c9d4SSatish Balay coef = user->coef; 467e77315bcSPatrick Sanan 468e77315bcSPatrick Sanan /* Get ghost points */ 4699566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, localX)); 4709566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, localX)); 4719566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm)); 4729566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, localX, &x)); 473e77315bcSPatrick Sanan 474e77315bcSPatrick Sanan /* Evaluate Jacobian of function */ 475e77315bcSPatrick Sanan for (k = zs; k < zs + zm; k++) { 476e77315bcSPatrick Sanan for (j = ys; j < ys + ym; j++) { 477e77315bcSPatrick Sanan for (i = xs; i < xs + xm; i++) { 478e77315bcSPatrick Sanan t0 = x[k][j][i]; 4799371c9d4SSatish Balay row.k = k; 4809371c9d4SSatish Balay row.j = j; 4819371c9d4SSatish Balay row.i = i; 482e77315bcSPatrick Sanan if (i > 0 && i < mx - 1 && j > 0 && j < my - 1 && k > 0 && k < mz - 1) { 483e77315bcSPatrick Sanan /* general interior volume */ 484e77315bcSPatrick Sanan 485e77315bcSPatrick Sanan tw = x[k][j][i - 1]; 486e77315bcSPatrick Sanan aw = 0.5 * (t0 + tw); 487e77315bcSPatrick Sanan bw = PetscPowScalar(aw, bm1); 488e77315bcSPatrick Sanan /* dw = bw * aw */ 489e77315bcSPatrick Sanan dw = PetscPowScalar(aw, beta); 490e77315bcSPatrick Sanan gw = coef * bw * (t0 - tw); 491e77315bcSPatrick Sanan 492e77315bcSPatrick Sanan te = x[k][j][i + 1]; 493e77315bcSPatrick Sanan ae = 0.5 * (t0 + te); 494e77315bcSPatrick Sanan be = PetscPowScalar(ae, bm1); 495e77315bcSPatrick Sanan /* de = be * ae; */ 496e77315bcSPatrick Sanan de = PetscPowScalar(ae, beta); 497e77315bcSPatrick Sanan ge = coef * be * (te - t0); 498e77315bcSPatrick Sanan 499e77315bcSPatrick Sanan ts = x[k][j - 1][i]; 500e77315bcSPatrick Sanan as = 0.5 * (t0 + ts); 501e77315bcSPatrick Sanan bs = PetscPowScalar(as, bm1); 502e77315bcSPatrick Sanan /* ds = bs * as; */ 503e77315bcSPatrick Sanan ds = PetscPowScalar(as, beta); 504e77315bcSPatrick Sanan gs = coef * bs * (t0 - ts); 505e77315bcSPatrick Sanan 506e77315bcSPatrick Sanan tn = x[k][j + 1][i]; 507e77315bcSPatrick Sanan an = 0.5 * (t0 + tn); 508e77315bcSPatrick Sanan bn = PetscPowScalar(an, bm1); 509e77315bcSPatrick Sanan /* dn = bn * an; */ 510e77315bcSPatrick Sanan dn = PetscPowScalar(an, beta); 511e77315bcSPatrick Sanan gn = coef * bn * (tn - t0); 512e77315bcSPatrick Sanan 513e77315bcSPatrick Sanan td = x[k - 1][j][i]; 514e77315bcSPatrick Sanan ad = 0.5 * (t0 + td); 515e77315bcSPatrick Sanan bd = PetscPowScalar(ad, bm1); 516e77315bcSPatrick Sanan /* dd = bd * ad; */ 517e77315bcSPatrick Sanan dd = PetscPowScalar(ad, beta); 518e77315bcSPatrick Sanan gd = coef * bd * (t0 - td); 519e77315bcSPatrick Sanan 520e77315bcSPatrick Sanan tu = x[k + 1][j][i]; 521e77315bcSPatrick Sanan au = 0.5 * (t0 + tu); 522e77315bcSPatrick Sanan bu = PetscPowScalar(au, bm1); 523e77315bcSPatrick Sanan /* du = bu * au; */ 524e77315bcSPatrick Sanan du = PetscPowScalar(au, beta); 525e77315bcSPatrick Sanan gu = coef * bu * (tu - t0); 526e77315bcSPatrick Sanan 5279371c9d4SSatish Balay c[0].k = k - 1; 5289371c9d4SSatish Balay c[0].j = j; 5299371c9d4SSatish Balay c[0].i = i; 5309371c9d4SSatish Balay v[0] = -hxhydhz * (dd - gd); 5319371c9d4SSatish Balay c[1].k = k; 5329371c9d4SSatish Balay c[1].j = j - 1; 5339371c9d4SSatish Balay c[1].i = i; 534e77315bcSPatrick Sanan v[1] = -hzhxdhy * (ds - gs); 5359371c9d4SSatish Balay c[2].k = k; 5369371c9d4SSatish Balay c[2].j = j; 5379371c9d4SSatish Balay c[2].i = i - 1; 538e77315bcSPatrick Sanan v[2] = -hyhzdhx * (dw - gw); 5399371c9d4SSatish Balay c[3].k = k; 5409371c9d4SSatish Balay c[3].j = j; 5419371c9d4SSatish Balay c[3].i = i; 542e77315bcSPatrick Sanan v[3] = hzhxdhy * (ds + dn + gs - gn) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (dd + du + gd - gu); 5439371c9d4SSatish Balay c[4].k = k; 5449371c9d4SSatish Balay c[4].j = j; 5459371c9d4SSatish Balay c[4].i = i + 1; 546e77315bcSPatrick Sanan v[4] = -hyhzdhx * (de + ge); 5479371c9d4SSatish Balay c[5].k = k; 5489371c9d4SSatish Balay c[5].j = j + 1; 5499371c9d4SSatish Balay c[5].i = i; 550e77315bcSPatrick Sanan v[5] = -hzhxdhy * (dn + gn); 5519371c9d4SSatish Balay c[6].k = k + 1; 5529371c9d4SSatish Balay c[6].j = j; 5539371c9d4SSatish Balay c[6].i = i; 554e77315bcSPatrick Sanan v[6] = -hxhydhz * (du + gu); 5559566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac, 1, &row, 7, c, v, INSERT_VALUES)); 556e77315bcSPatrick Sanan 557e77315bcSPatrick Sanan } else if (i == 0) { 558e77315bcSPatrick Sanan /* left-hand plane boundary */ 559e77315bcSPatrick Sanan tw = tleft; 560e77315bcSPatrick Sanan aw = 0.5 * (t0 + tw); 561e77315bcSPatrick Sanan bw = PetscPowScalar(aw, bm1); 562e77315bcSPatrick Sanan /* dw = bw * aw */ 563e77315bcSPatrick Sanan dw = PetscPowScalar(aw, beta); 564e77315bcSPatrick Sanan gw = coef * bw * (t0 - tw); 565e77315bcSPatrick Sanan 566e77315bcSPatrick Sanan te = x[k][j][i + 1]; 567e77315bcSPatrick Sanan ae = 0.5 * (t0 + te); 568e77315bcSPatrick Sanan be = PetscPowScalar(ae, bm1); 569e77315bcSPatrick Sanan /* de = be * ae; */ 570e77315bcSPatrick Sanan de = PetscPowScalar(ae, beta); 571e77315bcSPatrick Sanan ge = coef * be * (te - t0); 572e77315bcSPatrick Sanan 573e77315bcSPatrick Sanan /* left-hand bottom edge */ 574e77315bcSPatrick Sanan if (j == 0) { 575e77315bcSPatrick Sanan tn = x[k][j + 1][i]; 576e77315bcSPatrick Sanan an = 0.5 * (t0 + tn); 577e77315bcSPatrick Sanan bn = PetscPowScalar(an, bm1); 578e77315bcSPatrick Sanan /* dn = bn * an; */ 579e77315bcSPatrick Sanan dn = PetscPowScalar(an, beta); 580e77315bcSPatrick Sanan gn = coef * bn * (tn - t0); 581e77315bcSPatrick Sanan 582e77315bcSPatrick Sanan /* left-hand bottom down corner */ 583e77315bcSPatrick Sanan if (k == 0) { 584e77315bcSPatrick Sanan tu = x[k + 1][j][i]; 585e77315bcSPatrick Sanan au = 0.5 * (t0 + tu); 586e77315bcSPatrick Sanan bu = PetscPowScalar(au, bm1); 587e77315bcSPatrick Sanan /* du = bu * au; */ 588e77315bcSPatrick Sanan du = PetscPowScalar(au, beta); 589e77315bcSPatrick Sanan gu = coef * bu * (tu - t0); 590e77315bcSPatrick Sanan 5919371c9d4SSatish Balay c[0].k = k; 5929371c9d4SSatish Balay c[0].j = j; 5939371c9d4SSatish Balay c[0].i = i; 594e77315bcSPatrick Sanan v[0] = hzhxdhy * (dn - gn) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (du - gu); 5959371c9d4SSatish Balay c[1].k = k; 5969371c9d4SSatish Balay c[1].j = j; 5979371c9d4SSatish Balay c[1].i = i + 1; 598e77315bcSPatrick Sanan v[1] = -hyhzdhx * (de + ge); 5999371c9d4SSatish Balay c[2].k = k; 6009371c9d4SSatish Balay c[2].j = j + 1; 6019371c9d4SSatish Balay c[2].i = i; 602e77315bcSPatrick Sanan v[2] = -hzhxdhy * (dn + gn); 6039371c9d4SSatish Balay c[3].k = k + 1; 6049371c9d4SSatish Balay c[3].j = j; 6059371c9d4SSatish Balay c[3].i = i; 606e77315bcSPatrick Sanan v[3] = -hxhydhz * (du + gu); 6079566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac, 1, &row, 4, c, v, INSERT_VALUES)); 608e77315bcSPatrick Sanan 609e77315bcSPatrick Sanan /* left-hand bottom interior edge */ 610e77315bcSPatrick Sanan } else if (k < mz - 1) { 611e77315bcSPatrick Sanan tu = x[k + 1][j][i]; 612e77315bcSPatrick Sanan au = 0.5 * (t0 + tu); 613e77315bcSPatrick Sanan bu = PetscPowScalar(au, bm1); 614e77315bcSPatrick Sanan /* du = bu * au; */ 615e77315bcSPatrick Sanan du = PetscPowScalar(au, beta); 616e77315bcSPatrick Sanan gu = coef * bu * (tu - t0); 617e77315bcSPatrick Sanan 618e77315bcSPatrick Sanan td = x[k - 1][j][i]; 619e77315bcSPatrick Sanan ad = 0.5 * (t0 + td); 620e77315bcSPatrick Sanan bd = PetscPowScalar(ad, bm1); 621e77315bcSPatrick Sanan /* dd = bd * ad; */ 622e77315bcSPatrick Sanan dd = PetscPowScalar(ad, beta); 623e77315bcSPatrick Sanan gd = coef * bd * (td - t0); 624e77315bcSPatrick Sanan 6259371c9d4SSatish Balay c[0].k = k - 1; 6269371c9d4SSatish Balay c[0].j = j; 6279371c9d4SSatish Balay c[0].i = i; 628e77315bcSPatrick Sanan v[0] = -hxhydhz * (dd - gd); 6299371c9d4SSatish Balay c[1].k = k; 6309371c9d4SSatish Balay c[1].j = j; 6319371c9d4SSatish Balay c[1].i = i; 632e77315bcSPatrick Sanan v[1] = hzhxdhy * (dn - gn) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (dd + du + gd - gu); 6339371c9d4SSatish Balay c[2].k = k; 6349371c9d4SSatish Balay c[2].j = j; 6359371c9d4SSatish Balay c[2].i = i + 1; 636e77315bcSPatrick Sanan v[2] = -hyhzdhx * (de + ge); 6379371c9d4SSatish Balay c[3].k = k; 6389371c9d4SSatish Balay c[3].j = j + 1; 6399371c9d4SSatish Balay c[3].i = i; 640e77315bcSPatrick Sanan v[3] = -hzhxdhy * (dn + gn); 6419371c9d4SSatish Balay c[4].k = k + 1; 6429371c9d4SSatish Balay c[4].j = j; 6439371c9d4SSatish Balay c[4].i = i; 644e77315bcSPatrick Sanan v[4] = -hxhydhz * (du + gu); 6459566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac, 1, &row, 5, c, v, INSERT_VALUES)); 646e77315bcSPatrick Sanan 647e77315bcSPatrick Sanan /* left-hand bottom up corner */ 648e77315bcSPatrick Sanan } else { 649e77315bcSPatrick Sanan td = x[k - 1][j][i]; 650e77315bcSPatrick Sanan ad = 0.5 * (t0 + td); 651e77315bcSPatrick Sanan bd = PetscPowScalar(ad, bm1); 652e77315bcSPatrick Sanan /* dd = bd * ad; */ 653e77315bcSPatrick Sanan dd = PetscPowScalar(ad, beta); 654e77315bcSPatrick Sanan gd = coef * bd * (td - t0); 655e77315bcSPatrick Sanan 6569371c9d4SSatish Balay c[0].k = k - 1; 6579371c9d4SSatish Balay c[0].j = j; 6589371c9d4SSatish Balay c[0].i = i; 659e77315bcSPatrick Sanan v[0] = -hxhydhz * (dd - gd); 6609371c9d4SSatish Balay c[1].k = k; 6619371c9d4SSatish Balay c[1].j = j; 6629371c9d4SSatish Balay c[1].i = i; 663e77315bcSPatrick Sanan v[1] = hzhxdhy * (dn - gn) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (dd + gd); 6649371c9d4SSatish Balay c[2].k = k; 6659371c9d4SSatish Balay c[2].j = j; 6669371c9d4SSatish Balay c[2].i = i + 1; 667e77315bcSPatrick Sanan v[2] = -hyhzdhx * (de + ge); 6689371c9d4SSatish Balay c[3].k = k; 6699371c9d4SSatish Balay c[3].j = j + 1; 6709371c9d4SSatish Balay c[3].i = i; 671e77315bcSPatrick Sanan v[3] = -hzhxdhy * (dn + gn); 6729566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac, 1, &row, 4, c, v, INSERT_VALUES)); 673e77315bcSPatrick Sanan } 674e77315bcSPatrick Sanan 675e77315bcSPatrick Sanan /* left-hand top edge */ 676e77315bcSPatrick Sanan } else if (j == my - 1) { 677e77315bcSPatrick Sanan ts = x[k][j - 1][i]; 678e77315bcSPatrick Sanan as = 0.5 * (t0 + ts); 679e77315bcSPatrick Sanan bs = PetscPowScalar(as, bm1); 680e77315bcSPatrick Sanan /* ds = bs * as; */ 681e77315bcSPatrick Sanan ds = PetscPowScalar(as, beta); 682e77315bcSPatrick Sanan gs = coef * bs * (ts - t0); 683e77315bcSPatrick Sanan 684e77315bcSPatrick Sanan /* left-hand top down corner */ 685e77315bcSPatrick Sanan if (k == 0) { 686e77315bcSPatrick Sanan tu = x[k + 1][j][i]; 687e77315bcSPatrick Sanan au = 0.5 * (t0 + tu); 688e77315bcSPatrick Sanan bu = PetscPowScalar(au, bm1); 689e77315bcSPatrick Sanan /* du = bu * au; */ 690e77315bcSPatrick Sanan du = PetscPowScalar(au, beta); 691e77315bcSPatrick Sanan gu = coef * bu * (tu - t0); 692e77315bcSPatrick Sanan 6939371c9d4SSatish Balay c[0].k = k; 6949371c9d4SSatish Balay c[0].j = j - 1; 6959371c9d4SSatish Balay c[0].i = i; 696e77315bcSPatrick Sanan v[0] = -hzhxdhy * (ds - gs); 6979371c9d4SSatish Balay c[1].k = k; 6989371c9d4SSatish Balay c[1].j = j; 6999371c9d4SSatish Balay c[1].i = i; 700e77315bcSPatrick Sanan v[1] = hzhxdhy * (ds + gs) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (du - gu); 7019371c9d4SSatish Balay c[2].k = k; 7029371c9d4SSatish Balay c[2].j = j; 7039371c9d4SSatish Balay c[2].i = i + 1; 704e77315bcSPatrick Sanan v[2] = -hyhzdhx * (de + ge); 7059371c9d4SSatish Balay c[3].k = k + 1; 7069371c9d4SSatish Balay c[3].j = j; 7079371c9d4SSatish Balay c[3].i = i; 708e77315bcSPatrick Sanan v[3] = -hxhydhz * (du + gu); 7099566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac, 1, &row, 4, c, v, INSERT_VALUES)); 710e77315bcSPatrick Sanan 711e77315bcSPatrick Sanan /* left-hand top interior edge */ 712e77315bcSPatrick Sanan } else if (k < mz - 1) { 713e77315bcSPatrick Sanan tu = x[k + 1][j][i]; 714e77315bcSPatrick Sanan au = 0.5 * (t0 + tu); 715e77315bcSPatrick Sanan bu = PetscPowScalar(au, bm1); 716e77315bcSPatrick Sanan /* du = bu * au; */ 717e77315bcSPatrick Sanan du = PetscPowScalar(au, beta); 718e77315bcSPatrick Sanan gu = coef * bu * (tu - t0); 719e77315bcSPatrick Sanan 720e77315bcSPatrick Sanan td = x[k - 1][j][i]; 721e77315bcSPatrick Sanan ad = 0.5 * (t0 + td); 722e77315bcSPatrick Sanan bd = PetscPowScalar(ad, bm1); 723e77315bcSPatrick Sanan /* dd = bd * ad; */ 724e77315bcSPatrick Sanan dd = PetscPowScalar(ad, beta); 725e77315bcSPatrick Sanan gd = coef * bd * (td - t0); 726e77315bcSPatrick Sanan 7279371c9d4SSatish Balay c[0].k = k - 1; 7289371c9d4SSatish Balay c[0].j = j; 7299371c9d4SSatish Balay c[0].i = i; 730e77315bcSPatrick Sanan v[0] = -hxhydhz * (dd - gd); 7319371c9d4SSatish Balay c[1].k = k; 7329371c9d4SSatish Balay c[1].j = j - 1; 7339371c9d4SSatish Balay c[1].i = i; 734e77315bcSPatrick Sanan v[1] = -hzhxdhy * (ds - gs); 7359371c9d4SSatish Balay c[2].k = k; 7369371c9d4SSatish Balay c[2].j = j; 7379371c9d4SSatish Balay c[2].i = i; 738e77315bcSPatrick Sanan v[2] = hzhxdhy * (ds + gs) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (dd + du + gd - gu); 7399371c9d4SSatish Balay c[3].k = k; 7409371c9d4SSatish Balay c[3].j = j; 7419371c9d4SSatish Balay c[3].i = i + 1; 742e77315bcSPatrick Sanan v[3] = -hyhzdhx * (de + ge); 7439371c9d4SSatish Balay c[4].k = k + 1; 7449371c9d4SSatish Balay c[4].j = j; 7459371c9d4SSatish Balay c[4].i = i; 746e77315bcSPatrick Sanan v[4] = -hxhydhz * (du + gu); 7479566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac, 1, &row, 5, c, v, INSERT_VALUES)); 748e77315bcSPatrick Sanan 749e77315bcSPatrick Sanan /* left-hand top up corner */ 750e77315bcSPatrick Sanan } else { 751e77315bcSPatrick Sanan td = x[k - 1][j][i]; 752e77315bcSPatrick Sanan ad = 0.5 * (t0 + td); 753e77315bcSPatrick Sanan bd = PetscPowScalar(ad, bm1); 754e77315bcSPatrick Sanan /* dd = bd * ad; */ 755e77315bcSPatrick Sanan dd = PetscPowScalar(ad, beta); 756e77315bcSPatrick Sanan gd = coef * bd * (td - t0); 757e77315bcSPatrick Sanan 7589371c9d4SSatish Balay c[0].k = k - 1; 7599371c9d4SSatish Balay c[0].j = j; 7609371c9d4SSatish Balay c[0].i = i; 761e77315bcSPatrick Sanan v[0] = -hxhydhz * (dd - gd); 7629371c9d4SSatish Balay c[1].k = k; 7639371c9d4SSatish Balay c[1].j = j - 1; 7649371c9d4SSatish Balay c[1].i = i; 765e77315bcSPatrick Sanan v[1] = -hzhxdhy * (ds - gs); 7669371c9d4SSatish Balay c[2].k = k; 7679371c9d4SSatish Balay c[2].j = j; 7689371c9d4SSatish Balay c[2].i = i; 769e77315bcSPatrick Sanan v[2] = hzhxdhy * (ds + gs) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (dd + gd); 7709371c9d4SSatish Balay c[3].k = k; 7719371c9d4SSatish Balay c[3].j = j; 7729371c9d4SSatish Balay c[3].i = i + 1; 773e77315bcSPatrick Sanan v[3] = -hyhzdhx * (de + ge); 7749566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac, 1, &row, 4, c, v, INSERT_VALUES)); 775e77315bcSPatrick Sanan } 776e77315bcSPatrick Sanan 777e77315bcSPatrick Sanan } else { 778e77315bcSPatrick Sanan ts = x[k][j - 1][i]; 779e77315bcSPatrick Sanan as = 0.5 * (t0 + ts); 780e77315bcSPatrick Sanan bs = PetscPowScalar(as, bm1); 781e77315bcSPatrick Sanan /* ds = bs * as; */ 782e77315bcSPatrick Sanan ds = PetscPowScalar(as, beta); 783e77315bcSPatrick Sanan gs = coef * bs * (t0 - ts); 784e77315bcSPatrick Sanan 785e77315bcSPatrick Sanan tn = x[k][j + 1][i]; 786e77315bcSPatrick Sanan an = 0.5 * (t0 + tn); 787e77315bcSPatrick Sanan bn = PetscPowScalar(an, bm1); 788e77315bcSPatrick Sanan /* dn = bn * an; */ 789e77315bcSPatrick Sanan dn = PetscPowScalar(an, beta); 790e77315bcSPatrick Sanan gn = coef * bn * (tn - t0); 791e77315bcSPatrick Sanan 792e77315bcSPatrick Sanan /* left-hand down interior edge */ 793e77315bcSPatrick Sanan if (k == 0) { 794e77315bcSPatrick Sanan tu = x[k + 1][j][i]; 795e77315bcSPatrick Sanan au = 0.5 * (t0 + tu); 796e77315bcSPatrick Sanan bu = PetscPowScalar(au, bm1); 797e77315bcSPatrick Sanan /* du = bu * au; */ 798e77315bcSPatrick Sanan du = PetscPowScalar(au, beta); 799e77315bcSPatrick Sanan gu = coef * bu * (tu - t0); 800e77315bcSPatrick Sanan 8019371c9d4SSatish Balay c[0].k = k; 8029371c9d4SSatish Balay c[0].j = j - 1; 8039371c9d4SSatish Balay c[0].i = i; 804e77315bcSPatrick Sanan v[0] = -hzhxdhy * (ds - gs); 8059371c9d4SSatish Balay c[1].k = k; 8069371c9d4SSatish Balay c[1].j = j; 8079371c9d4SSatish Balay c[1].i = i; 808e77315bcSPatrick Sanan v[1] = hzhxdhy * (ds + dn + gs - gn) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (du - gu); 8099371c9d4SSatish Balay c[2].k = k; 8109371c9d4SSatish Balay c[2].j = j; 8119371c9d4SSatish Balay c[2].i = i + 1; 812e77315bcSPatrick Sanan v[2] = -hyhzdhx * (de + ge); 8139371c9d4SSatish Balay c[3].k = k; 8149371c9d4SSatish Balay c[3].j = j + 1; 8159371c9d4SSatish Balay c[3].i = i; 816e77315bcSPatrick Sanan v[3] = -hzhxdhy * (dn + gn); 8179371c9d4SSatish Balay c[4].k = k + 1; 8189371c9d4SSatish Balay c[4].j = j; 8199371c9d4SSatish Balay c[4].i = i; 820e77315bcSPatrick Sanan v[4] = -hxhydhz * (du + gu); 8219566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac, 1, &row, 5, c, v, INSERT_VALUES)); 822e77315bcSPatrick Sanan 823e77315bcSPatrick Sanan } else if (k == mz - 1) { /* left-hand up interior edge */ 824e77315bcSPatrick Sanan 825e77315bcSPatrick Sanan td = x[k - 1][j][i]; 826e77315bcSPatrick Sanan ad = 0.5 * (t0 + td); 827e77315bcSPatrick Sanan bd = PetscPowScalar(ad, bm1); 828e77315bcSPatrick Sanan /* dd = bd * ad; */ 829e77315bcSPatrick Sanan dd = PetscPowScalar(ad, beta); 830e77315bcSPatrick Sanan gd = coef * bd * (t0 - td); 831e77315bcSPatrick Sanan 8329371c9d4SSatish Balay c[0].k = k - 1; 8339371c9d4SSatish Balay c[0].j = j; 8349371c9d4SSatish Balay c[0].i = i; 835e77315bcSPatrick Sanan v[0] = -hxhydhz * (dd - gd); 8369371c9d4SSatish Balay c[1].k = k; 8379371c9d4SSatish Balay c[1].j = j - 1; 8389371c9d4SSatish Balay c[1].i = i; 839e77315bcSPatrick Sanan v[1] = -hzhxdhy * (ds - gs); 8409371c9d4SSatish Balay c[2].k = k; 8419371c9d4SSatish Balay c[2].j = j; 8429371c9d4SSatish Balay c[2].i = i; 843e77315bcSPatrick Sanan v[2] = hzhxdhy * (ds + dn + gs - gn) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (dd + gd); 8449371c9d4SSatish Balay c[3].k = k; 8459371c9d4SSatish Balay c[3].j = j; 8469371c9d4SSatish Balay c[3].i = i + 1; 847e77315bcSPatrick Sanan v[3] = -hyhzdhx * (de + ge); 8489371c9d4SSatish Balay c[4].k = k; 8499371c9d4SSatish Balay c[4].j = j + 1; 8509371c9d4SSatish Balay c[4].i = i; 851e77315bcSPatrick Sanan v[4] = -hzhxdhy * (dn + gn); 8529566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac, 1, &row, 5, c, v, INSERT_VALUES)); 853e77315bcSPatrick Sanan } else { /* left-hand interior plane */ 854e77315bcSPatrick Sanan 855e77315bcSPatrick Sanan td = x[k - 1][j][i]; 856e77315bcSPatrick Sanan ad = 0.5 * (t0 + td); 857e77315bcSPatrick Sanan bd = PetscPowScalar(ad, bm1); 858e77315bcSPatrick Sanan /* dd = bd * ad; */ 859e77315bcSPatrick Sanan dd = PetscPowScalar(ad, beta); 860e77315bcSPatrick Sanan gd = coef * bd * (t0 - td); 861e77315bcSPatrick Sanan 862e77315bcSPatrick Sanan tu = x[k + 1][j][i]; 863e77315bcSPatrick Sanan au = 0.5 * (t0 + tu); 864e77315bcSPatrick Sanan bu = PetscPowScalar(au, bm1); 865e77315bcSPatrick Sanan /* du = bu * au; */ 866e77315bcSPatrick Sanan du = PetscPowScalar(au, beta); 867e77315bcSPatrick Sanan gu = coef * bu * (tu - t0); 868e77315bcSPatrick Sanan 8699371c9d4SSatish Balay c[0].k = k - 1; 8709371c9d4SSatish Balay c[0].j = j; 8719371c9d4SSatish Balay c[0].i = i; 872e77315bcSPatrick Sanan v[0] = -hxhydhz * (dd - gd); 8739371c9d4SSatish Balay c[1].k = k; 8749371c9d4SSatish Balay c[1].j = j - 1; 8759371c9d4SSatish Balay c[1].i = i; 876e77315bcSPatrick Sanan v[1] = -hzhxdhy * (ds - gs); 8779371c9d4SSatish Balay c[2].k = k; 8789371c9d4SSatish Balay c[2].j = j; 8799371c9d4SSatish Balay c[2].i = i; 880e77315bcSPatrick Sanan v[2] = hzhxdhy * (ds + dn + gs - gn) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (dd + du + gd - gu); 8819371c9d4SSatish Balay c[3].k = k; 8829371c9d4SSatish Balay c[3].j = j; 8839371c9d4SSatish Balay c[3].i = i + 1; 884e77315bcSPatrick Sanan v[3] = -hyhzdhx * (de + ge); 8859371c9d4SSatish Balay c[4].k = k; 8869371c9d4SSatish Balay c[4].j = j + 1; 8879371c9d4SSatish Balay c[4].i = i; 888e77315bcSPatrick Sanan v[4] = -hzhxdhy * (dn + gn); 8899371c9d4SSatish Balay c[5].k = k + 1; 8909371c9d4SSatish Balay c[5].j = j; 8919371c9d4SSatish Balay c[5].i = i; 892e77315bcSPatrick Sanan v[5] = -hxhydhz * (du + gu); 8939566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac, 1, &row, 6, c, v, INSERT_VALUES)); 894e77315bcSPatrick Sanan } 895e77315bcSPatrick Sanan } 896e77315bcSPatrick Sanan 897e77315bcSPatrick Sanan } else if (i == mx - 1) { 898e77315bcSPatrick Sanan /* right-hand plane boundary */ 899e77315bcSPatrick Sanan tw = x[k][j][i - 1]; 900e77315bcSPatrick Sanan aw = 0.5 * (t0 + tw); 901e77315bcSPatrick Sanan bw = PetscPowScalar(aw, bm1); 902e77315bcSPatrick Sanan /* dw = bw * aw */ 903e77315bcSPatrick Sanan dw = PetscPowScalar(aw, beta); 904e77315bcSPatrick Sanan gw = coef * bw * (t0 - tw); 905e77315bcSPatrick Sanan 906e77315bcSPatrick Sanan te = tright; 907e77315bcSPatrick Sanan ae = 0.5 * (t0 + te); 908e77315bcSPatrick Sanan be = PetscPowScalar(ae, bm1); 909e77315bcSPatrick Sanan /* de = be * ae; */ 910e77315bcSPatrick Sanan de = PetscPowScalar(ae, beta); 911e77315bcSPatrick Sanan ge = coef * be * (te - t0); 912e77315bcSPatrick Sanan 913e77315bcSPatrick Sanan /* right-hand bottom edge */ 914e77315bcSPatrick Sanan if (j == 0) { 915e77315bcSPatrick Sanan tn = x[k][j + 1][i]; 916e77315bcSPatrick Sanan an = 0.5 * (t0 + tn); 917e77315bcSPatrick Sanan bn = PetscPowScalar(an, bm1); 918e77315bcSPatrick Sanan /* dn = bn * an; */ 919e77315bcSPatrick Sanan dn = PetscPowScalar(an, beta); 920e77315bcSPatrick Sanan gn = coef * bn * (tn - t0); 921e77315bcSPatrick Sanan 922e77315bcSPatrick Sanan /* right-hand bottom down corner */ 923e77315bcSPatrick Sanan if (k == 0) { 924e77315bcSPatrick Sanan tu = x[k + 1][j][i]; 925e77315bcSPatrick Sanan au = 0.5 * (t0 + tu); 926e77315bcSPatrick Sanan bu = PetscPowScalar(au, bm1); 927e77315bcSPatrick Sanan /* du = bu * au; */ 928e77315bcSPatrick Sanan du = PetscPowScalar(au, beta); 929e77315bcSPatrick Sanan gu = coef * bu * (tu - t0); 930e77315bcSPatrick Sanan 9319371c9d4SSatish Balay c[0].k = k; 9329371c9d4SSatish Balay c[0].j = j; 9339371c9d4SSatish Balay c[0].i = i - 1; 934e77315bcSPatrick Sanan v[0] = -hyhzdhx * (dw - gw); 9359371c9d4SSatish Balay c[1].k = k; 9369371c9d4SSatish Balay c[1].j = j; 9379371c9d4SSatish Balay c[1].i = i; 938e77315bcSPatrick Sanan v[1] = hzhxdhy * (dn - gn) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (du - gu); 9399371c9d4SSatish Balay c[2].k = k; 9409371c9d4SSatish Balay c[2].j = j + 1; 9419371c9d4SSatish Balay c[2].i = i; 942e77315bcSPatrick Sanan v[2] = -hzhxdhy * (dn + gn); 9439371c9d4SSatish Balay c[3].k = k + 1; 9449371c9d4SSatish Balay c[3].j = j; 9459371c9d4SSatish Balay c[3].i = i; 946e77315bcSPatrick Sanan v[3] = -hxhydhz * (du + gu); 9479566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac, 1, &row, 4, c, v, INSERT_VALUES)); 948e77315bcSPatrick Sanan 949e77315bcSPatrick Sanan /* right-hand bottom interior edge */ 950e77315bcSPatrick Sanan } else if (k < mz - 1) { 951e77315bcSPatrick Sanan tu = x[k + 1][j][i]; 952e77315bcSPatrick Sanan au = 0.5 * (t0 + tu); 953e77315bcSPatrick Sanan bu = PetscPowScalar(au, bm1); 954e77315bcSPatrick Sanan /* du = bu * au; */ 955e77315bcSPatrick Sanan du = PetscPowScalar(au, beta); 956e77315bcSPatrick Sanan gu = coef * bu * (tu - t0); 957e77315bcSPatrick Sanan 958e77315bcSPatrick Sanan td = x[k - 1][j][i]; 959e77315bcSPatrick Sanan ad = 0.5 * (t0 + td); 960e77315bcSPatrick Sanan bd = PetscPowScalar(ad, bm1); 961e77315bcSPatrick Sanan /* dd = bd * ad; */ 962e77315bcSPatrick Sanan dd = PetscPowScalar(ad, beta); 963e77315bcSPatrick Sanan gd = coef * bd * (td - t0); 964e77315bcSPatrick Sanan 9659371c9d4SSatish Balay c[0].k = k - 1; 9669371c9d4SSatish Balay c[0].j = j; 9679371c9d4SSatish Balay c[0].i = i; 968e77315bcSPatrick Sanan v[0] = -hxhydhz * (dd - gd); 9699371c9d4SSatish Balay c[1].k = k; 9709371c9d4SSatish Balay c[1].j = j; 9719371c9d4SSatish Balay c[1].i = i - 1; 972e77315bcSPatrick Sanan v[1] = -hyhzdhx * (dw - gw); 9739371c9d4SSatish Balay c[2].k = k; 9749371c9d4SSatish Balay c[2].j = j; 9759371c9d4SSatish Balay c[2].i = i; 976e77315bcSPatrick Sanan v[2] = hzhxdhy * (dn - gn) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (dd + du + gd - gu); 9779371c9d4SSatish Balay c[3].k = k; 9789371c9d4SSatish Balay c[3].j = j + 1; 9799371c9d4SSatish Balay c[3].i = i; 980e77315bcSPatrick Sanan v[3] = -hzhxdhy * (dn + gn); 9819371c9d4SSatish Balay c[4].k = k + 1; 9829371c9d4SSatish Balay c[4].j = j; 9839371c9d4SSatish Balay c[4].i = i; 984e77315bcSPatrick Sanan v[4] = -hxhydhz * (du + gu); 9859566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac, 1, &row, 5, c, v, INSERT_VALUES)); 986e77315bcSPatrick Sanan 987e77315bcSPatrick Sanan /* right-hand bottom up corner */ 988e77315bcSPatrick Sanan } else { 989e77315bcSPatrick Sanan td = x[k - 1][j][i]; 990e77315bcSPatrick Sanan ad = 0.5 * (t0 + td); 991e77315bcSPatrick Sanan bd = PetscPowScalar(ad, bm1); 992e77315bcSPatrick Sanan /* dd = bd * ad; */ 993e77315bcSPatrick Sanan dd = PetscPowScalar(ad, beta); 994e77315bcSPatrick Sanan gd = coef * bd * (td - t0); 995e77315bcSPatrick Sanan 9969371c9d4SSatish Balay c[0].k = k - 1; 9979371c9d4SSatish Balay c[0].j = j; 9989371c9d4SSatish Balay c[0].i = i; 999e77315bcSPatrick Sanan v[0] = -hxhydhz * (dd - gd); 10009371c9d4SSatish Balay c[1].k = k; 10019371c9d4SSatish Balay c[1].j = j; 10029371c9d4SSatish Balay c[1].i = i - 1; 1003e77315bcSPatrick Sanan v[1] = -hyhzdhx * (dw - gw); 10049371c9d4SSatish Balay c[2].k = k; 10059371c9d4SSatish Balay c[2].j = j; 10069371c9d4SSatish Balay c[2].i = i; 1007e77315bcSPatrick Sanan v[2] = hzhxdhy * (dn - gn) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (dd + gd); 10089371c9d4SSatish Balay c[3].k = k; 10099371c9d4SSatish Balay c[3].j = j + 1; 10109371c9d4SSatish Balay c[3].i = i; 1011e77315bcSPatrick Sanan v[3] = -hzhxdhy * (dn + gn); 10129566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac, 1, &row, 4, c, v, INSERT_VALUES)); 1013e77315bcSPatrick Sanan } 1014e77315bcSPatrick Sanan 1015e77315bcSPatrick Sanan /* right-hand top edge */ 1016e77315bcSPatrick Sanan } else if (j == my - 1) { 1017e77315bcSPatrick Sanan ts = x[k][j - 1][i]; 1018e77315bcSPatrick Sanan as = 0.5 * (t0 + ts); 1019e77315bcSPatrick Sanan bs = PetscPowScalar(as, bm1); 1020e77315bcSPatrick Sanan /* ds = bs * as; */ 1021e77315bcSPatrick Sanan ds = PetscPowScalar(as, beta); 1022e77315bcSPatrick Sanan gs = coef * bs * (ts - t0); 1023e77315bcSPatrick Sanan 1024e77315bcSPatrick Sanan /* right-hand top down corner */ 1025e77315bcSPatrick Sanan if (k == 0) { 1026e77315bcSPatrick Sanan tu = x[k + 1][j][i]; 1027e77315bcSPatrick Sanan au = 0.5 * (t0 + tu); 1028e77315bcSPatrick Sanan bu = PetscPowScalar(au, bm1); 1029e77315bcSPatrick Sanan /* du = bu * au; */ 1030e77315bcSPatrick Sanan du = PetscPowScalar(au, beta); 1031e77315bcSPatrick Sanan gu = coef * bu * (tu - t0); 1032e77315bcSPatrick Sanan 10339371c9d4SSatish Balay c[0].k = k; 10349371c9d4SSatish Balay c[0].j = j - 1; 10359371c9d4SSatish Balay c[0].i = i; 1036e77315bcSPatrick Sanan v[0] = -hzhxdhy * (ds - gs); 10379371c9d4SSatish Balay c[1].k = k; 10389371c9d4SSatish Balay c[1].j = j; 10399371c9d4SSatish Balay c[1].i = i - 1; 1040e77315bcSPatrick Sanan v[1] = -hyhzdhx * (dw - gw); 10419371c9d4SSatish Balay c[2].k = k; 10429371c9d4SSatish Balay c[2].j = j; 10439371c9d4SSatish Balay c[2].i = i; 1044e77315bcSPatrick Sanan v[2] = hzhxdhy * (ds + gs) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (du - gu); 10459371c9d4SSatish Balay c[3].k = k + 1; 10469371c9d4SSatish Balay c[3].j = j; 10479371c9d4SSatish Balay c[3].i = i; 1048e77315bcSPatrick Sanan v[3] = -hxhydhz * (du + gu); 10499566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac, 1, &row, 4, c, v, INSERT_VALUES)); 1050e77315bcSPatrick Sanan 1051e77315bcSPatrick Sanan /* right-hand top interior edge */ 1052e77315bcSPatrick Sanan } else if (k < mz - 1) { 1053e77315bcSPatrick Sanan tu = x[k + 1][j][i]; 1054e77315bcSPatrick Sanan au = 0.5 * (t0 + tu); 1055e77315bcSPatrick Sanan bu = PetscPowScalar(au, bm1); 1056e77315bcSPatrick Sanan /* du = bu * au; */ 1057e77315bcSPatrick Sanan du = PetscPowScalar(au, beta); 1058e77315bcSPatrick Sanan gu = coef * bu * (tu - t0); 1059e77315bcSPatrick Sanan 1060e77315bcSPatrick Sanan td = x[k - 1][j][i]; 1061e77315bcSPatrick Sanan ad = 0.5 * (t0 + td); 1062e77315bcSPatrick Sanan bd = PetscPowScalar(ad, bm1); 1063e77315bcSPatrick Sanan /* dd = bd * ad; */ 1064e77315bcSPatrick Sanan dd = PetscPowScalar(ad, beta); 1065e77315bcSPatrick Sanan gd = coef * bd * (td - t0); 1066e77315bcSPatrick Sanan 10679371c9d4SSatish Balay c[0].k = k - 1; 10689371c9d4SSatish Balay c[0].j = j; 10699371c9d4SSatish Balay c[0].i = i; 1070e77315bcSPatrick Sanan v[0] = -hxhydhz * (dd - gd); 10719371c9d4SSatish Balay c[1].k = k; 10729371c9d4SSatish Balay c[1].j = j - 1; 10739371c9d4SSatish Balay c[1].i = i; 1074e77315bcSPatrick Sanan v[1] = -hzhxdhy * (ds - gs); 10759371c9d4SSatish Balay c[2].k = k; 10769371c9d4SSatish Balay c[2].j = j; 10779371c9d4SSatish Balay c[2].i = i - 1; 1078e77315bcSPatrick Sanan v[2] = -hyhzdhx * (dw - gw); 10799371c9d4SSatish Balay c[3].k = k; 10809371c9d4SSatish Balay c[3].j = j; 10819371c9d4SSatish Balay c[3].i = i; 1082e77315bcSPatrick Sanan v[3] = hzhxdhy * (ds + gs) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (dd + du + gd - gu); 10839371c9d4SSatish Balay c[4].k = k + 1; 10849371c9d4SSatish Balay c[4].j = j; 10859371c9d4SSatish Balay c[4].i = i; 1086e77315bcSPatrick Sanan v[4] = -hxhydhz * (du + gu); 10879566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac, 1, &row, 5, c, v, INSERT_VALUES)); 1088e77315bcSPatrick Sanan 1089e77315bcSPatrick Sanan /* right-hand top up corner */ 1090e77315bcSPatrick Sanan } else { 1091e77315bcSPatrick Sanan td = x[k - 1][j][i]; 1092e77315bcSPatrick Sanan ad = 0.5 * (t0 + td); 1093e77315bcSPatrick Sanan bd = PetscPowScalar(ad, bm1); 1094e77315bcSPatrick Sanan /* dd = bd * ad; */ 1095e77315bcSPatrick Sanan dd = PetscPowScalar(ad, beta); 1096e77315bcSPatrick Sanan gd = coef * bd * (td - t0); 1097e77315bcSPatrick Sanan 10989371c9d4SSatish Balay c[0].k = k - 1; 10999371c9d4SSatish Balay c[0].j = j; 11009371c9d4SSatish Balay c[0].i = i; 1101e77315bcSPatrick Sanan v[0] = -hxhydhz * (dd - gd); 11029371c9d4SSatish Balay c[1].k = k; 11039371c9d4SSatish Balay c[1].j = j - 1; 11049371c9d4SSatish Balay c[1].i = i; 1105e77315bcSPatrick Sanan v[1] = -hzhxdhy * (ds - gs); 11069371c9d4SSatish Balay c[2].k = k; 11079371c9d4SSatish Balay c[2].j = j; 11089371c9d4SSatish Balay c[2].i = i - 1; 1109e77315bcSPatrick Sanan v[2] = -hyhzdhx * (dw - gw); 11109371c9d4SSatish Balay c[3].k = k; 11119371c9d4SSatish Balay c[3].j = j; 11129371c9d4SSatish Balay c[3].i = i; 1113e77315bcSPatrick Sanan v[3] = hzhxdhy * (ds + gs) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (dd + gd); 11149566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac, 1, &row, 4, c, v, INSERT_VALUES)); 1115e77315bcSPatrick Sanan } 1116e77315bcSPatrick Sanan 1117e77315bcSPatrick Sanan } else { 1118e77315bcSPatrick Sanan ts = x[k][j - 1][i]; 1119e77315bcSPatrick Sanan as = 0.5 * (t0 + ts); 1120e77315bcSPatrick Sanan bs = PetscPowScalar(as, bm1); 1121e77315bcSPatrick Sanan /* ds = bs * as; */ 1122e77315bcSPatrick Sanan ds = PetscPowScalar(as, beta); 1123e77315bcSPatrick Sanan gs = coef * bs * (t0 - ts); 1124e77315bcSPatrick Sanan 1125e77315bcSPatrick Sanan tn = x[k][j + 1][i]; 1126e77315bcSPatrick Sanan an = 0.5 * (t0 + tn); 1127e77315bcSPatrick Sanan bn = PetscPowScalar(an, bm1); 1128e77315bcSPatrick Sanan /* dn = bn * an; */ 1129e77315bcSPatrick Sanan dn = PetscPowScalar(an, beta); 1130e77315bcSPatrick Sanan gn = coef * bn * (tn - t0); 1131e77315bcSPatrick Sanan 1132e77315bcSPatrick Sanan /* right-hand down interior edge */ 1133e77315bcSPatrick Sanan if (k == 0) { 1134e77315bcSPatrick Sanan tu = x[k + 1][j][i]; 1135e77315bcSPatrick Sanan au = 0.5 * (t0 + tu); 1136e77315bcSPatrick Sanan bu = PetscPowScalar(au, bm1); 1137e77315bcSPatrick Sanan /* du = bu * au; */ 1138e77315bcSPatrick Sanan du = PetscPowScalar(au, beta); 1139e77315bcSPatrick Sanan gu = coef * bu * (tu - t0); 1140e77315bcSPatrick Sanan 11419371c9d4SSatish Balay c[0].k = k; 11429371c9d4SSatish Balay c[0].j = j - 1; 11439371c9d4SSatish Balay c[0].i = i; 1144e77315bcSPatrick Sanan v[0] = -hzhxdhy * (ds - gs); 11459371c9d4SSatish Balay c[1].k = k; 11469371c9d4SSatish Balay c[1].j = j; 11479371c9d4SSatish Balay c[1].i = i - 1; 1148e77315bcSPatrick Sanan v[1] = -hyhzdhx * (dw - gw); 11499371c9d4SSatish Balay c[2].k = k; 11509371c9d4SSatish Balay c[2].j = j; 11519371c9d4SSatish Balay c[2].i = i; 1152e77315bcSPatrick Sanan v[2] = hzhxdhy * (ds + dn + gs - gn) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (du - gu); 11539371c9d4SSatish Balay c[3].k = k; 11549371c9d4SSatish Balay c[3].j = j + 1; 11559371c9d4SSatish Balay c[3].i = i; 1156e77315bcSPatrick Sanan v[3] = -hzhxdhy * (dn + gn); 11579371c9d4SSatish Balay c[4].k = k + 1; 11589371c9d4SSatish Balay c[4].j = j; 11599371c9d4SSatish Balay c[4].i = i; 1160e77315bcSPatrick Sanan v[4] = -hxhydhz * (du + gu); 11619566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac, 1, &row, 5, c, v, INSERT_VALUES)); 1162e77315bcSPatrick Sanan 1163e77315bcSPatrick Sanan } else if (k == mz - 1) { /* right-hand up interior edge */ 1164e77315bcSPatrick Sanan 1165e77315bcSPatrick Sanan td = x[k - 1][j][i]; 1166e77315bcSPatrick Sanan ad = 0.5 * (t0 + td); 1167e77315bcSPatrick Sanan bd = PetscPowScalar(ad, bm1); 1168e77315bcSPatrick Sanan /* dd = bd * ad; */ 1169e77315bcSPatrick Sanan dd = PetscPowScalar(ad, beta); 1170e77315bcSPatrick Sanan gd = coef * bd * (t0 - td); 1171e77315bcSPatrick Sanan 11729371c9d4SSatish Balay c[0].k = k - 1; 11739371c9d4SSatish Balay c[0].j = j; 11749371c9d4SSatish Balay c[0].i = i; 1175e77315bcSPatrick Sanan v[0] = -hxhydhz * (dd - gd); 11769371c9d4SSatish Balay c[1].k = k; 11779371c9d4SSatish Balay c[1].j = j - 1; 11789371c9d4SSatish Balay c[1].i = i; 1179e77315bcSPatrick Sanan v[1] = -hzhxdhy * (ds - gs); 11809371c9d4SSatish Balay c[2].k = k; 11819371c9d4SSatish Balay c[2].j = j; 11829371c9d4SSatish Balay c[2].i = i - 1; 1183e77315bcSPatrick Sanan v[2] = -hyhzdhx * (dw - gw); 11849371c9d4SSatish Balay c[3].k = k; 11859371c9d4SSatish Balay c[3].j = j; 11869371c9d4SSatish Balay c[3].i = i; 1187e77315bcSPatrick Sanan v[3] = hzhxdhy * (ds + dn + gs - gn) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (dd + gd); 11889371c9d4SSatish Balay c[4].k = k; 11899371c9d4SSatish Balay c[4].j = j + 1; 11909371c9d4SSatish Balay c[4].i = i; 1191e77315bcSPatrick Sanan v[4] = -hzhxdhy * (dn + gn); 11929566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac, 1, &row, 5, c, v, INSERT_VALUES)); 1193e77315bcSPatrick Sanan 1194e77315bcSPatrick Sanan } else { /* right-hand interior plane */ 1195e77315bcSPatrick Sanan 1196e77315bcSPatrick Sanan td = x[k - 1][j][i]; 1197e77315bcSPatrick Sanan ad = 0.5 * (t0 + td); 1198e77315bcSPatrick Sanan bd = PetscPowScalar(ad, bm1); 1199e77315bcSPatrick Sanan /* dd = bd * ad; */ 1200e77315bcSPatrick Sanan dd = PetscPowScalar(ad, beta); 1201e77315bcSPatrick Sanan gd = coef * bd * (t0 - td); 1202e77315bcSPatrick Sanan 1203e77315bcSPatrick Sanan tu = x[k + 1][j][i]; 1204e77315bcSPatrick Sanan au = 0.5 * (t0 + tu); 1205e77315bcSPatrick Sanan bu = PetscPowScalar(au, bm1); 1206e77315bcSPatrick Sanan /* du = bu * au; */ 1207e77315bcSPatrick Sanan du = PetscPowScalar(au, beta); 1208e77315bcSPatrick Sanan gu = coef * bu * (tu - t0); 1209e77315bcSPatrick Sanan 12109371c9d4SSatish Balay c[0].k = k - 1; 12119371c9d4SSatish Balay c[0].j = j; 12129371c9d4SSatish Balay c[0].i = i; 1213e77315bcSPatrick Sanan v[0] = -hxhydhz * (dd - gd); 12149371c9d4SSatish Balay c[1].k = k; 12159371c9d4SSatish Balay c[1].j = j - 1; 12169371c9d4SSatish Balay c[1].i = i; 1217e77315bcSPatrick Sanan v[1] = -hzhxdhy * (ds - gs); 12189371c9d4SSatish Balay c[2].k = k; 12199371c9d4SSatish Balay c[2].j = j; 12209371c9d4SSatish Balay c[2].i = i - 1; 1221e77315bcSPatrick Sanan v[2] = -hyhzdhx * (dw - gw); 12229371c9d4SSatish Balay c[3].k = k; 12239371c9d4SSatish Balay c[3].j = j; 12249371c9d4SSatish Balay c[3].i = i; 1225e77315bcSPatrick Sanan v[3] = hzhxdhy * (ds + dn + gs - gn) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (dd + du + gd - gu); 12269371c9d4SSatish Balay c[4].k = k; 12279371c9d4SSatish Balay c[4].j = j + 1; 12289371c9d4SSatish Balay c[4].i = i; 1229e77315bcSPatrick Sanan v[4] = -hzhxdhy * (dn + gn); 12309371c9d4SSatish Balay c[5].k = k + 1; 12319371c9d4SSatish Balay c[5].j = j; 12329371c9d4SSatish Balay c[5].i = i; 1233e77315bcSPatrick Sanan v[5] = -hxhydhz * (du + gu); 12349566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac, 1, &row, 6, c, v, INSERT_VALUES)); 1235e77315bcSPatrick Sanan } 1236e77315bcSPatrick Sanan } 1237e77315bcSPatrick Sanan 1238e77315bcSPatrick Sanan } else if (j == 0) { 1239e77315bcSPatrick Sanan tw = x[k][j][i - 1]; 1240e77315bcSPatrick Sanan aw = 0.5 * (t0 + tw); 1241e77315bcSPatrick Sanan bw = PetscPowScalar(aw, bm1); 1242e77315bcSPatrick Sanan /* dw = bw * aw */ 1243e77315bcSPatrick Sanan dw = PetscPowScalar(aw, beta); 1244e77315bcSPatrick Sanan gw = coef * bw * (t0 - tw); 1245e77315bcSPatrick Sanan 1246e77315bcSPatrick Sanan te = x[k][j][i + 1]; 1247e77315bcSPatrick Sanan ae = 0.5 * (t0 + te); 1248e77315bcSPatrick Sanan be = PetscPowScalar(ae, bm1); 1249e77315bcSPatrick Sanan /* de = be * ae; */ 1250e77315bcSPatrick Sanan de = PetscPowScalar(ae, beta); 1251e77315bcSPatrick Sanan ge = coef * be * (te - t0); 1252e77315bcSPatrick Sanan 1253e77315bcSPatrick Sanan tn = x[k][j + 1][i]; 1254e77315bcSPatrick Sanan an = 0.5 * (t0 + tn); 1255e77315bcSPatrick Sanan bn = PetscPowScalar(an, bm1); 1256e77315bcSPatrick Sanan /* dn = bn * an; */ 1257e77315bcSPatrick Sanan dn = PetscPowScalar(an, beta); 1258e77315bcSPatrick Sanan gn = coef * bn * (tn - t0); 1259e77315bcSPatrick Sanan 1260e77315bcSPatrick Sanan /* bottom down interior edge */ 1261e77315bcSPatrick Sanan if (k == 0) { 1262e77315bcSPatrick Sanan tu = x[k + 1][j][i]; 1263e77315bcSPatrick Sanan au = 0.5 * (t0 + tu); 1264e77315bcSPatrick Sanan bu = PetscPowScalar(au, bm1); 1265e77315bcSPatrick Sanan /* du = bu * au; */ 1266e77315bcSPatrick Sanan du = PetscPowScalar(au, beta); 1267e77315bcSPatrick Sanan gu = coef * bu * (tu - t0); 1268e77315bcSPatrick Sanan 12699371c9d4SSatish Balay c[0].k = k; 12709371c9d4SSatish Balay c[0].j = j; 12719371c9d4SSatish Balay c[0].i = i - 1; 1272e77315bcSPatrick Sanan v[0] = -hyhzdhx * (dw - gw); 12739371c9d4SSatish Balay c[1].k = k; 12749371c9d4SSatish Balay c[1].j = j; 12759371c9d4SSatish Balay c[1].i = i; 1276e77315bcSPatrick Sanan v[1] = hzhxdhy * (dn - gn) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (du - gu); 12779371c9d4SSatish Balay c[2].k = k; 12789371c9d4SSatish Balay c[2].j = j; 12799371c9d4SSatish Balay c[2].i = i + 1; 1280e77315bcSPatrick Sanan v[2] = -hyhzdhx * (de + ge); 12819371c9d4SSatish Balay c[3].k = k; 12829371c9d4SSatish Balay c[3].j = j + 1; 12839371c9d4SSatish Balay c[3].i = i; 1284e77315bcSPatrick Sanan v[3] = -hzhxdhy * (dn + gn); 12859371c9d4SSatish Balay c[4].k = k + 1; 12869371c9d4SSatish Balay c[4].j = j; 12879371c9d4SSatish Balay c[4].i = i; 1288e77315bcSPatrick Sanan v[4] = -hxhydhz * (du + gu); 12899566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac, 1, &row, 5, c, v, INSERT_VALUES)); 1290e77315bcSPatrick Sanan 1291e77315bcSPatrick Sanan } else if (k == mz - 1) { /* bottom up interior edge */ 1292e77315bcSPatrick Sanan 1293e77315bcSPatrick Sanan td = x[k - 1][j][i]; 1294e77315bcSPatrick Sanan ad = 0.5 * (t0 + td); 1295e77315bcSPatrick Sanan bd = PetscPowScalar(ad, bm1); 1296e77315bcSPatrick Sanan /* dd = bd * ad; */ 1297e77315bcSPatrick Sanan dd = PetscPowScalar(ad, beta); 1298e77315bcSPatrick Sanan gd = coef * bd * (td - t0); 1299e77315bcSPatrick Sanan 13009371c9d4SSatish Balay c[0].k = k - 1; 13019371c9d4SSatish Balay c[0].j = j; 13029371c9d4SSatish Balay c[0].i = i; 1303e77315bcSPatrick Sanan v[0] = -hxhydhz * (dd - gd); 13049371c9d4SSatish Balay c[1].k = k; 13059371c9d4SSatish Balay c[1].j = j; 13069371c9d4SSatish Balay c[1].i = i - 1; 1307e77315bcSPatrick Sanan v[1] = -hyhzdhx * (dw - gw); 13089371c9d4SSatish Balay c[2].k = k; 13099371c9d4SSatish Balay c[2].j = j; 13109371c9d4SSatish Balay c[2].i = i; 1311e77315bcSPatrick Sanan v[2] = hzhxdhy * (dn - gn) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (dd + gd); 13129371c9d4SSatish Balay c[3].k = k; 13139371c9d4SSatish Balay c[3].j = j; 13149371c9d4SSatish Balay c[3].i = i + 1; 1315e77315bcSPatrick Sanan v[3] = -hyhzdhx * (de + ge); 13169371c9d4SSatish Balay c[4].k = k; 13179371c9d4SSatish Balay c[4].j = j + 1; 13189371c9d4SSatish Balay c[4].i = i; 1319e77315bcSPatrick Sanan v[4] = -hzhxdhy * (dn + gn); 13209566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac, 1, &row, 5, c, v, INSERT_VALUES)); 1321e77315bcSPatrick Sanan 1322e77315bcSPatrick Sanan } else { /* bottom interior plane */ 1323e77315bcSPatrick Sanan 1324e77315bcSPatrick Sanan tu = x[k + 1][j][i]; 1325e77315bcSPatrick Sanan au = 0.5 * (t0 + tu); 1326e77315bcSPatrick Sanan bu = PetscPowScalar(au, bm1); 1327e77315bcSPatrick Sanan /* du = bu * au; */ 1328e77315bcSPatrick Sanan du = PetscPowScalar(au, beta); 1329e77315bcSPatrick Sanan gu = coef * bu * (tu - t0); 1330e77315bcSPatrick Sanan 1331e77315bcSPatrick Sanan td = x[k - 1][j][i]; 1332e77315bcSPatrick Sanan ad = 0.5 * (t0 + td); 1333e77315bcSPatrick Sanan bd = PetscPowScalar(ad, bm1); 1334e77315bcSPatrick Sanan /* dd = bd * ad; */ 1335e77315bcSPatrick Sanan dd = PetscPowScalar(ad, beta); 1336e77315bcSPatrick Sanan gd = coef * bd * (td - t0); 1337e77315bcSPatrick Sanan 13389371c9d4SSatish Balay c[0].k = k - 1; 13399371c9d4SSatish Balay c[0].j = j; 13409371c9d4SSatish Balay c[0].i = i; 1341e77315bcSPatrick Sanan v[0] = -hxhydhz * (dd - gd); 13429371c9d4SSatish Balay c[1].k = k; 13439371c9d4SSatish Balay c[1].j = j; 13449371c9d4SSatish Balay c[1].i = i - 1; 1345e77315bcSPatrick Sanan v[1] = -hyhzdhx * (dw - gw); 13469371c9d4SSatish Balay c[2].k = k; 13479371c9d4SSatish Balay c[2].j = j; 13489371c9d4SSatish Balay c[2].i = i; 1349e77315bcSPatrick Sanan v[2] = hzhxdhy * (dn - gn) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (dd + du + gd - gu); 13509371c9d4SSatish Balay c[3].k = k; 13519371c9d4SSatish Balay c[3].j = j; 13529371c9d4SSatish Balay c[3].i = i + 1; 1353e77315bcSPatrick Sanan v[3] = -hyhzdhx * (de + ge); 13549371c9d4SSatish Balay c[4].k = k; 13559371c9d4SSatish Balay c[4].j = j + 1; 13569371c9d4SSatish Balay c[4].i = i; 1357e77315bcSPatrick Sanan v[4] = -hzhxdhy * (dn + gn); 13589371c9d4SSatish Balay c[5].k = k + 1; 13599371c9d4SSatish Balay c[5].j = j; 13609371c9d4SSatish Balay c[5].i = i; 1361e77315bcSPatrick Sanan v[5] = -hxhydhz * (du + gu); 13629566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac, 1, &row, 6, c, v, INSERT_VALUES)); 1363e77315bcSPatrick Sanan } 1364e77315bcSPatrick Sanan 1365e77315bcSPatrick Sanan } else if (j == my - 1) { 1366e77315bcSPatrick Sanan tw = x[k][j][i - 1]; 1367e77315bcSPatrick Sanan aw = 0.5 * (t0 + tw); 1368e77315bcSPatrick Sanan bw = PetscPowScalar(aw, bm1); 1369e77315bcSPatrick Sanan /* dw = bw * aw */ 1370e77315bcSPatrick Sanan dw = PetscPowScalar(aw, beta); 1371e77315bcSPatrick Sanan gw = coef * bw * (t0 - tw); 1372e77315bcSPatrick Sanan 1373e77315bcSPatrick Sanan te = x[k][j][i + 1]; 1374e77315bcSPatrick Sanan ae = 0.5 * (t0 + te); 1375e77315bcSPatrick Sanan be = PetscPowScalar(ae, bm1); 1376e77315bcSPatrick Sanan /* de = be * ae; */ 1377e77315bcSPatrick Sanan de = PetscPowScalar(ae, beta); 1378e77315bcSPatrick Sanan ge = coef * be * (te - t0); 1379e77315bcSPatrick Sanan 1380e77315bcSPatrick Sanan ts = x[k][j - 1][i]; 1381e77315bcSPatrick Sanan as = 0.5 * (t0 + ts); 1382e77315bcSPatrick Sanan bs = PetscPowScalar(as, bm1); 1383e77315bcSPatrick Sanan /* ds = bs * as; */ 1384e77315bcSPatrick Sanan ds = PetscPowScalar(as, beta); 1385e77315bcSPatrick Sanan gs = coef * bs * (t0 - ts); 1386e77315bcSPatrick Sanan 1387e77315bcSPatrick Sanan /* top down interior edge */ 1388e77315bcSPatrick Sanan if (k == 0) { 1389e77315bcSPatrick Sanan tu = x[k + 1][j][i]; 1390e77315bcSPatrick Sanan au = 0.5 * (t0 + tu); 1391e77315bcSPatrick Sanan bu = PetscPowScalar(au, bm1); 1392e77315bcSPatrick Sanan /* du = bu * au; */ 1393e77315bcSPatrick Sanan du = PetscPowScalar(au, beta); 1394e77315bcSPatrick Sanan gu = coef * bu * (tu - t0); 1395e77315bcSPatrick Sanan 13969371c9d4SSatish Balay c[0].k = k; 13979371c9d4SSatish Balay c[0].j = j - 1; 13989371c9d4SSatish Balay c[0].i = i; 1399e77315bcSPatrick Sanan v[0] = -hzhxdhy * (ds - gs); 14009371c9d4SSatish Balay c[1].k = k; 14019371c9d4SSatish Balay c[1].j = j; 14029371c9d4SSatish Balay c[1].i = i - 1; 1403e77315bcSPatrick Sanan v[1] = -hyhzdhx * (dw - gw); 14049371c9d4SSatish Balay c[2].k = k; 14059371c9d4SSatish Balay c[2].j = j; 14069371c9d4SSatish Balay c[2].i = i; 1407e77315bcSPatrick Sanan v[2] = hzhxdhy * (ds + gs) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (du - gu); 14089371c9d4SSatish Balay c[3].k = k; 14099371c9d4SSatish Balay c[3].j = j; 14109371c9d4SSatish Balay c[3].i = i + 1; 1411e77315bcSPatrick Sanan v[3] = -hyhzdhx * (de + ge); 14129371c9d4SSatish Balay c[4].k = k + 1; 14139371c9d4SSatish Balay c[4].j = j; 14149371c9d4SSatish Balay c[4].i = i; 1415e77315bcSPatrick Sanan v[4] = -hxhydhz * (du + gu); 14169566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac, 1, &row, 5, c, v, INSERT_VALUES)); 1417e77315bcSPatrick Sanan 1418e77315bcSPatrick Sanan } else if (k == mz - 1) { /* top up interior edge */ 1419e77315bcSPatrick Sanan 1420e77315bcSPatrick Sanan td = x[k - 1][j][i]; 1421e77315bcSPatrick Sanan ad = 0.5 * (t0 + td); 1422e77315bcSPatrick Sanan bd = PetscPowScalar(ad, bm1); 1423e77315bcSPatrick Sanan /* dd = bd * ad; */ 1424e77315bcSPatrick Sanan dd = PetscPowScalar(ad, beta); 1425e77315bcSPatrick Sanan gd = coef * bd * (td - t0); 1426e77315bcSPatrick Sanan 14279371c9d4SSatish Balay c[0].k = k - 1; 14289371c9d4SSatish Balay c[0].j = j; 14299371c9d4SSatish Balay c[0].i = i; 1430e77315bcSPatrick Sanan v[0] = -hxhydhz * (dd - gd); 14319371c9d4SSatish Balay c[1].k = k; 14329371c9d4SSatish Balay c[1].j = j - 1; 14339371c9d4SSatish Balay c[1].i = i; 1434e77315bcSPatrick Sanan v[1] = -hzhxdhy * (ds - gs); 14359371c9d4SSatish Balay c[2].k = k; 14369371c9d4SSatish Balay c[2].j = j; 14379371c9d4SSatish Balay c[2].i = i - 1; 1438e77315bcSPatrick Sanan v[2] = -hyhzdhx * (dw - gw); 14399371c9d4SSatish Balay c[3].k = k; 14409371c9d4SSatish Balay c[3].j = j; 14419371c9d4SSatish Balay c[3].i = i; 1442e77315bcSPatrick Sanan v[3] = hzhxdhy * (ds + gs) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (dd + gd); 14439371c9d4SSatish Balay c[4].k = k; 14449371c9d4SSatish Balay c[4].j = j; 14459371c9d4SSatish Balay c[4].i = i + 1; 1446e77315bcSPatrick Sanan v[4] = -hyhzdhx * (de + ge); 14479566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac, 1, &row, 5, c, v, INSERT_VALUES)); 1448e77315bcSPatrick Sanan 1449e77315bcSPatrick Sanan } else { /* top interior plane */ 1450e77315bcSPatrick Sanan 1451e77315bcSPatrick Sanan tu = x[k + 1][j][i]; 1452e77315bcSPatrick Sanan au = 0.5 * (t0 + tu); 1453e77315bcSPatrick Sanan bu = PetscPowScalar(au, bm1); 1454e77315bcSPatrick Sanan /* du = bu * au; */ 1455e77315bcSPatrick Sanan du = PetscPowScalar(au, beta); 1456e77315bcSPatrick Sanan gu = coef * bu * (tu - t0); 1457e77315bcSPatrick Sanan 1458e77315bcSPatrick Sanan td = x[k - 1][j][i]; 1459e77315bcSPatrick Sanan ad = 0.5 * (t0 + td); 1460e77315bcSPatrick Sanan bd = PetscPowScalar(ad, bm1); 1461e77315bcSPatrick Sanan /* dd = bd * ad; */ 1462e77315bcSPatrick Sanan dd = PetscPowScalar(ad, beta); 1463e77315bcSPatrick Sanan gd = coef * bd * (td - t0); 1464e77315bcSPatrick Sanan 14659371c9d4SSatish Balay c[0].k = k - 1; 14669371c9d4SSatish Balay c[0].j = j; 14679371c9d4SSatish Balay c[0].i = i; 1468e77315bcSPatrick Sanan v[0] = -hxhydhz * (dd - gd); 14699371c9d4SSatish Balay c[1].k = k; 14709371c9d4SSatish Balay c[1].j = j - 1; 14719371c9d4SSatish Balay c[1].i = i; 1472e77315bcSPatrick Sanan v[1] = -hzhxdhy * (ds - gs); 14739371c9d4SSatish Balay c[2].k = k; 14749371c9d4SSatish Balay c[2].j = j; 14759371c9d4SSatish Balay c[2].i = i - 1; 1476e77315bcSPatrick Sanan v[2] = -hyhzdhx * (dw - gw); 14779371c9d4SSatish Balay c[3].k = k; 14789371c9d4SSatish Balay c[3].j = j; 14799371c9d4SSatish Balay c[3].i = i; 1480e77315bcSPatrick Sanan v[3] = hzhxdhy * (ds + gs) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (dd + du + gd - gu); 14819371c9d4SSatish Balay c[4].k = k; 14829371c9d4SSatish Balay c[4].j = j; 14839371c9d4SSatish Balay c[4].i = i + 1; 1484e77315bcSPatrick Sanan v[4] = -hyhzdhx * (de + ge); 14859371c9d4SSatish Balay c[5].k = k + 1; 14869371c9d4SSatish Balay c[5].j = j; 14879371c9d4SSatish Balay c[5].i = i; 1488e77315bcSPatrick Sanan v[5] = -hxhydhz * (du + gu); 14899566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac, 1, &row, 6, c, v, INSERT_VALUES)); 1490e77315bcSPatrick Sanan } 1491e77315bcSPatrick Sanan 1492e77315bcSPatrick Sanan } else if (k == 0) { 1493e77315bcSPatrick Sanan /* down interior plane */ 1494e77315bcSPatrick Sanan 1495e77315bcSPatrick Sanan tw = x[k][j][i - 1]; 1496e77315bcSPatrick Sanan aw = 0.5 * (t0 + tw); 1497e77315bcSPatrick Sanan bw = PetscPowScalar(aw, bm1); 1498e77315bcSPatrick Sanan /* dw = bw * aw */ 1499e77315bcSPatrick Sanan dw = PetscPowScalar(aw, beta); 1500e77315bcSPatrick Sanan gw = coef * bw * (t0 - tw); 1501e77315bcSPatrick Sanan 1502e77315bcSPatrick Sanan te = x[k][j][i + 1]; 1503e77315bcSPatrick Sanan ae = 0.5 * (t0 + te); 1504e77315bcSPatrick Sanan be = PetscPowScalar(ae, bm1); 1505e77315bcSPatrick Sanan /* de = be * ae; */ 1506e77315bcSPatrick Sanan de = PetscPowScalar(ae, beta); 1507e77315bcSPatrick Sanan ge = coef * be * (te - t0); 1508e77315bcSPatrick Sanan 1509e77315bcSPatrick Sanan ts = x[k][j - 1][i]; 1510e77315bcSPatrick Sanan as = 0.5 * (t0 + ts); 1511e77315bcSPatrick Sanan bs = PetscPowScalar(as, bm1); 1512e77315bcSPatrick Sanan /* ds = bs * as; */ 1513e77315bcSPatrick Sanan ds = PetscPowScalar(as, beta); 1514e77315bcSPatrick Sanan gs = coef * bs * (t0 - ts); 1515e77315bcSPatrick Sanan 1516e77315bcSPatrick Sanan tn = x[k][j + 1][i]; 1517e77315bcSPatrick Sanan an = 0.5 * (t0 + tn); 1518e77315bcSPatrick Sanan bn = PetscPowScalar(an, bm1); 1519e77315bcSPatrick Sanan /* dn = bn * an; */ 1520e77315bcSPatrick Sanan dn = PetscPowScalar(an, beta); 1521e77315bcSPatrick Sanan gn = coef * bn * (tn - t0); 1522e77315bcSPatrick Sanan 1523e77315bcSPatrick Sanan tu = x[k + 1][j][i]; 1524e77315bcSPatrick Sanan au = 0.5 * (t0 + tu); 1525e77315bcSPatrick Sanan bu = PetscPowScalar(au, bm1); 1526e77315bcSPatrick Sanan /* du = bu * au; */ 1527e77315bcSPatrick Sanan du = PetscPowScalar(au, beta); 1528e77315bcSPatrick Sanan gu = coef * bu * (tu - t0); 1529e77315bcSPatrick Sanan 15309371c9d4SSatish Balay c[0].k = k; 15319371c9d4SSatish Balay c[0].j = j - 1; 15329371c9d4SSatish Balay c[0].i = i; 1533e77315bcSPatrick Sanan v[0] = -hzhxdhy * (ds - gs); 15349371c9d4SSatish Balay c[1].k = k; 15359371c9d4SSatish Balay c[1].j = j; 15369371c9d4SSatish Balay c[1].i = i - 1; 1537e77315bcSPatrick Sanan v[1] = -hyhzdhx * (dw - gw); 15389371c9d4SSatish Balay c[2].k = k; 15399371c9d4SSatish Balay c[2].j = j; 15409371c9d4SSatish Balay c[2].i = i; 1541e77315bcSPatrick Sanan v[2] = hzhxdhy * (ds + dn + gs - gn) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (du - gu); 15429371c9d4SSatish Balay c[3].k = k; 15439371c9d4SSatish Balay c[3].j = j; 15449371c9d4SSatish Balay c[3].i = i + 1; 1545e77315bcSPatrick Sanan v[3] = -hyhzdhx * (de + ge); 15469371c9d4SSatish Balay c[4].k = k; 15479371c9d4SSatish Balay c[4].j = j + 1; 15489371c9d4SSatish Balay c[4].i = i; 1549e77315bcSPatrick Sanan v[4] = -hzhxdhy * (dn + gn); 15509371c9d4SSatish Balay c[5].k = k + 1; 15519371c9d4SSatish Balay c[5].j = j; 15529371c9d4SSatish Balay c[5].i = i; 1553e77315bcSPatrick Sanan v[5] = -hxhydhz * (du + gu); 15549566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac, 1, &row, 6, c, v, INSERT_VALUES)); 1555e77315bcSPatrick Sanan 1556e77315bcSPatrick Sanan } else if (k == mz - 1) { 1557e77315bcSPatrick Sanan /* up interior plane */ 1558e77315bcSPatrick Sanan 1559e77315bcSPatrick Sanan tw = x[k][j][i - 1]; 1560e77315bcSPatrick Sanan aw = 0.5 * (t0 + tw); 1561e77315bcSPatrick Sanan bw = PetscPowScalar(aw, bm1); 1562e77315bcSPatrick Sanan /* dw = bw * aw */ 1563e77315bcSPatrick Sanan dw = PetscPowScalar(aw, beta); 1564e77315bcSPatrick Sanan gw = coef * bw * (t0 - tw); 1565e77315bcSPatrick Sanan 1566e77315bcSPatrick Sanan te = x[k][j][i + 1]; 1567e77315bcSPatrick Sanan ae = 0.5 * (t0 + te); 1568e77315bcSPatrick Sanan be = PetscPowScalar(ae, bm1); 1569e77315bcSPatrick Sanan /* de = be * ae; */ 1570e77315bcSPatrick Sanan de = PetscPowScalar(ae, beta); 1571e77315bcSPatrick Sanan ge = coef * be * (te - t0); 1572e77315bcSPatrick Sanan 1573e77315bcSPatrick Sanan ts = x[k][j - 1][i]; 1574e77315bcSPatrick Sanan as = 0.5 * (t0 + ts); 1575e77315bcSPatrick Sanan bs = PetscPowScalar(as, bm1); 1576e77315bcSPatrick Sanan /* ds = bs * as; */ 1577e77315bcSPatrick Sanan ds = PetscPowScalar(as, beta); 1578e77315bcSPatrick Sanan gs = coef * bs * (t0 - ts); 1579e77315bcSPatrick Sanan 1580e77315bcSPatrick Sanan tn = x[k][j + 1][i]; 1581e77315bcSPatrick Sanan an = 0.5 * (t0 + tn); 1582e77315bcSPatrick Sanan bn = PetscPowScalar(an, bm1); 1583e77315bcSPatrick Sanan /* dn = bn * an; */ 1584e77315bcSPatrick Sanan dn = PetscPowScalar(an, beta); 1585e77315bcSPatrick Sanan gn = coef * bn * (tn - t0); 1586e77315bcSPatrick Sanan 1587e77315bcSPatrick Sanan td = x[k - 1][j][i]; 1588e77315bcSPatrick Sanan ad = 0.5 * (t0 + td); 1589e77315bcSPatrick Sanan bd = PetscPowScalar(ad, bm1); 1590e77315bcSPatrick Sanan /* dd = bd * ad; */ 1591e77315bcSPatrick Sanan dd = PetscPowScalar(ad, beta); 1592e77315bcSPatrick Sanan gd = coef * bd * (t0 - td); 1593e77315bcSPatrick Sanan 15949371c9d4SSatish Balay c[0].k = k - 1; 15959371c9d4SSatish Balay c[0].j = j; 15969371c9d4SSatish Balay c[0].i = i; 1597e77315bcSPatrick Sanan v[0] = -hxhydhz * (dd - gd); 15989371c9d4SSatish Balay c[1].k = k; 15999371c9d4SSatish Balay c[1].j = j - 1; 16009371c9d4SSatish Balay c[1].i = i; 1601e77315bcSPatrick Sanan v[1] = -hzhxdhy * (ds - gs); 16029371c9d4SSatish Balay c[2].k = k; 16039371c9d4SSatish Balay c[2].j = j; 16049371c9d4SSatish Balay c[2].i = i - 1; 1605e77315bcSPatrick Sanan v[2] = -hyhzdhx * (dw - gw); 16069371c9d4SSatish Balay c[3].k = k; 16079371c9d4SSatish Balay c[3].j = j; 16089371c9d4SSatish Balay c[3].i = i; 1609e77315bcSPatrick Sanan v[3] = hzhxdhy * (ds + dn + gs - gn) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (dd + gd); 16109371c9d4SSatish Balay c[4].k = k; 16119371c9d4SSatish Balay c[4].j = j; 16129371c9d4SSatish Balay c[4].i = i + 1; 1613e77315bcSPatrick Sanan v[4] = -hyhzdhx * (de + ge); 16149371c9d4SSatish Balay c[5].k = k; 16159371c9d4SSatish Balay c[5].j = j + 1; 16169371c9d4SSatish Balay c[5].i = i; 1617e77315bcSPatrick Sanan v[5] = -hzhxdhy * (dn + gn); 16189566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac, 1, &row, 6, c, v, INSERT_VALUES)); 1619e77315bcSPatrick Sanan } 1620e77315bcSPatrick Sanan } 1621e77315bcSPatrick Sanan } 1622e77315bcSPatrick Sanan } 16239566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY)); 16249566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY)); 1625e77315bcSPatrick Sanan if (jac != J) { 16269566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 16279566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 1628e77315bcSPatrick Sanan } 16299566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, localX, &x)); 16309566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &localX)); 1631e77315bcSPatrick Sanan 16329566063dSJacob Faibussowitsch PetscCall(PetscLogFlops((41.0 + 8.0 * POWFLOP) * xm * ym)); 1633*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1634e77315bcSPatrick Sanan } 1635e77315bcSPatrick Sanan 1636e77315bcSPatrick Sanan /*TEST 1637e77315bcSPatrick Sanan 1638e77315bcSPatrick Sanan test: 1639e77315bcSPatrick Sanan nsize: 4 1640e77315bcSPatrick 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 1641e77315bcSPatrick Sanan requires: !single 1642e77315bcSPatrick Sanan 1643e77315bcSPatrick Sanan TEST*/ 1644