xref: /petsc/src/ts/tutorials/ex10.c (revision 97a1619fdac2806951a9b35e227115e7d92f3e41)
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 
829371c9d4SSatish Balay static PetscErrorCode RDDestroy(RD *rd) {
83c4762a1bSJed Brown   PetscFunctionBeginUser;
849566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&(*rd)->da));
859566063dSJacob Faibussowitsch   PetscCall(PetscFree(*rd));
86c4762a1bSJed Brown   PetscFunctionReturn(0);
87c4762a1bSJed Brown }
88c4762a1bSJed Brown 
89c4762a1bSJed Brown /* The paper has a time derivative for material energy (Eq 2) which is a dependent variable (computable from temperature
90c4762a1bSJed Brown  * and density through an uninvertible relation).  Computing this derivative is trivial for trapezoid rule (used in the
91c4762a1bSJed Brown  * paper), but does not generalize nicely to higher order integrators.  Here we use the implicit form which provides
92c4762a1bSJed Brown  * time derivatives of the independent variables (radiation energy and temperature), so we must compute the time
93c4762a1bSJed Brown  * derivative of material energy ourselves (could be done using AD).
94c4762a1bSJed Brown  *
95c4762a1bSJed Brown  * There are multiple ionization models, this interface dispatches to the one currently in use.
96c4762a1bSJed Brown  */
979371c9d4SSatish Balay static void RDMaterialEnergy(RD rd, const RDNode *n, PetscScalar *Em, RDNode *dEm) {
989371c9d4SSatish Balay   rd->MaterialEnergy(rd, n, Em, dEm);
999371c9d4SSatish Balay }
100c4762a1bSJed Brown 
101c4762a1bSJed Brown /* Solves a quadratic equation while propagating tangents */
1029371c9d4SSatish Balay static void QuadraticSolve(PetscScalar a, PetscScalar a_t, PetscScalar b, PetscScalar b_t, PetscScalar c, PetscScalar c_t, PetscScalar *x, PetscScalar *x_t) {
1039371c9d4SSatish 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 */
1049371c9d4SSatish Balay     num_t = -b_t + 0.5 / PetscSqrtScalar(disc) * disc_t, den = 2. * a, den_t = 2. * a_t;
105c4762a1bSJed Brown   *x   = num / den;
106c4762a1bSJed Brown   *x_t = (num_t * den - num * den_t) / PetscSqr(den);
107c4762a1bSJed Brown }
108c4762a1bSJed Brown 
109c4762a1bSJed Brown /* The primary model presented in the paper */
1109371c9d4SSatish Balay static void RDMaterialEnergy_Saha(RD rd, const RDNode *n, PetscScalar *inEm, RDNode *dEm) {
1119371c9d4SSatish 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 */
1129371c9d4SSatish Balay     b_t = -b * chi_t + 1.5 * b / chi * chi_t, c = -b, c_t = -b_t;
113c4762a1bSJed Brown   QuadraticSolve(a, a_t, b, b_t, c, c_t, &alpha, &alpha_t);      /* Solve Eq 7 for alpha */
114c4762a1bSJed Brown   Em = rd->k * T / rd->m_p * (1.5 * (1. + alpha) + alpha * chi); /* Eq 6 */
115c4762a1bSJed Brown   if (inEm) *inEm = Em;
116c4762a1bSJed Brown   if (dEm) {
117c4762a1bSJed Brown     dEm->E = 0;
118c4762a1bSJed Brown     dEm->T = Em / T * T_t + rd->k * T / rd->m_p * (1.5 * alpha_t + alpha_t * chi + alpha * chi_t);
119c4762a1bSJed Brown   }
120c4762a1bSJed Brown }
121c4762a1bSJed Brown /* Reduced ionization model, Eq 30 */
1229371c9d4SSatish Balay static void RDMaterialEnergy_Reduced(RD rd, const RDNode *n, PetscScalar *Em, RDNode *dEm) {
1239371c9d4SSatish 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;
124c4762a1bSJed Brown   QuadraticSolve(a, a_t, b, b_t, c, c_t, &alpha, &alpha_t);
125c4762a1bSJed Brown   if (Em) *Em = (1. + alpha) * T + 0.3 * alpha;
126c4762a1bSJed Brown   if (dEm) {
127c4762a1bSJed Brown     dEm->E = 0;
128c4762a1bSJed Brown     dEm->T = alpha_t * T + (1. + alpha) * T_t + 0.3 * alpha_t;
129c4762a1bSJed Brown   }
130c4762a1bSJed Brown }
131c4762a1bSJed Brown 
132c4762a1bSJed Brown /* Eq 5 */
1339371c9d4SSatish Balay static void RDSigma_R(RD rd, RDNode *n, PetscScalar *sigma_R, RDNode *dsigma_R) {
134c4762a1bSJed Brown   *sigma_R    = rd->K_R * rd->rho * PetscPowScalar(n->T, -rd->gamma);
135c4762a1bSJed Brown   dsigma_R->E = 0;
136c4762a1bSJed Brown   dsigma_R->T = -rd->gamma * (*sigma_R) / n->T;
137c4762a1bSJed Brown }
138c4762a1bSJed Brown 
139c4762a1bSJed Brown /* Eq 4 */
1409371c9d4SSatish Balay static void RDDiffusionCoefficient(RD rd, PetscBool limit, RDNode *n, RDNode *nx, PetscScalar *D_R, RDNode *dD_R, RDNode *dxD_R) {
141c4762a1bSJed Brown   PetscScalar sigma_R, denom;
142c4762a1bSJed Brown   RDNode      dsigma_R, ddenom, dxdenom;
143c4762a1bSJed Brown 
144c4762a1bSJed Brown   RDSigma_R(rd, n, &sigma_R, &dsigma_R);
145c4762a1bSJed Brown   denom     = 3. * rd->rho * sigma_R + (int)limit * PetscAbsScalar(nx->E) / n->E;
146c4762a1bSJed Brown   ddenom.E  = -(int)limit * PetscAbsScalar(nx->E) / PetscSqr(n->E);
147c4762a1bSJed Brown   ddenom.T  = 3. * rd->rho * dsigma_R.T;
148c4762a1bSJed Brown   dxdenom.E = (int)limit * (PetscRealPart(nx->E) < 0 ? -1. : 1.) / n->E;
149c4762a1bSJed Brown   dxdenom.T = 0;
150c4762a1bSJed Brown   *D_R      = rd->c / denom;
151c4762a1bSJed Brown   if (dD_R) {
152c4762a1bSJed Brown     dD_R->E = -rd->c / PetscSqr(denom) * ddenom.E;
153c4762a1bSJed Brown     dD_R->T = -rd->c / PetscSqr(denom) * ddenom.T;
154c4762a1bSJed Brown   }
155c4762a1bSJed Brown   if (dxD_R) {
156c4762a1bSJed Brown     dxD_R->E = -rd->c / PetscSqr(denom) * dxdenom.E;
157c4762a1bSJed Brown     dxD_R->T = -rd->c / PetscSqr(denom) * dxdenom.T;
158c4762a1bSJed Brown   }
159c4762a1bSJed Brown }
160c4762a1bSJed Brown 
1619371c9d4SSatish Balay static PetscErrorCode RDStateView(RD rd, Vec X, Vec Xdot, Vec F) {
162c4762a1bSJed Brown   DMDALocalInfo info;
163c4762a1bSJed Brown   PetscInt      i;
164c4762a1bSJed Brown   const RDNode *x, *xdot, *f;
165c4762a1bSJed Brown   MPI_Comm      comm;
166c4762a1bSJed Brown 
167c4762a1bSJed Brown   PetscFunctionBeginUser;
1689566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)rd->da, &comm));
1699566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(rd->da, &info));
1709566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(rd->da, X, (void *)&x));
1719566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(rd->da, Xdot, (void *)&xdot));
1729566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(rd->da, F, (void *)&f));
173c4762a1bSJed Brown   for (i = info.xs; i < info.xs + info.xm; i++) {
1749371c9d4SSatish 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),
1759371c9d4SSatish Balay                                       (double)PetscRealPart(f[i].E), (double)PetscRealPart(f[i].T)));
176c4762a1bSJed Brown   }
1779566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(rd->da, X, (void *)&x));
1789566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(rd->da, Xdot, (void *)&xdot));
1799566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(rd->da, F, (void *)&f));
1809566063dSJacob Faibussowitsch   PetscCall(PetscSynchronizedFlush(comm, PETSC_STDOUT));
181c4762a1bSJed Brown   PetscFunctionReturn(0);
182c4762a1bSJed Brown }
183c4762a1bSJed Brown 
1849371c9d4SSatish Balay static PetscScalar RDRadiation(RD rd, const RDNode *n, RDNode *dn) {
1859371c9d4SSatish 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);
186c4762a1bSJed Brown   if (dn) {
187c4762a1bSJed Brown     dn->E = rad_E;
188c4762a1bSJed Brown     dn->T = rad_T;
189c4762a1bSJed Brown   }
190c4762a1bSJed Brown   return rad;
191c4762a1bSJed Brown }
192c4762a1bSJed Brown 
1939371c9d4SSatish Balay static PetscScalar RDDiffusion(RD rd, PetscReal hx, const RDNode x[], PetscInt i, RDNode d[]) {
194c4762a1bSJed Brown   PetscReal   ihx = 1. / hx;
195c4762a1bSJed Brown   RDNode      n_L, nx_L, n_R, nx_R, dD_L, dxD_L, dD_R, dxD_R, dfluxL[2], dfluxR[2];
196c4762a1bSJed Brown   PetscScalar D_L, D_R, fluxL, fluxR;
197c4762a1bSJed Brown 
198c4762a1bSJed Brown   n_L.E  = 0.5 * (x[i - 1].E + x[i].E);
199c4762a1bSJed Brown   n_L.T  = 0.5 * (x[i - 1].T + x[i].T);
200c4762a1bSJed Brown   nx_L.E = (x[i].E - x[i - 1].E) / hx;
201c4762a1bSJed Brown   nx_L.T = (x[i].T - x[i - 1].T) / hx;
202c4762a1bSJed Brown   RDDiffusionCoefficient(rd, PETSC_TRUE, &n_L, &nx_L, &D_L, &dD_L, &dxD_L);
203c4762a1bSJed Brown   fluxL       = D_L * nx_L.E;
204c4762a1bSJed Brown   dfluxL[0].E = -ihx * D_L + (0.5 * dD_L.E - ihx * dxD_L.E) * nx_L.E;
205c4762a1bSJed Brown   dfluxL[1].E = +ihx * D_L + (0.5 * dD_L.E + ihx * dxD_L.E) * nx_L.E;
206c4762a1bSJed Brown   dfluxL[0].T = (0.5 * dD_L.T - ihx * dxD_L.T) * nx_L.E;
207c4762a1bSJed Brown   dfluxL[1].T = (0.5 * dD_L.T + ihx * dxD_L.T) * nx_L.E;
208c4762a1bSJed Brown 
209c4762a1bSJed Brown   n_R.E  = 0.5 * (x[i].E + x[i + 1].E);
210c4762a1bSJed Brown   n_R.T  = 0.5 * (x[i].T + x[i + 1].T);
211c4762a1bSJed Brown   nx_R.E = (x[i + 1].E - x[i].E) / hx;
212c4762a1bSJed Brown   nx_R.T = (x[i + 1].T - x[i].T) / hx;
213c4762a1bSJed Brown   RDDiffusionCoefficient(rd, PETSC_TRUE, &n_R, &nx_R, &D_R, &dD_R, &dxD_R);
214c4762a1bSJed Brown   fluxR       = D_R * nx_R.E;
215c4762a1bSJed Brown   dfluxR[0].E = -ihx * D_R + (0.5 * dD_R.E - ihx * dxD_R.E) * nx_R.E;
216c4762a1bSJed Brown   dfluxR[1].E = +ihx * D_R + (0.5 * dD_R.E + ihx * dxD_R.E) * nx_R.E;
217c4762a1bSJed Brown   dfluxR[0].T = (0.5 * dD_R.T - ihx * dxD_R.T) * nx_R.E;
218c4762a1bSJed Brown   dfluxR[1].T = (0.5 * dD_R.T + ihx * dxD_R.T) * nx_R.E;
219c4762a1bSJed Brown 
220c4762a1bSJed Brown   if (d) {
221c4762a1bSJed Brown     d[0].E = -ihx * dfluxL[0].E;
222c4762a1bSJed Brown     d[0].T = -ihx * dfluxL[0].T;
223c4762a1bSJed Brown     d[1].E = ihx * (dfluxR[0].E - dfluxL[1].E);
224c4762a1bSJed Brown     d[1].T = ihx * (dfluxR[0].T - dfluxL[1].T);
225c4762a1bSJed Brown     d[2].E = ihx * dfluxR[1].E;
226c4762a1bSJed Brown     d[2].T = ihx * dfluxR[1].T;
227c4762a1bSJed Brown   }
228c4762a1bSJed Brown   return ihx * (fluxR - fluxL);
229c4762a1bSJed Brown }
230c4762a1bSJed Brown 
2319371c9d4SSatish Balay 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) {
232c4762a1bSJed Brown   PetscBool istheta;
233c4762a1bSJed Brown 
234c4762a1bSJed Brown   PetscFunctionBeginUser;
2359566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(rd->da, X0loc));
2369566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(rd->da, Xloc));
2379566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(rd->da, Xloc_t));
238c4762a1bSJed Brown 
2399566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(rd->da, X, INSERT_VALUES, *Xloc));
2409566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(rd->da, X, INSERT_VALUES, *Xloc));
2419566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(rd->da, Xdot, INSERT_VALUES, *Xloc_t));
2429566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(rd->da, Xdot, INSERT_VALUES, *Xloc_t));
243c4762a1bSJed Brown 
244c4762a1bSJed Brown   /*
245c4762a1bSJed Brown     The following is a hack to subvert TSTHETA which is like an implicit midpoint method to behave more like a trapezoid
246c4762a1bSJed Brown     rule.  These methods have equivalent linear stability, but the nonlinear stability is somewhat different.  The
247c4762a1bSJed Brown     radiation system is inconvenient to write in explicit form because the ionization model is "on the left".
248c4762a1bSJed Brown    */
2499566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)ts, TSTHETA, &istheta));
250c4762a1bSJed Brown   if (istheta && rd->endpoint) {
2519566063dSJacob Faibussowitsch     PetscCall(TSThetaGetTheta(ts, Theta));
252c4762a1bSJed Brown   } else *Theta = 1.;
253c4762a1bSJed Brown 
2549566063dSJacob Faibussowitsch   PetscCall(TSGetTimeStep(ts, dt));
2559566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(*X0loc, -(*Theta) * (*dt), *Xloc_t, *Xloc)); /* back out the value at the start of this step */
2569371c9d4SSatish Balay   if (rd->endpoint) { PetscCall(VecWAXPY(*Xloc, *dt, *Xloc_t, *X0loc)); /* move the abscissa to the end of the step */ }
257c4762a1bSJed Brown 
2589566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(rd->da, *X0loc, x0));
2599566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(rd->da, *Xloc, x));
2609566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(rd->da, *Xloc_t, xdot));
261c4762a1bSJed Brown   PetscFunctionReturn(0);
262c4762a1bSJed Brown }
263c4762a1bSJed Brown 
2649371c9d4SSatish Balay static PetscErrorCode RDRestoreLocalArrays(RD rd, Vec *X0loc, RDNode **x0, Vec *Xloc, RDNode **x, Vec *Xloc_t, RDNode **xdot) {
265c4762a1bSJed Brown   PetscFunctionBeginUser;
2669566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(rd->da, *X0loc, x0));
2679566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(rd->da, *Xloc, x));
2689566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(rd->da, *Xloc_t, xdot));
2699566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(rd->da, X0loc));
2709566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(rd->da, Xloc));
2719566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(rd->da, Xloc_t));
272c4762a1bSJed Brown   PetscFunctionReturn(0);
273c4762a1bSJed Brown }
274c4762a1bSJed Brown 
2759371c9d4SSatish Balay static PetscErrorCode PETSC_UNUSED RDCheckDomain_Private(RD rd, TS ts, Vec X, PetscBool *in) {
276c4762a1bSJed Brown   PetscInt  minloc;
277c4762a1bSJed Brown   PetscReal min;
278c4762a1bSJed Brown 
279c4762a1bSJed Brown   PetscFunctionBeginUser;
2809566063dSJacob Faibussowitsch   PetscCall(VecMin(X, &minloc, &min));
281c4762a1bSJed Brown   if (min < 0) {
282c4762a1bSJed Brown     SNES snes;
283c4762a1bSJed Brown     *in = PETSC_FALSE;
2849566063dSJacob Faibussowitsch     PetscCall(TSGetSNES(ts, &snes));
2859566063dSJacob Faibussowitsch     PetscCall(SNESSetFunctionDomainError(snes));
28663a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(ts, "Domain violation at %" PetscInt_FMT " field %" PetscInt_FMT " value %g\n", minloc / 2, minloc % 2, (double)min));
287c4762a1bSJed Brown   } else *in = PETSC_TRUE;
288c4762a1bSJed Brown   PetscFunctionReturn(0);
289c4762a1bSJed Brown }
290c4762a1bSJed Brown 
291c4762a1bSJed Brown /* Energy and temperature must remain positive */
2929371c9d4SSatish Balay #define RDCheckDomain(rd, ts, X) \
2939371c9d4SSatish Balay   do { \
294c4762a1bSJed Brown     PetscBool _in; \
2959566063dSJacob Faibussowitsch     PetscCall(RDCheckDomain_Private(rd, ts, X, &_in)); \
296c4762a1bSJed Brown     if (!_in) PetscFunctionReturn(0); \
297c4762a1bSJed Brown   } while (0)
298c4762a1bSJed Brown 
2999371c9d4SSatish Balay static PetscErrorCode RDIFunction_FD(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, void *ctx) {
300c4762a1bSJed Brown   RD            rd = (RD)ctx;
301c4762a1bSJed Brown   RDNode       *x, *x0, *xdot, *f;
302c4762a1bSJed Brown   Vec           X0loc, Xloc, Xloc_t;
303c4762a1bSJed Brown   PetscReal     hx, Theta, dt;
304c4762a1bSJed Brown   DMDALocalInfo info;
305c4762a1bSJed Brown   PetscInt      i;
306c4762a1bSJed Brown 
307c4762a1bSJed Brown   PetscFunctionBeginUser;
3089566063dSJacob Faibussowitsch   PetscCall(RDGetLocalArrays(rd, ts, X, Xdot, &Theta, &dt, &X0loc, &x0, &Xloc, &x, &Xloc_t, &xdot));
3099566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(rd->da, F, &f));
3109566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(rd->da, &info));
311c4762a1bSJed Brown   VecZeroEntries(F);
312c4762a1bSJed Brown 
313c4762a1bSJed Brown   hx = rd->L / (info.mx - 1);
314c4762a1bSJed Brown 
315c4762a1bSJed Brown   for (i = info.xs; i < info.xs + info.xm; i++) {
316c4762a1bSJed Brown     PetscReal   rho = rd->rho;
317c4762a1bSJed Brown     PetscScalar Em_t, rad;
318c4762a1bSJed Brown 
319c4762a1bSJed Brown     rad = (1. - Theta) * RDRadiation(rd, &x0[i], 0) + Theta * RDRadiation(rd, &x[i], 0);
320c4762a1bSJed Brown     if (rd->endpoint) {
321c4762a1bSJed Brown       PetscScalar Em0, Em1;
322c4762a1bSJed Brown       RDMaterialEnergy(rd, &x0[i], &Em0, NULL);
323c4762a1bSJed Brown       RDMaterialEnergy(rd, &x[i], &Em1, NULL);
324c4762a1bSJed Brown       Em_t = (Em1 - Em0) / dt;
325c4762a1bSJed Brown     } else {
326c4762a1bSJed Brown       RDNode dEm;
327c4762a1bSJed Brown       RDMaterialEnergy(rd, &x[i], NULL, &dEm);
328c4762a1bSJed Brown       Em_t = dEm.E * xdot[i].E + dEm.T * xdot[i].T;
329c4762a1bSJed Brown     }
330c4762a1bSJed Brown     /* Residuals are multiplied by the volume element (hx).  */
331c4762a1bSJed Brown     /* The temperature equation does not have boundary conditions */
332c4762a1bSJed Brown     f[i].T = hx * (rho * Em_t + rad);
333c4762a1bSJed Brown 
334c4762a1bSJed Brown     if (i == 0) { /* Left boundary condition */
335c4762a1bSJed Brown       PetscScalar D_R, bcTheta = rd->bcmidpoint ? Theta : 1.;
336c4762a1bSJed Brown       RDNode      n, nx;
337c4762a1bSJed Brown 
338c4762a1bSJed Brown       n.E  = (1. - bcTheta) * x0[0].E + bcTheta * x[0].E;
339c4762a1bSJed Brown       n.T  = (1. - bcTheta) * x0[0].T + bcTheta * x[0].T;
340c4762a1bSJed Brown       nx.E = ((1. - bcTheta) * (x0[1].E - x0[0].E) + bcTheta * (x[1].E - x[0].E)) / hx;
341c4762a1bSJed Brown       nx.T = ((1. - bcTheta) * (x0[1].T - x0[0].T) + bcTheta * (x[1].T - x[0].T)) / hx;
342c4762a1bSJed Brown       switch (rd->leftbc) {
343c4762a1bSJed Brown       case BC_ROBIN:
344c4762a1bSJed Brown         RDDiffusionCoefficient(rd, rd->bclimit, &n, &nx, &D_R, 0, 0);
345c4762a1bSJed Brown         f[0].E = hx * (n.E - 2. * D_R * nx.E - rd->Eapplied);
346c4762a1bSJed Brown         break;
3479371c9d4SSatish Balay       case BC_NEUMANN: f[0].E = x[1].E - x[0].E; break;
34863a3b9bcSJacob Faibussowitsch       default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Case %" PetscInt_FMT, rd->initial);
349c4762a1bSJed Brown       }
350c4762a1bSJed Brown     } else if (i == info.mx - 1) {  /* Right boundary */
351c4762a1bSJed Brown       f[i].E = x[i].E - x[i - 1].E; /* Homogeneous Neumann */
352c4762a1bSJed Brown     } else {
353c4762a1bSJed Brown       PetscScalar diff = (1. - Theta) * RDDiffusion(rd, hx, x0, i, 0) + Theta * RDDiffusion(rd, hx, x, i, 0);
354c4762a1bSJed Brown       f[i].E           = hx * (xdot[i].E - diff - rad);
355c4762a1bSJed Brown     }
356c4762a1bSJed Brown   }
3579566063dSJacob Faibussowitsch   PetscCall(RDRestoreLocalArrays(rd, &X0loc, &x0, &Xloc, &x, &Xloc_t, &xdot));
3589566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(rd->da, F, &f));
3599566063dSJacob Faibussowitsch   if (rd->monitor_residual) PetscCall(RDStateView(rd, X, Xdot, F));
360c4762a1bSJed Brown   PetscFunctionReturn(0);
361c4762a1bSJed Brown }
362c4762a1bSJed Brown 
3639371c9d4SSatish Balay static PetscErrorCode RDIJacobian_FD(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal a, Mat A, Mat B, void *ctx) {
364c4762a1bSJed Brown   RD            rd = (RD)ctx;
365c4762a1bSJed Brown   RDNode       *x, *x0, *xdot;
366c4762a1bSJed Brown   Vec           X0loc, Xloc, Xloc_t;
367c4762a1bSJed Brown   PetscReal     hx, Theta, dt;
368c4762a1bSJed Brown   DMDALocalInfo info;
369c4762a1bSJed Brown   PetscInt      i;
370c4762a1bSJed Brown 
371c4762a1bSJed Brown   PetscFunctionBeginUser;
3729566063dSJacob Faibussowitsch   PetscCall(RDGetLocalArrays(rd, ts, X, Xdot, &Theta, &dt, &X0loc, &x0, &Xloc, &x, &Xloc_t, &xdot));
3739566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(rd->da, &info));
374c4762a1bSJed Brown   hx = rd->L / (info.mx - 1);
3759566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(B));
376c4762a1bSJed Brown 
377c4762a1bSJed Brown   for (i = info.xs; i < info.xs + info.xm; i++) {
378c4762a1bSJed Brown     PetscInt                  col[3];
379c4762a1bSJed Brown     PetscReal                 rho = rd->rho;
380c4762a1bSJed Brown     PetscScalar /*Em_t,rad,*/ K[2][6];
381c4762a1bSJed Brown     RDNode                    dEm_t, drad;
382c4762a1bSJed Brown 
3839371c9d4SSatish Balay     /*rad = (1.-Theta)* */ RDRadiation(rd, &x0[i], 0); /* + Theta* */
3849371c9d4SSatish Balay     RDRadiation(rd, &x[i], &drad);
385c4762a1bSJed Brown 
386c4762a1bSJed Brown     if (rd->endpoint) {
387c4762a1bSJed Brown       PetscScalar Em0, Em1;
388c4762a1bSJed Brown       RDNode      dEm1;
389c4762a1bSJed Brown       RDMaterialEnergy(rd, &x0[i], &Em0, NULL);
390c4762a1bSJed Brown       RDMaterialEnergy(rd, &x[i], &Em1, &dEm1);
391c4762a1bSJed Brown       /*Em_t = (Em1 - Em0) / (Theta*dt);*/
392c4762a1bSJed Brown       dEm_t.E = dEm1.E / (Theta * dt);
393c4762a1bSJed Brown       dEm_t.T = dEm1.T / (Theta * dt);
394c4762a1bSJed Brown     } else {
395c4762a1bSJed Brown       const PetscScalar epsilon = x[i].T * PETSC_SQRT_MACHINE_EPSILON;
396c4762a1bSJed Brown       RDNode            n1;
397c4762a1bSJed Brown       RDNode            dEm, dEm1;
398c4762a1bSJed Brown       PetscScalar       Em_TT;
399c4762a1bSJed Brown 
400c4762a1bSJed Brown       n1.E = x[i].E;
401c4762a1bSJed Brown       n1.T = x[i].T + epsilon;
402c4762a1bSJed Brown       RDMaterialEnergy(rd, &x[i], NULL, &dEm);
403c4762a1bSJed Brown       RDMaterialEnergy(rd, &n1, NULL, &dEm1);
404c4762a1bSJed Brown       /* The Jacobian needs another derivative.  We finite difference here instead of
405c4762a1bSJed Brown        * propagating second derivatives through the ionization model. */
406c4762a1bSJed Brown       Em_TT   = (dEm1.T - dEm.T) / epsilon;
407c4762a1bSJed Brown       /*Em_t = dEm.E * xdot[i].E + dEm.T * xdot[i].T;*/
408c4762a1bSJed Brown       dEm_t.E = dEm.E * a;
409c4762a1bSJed Brown       dEm_t.T = dEm.T * a + Em_TT * xdot[i].T;
410c4762a1bSJed Brown     }
411c4762a1bSJed Brown 
4129566063dSJacob Faibussowitsch     PetscCall(PetscMemzero(K, sizeof(K)));
413c4762a1bSJed Brown     /* Residuals are multiplied by the volume element (hx).  */
414c4762a1bSJed Brown     if (i == 0) {
415c4762a1bSJed Brown       PetscScalar D, bcTheta = rd->bcmidpoint ? Theta : 1.;
416c4762a1bSJed Brown       RDNode      n, nx;
417c4762a1bSJed Brown       RDNode      dD, dxD;
418c4762a1bSJed Brown 
419c4762a1bSJed Brown       n.E  = (1. - bcTheta) * x0[0].E + bcTheta * x[0].E;
420c4762a1bSJed Brown       n.T  = (1. - bcTheta) * x0[0].T + bcTheta * x[0].T;
421c4762a1bSJed Brown       nx.E = ((1. - bcTheta) * (x0[1].E - x0[0].E) + bcTheta * (x[1].E - x[0].E)) / hx;
422c4762a1bSJed Brown       nx.T = ((1. - bcTheta) * (x0[1].T - x0[0].T) + bcTheta * (x[1].T - x[0].T)) / hx;
423c4762a1bSJed Brown       switch (rd->leftbc) {
424c4762a1bSJed Brown       case BC_ROBIN:
425c4762a1bSJed Brown         RDDiffusionCoefficient(rd, rd->bclimit, &n, &nx, &D, &dD, &dxD);
426c4762a1bSJed 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);
427c4762a1bSJed Brown         K[0][1 * 2 + 1] = (bcTheta / Theta) * hx * (-2. * nx.E * dD.T);
428c4762a1bSJed Brown         K[0][2 * 2 + 0] = (bcTheta / Theta) * hx * (-2. * D * (1. / hx) - 2. * nx.E * dD.E - 2. * nx.E * dxD.E / hx);
429c4762a1bSJed Brown         break;
430c4762a1bSJed Brown       case BC_NEUMANN:
431c4762a1bSJed Brown         K[0][1 * 2 + 0] = -1. / Theta;
432c4762a1bSJed Brown         K[0][2 * 2 + 0] = 1. / Theta;
433c4762a1bSJed Brown         break;
43463a3b9bcSJacob Faibussowitsch       default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Case %" PetscInt_FMT, rd->initial);
435c4762a1bSJed Brown       }
436c4762a1bSJed Brown     } else if (i == info.mx - 1) {
437c4762a1bSJed Brown       K[0][0 * 2 + 0] = -1. / Theta;
438c4762a1bSJed Brown       K[0][1 * 2 + 0] = 1. / Theta;
439c4762a1bSJed Brown     } else {
440c4762a1bSJed Brown       /*PetscScalar diff;*/
441c4762a1bSJed Brown       RDNode ddiff[3];
442c4762a1bSJed Brown       /*diff = (1.-Theta)*RDDiffusion(rd,hx,x0,i,0) + Theta* */ RDDiffusion(rd, hx, x, i, ddiff);
443c4762a1bSJed Brown       K[0][0 * 2 + 0] = -hx * ddiff[0].E;
444c4762a1bSJed Brown       K[0][0 * 2 + 1] = -hx * ddiff[0].T;
445c4762a1bSJed Brown       K[0][1 * 2 + 0] = hx * (a - ddiff[1].E - drad.E);
446c4762a1bSJed Brown       K[0][1 * 2 + 1] = hx * (-ddiff[1].T - drad.T);
447c4762a1bSJed Brown       K[0][2 * 2 + 0] = -hx * ddiff[2].E;
448c4762a1bSJed Brown       K[0][2 * 2 + 1] = -hx * ddiff[2].T;
449c4762a1bSJed Brown     }
450c4762a1bSJed Brown 
451c4762a1bSJed Brown     K[1][1 * 2 + 0] = hx * (rho * dEm_t.E + drad.E);
452c4762a1bSJed Brown     K[1][1 * 2 + 1] = hx * (rho * dEm_t.T + drad.T);
453c4762a1bSJed Brown 
454c4762a1bSJed Brown     col[0] = i - 1;
455c4762a1bSJed Brown     col[1] = i;
456c4762a1bSJed Brown     col[2] = i + 1 < info.mx ? i + 1 : -1;
4579566063dSJacob Faibussowitsch     PetscCall(MatSetValuesBlocked(B, 1, &i, 3, col, &K[0][0], INSERT_VALUES));
458c4762a1bSJed Brown   }
4599566063dSJacob Faibussowitsch   PetscCall(RDRestoreLocalArrays(rd, &X0loc, &x0, &Xloc, &x, &Xloc_t, &xdot));
4609566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
4619566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
462c4762a1bSJed Brown   if (A != B) {
4639566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
4649566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
465c4762a1bSJed Brown   }
466c4762a1bSJed Brown   PetscFunctionReturn(0);
467c4762a1bSJed Brown }
468c4762a1bSJed Brown 
469c4762a1bSJed Brown /* Evaluate interpolants and derivatives at a select quadrature point */
4709371c9d4SSatish Balay static void RDEvaluate(PetscReal interp[][2], PetscReal deriv[][2], PetscInt q, const RDNode x[], PetscInt i, RDNode *n, RDNode *nx) {
471c4762a1bSJed Brown   PetscInt j;
4729371c9d4SSatish Balay   n->E  = 0;
4739371c9d4SSatish Balay   n->T  = 0;
4749371c9d4SSatish Balay   nx->E = 0;
4759371c9d4SSatish Balay   nx->T = 0;
476c4762a1bSJed Brown   for (j = 0; j < 2; j++) {
477c4762a1bSJed Brown     n->E += interp[q][j] * x[i + j].E;
478c4762a1bSJed Brown     n->T += interp[q][j] * x[i + j].T;
479c4762a1bSJed Brown     nx->E += deriv[q][j] * x[i + j].E;
480c4762a1bSJed Brown     nx->T += deriv[q][j] * x[i + j].T;
481c4762a1bSJed Brown   }
482c4762a1bSJed Brown }
483c4762a1bSJed Brown 
484c4762a1bSJed Brown /*
485c4762a1bSJed Brown  Various quadrature rules.  The nonlinear terms are non-polynomial so no standard quadrature will be exact.
486c4762a1bSJed Brown */
4879371c9d4SSatish Balay static PetscErrorCode RDGetQuadrature(RD rd, PetscReal hx, PetscInt *nq, PetscReal weight[], PetscReal interp[][2], PetscReal deriv[][2]) {
488c4762a1bSJed Brown   PetscInt         q, j;
489c4762a1bSJed Brown   const PetscReal *refweight, (*refinterp)[2], (*refderiv)[2];
490c4762a1bSJed Brown 
491c4762a1bSJed Brown   PetscFunctionBeginUser;
492c4762a1bSJed Brown   switch (rd->quadrature) {
493c4762a1bSJed Brown   case QUADRATURE_GAUSS1: {
494*97a1619fSSatish Balay     static const PetscReal ww[1]    = {1.};
495*97a1619fSSatish Balay     static const PetscReal ii[1][2] = {
496*97a1619fSSatish Balay       {0.5, 0.5}
497*97a1619fSSatish Balay     };
498*97a1619fSSatish Balay     static const PetscReal dd[1][2] = {
499*97a1619fSSatish Balay       {-1., 1.}
500*97a1619fSSatish Balay     };
5019371c9d4SSatish Balay     *nq       = 1;
5029371c9d4SSatish Balay     refweight = ww;
5039371c9d4SSatish Balay     refinterp = ii;
5049371c9d4SSatish Balay     refderiv  = dd;
505c4762a1bSJed Brown   } break;
506c4762a1bSJed Brown   case QUADRATURE_GAUSS2: {
507*97a1619fSSatish Balay     static const PetscReal ii[2][2] = {
5089371c9d4SSatish Balay       {0.78867513459481287, 0.21132486540518713},
5099371c9d4SSatish Balay       {0.21132486540518713, 0.78867513459481287}
510*97a1619fSSatish Balay     };
511*97a1619fSSatish Balay     static const PetscReal dd[2][2] = {
512*97a1619fSSatish Balay       {-1., 1.},
513*97a1619fSSatish Balay       {-1., 1.}
514*97a1619fSSatish Balay     };
515*97a1619fSSatish Balay     static const PetscReal ww[2] = {0.5, 0.5};
5169371c9d4SSatish Balay     *nq                          = 2;
5179371c9d4SSatish Balay     refweight                    = ww;
5189371c9d4SSatish Balay     refinterp                    = ii;
5199371c9d4SSatish Balay     refderiv                     = dd;
520c4762a1bSJed Brown   } break;
521c4762a1bSJed Brown   case QUADRATURE_GAUSS3: {
522*97a1619fSSatish Balay     static const PetscReal ii[3][2] = {
5239371c9d4SSatish Balay       {0.8872983346207417, 0.1127016653792583},
5249371c9d4SSatish Balay       {0.5,                0.5               },
5259371c9d4SSatish Balay       {0.1127016653792583, 0.8872983346207417}
526*97a1619fSSatish Balay     };
527*97a1619fSSatish Balay     static const PetscReal dd[3][2] = {
528*97a1619fSSatish Balay       {-1, 1},
529*97a1619fSSatish Balay       {-1, 1},
530*97a1619fSSatish Balay       {-1, 1}
531*97a1619fSSatish Balay     };
532*97a1619fSSatish Balay     static const PetscReal ww[3] = {5. / 18, 8. / 18, 5. / 18};
5339371c9d4SSatish Balay     *nq                          = 3;
5349371c9d4SSatish Balay     refweight                    = ww;
5359371c9d4SSatish Balay     refinterp                    = ii;
5369371c9d4SSatish Balay     refderiv                     = dd;
537c4762a1bSJed Brown   } break;
538c4762a1bSJed Brown   case QUADRATURE_GAUSS4: {
539*97a1619fSSatish Balay     static const PetscReal ii[][2] = {
5409371c9d4SSatish Balay       {0.93056815579702623,  0.069431844202973658},
541c4762a1bSJed Brown       {0.66999052179242813,  0.33000947820757187 },
542c4762a1bSJed Brown       {0.33000947820757187,  0.66999052179242813 },
5439371c9d4SSatish Balay       {0.069431844202973658, 0.93056815579702623 }
544*97a1619fSSatish Balay     };
545*97a1619fSSatish Balay     static const PetscReal dd[][2] = {
546*97a1619fSSatish Balay       {-1, 1},
547*97a1619fSSatish Balay       {-1, 1},
548*97a1619fSSatish Balay       {-1, 1},
549*97a1619fSSatish Balay       {-1, 1}
550*97a1619fSSatish Balay     };
551*97a1619fSSatish Balay     static const PetscReal ww[] = {0.17392742256872692, 0.3260725774312731, 0.3260725774312731, 0.17392742256872692};
552c4762a1bSJed Brown 
5539371c9d4SSatish Balay     *nq       = 4;
5549371c9d4SSatish Balay     refweight = ww;
5559371c9d4SSatish Balay     refinterp = ii;
5569371c9d4SSatish Balay     refderiv  = dd;
557c4762a1bSJed Brown   } break;
558c4762a1bSJed Brown   case QUADRATURE_LOBATTO2: {
559*97a1619fSSatish Balay     static const PetscReal ii[2][2] = {
5609371c9d4SSatish Balay       {1., 0.},
5619371c9d4SSatish Balay       {0., 1.}
562*97a1619fSSatish Balay     };
563*97a1619fSSatish Balay     static const PetscReal dd[2][2] = {
564*97a1619fSSatish Balay       {-1., 1.},
565*97a1619fSSatish Balay       {-1., 1.}
566*97a1619fSSatish Balay     };
567*97a1619fSSatish Balay     static const PetscReal ww[2] = {0.5, 0.5};
5689371c9d4SSatish Balay     *nq                          = 2;
5699371c9d4SSatish Balay     refweight                    = ww;
5709371c9d4SSatish Balay     refinterp                    = ii;
5719371c9d4SSatish Balay     refderiv                     = dd;
572c4762a1bSJed Brown   } break;
573c4762a1bSJed Brown   case QUADRATURE_LOBATTO3: {
574*97a1619fSSatish Balay     static const PetscReal ii[3][2] = {
5759371c9d4SSatish Balay       {1,   0  },
5769371c9d4SSatish Balay       {0.5, 0.5},
5779371c9d4SSatish Balay       {0,   1  }
578*97a1619fSSatish Balay     };
579*97a1619fSSatish Balay     static const PetscReal dd[3][2] = {
580*97a1619fSSatish Balay       {-1, 1},
581*97a1619fSSatish Balay       {-1, 1},
582*97a1619fSSatish Balay       {-1, 1}
583*97a1619fSSatish Balay     };
584*97a1619fSSatish Balay     static const PetscReal ww[3] = {1. / 6, 4. / 6, 1. / 6};
5859371c9d4SSatish Balay     *nq                          = 3;
5869371c9d4SSatish Balay     refweight                    = ww;
5879371c9d4SSatish Balay     refinterp                    = ii;
5889371c9d4SSatish Balay     refderiv                     = dd;
589c4762a1bSJed Brown   } break;
59098921bdaSJacob Faibussowitsch   default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unknown quadrature %d", (int)rd->quadrature);
591c4762a1bSJed Brown   }
592c4762a1bSJed Brown 
593c4762a1bSJed Brown   for (q = 0; q < *nq; q++) {
594c4762a1bSJed Brown     weight[q] = refweight[q] * hx;
595c4762a1bSJed Brown     for (j = 0; j < 2; j++) {
596c4762a1bSJed Brown       interp[q][j] = refinterp[q][j];
597c4762a1bSJed Brown       deriv[q][j]  = refderiv[q][j] / hx;
598c4762a1bSJed Brown     }
599c4762a1bSJed Brown   }
600c4762a1bSJed Brown   PetscFunctionReturn(0);
601c4762a1bSJed Brown }
602c4762a1bSJed Brown 
603c4762a1bSJed Brown /*
604c4762a1bSJed Brown  Finite element version
605c4762a1bSJed Brown */
6069371c9d4SSatish Balay static PetscErrorCode RDIFunction_FE(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, void *ctx) {
607c4762a1bSJed Brown   RD            rd = (RD)ctx;
608c4762a1bSJed Brown   RDNode       *x, *x0, *xdot, *f;
609c4762a1bSJed Brown   Vec           X0loc, Xloc, Xloc_t, Floc;
610c4762a1bSJed Brown   PetscReal     hx, Theta, dt, weight[5], interp[5][2], deriv[5][2];
611c4762a1bSJed Brown   DMDALocalInfo info;
612c4762a1bSJed Brown   PetscInt      i, j, q, nq;
613c4762a1bSJed Brown 
614c4762a1bSJed Brown   PetscFunctionBeginUser;
6159566063dSJacob Faibussowitsch   PetscCall(RDGetLocalArrays(rd, ts, X, Xdot, &Theta, &dt, &X0loc, &x0, &Xloc, &x, &Xloc_t, &xdot));
616c4762a1bSJed Brown 
6179566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(rd->da, &Floc));
6189566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(Floc));
6199566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(rd->da, Floc, &f));
6209566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(rd->da, &info));
621c4762a1bSJed Brown 
622c4762a1bSJed Brown   /* Set up shape functions and quadrature for elements (assumes a uniform grid) */
623c4762a1bSJed Brown   hx = rd->L / (info.mx - 1);
6249566063dSJacob Faibussowitsch   PetscCall(RDGetQuadrature(rd, hx, &nq, weight, interp, deriv));
625c4762a1bSJed Brown 
626c4762a1bSJed Brown   for (i = info.xs; i < PetscMin(info.xs + info.xm, info.mx - 1); i++) {
627c4762a1bSJed Brown     for (q = 0; q < nq; q++) {
628c4762a1bSJed Brown       PetscReal   rho = rd->rho;
629c4762a1bSJed Brown       PetscScalar Em_t, rad, D_R, D0_R;
630c4762a1bSJed Brown       RDNode      n, n0, nx, n0x, nt, ntx;
631c4762a1bSJed Brown       RDEvaluate(interp, deriv, q, x, i, &n, &nx);
632c4762a1bSJed Brown       RDEvaluate(interp, deriv, q, x0, i, &n0, &n0x);
633c4762a1bSJed Brown       RDEvaluate(interp, deriv, q, xdot, i, &nt, &ntx);
634c4762a1bSJed Brown 
635c4762a1bSJed Brown       rad = (1. - Theta) * RDRadiation(rd, &n0, 0) + Theta * RDRadiation(rd, &n, 0);
636c4762a1bSJed Brown       if (rd->endpoint) {
637c4762a1bSJed Brown         PetscScalar Em0, Em1;
638c4762a1bSJed Brown         RDMaterialEnergy(rd, &n0, &Em0, NULL);
639c4762a1bSJed Brown         RDMaterialEnergy(rd, &n, &Em1, NULL);
640c4762a1bSJed Brown         Em_t = (Em1 - Em0) / dt;
641c4762a1bSJed Brown       } else {
642c4762a1bSJed Brown         RDNode dEm;
643c4762a1bSJed Brown         RDMaterialEnergy(rd, &n, NULL, &dEm);
644c4762a1bSJed Brown         Em_t = dEm.E * nt.E + dEm.T * nt.T;
645c4762a1bSJed Brown       }
646c4762a1bSJed Brown       RDDiffusionCoefficient(rd, PETSC_TRUE, &n0, &n0x, &D0_R, 0, 0);
647c4762a1bSJed Brown       RDDiffusionCoefficient(rd, PETSC_TRUE, &n, &nx, &D_R, 0, 0);
648c4762a1bSJed Brown       for (j = 0; j < 2; j++) {
6499371c9d4SSatish 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));
650c4762a1bSJed Brown         f[i + j].T += interp[q][j] * weight[q] * (rho * Em_t + rad);
651c4762a1bSJed Brown       }
652c4762a1bSJed Brown     }
653c4762a1bSJed Brown   }
654c4762a1bSJed Brown   if (info.xs == 0) {
655c4762a1bSJed Brown     switch (rd->leftbc) {
656c4762a1bSJed Brown     case BC_ROBIN: {
657c4762a1bSJed Brown       PetscScalar D_R, D_R_bc;
658c4762a1bSJed Brown       PetscReal   ratio, bcTheta = rd->bcmidpoint ? Theta : 1.;
659c4762a1bSJed Brown       RDNode      n, nx;
660c4762a1bSJed Brown 
661c4762a1bSJed Brown       n.E  = (1 - bcTheta) * x0[0].E + bcTheta * x[0].E;
662c4762a1bSJed Brown       n.T  = (1 - bcTheta) * x0[0].T + bcTheta * x[0].T;
663c4762a1bSJed Brown       nx.E = (x[1].E - x[0].E) / hx;
664c4762a1bSJed Brown       nx.T = (x[1].T - x[0].T) / hx;
665c4762a1bSJed Brown       RDDiffusionCoefficient(rd, PETSC_TRUE, &n, &nx, &D_R, 0, 0);
666c4762a1bSJed Brown       RDDiffusionCoefficient(rd, rd->bclimit, &n, &nx, &D_R_bc, 0, 0);
667c4762a1bSJed Brown       ratio = PetscRealPart(D_R / D_R_bc);
6683c633725SBarry Smith       PetscCheck(ratio <= 1., PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Limited diffusivity is greater than unlimited");
6693c633725SBarry Smith       PetscCheck(ratio >= 1e-3, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Heavily limited diffusivity");
670c4762a1bSJed Brown       f[0].E += -ratio * 0.5 * (rd->Eapplied - n.E);
671c4762a1bSJed Brown     } break;
672c4762a1bSJed Brown     case BC_NEUMANN:
673c4762a1bSJed Brown       /* homogeneous Neumann is the natural condition */
674c4762a1bSJed Brown       break;
67563a3b9bcSJacob Faibussowitsch     default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Case %" PetscInt_FMT, rd->initial);
676c4762a1bSJed Brown     }
677c4762a1bSJed Brown   }
678c4762a1bSJed Brown 
6799566063dSJacob Faibussowitsch   PetscCall(RDRestoreLocalArrays(rd, &X0loc, &x0, &Xloc, &x, &Xloc_t, &xdot));
6809566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(rd->da, Floc, &f));
6819566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(F));
6829566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalBegin(rd->da, Floc, ADD_VALUES, F));
6839566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalEnd(rd->da, Floc, ADD_VALUES, F));
6849566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(rd->da, &Floc));
685c4762a1bSJed Brown 
6869566063dSJacob Faibussowitsch   if (rd->monitor_residual) PetscCall(RDStateView(rd, X, Xdot, F));
687c4762a1bSJed Brown   PetscFunctionReturn(0);
688c4762a1bSJed Brown }
689c4762a1bSJed Brown 
6909371c9d4SSatish Balay static PetscErrorCode RDIJacobian_FE(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal a, Mat A, Mat B, void *ctx) {
691c4762a1bSJed Brown   RD            rd = (RD)ctx;
692c4762a1bSJed Brown   RDNode       *x, *x0, *xdot;
693c4762a1bSJed Brown   Vec           X0loc, Xloc, Xloc_t;
694c4762a1bSJed Brown   PetscReal     hx, Theta, dt, weight[5], interp[5][2], deriv[5][2];
695c4762a1bSJed Brown   DMDALocalInfo info;
696c4762a1bSJed Brown   PetscInt      i, j, k, q, nq;
697c4762a1bSJed Brown   PetscScalar   K[4][4];
698c4762a1bSJed Brown 
699c4762a1bSJed Brown   PetscFunctionBeginUser;
7009566063dSJacob Faibussowitsch   PetscCall(RDGetLocalArrays(rd, ts, X, Xdot, &Theta, &dt, &X0loc, &x0, &Xloc, &x, &Xloc_t, &xdot));
7019566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(rd->da, &info));
702c4762a1bSJed Brown   hx = rd->L / (info.mx - 1);
7039566063dSJacob Faibussowitsch   PetscCall(RDGetQuadrature(rd, hx, &nq, weight, interp, deriv));
7049566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(B));
705c4762a1bSJed Brown   for (i = info.xs; i < PetscMin(info.xs + info.xm, info.mx - 1); i++) {
706c4762a1bSJed Brown     PetscInt rc[2];
707c4762a1bSJed Brown 
7089371c9d4SSatish Balay     rc[0] = i;
7099371c9d4SSatish Balay     rc[1] = i + 1;
7109566063dSJacob Faibussowitsch     PetscCall(PetscMemzero(K, sizeof(K)));
711c4762a1bSJed Brown     for (q = 0; q < nq; q++) {
712c4762a1bSJed Brown       PetscScalar              D_R;
713c4762a1bSJed Brown       PETSC_UNUSED PetscScalar rad;
714c4762a1bSJed Brown       RDNode                   n, nx, nt, ntx, drad, dD_R, dxD_R, dEm;
715c4762a1bSJed Brown       RDEvaluate(interp, deriv, q, x, i, &n, &nx);
716c4762a1bSJed Brown       RDEvaluate(interp, deriv, q, xdot, i, &nt, &ntx);
717c4762a1bSJed Brown       rad = RDRadiation(rd, &n, &drad);
718c4762a1bSJed Brown       RDDiffusionCoefficient(rd, PETSC_TRUE, &n, &nx, &D_R, &dD_R, &dxD_R);
719c4762a1bSJed Brown       RDMaterialEnergy(rd, &n, NULL, &dEm);
720c4762a1bSJed Brown       for (j = 0; j < 2; j++) {
721c4762a1bSJed Brown         for (k = 0; k < 2; k++) {
7229371c9d4SSatish 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]));
7239371c9d4SSatish 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);
724c4762a1bSJed Brown           K[j * 2 + 1][k * 2 + 0] += interp[q][j] * weight[q] * drad.E * interp[q][k];
725c4762a1bSJed Brown           K[j * 2 + 1][k * 2 + 1] += interp[q][j] * weight[q] * (a * rd->rho * dEm.T + drad.T) * interp[q][k];
726c4762a1bSJed Brown         }
727c4762a1bSJed Brown       }
728c4762a1bSJed Brown     }
7299566063dSJacob Faibussowitsch     PetscCall(MatSetValuesBlocked(B, 2, rc, 2, rc, &K[0][0], ADD_VALUES));
730c4762a1bSJed Brown   }
731c4762a1bSJed Brown   if (info.xs == 0) {
732c4762a1bSJed Brown     switch (rd->leftbc) {
733c4762a1bSJed Brown     case BC_ROBIN: {
734c4762a1bSJed Brown       PetscScalar D_R, D_R_bc;
735c4762a1bSJed Brown       PetscReal   ratio;
736c4762a1bSJed Brown       RDNode      n, nx;
737c4762a1bSJed Brown 
738c4762a1bSJed Brown       n.E  = (1 - Theta) * x0[0].E + Theta * x[0].E;
739c4762a1bSJed Brown       n.T  = (1 - Theta) * x0[0].T + Theta * x[0].T;
740c4762a1bSJed Brown       nx.E = (x[1].E - x[0].E) / hx;
741c4762a1bSJed Brown       nx.T = (x[1].T - x[0].T) / hx;
742c4762a1bSJed Brown       RDDiffusionCoefficient(rd, PETSC_TRUE, &n, &nx, &D_R, 0, 0);
743c4762a1bSJed Brown       RDDiffusionCoefficient(rd, rd->bclimit, &n, &nx, &D_R_bc, 0, 0);
744c4762a1bSJed Brown       ratio = PetscRealPart(D_R / D_R_bc);
7459566063dSJacob Faibussowitsch       PetscCall(MatSetValue(B, 0, 0, ratio * 0.5, ADD_VALUES));
746c4762a1bSJed Brown     } break;
747c4762a1bSJed Brown     case BC_NEUMANN:
748c4762a1bSJed Brown       /* homogeneous Neumann is the natural condition */
749c4762a1bSJed Brown       break;
75063a3b9bcSJacob Faibussowitsch     default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Case %" PetscInt_FMT, rd->initial);
751c4762a1bSJed Brown     }
752c4762a1bSJed Brown   }
753c4762a1bSJed Brown 
7549566063dSJacob Faibussowitsch   PetscCall(RDRestoreLocalArrays(rd, &X0loc, &x0, &Xloc, &x, &Xloc_t, &xdot));
7559566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
7569566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
757c4762a1bSJed Brown   if (A != B) {
7589566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
7599566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
760c4762a1bSJed Brown   }
761c4762a1bSJed Brown   PetscFunctionReturn(0);
762c4762a1bSJed Brown }
763c4762a1bSJed Brown 
764c4762a1bSJed Brown /* Temperature that is in equilibrium with the radiation density */
7659371c9d4SSatish Balay static PetscScalar RDRadiationTemperature(RD rd, PetscScalar E) {
7669371c9d4SSatish Balay   return PetscPowScalar(E * rd->c / (4. * rd->sigma_b), 0.25);
7679371c9d4SSatish Balay }
768c4762a1bSJed Brown 
7699371c9d4SSatish Balay static PetscErrorCode RDInitialState(RD rd, Vec X) {
770c4762a1bSJed Brown   DMDALocalInfo info;
771c4762a1bSJed Brown   PetscInt      i;
772c4762a1bSJed Brown   RDNode       *x;
773c4762a1bSJed Brown 
774c4762a1bSJed Brown   PetscFunctionBeginUser;
7759566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(rd->da, &info));
7769566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(rd->da, X, &x));
777c4762a1bSJed Brown   for (i = info.xs; i < info.xs + info.xm; i++) {
778c4762a1bSJed Brown     PetscReal coord = i * rd->L / (info.mx - 1);
779c4762a1bSJed Brown     switch (rd->initial) {
780c4762a1bSJed Brown     case 1:
781c4762a1bSJed Brown       x[i].E = 0.001;
782c4762a1bSJed Brown       x[i].T = RDRadiationTemperature(rd, x[i].E);
783c4762a1bSJed Brown       break;
784c4762a1bSJed Brown     case 2:
785c4762a1bSJed Brown       x[i].E = 0.001 + 100. * PetscExpReal(-PetscSqr(coord / 0.1));
786c4762a1bSJed Brown       x[i].T = RDRadiationTemperature(rd, x[i].E);
787c4762a1bSJed Brown       break;
788c4762a1bSJed Brown     case 3:
789c4762a1bSJed Brown       x[i].E = 7.56e-2 * rd->unit.Joule / PetscPowScalarInt(rd->unit.meter, 3);
790c4762a1bSJed Brown       x[i].T = RDRadiationTemperature(rd, x[i].E);
791c4762a1bSJed Brown       break;
79263a3b9bcSJacob Faibussowitsch     default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "No initial state %" PetscInt_FMT, rd->initial);
793c4762a1bSJed Brown     }
794c4762a1bSJed Brown   }
7959566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(rd->da, X, &x));
796c4762a1bSJed Brown   PetscFunctionReturn(0);
797c4762a1bSJed Brown }
798c4762a1bSJed Brown 
7999371c9d4SSatish Balay static PetscErrorCode RDView(RD rd, Vec X, PetscViewer viewer) {
800c4762a1bSJed Brown   Vec             Y;
801c4762a1bSJed Brown   const RDNode   *x;
802c4762a1bSJed Brown   PetscScalar    *y;
803c4762a1bSJed Brown   PetscInt        i, m, M;
804c4762a1bSJed Brown   const PetscInt *lx;
805c4762a1bSJed Brown   DM              da;
806c4762a1bSJed Brown   MPI_Comm        comm;
807c4762a1bSJed Brown 
808c4762a1bSJed Brown   PetscFunctionBeginUser;
809c4762a1bSJed Brown   /*
810c4762a1bSJed Brown     Create a DMDA (one dof per node, zero stencil width, same layout) to hold Trad
811c4762a1bSJed Brown     (radiation temperature).  It is not necessary to create a DMDA for this, but this way
812c4762a1bSJed Brown     output and visualization will have meaningful variable names and correct scales.
813c4762a1bSJed Brown   */
8149566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(rd->da, 0, &M, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0));
8159566063dSJacob Faibussowitsch   PetscCall(DMDAGetOwnershipRanges(rd->da, &lx, 0, 0));
8169566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)rd->da, &comm));
8179566063dSJacob Faibussowitsch   PetscCall(DMDACreate1d(comm, DM_BOUNDARY_NONE, M, 1, 0, lx, &da));
8189566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
8199566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
8209566063dSJacob Faibussowitsch   PetscCall(DMDASetUniformCoordinates(da, 0., rd->L, 0., 0., 0., 0.));
8219566063dSJacob Faibussowitsch   PetscCall(DMDASetFieldName(da, 0, "T_rad"));
8229566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(da, &Y));
823c4762a1bSJed Brown 
824c4762a1bSJed Brown   /* Compute the radiation temperature from the solution at each node */
8259566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(Y, &m));
8269566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, (const PetscScalar **)&x));
8279566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Y, &y));
828c4762a1bSJed Brown   for (i = 0; i < m; i++) y[i] = RDRadiationTemperature(rd, x[i].E);
8299566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, (const PetscScalar **)&x));
8309566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Y, &y));
831c4762a1bSJed Brown 
8329566063dSJacob Faibussowitsch   PetscCall(VecView(Y, viewer));
8339566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Y));
8349566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
835c4762a1bSJed Brown   PetscFunctionReturn(0);
836c4762a1bSJed Brown }
837c4762a1bSJed Brown 
8389371c9d4SSatish Balay static PetscErrorCode RDTestDifferentiation(RD rd) {
839c4762a1bSJed Brown   MPI_Comm    comm;
840c4762a1bSJed Brown   RDNode      n, nx;
841c4762a1bSJed Brown   PetscScalar epsilon;
842c4762a1bSJed Brown 
843c4762a1bSJed Brown   PetscFunctionBeginUser;
8449566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)rd->da, &comm));
845c4762a1bSJed Brown   epsilon = 1e-8;
846c4762a1bSJed Brown   {
847c4762a1bSJed Brown     RDNode      dEm, fdEm;
848c4762a1bSJed Brown     PetscScalar T0 = 1000., T1 = T0 * (1. + epsilon), Em0, Em1;
849c4762a1bSJed Brown     n.E = 1.;
850c4762a1bSJed Brown     n.T = T0;
851c4762a1bSJed Brown     rd->MaterialEnergy(rd, &n, &Em0, &dEm);
852c4762a1bSJed Brown     n.E = 1. + epsilon;
853c4762a1bSJed Brown     n.T = T0;
854c4762a1bSJed Brown     rd->MaterialEnergy(rd, &n, &Em1, 0);
855c4762a1bSJed Brown     fdEm.E = (Em1 - Em0) / epsilon;
856c4762a1bSJed Brown     n.E    = 1.;
857c4762a1bSJed Brown     n.T    = T1;
858c4762a1bSJed Brown     rd->MaterialEnergy(rd, &n, &Em1, 0);
859c4762a1bSJed Brown     fdEm.T = (Em1 - Em0) / (T0 * epsilon);
8609371c9d4SSatish 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),
8619371c9d4SSatish Balay                           (double)PetscRealPart(dEm.T - fdEm.T)));
862c4762a1bSJed Brown   }
863c4762a1bSJed Brown   {
864c4762a1bSJed Brown     PetscScalar D0, D;
865c4762a1bSJed Brown     RDNode      dD, dxD, fdD, fdxD;
8669371c9d4SSatish Balay     n.E  = 1.;
8679371c9d4SSatish Balay     n.T  = 1.;
8689371c9d4SSatish Balay     nx.E = 1.;
8699371c9d4SSatish Balay     n.T  = 1.;
870c4762a1bSJed Brown     RDDiffusionCoefficient(rd, rd->bclimit, &n, &nx, &D0, &dD, &dxD);
8719371c9d4SSatish Balay     n.E  = 1. + epsilon;
8729371c9d4SSatish Balay     n.T  = 1.;
8739371c9d4SSatish Balay     nx.E = 1.;
8749371c9d4SSatish Balay     n.T  = 1.;
8759371c9d4SSatish Balay     RDDiffusionCoefficient(rd, rd->bclimit, &n, &nx, &D, 0, 0);
8769371c9d4SSatish Balay     fdD.E = (D - D0) / epsilon;
8779371c9d4SSatish Balay     n.E   = 1;
8789371c9d4SSatish Balay     n.T   = 1. + epsilon;
8799371c9d4SSatish Balay     nx.E  = 1.;
8809371c9d4SSatish Balay     n.T   = 1.;
8819371c9d4SSatish Balay     RDDiffusionCoefficient(rd, rd->bclimit, &n, &nx, &D, 0, 0);
8829371c9d4SSatish Balay     fdD.T = (D - D0) / epsilon;
8839371c9d4SSatish Balay     n.E   = 1;
8849371c9d4SSatish Balay     n.T   = 1.;
8859371c9d4SSatish Balay     nx.E  = 1. + epsilon;
8869371c9d4SSatish Balay     n.T   = 1.;
8879371c9d4SSatish Balay     RDDiffusionCoefficient(rd, rd->bclimit, &n, &nx, &D, 0, 0);
8889371c9d4SSatish Balay     fdxD.E = (D - D0) / epsilon;
8899371c9d4SSatish Balay     n.E    = 1;
8909371c9d4SSatish Balay     n.T    = 1.;
8919371c9d4SSatish Balay     nx.E   = 1.;
8929371c9d4SSatish Balay     n.T    = 1. + epsilon;
8939371c9d4SSatish Balay     RDDiffusionCoefficient(rd, rd->bclimit, &n, &nx, &D, 0, 0);
8949371c9d4SSatish Balay     fdxD.T = (D - D0) / epsilon;
8959371c9d4SSatish 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),
8969371c9d4SSatish Balay                           (double)PetscRealPart(dD.T - fdD.T)));
8979371c9d4SSatish 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),
8989371c9d4SSatish Balay                           (double)PetscRealPart(dxD.T - fdxD.T)));
899c4762a1bSJed Brown   }
900c4762a1bSJed Brown   {
901c4762a1bSJed Brown     PetscInt    i;
902c4762a1bSJed Brown     PetscReal   hx = 1.;
903c4762a1bSJed Brown     PetscScalar a0;
904c4762a1bSJed Brown     RDNode      n0[3], n1[3], d[3], fd[3];
905c4762a1bSJed Brown 
906c4762a1bSJed Brown     n0[0].E = 1.;
907c4762a1bSJed Brown     n0[0].T = 1.;
908c4762a1bSJed Brown     n0[1].E = 5.;
909c4762a1bSJed Brown     n0[1].T = 3.;
910c4762a1bSJed Brown     n0[2].E = 4.;
911c4762a1bSJed Brown     n0[2].T = 2.;
912c4762a1bSJed Brown     a0      = RDDiffusion(rd, hx, n0, 1, d);
913c4762a1bSJed Brown     for (i = 0; i < 3; i++) {
9149371c9d4SSatish Balay       PetscCall(PetscMemcpy(n1, n0, sizeof(n0)));
9159371c9d4SSatish Balay       n1[i].E += epsilon;
916c4762a1bSJed Brown       fd[i].E = (RDDiffusion(rd, hx, n1, 1, 0) - a0) / epsilon;
9179371c9d4SSatish Balay       PetscCall(PetscMemcpy(n1, n0, sizeof(n0)));
9189371c9d4SSatish Balay       n1[i].T += epsilon;
919c4762a1bSJed Brown       fd[i].T = (RDDiffusion(rd, hx, n1, 1, 0) - a0) / epsilon;
9209371c9d4SSatish 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),
9219371c9d4SSatish Balay                             (double)PetscRealPart(d[i].E - fd[i].E), (double)PetscRealPart(d[i].T - fd[i].T)));
922c4762a1bSJed Brown     }
923c4762a1bSJed Brown   }
924c4762a1bSJed Brown   {
925c4762a1bSJed Brown     PetscScalar rad0, rad;
926c4762a1bSJed Brown     RDNode      drad, fdrad;
9279371c9d4SSatish Balay     n.E     = 1.;
9289371c9d4SSatish Balay     n.T     = 1.;
929c4762a1bSJed Brown     rad0    = RDRadiation(rd, &n, &drad);
9309371c9d4SSatish Balay     n.E     = 1. + epsilon;
9319371c9d4SSatish Balay     n.T     = 1.;
9329371c9d4SSatish Balay     rad     = RDRadiation(rd, &n, 0);
9339371c9d4SSatish Balay     fdrad.E = (rad - rad0) / epsilon;
9349371c9d4SSatish Balay     n.E     = 1.;
9359371c9d4SSatish Balay     n.T     = 1. + epsilon;
9369371c9d4SSatish Balay     rad     = RDRadiation(rd, &n, 0);
9379371c9d4SSatish Balay     fdrad.T = (rad - rad0) / epsilon;
9389371c9d4SSatish 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),
9399371c9d4SSatish Balay                           (double)PetscRealPart(drad.T - fdrad.T)));
940c4762a1bSJed Brown   }
941c4762a1bSJed Brown   PetscFunctionReturn(0);
942c4762a1bSJed Brown }
943c4762a1bSJed Brown 
9449371c9d4SSatish Balay static PetscErrorCode RDCreate(MPI_Comm comm, RD *inrd) {
945c4762a1bSJed Brown   RD        rd;
946c4762a1bSJed Brown   PetscReal meter = 0, kilogram = 0, second = 0, Kelvin = 0, Joule = 0, Watt = 0;
947c4762a1bSJed Brown 
948c4762a1bSJed Brown   PetscFunctionBeginUser;
949c4762a1bSJed Brown   *inrd = 0;
9509566063dSJacob Faibussowitsch   PetscCall(PetscNew(&rd));
951c4762a1bSJed Brown 
952d0609cedSBarry Smith   PetscOptionsBegin(comm, NULL, "Options for nonequilibrium radiation-diffusion with RD ionization", NULL);
953c4762a1bSJed Brown   {
954c4762a1bSJed Brown     rd->initial = 1;
9559566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-rd_initial", "Initial condition (1=Marshak, 2=Blast, 3=Marshak+)", "", rd->initial, &rd->initial, 0));
956c4762a1bSJed Brown     switch (rd->initial) {
957c4762a1bSJed Brown     case 1:
958c4762a1bSJed Brown     case 2:
959c4762a1bSJed Brown       rd->unit.kilogram = 1.;
960c4762a1bSJed Brown       rd->unit.meter    = 1.;
961c4762a1bSJed Brown       rd->unit.second   = 1.;
962c4762a1bSJed Brown       rd->unit.Kelvin   = 1.;
963c4762a1bSJed Brown       break;
964c4762a1bSJed Brown     case 3:
965c4762a1bSJed Brown       rd->unit.kilogram = 1.e12;
966c4762a1bSJed Brown       rd->unit.meter    = 1.;
967c4762a1bSJed Brown       rd->unit.second   = 1.e9;
968c4762a1bSJed Brown       rd->unit.Kelvin   = 1.;
969c4762a1bSJed Brown       break;
97063a3b9bcSJacob Faibussowitsch     default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unknown initial condition %" PetscInt_FMT, rd->initial);
971c4762a1bSJed Brown     }
972c4762a1bSJed Brown     /* Fundamental units */
9739566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-rd_unit_meter", "Length of 1 meter in nondimensional units", "", rd->unit.meter, &rd->unit.meter, 0));
9749566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-rd_unit_kilogram", "Mass of 1 kilogram in nondimensional units", "", rd->unit.kilogram, &rd->unit.kilogram, 0));
9759566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-rd_unit_second", "Time of a second in nondimensional units", "", rd->unit.second, &rd->unit.second, 0));
9769566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-rd_unit_Kelvin", "Temperature of a Kelvin in nondimensional units", "", rd->unit.Kelvin, &rd->unit.Kelvin, 0));
977c4762a1bSJed Brown     /* Derived units */
978c4762a1bSJed Brown     rd->unit.Joule = rd->unit.kilogram * PetscSqr(rd->unit.meter / rd->unit.second);
979c4762a1bSJed Brown     rd->unit.Watt  = rd->unit.Joule / rd->unit.second;
980c4762a1bSJed Brown     /* Local aliases */
981c4762a1bSJed Brown     meter          = rd->unit.meter;
982c4762a1bSJed Brown     kilogram       = rd->unit.kilogram;
983c4762a1bSJed Brown     second         = rd->unit.second;
984c4762a1bSJed Brown     Kelvin         = rd->unit.Kelvin;
985c4762a1bSJed Brown     Joule          = rd->unit.Joule;
986c4762a1bSJed Brown     Watt           = rd->unit.Watt;
987c4762a1bSJed Brown 
9889566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-rd_monitor_residual", "Display residuals every time they are evaluated", "", rd->monitor_residual, &rd->monitor_residual, NULL));
9899566063dSJacob Faibussowitsch     PetscCall(PetscOptionsEnum("-rd_discretization", "Discretization type", "", DiscretizationTypes, (PetscEnum)rd->discretization, (PetscEnum *)&rd->discretization, NULL));
990c4762a1bSJed Brown     if (rd->discretization == DISCRETIZATION_FE) {
991c4762a1bSJed Brown       rd->quadrature = QUADRATURE_GAUSS2;
9929566063dSJacob Faibussowitsch       PetscCall(PetscOptionsEnum("-rd_quadrature", "Finite element quadrature", "", QuadratureTypes, (PetscEnum)rd->quadrature, (PetscEnum *)&rd->quadrature, NULL));
993c4762a1bSJed Brown     }
9949566063dSJacob Faibussowitsch     PetscCall(PetscOptionsEnum("-rd_jacobian", "Type of finite difference Jacobian", "", JacobianTypes, (PetscEnum)rd->jacobian, (PetscEnum *)&rd->jacobian, NULL));
995c4762a1bSJed Brown     switch (rd->initial) {
996c4762a1bSJed Brown     case 1:
997c4762a1bSJed Brown       rd->leftbc     = BC_ROBIN;
998c4762a1bSJed Brown       rd->Eapplied   = 4 * rd->unit.Joule / PetscPowRealInt(rd->unit.meter, 3);
999c4762a1bSJed Brown       rd->L          = 1. * rd->unit.meter;
1000c4762a1bSJed Brown       rd->beta       = 3.0;
1001c4762a1bSJed Brown       rd->gamma      = 3.0;
1002c4762a1bSJed Brown       rd->final_time = 3 * second;
1003c4762a1bSJed Brown       break;
1004c4762a1bSJed Brown     case 2:
1005c4762a1bSJed Brown       rd->leftbc     = BC_NEUMANN;
1006c4762a1bSJed Brown       rd->Eapplied   = 0.;
1007c4762a1bSJed Brown       rd->L          = 1. * rd->unit.meter;
1008c4762a1bSJed Brown       rd->beta       = 3.0;
1009c4762a1bSJed Brown       rd->gamma      = 3.0;
1010c4762a1bSJed Brown       rd->final_time = 1 * second;
1011c4762a1bSJed Brown       break;
1012c4762a1bSJed Brown     case 3:
1013c4762a1bSJed Brown       rd->leftbc     = BC_ROBIN;
1014c4762a1bSJed Brown       rd->Eapplied   = 7.503e6 * rd->unit.Joule / PetscPowRealInt(rd->unit.meter, 3);
1015c4762a1bSJed Brown       rd->L          = 5. * rd->unit.meter;
1016c4762a1bSJed Brown       rd->beta       = 3.5;
1017c4762a1bSJed Brown       rd->gamma      = 3.5;
1018c4762a1bSJed Brown       rd->final_time = 20e-9 * second;
1019c4762a1bSJed Brown       break;
102063a3b9bcSJacob Faibussowitsch     default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Initial %" PetscInt_FMT, rd->initial);
1021c4762a1bSJed Brown     }
10229566063dSJacob Faibussowitsch     PetscCall(PetscOptionsEnum("-rd_leftbc", "Left boundary condition", "", BCTypes, (PetscEnum)rd->leftbc, (PetscEnum *)&rd->leftbc, NULL));
10239566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-rd_E_applied", "Radiation flux at left end of domain", "", rd->Eapplied, &rd->Eapplied, NULL));
10249566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-rd_beta", "Thermal exponent for photon absorption", "", rd->beta, &rd->beta, NULL));
10259566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-rd_gamma", "Thermal exponent for diffusion coefficient", "", rd->gamma, &rd->gamma, NULL));
10269566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-rd_view_draw", "Draw final solution", "", rd->view_draw, &rd->view_draw, NULL));
10279566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-rd_endpoint", "Discretize using endpoints (like trapezoid rule) instead of midpoint", "", rd->endpoint, &rd->endpoint, NULL));
10289566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-rd_bcmidpoint", "Impose the boundary condition at the midpoint (Theta) of the interval", "", rd->bcmidpoint, &rd->bcmidpoint, NULL));
10299566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-rd_bclimit", "Limit diffusion coefficient in definition of Robin boundary condition", "", rd->bclimit, &rd->bclimit, NULL));
10309566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-rd_test_diff", "Test differentiation in constitutive relations", "", rd->test_diff, &rd->test_diff, NULL));
10319566063dSJacob Faibussowitsch     PetscCall(PetscOptionsString("-rd_view_binary", "File name to hold final solution", "", rd->view_binary, rd->view_binary, sizeof(rd->view_binary), NULL));
1032c4762a1bSJed Brown   }
1033d0609cedSBarry Smith   PetscOptionsEnd();
1034c4762a1bSJed Brown 
1035c4762a1bSJed Brown   switch (rd->initial) {
1036c4762a1bSJed Brown   case 1:
1037c4762a1bSJed Brown   case 2:
1038c4762a1bSJed Brown     rd->rho            = 1.;
1039c4762a1bSJed Brown     rd->c              = 1.;
1040c4762a1bSJed Brown     rd->K_R            = 1.;
1041c4762a1bSJed Brown     rd->K_p            = 1.;
1042c4762a1bSJed Brown     rd->sigma_b        = 0.25;
1043c4762a1bSJed Brown     rd->MaterialEnergy = RDMaterialEnergy_Reduced;
1044c4762a1bSJed Brown     break;
1045c4762a1bSJed Brown   case 3:
1046c4762a1bSJed Brown     /* Table 2 */
1047c4762a1bSJed Brown     rd->rho            = 1.17e-3 * kilogram / (meter * meter * meter);                                                    /* density */
1048c4762a1bSJed Brown     rd->K_R            = 7.44e18 * PetscPowRealInt(meter, 5) * PetscPowReal(Kelvin, 3.5) * PetscPowRealInt(kilogram, -2); /*  */
1049c4762a1bSJed Brown     rd->K_p            = 2.33e20 * PetscPowRealInt(meter, 5) * PetscPowReal(Kelvin, 3.5) * PetscPowRealInt(kilogram, -2); /*  */
1050c4762a1bSJed Brown     rd->I_H            = 2.179e-18 * Joule;                                                                               /* Hydrogen ionization potential */
1051c4762a1bSJed Brown     rd->m_p            = 1.673e-27 * kilogram;                                                                            /* proton mass */
1052c4762a1bSJed Brown     rd->m_e            = 9.109e-31 * kilogram;                                                                            /* electron mass */
1053c4762a1bSJed Brown     rd->h              = 6.626e-34 * Joule * second;                                                                      /* Planck's constant */
1054c4762a1bSJed Brown     rd->k              = 1.381e-23 * Joule / Kelvin;                                                                      /* Boltzman constant */
1055c4762a1bSJed Brown     rd->c              = 3.00e8 * meter / second;                                                                         /* speed of light */
1056c4762a1bSJed Brown     rd->sigma_b        = 5.67e-8 * Watt * PetscPowRealInt(meter, -2) * PetscPowRealInt(Kelvin, -4);                       /* Stefan-Boltzman constant */
1057c4762a1bSJed Brown     rd->MaterialEnergy = RDMaterialEnergy_Saha;
1058c4762a1bSJed Brown     break;
1059c4762a1bSJed Brown   }
1060c4762a1bSJed Brown 
10619566063dSJacob Faibussowitsch   PetscCall(DMDACreate1d(comm, DM_BOUNDARY_NONE, 20, sizeof(RDNode) / sizeof(PetscScalar), 1, NULL, &rd->da));
10629566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(rd->da));
10639566063dSJacob Faibussowitsch   PetscCall(DMSetUp(rd->da));
10649566063dSJacob Faibussowitsch   PetscCall(DMDASetFieldName(rd->da, 0, "E"));
10659566063dSJacob Faibussowitsch   PetscCall(DMDASetFieldName(rd->da, 1, "T"));
10669566063dSJacob Faibussowitsch   PetscCall(DMDASetUniformCoordinates(rd->da, 0., 1., 0., 0., 0., 0.));
1067c4762a1bSJed Brown 
1068c4762a1bSJed Brown   *inrd = rd;
1069c4762a1bSJed Brown   PetscFunctionReturn(0);
1070c4762a1bSJed Brown }
1071c4762a1bSJed Brown 
10729371c9d4SSatish Balay int main(int argc, char *argv[]) {
1073c4762a1bSJed Brown   RD        rd;
1074c4762a1bSJed Brown   TS        ts;
1075c4762a1bSJed Brown   SNES      snes;
1076c4762a1bSJed Brown   Vec       X;
1077c4762a1bSJed Brown   Mat       A, B;
1078c4762a1bSJed Brown   PetscInt  steps;
1079c4762a1bSJed Brown   PetscReal ftime;
1080c4762a1bSJed Brown 
1081327415f7SBarry Smith   PetscFunctionBeginUser;
10829ded082cSBarry Smith   PetscCall(PetscInitialize(&argc, &argv, 0, NULL));
10839566063dSJacob Faibussowitsch   PetscCall(RDCreate(PETSC_COMM_WORLD, &rd));
10849566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(rd->da, &X));
10859566063dSJacob Faibussowitsch   PetscCall(DMSetMatType(rd->da, MATAIJ));
10869566063dSJacob Faibussowitsch   PetscCall(DMCreateMatrix(rd->da, &B));
10879566063dSJacob Faibussowitsch   PetscCall(RDInitialState(rd, X));
1088c4762a1bSJed Brown 
10899566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
10909566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
10919566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts, TSTHETA));
10929566063dSJacob Faibussowitsch   PetscCall(TSSetDM(ts, rd->da));
1093c4762a1bSJed Brown   switch (rd->discretization) {
1094c4762a1bSJed Brown   case DISCRETIZATION_FD:
10959566063dSJacob Faibussowitsch     PetscCall(TSSetIFunction(ts, NULL, RDIFunction_FD, rd));
10969566063dSJacob Faibussowitsch     if (rd->jacobian == JACOBIAN_ANALYTIC) PetscCall(TSSetIJacobian(ts, B, B, RDIJacobian_FD, rd));
1097c4762a1bSJed Brown     break;
1098c4762a1bSJed Brown   case DISCRETIZATION_FE:
10999566063dSJacob Faibussowitsch     PetscCall(TSSetIFunction(ts, NULL, RDIFunction_FE, rd));
11009566063dSJacob Faibussowitsch     if (rd->jacobian == JACOBIAN_ANALYTIC) PetscCall(TSSetIJacobian(ts, B, B, RDIJacobian_FE, rd));
1101c4762a1bSJed Brown     break;
1102c4762a1bSJed Brown   }
11039566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, rd->final_time));
11049566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, 1e-3));
11059566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
11069566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
1107c4762a1bSJed Brown 
1108c4762a1bSJed Brown   A = B;
11099566063dSJacob Faibussowitsch   PetscCall(TSGetSNES(ts, &snes));
1110c4762a1bSJed Brown   switch (rd->jacobian) {
11119371c9d4SSatish Balay   case JACOBIAN_ANALYTIC: break;
11129371c9d4SSatish Balay   case JACOBIAN_MATRIXFREE: break;
1113c4762a1bSJed Brown   case JACOBIAN_FD_COLORING: {
11149566063dSJacob Faibussowitsch     PetscCall(SNESSetJacobian(snes, A, B, SNESComputeJacobianDefaultColor, 0));
1115c4762a1bSJed Brown   } break;
11169371c9d4SSatish Balay   case JACOBIAN_FD_FULL: PetscCall(SNESSetJacobian(snes, A, B, SNESComputeJacobianDefault, ts)); break;
1117c4762a1bSJed Brown   }
1118c4762a1bSJed Brown 
11191baa6e33SBarry Smith   if (rd->test_diff) PetscCall(RDTestDifferentiation(rd));
11209566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, X));
11219566063dSJacob Faibussowitsch   PetscCall(TSGetSolveTime(ts, &ftime));
11229566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts, &steps));
112363a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Steps %" PetscInt_FMT "  final time %g\n", steps, (double)ftime));
11241baa6e33SBarry Smith   if (rd->view_draw) PetscCall(RDView(rd, X, PETSC_VIEWER_DRAW_WORLD));
1125c4762a1bSJed Brown   if (rd->view_binary[0]) {
1126c4762a1bSJed Brown     PetscViewer viewer;
11279566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, rd->view_binary, FILE_MODE_WRITE, &viewer));
11289566063dSJacob Faibussowitsch     PetscCall(RDView(rd, X, viewer));
11299566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewer));
1130c4762a1bSJed Brown   }
11319566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&X));
11329566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
11339566063dSJacob Faibussowitsch   PetscCall(RDDestroy(&rd));
11349566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
11359566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
1136b122ec5aSJacob Faibussowitsch   return 0;
1137c4762a1bSJed Brown }
1138c4762a1bSJed Brown /*TEST
1139c4762a1bSJed Brown 
1140c4762a1bSJed Brown     test:
1141c4762a1bSJed 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
1142c4762a1bSJed Brown       requires: !single
1143c4762a1bSJed Brown 
1144c4762a1bSJed Brown     test:
1145c4762a1bSJed Brown       suffix: 2
1146c4762a1bSJed 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
1147c4762a1bSJed Brown       requires: !single
1148c4762a1bSJed Brown 
1149c4762a1bSJed Brown     test:
1150c4762a1bSJed Brown       suffix: 3
1151c4762a1bSJed Brown       nsize: 2
1152c4762a1bSJed 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
1153c4762a1bSJed Brown       requires: !single
1154c4762a1bSJed Brown 
1155c4762a1bSJed Brown TEST*/
1156