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