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