xref: /petsc/src/ts/tutorials/ex10.c (revision 3f02e49b19195914bf17f317a25cb39636853415)
1 /*
2   This example implements the model described in
3 
4     Rauenzahn, Mousseau, Knoll. "Temporal accuracy of the nonequilibrium radiation diffusion
5     equations employing a Saha ionization model" 2005.
6 
7   The paper discusses three examples, the first two are nondimensional with a simple
8   ionization model.  The third example is fully dimensional and uses the Saha ionization
9   model with realistic parameters.
10 */
11 
12 #include <petscts.h>
13 #include <petscdm.h>
14 #include <petscdmda.h>
15 
16 typedef enum {
17   BC_DIRICHLET,
18   BC_NEUMANN,
19   BC_ROBIN
20 } BCType;
21 static const char *const BCTypes[] = {"DIRICHLET", "NEUMANN", "ROBIN", "BCType", "BC_", 0};
22 typedef enum {
23   JACOBIAN_ANALYTIC,
24   JACOBIAN_MATRIXFREE,
25   JACOBIAN_FD_COLORING,
26   JACOBIAN_FD_FULL
27 } JacobianType;
28 static const char *const JacobianTypes[] = {"ANALYTIC", "MATRIXFREE", "FD_COLORING", "FD_FULL", "JacobianType", "FD_", 0};
29 typedef enum {
30   DISCRETIZATION_FD,
31   DISCRETIZATION_FE
32 } DiscretizationType;
33 static const char *const DiscretizationTypes[] = {"FD", "FE", "DiscretizationType", "DISCRETIZATION_", 0};
34 typedef enum {
35   QUADRATURE_GAUSS1,
36   QUADRATURE_GAUSS2,
37   QUADRATURE_GAUSS3,
38   QUADRATURE_GAUSS4,
39   QUADRATURE_LOBATTO2,
40   QUADRATURE_LOBATTO3
41 } QuadratureType;
42 static const char *const QuadratureTypes[] = {"GAUSS1", "GAUSS2", "GAUSS3", "GAUSS4", "LOBATTO2", "LOBATTO3", "QuadratureType", "QUADRATURE_", 0};
43 
44 typedef struct {
45   PetscScalar E; /* radiation energy */
46   PetscScalar T; /* material temperature */
47 } RDNode;
48 
49 typedef struct {
50   PetscReal meter, kilogram, second, Kelvin; /* Fundamental units */
51   PetscReal Joule, Watt;                     /* Derived units */
52 } RDUnit;
53 
54 typedef struct _n_RD *RD;
55 
56 struct _n_RD {
57   void (*MaterialEnergy)(RD, const RDNode *, PetscScalar *, RDNode *);
58   DM                 da;
59   PetscBool          monitor_residual;
60   DiscretizationType discretization;
61   QuadratureType     quadrature;
62   JacobianType       jacobian;
63   PetscInt           initial;
64   BCType             leftbc;
65   PetscBool          view_draw;
66   char               view_binary[PETSC_MAX_PATH_LEN];
67   PetscBool          test_diff;
68   PetscBool          endpoint;
69   PetscBool          bclimit;
70   PetscBool          bcmidpoint;
71   RDUnit             unit;
72 
73   /* model constants, see Table 2 and RDCreate() */
74   PetscReal rho, K_R, K_p, I_H, m_p, m_e, h, k, c, sigma_b, beta, gamma;
75 
76   /* Domain and boundary conditions */
77   PetscReal Eapplied; /* Radiation flux from the left */
78   PetscReal L;        /* Length of domain */
79   PetscReal final_time;
80 };
81 
82 static PetscErrorCode RDDestroy(RD *rd)
83 {
84   PetscFunctionBeginUser;
85   PetscCall(DMDestroy(&(*rd)->da));
86   PetscCall(PetscFree(*rd));
87   PetscFunctionReturn(PETSC_SUCCESS);
88 }
89 
90 /* The paper has a time derivative for material energy (Eq 2) which is a dependent variable (computable from temperature
91  * and density through an uninvertible relation).  Computing this derivative is trivial for trapezoid rule (used in the
92  * paper), but does not generalize nicely to higher order integrators.  Here we use the implicit form which provides
93  * time derivatives of the independent variables (radiation energy and temperature), so we must compute the time
94  * derivative of material energy ourselves (could be done using AD).
95  *
96  * There are multiple ionization models, this interface dispatches to the one currently in use.
97  */
98 static void RDMaterialEnergy(RD rd, const RDNode *n, PetscScalar *Em, RDNode *dEm)
99 {
100   rd->MaterialEnergy(rd, n, Em, dEm);
101 }
102 
103 /* Solves a quadratic equation while propagating tangents */
104 static void QuadraticSolve(PetscScalar a, PetscScalar a_t, PetscScalar b, PetscScalar b_t, PetscScalar c, PetscScalar c_t, PetscScalar *x, PetscScalar *x_t)
105 {
106   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 */
107     num_t = -b_t + 0.5 / PetscSqrtScalar(disc) * disc_t, den = 2. * a, den_t = 2. * a_t;
108   *x   = num / den;
109   *x_t = (num_t * den - num * den_t) / PetscSqr(den);
110 }
111 
112 /* The primary model presented in the paper */
113 static void RDMaterialEnergy_Saha(RD rd, const RDNode *n, PetscScalar *inEm, RDNode *dEm)
114 {
115   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 */
116     b_t = -b * chi_t + 1.5 * b / chi * chi_t, c = -b, c_t = -b_t;
117   QuadraticSolve(a, a_t, b, b_t, c, c_t, &alpha, &alpha_t);      /* Solve Eq 7 for alpha */
118   Em = rd->k * T / rd->m_p * (1.5 * (1. + alpha) + alpha * chi); /* Eq 6 */
119   if (inEm) *inEm = Em;
120   if (dEm) {
121     dEm->E = 0;
122     dEm->T = Em / T * T_t + rd->k * T / rd->m_p * (1.5 * alpha_t + alpha_t * chi + alpha * chi_t);
123   }
124 }
125 /* Reduced ionization model, Eq 30 */
126 static void RDMaterialEnergy_Reduced(RD rd, const RDNode *n, PetscScalar *Em, RDNode *dEm)
127 {
128   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;
129   QuadraticSolve(a, a_t, b, b_t, c, c_t, &alpha, &alpha_t);
130   if (Em) *Em = (1. + alpha) * T + 0.3 * alpha;
131   if (dEm) {
132     dEm->E = 0;
133     dEm->T = alpha_t * T + (1. + alpha) * T_t + 0.3 * alpha_t;
134   }
135 }
136 
137 /* Eq 5 */
138 static void RDSigma_R(RD rd, RDNode *n, PetscScalar *sigma_R, RDNode *dsigma_R)
139 {
140   *sigma_R    = rd->K_R * rd->rho * PetscPowScalar(n->T, -rd->gamma);
141   dsigma_R->E = 0;
142   dsigma_R->T = -rd->gamma * (*sigma_R) / n->T;
143 }
144 
145 /* Eq 4 */
146 static void RDDiffusionCoefficient(RD rd, PetscBool limit, RDNode *n, RDNode *nx, PetscScalar *D_R, RDNode *dD_R, RDNode *dxD_R)
147 {
148   PetscScalar sigma_R, denom;
149   RDNode      dsigma_R, ddenom, dxdenom;
150 
151   RDSigma_R(rd, n, &sigma_R, &dsigma_R);
152   denom     = 3. * rd->rho * sigma_R + (int)limit * PetscAbsScalar(nx->E) / n->E;
153   ddenom.E  = -(int)limit * PetscAbsScalar(nx->E) / PetscSqr(n->E);
154   ddenom.T  = 3. * rd->rho * dsigma_R.T;
155   dxdenom.E = (int)limit * (PetscRealPart(nx->E) < 0 ? -1. : 1.) / n->E;
156   dxdenom.T = 0;
157   *D_R      = rd->c / denom;
158   if (dD_R) {
159     dD_R->E = -rd->c / PetscSqr(denom) * ddenom.E;
160     dD_R->T = -rd->c / PetscSqr(denom) * ddenom.T;
161   }
162   if (dxD_R) {
163     dxD_R->E = -rd->c / PetscSqr(denom) * dxdenom.E;
164     dxD_R->T = -rd->c / PetscSqr(denom) * dxdenom.T;
165   }
166 }
167 
168 static PetscErrorCode RDStateView(RD rd, Vec X, Vec Xdot, Vec F)
169 {
170   DMDALocalInfo info;
171   PetscInt      i;
172   const RDNode *x, *xdot, *f;
173   MPI_Comm      comm;
174 
175   PetscFunctionBeginUser;
176   PetscCall(PetscObjectGetComm((PetscObject)rd->da, &comm));
177   PetscCall(DMDAGetLocalInfo(rd->da, &info));
178   PetscCall(DMDAVecGetArrayRead(rd->da, X, (void *)&x));
179   PetscCall(DMDAVecGetArrayRead(rd->da, Xdot, (void *)&xdot));
180   PetscCall(DMDAVecGetArrayRead(rd->da, F, (void *)&f));
181   for (i = info.xs; i < info.xs + info.xm; i++) {
182     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),
183                                       (double)PetscRealPart(f[i].E), (double)PetscRealPart(f[i].T)));
184   }
185   PetscCall(DMDAVecRestoreArrayRead(rd->da, X, (void *)&x));
186   PetscCall(DMDAVecRestoreArrayRead(rd->da, Xdot, (void *)&xdot));
187   PetscCall(DMDAVecRestoreArrayRead(rd->da, F, (void *)&f));
188   PetscCall(PetscSynchronizedFlush(comm, PETSC_STDOUT));
189   PetscFunctionReturn(PETSC_SUCCESS);
190 }
191 
192 static PetscScalar RDRadiation(RD rd, const RDNode *n, RDNode *dn)
193 {
194   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);
195   if (dn) {
196     dn->E = rad_E;
197     dn->T = rad_T;
198   }
199   return rad;
200 }
201 
202 static PetscScalar RDDiffusion(RD rd, PetscReal hx, const RDNode x[], PetscInt i, RDNode d[])
203 {
204   PetscReal   ihx = 1. / hx;
205   RDNode      n_L, nx_L, n_R, nx_R, dD_L, dxD_L, dD_R, dxD_R, dfluxL[2], dfluxR[2];
206   PetscScalar D_L, D_R, fluxL, fluxR;
207 
208   n_L.E  = 0.5 * (x[i - 1].E + x[i].E);
209   n_L.T  = 0.5 * (x[i - 1].T + x[i].T);
210   nx_L.E = (x[i].E - x[i - 1].E) / hx;
211   nx_L.T = (x[i].T - x[i - 1].T) / hx;
212   RDDiffusionCoefficient(rd, PETSC_TRUE, &n_L, &nx_L, &D_L, &dD_L, &dxD_L);
213   fluxL       = D_L * nx_L.E;
214   dfluxL[0].E = -ihx * D_L + (0.5 * dD_L.E - ihx * dxD_L.E) * nx_L.E;
215   dfluxL[1].E = +ihx * D_L + (0.5 * dD_L.E + ihx * dxD_L.E) * nx_L.E;
216   dfluxL[0].T = (0.5 * dD_L.T - ihx * dxD_L.T) * nx_L.E;
217   dfluxL[1].T = (0.5 * dD_L.T + ihx * dxD_L.T) * nx_L.E;
218 
219   n_R.E  = 0.5 * (x[i].E + x[i + 1].E);
220   n_R.T  = 0.5 * (x[i].T + x[i + 1].T);
221   nx_R.E = (x[i + 1].E - x[i].E) / hx;
222   nx_R.T = (x[i + 1].T - x[i].T) / hx;
223   RDDiffusionCoefficient(rd, PETSC_TRUE, &n_R, &nx_R, &D_R, &dD_R, &dxD_R);
224   fluxR       = D_R * nx_R.E;
225   dfluxR[0].E = -ihx * D_R + (0.5 * dD_R.E - ihx * dxD_R.E) * nx_R.E;
226   dfluxR[1].E = +ihx * D_R + (0.5 * dD_R.E + ihx * dxD_R.E) * nx_R.E;
227   dfluxR[0].T = (0.5 * dD_R.T - ihx * dxD_R.T) * nx_R.E;
228   dfluxR[1].T = (0.5 * dD_R.T + ihx * dxD_R.T) * nx_R.E;
229 
230   if (d) {
231     d[0].E = -ihx * dfluxL[0].E;
232     d[0].T = -ihx * dfluxL[0].T;
233     d[1].E = ihx * (dfluxR[0].E - dfluxL[1].E);
234     d[1].T = ihx * (dfluxR[0].T - dfluxL[1].T);
235     d[2].E = ihx * dfluxR[1].E;
236     d[2].T = ihx * dfluxR[1].T;
237   }
238   return ihx * (fluxR - fluxL);
239 }
240 
241 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)
242 {
243   PetscBool istheta;
244 
245   PetscFunctionBeginUser;
246   PetscCall(DMGetLocalVector(rd->da, X0loc));
247   PetscCall(DMGetLocalVector(rd->da, Xloc));
248   PetscCall(DMGetLocalVector(rd->da, Xloc_t));
249 
250   PetscCall(DMGlobalToLocalBegin(rd->da, X, INSERT_VALUES, *Xloc));
251   PetscCall(DMGlobalToLocalEnd(rd->da, X, INSERT_VALUES, *Xloc));
252   PetscCall(DMGlobalToLocalBegin(rd->da, Xdot, INSERT_VALUES, *Xloc_t));
253   PetscCall(DMGlobalToLocalEnd(rd->da, Xdot, INSERT_VALUES, *Xloc_t));
254 
255   /*
256     The following is a hack to subvert TSTHETA which is like an implicit midpoint method to behave more like a trapezoid
257     rule.  These methods have equivalent linear stability, but the nonlinear stability is somewhat different.  The
258     radiation system is inconvenient to write in explicit form because the ionization model is "on the left".
259    */
260   PetscCall(PetscObjectTypeCompare((PetscObject)ts, TSTHETA, &istheta));
261   if (istheta && rd->endpoint) PetscCall(TSThetaGetTheta(ts, Theta));
262   else *Theta = 1.;
263 
264   PetscCall(TSGetTimeStep(ts, dt));
265   PetscCall(VecWAXPY(*X0loc, -(*Theta) * (*dt), *Xloc_t, *Xloc));     /* back out the value at the start of this step */
266   if (rd->endpoint) PetscCall(VecWAXPY(*Xloc, *dt, *Xloc_t, *X0loc)); /* move the abscissa to the end of the step */
267 
268   PetscCall(DMDAVecGetArray(rd->da, *X0loc, x0));
269   PetscCall(DMDAVecGetArray(rd->da, *Xloc, x));
270   PetscCall(DMDAVecGetArray(rd->da, *Xloc_t, xdot));
271   PetscFunctionReturn(PETSC_SUCCESS);
272 }
273 
274 static PetscErrorCode RDRestoreLocalArrays(RD rd, Vec *X0loc, RDNode **x0, Vec *Xloc, RDNode **x, Vec *Xloc_t, RDNode **xdot)
275 {
276   PetscFunctionBeginUser;
277   PetscCall(DMDAVecRestoreArray(rd->da, *X0loc, x0));
278   PetscCall(DMDAVecRestoreArray(rd->da, *Xloc, x));
279   PetscCall(DMDAVecRestoreArray(rd->da, *Xloc_t, xdot));
280   PetscCall(DMRestoreLocalVector(rd->da, X0loc));
281   PetscCall(DMRestoreLocalVector(rd->da, Xloc));
282   PetscCall(DMRestoreLocalVector(rd->da, Xloc_t));
283   PetscFunctionReturn(PETSC_SUCCESS);
284 }
285 
286 static PetscErrorCode PETSC_UNUSED RDCheckDomain_Private(RD rd, TS ts, Vec X, PetscBool *in)
287 {
288   PetscInt  minloc;
289   PetscReal min;
290 
291   PetscFunctionBeginUser;
292   PetscCall(VecMin(X, &minloc, &min));
293   if (min < 0) {
294     SNES snes;
295     *in = PETSC_FALSE;
296     PetscCall(TSGetSNES(ts, &snes));
297     PetscCall(SNESSetFunctionDomainError(snes));
298     PetscCall(PetscInfo(ts, "Domain violation at %" PetscInt_FMT " field %" PetscInt_FMT " value %g\n", minloc / 2, minloc % 2, (double)min));
299   } else *in = PETSC_TRUE;
300   PetscFunctionReturn(PETSC_SUCCESS);
301 }
302 
303 /* Energy and temperature must remain positive */
304 #define RDCheckDomain(rd, ts, X) \
305   do { \
306     PetscBool _in; \
307     PetscCall(RDCheckDomain_Private(rd, ts, X, &_in)); \
308     if (!_in) PetscFunctionReturn(PETSC_SUCCESS); \
309   } while (0)
310 
311 static PetscErrorCode RDIFunction_FD(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, void *ctx)
312 {
313   RD            rd = (RD)ctx;
314   RDNode       *x, *x0, *xdot, *f;
315   Vec           X0loc, Xloc, Xloc_t;
316   PetscReal     hx, Theta, dt;
317   DMDALocalInfo info;
318   PetscInt      i;
319 
320   PetscFunctionBeginUser;
321   PetscCall(RDGetLocalArrays(rd, ts, X, Xdot, &Theta, &dt, &X0loc, &x0, &Xloc, &x, &Xloc_t, &xdot));
322   PetscCall(DMDAVecGetArray(rd->da, F, &f));
323   PetscCall(DMDAGetLocalInfo(rd->da, &info));
324   PetscCall(VecZeroEntries(F));
325 
326   hx = rd->L / (info.mx - 1);
327 
328   for (i = info.xs; i < info.xs + info.xm; i++) {
329     PetscReal   rho = rd->rho;
330     PetscScalar Em_t, rad;
331 
332     rad = (1. - Theta) * RDRadiation(rd, &x0[i], 0) + Theta * RDRadiation(rd, &x[i], 0);
333     if (rd->endpoint) {
334       PetscScalar Em0, Em1;
335       RDMaterialEnergy(rd, &x0[i], &Em0, NULL);
336       RDMaterialEnergy(rd, &x[i], &Em1, NULL);
337       Em_t = (Em1 - Em0) / dt;
338     } else {
339       RDNode dEm;
340       RDMaterialEnergy(rd, &x[i], NULL, &dEm);
341       Em_t = dEm.E * xdot[i].E + dEm.T * xdot[i].T;
342     }
343     /* Residuals are multiplied by the volume element (hx).  */
344     /* The temperature equation does not have boundary conditions */
345     f[i].T = hx * (rho * Em_t + rad);
346 
347     if (i == 0) { /* Left boundary condition */
348       PetscScalar D_R, bcTheta = rd->bcmidpoint ? Theta : 1.;
349       RDNode      n, nx;
350 
351       n.E  = (1. - bcTheta) * x0[0].E + bcTheta * x[0].E;
352       n.T  = (1. - bcTheta) * x0[0].T + bcTheta * x[0].T;
353       nx.E = ((1. - bcTheta) * (x0[1].E - x0[0].E) + bcTheta * (x[1].E - x[0].E)) / hx;
354       nx.T = ((1. - bcTheta) * (x0[1].T - x0[0].T) + bcTheta * (x[1].T - x[0].T)) / hx;
355       switch (rd->leftbc) {
356       case BC_ROBIN:
357         RDDiffusionCoefficient(rd, rd->bclimit, &n, &nx, &D_R, 0, 0);
358         f[0].E = hx * (n.E - 2. * D_R * nx.E - rd->Eapplied);
359         break;
360       case BC_NEUMANN:
361         f[0].E = x[1].E - x[0].E;
362         break;
363       default:
364         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Case %" PetscInt_FMT, rd->initial);
365       }
366     } else if (i == info.mx - 1) {  /* Right boundary */
367       f[i].E = x[i].E - x[i - 1].E; /* Homogeneous Neumann */
368     } else {
369       PetscScalar diff = (1. - Theta) * RDDiffusion(rd, hx, x0, i, 0) + Theta * RDDiffusion(rd, hx, x, i, 0);
370       f[i].E           = hx * (xdot[i].E - diff - rad);
371     }
372   }
373   PetscCall(RDRestoreLocalArrays(rd, &X0loc, &x0, &Xloc, &x, &Xloc_t, &xdot));
374   PetscCall(DMDAVecRestoreArray(rd->da, F, &f));
375   if (rd->monitor_residual) PetscCall(RDStateView(rd, X, Xdot, F));
376   PetscFunctionReturn(PETSC_SUCCESS);
377 }
378 
379 static PetscErrorCode RDIJacobian_FD(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal a, Mat A, Mat B, void *ctx)
380 {
381   RD            rd = (RD)ctx;
382   RDNode       *x, *x0, *xdot;
383   Vec           X0loc, Xloc, Xloc_t;
384   PetscReal     hx, Theta, dt;
385   DMDALocalInfo info;
386   PetscInt      i;
387 
388   PetscFunctionBeginUser;
389   PetscCall(RDGetLocalArrays(rd, ts, X, Xdot, &Theta, &dt, &X0loc, &x0, &Xloc, &x, &Xloc_t, &xdot));
390   PetscCall(DMDAGetLocalInfo(rd->da, &info));
391   hx = rd->L / (info.mx - 1);
392   PetscCall(MatZeroEntries(B));
393 
394   for (i = info.xs; i < info.xs + info.xm; i++) {
395     PetscInt                  col[3];
396     PetscReal                 rho = rd->rho;
397     PetscScalar /*Em_t,rad,*/ K[2][6];
398     RDNode                    dEm_t, drad;
399 
400     /*rad = (1.-Theta)* */ RDRadiation(rd, &x0[i], 0); /* + Theta* */
401     RDRadiation(rd, &x[i], &drad);
402 
403     if (rd->endpoint) {
404       PetscScalar Em0, Em1;
405       RDNode      dEm1;
406       RDMaterialEnergy(rd, &x0[i], &Em0, NULL);
407       RDMaterialEnergy(rd, &x[i], &Em1, &dEm1);
408       /*Em_t = (Em1 - Em0) / (Theta*dt);*/
409       dEm_t.E = dEm1.E / (Theta * dt);
410       dEm_t.T = dEm1.T / (Theta * dt);
411     } else {
412       const PetscScalar epsilon = x[i].T * PETSC_SQRT_MACHINE_EPSILON;
413       RDNode            n1;
414       RDNode            dEm, dEm1;
415       PetscScalar       Em_TT;
416 
417       n1.E = x[i].E;
418       n1.T = x[i].T + epsilon;
419       RDMaterialEnergy(rd, &x[i], NULL, &dEm);
420       RDMaterialEnergy(rd, &n1, NULL, &dEm1);
421       /* The Jacobian needs another derivative.  We finite difference here instead of
422        * propagating second derivatives through the ionization model. */
423       Em_TT = (dEm1.T - dEm.T) / epsilon;
424       /*Em_t = dEm.E * xdot[i].E + dEm.T * xdot[i].T;*/
425       dEm_t.E = dEm.E * a;
426       dEm_t.T = dEm.T * a + Em_TT * xdot[i].T;
427     }
428 
429     PetscCall(PetscMemzero(K, sizeof(K)));
430     /* Residuals are multiplied by the volume element (hx).  */
431     if (i == 0) {
432       PetscScalar D, bcTheta = rd->bcmidpoint ? Theta : 1.;
433       RDNode      n, nx;
434       RDNode      dD, dxD;
435 
436       n.E  = (1. - bcTheta) * x0[0].E + bcTheta * x[0].E;
437       n.T  = (1. - bcTheta) * x0[0].T + bcTheta * x[0].T;
438       nx.E = ((1. - bcTheta) * (x0[1].E - x0[0].E) + bcTheta * (x[1].E - x[0].E)) / hx;
439       nx.T = ((1. - bcTheta) * (x0[1].T - x0[0].T) + bcTheta * (x[1].T - x[0].T)) / hx;
440       switch (rd->leftbc) {
441       case BC_ROBIN:
442         RDDiffusionCoefficient(rd, rd->bclimit, &n, &nx, &D, &dD, &dxD);
443         K[0][1 * 2 + 0] = (bcTheta / Theta) * hx * (1. - 2. * D * (-1. / hx) - 2. * nx.E * dD.E + 2. * nx.E * dxD.E / hx);
444         K[0][1 * 2 + 1] = (bcTheta / Theta) * hx * (-2. * nx.E * dD.T);
445         K[0][2 * 2 + 0] = (bcTheta / Theta) * hx * (-2. * D * (1. / hx) - 2. * nx.E * dD.E - 2. * nx.E * dxD.E / hx);
446         break;
447       case BC_NEUMANN:
448         K[0][1 * 2 + 0] = -1. / Theta;
449         K[0][2 * 2 + 0] = 1. / Theta;
450         break;
451       default:
452         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Case %" PetscInt_FMT, rd->initial);
453       }
454     } else if (i == info.mx - 1) {
455       K[0][0 * 2 + 0] = -1. / Theta;
456       K[0][1 * 2 + 0] = 1. / Theta;
457     } else {
458       /*PetscScalar diff;*/
459       RDNode ddiff[3];
460       /*diff = (1.-Theta)*RDDiffusion(rd,hx,x0,i,0) + Theta* */ RDDiffusion(rd, hx, x, i, ddiff);
461       K[0][0 * 2 + 0] = -hx * ddiff[0].E;
462       K[0][0 * 2 + 1] = -hx * ddiff[0].T;
463       K[0][1 * 2 + 0] = hx * (a - ddiff[1].E - drad.E);
464       K[0][1 * 2 + 1] = hx * (-ddiff[1].T - drad.T);
465       K[0][2 * 2 + 0] = -hx * ddiff[2].E;
466       K[0][2 * 2 + 1] = -hx * ddiff[2].T;
467     }
468 
469     K[1][1 * 2 + 0] = hx * (rho * dEm_t.E + drad.E);
470     K[1][1 * 2 + 1] = hx * (rho * dEm_t.T + drad.T);
471 
472     col[0] = i - 1;
473     col[1] = i;
474     col[2] = i + 1 < info.mx ? i + 1 : -1;
475     PetscCall(MatSetValuesBlocked(B, 1, &i, 3, col, &K[0][0], INSERT_VALUES));
476   }
477   PetscCall(RDRestoreLocalArrays(rd, &X0loc, &x0, &Xloc, &x, &Xloc_t, &xdot));
478   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
479   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
480   if (A != B) {
481     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
482     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
483   }
484   PetscFunctionReturn(PETSC_SUCCESS);
485 }
486 
487 /* Evaluate interpolants and derivatives at a select quadrature point */
488 static void RDEvaluate(PetscReal interp[][2], PetscReal deriv[][2], PetscInt q, const RDNode x[], PetscInt i, RDNode *n, RDNode *nx)
489 {
490   PetscInt j;
491   n->E  = 0;
492   n->T  = 0;
493   nx->E = 0;
494   nx->T = 0;
495   for (j = 0; j < 2; j++) {
496     n->E += interp[q][j] * x[i + j].E;
497     n->T += interp[q][j] * x[i + j].T;
498     nx->E += deriv[q][j] * x[i + j].E;
499     nx->T += deriv[q][j] * x[i + j].T;
500   }
501 }
502 
503 /*
504  Various quadrature rules.  The nonlinear terms are non-polynomial so no standard quadrature will be exact.
505 */
506 static PetscErrorCode RDGetQuadrature(RD rd, PetscReal hx, PetscInt *nq, PetscReal weight[], PetscReal interp[][2], PetscReal deriv[][2])
507 {
508   PetscInt         q, j;
509   const PetscReal *refweight, (*refinterp)[2], (*refderiv)[2];
510 
511   PetscFunctionBeginUser;
512   switch (rd->quadrature) {
513   case QUADRATURE_GAUSS1: {
514     static const PetscReal ww[1]    = {1.};
515     static const PetscReal ii[1][2] = {
516       {0.5, 0.5}
517     };
518     static const PetscReal dd[1][2] = {
519       {-1., 1.}
520     };
521     *nq       = 1;
522     refweight = ww;
523     refinterp = ii;
524     refderiv  = dd;
525   } break;
526   case QUADRATURE_GAUSS2: {
527     static const PetscReal ii[2][2] = {
528       {0.78867513459481287, 0.21132486540518713},
529       {0.21132486540518713, 0.78867513459481287}
530     };
531     static const PetscReal dd[2][2] = {
532       {-1., 1.},
533       {-1., 1.}
534     };
535     static const PetscReal ww[2] = {0.5, 0.5};
536     *nq                          = 2;
537     refweight                    = ww;
538     refinterp                    = ii;
539     refderiv                     = dd;
540   } break;
541   case QUADRATURE_GAUSS3: {
542     static const PetscReal ii[3][2] = {
543       {0.8872983346207417, 0.1127016653792583},
544       {0.5,                0.5               },
545       {0.1127016653792583, 0.8872983346207417}
546     };
547     static const PetscReal dd[3][2] = {
548       {-1, 1},
549       {-1, 1},
550       {-1, 1}
551     };
552     static const PetscReal ww[3] = {5. / 18, 8. / 18, 5. / 18};
553     *nq                          = 3;
554     refweight                    = ww;
555     refinterp                    = ii;
556     refderiv                     = dd;
557   } break;
558   case QUADRATURE_GAUSS4: {
559     static const PetscReal ii[][2] = {
560       {0.93056815579702623,  0.069431844202973658},
561       {0.66999052179242813,  0.33000947820757187 },
562       {0.33000947820757187,  0.66999052179242813 },
563       {0.069431844202973658, 0.93056815579702623 }
564     };
565     static const PetscReal dd[][2] = {
566       {-1, 1},
567       {-1, 1},
568       {-1, 1},
569       {-1, 1}
570     };
571     static const PetscReal ww[] = {0.17392742256872692, 0.3260725774312731, 0.3260725774312731, 0.17392742256872692};
572 
573     *nq       = 4;
574     refweight = ww;
575     refinterp = ii;
576     refderiv  = dd;
577   } break;
578   case QUADRATURE_LOBATTO2: {
579     static const PetscReal ii[2][2] = {
580       {1., 0.},
581       {0., 1.}
582     };
583     static const PetscReal dd[2][2] = {
584       {-1., 1.},
585       {-1., 1.}
586     };
587     static const PetscReal ww[2] = {0.5, 0.5};
588     *nq                          = 2;
589     refweight                    = ww;
590     refinterp                    = ii;
591     refderiv                     = dd;
592   } break;
593   case QUADRATURE_LOBATTO3: {
594     static const PetscReal ii[3][2] = {
595       {1,   0  },
596       {0.5, 0.5},
597       {0,   1  }
598     };
599     static const PetscReal dd[3][2] = {
600       {-1, 1},
601       {-1, 1},
602       {-1, 1}
603     };
604     static const PetscReal ww[3] = {1. / 6, 4. / 6, 1. / 6};
605     *nq                          = 3;
606     refweight                    = ww;
607     refinterp                    = ii;
608     refderiv                     = dd;
609   } break;
610   default:
611     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unknown quadrature %d", (int)rd->quadrature);
612   }
613 
614   for (q = 0; q < *nq; q++) {
615     weight[q] = refweight[q] * hx;
616     for (j = 0; j < 2; j++) {
617       interp[q][j] = refinterp[q][j];
618       deriv[q][j]  = refderiv[q][j] / hx;
619     }
620   }
621   PetscFunctionReturn(PETSC_SUCCESS);
622 }
623 
624 /*
625  Finite element version
626 */
627 static PetscErrorCode RDIFunction_FE(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, void *ctx)
628 {
629   RD            rd = (RD)ctx;
630   RDNode       *x, *x0, *xdot, *f;
631   Vec           X0loc, Xloc, Xloc_t, Floc;
632   PetscReal     hx, Theta, dt, weight[5], interp[5][2], deriv[5][2];
633   DMDALocalInfo info;
634   PetscInt      i, j, q, nq;
635 
636   PetscFunctionBeginUser;
637   PetscCall(RDGetLocalArrays(rd, ts, X, Xdot, &Theta, &dt, &X0loc, &x0, &Xloc, &x, &Xloc_t, &xdot));
638 
639   PetscCall(DMGetLocalVector(rd->da, &Floc));
640   PetscCall(VecZeroEntries(Floc));
641   PetscCall(DMDAVecGetArray(rd->da, Floc, &f));
642   PetscCall(DMDAGetLocalInfo(rd->da, &info));
643 
644   /* Set up shape functions and quadrature for elements (assumes a uniform grid) */
645   hx = rd->L / (info.mx - 1);
646   PetscCall(RDGetQuadrature(rd, hx, &nq, weight, interp, deriv));
647 
648   for (i = info.xs; i < PetscMin(info.xs + info.xm, info.mx - 1); i++) {
649     for (q = 0; q < nq; q++) {
650       PetscReal   rho = rd->rho;
651       PetscScalar Em_t, rad, D_R, D0_R;
652       RDNode      n, n0, nx, n0x, nt, ntx;
653       RDEvaluate(interp, deriv, q, x, i, &n, &nx);
654       RDEvaluate(interp, deriv, q, x0, i, &n0, &n0x);
655       RDEvaluate(interp, deriv, q, xdot, i, &nt, &ntx);
656 
657       rad = (1. - Theta) * RDRadiation(rd, &n0, 0) + Theta * RDRadiation(rd, &n, 0);
658       if (rd->endpoint) {
659         PetscScalar Em0, Em1;
660         RDMaterialEnergy(rd, &n0, &Em0, NULL);
661         RDMaterialEnergy(rd, &n, &Em1, NULL);
662         Em_t = (Em1 - Em0) / dt;
663       } else {
664         RDNode dEm;
665         RDMaterialEnergy(rd, &n, NULL, &dEm);
666         Em_t = dEm.E * nt.E + dEm.T * nt.T;
667       }
668       RDDiffusionCoefficient(rd, PETSC_TRUE, &n0, &n0x, &D0_R, 0, 0);
669       RDDiffusionCoefficient(rd, PETSC_TRUE, &n, &nx, &D_R, 0, 0);
670       for (j = 0; j < 2; j++) {
671         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));
672         f[i + j].T += interp[q][j] * weight[q] * (rho * Em_t + rad);
673       }
674     }
675   }
676   if (info.xs == 0) {
677     switch (rd->leftbc) {
678     case BC_ROBIN: {
679       PetscScalar D_R, D_R_bc;
680       PetscReal   ratio, bcTheta = rd->bcmidpoint ? Theta : 1.;
681       RDNode      n, nx;
682 
683       n.E  = (1 - bcTheta) * x0[0].E + bcTheta * x[0].E;
684       n.T  = (1 - bcTheta) * x0[0].T + bcTheta * x[0].T;
685       nx.E = (x[1].E - x[0].E) / hx;
686       nx.T = (x[1].T - x[0].T) / hx;
687       RDDiffusionCoefficient(rd, PETSC_TRUE, &n, &nx, &D_R, 0, 0);
688       RDDiffusionCoefficient(rd, rd->bclimit, &n, &nx, &D_R_bc, 0, 0);
689       ratio = PetscRealPart(D_R / D_R_bc);
690       PetscCheck(ratio <= 1., PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Limited diffusivity is greater than unlimited");
691       PetscCheck(ratio >= 1e-3, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Heavily limited diffusivity");
692       f[0].E += -ratio * 0.5 * (rd->Eapplied - n.E);
693     } break;
694     case BC_NEUMANN:
695       /* homogeneous Neumann is the natural condition */
696       break;
697     default:
698       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Case %" PetscInt_FMT, rd->initial);
699     }
700   }
701 
702   PetscCall(RDRestoreLocalArrays(rd, &X0loc, &x0, &Xloc, &x, &Xloc_t, &xdot));
703   PetscCall(DMDAVecRestoreArray(rd->da, Floc, &f));
704   PetscCall(VecZeroEntries(F));
705   PetscCall(DMLocalToGlobalBegin(rd->da, Floc, ADD_VALUES, F));
706   PetscCall(DMLocalToGlobalEnd(rd->da, Floc, ADD_VALUES, F));
707   PetscCall(DMRestoreLocalVector(rd->da, &Floc));
708 
709   if (rd->monitor_residual) PetscCall(RDStateView(rd, X, Xdot, F));
710   PetscFunctionReturn(PETSC_SUCCESS);
711 }
712 
713 static PetscErrorCode RDIJacobian_FE(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal a, Mat A, Mat B, void *ctx)
714 {
715   RD            rd = (RD)ctx;
716   RDNode       *x, *x0, *xdot;
717   Vec           X0loc, Xloc, Xloc_t;
718   PetscReal     hx, Theta, dt, weight[5], interp[5][2], deriv[5][2];
719   DMDALocalInfo info;
720   PetscInt      i, j, k, q, nq;
721   PetscScalar   K[4][4];
722 
723   PetscFunctionBeginUser;
724   PetscCall(RDGetLocalArrays(rd, ts, X, Xdot, &Theta, &dt, &X0loc, &x0, &Xloc, &x, &Xloc_t, &xdot));
725   PetscCall(DMDAGetLocalInfo(rd->da, &info));
726   hx = rd->L / (info.mx - 1);
727   PetscCall(RDGetQuadrature(rd, hx, &nq, weight, interp, deriv));
728   PetscCall(MatZeroEntries(B));
729   for (i = info.xs; i < PetscMin(info.xs + info.xm, info.mx - 1); i++) {
730     PetscInt rc[2];
731 
732     rc[0] = i;
733     rc[1] = i + 1;
734     PetscCall(PetscMemzero(K, sizeof(K)));
735     for (q = 0; q < nq; q++) {
736       PetscScalar              D_R;
737       PETSC_UNUSED PetscScalar rad;
738       RDNode                   n, nx, nt, ntx, drad, dD_R, dxD_R, dEm;
739       RDEvaluate(interp, deriv, q, x, i, &n, &nx);
740       RDEvaluate(interp, deriv, q, xdot, i, &nt, &ntx);
741       rad = RDRadiation(rd, &n, &drad);
742       RDDiffusionCoefficient(rd, PETSC_TRUE, &n, &nx, &D_R, &dD_R, &dxD_R);
743       RDMaterialEnergy(rd, &n, NULL, &dEm);
744       for (j = 0; j < 2; j++) {
745         for (k = 0; k < 2; k++) {
746           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]));
747           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);
748           K[j * 2 + 1][k * 2 + 0] += interp[q][j] * weight[q] * drad.E * interp[q][k];
749           K[j * 2 + 1][k * 2 + 1] += interp[q][j] * weight[q] * (a * rd->rho * dEm.T + drad.T) * interp[q][k];
750         }
751       }
752     }
753     PetscCall(MatSetValuesBlocked(B, 2, rc, 2, rc, &K[0][0], ADD_VALUES));
754   }
755   if (info.xs == 0) {
756     switch (rd->leftbc) {
757     case BC_ROBIN: {
758       PetscScalar D_R, D_R_bc;
759       PetscReal   ratio;
760       RDNode      n, nx;
761 
762       n.E  = (1 - Theta) * x0[0].E + Theta * x[0].E;
763       n.T  = (1 - Theta) * x0[0].T + Theta * x[0].T;
764       nx.E = (x[1].E - x[0].E) / hx;
765       nx.T = (x[1].T - x[0].T) / hx;
766       RDDiffusionCoefficient(rd, PETSC_TRUE, &n, &nx, &D_R, 0, 0);
767       RDDiffusionCoefficient(rd, rd->bclimit, &n, &nx, &D_R_bc, 0, 0);
768       ratio = PetscRealPart(D_R / D_R_bc);
769       PetscCall(MatSetValue(B, 0, 0, ratio * 0.5, ADD_VALUES));
770     } break;
771     case BC_NEUMANN:
772       /* homogeneous Neumann is the natural condition */
773       break;
774     default:
775       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Case %" PetscInt_FMT, rd->initial);
776     }
777   }
778 
779   PetscCall(RDRestoreLocalArrays(rd, &X0loc, &x0, &Xloc, &x, &Xloc_t, &xdot));
780   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
781   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
782   if (A != B) {
783     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
784     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
785   }
786   PetscFunctionReturn(PETSC_SUCCESS);
787 }
788 
789 /* Temperature that is in equilibrium with the radiation density */
790 static PetscScalar RDRadiationTemperature(RD rd, PetscScalar E)
791 {
792   return PetscPowScalar(E * rd->c / (4. * rd->sigma_b), 0.25);
793 }
794 
795 static PetscErrorCode RDInitialState(RD rd, Vec X)
796 {
797   DMDALocalInfo info;
798   PetscInt      i;
799   RDNode       *x;
800 
801   PetscFunctionBeginUser;
802   PetscCall(DMDAGetLocalInfo(rd->da, &info));
803   PetscCall(DMDAVecGetArray(rd->da, X, &x));
804   for (i = info.xs; i < info.xs + info.xm; i++) {
805     PetscReal coord = i * rd->L / (info.mx - 1);
806     switch (rd->initial) {
807     case 1:
808       x[i].E = 0.001;
809       x[i].T = RDRadiationTemperature(rd, x[i].E);
810       break;
811     case 2:
812       x[i].E = 0.001 + 100. * PetscExpReal(-PetscSqr(coord / 0.1));
813       x[i].T = RDRadiationTemperature(rd, x[i].E);
814       break;
815     case 3:
816       x[i].E = 7.56e-2 * rd->unit.Joule / PetscPowScalarInt(rd->unit.meter, 3);
817       x[i].T = RDRadiationTemperature(rd, x[i].E);
818       break;
819     default:
820       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "No initial state %" PetscInt_FMT, rd->initial);
821     }
822   }
823   PetscCall(DMDAVecRestoreArray(rd->da, X, &x));
824   PetscFunctionReturn(PETSC_SUCCESS);
825 }
826 
827 static PetscErrorCode RDView(RD rd, Vec X, PetscViewer viewer)
828 {
829   Vec             Y;
830   const RDNode   *x;
831   PetscScalar    *y;
832   PetscInt        i, m, M;
833   const PetscInt *lx;
834   DM              da;
835   MPI_Comm        comm;
836 
837   PetscFunctionBeginUser;
838   /*
839     Create a DMDA (one dof per node, zero stencil width, same layout) to hold Trad
840     (radiation temperature).  It is not necessary to create a DMDA for this, but this way
841     output and visualization will have meaningful variable names and correct scales.
842   */
843   PetscCall(DMDAGetInfo(rd->da, 0, &M, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0));
844   PetscCall(DMDAGetOwnershipRanges(rd->da, &lx, 0, 0));
845   PetscCall(PetscObjectGetComm((PetscObject)rd->da, &comm));
846   PetscCall(DMDACreate1d(comm, DM_BOUNDARY_NONE, M, 1, 0, lx, &da));
847   PetscCall(DMSetFromOptions(da));
848   PetscCall(DMSetUp(da));
849   PetscCall(DMDASetUniformCoordinates(da, 0., rd->L, 0., 0., 0., 0.));
850   PetscCall(DMDASetFieldName(da, 0, "T_rad"));
851   PetscCall(DMCreateGlobalVector(da, &Y));
852 
853   /* Compute the radiation temperature from the solution at each node */
854   PetscCall(VecGetLocalSize(Y, &m));
855   PetscCall(VecGetArrayRead(X, (const PetscScalar **)&x));
856   PetscCall(VecGetArray(Y, &y));
857   for (i = 0; i < m; i++) y[i] = RDRadiationTemperature(rd, x[i].E);
858   PetscCall(VecRestoreArrayRead(X, (const PetscScalar **)&x));
859   PetscCall(VecRestoreArray(Y, &y));
860 
861   PetscCall(VecView(Y, viewer));
862   PetscCall(VecDestroy(&Y));
863   PetscCall(DMDestroy(&da));
864   PetscFunctionReturn(PETSC_SUCCESS);
865 }
866 
867 static PetscErrorCode RDTestDifferentiation(RD rd)
868 {
869   MPI_Comm    comm;
870   RDNode      n, nx;
871   PetscScalar epsilon;
872 
873   PetscFunctionBeginUser;
874   PetscCall(PetscObjectGetComm((PetscObject)rd->da, &comm));
875   epsilon = 1e-8;
876   {
877     RDNode      dEm, fdEm;
878     PetscScalar T0 = 1000., T1 = T0 * (1. + epsilon), Em0, Em1;
879     n.E = 1.;
880     n.T = T0;
881     rd->MaterialEnergy(rd, &n, &Em0, &dEm);
882     n.E = 1. + epsilon;
883     n.T = T0;
884     rd->MaterialEnergy(rd, &n, &Em1, 0);
885     fdEm.E = (Em1 - Em0) / epsilon;
886     n.E    = 1.;
887     n.T    = T1;
888     rd->MaterialEnergy(rd, &n, &Em1, 0);
889     fdEm.T = (Em1 - Em0) / (T0 * epsilon);
890     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),
891                           (double)PetscRealPart(dEm.T - fdEm.T)));
892   }
893   {
894     PetscScalar D0, D;
895     RDNode      dD, dxD, fdD, fdxD;
896     n.E  = 1.;
897     n.T  = 1.;
898     nx.E = 1.;
899     n.T  = 1.;
900     RDDiffusionCoefficient(rd, rd->bclimit, &n, &nx, &D0, &dD, &dxD);
901     n.E  = 1. + epsilon;
902     n.T  = 1.;
903     nx.E = 1.;
904     n.T  = 1.;
905     RDDiffusionCoefficient(rd, rd->bclimit, &n, &nx, &D, 0, 0);
906     fdD.E = (D - D0) / epsilon;
907     n.E   = 1;
908     n.T   = 1. + epsilon;
909     nx.E  = 1.;
910     n.T   = 1.;
911     RDDiffusionCoefficient(rd, rd->bclimit, &n, &nx, &D, 0, 0);
912     fdD.T = (D - D0) / epsilon;
913     n.E   = 1;
914     n.T   = 1.;
915     nx.E  = 1. + epsilon;
916     n.T   = 1.;
917     RDDiffusionCoefficient(rd, rd->bclimit, &n, &nx, &D, 0, 0);
918     fdxD.E = (D - D0) / epsilon;
919     n.E    = 1;
920     n.T    = 1.;
921     nx.E   = 1.;
922     n.T    = 1. + epsilon;
923     RDDiffusionCoefficient(rd, rd->bclimit, &n, &nx, &D, 0, 0);
924     fdxD.T = (D - D0) / epsilon;
925     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),
926                           (double)PetscRealPart(dD.T - fdD.T)));
927     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),
928                           (double)PetscRealPart(dxD.T - fdxD.T)));
929   }
930   {
931     PetscInt    i;
932     PetscReal   hx = 1.;
933     PetscScalar a0;
934     RDNode      n0[3], n1[3], d[3], fd[3];
935 
936     n0[0].E = 1.;
937     n0[0].T = 1.;
938     n0[1].E = 5.;
939     n0[1].T = 3.;
940     n0[2].E = 4.;
941     n0[2].T = 2.;
942     a0      = RDDiffusion(rd, hx, n0, 1, d);
943     for (i = 0; i < 3; i++) {
944       PetscCall(PetscMemcpy(n1, n0, sizeof(n0)));
945       n1[i].E += epsilon;
946       fd[i].E = (RDDiffusion(rd, hx, n1, 1, 0) - a0) / epsilon;
947       PetscCall(PetscMemcpy(n1, n0, sizeof(n0)));
948       n1[i].T += epsilon;
949       fd[i].T = (RDDiffusion(rd, hx, n1, 1, 0) - a0) / epsilon;
950       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),
951                             (double)PetscRealPart(d[i].E - fd[i].E), (double)PetscRealPart(d[i].T - fd[i].T)));
952     }
953   }
954   {
955     PetscScalar rad0, rad;
956     RDNode      drad, fdrad;
957     n.E     = 1.;
958     n.T     = 1.;
959     rad0    = RDRadiation(rd, &n, &drad);
960     n.E     = 1. + epsilon;
961     n.T     = 1.;
962     rad     = RDRadiation(rd, &n, 0);
963     fdrad.E = (rad - rad0) / epsilon;
964     n.E     = 1.;
965     n.T     = 1. + epsilon;
966     rad     = RDRadiation(rd, &n, 0);
967     fdrad.T = (rad - rad0) / epsilon;
968     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),
969                           (double)PetscRealPart(drad.T - fdrad.T)));
970   }
971   PetscFunctionReturn(PETSC_SUCCESS);
972 }
973 
974 static PetscErrorCode RDCreate(MPI_Comm comm, RD *inrd)
975 {
976   RD        rd;
977   PetscReal meter = 0, kilogram = 0, second = 0, Kelvin = 0, Joule = 0, Watt = 0;
978 
979   PetscFunctionBeginUser;
980   *inrd = 0;
981   PetscCall(PetscNew(&rd));
982 
983   PetscOptionsBegin(comm, NULL, "Options for nonequilibrium radiation-diffusion with RD ionization", NULL);
984   {
985     rd->initial = 1;
986     PetscCall(PetscOptionsInt("-rd_initial", "Initial condition (1=Marshak, 2=Blast, 3=Marshak+)", "", rd->initial, &rd->initial, 0));
987     switch (rd->initial) {
988     case 1:
989     case 2:
990       rd->unit.kilogram = 1.;
991       rd->unit.meter    = 1.;
992       rd->unit.second   = 1.;
993       rd->unit.Kelvin   = 1.;
994       break;
995     case 3:
996       rd->unit.kilogram = 1.e12;
997       rd->unit.meter    = 1.;
998       rd->unit.second   = 1.e9;
999       rd->unit.Kelvin   = 1.;
1000       break;
1001     default:
1002       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unknown initial condition %" PetscInt_FMT, rd->initial);
1003     }
1004     /* Fundamental units */
1005     PetscCall(PetscOptionsReal("-rd_unit_meter", "Length of 1 meter in nondimensional units", "", rd->unit.meter, &rd->unit.meter, 0));
1006     PetscCall(PetscOptionsReal("-rd_unit_kilogram", "Mass of 1 kilogram in nondimensional units", "", rd->unit.kilogram, &rd->unit.kilogram, 0));
1007     PetscCall(PetscOptionsReal("-rd_unit_second", "Time of a second in nondimensional units", "", rd->unit.second, &rd->unit.second, 0));
1008     PetscCall(PetscOptionsReal("-rd_unit_Kelvin", "Temperature of a Kelvin in nondimensional units", "", rd->unit.Kelvin, &rd->unit.Kelvin, 0));
1009     /* Derived units */
1010     rd->unit.Joule = rd->unit.kilogram * PetscSqr(rd->unit.meter / rd->unit.second);
1011     rd->unit.Watt  = rd->unit.Joule / rd->unit.second;
1012     /* Local aliases */
1013     meter    = rd->unit.meter;
1014     kilogram = rd->unit.kilogram;
1015     second   = rd->unit.second;
1016     Kelvin   = rd->unit.Kelvin;
1017     Joule    = rd->unit.Joule;
1018     Watt     = rd->unit.Watt;
1019 
1020     PetscCall(PetscOptionsBool("-rd_monitor_residual", "Display residuals every time they are evaluated", "", rd->monitor_residual, &rd->monitor_residual, NULL));
1021     PetscCall(PetscOptionsEnum("-rd_discretization", "Discretization type", "", DiscretizationTypes, (PetscEnum)rd->discretization, (PetscEnum *)&rd->discretization, NULL));
1022     if (rd->discretization == DISCRETIZATION_FE) {
1023       rd->quadrature = QUADRATURE_GAUSS2;
1024       PetscCall(PetscOptionsEnum("-rd_quadrature", "Finite element quadrature", "", QuadratureTypes, (PetscEnum)rd->quadrature, (PetscEnum *)&rd->quadrature, NULL));
1025     }
1026     PetscCall(PetscOptionsEnum("-rd_jacobian", "Type of finite difference Jacobian", "", JacobianTypes, (PetscEnum)rd->jacobian, (PetscEnum *)&rd->jacobian, NULL));
1027     switch (rd->initial) {
1028     case 1:
1029       rd->leftbc     = BC_ROBIN;
1030       rd->Eapplied   = 4 * rd->unit.Joule / PetscPowRealInt(rd->unit.meter, 3);
1031       rd->L          = 1. * rd->unit.meter;
1032       rd->beta       = 3.0;
1033       rd->gamma      = 3.0;
1034       rd->final_time = 3 * second;
1035       break;
1036     case 2:
1037       rd->leftbc     = BC_NEUMANN;
1038       rd->Eapplied   = 0.;
1039       rd->L          = 1. * rd->unit.meter;
1040       rd->beta       = 3.0;
1041       rd->gamma      = 3.0;
1042       rd->final_time = 1 * second;
1043       break;
1044     case 3:
1045       rd->leftbc     = BC_ROBIN;
1046       rd->Eapplied   = 7.503e6 * rd->unit.Joule / PetscPowRealInt(rd->unit.meter, 3);
1047       rd->L          = 5. * rd->unit.meter;
1048       rd->beta       = 3.5;
1049       rd->gamma      = 3.5;
1050       rd->final_time = 20e-9 * second;
1051       break;
1052     default:
1053       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Initial %" PetscInt_FMT, rd->initial);
1054     }
1055     PetscCall(PetscOptionsEnum("-rd_leftbc", "Left boundary condition", "", BCTypes, (PetscEnum)rd->leftbc, (PetscEnum *)&rd->leftbc, NULL));
1056     PetscCall(PetscOptionsReal("-rd_E_applied", "Radiation flux at left end of domain", "", rd->Eapplied, &rd->Eapplied, NULL));
1057     PetscCall(PetscOptionsReal("-rd_beta", "Thermal exponent for photon absorption", "", rd->beta, &rd->beta, NULL));
1058     PetscCall(PetscOptionsReal("-rd_gamma", "Thermal exponent for diffusion coefficient", "", rd->gamma, &rd->gamma, NULL));
1059     PetscCall(PetscOptionsBool("-rd_view_draw", "Draw final solution", "", rd->view_draw, &rd->view_draw, NULL));
1060     PetscCall(PetscOptionsBool("-rd_endpoint", "Discretize using endpoints (like trapezoid rule) instead of midpoint", "", rd->endpoint, &rd->endpoint, NULL));
1061     PetscCall(PetscOptionsBool("-rd_bcmidpoint", "Impose the boundary condition at the midpoint (Theta) of the interval", "", rd->bcmidpoint, &rd->bcmidpoint, NULL));
1062     PetscCall(PetscOptionsBool("-rd_bclimit", "Limit diffusion coefficient in definition of Robin boundary condition", "", rd->bclimit, &rd->bclimit, NULL));
1063     PetscCall(PetscOptionsBool("-rd_test_diff", "Test differentiation in constitutive relations", "", rd->test_diff, &rd->test_diff, NULL));
1064     PetscCall(PetscOptionsString("-rd_view_binary", "File name to hold final solution", "", rd->view_binary, rd->view_binary, sizeof(rd->view_binary), NULL));
1065   }
1066   PetscOptionsEnd();
1067 
1068   switch (rd->initial) {
1069   case 1:
1070   case 2:
1071     rd->rho            = 1.;
1072     rd->c              = 1.;
1073     rd->K_R            = 1.;
1074     rd->K_p            = 1.;
1075     rd->sigma_b        = 0.25;
1076     rd->MaterialEnergy = RDMaterialEnergy_Reduced;
1077     break;
1078   case 3:
1079     /* Table 2 */
1080     rd->rho            = 1.17e-3 * kilogram / (meter * meter * meter); /* density */
1081     rd->K_R            = 7.44e18 * PetscPowRealInt(meter, 5) * PetscPowReal(Kelvin, 3.5) * PetscPowRealInt(kilogram, -2);
1082     rd->K_p            = 2.33e20 * PetscPowRealInt(meter, 5) * PetscPowReal(Kelvin, 3.5) * PetscPowRealInt(kilogram, -2);
1083     rd->I_H            = 2.179e-18 * Joule;                                                         /* Hydrogen ionization potential */
1084     rd->m_p            = 1.673e-27 * kilogram;                                                      /* proton mass */
1085     rd->m_e            = 9.109e-31 * kilogram;                                                      /* electron mass */
1086     rd->h              = 6.626e-34 * Joule * second;                                                /* Planck's constant */
1087     rd->k              = 1.381e-23 * Joule / Kelvin;                                                /* Boltzman constant */
1088     rd->c              = 3.00e8 * meter / second;                                                   /* speed of light */
1089     rd->sigma_b        = 5.67e-8 * Watt * PetscPowRealInt(meter, -2) * PetscPowRealInt(Kelvin, -4); /* Stefan-Boltzman constant */
1090     rd->MaterialEnergy = RDMaterialEnergy_Saha;
1091     break;
1092   }
1093 
1094   PetscCall(DMDACreate1d(comm, DM_BOUNDARY_NONE, 20, sizeof(RDNode) / sizeof(PetscScalar), 1, NULL, &rd->da));
1095   PetscCall(DMSetFromOptions(rd->da));
1096   PetscCall(DMSetUp(rd->da));
1097   PetscCall(DMDASetFieldName(rd->da, 0, "E"));
1098   PetscCall(DMDASetFieldName(rd->da, 1, "T"));
1099   PetscCall(DMDASetUniformCoordinates(rd->da, 0., 1., 0., 0., 0., 0.));
1100 
1101   *inrd = rd;
1102   PetscFunctionReturn(PETSC_SUCCESS);
1103 }
1104 
1105 int main(int argc, char *argv[])
1106 {
1107   RD        rd;
1108   TS        ts;
1109   SNES      snes;
1110   Vec       X;
1111   Mat       A, B;
1112   PetscInt  steps;
1113   PetscReal ftime;
1114 
1115   PetscFunctionBeginUser;
1116   PetscCall(PetscInitialize(&argc, &argv, 0, NULL));
1117   PetscCall(RDCreate(PETSC_COMM_WORLD, &rd));
1118   PetscCall(DMCreateGlobalVector(rd->da, &X));
1119   PetscCall(DMSetMatType(rd->da, MATAIJ));
1120   PetscCall(DMCreateMatrix(rd->da, &B));
1121   PetscCall(RDInitialState(rd, X));
1122 
1123   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
1124   PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
1125   PetscCall(TSSetType(ts, TSTHETA));
1126   PetscCall(TSSetDM(ts, rd->da));
1127   switch (rd->discretization) {
1128   case DISCRETIZATION_FD:
1129     PetscCall(TSSetIFunction(ts, NULL, RDIFunction_FD, rd));
1130     if (rd->jacobian == JACOBIAN_ANALYTIC) PetscCall(TSSetIJacobian(ts, B, B, RDIJacobian_FD, rd));
1131     break;
1132   case DISCRETIZATION_FE:
1133     PetscCall(TSSetIFunction(ts, NULL, RDIFunction_FE, rd));
1134     if (rd->jacobian == JACOBIAN_ANALYTIC) PetscCall(TSSetIJacobian(ts, B, B, RDIJacobian_FE, rd));
1135     break;
1136   }
1137   PetscCall(TSSetMaxTime(ts, rd->final_time));
1138   PetscCall(TSSetTimeStep(ts, 1e-3));
1139   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
1140   PetscCall(TSSetFromOptions(ts));
1141 
1142   A = B;
1143   PetscCall(TSGetSNES(ts, &snes));
1144   switch (rd->jacobian) {
1145   case JACOBIAN_ANALYTIC:
1146     break;
1147   case JACOBIAN_MATRIXFREE:
1148     break;
1149   case JACOBIAN_FD_COLORING: {
1150     PetscCall(SNESSetJacobian(snes, A, B, SNESComputeJacobianDefaultColor, 0));
1151   } break;
1152   case JACOBIAN_FD_FULL:
1153     PetscCall(SNESSetJacobian(snes, A, B, SNESComputeJacobianDefault, ts));
1154     break;
1155   }
1156 
1157   if (rd->test_diff) PetscCall(RDTestDifferentiation(rd));
1158   PetscCall(TSSolve(ts, X));
1159   PetscCall(TSGetSolveTime(ts, &ftime));
1160   PetscCall(TSGetStepNumber(ts, &steps));
1161   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Steps %" PetscInt_FMT "  final time %g\n", steps, (double)ftime));
1162   if (rd->view_draw) PetscCall(RDView(rd, X, PETSC_VIEWER_DRAW_WORLD));
1163   if (rd->view_binary[0]) {
1164     PetscViewer viewer;
1165     PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, rd->view_binary, FILE_MODE_WRITE, &viewer));
1166     PetscCall(RDView(rd, X, viewer));
1167     PetscCall(PetscViewerDestroy(&viewer));
1168   }
1169   PetscCall(VecDestroy(&X));
1170   PetscCall(MatDestroy(&B));
1171   PetscCall(RDDestroy(&rd));
1172   PetscCall(TSDestroy(&ts));
1173   PetscCall(PetscFinalize());
1174   return 0;
1175 }
1176 /*TEST
1177 
1178     test:
1179       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
1180       requires: !single
1181 
1182     test:
1183       suffix: 2
1184       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
1185       requires: !single
1186 
1187     test:
1188       suffix: 3
1189       nsize: 2
1190       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
1191       requires: !single
1192 
1193 TEST*/
1194