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