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