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