1c4762a1bSJed Brown /* 2c4762a1bSJed Brown This example implements the model described in 3c4762a1bSJed Brown 4c4762a1bSJed Brown Rauenzahn, Mousseau, Knoll. "Temporal accuracy of the nonequilibrium radiation diffusion 5c4762a1bSJed Brown equations employing a Saha ionization model" 2005. 6c4762a1bSJed Brown 7c4762a1bSJed Brown The paper discusses three examples, the first two are nondimensional with a simple 8c4762a1bSJed Brown ionization model. The third example is fully dimensional and uses the Saha ionization 9c4762a1bSJed Brown model with realistic parameters. 10c4762a1bSJed Brown */ 11c4762a1bSJed Brown 12c4762a1bSJed Brown #include <petscts.h> 13c4762a1bSJed Brown #include <petscdm.h> 14c4762a1bSJed Brown #include <petscdmda.h> 15c4762a1bSJed Brown 169371c9d4SSatish Balay typedef enum { 179371c9d4SSatish Balay BC_DIRICHLET, 189371c9d4SSatish Balay BC_NEUMANN, 199371c9d4SSatish Balay BC_ROBIN 209371c9d4SSatish Balay } BCType; 21c4762a1bSJed Brown static const char *const BCTypes[] = {"DIRICHLET", "NEUMANN", "ROBIN", "BCType", "BC_", 0}; 229371c9d4SSatish Balay typedef enum { 239371c9d4SSatish Balay JACOBIAN_ANALYTIC, 249371c9d4SSatish Balay JACOBIAN_MATRIXFREE, 259371c9d4SSatish Balay JACOBIAN_FD_COLORING, 269371c9d4SSatish Balay JACOBIAN_FD_FULL 279371c9d4SSatish Balay } JacobianType; 28c4762a1bSJed Brown static const char *const JacobianTypes[] = {"ANALYTIC", "MATRIXFREE", "FD_COLORING", "FD_FULL", "JacobianType", "FD_", 0}; 299371c9d4SSatish Balay typedef enum { 309371c9d4SSatish Balay DISCRETIZATION_FD, 319371c9d4SSatish Balay DISCRETIZATION_FE 329371c9d4SSatish Balay } DiscretizationType; 33c4762a1bSJed Brown static const char *const DiscretizationTypes[] = {"FD", "FE", "DiscretizationType", "DISCRETIZATION_", 0}; 349371c9d4SSatish Balay typedef enum { 359371c9d4SSatish Balay QUADRATURE_GAUSS1, 369371c9d4SSatish Balay QUADRATURE_GAUSS2, 379371c9d4SSatish Balay QUADRATURE_GAUSS3, 389371c9d4SSatish Balay QUADRATURE_GAUSS4, 399371c9d4SSatish Balay QUADRATURE_LOBATTO2, 409371c9d4SSatish Balay QUADRATURE_LOBATTO3 419371c9d4SSatish Balay } QuadratureType; 42c4762a1bSJed Brown static const char *const QuadratureTypes[] = {"GAUSS1", "GAUSS2", "GAUSS3", "GAUSS4", "LOBATTO2", "LOBATTO3", "QuadratureType", "QUADRATURE_", 0}; 43c4762a1bSJed Brown 44c4762a1bSJed Brown typedef struct { 45c4762a1bSJed Brown PetscScalar E; /* radiation energy */ 46c4762a1bSJed Brown PetscScalar T; /* material temperature */ 47c4762a1bSJed Brown } RDNode; 48c4762a1bSJed Brown 49c4762a1bSJed Brown typedef struct { 50c4762a1bSJed Brown PetscReal meter, kilogram, second, Kelvin; /* Fundamental units */ 51c4762a1bSJed Brown PetscReal Joule, Watt; /* Derived units */ 52c4762a1bSJed Brown } RDUnit; 53c4762a1bSJed Brown 54c4762a1bSJed Brown typedef struct _n_RD *RD; 55c4762a1bSJed Brown 56c4762a1bSJed Brown struct _n_RD { 57c4762a1bSJed Brown void (*MaterialEnergy)(RD, const RDNode *, PetscScalar *, RDNode *); 58c4762a1bSJed Brown DM da; 59c4762a1bSJed Brown PetscBool monitor_residual; 60c4762a1bSJed Brown DiscretizationType discretization; 61c4762a1bSJed Brown QuadratureType quadrature; 62c4762a1bSJed Brown JacobianType jacobian; 63c4762a1bSJed Brown PetscInt initial; 64c4762a1bSJed Brown BCType leftbc; 65c4762a1bSJed Brown PetscBool view_draw; 66c4762a1bSJed Brown char view_binary[PETSC_MAX_PATH_LEN]; 67c4762a1bSJed Brown PetscBool test_diff; 68c4762a1bSJed Brown PetscBool endpoint; 69c4762a1bSJed Brown PetscBool bclimit; 70c4762a1bSJed Brown PetscBool bcmidpoint; 71c4762a1bSJed Brown RDUnit unit; 72c4762a1bSJed Brown 73c4762a1bSJed Brown /* model constants, see Table 2 and RDCreate() */ 74c4762a1bSJed Brown PetscReal rho, K_R, K_p, I_H, m_p, m_e, h, k, c, sigma_b, beta, gamma; 75c4762a1bSJed Brown 76c4762a1bSJed Brown /* Domain and boundary conditions */ 77c4762a1bSJed Brown PetscReal Eapplied; /* Radiation flux from the left */ 78c4762a1bSJed Brown PetscReal L; /* Length of domain */ 79c4762a1bSJed Brown PetscReal final_time; 80c4762a1bSJed Brown }; 81c4762a1bSJed Brown 82d71ae5a4SJacob Faibussowitsch static PetscErrorCode RDDestroy(RD *rd) 83d71ae5a4SJacob Faibussowitsch { 84c4762a1bSJed Brown PetscFunctionBeginUser; 859566063dSJacob Faibussowitsch PetscCall(DMDestroy(&(*rd)->da)); 869566063dSJacob Faibussowitsch PetscCall(PetscFree(*rd)); 87*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 88c4762a1bSJed Brown } 89c4762a1bSJed Brown 90c4762a1bSJed Brown /* The paper has a time derivative for material energy (Eq 2) which is a dependent variable (computable from temperature 91c4762a1bSJed Brown * and density through an uninvertible relation). Computing this derivative is trivial for trapezoid rule (used in the 92c4762a1bSJed Brown * paper), but does not generalize nicely to higher order integrators. Here we use the implicit form which provides 93c4762a1bSJed Brown * time derivatives of the independent variables (radiation energy and temperature), so we must compute the time 94c4762a1bSJed Brown * derivative of material energy ourselves (could be done using AD). 95c4762a1bSJed Brown * 96c4762a1bSJed Brown * There are multiple ionization models, this interface dispatches to the one currently in use. 97c4762a1bSJed Brown */ 98d71ae5a4SJacob Faibussowitsch static void RDMaterialEnergy(RD rd, const RDNode *n, PetscScalar *Em, RDNode *dEm) 99d71ae5a4SJacob Faibussowitsch { 1009371c9d4SSatish Balay rd->MaterialEnergy(rd, n, Em, dEm); 1019371c9d4SSatish Balay } 102c4762a1bSJed Brown 103c4762a1bSJed Brown /* Solves a quadratic equation while propagating tangents */ 104d71ae5a4SJacob Faibussowitsch static void QuadraticSolve(PetscScalar a, PetscScalar a_t, PetscScalar b, PetscScalar b_t, PetscScalar c, PetscScalar c_t, PetscScalar *x, PetscScalar *x_t) 105d71ae5a4SJacob Faibussowitsch { 1069371c9d4SSatish Balay PetscScalar disc = b * b - 4. * a * c, disc_t = 2. * b * b_t - 4. * a_t * c - 4. * a * c_t, num = -b + PetscSqrtScalar(disc), /* choose positive sign */ 1079371c9d4SSatish Balay num_t = -b_t + 0.5 / PetscSqrtScalar(disc) * disc_t, den = 2. * a, den_t = 2. * a_t; 108c4762a1bSJed Brown *x = num / den; 109c4762a1bSJed Brown *x_t = (num_t * den - num * den_t) / PetscSqr(den); 110c4762a1bSJed Brown } 111c4762a1bSJed Brown 112c4762a1bSJed Brown /* The primary model presented in the paper */ 113d71ae5a4SJacob Faibussowitsch static void RDMaterialEnergy_Saha(RD rd, const RDNode *n, PetscScalar *inEm, RDNode *dEm) 114d71ae5a4SJacob Faibussowitsch { 1159371c9d4SSatish Balay PetscScalar Em, alpha, alpha_t, T = n->T, T_t = 1., chi = rd->I_H / (rd->k * T), chi_t = -chi / T * T_t, a = 1., a_t = 0, b = 4. * rd->m_p / rd->rho * PetscPowScalarReal(2. * PETSC_PI * rd->m_e * rd->I_H / PetscSqr(rd->h), 1.5) * PetscExpScalar(-chi) * PetscPowScalarReal(chi, 1.5), /* Eq 7 */ 1169371c9d4SSatish Balay b_t = -b * chi_t + 1.5 * b / chi * chi_t, c = -b, c_t = -b_t; 117c4762a1bSJed Brown QuadraticSolve(a, a_t, b, b_t, c, c_t, &alpha, &alpha_t); /* Solve Eq 7 for alpha */ 118c4762a1bSJed Brown Em = rd->k * T / rd->m_p * (1.5 * (1. + alpha) + alpha * chi); /* Eq 6 */ 119c4762a1bSJed Brown if (inEm) *inEm = Em; 120c4762a1bSJed Brown if (dEm) { 121c4762a1bSJed Brown dEm->E = 0; 122c4762a1bSJed Brown dEm->T = Em / T * T_t + rd->k * T / rd->m_p * (1.5 * alpha_t + alpha_t * chi + alpha * chi_t); 123c4762a1bSJed Brown } 124c4762a1bSJed Brown } 125c4762a1bSJed Brown /* Reduced ionization model, Eq 30 */ 126d71ae5a4SJacob Faibussowitsch static void RDMaterialEnergy_Reduced(RD rd, const RDNode *n, PetscScalar *Em, RDNode *dEm) 127d71ae5a4SJacob Faibussowitsch { 1289371c9d4SSatish Balay PetscScalar alpha, alpha_t, T = n->T, T_t = 1., chi = -0.3 / T, chi_t = -chi / T * T_t, a = 1., a_t = 0., b = PetscExpScalar(chi), b_t = b * chi_t, c = -b, c_t = -b_t; 129c4762a1bSJed Brown QuadraticSolve(a, a_t, b, b_t, c, c_t, &alpha, &alpha_t); 130c4762a1bSJed Brown if (Em) *Em = (1. + alpha) * T + 0.3 * alpha; 131c4762a1bSJed Brown if (dEm) { 132c4762a1bSJed Brown dEm->E = 0; 133c4762a1bSJed Brown dEm->T = alpha_t * T + (1. + alpha) * T_t + 0.3 * alpha_t; 134c4762a1bSJed Brown } 135c4762a1bSJed Brown } 136c4762a1bSJed Brown 137c4762a1bSJed Brown /* Eq 5 */ 138d71ae5a4SJacob Faibussowitsch static void RDSigma_R(RD rd, RDNode *n, PetscScalar *sigma_R, RDNode *dsigma_R) 139d71ae5a4SJacob Faibussowitsch { 140c4762a1bSJed Brown *sigma_R = rd->K_R * rd->rho * PetscPowScalar(n->T, -rd->gamma); 141c4762a1bSJed Brown dsigma_R->E = 0; 142c4762a1bSJed Brown dsigma_R->T = -rd->gamma * (*sigma_R) / n->T; 143c4762a1bSJed Brown } 144c4762a1bSJed Brown 145c4762a1bSJed Brown /* Eq 4 */ 146d71ae5a4SJacob Faibussowitsch static void RDDiffusionCoefficient(RD rd, PetscBool limit, RDNode *n, RDNode *nx, PetscScalar *D_R, RDNode *dD_R, RDNode *dxD_R) 147d71ae5a4SJacob Faibussowitsch { 148c4762a1bSJed Brown PetscScalar sigma_R, denom; 149c4762a1bSJed Brown RDNode dsigma_R, ddenom, dxdenom; 150c4762a1bSJed Brown 151c4762a1bSJed Brown RDSigma_R(rd, n, &sigma_R, &dsigma_R); 152c4762a1bSJed Brown denom = 3. * rd->rho * sigma_R + (int)limit * PetscAbsScalar(nx->E) / n->E; 153c4762a1bSJed Brown ddenom.E = -(int)limit * PetscAbsScalar(nx->E) / PetscSqr(n->E); 154c4762a1bSJed Brown ddenom.T = 3. * rd->rho * dsigma_R.T; 155c4762a1bSJed Brown dxdenom.E = (int)limit * (PetscRealPart(nx->E) < 0 ? -1. : 1.) / n->E; 156c4762a1bSJed Brown dxdenom.T = 0; 157c4762a1bSJed Brown *D_R = rd->c / denom; 158c4762a1bSJed Brown if (dD_R) { 159c4762a1bSJed Brown dD_R->E = -rd->c / PetscSqr(denom) * ddenom.E; 160c4762a1bSJed Brown dD_R->T = -rd->c / PetscSqr(denom) * ddenom.T; 161c4762a1bSJed Brown } 162c4762a1bSJed Brown if (dxD_R) { 163c4762a1bSJed Brown dxD_R->E = -rd->c / PetscSqr(denom) * dxdenom.E; 164c4762a1bSJed Brown dxD_R->T = -rd->c / PetscSqr(denom) * dxdenom.T; 165c4762a1bSJed Brown } 166c4762a1bSJed Brown } 167c4762a1bSJed Brown 168d71ae5a4SJacob Faibussowitsch static PetscErrorCode RDStateView(RD rd, Vec X, Vec Xdot, Vec F) 169d71ae5a4SJacob Faibussowitsch { 170c4762a1bSJed Brown DMDALocalInfo info; 171c4762a1bSJed Brown PetscInt i; 172c4762a1bSJed Brown const RDNode *x, *xdot, *f; 173c4762a1bSJed Brown MPI_Comm comm; 174c4762a1bSJed Brown 175c4762a1bSJed Brown PetscFunctionBeginUser; 1769566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)rd->da, &comm)); 1779566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(rd->da, &info)); 1789566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(rd->da, X, (void *)&x)); 1799566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(rd->da, Xdot, (void *)&xdot)); 1809566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(rd->da, F, (void *)&f)); 181c4762a1bSJed Brown for (i = info.xs; i < info.xs + info.xm; i++) { 1829371c9d4SSatish Balay PetscCall(PetscSynchronizedPrintf(comm, "x[%" PetscInt_FMT "] (%10.2G,%10.2G) (%10.2G,%10.2G) (%10.2G,%10.2G)\n", i, (double)PetscRealPart(x[i].E), (double)PetscRealPart(x[i].T), (double)PetscRealPart(xdot[i].E), (double)PetscRealPart(xdot[i].T), 1839371c9d4SSatish Balay (double)PetscRealPart(f[i].E), (double)PetscRealPart(f[i].T))); 184c4762a1bSJed Brown } 1859566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(rd->da, X, (void *)&x)); 1869566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(rd->da, Xdot, (void *)&xdot)); 1879566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(rd->da, F, (void *)&f)); 1889566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(comm, PETSC_STDOUT)); 189*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 190c4762a1bSJed Brown } 191c4762a1bSJed Brown 192d71ae5a4SJacob Faibussowitsch static PetscScalar RDRadiation(RD rd, const RDNode *n, RDNode *dn) 193d71ae5a4SJacob Faibussowitsch { 1949371c9d4SSatish Balay PetscScalar sigma_p = rd->K_p * rd->rho * PetscPowScalar(n->T, -rd->beta), sigma_p_T = -rd->beta * sigma_p / n->T, tmp = 4. * rd->sigma_b * PetscSqr(PetscSqr(n->T)) / rd->c - n->E, tmp_E = -1., tmp_T = 4. * rd->sigma_b * 4 * n->T * (PetscSqr(n->T)) / rd->c, rad = sigma_p * rd->c * rd->rho * tmp, rad_E = sigma_p * rd->c * rd->rho * tmp_E, rad_T = rd->c * rd->rho * (sigma_p_T * tmp + sigma_p * tmp_T); 195c4762a1bSJed Brown if (dn) { 196c4762a1bSJed Brown dn->E = rad_E; 197c4762a1bSJed Brown dn->T = rad_T; 198c4762a1bSJed Brown } 199c4762a1bSJed Brown return rad; 200c4762a1bSJed Brown } 201c4762a1bSJed Brown 202d71ae5a4SJacob Faibussowitsch static PetscScalar RDDiffusion(RD rd, PetscReal hx, const RDNode x[], PetscInt i, RDNode d[]) 203d71ae5a4SJacob Faibussowitsch { 204c4762a1bSJed Brown PetscReal ihx = 1. / hx; 205c4762a1bSJed Brown RDNode n_L, nx_L, n_R, nx_R, dD_L, dxD_L, dD_R, dxD_R, dfluxL[2], dfluxR[2]; 206c4762a1bSJed Brown PetscScalar D_L, D_R, fluxL, fluxR; 207c4762a1bSJed Brown 208c4762a1bSJed Brown n_L.E = 0.5 * (x[i - 1].E + x[i].E); 209c4762a1bSJed Brown n_L.T = 0.5 * (x[i - 1].T + x[i].T); 210c4762a1bSJed Brown nx_L.E = (x[i].E - x[i - 1].E) / hx; 211c4762a1bSJed Brown nx_L.T = (x[i].T - x[i - 1].T) / hx; 212c4762a1bSJed Brown RDDiffusionCoefficient(rd, PETSC_TRUE, &n_L, &nx_L, &D_L, &dD_L, &dxD_L); 213c4762a1bSJed Brown fluxL = D_L * nx_L.E; 214c4762a1bSJed Brown dfluxL[0].E = -ihx * D_L + (0.5 * dD_L.E - ihx * dxD_L.E) * nx_L.E; 215c4762a1bSJed Brown dfluxL[1].E = +ihx * D_L + (0.5 * dD_L.E + ihx * dxD_L.E) * nx_L.E; 216c4762a1bSJed Brown dfluxL[0].T = (0.5 * dD_L.T - ihx * dxD_L.T) * nx_L.E; 217c4762a1bSJed Brown dfluxL[1].T = (0.5 * dD_L.T + ihx * dxD_L.T) * nx_L.E; 218c4762a1bSJed Brown 219c4762a1bSJed Brown n_R.E = 0.5 * (x[i].E + x[i + 1].E); 220c4762a1bSJed Brown n_R.T = 0.5 * (x[i].T + x[i + 1].T); 221c4762a1bSJed Brown nx_R.E = (x[i + 1].E - x[i].E) / hx; 222c4762a1bSJed Brown nx_R.T = (x[i + 1].T - x[i].T) / hx; 223c4762a1bSJed Brown RDDiffusionCoefficient(rd, PETSC_TRUE, &n_R, &nx_R, &D_R, &dD_R, &dxD_R); 224c4762a1bSJed Brown fluxR = D_R * nx_R.E; 225c4762a1bSJed Brown dfluxR[0].E = -ihx * D_R + (0.5 * dD_R.E - ihx * dxD_R.E) * nx_R.E; 226c4762a1bSJed Brown dfluxR[1].E = +ihx * D_R + (0.5 * dD_R.E + ihx * dxD_R.E) * nx_R.E; 227c4762a1bSJed Brown dfluxR[0].T = (0.5 * dD_R.T - ihx * dxD_R.T) * nx_R.E; 228c4762a1bSJed Brown dfluxR[1].T = (0.5 * dD_R.T + ihx * dxD_R.T) * nx_R.E; 229c4762a1bSJed Brown 230c4762a1bSJed Brown if (d) { 231c4762a1bSJed Brown d[0].E = -ihx * dfluxL[0].E; 232c4762a1bSJed Brown d[0].T = -ihx * dfluxL[0].T; 233c4762a1bSJed Brown d[1].E = ihx * (dfluxR[0].E - dfluxL[1].E); 234c4762a1bSJed Brown d[1].T = ihx * (dfluxR[0].T - dfluxL[1].T); 235c4762a1bSJed Brown d[2].E = ihx * dfluxR[1].E; 236c4762a1bSJed Brown d[2].T = ihx * dfluxR[1].T; 237c4762a1bSJed Brown } 238c4762a1bSJed Brown return ihx * (fluxR - fluxL); 239c4762a1bSJed Brown } 240c4762a1bSJed Brown 241d71ae5a4SJacob Faibussowitsch static PetscErrorCode RDGetLocalArrays(RD rd, TS ts, Vec X, Vec Xdot, PetscReal *Theta, PetscReal *dt, Vec *X0loc, RDNode **x0, Vec *Xloc, RDNode **x, Vec *Xloc_t, RDNode **xdot) 242d71ae5a4SJacob Faibussowitsch { 243c4762a1bSJed Brown PetscBool istheta; 244c4762a1bSJed Brown 245c4762a1bSJed Brown PetscFunctionBeginUser; 2469566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(rd->da, X0loc)); 2479566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(rd->da, Xloc)); 2489566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(rd->da, Xloc_t)); 249c4762a1bSJed Brown 2509566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(rd->da, X, INSERT_VALUES, *Xloc)); 2519566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(rd->da, X, INSERT_VALUES, *Xloc)); 2529566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(rd->da, Xdot, INSERT_VALUES, *Xloc_t)); 2539566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(rd->da, Xdot, INSERT_VALUES, *Xloc_t)); 254c4762a1bSJed Brown 255c4762a1bSJed Brown /* 256c4762a1bSJed Brown The following is a hack to subvert TSTHETA which is like an implicit midpoint method to behave more like a trapezoid 257c4762a1bSJed Brown rule. These methods have equivalent linear stability, but the nonlinear stability is somewhat different. The 258c4762a1bSJed Brown radiation system is inconvenient to write in explicit form because the ionization model is "on the left". 259c4762a1bSJed Brown */ 2609566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)ts, TSTHETA, &istheta)); 261c4762a1bSJed Brown if (istheta && rd->endpoint) { 2629566063dSJacob Faibussowitsch PetscCall(TSThetaGetTheta(ts, Theta)); 263c4762a1bSJed Brown } else *Theta = 1.; 264c4762a1bSJed Brown 2659566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ts, dt)); 2669566063dSJacob Faibussowitsch PetscCall(VecWAXPY(*X0loc, -(*Theta) * (*dt), *Xloc_t, *Xloc)); /* back out the value at the start of this step */ 2679371c9d4SSatish Balay if (rd->endpoint) { PetscCall(VecWAXPY(*Xloc, *dt, *Xloc_t, *X0loc)); /* move the abscissa to the end of the step */ } 268c4762a1bSJed Brown 2699566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(rd->da, *X0loc, x0)); 2709566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(rd->da, *Xloc, x)); 2719566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(rd->da, *Xloc_t, xdot)); 272*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 273c4762a1bSJed Brown } 274c4762a1bSJed Brown 275d71ae5a4SJacob Faibussowitsch static PetscErrorCode RDRestoreLocalArrays(RD rd, Vec *X0loc, RDNode **x0, Vec *Xloc, RDNode **x, Vec *Xloc_t, RDNode **xdot) 276d71ae5a4SJacob Faibussowitsch { 277c4762a1bSJed Brown PetscFunctionBeginUser; 2789566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(rd->da, *X0loc, x0)); 2799566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(rd->da, *Xloc, x)); 2809566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(rd->da, *Xloc_t, xdot)); 2819566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(rd->da, X0loc)); 2829566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(rd->da, Xloc)); 2839566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(rd->da, Xloc_t)); 284*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 285c4762a1bSJed Brown } 286c4762a1bSJed Brown 287d71ae5a4SJacob Faibussowitsch static PetscErrorCode PETSC_UNUSED RDCheckDomain_Private(RD rd, TS ts, Vec X, PetscBool *in) 288d71ae5a4SJacob Faibussowitsch { 289c4762a1bSJed Brown PetscInt minloc; 290c4762a1bSJed Brown PetscReal min; 291c4762a1bSJed Brown 292c4762a1bSJed Brown PetscFunctionBeginUser; 2939566063dSJacob Faibussowitsch PetscCall(VecMin(X, &minloc, &min)); 294c4762a1bSJed Brown if (min < 0) { 295c4762a1bSJed Brown SNES snes; 296c4762a1bSJed Brown *in = PETSC_FALSE; 2979566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &snes)); 2989566063dSJacob Faibussowitsch PetscCall(SNESSetFunctionDomainError(snes)); 29963a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(ts, "Domain violation at %" PetscInt_FMT " field %" PetscInt_FMT " value %g\n", minloc / 2, minloc % 2, (double)min)); 300c4762a1bSJed Brown } else *in = PETSC_TRUE; 301*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 302c4762a1bSJed Brown } 303c4762a1bSJed Brown 304c4762a1bSJed Brown /* Energy and temperature must remain positive */ 3059371c9d4SSatish Balay #define RDCheckDomain(rd, ts, X) \ 3069371c9d4SSatish Balay do { \ 307c4762a1bSJed Brown PetscBool _in; \ 3089566063dSJacob Faibussowitsch PetscCall(RDCheckDomain_Private(rd, ts, X, &_in)); \ 309*3ba16761SJacob Faibussowitsch if (!_in) PetscFunctionReturn(PETSC_SUCCESS); \ 310c4762a1bSJed Brown } while (0) 311c4762a1bSJed Brown 312d71ae5a4SJacob Faibussowitsch static PetscErrorCode RDIFunction_FD(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, void *ctx) 313d71ae5a4SJacob Faibussowitsch { 314c4762a1bSJed Brown RD rd = (RD)ctx; 315c4762a1bSJed Brown RDNode *x, *x0, *xdot, *f; 316c4762a1bSJed Brown Vec X0loc, Xloc, Xloc_t; 317c4762a1bSJed Brown PetscReal hx, Theta, dt; 318c4762a1bSJed Brown DMDALocalInfo info; 319c4762a1bSJed Brown PetscInt i; 320c4762a1bSJed Brown 321c4762a1bSJed Brown PetscFunctionBeginUser; 3229566063dSJacob Faibussowitsch PetscCall(RDGetLocalArrays(rd, ts, X, Xdot, &Theta, &dt, &X0loc, &x0, &Xloc, &x, &Xloc_t, &xdot)); 3239566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(rd->da, F, &f)); 3249566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(rd->da, &info)); 325*3ba16761SJacob Faibussowitsch PetscCall(VecZeroEntries(F)); 326c4762a1bSJed Brown 327c4762a1bSJed Brown hx = rd->L / (info.mx - 1); 328c4762a1bSJed Brown 329c4762a1bSJed Brown for (i = info.xs; i < info.xs + info.xm; i++) { 330c4762a1bSJed Brown PetscReal rho = rd->rho; 331c4762a1bSJed Brown PetscScalar Em_t, rad; 332c4762a1bSJed Brown 333c4762a1bSJed Brown rad = (1. - Theta) * RDRadiation(rd, &x0[i], 0) + Theta * RDRadiation(rd, &x[i], 0); 334c4762a1bSJed Brown if (rd->endpoint) { 335c4762a1bSJed Brown PetscScalar Em0, Em1; 336c4762a1bSJed Brown RDMaterialEnergy(rd, &x0[i], &Em0, NULL); 337c4762a1bSJed Brown RDMaterialEnergy(rd, &x[i], &Em1, NULL); 338c4762a1bSJed Brown Em_t = (Em1 - Em0) / dt; 339c4762a1bSJed Brown } else { 340c4762a1bSJed Brown RDNode dEm; 341c4762a1bSJed Brown RDMaterialEnergy(rd, &x[i], NULL, &dEm); 342c4762a1bSJed Brown Em_t = dEm.E * xdot[i].E + dEm.T * xdot[i].T; 343c4762a1bSJed Brown } 344c4762a1bSJed Brown /* Residuals are multiplied by the volume element (hx). */ 345c4762a1bSJed Brown /* The temperature equation does not have boundary conditions */ 346c4762a1bSJed Brown f[i].T = hx * (rho * Em_t + rad); 347c4762a1bSJed Brown 348c4762a1bSJed Brown if (i == 0) { /* Left boundary condition */ 349c4762a1bSJed Brown PetscScalar D_R, bcTheta = rd->bcmidpoint ? Theta : 1.; 350c4762a1bSJed Brown RDNode n, nx; 351c4762a1bSJed Brown 352c4762a1bSJed Brown n.E = (1. - bcTheta) * x0[0].E + bcTheta * x[0].E; 353c4762a1bSJed Brown n.T = (1. - bcTheta) * x0[0].T + bcTheta * x[0].T; 354c4762a1bSJed Brown nx.E = ((1. - bcTheta) * (x0[1].E - x0[0].E) + bcTheta * (x[1].E - x[0].E)) / hx; 355c4762a1bSJed Brown nx.T = ((1. - bcTheta) * (x0[1].T - x0[0].T) + bcTheta * (x[1].T - x[0].T)) / hx; 356c4762a1bSJed Brown switch (rd->leftbc) { 357c4762a1bSJed Brown case BC_ROBIN: 358c4762a1bSJed Brown RDDiffusionCoefficient(rd, rd->bclimit, &n, &nx, &D_R, 0, 0); 359c4762a1bSJed Brown f[0].E = hx * (n.E - 2. * D_R * nx.E - rd->Eapplied); 360c4762a1bSJed Brown break; 361d71ae5a4SJacob Faibussowitsch case BC_NEUMANN: 362d71ae5a4SJacob Faibussowitsch f[0].E = x[1].E - x[0].E; 363d71ae5a4SJacob Faibussowitsch break; 364d71ae5a4SJacob Faibussowitsch default: 365d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Case %" PetscInt_FMT, rd->initial); 366c4762a1bSJed Brown } 367c4762a1bSJed Brown } else if (i == info.mx - 1) { /* Right boundary */ 368c4762a1bSJed Brown f[i].E = x[i].E - x[i - 1].E; /* Homogeneous Neumann */ 369c4762a1bSJed Brown } else { 370c4762a1bSJed Brown PetscScalar diff = (1. - Theta) * RDDiffusion(rd, hx, x0, i, 0) + Theta * RDDiffusion(rd, hx, x, i, 0); 371c4762a1bSJed Brown f[i].E = hx * (xdot[i].E - diff - rad); 372c4762a1bSJed Brown } 373c4762a1bSJed Brown } 3749566063dSJacob Faibussowitsch PetscCall(RDRestoreLocalArrays(rd, &X0loc, &x0, &Xloc, &x, &Xloc_t, &xdot)); 3759566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(rd->da, F, &f)); 3769566063dSJacob Faibussowitsch if (rd->monitor_residual) PetscCall(RDStateView(rd, X, Xdot, F)); 377*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 378c4762a1bSJed Brown } 379c4762a1bSJed Brown 380d71ae5a4SJacob Faibussowitsch static PetscErrorCode RDIJacobian_FD(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal a, Mat A, Mat B, void *ctx) 381d71ae5a4SJacob Faibussowitsch { 382c4762a1bSJed Brown RD rd = (RD)ctx; 383c4762a1bSJed Brown RDNode *x, *x0, *xdot; 384c4762a1bSJed Brown Vec X0loc, Xloc, Xloc_t; 385c4762a1bSJed Brown PetscReal hx, Theta, dt; 386c4762a1bSJed Brown DMDALocalInfo info; 387c4762a1bSJed Brown PetscInt i; 388c4762a1bSJed Brown 389c4762a1bSJed Brown PetscFunctionBeginUser; 3909566063dSJacob Faibussowitsch PetscCall(RDGetLocalArrays(rd, ts, X, Xdot, &Theta, &dt, &X0loc, &x0, &Xloc, &x, &Xloc_t, &xdot)); 3919566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(rd->da, &info)); 392c4762a1bSJed Brown hx = rd->L / (info.mx - 1); 3939566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(B)); 394c4762a1bSJed Brown 395c4762a1bSJed Brown for (i = info.xs; i < info.xs + info.xm; i++) { 396c4762a1bSJed Brown PetscInt col[3]; 397c4762a1bSJed Brown PetscReal rho = rd->rho; 398c4762a1bSJed Brown PetscScalar /*Em_t,rad,*/ K[2][6]; 399c4762a1bSJed Brown RDNode dEm_t, drad; 400c4762a1bSJed Brown 4019371c9d4SSatish Balay /*rad = (1.-Theta)* */ RDRadiation(rd, &x0[i], 0); /* + Theta* */ 4029371c9d4SSatish Balay RDRadiation(rd, &x[i], &drad); 403c4762a1bSJed Brown 404c4762a1bSJed Brown if (rd->endpoint) { 405c4762a1bSJed Brown PetscScalar Em0, Em1; 406c4762a1bSJed Brown RDNode dEm1; 407c4762a1bSJed Brown RDMaterialEnergy(rd, &x0[i], &Em0, NULL); 408c4762a1bSJed Brown RDMaterialEnergy(rd, &x[i], &Em1, &dEm1); 409c4762a1bSJed Brown /*Em_t = (Em1 - Em0) / (Theta*dt);*/ 410c4762a1bSJed Brown dEm_t.E = dEm1.E / (Theta * dt); 411c4762a1bSJed Brown dEm_t.T = dEm1.T / (Theta * dt); 412c4762a1bSJed Brown } else { 413c4762a1bSJed Brown const PetscScalar epsilon = x[i].T * PETSC_SQRT_MACHINE_EPSILON; 414c4762a1bSJed Brown RDNode n1; 415c4762a1bSJed Brown RDNode dEm, dEm1; 416c4762a1bSJed Brown PetscScalar Em_TT; 417c4762a1bSJed Brown 418c4762a1bSJed Brown n1.E = x[i].E; 419c4762a1bSJed Brown n1.T = x[i].T + epsilon; 420c4762a1bSJed Brown RDMaterialEnergy(rd, &x[i], NULL, &dEm); 421c4762a1bSJed Brown RDMaterialEnergy(rd, &n1, NULL, &dEm1); 422c4762a1bSJed Brown /* The Jacobian needs another derivative. We finite difference here instead of 423c4762a1bSJed Brown * propagating second derivatives through the ionization model. */ 424c4762a1bSJed Brown Em_TT = (dEm1.T - dEm.T) / epsilon; 425c4762a1bSJed Brown /*Em_t = dEm.E * xdot[i].E + dEm.T * xdot[i].T;*/ 426c4762a1bSJed Brown dEm_t.E = dEm.E * a; 427c4762a1bSJed Brown dEm_t.T = dEm.T * a + Em_TT * xdot[i].T; 428c4762a1bSJed Brown } 429c4762a1bSJed Brown 4309566063dSJacob Faibussowitsch PetscCall(PetscMemzero(K, sizeof(K))); 431c4762a1bSJed Brown /* Residuals are multiplied by the volume element (hx). */ 432c4762a1bSJed Brown if (i == 0) { 433c4762a1bSJed Brown PetscScalar D, bcTheta = rd->bcmidpoint ? Theta : 1.; 434c4762a1bSJed Brown RDNode n, nx; 435c4762a1bSJed Brown RDNode dD, dxD; 436c4762a1bSJed Brown 437c4762a1bSJed Brown n.E = (1. - bcTheta) * x0[0].E + bcTheta * x[0].E; 438c4762a1bSJed Brown n.T = (1. - bcTheta) * x0[0].T + bcTheta * x[0].T; 439c4762a1bSJed Brown nx.E = ((1. - bcTheta) * (x0[1].E - x0[0].E) + bcTheta * (x[1].E - x[0].E)) / hx; 440c4762a1bSJed Brown nx.T = ((1. - bcTheta) * (x0[1].T - x0[0].T) + bcTheta * (x[1].T - x[0].T)) / hx; 441c4762a1bSJed Brown switch (rd->leftbc) { 442c4762a1bSJed Brown case BC_ROBIN: 443c4762a1bSJed Brown RDDiffusionCoefficient(rd, rd->bclimit, &n, &nx, &D, &dD, &dxD); 444c4762a1bSJed Brown K[0][1 * 2 + 0] = (bcTheta / Theta) * hx * (1. - 2. * D * (-1. / hx) - 2. * nx.E * dD.E + 2. * nx.E * dxD.E / hx); 445c4762a1bSJed Brown K[0][1 * 2 + 1] = (bcTheta / Theta) * hx * (-2. * nx.E * dD.T); 446c4762a1bSJed Brown K[0][2 * 2 + 0] = (bcTheta / Theta) * hx * (-2. * D * (1. / hx) - 2. * nx.E * dD.E - 2. * nx.E * dxD.E / hx); 447c4762a1bSJed Brown break; 448c4762a1bSJed Brown case BC_NEUMANN: 449c4762a1bSJed Brown K[0][1 * 2 + 0] = -1. / Theta; 450c4762a1bSJed Brown K[0][2 * 2 + 0] = 1. / Theta; 451c4762a1bSJed Brown break; 452d71ae5a4SJacob Faibussowitsch default: 453d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Case %" PetscInt_FMT, rd->initial); 454c4762a1bSJed Brown } 455c4762a1bSJed Brown } else if (i == info.mx - 1) { 456c4762a1bSJed Brown K[0][0 * 2 + 0] = -1. / Theta; 457c4762a1bSJed Brown K[0][1 * 2 + 0] = 1. / Theta; 458c4762a1bSJed Brown } else { 459c4762a1bSJed Brown /*PetscScalar diff;*/ 460c4762a1bSJed Brown RDNode ddiff[3]; 461c4762a1bSJed Brown /*diff = (1.-Theta)*RDDiffusion(rd,hx,x0,i,0) + Theta* */ RDDiffusion(rd, hx, x, i, ddiff); 462c4762a1bSJed Brown K[0][0 * 2 + 0] = -hx * ddiff[0].E; 463c4762a1bSJed Brown K[0][0 * 2 + 1] = -hx * ddiff[0].T; 464c4762a1bSJed Brown K[0][1 * 2 + 0] = hx * (a - ddiff[1].E - drad.E); 465c4762a1bSJed Brown K[0][1 * 2 + 1] = hx * (-ddiff[1].T - drad.T); 466c4762a1bSJed Brown K[0][2 * 2 + 0] = -hx * ddiff[2].E; 467c4762a1bSJed Brown K[0][2 * 2 + 1] = -hx * ddiff[2].T; 468c4762a1bSJed Brown } 469c4762a1bSJed Brown 470c4762a1bSJed Brown K[1][1 * 2 + 0] = hx * (rho * dEm_t.E + drad.E); 471c4762a1bSJed Brown K[1][1 * 2 + 1] = hx * (rho * dEm_t.T + drad.T); 472c4762a1bSJed Brown 473c4762a1bSJed Brown col[0] = i - 1; 474c4762a1bSJed Brown col[1] = i; 475c4762a1bSJed Brown col[2] = i + 1 < info.mx ? i + 1 : -1; 4769566063dSJacob Faibussowitsch PetscCall(MatSetValuesBlocked(B, 1, &i, 3, col, &K[0][0], INSERT_VALUES)); 477c4762a1bSJed Brown } 4789566063dSJacob Faibussowitsch PetscCall(RDRestoreLocalArrays(rd, &X0loc, &x0, &Xloc, &x, &Xloc_t, &xdot)); 4799566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 4809566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 481c4762a1bSJed Brown if (A != B) { 4829566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 4839566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 484c4762a1bSJed Brown } 485*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 486c4762a1bSJed Brown } 487c4762a1bSJed Brown 488c4762a1bSJed Brown /* Evaluate interpolants and derivatives at a select quadrature point */ 489d71ae5a4SJacob Faibussowitsch static void RDEvaluate(PetscReal interp[][2], PetscReal deriv[][2], PetscInt q, const RDNode x[], PetscInt i, RDNode *n, RDNode *nx) 490d71ae5a4SJacob Faibussowitsch { 491c4762a1bSJed Brown PetscInt j; 4929371c9d4SSatish Balay n->E = 0; 4939371c9d4SSatish Balay n->T = 0; 4949371c9d4SSatish Balay nx->E = 0; 4959371c9d4SSatish Balay nx->T = 0; 496c4762a1bSJed Brown for (j = 0; j < 2; j++) { 497c4762a1bSJed Brown n->E += interp[q][j] * x[i + j].E; 498c4762a1bSJed Brown n->T += interp[q][j] * x[i + j].T; 499c4762a1bSJed Brown nx->E += deriv[q][j] * x[i + j].E; 500c4762a1bSJed Brown nx->T += deriv[q][j] * x[i + j].T; 501c4762a1bSJed Brown } 502c4762a1bSJed Brown } 503c4762a1bSJed Brown 504c4762a1bSJed Brown /* 505c4762a1bSJed Brown Various quadrature rules. The nonlinear terms are non-polynomial so no standard quadrature will be exact. 506c4762a1bSJed Brown */ 507d71ae5a4SJacob Faibussowitsch static PetscErrorCode RDGetQuadrature(RD rd, PetscReal hx, PetscInt *nq, PetscReal weight[], PetscReal interp[][2], PetscReal deriv[][2]) 508d71ae5a4SJacob Faibussowitsch { 509c4762a1bSJed Brown PetscInt q, j; 510c4762a1bSJed Brown const PetscReal *refweight, (*refinterp)[2], (*refderiv)[2]; 511c4762a1bSJed Brown 512c4762a1bSJed Brown PetscFunctionBeginUser; 513c4762a1bSJed Brown switch (rd->quadrature) { 514c4762a1bSJed Brown case QUADRATURE_GAUSS1: { 51597a1619fSSatish Balay static const PetscReal ww[1] = {1.}; 51697a1619fSSatish Balay static const PetscReal ii[1][2] = { 51797a1619fSSatish Balay {0.5, 0.5} 51897a1619fSSatish Balay }; 51997a1619fSSatish Balay static const PetscReal dd[1][2] = { 52097a1619fSSatish Balay {-1., 1.} 52197a1619fSSatish Balay }; 5229371c9d4SSatish Balay *nq = 1; 5239371c9d4SSatish Balay refweight = ww; 5249371c9d4SSatish Balay refinterp = ii; 5259371c9d4SSatish Balay refderiv = dd; 526c4762a1bSJed Brown } break; 527c4762a1bSJed Brown case QUADRATURE_GAUSS2: { 52897a1619fSSatish Balay static const PetscReal ii[2][2] = { 5299371c9d4SSatish Balay {0.78867513459481287, 0.21132486540518713}, 5309371c9d4SSatish Balay {0.21132486540518713, 0.78867513459481287} 53197a1619fSSatish Balay }; 53297a1619fSSatish Balay static const PetscReal dd[2][2] = { 53397a1619fSSatish Balay {-1., 1.}, 53497a1619fSSatish Balay {-1., 1.} 53597a1619fSSatish Balay }; 53697a1619fSSatish Balay static const PetscReal ww[2] = {0.5, 0.5}; 5379371c9d4SSatish Balay *nq = 2; 5389371c9d4SSatish Balay refweight = ww; 5399371c9d4SSatish Balay refinterp = ii; 5409371c9d4SSatish Balay refderiv = dd; 541c4762a1bSJed Brown } break; 542c4762a1bSJed Brown case QUADRATURE_GAUSS3: { 54397a1619fSSatish Balay static const PetscReal ii[3][2] = { 5449371c9d4SSatish Balay {0.8872983346207417, 0.1127016653792583}, 5459371c9d4SSatish Balay {0.5, 0.5 }, 5469371c9d4SSatish Balay {0.1127016653792583, 0.8872983346207417} 54797a1619fSSatish Balay }; 54897a1619fSSatish Balay static const PetscReal dd[3][2] = { 54997a1619fSSatish Balay {-1, 1}, 55097a1619fSSatish Balay {-1, 1}, 55197a1619fSSatish Balay {-1, 1} 55297a1619fSSatish Balay }; 55397a1619fSSatish Balay static const PetscReal ww[3] = {5. / 18, 8. / 18, 5. / 18}; 5549371c9d4SSatish Balay *nq = 3; 5559371c9d4SSatish Balay refweight = ww; 5569371c9d4SSatish Balay refinterp = ii; 5579371c9d4SSatish Balay refderiv = dd; 558c4762a1bSJed Brown } break; 559c4762a1bSJed Brown case QUADRATURE_GAUSS4: { 56097a1619fSSatish Balay static const PetscReal ii[][2] = { 5619371c9d4SSatish Balay {0.93056815579702623, 0.069431844202973658}, 562c4762a1bSJed Brown {0.66999052179242813, 0.33000947820757187 }, 563c4762a1bSJed Brown {0.33000947820757187, 0.66999052179242813 }, 5649371c9d4SSatish Balay {0.069431844202973658, 0.93056815579702623 } 56597a1619fSSatish Balay }; 56697a1619fSSatish Balay static const PetscReal dd[][2] = { 56797a1619fSSatish Balay {-1, 1}, 56897a1619fSSatish Balay {-1, 1}, 56997a1619fSSatish Balay {-1, 1}, 57097a1619fSSatish Balay {-1, 1} 57197a1619fSSatish Balay }; 57297a1619fSSatish Balay static const PetscReal ww[] = {0.17392742256872692, 0.3260725774312731, 0.3260725774312731, 0.17392742256872692}; 573c4762a1bSJed Brown 5749371c9d4SSatish Balay *nq = 4; 5759371c9d4SSatish Balay refweight = ww; 5769371c9d4SSatish Balay refinterp = ii; 5779371c9d4SSatish Balay refderiv = dd; 578c4762a1bSJed Brown } break; 579c4762a1bSJed Brown case QUADRATURE_LOBATTO2: { 58097a1619fSSatish Balay static const PetscReal ii[2][2] = { 5819371c9d4SSatish Balay {1., 0.}, 5829371c9d4SSatish Balay {0., 1.} 58397a1619fSSatish Balay }; 58497a1619fSSatish Balay static const PetscReal dd[2][2] = { 58597a1619fSSatish Balay {-1., 1.}, 58697a1619fSSatish Balay {-1., 1.} 58797a1619fSSatish Balay }; 58897a1619fSSatish Balay static const PetscReal ww[2] = {0.5, 0.5}; 5899371c9d4SSatish Balay *nq = 2; 5909371c9d4SSatish Balay refweight = ww; 5919371c9d4SSatish Balay refinterp = ii; 5929371c9d4SSatish Balay refderiv = dd; 593c4762a1bSJed Brown } break; 594c4762a1bSJed Brown case QUADRATURE_LOBATTO3: { 59597a1619fSSatish Balay static const PetscReal ii[3][2] = { 5969371c9d4SSatish Balay {1, 0 }, 5979371c9d4SSatish Balay {0.5, 0.5}, 5989371c9d4SSatish Balay {0, 1 } 59997a1619fSSatish Balay }; 60097a1619fSSatish Balay static const PetscReal dd[3][2] = { 60197a1619fSSatish Balay {-1, 1}, 60297a1619fSSatish Balay {-1, 1}, 60397a1619fSSatish Balay {-1, 1} 60497a1619fSSatish Balay }; 60597a1619fSSatish Balay static const PetscReal ww[3] = {1. / 6, 4. / 6, 1. / 6}; 6069371c9d4SSatish Balay *nq = 3; 6079371c9d4SSatish Balay refweight = ww; 6089371c9d4SSatish Balay refinterp = ii; 6099371c9d4SSatish Balay refderiv = dd; 610c4762a1bSJed Brown } break; 611d71ae5a4SJacob Faibussowitsch default: 612d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unknown quadrature %d", (int)rd->quadrature); 613c4762a1bSJed Brown } 614c4762a1bSJed Brown 615c4762a1bSJed Brown for (q = 0; q < *nq; q++) { 616c4762a1bSJed Brown weight[q] = refweight[q] * hx; 617c4762a1bSJed Brown for (j = 0; j < 2; j++) { 618c4762a1bSJed Brown interp[q][j] = refinterp[q][j]; 619c4762a1bSJed Brown deriv[q][j] = refderiv[q][j] / hx; 620c4762a1bSJed Brown } 621c4762a1bSJed Brown } 622*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 623c4762a1bSJed Brown } 624c4762a1bSJed Brown 625c4762a1bSJed Brown /* 626c4762a1bSJed Brown Finite element version 627c4762a1bSJed Brown */ 628d71ae5a4SJacob Faibussowitsch static PetscErrorCode RDIFunction_FE(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, void *ctx) 629d71ae5a4SJacob Faibussowitsch { 630c4762a1bSJed Brown RD rd = (RD)ctx; 631c4762a1bSJed Brown RDNode *x, *x0, *xdot, *f; 632c4762a1bSJed Brown Vec X0loc, Xloc, Xloc_t, Floc; 633c4762a1bSJed Brown PetscReal hx, Theta, dt, weight[5], interp[5][2], deriv[5][2]; 634c4762a1bSJed Brown DMDALocalInfo info; 635c4762a1bSJed Brown PetscInt i, j, q, nq; 636c4762a1bSJed Brown 637c4762a1bSJed Brown PetscFunctionBeginUser; 6389566063dSJacob Faibussowitsch PetscCall(RDGetLocalArrays(rd, ts, X, Xdot, &Theta, &dt, &X0loc, &x0, &Xloc, &x, &Xloc_t, &xdot)); 639c4762a1bSJed Brown 6409566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(rd->da, &Floc)); 6419566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(Floc)); 6429566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(rd->da, Floc, &f)); 6439566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(rd->da, &info)); 644c4762a1bSJed Brown 645c4762a1bSJed Brown /* Set up shape functions and quadrature for elements (assumes a uniform grid) */ 646c4762a1bSJed Brown hx = rd->L / (info.mx - 1); 6479566063dSJacob Faibussowitsch PetscCall(RDGetQuadrature(rd, hx, &nq, weight, interp, deriv)); 648c4762a1bSJed Brown 649c4762a1bSJed Brown for (i = info.xs; i < PetscMin(info.xs + info.xm, info.mx - 1); i++) { 650c4762a1bSJed Brown for (q = 0; q < nq; q++) { 651c4762a1bSJed Brown PetscReal rho = rd->rho; 652c4762a1bSJed Brown PetscScalar Em_t, rad, D_R, D0_R; 653c4762a1bSJed Brown RDNode n, n0, nx, n0x, nt, ntx; 654c4762a1bSJed Brown RDEvaluate(interp, deriv, q, x, i, &n, &nx); 655c4762a1bSJed Brown RDEvaluate(interp, deriv, q, x0, i, &n0, &n0x); 656c4762a1bSJed Brown RDEvaluate(interp, deriv, q, xdot, i, &nt, &ntx); 657c4762a1bSJed Brown 658c4762a1bSJed Brown rad = (1. - Theta) * RDRadiation(rd, &n0, 0) + Theta * RDRadiation(rd, &n, 0); 659c4762a1bSJed Brown if (rd->endpoint) { 660c4762a1bSJed Brown PetscScalar Em0, Em1; 661c4762a1bSJed Brown RDMaterialEnergy(rd, &n0, &Em0, NULL); 662c4762a1bSJed Brown RDMaterialEnergy(rd, &n, &Em1, NULL); 663c4762a1bSJed Brown Em_t = (Em1 - Em0) / dt; 664c4762a1bSJed Brown } else { 665c4762a1bSJed Brown RDNode dEm; 666c4762a1bSJed Brown RDMaterialEnergy(rd, &n, NULL, &dEm); 667c4762a1bSJed Brown Em_t = dEm.E * nt.E + dEm.T * nt.T; 668c4762a1bSJed Brown } 669c4762a1bSJed Brown RDDiffusionCoefficient(rd, PETSC_TRUE, &n0, &n0x, &D0_R, 0, 0); 670c4762a1bSJed Brown RDDiffusionCoefficient(rd, PETSC_TRUE, &n, &nx, &D_R, 0, 0); 671c4762a1bSJed Brown for (j = 0; j < 2; j++) { 6729371c9d4SSatish Balay f[i + j].E += (deriv[q][j] * weight[q] * ((1. - Theta) * D0_R * n0x.E + Theta * D_R * nx.E) + interp[q][j] * weight[q] * (nt.E - rad)); 673c4762a1bSJed Brown f[i + j].T += interp[q][j] * weight[q] * (rho * Em_t + rad); 674c4762a1bSJed Brown } 675c4762a1bSJed Brown } 676c4762a1bSJed Brown } 677c4762a1bSJed Brown if (info.xs == 0) { 678c4762a1bSJed Brown switch (rd->leftbc) { 679c4762a1bSJed Brown case BC_ROBIN: { 680c4762a1bSJed Brown PetscScalar D_R, D_R_bc; 681c4762a1bSJed Brown PetscReal ratio, bcTheta = rd->bcmidpoint ? Theta : 1.; 682c4762a1bSJed Brown RDNode n, nx; 683c4762a1bSJed Brown 684c4762a1bSJed Brown n.E = (1 - bcTheta) * x0[0].E + bcTheta * x[0].E; 685c4762a1bSJed Brown n.T = (1 - bcTheta) * x0[0].T + bcTheta * x[0].T; 686c4762a1bSJed Brown nx.E = (x[1].E - x[0].E) / hx; 687c4762a1bSJed Brown nx.T = (x[1].T - x[0].T) / hx; 688c4762a1bSJed Brown RDDiffusionCoefficient(rd, PETSC_TRUE, &n, &nx, &D_R, 0, 0); 689c4762a1bSJed Brown RDDiffusionCoefficient(rd, rd->bclimit, &n, &nx, &D_R_bc, 0, 0); 690c4762a1bSJed Brown ratio = PetscRealPart(D_R / D_R_bc); 6913c633725SBarry Smith PetscCheck(ratio <= 1., PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Limited diffusivity is greater than unlimited"); 6923c633725SBarry Smith PetscCheck(ratio >= 1e-3, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Heavily limited diffusivity"); 693c4762a1bSJed Brown f[0].E += -ratio * 0.5 * (rd->Eapplied - n.E); 694c4762a1bSJed Brown } break; 695c4762a1bSJed Brown case BC_NEUMANN: 696c4762a1bSJed Brown /* homogeneous Neumann is the natural condition */ 697c4762a1bSJed Brown break; 698d71ae5a4SJacob Faibussowitsch default: 699d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Case %" PetscInt_FMT, rd->initial); 700c4762a1bSJed Brown } 701c4762a1bSJed Brown } 702c4762a1bSJed Brown 7039566063dSJacob Faibussowitsch PetscCall(RDRestoreLocalArrays(rd, &X0loc, &x0, &Xloc, &x, &Xloc_t, &xdot)); 7049566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(rd->da, Floc, &f)); 7059566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(F)); 7069566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(rd->da, Floc, ADD_VALUES, F)); 7079566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(rd->da, Floc, ADD_VALUES, F)); 7089566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(rd->da, &Floc)); 709c4762a1bSJed Brown 7109566063dSJacob Faibussowitsch if (rd->monitor_residual) PetscCall(RDStateView(rd, X, Xdot, F)); 711*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 712c4762a1bSJed Brown } 713c4762a1bSJed Brown 714d71ae5a4SJacob Faibussowitsch static PetscErrorCode RDIJacobian_FE(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal a, Mat A, Mat B, void *ctx) 715d71ae5a4SJacob Faibussowitsch { 716c4762a1bSJed Brown RD rd = (RD)ctx; 717c4762a1bSJed Brown RDNode *x, *x0, *xdot; 718c4762a1bSJed Brown Vec X0loc, Xloc, Xloc_t; 719c4762a1bSJed Brown PetscReal hx, Theta, dt, weight[5], interp[5][2], deriv[5][2]; 720c4762a1bSJed Brown DMDALocalInfo info; 721c4762a1bSJed Brown PetscInt i, j, k, q, nq; 722c4762a1bSJed Brown PetscScalar K[4][4]; 723c4762a1bSJed Brown 724c4762a1bSJed Brown PetscFunctionBeginUser; 7259566063dSJacob Faibussowitsch PetscCall(RDGetLocalArrays(rd, ts, X, Xdot, &Theta, &dt, &X0loc, &x0, &Xloc, &x, &Xloc_t, &xdot)); 7269566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(rd->da, &info)); 727c4762a1bSJed Brown hx = rd->L / (info.mx - 1); 7289566063dSJacob Faibussowitsch PetscCall(RDGetQuadrature(rd, hx, &nq, weight, interp, deriv)); 7299566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(B)); 730c4762a1bSJed Brown for (i = info.xs; i < PetscMin(info.xs + info.xm, info.mx - 1); i++) { 731c4762a1bSJed Brown PetscInt rc[2]; 732c4762a1bSJed Brown 7339371c9d4SSatish Balay rc[0] = i; 7349371c9d4SSatish Balay rc[1] = i + 1; 7359566063dSJacob Faibussowitsch PetscCall(PetscMemzero(K, sizeof(K))); 736c4762a1bSJed Brown for (q = 0; q < nq; q++) { 737c4762a1bSJed Brown PetscScalar D_R; 738c4762a1bSJed Brown PETSC_UNUSED PetscScalar rad; 739c4762a1bSJed Brown RDNode n, nx, nt, ntx, drad, dD_R, dxD_R, dEm; 740c4762a1bSJed Brown RDEvaluate(interp, deriv, q, x, i, &n, &nx); 741c4762a1bSJed Brown RDEvaluate(interp, deriv, q, xdot, i, &nt, &ntx); 742c4762a1bSJed Brown rad = RDRadiation(rd, &n, &drad); 743c4762a1bSJed Brown RDDiffusionCoefficient(rd, PETSC_TRUE, &n, &nx, &D_R, &dD_R, &dxD_R); 744c4762a1bSJed Brown RDMaterialEnergy(rd, &n, NULL, &dEm); 745c4762a1bSJed Brown for (j = 0; j < 2; j++) { 746c4762a1bSJed Brown for (k = 0; k < 2; k++) { 7479371c9d4SSatish Balay K[j * 2 + 0][k * 2 + 0] += (+interp[q][j] * weight[q] * (a - drad.E) * interp[q][k] + deriv[q][j] * weight[q] * ((D_R + dxD_R.E * nx.E) * deriv[q][k] + dD_R.E * nx.E * interp[q][k])); 7489371c9d4SSatish Balay K[j * 2 + 0][k * 2 + 1] += (+interp[q][j] * weight[q] * (-drad.T * interp[q][k]) + deriv[q][j] * weight[q] * (dxD_R.T * deriv[q][k] + dD_R.T * interp[q][k]) * nx.E); 749c4762a1bSJed Brown K[j * 2 + 1][k * 2 + 0] += interp[q][j] * weight[q] * drad.E * interp[q][k]; 750c4762a1bSJed Brown K[j * 2 + 1][k * 2 + 1] += interp[q][j] * weight[q] * (a * rd->rho * dEm.T + drad.T) * interp[q][k]; 751c4762a1bSJed Brown } 752c4762a1bSJed Brown } 753c4762a1bSJed Brown } 7549566063dSJacob Faibussowitsch PetscCall(MatSetValuesBlocked(B, 2, rc, 2, rc, &K[0][0], ADD_VALUES)); 755c4762a1bSJed Brown } 756c4762a1bSJed Brown if (info.xs == 0) { 757c4762a1bSJed Brown switch (rd->leftbc) { 758c4762a1bSJed Brown case BC_ROBIN: { 759c4762a1bSJed Brown PetscScalar D_R, D_R_bc; 760c4762a1bSJed Brown PetscReal ratio; 761c4762a1bSJed Brown RDNode n, nx; 762c4762a1bSJed Brown 763c4762a1bSJed Brown n.E = (1 - Theta) * x0[0].E + Theta * x[0].E; 764c4762a1bSJed Brown n.T = (1 - Theta) * x0[0].T + Theta * x[0].T; 765c4762a1bSJed Brown nx.E = (x[1].E - x[0].E) / hx; 766c4762a1bSJed Brown nx.T = (x[1].T - x[0].T) / hx; 767c4762a1bSJed Brown RDDiffusionCoefficient(rd, PETSC_TRUE, &n, &nx, &D_R, 0, 0); 768c4762a1bSJed Brown RDDiffusionCoefficient(rd, rd->bclimit, &n, &nx, &D_R_bc, 0, 0); 769c4762a1bSJed Brown ratio = PetscRealPart(D_R / D_R_bc); 7709566063dSJacob Faibussowitsch PetscCall(MatSetValue(B, 0, 0, ratio * 0.5, ADD_VALUES)); 771c4762a1bSJed Brown } break; 772c4762a1bSJed Brown case BC_NEUMANN: 773c4762a1bSJed Brown /* homogeneous Neumann is the natural condition */ 774c4762a1bSJed Brown break; 775d71ae5a4SJacob Faibussowitsch default: 776d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Case %" PetscInt_FMT, rd->initial); 777c4762a1bSJed Brown } 778c4762a1bSJed Brown } 779c4762a1bSJed Brown 7809566063dSJacob Faibussowitsch PetscCall(RDRestoreLocalArrays(rd, &X0loc, &x0, &Xloc, &x, &Xloc_t, &xdot)); 7819566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 7829566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 783c4762a1bSJed Brown if (A != B) { 7849566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 7859566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 786c4762a1bSJed Brown } 787*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 788c4762a1bSJed Brown } 789c4762a1bSJed Brown 790c4762a1bSJed Brown /* Temperature that is in equilibrium with the radiation density */ 791d71ae5a4SJacob Faibussowitsch static PetscScalar RDRadiationTemperature(RD rd, PetscScalar E) 792d71ae5a4SJacob Faibussowitsch { 7939371c9d4SSatish Balay return PetscPowScalar(E * rd->c / (4. * rd->sigma_b), 0.25); 7949371c9d4SSatish Balay } 795c4762a1bSJed Brown 796d71ae5a4SJacob Faibussowitsch static PetscErrorCode RDInitialState(RD rd, Vec X) 797d71ae5a4SJacob Faibussowitsch { 798c4762a1bSJed Brown DMDALocalInfo info; 799c4762a1bSJed Brown PetscInt i; 800c4762a1bSJed Brown RDNode *x; 801c4762a1bSJed Brown 802c4762a1bSJed Brown PetscFunctionBeginUser; 8039566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(rd->da, &info)); 8049566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(rd->da, X, &x)); 805c4762a1bSJed Brown for (i = info.xs; i < info.xs + info.xm; i++) { 806c4762a1bSJed Brown PetscReal coord = i * rd->L / (info.mx - 1); 807c4762a1bSJed Brown switch (rd->initial) { 808c4762a1bSJed Brown case 1: 809c4762a1bSJed Brown x[i].E = 0.001; 810c4762a1bSJed Brown x[i].T = RDRadiationTemperature(rd, x[i].E); 811c4762a1bSJed Brown break; 812c4762a1bSJed Brown case 2: 813c4762a1bSJed Brown x[i].E = 0.001 + 100. * PetscExpReal(-PetscSqr(coord / 0.1)); 814c4762a1bSJed Brown x[i].T = RDRadiationTemperature(rd, x[i].E); 815c4762a1bSJed Brown break; 816c4762a1bSJed Brown case 3: 817c4762a1bSJed Brown x[i].E = 7.56e-2 * rd->unit.Joule / PetscPowScalarInt(rd->unit.meter, 3); 818c4762a1bSJed Brown x[i].T = RDRadiationTemperature(rd, x[i].E); 819c4762a1bSJed Brown break; 820d71ae5a4SJacob Faibussowitsch default: 821d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "No initial state %" PetscInt_FMT, rd->initial); 822c4762a1bSJed Brown } 823c4762a1bSJed Brown } 8249566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(rd->da, X, &x)); 825*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 826c4762a1bSJed Brown } 827c4762a1bSJed Brown 828d71ae5a4SJacob Faibussowitsch static PetscErrorCode RDView(RD rd, Vec X, PetscViewer viewer) 829d71ae5a4SJacob Faibussowitsch { 830c4762a1bSJed Brown Vec Y; 831c4762a1bSJed Brown const RDNode *x; 832c4762a1bSJed Brown PetscScalar *y; 833c4762a1bSJed Brown PetscInt i, m, M; 834c4762a1bSJed Brown const PetscInt *lx; 835c4762a1bSJed Brown DM da; 836c4762a1bSJed Brown MPI_Comm comm; 837c4762a1bSJed Brown 838c4762a1bSJed Brown PetscFunctionBeginUser; 839c4762a1bSJed Brown /* 840c4762a1bSJed Brown Create a DMDA (one dof per node, zero stencil width, same layout) to hold Trad 841c4762a1bSJed Brown (radiation temperature). It is not necessary to create a DMDA for this, but this way 842c4762a1bSJed Brown output and visualization will have meaningful variable names and correct scales. 843c4762a1bSJed Brown */ 8449566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(rd->da, 0, &M, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0)); 8459566063dSJacob Faibussowitsch PetscCall(DMDAGetOwnershipRanges(rd->da, &lx, 0, 0)); 8469566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)rd->da, &comm)); 8479566063dSJacob Faibussowitsch PetscCall(DMDACreate1d(comm, DM_BOUNDARY_NONE, M, 1, 0, lx, &da)); 8489566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da)); 8499566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 8509566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(da, 0., rd->L, 0., 0., 0., 0.)); 8519566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(da, 0, "T_rad")); 8529566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da, &Y)); 853c4762a1bSJed Brown 854c4762a1bSJed Brown /* Compute the radiation temperature from the solution at each node */ 8559566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(Y, &m)); 8569566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, (const PetscScalar **)&x)); 8579566063dSJacob Faibussowitsch PetscCall(VecGetArray(Y, &y)); 858c4762a1bSJed Brown for (i = 0; i < m; i++) y[i] = RDRadiationTemperature(rd, x[i].E); 8599566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, (const PetscScalar **)&x)); 8609566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Y, &y)); 861c4762a1bSJed Brown 8629566063dSJacob Faibussowitsch PetscCall(VecView(Y, viewer)); 8639566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Y)); 8649566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da)); 865*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 866c4762a1bSJed Brown } 867c4762a1bSJed Brown 868d71ae5a4SJacob Faibussowitsch static PetscErrorCode RDTestDifferentiation(RD rd) 869d71ae5a4SJacob Faibussowitsch { 870c4762a1bSJed Brown MPI_Comm comm; 871c4762a1bSJed Brown RDNode n, nx; 872c4762a1bSJed Brown PetscScalar epsilon; 873c4762a1bSJed Brown 874c4762a1bSJed Brown PetscFunctionBeginUser; 8759566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)rd->da, &comm)); 876c4762a1bSJed Brown epsilon = 1e-8; 877c4762a1bSJed Brown { 878c4762a1bSJed Brown RDNode dEm, fdEm; 879c4762a1bSJed Brown PetscScalar T0 = 1000., T1 = T0 * (1. + epsilon), Em0, Em1; 880c4762a1bSJed Brown n.E = 1.; 881c4762a1bSJed Brown n.T = T0; 882c4762a1bSJed Brown rd->MaterialEnergy(rd, &n, &Em0, &dEm); 883c4762a1bSJed Brown n.E = 1. + epsilon; 884c4762a1bSJed Brown n.T = T0; 885c4762a1bSJed Brown rd->MaterialEnergy(rd, &n, &Em1, 0); 886c4762a1bSJed Brown fdEm.E = (Em1 - Em0) / epsilon; 887c4762a1bSJed Brown n.E = 1.; 888c4762a1bSJed Brown n.T = T1; 889c4762a1bSJed Brown rd->MaterialEnergy(rd, &n, &Em1, 0); 890c4762a1bSJed Brown fdEm.T = (Em1 - Em0) / (T0 * epsilon); 8919371c9d4SSatish Balay PetscCall(PetscPrintf(comm, "dEm {%g,%g}, fdEm {%g,%g}, diff {%g,%g}\n", (double)PetscRealPart(dEm.E), (double)PetscRealPart(dEm.T), (double)PetscRealPart(fdEm.E), (double)PetscRealPart(fdEm.T), (double)PetscRealPart(dEm.E - fdEm.E), 8929371c9d4SSatish Balay (double)PetscRealPart(dEm.T - fdEm.T))); 893c4762a1bSJed Brown } 894c4762a1bSJed Brown { 895c4762a1bSJed Brown PetscScalar D0, D; 896c4762a1bSJed Brown RDNode dD, dxD, fdD, fdxD; 8979371c9d4SSatish Balay n.E = 1.; 8989371c9d4SSatish Balay n.T = 1.; 8999371c9d4SSatish Balay nx.E = 1.; 9009371c9d4SSatish Balay n.T = 1.; 901c4762a1bSJed Brown RDDiffusionCoefficient(rd, rd->bclimit, &n, &nx, &D0, &dD, &dxD); 9029371c9d4SSatish Balay n.E = 1. + epsilon; 9039371c9d4SSatish Balay n.T = 1.; 9049371c9d4SSatish Balay nx.E = 1.; 9059371c9d4SSatish Balay n.T = 1.; 9069371c9d4SSatish Balay RDDiffusionCoefficient(rd, rd->bclimit, &n, &nx, &D, 0, 0); 9079371c9d4SSatish Balay fdD.E = (D - D0) / epsilon; 9089371c9d4SSatish Balay n.E = 1; 9099371c9d4SSatish Balay n.T = 1. + epsilon; 9109371c9d4SSatish Balay nx.E = 1.; 9119371c9d4SSatish Balay n.T = 1.; 9129371c9d4SSatish Balay RDDiffusionCoefficient(rd, rd->bclimit, &n, &nx, &D, 0, 0); 9139371c9d4SSatish Balay fdD.T = (D - D0) / epsilon; 9149371c9d4SSatish Balay n.E = 1; 9159371c9d4SSatish Balay n.T = 1.; 9169371c9d4SSatish Balay nx.E = 1. + epsilon; 9179371c9d4SSatish Balay n.T = 1.; 9189371c9d4SSatish Balay RDDiffusionCoefficient(rd, rd->bclimit, &n, &nx, &D, 0, 0); 9199371c9d4SSatish Balay fdxD.E = (D - D0) / epsilon; 9209371c9d4SSatish Balay n.E = 1; 9219371c9d4SSatish Balay n.T = 1.; 9229371c9d4SSatish Balay nx.E = 1.; 9239371c9d4SSatish Balay n.T = 1. + epsilon; 9249371c9d4SSatish Balay RDDiffusionCoefficient(rd, rd->bclimit, &n, &nx, &D, 0, 0); 9259371c9d4SSatish Balay fdxD.T = (D - D0) / epsilon; 9269371c9d4SSatish Balay PetscCall(PetscPrintf(comm, "dD {%g,%g}, fdD {%g,%g}, diff {%g,%g}\n", (double)PetscRealPart(dD.E), (double)PetscRealPart(dD.T), (double)PetscRealPart(fdD.E), (double)PetscRealPart(fdD.T), (double)PetscRealPart(dD.E - fdD.E), 9279371c9d4SSatish Balay (double)PetscRealPart(dD.T - fdD.T))); 9289371c9d4SSatish Balay PetscCall(PetscPrintf(comm, "dxD {%g,%g}, fdxD {%g,%g}, diffx {%g,%g}\n", (double)PetscRealPart(dxD.E), (double)PetscRealPart(dxD.T), (double)PetscRealPart(fdxD.E), (double)PetscRealPart(fdxD.T), (double)PetscRealPart(dxD.E - fdxD.E), 9299371c9d4SSatish Balay (double)PetscRealPart(dxD.T - fdxD.T))); 930c4762a1bSJed Brown } 931c4762a1bSJed Brown { 932c4762a1bSJed Brown PetscInt i; 933c4762a1bSJed Brown PetscReal hx = 1.; 934c4762a1bSJed Brown PetscScalar a0; 935c4762a1bSJed Brown RDNode n0[3], n1[3], d[3], fd[3]; 936c4762a1bSJed Brown 937c4762a1bSJed Brown n0[0].E = 1.; 938c4762a1bSJed Brown n0[0].T = 1.; 939c4762a1bSJed Brown n0[1].E = 5.; 940c4762a1bSJed Brown n0[1].T = 3.; 941c4762a1bSJed Brown n0[2].E = 4.; 942c4762a1bSJed Brown n0[2].T = 2.; 943c4762a1bSJed Brown a0 = RDDiffusion(rd, hx, n0, 1, d); 944c4762a1bSJed Brown for (i = 0; i < 3; i++) { 9459371c9d4SSatish Balay PetscCall(PetscMemcpy(n1, n0, sizeof(n0))); 9469371c9d4SSatish Balay n1[i].E += epsilon; 947c4762a1bSJed Brown fd[i].E = (RDDiffusion(rd, hx, n1, 1, 0) - a0) / epsilon; 9489371c9d4SSatish Balay PetscCall(PetscMemcpy(n1, n0, sizeof(n0))); 9499371c9d4SSatish Balay n1[i].T += epsilon; 950c4762a1bSJed Brown fd[i].T = (RDDiffusion(rd, hx, n1, 1, 0) - a0) / epsilon; 9519371c9d4SSatish Balay PetscCall(PetscPrintf(comm, "ddiff[%" PetscInt_FMT "] {%g,%g}, fd {%g %g}, diff {%g,%g}\n", i, (double)PetscRealPart(d[i].E), (double)PetscRealPart(d[i].T), (double)PetscRealPart(fd[i].E), (double)PetscRealPart(fd[i].T), 9529371c9d4SSatish Balay (double)PetscRealPart(d[i].E - fd[i].E), (double)PetscRealPart(d[i].T - fd[i].T))); 953c4762a1bSJed Brown } 954c4762a1bSJed Brown } 955c4762a1bSJed Brown { 956c4762a1bSJed Brown PetscScalar rad0, rad; 957c4762a1bSJed Brown RDNode drad, fdrad; 9589371c9d4SSatish Balay n.E = 1.; 9599371c9d4SSatish Balay n.T = 1.; 960c4762a1bSJed Brown rad0 = RDRadiation(rd, &n, &drad); 9619371c9d4SSatish Balay n.E = 1. + epsilon; 9629371c9d4SSatish Balay n.T = 1.; 9639371c9d4SSatish Balay rad = RDRadiation(rd, &n, 0); 9649371c9d4SSatish Balay fdrad.E = (rad - rad0) / epsilon; 9659371c9d4SSatish Balay n.E = 1.; 9669371c9d4SSatish Balay n.T = 1. + epsilon; 9679371c9d4SSatish Balay rad = RDRadiation(rd, &n, 0); 9689371c9d4SSatish Balay fdrad.T = (rad - rad0) / epsilon; 9699371c9d4SSatish Balay PetscCall(PetscPrintf(comm, "drad {%g,%g}, fdrad {%g,%g}, diff {%g,%g}\n", (double)PetscRealPart(drad.E), (double)PetscRealPart(drad.T), (double)PetscRealPart(fdrad.E), (double)PetscRealPart(fdrad.T), (double)PetscRealPart(drad.E - drad.E), 9709371c9d4SSatish Balay (double)PetscRealPart(drad.T - fdrad.T))); 971c4762a1bSJed Brown } 972*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 973c4762a1bSJed Brown } 974c4762a1bSJed Brown 975d71ae5a4SJacob Faibussowitsch static PetscErrorCode RDCreate(MPI_Comm comm, RD *inrd) 976d71ae5a4SJacob Faibussowitsch { 977c4762a1bSJed Brown RD rd; 978c4762a1bSJed Brown PetscReal meter = 0, kilogram = 0, second = 0, Kelvin = 0, Joule = 0, Watt = 0; 979c4762a1bSJed Brown 980c4762a1bSJed Brown PetscFunctionBeginUser; 981c4762a1bSJed Brown *inrd = 0; 9829566063dSJacob Faibussowitsch PetscCall(PetscNew(&rd)); 983c4762a1bSJed Brown 984d0609cedSBarry Smith PetscOptionsBegin(comm, NULL, "Options for nonequilibrium radiation-diffusion with RD ionization", NULL); 985c4762a1bSJed Brown { 986c4762a1bSJed Brown rd->initial = 1; 9879566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-rd_initial", "Initial condition (1=Marshak, 2=Blast, 3=Marshak+)", "", rd->initial, &rd->initial, 0)); 988c4762a1bSJed Brown switch (rd->initial) { 989c4762a1bSJed Brown case 1: 990c4762a1bSJed Brown case 2: 991c4762a1bSJed Brown rd->unit.kilogram = 1.; 992c4762a1bSJed Brown rd->unit.meter = 1.; 993c4762a1bSJed Brown rd->unit.second = 1.; 994c4762a1bSJed Brown rd->unit.Kelvin = 1.; 995c4762a1bSJed Brown break; 996c4762a1bSJed Brown case 3: 997c4762a1bSJed Brown rd->unit.kilogram = 1.e12; 998c4762a1bSJed Brown rd->unit.meter = 1.; 999c4762a1bSJed Brown rd->unit.second = 1.e9; 1000c4762a1bSJed Brown rd->unit.Kelvin = 1.; 1001c4762a1bSJed Brown break; 1002d71ae5a4SJacob Faibussowitsch default: 1003d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unknown initial condition %" PetscInt_FMT, rd->initial); 1004c4762a1bSJed Brown } 1005c4762a1bSJed Brown /* Fundamental units */ 10069566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-rd_unit_meter", "Length of 1 meter in nondimensional units", "", rd->unit.meter, &rd->unit.meter, 0)); 10079566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-rd_unit_kilogram", "Mass of 1 kilogram in nondimensional units", "", rd->unit.kilogram, &rd->unit.kilogram, 0)); 10089566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-rd_unit_second", "Time of a second in nondimensional units", "", rd->unit.second, &rd->unit.second, 0)); 10099566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-rd_unit_Kelvin", "Temperature of a Kelvin in nondimensional units", "", rd->unit.Kelvin, &rd->unit.Kelvin, 0)); 1010c4762a1bSJed Brown /* Derived units */ 1011c4762a1bSJed Brown rd->unit.Joule = rd->unit.kilogram * PetscSqr(rd->unit.meter / rd->unit.second); 1012c4762a1bSJed Brown rd->unit.Watt = rd->unit.Joule / rd->unit.second; 1013c4762a1bSJed Brown /* Local aliases */ 1014c4762a1bSJed Brown meter = rd->unit.meter; 1015c4762a1bSJed Brown kilogram = rd->unit.kilogram; 1016c4762a1bSJed Brown second = rd->unit.second; 1017c4762a1bSJed Brown Kelvin = rd->unit.Kelvin; 1018c4762a1bSJed Brown Joule = rd->unit.Joule; 1019c4762a1bSJed Brown Watt = rd->unit.Watt; 1020c4762a1bSJed Brown 10219566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-rd_monitor_residual", "Display residuals every time they are evaluated", "", rd->monitor_residual, &rd->monitor_residual, NULL)); 10229566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-rd_discretization", "Discretization type", "", DiscretizationTypes, (PetscEnum)rd->discretization, (PetscEnum *)&rd->discretization, NULL)); 1023c4762a1bSJed Brown if (rd->discretization == DISCRETIZATION_FE) { 1024c4762a1bSJed Brown rd->quadrature = QUADRATURE_GAUSS2; 10259566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-rd_quadrature", "Finite element quadrature", "", QuadratureTypes, (PetscEnum)rd->quadrature, (PetscEnum *)&rd->quadrature, NULL)); 1026c4762a1bSJed Brown } 10279566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-rd_jacobian", "Type of finite difference Jacobian", "", JacobianTypes, (PetscEnum)rd->jacobian, (PetscEnum *)&rd->jacobian, NULL)); 1028c4762a1bSJed Brown switch (rd->initial) { 1029c4762a1bSJed Brown case 1: 1030c4762a1bSJed Brown rd->leftbc = BC_ROBIN; 1031c4762a1bSJed Brown rd->Eapplied = 4 * rd->unit.Joule / PetscPowRealInt(rd->unit.meter, 3); 1032c4762a1bSJed Brown rd->L = 1. * rd->unit.meter; 1033c4762a1bSJed Brown rd->beta = 3.0; 1034c4762a1bSJed Brown rd->gamma = 3.0; 1035c4762a1bSJed Brown rd->final_time = 3 * second; 1036c4762a1bSJed Brown break; 1037c4762a1bSJed Brown case 2: 1038c4762a1bSJed Brown rd->leftbc = BC_NEUMANN; 1039c4762a1bSJed Brown rd->Eapplied = 0.; 1040c4762a1bSJed Brown rd->L = 1. * rd->unit.meter; 1041c4762a1bSJed Brown rd->beta = 3.0; 1042c4762a1bSJed Brown rd->gamma = 3.0; 1043c4762a1bSJed Brown rd->final_time = 1 * second; 1044c4762a1bSJed Brown break; 1045c4762a1bSJed Brown case 3: 1046c4762a1bSJed Brown rd->leftbc = BC_ROBIN; 1047c4762a1bSJed Brown rd->Eapplied = 7.503e6 * rd->unit.Joule / PetscPowRealInt(rd->unit.meter, 3); 1048c4762a1bSJed Brown rd->L = 5. * rd->unit.meter; 1049c4762a1bSJed Brown rd->beta = 3.5; 1050c4762a1bSJed Brown rd->gamma = 3.5; 1051c4762a1bSJed Brown rd->final_time = 20e-9 * second; 1052c4762a1bSJed Brown break; 1053d71ae5a4SJacob Faibussowitsch default: 1054d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Initial %" PetscInt_FMT, rd->initial); 1055c4762a1bSJed Brown } 10569566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-rd_leftbc", "Left boundary condition", "", BCTypes, (PetscEnum)rd->leftbc, (PetscEnum *)&rd->leftbc, NULL)); 10579566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-rd_E_applied", "Radiation flux at left end of domain", "", rd->Eapplied, &rd->Eapplied, NULL)); 10589566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-rd_beta", "Thermal exponent for photon absorption", "", rd->beta, &rd->beta, NULL)); 10599566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-rd_gamma", "Thermal exponent for diffusion coefficient", "", rd->gamma, &rd->gamma, NULL)); 10609566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-rd_view_draw", "Draw final solution", "", rd->view_draw, &rd->view_draw, NULL)); 10619566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-rd_endpoint", "Discretize using endpoints (like trapezoid rule) instead of midpoint", "", rd->endpoint, &rd->endpoint, NULL)); 10629566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-rd_bcmidpoint", "Impose the boundary condition at the midpoint (Theta) of the interval", "", rd->bcmidpoint, &rd->bcmidpoint, NULL)); 10639566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-rd_bclimit", "Limit diffusion coefficient in definition of Robin boundary condition", "", rd->bclimit, &rd->bclimit, NULL)); 10649566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-rd_test_diff", "Test differentiation in constitutive relations", "", rd->test_diff, &rd->test_diff, NULL)); 10659566063dSJacob Faibussowitsch PetscCall(PetscOptionsString("-rd_view_binary", "File name to hold final solution", "", rd->view_binary, rd->view_binary, sizeof(rd->view_binary), NULL)); 1066c4762a1bSJed Brown } 1067d0609cedSBarry Smith PetscOptionsEnd(); 1068c4762a1bSJed Brown 1069c4762a1bSJed Brown switch (rd->initial) { 1070c4762a1bSJed Brown case 1: 1071c4762a1bSJed Brown case 2: 1072c4762a1bSJed Brown rd->rho = 1.; 1073c4762a1bSJed Brown rd->c = 1.; 1074c4762a1bSJed Brown rd->K_R = 1.; 1075c4762a1bSJed Brown rd->K_p = 1.; 1076c4762a1bSJed Brown rd->sigma_b = 0.25; 1077c4762a1bSJed Brown rd->MaterialEnergy = RDMaterialEnergy_Reduced; 1078c4762a1bSJed Brown break; 1079c4762a1bSJed Brown case 3: 1080c4762a1bSJed Brown /* Table 2 */ 1081c4762a1bSJed Brown rd->rho = 1.17e-3 * kilogram / (meter * meter * meter); /* density */ 1082c4762a1bSJed Brown rd->K_R = 7.44e18 * PetscPowRealInt(meter, 5) * PetscPowReal(Kelvin, 3.5) * PetscPowRealInt(kilogram, -2); /* */ 1083c4762a1bSJed Brown rd->K_p = 2.33e20 * PetscPowRealInt(meter, 5) * PetscPowReal(Kelvin, 3.5) * PetscPowRealInt(kilogram, -2); /* */ 1084c4762a1bSJed Brown rd->I_H = 2.179e-18 * Joule; /* Hydrogen ionization potential */ 1085c4762a1bSJed Brown rd->m_p = 1.673e-27 * kilogram; /* proton mass */ 1086c4762a1bSJed Brown rd->m_e = 9.109e-31 * kilogram; /* electron mass */ 1087c4762a1bSJed Brown rd->h = 6.626e-34 * Joule * second; /* Planck's constant */ 1088c4762a1bSJed Brown rd->k = 1.381e-23 * Joule / Kelvin; /* Boltzman constant */ 1089c4762a1bSJed Brown rd->c = 3.00e8 * meter / second; /* speed of light */ 1090c4762a1bSJed Brown rd->sigma_b = 5.67e-8 * Watt * PetscPowRealInt(meter, -2) * PetscPowRealInt(Kelvin, -4); /* Stefan-Boltzman constant */ 1091c4762a1bSJed Brown rd->MaterialEnergy = RDMaterialEnergy_Saha; 1092c4762a1bSJed Brown break; 1093c4762a1bSJed Brown } 1094c4762a1bSJed Brown 10959566063dSJacob Faibussowitsch PetscCall(DMDACreate1d(comm, DM_BOUNDARY_NONE, 20, sizeof(RDNode) / sizeof(PetscScalar), 1, NULL, &rd->da)); 10969566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(rd->da)); 10979566063dSJacob Faibussowitsch PetscCall(DMSetUp(rd->da)); 10989566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(rd->da, 0, "E")); 10999566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(rd->da, 1, "T")); 11009566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(rd->da, 0., 1., 0., 0., 0., 0.)); 1101c4762a1bSJed Brown 1102c4762a1bSJed Brown *inrd = rd; 1103*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1104c4762a1bSJed Brown } 1105c4762a1bSJed Brown 1106d71ae5a4SJacob Faibussowitsch int main(int argc, char *argv[]) 1107d71ae5a4SJacob Faibussowitsch { 1108c4762a1bSJed Brown RD rd; 1109c4762a1bSJed Brown TS ts; 1110c4762a1bSJed Brown SNES snes; 1111c4762a1bSJed Brown Vec X; 1112c4762a1bSJed Brown Mat A, B; 1113c4762a1bSJed Brown PetscInt steps; 1114c4762a1bSJed Brown PetscReal ftime; 1115c4762a1bSJed Brown 1116327415f7SBarry Smith PetscFunctionBeginUser; 11179ded082cSBarry Smith PetscCall(PetscInitialize(&argc, &argv, 0, NULL)); 11189566063dSJacob Faibussowitsch PetscCall(RDCreate(PETSC_COMM_WORLD, &rd)); 11199566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(rd->da, &X)); 11209566063dSJacob Faibussowitsch PetscCall(DMSetMatType(rd->da, MATAIJ)); 11219566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(rd->da, &B)); 11229566063dSJacob Faibussowitsch PetscCall(RDInitialState(rd, X)); 1123c4762a1bSJed Brown 11249566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 11259566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ts, TS_NONLINEAR)); 11269566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSTHETA)); 11279566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts, rd->da)); 1128c4762a1bSJed Brown switch (rd->discretization) { 1129c4762a1bSJed Brown case DISCRETIZATION_FD: 11309566063dSJacob Faibussowitsch PetscCall(TSSetIFunction(ts, NULL, RDIFunction_FD, rd)); 11319566063dSJacob Faibussowitsch if (rd->jacobian == JACOBIAN_ANALYTIC) PetscCall(TSSetIJacobian(ts, B, B, RDIJacobian_FD, rd)); 1132c4762a1bSJed Brown break; 1133c4762a1bSJed Brown case DISCRETIZATION_FE: 11349566063dSJacob Faibussowitsch PetscCall(TSSetIFunction(ts, NULL, RDIFunction_FE, rd)); 11359566063dSJacob Faibussowitsch if (rd->jacobian == JACOBIAN_ANALYTIC) PetscCall(TSSetIJacobian(ts, B, B, RDIJacobian_FE, rd)); 1136c4762a1bSJed Brown break; 1137c4762a1bSJed Brown } 11389566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, rd->final_time)); 11399566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, 1e-3)); 11409566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER)); 11419566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 1142c4762a1bSJed Brown 1143c4762a1bSJed Brown A = B; 11449566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &snes)); 1145c4762a1bSJed Brown switch (rd->jacobian) { 1146d71ae5a4SJacob Faibussowitsch case JACOBIAN_ANALYTIC: 1147d71ae5a4SJacob Faibussowitsch break; 1148d71ae5a4SJacob Faibussowitsch case JACOBIAN_MATRIXFREE: 1149d71ae5a4SJacob Faibussowitsch break; 1150c4762a1bSJed Brown case JACOBIAN_FD_COLORING: { 11519566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, A, B, SNESComputeJacobianDefaultColor, 0)); 1152c4762a1bSJed Brown } break; 1153d71ae5a4SJacob Faibussowitsch case JACOBIAN_FD_FULL: 1154d71ae5a4SJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, A, B, SNESComputeJacobianDefault, ts)); 1155d71ae5a4SJacob Faibussowitsch break; 1156c4762a1bSJed Brown } 1157c4762a1bSJed Brown 11581baa6e33SBarry Smith if (rd->test_diff) PetscCall(RDTestDifferentiation(rd)); 11599566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, X)); 11609566063dSJacob Faibussowitsch PetscCall(TSGetSolveTime(ts, &ftime)); 11619566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts, &steps)); 116263a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Steps %" PetscInt_FMT " final time %g\n", steps, (double)ftime)); 11631baa6e33SBarry Smith if (rd->view_draw) PetscCall(RDView(rd, X, PETSC_VIEWER_DRAW_WORLD)); 1164c4762a1bSJed Brown if (rd->view_binary[0]) { 1165c4762a1bSJed Brown PetscViewer viewer; 11669566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, rd->view_binary, FILE_MODE_WRITE, &viewer)); 11679566063dSJacob Faibussowitsch PetscCall(RDView(rd, X, viewer)); 11689566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 1169c4762a1bSJed Brown } 11709566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X)); 11719566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 11729566063dSJacob Faibussowitsch PetscCall(RDDestroy(&rd)); 11739566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 11749566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 1175b122ec5aSJacob Faibussowitsch return 0; 1176c4762a1bSJed Brown } 1177c4762a1bSJed Brown /*TEST 1178c4762a1bSJed Brown 1179c4762a1bSJed Brown test: 1180c4762a1bSJed Brown args: -da_grid_x 20 -rd_initial 1 -rd_discretization fd -rd_jacobian fd_coloring -rd_endpoint -ts_max_time 1 -ts_dt 2e-1 -ts_theta_initial_guess_extrapolate 0 -ts_monitor -snes_monitor_short -ksp_monitor_short 1181c4762a1bSJed Brown requires: !single 1182c4762a1bSJed Brown 1183c4762a1bSJed Brown test: 1184c4762a1bSJed Brown suffix: 2 1185c4762a1bSJed Brown args: -da_grid_x 20 -rd_initial 1 -rd_discretization fe -rd_quadrature lobatto2 -rd_jacobian fd_coloring -rd_endpoint -ts_max_time 1 -ts_dt 2e-1 -ts_theta_initial_guess_extrapolate 0 -ts_monitor -snes_monitor_short -ksp_monitor_short 1186c4762a1bSJed Brown requires: !single 1187c4762a1bSJed Brown 1188c4762a1bSJed Brown test: 1189c4762a1bSJed Brown suffix: 3 1190c4762a1bSJed Brown nsize: 2 1191c4762a1bSJed Brown args: -da_grid_x 20 -rd_initial 1 -rd_discretization fd -rd_jacobian analytic -rd_endpoint -ts_max_time 3 -ts_dt 1e-1 -ts_theta_initial_guess_extrapolate 0 -ts_monitor -snes_monitor_short -ksp_monitor_short 1192c4762a1bSJed Brown requires: !single 1193c4762a1bSJed Brown 1194c4762a1bSJed Brown TEST*/ 1195