xref: /petsc/src/ts/tutorials/ex10.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
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 
16*9371c9d4SSatish Balay typedef enum {
17*9371c9d4SSatish Balay   BC_DIRICHLET,
18*9371c9d4SSatish Balay   BC_NEUMANN,
19*9371c9d4SSatish Balay   BC_ROBIN
20*9371c9d4SSatish Balay } BCType;
21c4762a1bSJed Brown static const char *const BCTypes[] = {"DIRICHLET", "NEUMANN", "ROBIN", "BCType", "BC_", 0};
22*9371c9d4SSatish Balay typedef enum {
23*9371c9d4SSatish Balay   JACOBIAN_ANALYTIC,
24*9371c9d4SSatish Balay   JACOBIAN_MATRIXFREE,
25*9371c9d4SSatish Balay   JACOBIAN_FD_COLORING,
26*9371c9d4SSatish Balay   JACOBIAN_FD_FULL
27*9371c9d4SSatish Balay } JacobianType;
28c4762a1bSJed Brown static const char *const JacobianTypes[] = {"ANALYTIC", "MATRIXFREE", "FD_COLORING", "FD_FULL", "JacobianType", "FD_", 0};
29*9371c9d4SSatish Balay typedef enum {
30*9371c9d4SSatish Balay   DISCRETIZATION_FD,
31*9371c9d4SSatish Balay   DISCRETIZATION_FE
32*9371c9d4SSatish Balay } DiscretizationType;
33c4762a1bSJed Brown static const char *const DiscretizationTypes[] = {"FD", "FE", "DiscretizationType", "DISCRETIZATION_", 0};
34*9371c9d4SSatish Balay typedef enum {
35*9371c9d4SSatish Balay   QUADRATURE_GAUSS1,
36*9371c9d4SSatish Balay   QUADRATURE_GAUSS2,
37*9371c9d4SSatish Balay   QUADRATURE_GAUSS3,
38*9371c9d4SSatish Balay   QUADRATURE_GAUSS4,
39*9371c9d4SSatish Balay   QUADRATURE_LOBATTO2,
40*9371c9d4SSatish Balay   QUADRATURE_LOBATTO3
41*9371c9d4SSatish 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 
82*9371c9d4SSatish 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  */
97*9371c9d4SSatish Balay static void RDMaterialEnergy(RD rd, const RDNode *n, PetscScalar *Em, RDNode *dEm) {
98*9371c9d4SSatish Balay   rd->MaterialEnergy(rd, n, Em, dEm);
99*9371c9d4SSatish Balay }
100c4762a1bSJed Brown 
101c4762a1bSJed Brown /* Solves a quadratic equation while propagating tangents */
102*9371c9d4SSatish Balay static void QuadraticSolve(PetscScalar a, PetscScalar a_t, PetscScalar b, PetscScalar b_t, PetscScalar c, PetscScalar c_t, PetscScalar *x, PetscScalar *x_t) {
103*9371c9d4SSatish 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 */
104*9371c9d4SSatish 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 */
110*9371c9d4SSatish Balay static void RDMaterialEnergy_Saha(RD rd, const RDNode *n, PetscScalar *inEm, RDNode *dEm) {
111*9371c9d4SSatish 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 */
112*9371c9d4SSatish 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 */
122*9371c9d4SSatish Balay static void RDMaterialEnergy_Reduced(RD rd, const RDNode *n, PetscScalar *Em, RDNode *dEm) {
123*9371c9d4SSatish 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 */
133*9371c9d4SSatish 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 */
140*9371c9d4SSatish 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 
161*9371c9d4SSatish 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++) {
174*9371c9d4SSatish 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),
175*9371c9d4SSatish 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 
184*9371c9d4SSatish Balay static PetscScalar RDRadiation(RD rd, const RDNode *n, RDNode *dn) {
185*9371c9d4SSatish 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 
193*9371c9d4SSatish 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 
231*9371c9d4SSatish 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 */
256*9371c9d4SSatish 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 
264*9371c9d4SSatish 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 
275*9371c9d4SSatish 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 */
292*9371c9d4SSatish Balay #define RDCheckDomain(rd, ts, X) \
293*9371c9d4SSatish 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 
299*9371c9d4SSatish 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;
347*9371c9d4SSatish 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 
363*9371c9d4SSatish 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 
383*9371c9d4SSatish Balay     /*rad = (1.-Theta)* */ RDRadiation(rd, &x0[i], 0); /* + Theta* */
384*9371c9d4SSatish 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 */
470*9371c9d4SSatish 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;
472*9371c9d4SSatish Balay   n->E  = 0;
473*9371c9d4SSatish Balay   n->T  = 0;
474*9371c9d4SSatish Balay   nx->E = 0;
475*9371c9d4SSatish 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 */
487*9371c9d4SSatish 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*9371c9d4SSatish Balay     static const PetscReal ww[1] =
495*9371c9d4SSatish Balay       {
496*9371c9d4SSatish Balay         1.
497*9371c9d4SSatish Balay     },
498*9371c9d4SSatish Balay                            ii[1][2] = {{0.5, 0.5}}, dd[1][2] = {{-1., 1.}};
499*9371c9d4SSatish Balay     *nq       = 1;
500*9371c9d4SSatish Balay     refweight = ww;
501*9371c9d4SSatish Balay     refinterp = ii;
502*9371c9d4SSatish Balay     refderiv  = dd;
503c4762a1bSJed Brown   } break;
504c4762a1bSJed Brown   case QUADRATURE_GAUSS2: {
505*9371c9d4SSatish Balay     static const PetscReal ii[2][2] =
506*9371c9d4SSatish Balay       {
507*9371c9d4SSatish Balay         {0.78867513459481287, 0.21132486540518713},
508*9371c9d4SSatish Balay         {0.21132486540518713, 0.78867513459481287}
509*9371c9d4SSatish Balay     },
510*9371c9d4SSatish Balay                            dd[2][2] = {{-1., 1.}, {-1., 1.}}, ww[2] = {0.5, 0.5};
511*9371c9d4SSatish Balay     *nq       = 2;
512*9371c9d4SSatish Balay     refweight = ww;
513*9371c9d4SSatish Balay     refinterp = ii;
514*9371c9d4SSatish Balay     refderiv  = dd;
515c4762a1bSJed Brown   } break;
516c4762a1bSJed Brown   case QUADRATURE_GAUSS3: {
517*9371c9d4SSatish Balay     static const PetscReal ii[3][2] =
518*9371c9d4SSatish Balay       {
519*9371c9d4SSatish Balay         {0.8872983346207417, 0.1127016653792583},
520*9371c9d4SSatish Balay         {0.5,                0.5               },
521*9371c9d4SSatish Balay         {0.1127016653792583, 0.8872983346207417}
522*9371c9d4SSatish Balay     },
523c4762a1bSJed Brown                            dd[3][2] = {{-1, 1}, {-1, 1}, {-1, 1}}, ww[3] = {5. / 18, 8. / 18, 5. / 18};
524*9371c9d4SSatish Balay     *nq       = 3;
525*9371c9d4SSatish Balay     refweight = ww;
526*9371c9d4SSatish Balay     refinterp = ii;
527*9371c9d4SSatish Balay     refderiv  = dd;
528c4762a1bSJed Brown   } break;
529c4762a1bSJed Brown   case QUADRATURE_GAUSS4: {
530*9371c9d4SSatish Balay     static const PetscReal ii[][2] =
531*9371c9d4SSatish Balay       {
532*9371c9d4SSatish Balay         {0.93056815579702623,  0.069431844202973658},
533c4762a1bSJed Brown         {0.66999052179242813,  0.33000947820757187 },
534c4762a1bSJed Brown         {0.33000947820757187,  0.66999052179242813 },
535*9371c9d4SSatish Balay         {0.069431844202973658, 0.93056815579702623 }
536*9371c9d4SSatish Balay     },
537c4762a1bSJed Brown                            dd[][2] = {{-1, 1}, {-1, 1}, {-1, 1}, {-1, 1}}, ww[] = {0.17392742256872692, 0.3260725774312731, 0.3260725774312731, 0.17392742256872692};
538c4762a1bSJed Brown 
539*9371c9d4SSatish Balay     *nq       = 4;
540*9371c9d4SSatish Balay     refweight = ww;
541*9371c9d4SSatish Balay     refinterp = ii;
542*9371c9d4SSatish Balay     refderiv  = dd;
543c4762a1bSJed Brown   } break;
544c4762a1bSJed Brown   case QUADRATURE_LOBATTO2: {
545*9371c9d4SSatish Balay     static const PetscReal ii[2][2] =
546*9371c9d4SSatish Balay       {
547*9371c9d4SSatish Balay         {1., 0.},
548*9371c9d4SSatish Balay         {0., 1.}
549*9371c9d4SSatish Balay     },
550*9371c9d4SSatish Balay                            dd[2][2] = {{-1., 1.}, {-1., 1.}}, ww[2] = {0.5, 0.5};
551*9371c9d4SSatish Balay     *nq       = 2;
552*9371c9d4SSatish Balay     refweight = ww;
553*9371c9d4SSatish Balay     refinterp = ii;
554*9371c9d4SSatish Balay     refderiv  = dd;
555c4762a1bSJed Brown   } break;
556c4762a1bSJed Brown   case QUADRATURE_LOBATTO3: {
557*9371c9d4SSatish Balay     static const PetscReal ii[3][2] =
558*9371c9d4SSatish Balay       {
559*9371c9d4SSatish Balay         {1,   0  },
560*9371c9d4SSatish Balay         {0.5, 0.5},
561*9371c9d4SSatish Balay         {0,   1  }
562*9371c9d4SSatish Balay     },
563*9371c9d4SSatish Balay                            dd[3][2] = {{-1, 1}, {-1, 1}, {-1, 1}}, ww[3] = {1. / 6, 4. / 6, 1. / 6};
564*9371c9d4SSatish Balay     *nq       = 3;
565*9371c9d4SSatish Balay     refweight = ww;
566*9371c9d4SSatish Balay     refinterp = ii;
567*9371c9d4SSatish Balay     refderiv  = dd;
568c4762a1bSJed Brown   } break;
56998921bdaSJacob Faibussowitsch   default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unknown quadrature %d", (int)rd->quadrature);
570c4762a1bSJed Brown   }
571c4762a1bSJed Brown 
572c4762a1bSJed Brown   for (q = 0; q < *nq; q++) {
573c4762a1bSJed Brown     weight[q] = refweight[q] * hx;
574c4762a1bSJed Brown     for (j = 0; j < 2; j++) {
575c4762a1bSJed Brown       interp[q][j] = refinterp[q][j];
576c4762a1bSJed Brown       deriv[q][j]  = refderiv[q][j] / hx;
577c4762a1bSJed Brown     }
578c4762a1bSJed Brown   }
579c4762a1bSJed Brown   PetscFunctionReturn(0);
580c4762a1bSJed Brown }
581c4762a1bSJed Brown 
582c4762a1bSJed Brown /*
583c4762a1bSJed Brown  Finite element version
584c4762a1bSJed Brown */
585*9371c9d4SSatish Balay static PetscErrorCode RDIFunction_FE(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, void *ctx) {
586c4762a1bSJed Brown   RD            rd = (RD)ctx;
587c4762a1bSJed Brown   RDNode       *x, *x0, *xdot, *f;
588c4762a1bSJed Brown   Vec           X0loc, Xloc, Xloc_t, Floc;
589c4762a1bSJed Brown   PetscReal     hx, Theta, dt, weight[5], interp[5][2], deriv[5][2];
590c4762a1bSJed Brown   DMDALocalInfo info;
591c4762a1bSJed Brown   PetscInt      i, j, q, nq;
592c4762a1bSJed Brown 
593c4762a1bSJed Brown   PetscFunctionBeginUser;
5949566063dSJacob Faibussowitsch   PetscCall(RDGetLocalArrays(rd, ts, X, Xdot, &Theta, &dt, &X0loc, &x0, &Xloc, &x, &Xloc_t, &xdot));
595c4762a1bSJed Brown 
5969566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(rd->da, &Floc));
5979566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(Floc));
5989566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(rd->da, Floc, &f));
5999566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(rd->da, &info));
600c4762a1bSJed Brown 
601c4762a1bSJed Brown   /* Set up shape functions and quadrature for elements (assumes a uniform grid) */
602c4762a1bSJed Brown   hx = rd->L / (info.mx - 1);
6039566063dSJacob Faibussowitsch   PetscCall(RDGetQuadrature(rd, hx, &nq, weight, interp, deriv));
604c4762a1bSJed Brown 
605c4762a1bSJed Brown   for (i = info.xs; i < PetscMin(info.xs + info.xm, info.mx - 1); i++) {
606c4762a1bSJed Brown     for (q = 0; q < nq; q++) {
607c4762a1bSJed Brown       PetscReal   rho = rd->rho;
608c4762a1bSJed Brown       PetscScalar Em_t, rad, D_R, D0_R;
609c4762a1bSJed Brown       RDNode      n, n0, nx, n0x, nt, ntx;
610c4762a1bSJed Brown       RDEvaluate(interp, deriv, q, x, i, &n, &nx);
611c4762a1bSJed Brown       RDEvaluate(interp, deriv, q, x0, i, &n0, &n0x);
612c4762a1bSJed Brown       RDEvaluate(interp, deriv, q, xdot, i, &nt, &ntx);
613c4762a1bSJed Brown 
614c4762a1bSJed Brown       rad = (1. - Theta) * RDRadiation(rd, &n0, 0) + Theta * RDRadiation(rd, &n, 0);
615c4762a1bSJed Brown       if (rd->endpoint) {
616c4762a1bSJed Brown         PetscScalar Em0, Em1;
617c4762a1bSJed Brown         RDMaterialEnergy(rd, &n0, &Em0, NULL);
618c4762a1bSJed Brown         RDMaterialEnergy(rd, &n, &Em1, NULL);
619c4762a1bSJed Brown         Em_t = (Em1 - Em0) / dt;
620c4762a1bSJed Brown       } else {
621c4762a1bSJed Brown         RDNode dEm;
622c4762a1bSJed Brown         RDMaterialEnergy(rd, &n, NULL, &dEm);
623c4762a1bSJed Brown         Em_t = dEm.E * nt.E + dEm.T * nt.T;
624c4762a1bSJed Brown       }
625c4762a1bSJed Brown       RDDiffusionCoefficient(rd, PETSC_TRUE, &n0, &n0x, &D0_R, 0, 0);
626c4762a1bSJed Brown       RDDiffusionCoefficient(rd, PETSC_TRUE, &n, &nx, &D_R, 0, 0);
627c4762a1bSJed Brown       for (j = 0; j < 2; j++) {
628*9371c9d4SSatish 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));
629c4762a1bSJed Brown         f[i + j].T += interp[q][j] * weight[q] * (rho * Em_t + rad);
630c4762a1bSJed Brown       }
631c4762a1bSJed Brown     }
632c4762a1bSJed Brown   }
633c4762a1bSJed Brown   if (info.xs == 0) {
634c4762a1bSJed Brown     switch (rd->leftbc) {
635c4762a1bSJed Brown     case BC_ROBIN: {
636c4762a1bSJed Brown       PetscScalar D_R, D_R_bc;
637c4762a1bSJed Brown       PetscReal   ratio, bcTheta = rd->bcmidpoint ? Theta : 1.;
638c4762a1bSJed Brown       RDNode      n, nx;
639c4762a1bSJed Brown 
640c4762a1bSJed Brown       n.E  = (1 - bcTheta) * x0[0].E + bcTheta * x[0].E;
641c4762a1bSJed Brown       n.T  = (1 - bcTheta) * x0[0].T + bcTheta * x[0].T;
642c4762a1bSJed Brown       nx.E = (x[1].E - x[0].E) / hx;
643c4762a1bSJed Brown       nx.T = (x[1].T - x[0].T) / hx;
644c4762a1bSJed Brown       RDDiffusionCoefficient(rd, PETSC_TRUE, &n, &nx, &D_R, 0, 0);
645c4762a1bSJed Brown       RDDiffusionCoefficient(rd, rd->bclimit, &n, &nx, &D_R_bc, 0, 0);
646c4762a1bSJed Brown       ratio = PetscRealPart(D_R / D_R_bc);
6473c633725SBarry Smith       PetscCheck(ratio <= 1., PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Limited diffusivity is greater than unlimited");
6483c633725SBarry Smith       PetscCheck(ratio >= 1e-3, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Heavily limited diffusivity");
649c4762a1bSJed Brown       f[0].E += -ratio * 0.5 * (rd->Eapplied - n.E);
650c4762a1bSJed Brown     } break;
651c4762a1bSJed Brown     case BC_NEUMANN:
652c4762a1bSJed Brown       /* homogeneous Neumann is the natural condition */
653c4762a1bSJed Brown       break;
65463a3b9bcSJacob Faibussowitsch     default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Case %" PetscInt_FMT, rd->initial);
655c4762a1bSJed Brown     }
656c4762a1bSJed Brown   }
657c4762a1bSJed Brown 
6589566063dSJacob Faibussowitsch   PetscCall(RDRestoreLocalArrays(rd, &X0loc, &x0, &Xloc, &x, &Xloc_t, &xdot));
6599566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(rd->da, Floc, &f));
6609566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(F));
6619566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalBegin(rd->da, Floc, ADD_VALUES, F));
6629566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalEnd(rd->da, Floc, ADD_VALUES, F));
6639566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(rd->da, &Floc));
664c4762a1bSJed Brown 
6659566063dSJacob Faibussowitsch   if (rd->monitor_residual) PetscCall(RDStateView(rd, X, Xdot, F));
666c4762a1bSJed Brown   PetscFunctionReturn(0);
667c4762a1bSJed Brown }
668c4762a1bSJed Brown 
669*9371c9d4SSatish Balay static PetscErrorCode RDIJacobian_FE(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal a, Mat A, Mat B, void *ctx) {
670c4762a1bSJed Brown   RD            rd = (RD)ctx;
671c4762a1bSJed Brown   RDNode       *x, *x0, *xdot;
672c4762a1bSJed Brown   Vec           X0loc, Xloc, Xloc_t;
673c4762a1bSJed Brown   PetscReal     hx, Theta, dt, weight[5], interp[5][2], deriv[5][2];
674c4762a1bSJed Brown   DMDALocalInfo info;
675c4762a1bSJed Brown   PetscInt      i, j, k, q, nq;
676c4762a1bSJed Brown   PetscScalar   K[4][4];
677c4762a1bSJed Brown 
678c4762a1bSJed Brown   PetscFunctionBeginUser;
6799566063dSJacob Faibussowitsch   PetscCall(RDGetLocalArrays(rd, ts, X, Xdot, &Theta, &dt, &X0loc, &x0, &Xloc, &x, &Xloc_t, &xdot));
6809566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(rd->da, &info));
681c4762a1bSJed Brown   hx = rd->L / (info.mx - 1);
6829566063dSJacob Faibussowitsch   PetscCall(RDGetQuadrature(rd, hx, &nq, weight, interp, deriv));
6839566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(B));
684c4762a1bSJed Brown   for (i = info.xs; i < PetscMin(info.xs + info.xm, info.mx - 1); i++) {
685c4762a1bSJed Brown     PetscInt rc[2];
686c4762a1bSJed Brown 
687*9371c9d4SSatish Balay     rc[0] = i;
688*9371c9d4SSatish Balay     rc[1] = i + 1;
6899566063dSJacob Faibussowitsch     PetscCall(PetscMemzero(K, sizeof(K)));
690c4762a1bSJed Brown     for (q = 0; q < nq; q++) {
691c4762a1bSJed Brown       PetscScalar              D_R;
692c4762a1bSJed Brown       PETSC_UNUSED PetscScalar rad;
693c4762a1bSJed Brown       RDNode                   n, nx, nt, ntx, drad, dD_R, dxD_R, dEm;
694c4762a1bSJed Brown       RDEvaluate(interp, deriv, q, x, i, &n, &nx);
695c4762a1bSJed Brown       RDEvaluate(interp, deriv, q, xdot, i, &nt, &ntx);
696c4762a1bSJed Brown       rad = RDRadiation(rd, &n, &drad);
697c4762a1bSJed Brown       RDDiffusionCoefficient(rd, PETSC_TRUE, &n, &nx, &D_R, &dD_R, &dxD_R);
698c4762a1bSJed Brown       RDMaterialEnergy(rd, &n, NULL, &dEm);
699c4762a1bSJed Brown       for (j = 0; j < 2; j++) {
700c4762a1bSJed Brown         for (k = 0; k < 2; k++) {
701*9371c9d4SSatish 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]));
702*9371c9d4SSatish 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);
703c4762a1bSJed Brown           K[j * 2 + 1][k * 2 + 0] += interp[q][j] * weight[q] * drad.E * interp[q][k];
704c4762a1bSJed Brown           K[j * 2 + 1][k * 2 + 1] += interp[q][j] * weight[q] * (a * rd->rho * dEm.T + drad.T) * interp[q][k];
705c4762a1bSJed Brown         }
706c4762a1bSJed Brown       }
707c4762a1bSJed Brown     }
7089566063dSJacob Faibussowitsch     PetscCall(MatSetValuesBlocked(B, 2, rc, 2, rc, &K[0][0], ADD_VALUES));
709c4762a1bSJed Brown   }
710c4762a1bSJed Brown   if (info.xs == 0) {
711c4762a1bSJed Brown     switch (rd->leftbc) {
712c4762a1bSJed Brown     case BC_ROBIN: {
713c4762a1bSJed Brown       PetscScalar D_R, D_R_bc;
714c4762a1bSJed Brown       PetscReal   ratio;
715c4762a1bSJed Brown       RDNode      n, nx;
716c4762a1bSJed Brown 
717c4762a1bSJed Brown       n.E  = (1 - Theta) * x0[0].E + Theta * x[0].E;
718c4762a1bSJed Brown       n.T  = (1 - Theta) * x0[0].T + Theta * x[0].T;
719c4762a1bSJed Brown       nx.E = (x[1].E - x[0].E) / hx;
720c4762a1bSJed Brown       nx.T = (x[1].T - x[0].T) / hx;
721c4762a1bSJed Brown       RDDiffusionCoefficient(rd, PETSC_TRUE, &n, &nx, &D_R, 0, 0);
722c4762a1bSJed Brown       RDDiffusionCoefficient(rd, rd->bclimit, &n, &nx, &D_R_bc, 0, 0);
723c4762a1bSJed Brown       ratio = PetscRealPart(D_R / D_R_bc);
7249566063dSJacob Faibussowitsch       PetscCall(MatSetValue(B, 0, 0, ratio * 0.5, ADD_VALUES));
725c4762a1bSJed Brown     } break;
726c4762a1bSJed Brown     case BC_NEUMANN:
727c4762a1bSJed Brown       /* homogeneous Neumann is the natural condition */
728c4762a1bSJed Brown       break;
72963a3b9bcSJacob Faibussowitsch     default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Case %" PetscInt_FMT, rd->initial);
730c4762a1bSJed Brown     }
731c4762a1bSJed Brown   }
732c4762a1bSJed Brown 
7339566063dSJacob Faibussowitsch   PetscCall(RDRestoreLocalArrays(rd, &X0loc, &x0, &Xloc, &x, &Xloc_t, &xdot));
7349566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
7359566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
736c4762a1bSJed Brown   if (A != B) {
7379566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
7389566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
739c4762a1bSJed Brown   }
740c4762a1bSJed Brown   PetscFunctionReturn(0);
741c4762a1bSJed Brown }
742c4762a1bSJed Brown 
743c4762a1bSJed Brown /* Temperature that is in equilibrium with the radiation density */
744*9371c9d4SSatish Balay static PetscScalar RDRadiationTemperature(RD rd, PetscScalar E) {
745*9371c9d4SSatish Balay   return PetscPowScalar(E * rd->c / (4. * rd->sigma_b), 0.25);
746*9371c9d4SSatish Balay }
747c4762a1bSJed Brown 
748*9371c9d4SSatish Balay static PetscErrorCode RDInitialState(RD rd, Vec X) {
749c4762a1bSJed Brown   DMDALocalInfo info;
750c4762a1bSJed Brown   PetscInt      i;
751c4762a1bSJed Brown   RDNode       *x;
752c4762a1bSJed Brown 
753c4762a1bSJed Brown   PetscFunctionBeginUser;
7549566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(rd->da, &info));
7559566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(rd->da, X, &x));
756c4762a1bSJed Brown   for (i = info.xs; i < info.xs + info.xm; i++) {
757c4762a1bSJed Brown     PetscReal coord = i * rd->L / (info.mx - 1);
758c4762a1bSJed Brown     switch (rd->initial) {
759c4762a1bSJed Brown     case 1:
760c4762a1bSJed Brown       x[i].E = 0.001;
761c4762a1bSJed Brown       x[i].T = RDRadiationTemperature(rd, x[i].E);
762c4762a1bSJed Brown       break;
763c4762a1bSJed Brown     case 2:
764c4762a1bSJed Brown       x[i].E = 0.001 + 100. * PetscExpReal(-PetscSqr(coord / 0.1));
765c4762a1bSJed Brown       x[i].T = RDRadiationTemperature(rd, x[i].E);
766c4762a1bSJed Brown       break;
767c4762a1bSJed Brown     case 3:
768c4762a1bSJed Brown       x[i].E = 7.56e-2 * rd->unit.Joule / PetscPowScalarInt(rd->unit.meter, 3);
769c4762a1bSJed Brown       x[i].T = RDRadiationTemperature(rd, x[i].E);
770c4762a1bSJed Brown       break;
77163a3b9bcSJacob Faibussowitsch     default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "No initial state %" PetscInt_FMT, rd->initial);
772c4762a1bSJed Brown     }
773c4762a1bSJed Brown   }
7749566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(rd->da, X, &x));
775c4762a1bSJed Brown   PetscFunctionReturn(0);
776c4762a1bSJed Brown }
777c4762a1bSJed Brown 
778*9371c9d4SSatish Balay static PetscErrorCode RDView(RD rd, Vec X, PetscViewer viewer) {
779c4762a1bSJed Brown   Vec             Y;
780c4762a1bSJed Brown   const RDNode   *x;
781c4762a1bSJed Brown   PetscScalar    *y;
782c4762a1bSJed Brown   PetscInt        i, m, M;
783c4762a1bSJed Brown   const PetscInt *lx;
784c4762a1bSJed Brown   DM              da;
785c4762a1bSJed Brown   MPI_Comm        comm;
786c4762a1bSJed Brown 
787c4762a1bSJed Brown   PetscFunctionBeginUser;
788c4762a1bSJed Brown   /*
789c4762a1bSJed Brown     Create a DMDA (one dof per node, zero stencil width, same layout) to hold Trad
790c4762a1bSJed Brown     (radiation temperature).  It is not necessary to create a DMDA for this, but this way
791c4762a1bSJed Brown     output and visualization will have meaningful variable names and correct scales.
792c4762a1bSJed Brown   */
7939566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(rd->da, 0, &M, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0));
7949566063dSJacob Faibussowitsch   PetscCall(DMDAGetOwnershipRanges(rd->da, &lx, 0, 0));
7959566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)rd->da, &comm));
7969566063dSJacob Faibussowitsch   PetscCall(DMDACreate1d(comm, DM_BOUNDARY_NONE, M, 1, 0, lx, &da));
7979566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
7989566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
7999566063dSJacob Faibussowitsch   PetscCall(DMDASetUniformCoordinates(da, 0., rd->L, 0., 0., 0., 0.));
8009566063dSJacob Faibussowitsch   PetscCall(DMDASetFieldName(da, 0, "T_rad"));
8019566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(da, &Y));
802c4762a1bSJed Brown 
803c4762a1bSJed Brown   /* Compute the radiation temperature from the solution at each node */
8049566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(Y, &m));
8059566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, (const PetscScalar **)&x));
8069566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Y, &y));
807c4762a1bSJed Brown   for (i = 0; i < m; i++) y[i] = RDRadiationTemperature(rd, x[i].E);
8089566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, (const PetscScalar **)&x));
8099566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Y, &y));
810c4762a1bSJed Brown 
8119566063dSJacob Faibussowitsch   PetscCall(VecView(Y, viewer));
8129566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Y));
8139566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
814c4762a1bSJed Brown   PetscFunctionReturn(0);
815c4762a1bSJed Brown }
816c4762a1bSJed Brown 
817*9371c9d4SSatish Balay static PetscErrorCode RDTestDifferentiation(RD rd) {
818c4762a1bSJed Brown   MPI_Comm    comm;
819c4762a1bSJed Brown   RDNode      n, nx;
820c4762a1bSJed Brown   PetscScalar epsilon;
821c4762a1bSJed Brown 
822c4762a1bSJed Brown   PetscFunctionBeginUser;
8239566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)rd->da, &comm));
824c4762a1bSJed Brown   epsilon = 1e-8;
825c4762a1bSJed Brown   {
826c4762a1bSJed Brown     RDNode      dEm, fdEm;
827c4762a1bSJed Brown     PetscScalar T0 = 1000., T1 = T0 * (1. + epsilon), Em0, Em1;
828c4762a1bSJed Brown     n.E = 1.;
829c4762a1bSJed Brown     n.T = T0;
830c4762a1bSJed Brown     rd->MaterialEnergy(rd, &n, &Em0, &dEm);
831c4762a1bSJed Brown     n.E = 1. + epsilon;
832c4762a1bSJed Brown     n.T = T0;
833c4762a1bSJed Brown     rd->MaterialEnergy(rd, &n, &Em1, 0);
834c4762a1bSJed Brown     fdEm.E = (Em1 - Em0) / epsilon;
835c4762a1bSJed Brown     n.E    = 1.;
836c4762a1bSJed Brown     n.T    = T1;
837c4762a1bSJed Brown     rd->MaterialEnergy(rd, &n, &Em1, 0);
838c4762a1bSJed Brown     fdEm.T = (Em1 - Em0) / (T0 * epsilon);
839*9371c9d4SSatish 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),
840*9371c9d4SSatish Balay                           (double)PetscRealPart(dEm.T - fdEm.T)));
841c4762a1bSJed Brown   }
842c4762a1bSJed Brown   {
843c4762a1bSJed Brown     PetscScalar D0, D;
844c4762a1bSJed Brown     RDNode      dD, dxD, fdD, fdxD;
845*9371c9d4SSatish Balay     n.E  = 1.;
846*9371c9d4SSatish Balay     n.T  = 1.;
847*9371c9d4SSatish Balay     nx.E = 1.;
848*9371c9d4SSatish Balay     n.T  = 1.;
849c4762a1bSJed Brown     RDDiffusionCoefficient(rd, rd->bclimit, &n, &nx, &D0, &dD, &dxD);
850*9371c9d4SSatish Balay     n.E  = 1. + epsilon;
851*9371c9d4SSatish Balay     n.T  = 1.;
852*9371c9d4SSatish Balay     nx.E = 1.;
853*9371c9d4SSatish Balay     n.T  = 1.;
854*9371c9d4SSatish Balay     RDDiffusionCoefficient(rd, rd->bclimit, &n, &nx, &D, 0, 0);
855*9371c9d4SSatish Balay     fdD.E = (D - D0) / epsilon;
856*9371c9d4SSatish Balay     n.E   = 1;
857*9371c9d4SSatish Balay     n.T   = 1. + epsilon;
858*9371c9d4SSatish Balay     nx.E  = 1.;
859*9371c9d4SSatish Balay     n.T   = 1.;
860*9371c9d4SSatish Balay     RDDiffusionCoefficient(rd, rd->bclimit, &n, &nx, &D, 0, 0);
861*9371c9d4SSatish Balay     fdD.T = (D - D0) / epsilon;
862*9371c9d4SSatish Balay     n.E   = 1;
863*9371c9d4SSatish Balay     n.T   = 1.;
864*9371c9d4SSatish Balay     nx.E  = 1. + epsilon;
865*9371c9d4SSatish Balay     n.T   = 1.;
866*9371c9d4SSatish Balay     RDDiffusionCoefficient(rd, rd->bclimit, &n, &nx, &D, 0, 0);
867*9371c9d4SSatish Balay     fdxD.E = (D - D0) / epsilon;
868*9371c9d4SSatish Balay     n.E    = 1;
869*9371c9d4SSatish Balay     n.T    = 1.;
870*9371c9d4SSatish Balay     nx.E   = 1.;
871*9371c9d4SSatish Balay     n.T    = 1. + epsilon;
872*9371c9d4SSatish Balay     RDDiffusionCoefficient(rd, rd->bclimit, &n, &nx, &D, 0, 0);
873*9371c9d4SSatish Balay     fdxD.T = (D - D0) / epsilon;
874*9371c9d4SSatish 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),
875*9371c9d4SSatish Balay                           (double)PetscRealPart(dD.T - fdD.T)));
876*9371c9d4SSatish 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),
877*9371c9d4SSatish Balay                           (double)PetscRealPart(dxD.T - fdxD.T)));
878c4762a1bSJed Brown   }
879c4762a1bSJed Brown   {
880c4762a1bSJed Brown     PetscInt    i;
881c4762a1bSJed Brown     PetscReal   hx = 1.;
882c4762a1bSJed Brown     PetscScalar a0;
883c4762a1bSJed Brown     RDNode      n0[3], n1[3], d[3], fd[3];
884c4762a1bSJed Brown 
885c4762a1bSJed Brown     n0[0].E = 1.;
886c4762a1bSJed Brown     n0[0].T = 1.;
887c4762a1bSJed Brown     n0[1].E = 5.;
888c4762a1bSJed Brown     n0[1].T = 3.;
889c4762a1bSJed Brown     n0[2].E = 4.;
890c4762a1bSJed Brown     n0[2].T = 2.;
891c4762a1bSJed Brown     a0      = RDDiffusion(rd, hx, n0, 1, d);
892c4762a1bSJed Brown     for (i = 0; i < 3; i++) {
893*9371c9d4SSatish Balay       PetscCall(PetscMemcpy(n1, n0, sizeof(n0)));
894*9371c9d4SSatish Balay       n1[i].E += epsilon;
895c4762a1bSJed Brown       fd[i].E = (RDDiffusion(rd, hx, n1, 1, 0) - a0) / epsilon;
896*9371c9d4SSatish Balay       PetscCall(PetscMemcpy(n1, n0, sizeof(n0)));
897*9371c9d4SSatish Balay       n1[i].T += epsilon;
898c4762a1bSJed Brown       fd[i].T = (RDDiffusion(rd, hx, n1, 1, 0) - a0) / epsilon;
899*9371c9d4SSatish 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),
900*9371c9d4SSatish Balay                             (double)PetscRealPart(d[i].E - fd[i].E), (double)PetscRealPart(d[i].T - fd[i].T)));
901c4762a1bSJed Brown     }
902c4762a1bSJed Brown   }
903c4762a1bSJed Brown   {
904c4762a1bSJed Brown     PetscScalar rad0, rad;
905c4762a1bSJed Brown     RDNode      drad, fdrad;
906*9371c9d4SSatish Balay     n.E     = 1.;
907*9371c9d4SSatish Balay     n.T     = 1.;
908c4762a1bSJed Brown     rad0    = RDRadiation(rd, &n, &drad);
909*9371c9d4SSatish Balay     n.E     = 1. + epsilon;
910*9371c9d4SSatish Balay     n.T     = 1.;
911*9371c9d4SSatish Balay     rad     = RDRadiation(rd, &n, 0);
912*9371c9d4SSatish Balay     fdrad.E = (rad - rad0) / epsilon;
913*9371c9d4SSatish Balay     n.E     = 1.;
914*9371c9d4SSatish Balay     n.T     = 1. + epsilon;
915*9371c9d4SSatish Balay     rad     = RDRadiation(rd, &n, 0);
916*9371c9d4SSatish Balay     fdrad.T = (rad - rad0) / epsilon;
917*9371c9d4SSatish 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),
918*9371c9d4SSatish Balay                           (double)PetscRealPart(drad.T - fdrad.T)));
919c4762a1bSJed Brown   }
920c4762a1bSJed Brown   PetscFunctionReturn(0);
921c4762a1bSJed Brown }
922c4762a1bSJed Brown 
923*9371c9d4SSatish Balay static PetscErrorCode RDCreate(MPI_Comm comm, RD *inrd) {
924c4762a1bSJed Brown   RD        rd;
925c4762a1bSJed Brown   PetscReal meter = 0, kilogram = 0, second = 0, Kelvin = 0, Joule = 0, Watt = 0;
926c4762a1bSJed Brown 
927c4762a1bSJed Brown   PetscFunctionBeginUser;
928c4762a1bSJed Brown   *inrd = 0;
9299566063dSJacob Faibussowitsch   PetscCall(PetscNew(&rd));
930c4762a1bSJed Brown 
931d0609cedSBarry Smith   PetscOptionsBegin(comm, NULL, "Options for nonequilibrium radiation-diffusion with RD ionization", NULL);
932c4762a1bSJed Brown   {
933c4762a1bSJed Brown     rd->initial = 1;
9349566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-rd_initial", "Initial condition (1=Marshak, 2=Blast, 3=Marshak+)", "", rd->initial, &rd->initial, 0));
935c4762a1bSJed Brown     switch (rd->initial) {
936c4762a1bSJed Brown     case 1:
937c4762a1bSJed Brown     case 2:
938c4762a1bSJed Brown       rd->unit.kilogram = 1.;
939c4762a1bSJed Brown       rd->unit.meter    = 1.;
940c4762a1bSJed Brown       rd->unit.second   = 1.;
941c4762a1bSJed Brown       rd->unit.Kelvin   = 1.;
942c4762a1bSJed Brown       break;
943c4762a1bSJed Brown     case 3:
944c4762a1bSJed Brown       rd->unit.kilogram = 1.e12;
945c4762a1bSJed Brown       rd->unit.meter    = 1.;
946c4762a1bSJed Brown       rd->unit.second   = 1.e9;
947c4762a1bSJed Brown       rd->unit.Kelvin   = 1.;
948c4762a1bSJed Brown       break;
94963a3b9bcSJacob Faibussowitsch     default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unknown initial condition %" PetscInt_FMT, rd->initial);
950c4762a1bSJed Brown     }
951c4762a1bSJed Brown     /* Fundamental units */
9529566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-rd_unit_meter", "Length of 1 meter in nondimensional units", "", rd->unit.meter, &rd->unit.meter, 0));
9539566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-rd_unit_kilogram", "Mass of 1 kilogram in nondimensional units", "", rd->unit.kilogram, &rd->unit.kilogram, 0));
9549566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-rd_unit_second", "Time of a second in nondimensional units", "", rd->unit.second, &rd->unit.second, 0));
9559566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-rd_unit_Kelvin", "Temperature of a Kelvin in nondimensional units", "", rd->unit.Kelvin, &rd->unit.Kelvin, 0));
956c4762a1bSJed Brown     /* Derived units */
957c4762a1bSJed Brown     rd->unit.Joule = rd->unit.kilogram * PetscSqr(rd->unit.meter / rd->unit.second);
958c4762a1bSJed Brown     rd->unit.Watt  = rd->unit.Joule / rd->unit.second;
959c4762a1bSJed Brown     /* Local aliases */
960c4762a1bSJed Brown     meter          = rd->unit.meter;
961c4762a1bSJed Brown     kilogram       = rd->unit.kilogram;
962c4762a1bSJed Brown     second         = rd->unit.second;
963c4762a1bSJed Brown     Kelvin         = rd->unit.Kelvin;
964c4762a1bSJed Brown     Joule          = rd->unit.Joule;
965c4762a1bSJed Brown     Watt           = rd->unit.Watt;
966c4762a1bSJed Brown 
9679566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-rd_monitor_residual", "Display residuals every time they are evaluated", "", rd->monitor_residual, &rd->monitor_residual, NULL));
9689566063dSJacob Faibussowitsch     PetscCall(PetscOptionsEnum("-rd_discretization", "Discretization type", "", DiscretizationTypes, (PetscEnum)rd->discretization, (PetscEnum *)&rd->discretization, NULL));
969c4762a1bSJed Brown     if (rd->discretization == DISCRETIZATION_FE) {
970c4762a1bSJed Brown       rd->quadrature = QUADRATURE_GAUSS2;
9719566063dSJacob Faibussowitsch       PetscCall(PetscOptionsEnum("-rd_quadrature", "Finite element quadrature", "", QuadratureTypes, (PetscEnum)rd->quadrature, (PetscEnum *)&rd->quadrature, NULL));
972c4762a1bSJed Brown     }
9739566063dSJacob Faibussowitsch     PetscCall(PetscOptionsEnum("-rd_jacobian", "Type of finite difference Jacobian", "", JacobianTypes, (PetscEnum)rd->jacobian, (PetscEnum *)&rd->jacobian, NULL));
974c4762a1bSJed Brown     switch (rd->initial) {
975c4762a1bSJed Brown     case 1:
976c4762a1bSJed Brown       rd->leftbc     = BC_ROBIN;
977c4762a1bSJed Brown       rd->Eapplied   = 4 * rd->unit.Joule / PetscPowRealInt(rd->unit.meter, 3);
978c4762a1bSJed Brown       rd->L          = 1. * rd->unit.meter;
979c4762a1bSJed Brown       rd->beta       = 3.0;
980c4762a1bSJed Brown       rd->gamma      = 3.0;
981c4762a1bSJed Brown       rd->final_time = 3 * second;
982c4762a1bSJed Brown       break;
983c4762a1bSJed Brown     case 2:
984c4762a1bSJed Brown       rd->leftbc     = BC_NEUMANN;
985c4762a1bSJed Brown       rd->Eapplied   = 0.;
986c4762a1bSJed Brown       rd->L          = 1. * rd->unit.meter;
987c4762a1bSJed Brown       rd->beta       = 3.0;
988c4762a1bSJed Brown       rd->gamma      = 3.0;
989c4762a1bSJed Brown       rd->final_time = 1 * second;
990c4762a1bSJed Brown       break;
991c4762a1bSJed Brown     case 3:
992c4762a1bSJed Brown       rd->leftbc     = BC_ROBIN;
993c4762a1bSJed Brown       rd->Eapplied   = 7.503e6 * rd->unit.Joule / PetscPowRealInt(rd->unit.meter, 3);
994c4762a1bSJed Brown       rd->L          = 5. * rd->unit.meter;
995c4762a1bSJed Brown       rd->beta       = 3.5;
996c4762a1bSJed Brown       rd->gamma      = 3.5;
997c4762a1bSJed Brown       rd->final_time = 20e-9 * second;
998c4762a1bSJed Brown       break;
99963a3b9bcSJacob Faibussowitsch     default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Initial %" PetscInt_FMT, rd->initial);
1000c4762a1bSJed Brown     }
10019566063dSJacob Faibussowitsch     PetscCall(PetscOptionsEnum("-rd_leftbc", "Left boundary condition", "", BCTypes, (PetscEnum)rd->leftbc, (PetscEnum *)&rd->leftbc, NULL));
10029566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-rd_E_applied", "Radiation flux at left end of domain", "", rd->Eapplied, &rd->Eapplied, NULL));
10039566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-rd_beta", "Thermal exponent for photon absorption", "", rd->beta, &rd->beta, NULL));
10049566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-rd_gamma", "Thermal exponent for diffusion coefficient", "", rd->gamma, &rd->gamma, NULL));
10059566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-rd_view_draw", "Draw final solution", "", rd->view_draw, &rd->view_draw, NULL));
10069566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-rd_endpoint", "Discretize using endpoints (like trapezoid rule) instead of midpoint", "", rd->endpoint, &rd->endpoint, NULL));
10079566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-rd_bcmidpoint", "Impose the boundary condition at the midpoint (Theta) of the interval", "", rd->bcmidpoint, &rd->bcmidpoint, NULL));
10089566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-rd_bclimit", "Limit diffusion coefficient in definition of Robin boundary condition", "", rd->bclimit, &rd->bclimit, NULL));
10099566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-rd_test_diff", "Test differentiation in constitutive relations", "", rd->test_diff, &rd->test_diff, NULL));
10109566063dSJacob Faibussowitsch     PetscCall(PetscOptionsString("-rd_view_binary", "File name to hold final solution", "", rd->view_binary, rd->view_binary, sizeof(rd->view_binary), NULL));
1011c4762a1bSJed Brown   }
1012d0609cedSBarry Smith   PetscOptionsEnd();
1013c4762a1bSJed Brown 
1014c4762a1bSJed Brown   switch (rd->initial) {
1015c4762a1bSJed Brown   case 1:
1016c4762a1bSJed Brown   case 2:
1017c4762a1bSJed Brown     rd->rho            = 1.;
1018c4762a1bSJed Brown     rd->c              = 1.;
1019c4762a1bSJed Brown     rd->K_R            = 1.;
1020c4762a1bSJed Brown     rd->K_p            = 1.;
1021c4762a1bSJed Brown     rd->sigma_b        = 0.25;
1022c4762a1bSJed Brown     rd->MaterialEnergy = RDMaterialEnergy_Reduced;
1023c4762a1bSJed Brown     break;
1024c4762a1bSJed Brown   case 3:
1025c4762a1bSJed Brown     /* Table 2 */
1026c4762a1bSJed Brown     rd->rho            = 1.17e-3 * kilogram / (meter * meter * meter);                                                    /* density */
1027c4762a1bSJed Brown     rd->K_R            = 7.44e18 * PetscPowRealInt(meter, 5) * PetscPowReal(Kelvin, 3.5) * PetscPowRealInt(kilogram, -2); /*  */
1028c4762a1bSJed Brown     rd->K_p            = 2.33e20 * PetscPowRealInt(meter, 5) * PetscPowReal(Kelvin, 3.5) * PetscPowRealInt(kilogram, -2); /*  */
1029c4762a1bSJed Brown     rd->I_H            = 2.179e-18 * Joule;                                                                               /* Hydrogen ionization potential */
1030c4762a1bSJed Brown     rd->m_p            = 1.673e-27 * kilogram;                                                                            /* proton mass */
1031c4762a1bSJed Brown     rd->m_e            = 9.109e-31 * kilogram;                                                                            /* electron mass */
1032c4762a1bSJed Brown     rd->h              = 6.626e-34 * Joule * second;                                                                      /* Planck's constant */
1033c4762a1bSJed Brown     rd->k              = 1.381e-23 * Joule / Kelvin;                                                                      /* Boltzman constant */
1034c4762a1bSJed Brown     rd->c              = 3.00e8 * meter / second;                                                                         /* speed of light */
1035c4762a1bSJed Brown     rd->sigma_b        = 5.67e-8 * Watt * PetscPowRealInt(meter, -2) * PetscPowRealInt(Kelvin, -4);                       /* Stefan-Boltzman constant */
1036c4762a1bSJed Brown     rd->MaterialEnergy = RDMaterialEnergy_Saha;
1037c4762a1bSJed Brown     break;
1038c4762a1bSJed Brown   }
1039c4762a1bSJed Brown 
10409566063dSJacob Faibussowitsch   PetscCall(DMDACreate1d(comm, DM_BOUNDARY_NONE, 20, sizeof(RDNode) / sizeof(PetscScalar), 1, NULL, &rd->da));
10419566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(rd->da));
10429566063dSJacob Faibussowitsch   PetscCall(DMSetUp(rd->da));
10439566063dSJacob Faibussowitsch   PetscCall(DMDASetFieldName(rd->da, 0, "E"));
10449566063dSJacob Faibussowitsch   PetscCall(DMDASetFieldName(rd->da, 1, "T"));
10459566063dSJacob Faibussowitsch   PetscCall(DMDASetUniformCoordinates(rd->da, 0., 1., 0., 0., 0., 0.));
1046c4762a1bSJed Brown 
1047c4762a1bSJed Brown   *inrd = rd;
1048c4762a1bSJed Brown   PetscFunctionReturn(0);
1049c4762a1bSJed Brown }
1050c4762a1bSJed Brown 
1051*9371c9d4SSatish Balay int main(int argc, char *argv[]) {
1052c4762a1bSJed Brown   RD        rd;
1053c4762a1bSJed Brown   TS        ts;
1054c4762a1bSJed Brown   SNES      snes;
1055c4762a1bSJed Brown   Vec       X;
1056c4762a1bSJed Brown   Mat       A, B;
1057c4762a1bSJed Brown   PetscInt  steps;
1058c4762a1bSJed Brown   PetscReal ftime;
1059c4762a1bSJed Brown 
1060327415f7SBarry Smith   PetscFunctionBeginUser;
10619ded082cSBarry Smith   PetscCall(PetscInitialize(&argc, &argv, 0, NULL));
10629566063dSJacob Faibussowitsch   PetscCall(RDCreate(PETSC_COMM_WORLD, &rd));
10639566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(rd->da, &X));
10649566063dSJacob Faibussowitsch   PetscCall(DMSetMatType(rd->da, MATAIJ));
10659566063dSJacob Faibussowitsch   PetscCall(DMCreateMatrix(rd->da, &B));
10669566063dSJacob Faibussowitsch   PetscCall(RDInitialState(rd, X));
1067c4762a1bSJed Brown 
10689566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
10699566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
10709566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts, TSTHETA));
10719566063dSJacob Faibussowitsch   PetscCall(TSSetDM(ts, rd->da));
1072c4762a1bSJed Brown   switch (rd->discretization) {
1073c4762a1bSJed Brown   case DISCRETIZATION_FD:
10749566063dSJacob Faibussowitsch     PetscCall(TSSetIFunction(ts, NULL, RDIFunction_FD, rd));
10759566063dSJacob Faibussowitsch     if (rd->jacobian == JACOBIAN_ANALYTIC) PetscCall(TSSetIJacobian(ts, B, B, RDIJacobian_FD, rd));
1076c4762a1bSJed Brown     break;
1077c4762a1bSJed Brown   case DISCRETIZATION_FE:
10789566063dSJacob Faibussowitsch     PetscCall(TSSetIFunction(ts, NULL, RDIFunction_FE, rd));
10799566063dSJacob Faibussowitsch     if (rd->jacobian == JACOBIAN_ANALYTIC) PetscCall(TSSetIJacobian(ts, B, B, RDIJacobian_FE, rd));
1080c4762a1bSJed Brown     break;
1081c4762a1bSJed Brown   }
10829566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, rd->final_time));
10839566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, 1e-3));
10849566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
10859566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
1086c4762a1bSJed Brown 
1087c4762a1bSJed Brown   A = B;
10889566063dSJacob Faibussowitsch   PetscCall(TSGetSNES(ts, &snes));
1089c4762a1bSJed Brown   switch (rd->jacobian) {
1090*9371c9d4SSatish Balay   case JACOBIAN_ANALYTIC: break;
1091*9371c9d4SSatish Balay   case JACOBIAN_MATRIXFREE: break;
1092c4762a1bSJed Brown   case JACOBIAN_FD_COLORING: {
10939566063dSJacob Faibussowitsch     PetscCall(SNESSetJacobian(snes, A, B, SNESComputeJacobianDefaultColor, 0));
1094c4762a1bSJed Brown   } break;
1095*9371c9d4SSatish Balay   case JACOBIAN_FD_FULL: PetscCall(SNESSetJacobian(snes, A, B, SNESComputeJacobianDefault, ts)); break;
1096c4762a1bSJed Brown   }
1097c4762a1bSJed Brown 
10981baa6e33SBarry Smith   if (rd->test_diff) PetscCall(RDTestDifferentiation(rd));
10999566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, X));
11009566063dSJacob Faibussowitsch   PetscCall(TSGetSolveTime(ts, &ftime));
11019566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts, &steps));
110263a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Steps %" PetscInt_FMT "  final time %g\n", steps, (double)ftime));
11031baa6e33SBarry Smith   if (rd->view_draw) PetscCall(RDView(rd, X, PETSC_VIEWER_DRAW_WORLD));
1104c4762a1bSJed Brown   if (rd->view_binary[0]) {
1105c4762a1bSJed Brown     PetscViewer viewer;
11069566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, rd->view_binary, FILE_MODE_WRITE, &viewer));
11079566063dSJacob Faibussowitsch     PetscCall(RDView(rd, X, viewer));
11089566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewer));
1109c4762a1bSJed Brown   }
11109566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&X));
11119566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
11129566063dSJacob Faibussowitsch   PetscCall(RDDestroy(&rd));
11139566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
11149566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
1115b122ec5aSJacob Faibussowitsch   return 0;
1116c4762a1bSJed Brown }
1117c4762a1bSJed Brown /*TEST
1118c4762a1bSJed Brown 
1119c4762a1bSJed Brown     test:
1120c4762a1bSJed 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
1121c4762a1bSJed Brown       requires: !single
1122c4762a1bSJed Brown 
1123c4762a1bSJed Brown     test:
1124c4762a1bSJed Brown       suffix: 2
1125c4762a1bSJed 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
1126c4762a1bSJed Brown       requires: !single
1127c4762a1bSJed Brown 
1128c4762a1bSJed Brown     test:
1129c4762a1bSJed Brown       suffix: 3
1130c4762a1bSJed Brown       nsize: 2
1131c4762a1bSJed 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
1132c4762a1bSJed Brown       requires: !single
1133c4762a1bSJed Brown 
1134c4762a1bSJed Brown TEST*/
1135