13 template <
int dim, 
int nstate, 
typename real>
    16     const dealii::Point<dim,real> &pos,
    17     const std::array<real,nstate> &conservative_soln,
    18     const real current_time,
    19     const dealii::types::global_dof_index )
 const    21     return source_term(pos,conservative_soln,current_time);
    24 template <
int dim, 
int nstate, 
typename real>
    27     const dealii::Point<dim,real> &,
    28     const std::array<real,nstate> &,
    31     std::array<real,nstate> source_term;
    32     for (
int s=0; s<nstate; s++) {
    40 template <
int dim, 
int nstate, 
typename real>
    44     std::array<real, nstate> primitive_soln;
    46     real density = conservative_soln[0];
    47     dealii::Tensor<1,dim,real> vel = compute_velocities (conservative_soln);
    48     real pressure = compute_pressure (conservative_soln);
    50     primitive_soln[0] = density;
    51     for (
int d=0; d<dim; ++d) {
    52         primitive_soln[1+d] = vel[d];
    54     primitive_soln[nstate-1] = pressure;
    55     return primitive_soln;
    59 template <
int dim, 
int nstate, 
typename real>
    64     const real density = primitive_soln[0];
    65     const dealii::Tensor<1,dim,real> velocities = extract_velocities_from_primitive(primitive_soln);
    67     std::array<real, nstate> conservative_soln;
    68     conservative_soln[0] = density;
    69     for (
int d=0; d<dim; ++d) {
    70         conservative_soln[1+d] = density*velocities[d];
    72     conservative_soln[nstate-1] = compute_total_energy(primitive_soln);
    74     return conservative_soln;
    77 template <
int dim, 
int nstate, 
typename real>
    81     const real density = conservative_soln[0];
    82     dealii::Tensor<1,dim,real> vel;
    83     for (
int d=0; d<dim; ++d) { vel[d] = conservative_soln[1+d]/density; }
    87 template <
int dim, 
int nstate, 
typename real>
    92     for (
int d=0; d<dim; d++) { vel2 = vel2 + velocities[d]*velocities[d]; }
    96 template <
int dim, 
int nstate, 
typename real>
   100     dealii::Tensor<1,dim,real> velocities;
   101     for (
int d=0; d<dim; d++) { velocities[d] = primitive_soln[1+d]; }
   107 template <
int dim, 
int nstate, 
typename real>
   111     const real density = primitive_soln[0];
   112     const real pressure = primitive_soln[nstate-1];
   113     const dealii::Tensor<1,dim,real> velocities = extract_velocities_from_primitive(primitive_soln);
   114     const real vel2 = compute_velocity_squared(velocities);
   116     const real tot_energy = pressure / gamm1 + 0.5*density*vel2;
   121 template <
int dim, 
int nstate, 
typename real>
   125     const real density = conservative_soln[0];
   126     const real pressure = compute_pressure(conservative_soln);
   127     const real entropy_measure = pressure*pow(density,-gam);
   128     return entropy_measure;
   132 template <
int dim, 
int nstate, 
typename real>
   136     const real density = conservative_soln[0];
   137     const real total_energy = conservative_soln[nstate-1];
   138     const real specific_enthalpy = (total_energy+pressure)/density;
   139     return specific_enthalpy;
   143 template <
int dim, 
int nstate, 
typename real>
   147     const real density = primitive_soln[0];
   148     const real pressure = primitive_soln[nstate-1];
   149     const real temperature = gam*pressure/density;
   153 template <
int dim, 
int nstate, 
typename real>
   157     const real dimensional_temperature = compute_dimensional_temperature(primitive_soln);
   158     const real temperature = dimensional_temperature ;
   162 template <
int dim, 
int nstate, 
typename real>
   166     const real density = gam*pressure/temperature ;
   169 template <
int dim, 
int nstate, 
typename real>
   173     const real temperature = gam*pressure/density ;
   178 template <
int dim, 
int nstate, 
typename real>
   182     const real density = conservative_soln[0];
   185     const real tot_energy  = conservative_soln[nstate-1];
   188     const dealii::Tensor<1,dim,real> vel = compute_velocities(conservative_soln);
   193     const real vel2 = compute_velocity_squared(vel);
   195     real pressure = gamm1*(tot_energy - 0.5*density*vel2);
   198         std::cout<<
"Cannot compute pressure..."<<std::endl;
   199         std::cout<<
"density "<<density<<std::endl;
   200         for(
int d=0;d<dim;d++) std::cout<<
"vel"<<d<<
" "<<vel[d]<<std::endl;
   201         std::cout<<
"energy "<<tot_energy<<std::endl;
   203     assert(pressure>0.0);
   208 template <
int dim, 
int nstate, 
typename real>
   212     real density = conservative_soln[0];
   215         std::cout<<
"density"<<density<<std::endl;
   219     const real pressure = compute_pressure(conservative_soln);
   221     const real sound = sqrt(pressure*gam/density);
   226 template <
int dim, 
int nstate, 
typename real>
   231     const real sound = sqrt(pressure*gam/density);
   235 template <
int dim, 
int nstate, 
typename real>
   239     const dealii::Tensor<1,dim,real> vel = compute_velocities(conservative_soln);
   240     const real velocity = sqrt(compute_velocity_squared(vel));
   241     const real sound = compute_sound (conservative_soln);
   242     const real mach_number = velocity/sound;
   246 template <
int dim, 
int nstate, 
typename real>
   250     real magnetic_energy = 0;
   251     for (
int i = 1; i <= 3; ++i)
   252         magnetic_energy += 1./2. * (conservative_soln[nstate - i] * conservative_soln[nstate - i] );
   253     return magnetic_energy;
   259 template <
int dim, 
int nstate, 
typename real>
   262                           const std::array<real,nstate> &soln_loop)
 const   264     return (soln_const[0] + soln_loop[0])/2.;
   267 template <
int dim, 
int nstate, 
typename real>
   270                       const std::array<real,nstate> &soln_loop)
 const   272     real pressure_const = compute_pressure(soln_const);
   273     real pressure_loop = compute_pressure(soln_loop);
   274     return (pressure_const + pressure_loop)/2.;
   277 template <
int dim, 
int nstate, 
typename real>
   280                         const std::array<real,nstate> &soln_loop)
 const   282     dealii::Tensor<1,dim,real> vel_const = compute_velocities(soln_const);
   283     dealii::Tensor<1,dim,real> vel_loop = compute_velocities(soln_loop);
   285     dealii::Tensor<1,dim,real> mean_vel;
   286     for (
int d=0; d<0; ++d) {
   287         mean_vel[d] = (vel_const[d] + vel_loop[d]) * 0.5;
   292 template <
int dim, 
int nstate, 
typename real>
   295                              const std::array<real,nstate> &soln_loop)
 const   297     return ((soln_const[nstate-1]/soln_const[0]) + (soln_loop[nstate-1]/soln_loop[0]))/2.;
   301 template <
int dim, 
int nstate, 
typename real>
   305     std::array<dealii::Tensor<1,dim,real>,nstate> conv_flux;
   306     const real density = conservative_soln[0];
   307     const real pressure = compute_pressure (conservative_soln);
   308     const dealii::Tensor<1,dim,real> vel = compute_velocities(conservative_soln);
   309     const real specific_total_energy = conservative_soln[nstate-1]/conservative_soln[0];
   310     const real specific_total_enthalpy = specific_total_energy + pressure/density;
   311     const real magnetic_energy = compute_magnetic_energy(conservative_soln);
   313     for (
int flux_dim=0; flux_dim<dim; ++flux_dim) {
   315         conv_flux[0][flux_dim] = conservative_soln[1+flux_dim];
   317         for (
int velocity_dim=0; velocity_dim<dim; ++velocity_dim){
   318             conv_flux[1+velocity_dim][flux_dim] = density*vel[flux_dim]*vel[velocity_dim];
   320         conv_flux[1+flux_dim][flux_dim] += pressure + magnetic_energy; 
   322         conv_flux[nstate-4][flux_dim] = density*vel[flux_dim]*specific_total_enthalpy;
   327 template <
int dim, 
int nstate, 
typename real>
   330     const std::array<real,nstate> &conservative_soln)
 const   332     std::cout<<
"Entropy variables for MHD hasn't been done yet."<<std::endl;
   334     return conservative_soln;
   337 template <
int dim, 
int nstate, 
typename real>
   340     const std::array<real,nstate> &entropy_var)
 const   342     std::cout<<
"Entropy variables for MHD hasn't been done yet."<<std::endl;
   347 template <
int dim, 
int nstate, 
typename real>
   351     std::array<real, nstate> conv_normal_flux;
   352     const real density = conservative_soln[0];
   353     const real pressure = compute_pressure (conservative_soln);
   354     const dealii::Tensor<1,dim,real> vel = compute_velocities(conservative_soln);
   356     real normal_vel = 0.0;
   357     for (
int d=0; d<dim; ++d) {
   358         normal_vel += vel[d]*normal[d];
   360     const real total_energy = conservative_soln[nstate-1];
   361     const real specific_total_enthalpy = (total_energy + pressure) / density;
   363     const real rhoV = density*normal_vel;
   365     conv_normal_flux[0] = rhoV;
   367     for (
int velocity_dim=0; velocity_dim<dim; ++velocity_dim){
   368         conv_normal_flux[1+velocity_dim] = rhoV*vel[velocity_dim] + normal[velocity_dim] * pressure;
   371     conv_normal_flux[nstate-1] = rhoV*specific_total_enthalpy;
   372     return conv_normal_flux;
   375 template <
int dim, 
int nstate, 
typename real>
   378     const std::array<real,nstate> &conservative_soln,
   379     const dealii::Tensor<1,dim,real> &normal)
 const   382     const dealii::Tensor<1,dim,real> vel = compute_velocities(conservative_soln);
   383     real vel_normal = 0.0;
   384     for (
int d=0;d<dim;d++) { vel_normal += vel[d] * normal[d]; }
   386     const real vel2 = compute_velocity_squared(vel);
   387     const real phi = 0.5*gamm1 * vel2;
   389     const real density = conservative_soln[0];
   390     const real tot_energy = conservative_soln[nstate-1];
   391     const real E = tot_energy / density;
   392     const real a1 = gam*E-phi;
   393     const real a2 = gam-1.0;
   394     const real a3 = gam-2.0;
   396     dealii::Tensor<2,nstate,real> jacobian;
   397     for (
int d=0; d<dim; ++d) {
   398         jacobian[0][1+d] = normal[d];
   400     for (
int row_dim=0; row_dim<dim; ++row_dim) {
   401         jacobian[1+row_dim][0] = normal[row_dim]*phi - vel[row_dim] * vel_normal;
   402         for (
int col_dim=0; col_dim<dim; ++col_dim){
   403             if (row_dim == col_dim) {
   404                 jacobian[1+row_dim][1+col_dim] = vel_normal - a3*normal[row_dim]*vel[row_dim];
   406                 jacobian[1+row_dim][1+col_dim] = normal[col_dim]*vel[row_dim] - a2*normal[row_dim]*vel[col_dim];
   409         jacobian[1+row_dim][nstate-1] = normal[row_dim]*a2;
   411     jacobian[nstate-1][0] = vel_normal*(phi-a1);
   412     for (
int d=0; d<dim; ++d){
   413         jacobian[nstate-1][1+d] = normal[d]*a1 - a2*vel[d]*vel_normal;
   415     jacobian[nstate-1][nstate-1] = gam*vel_normal;
   420 template <
int dim, 
int nstate, 
typename real>
   423     const std::array<real,nstate> &conservative_soln,
   424     const dealii::Tensor<1,dim,real> &normal)
 const   426     const dealii::Tensor<1,dim,real> vel = compute_velocities(conservative_soln);
   427     std::array<real,nstate> eig;
   428     real vel_dot_n = 0.0;
   429     for (
int d=0;d<dim;++d) { vel_dot_n += vel[d]*normal[d]; };
   430     for (
int i=0; i<nstate; i++) {
   439 template <
int dim, 
int nstate, 
typename real>
   444     const dealii::Tensor<1,dim,real> vel = compute_velocities(conservative_soln);
   447     const real sound = compute_sound (conservative_soln);
   450      real vel2 = compute_velocity_squared(vel);
   456     const real max_eig = sqrt(vel2) + sound;
   462 template <
int dim, 
int nstate, 
typename real>
   469 template <
int dim, 
int nstate, 
typename real>
   472     const std::array<real,nstate> &conservative_soln,
   473     const std::array<dealii::Tensor<1,dim,real>,nstate> &solution_gradient,
   474     const dealii::types::global_dof_index )
 const   476     return dissipative_flux(conservative_soln,solution_gradient);
   479 template <
int dim, 
int nstate, 
typename real>
   482     const std::array<real,nstate> &,
   483     const std::array<dealii::Tensor<1,dim,real>,nstate> &)
 const   485     std::array<dealii::Tensor<1,dim,real>,nstate> diss_flux;
   487     for (
int i=0; i<nstate; i++) {
   493 template <
int dim, 
int nstate, 
typename real>
   497    const dealii::Point<dim, real> &pos,
   498    const dealii::Tensor<1,dim,real> &normal_int,
   499    const std::array<real,nstate> &soln_int,
   500    const std::array<dealii::Tensor<1,dim,real>,nstate> &soln_grad_int,
   501    std::array<real,nstate> &soln_bc,
   502    std::array<dealii::Tensor<1,dim,real>,nstate> &soln_grad_bc)
 const   504     std::array<real,nstate> boundary_values;
   505     std::array<dealii::Tensor<1,dim,real>,nstate> boundary_gradients;
   506     for (
int s=0; s<nstate; s++) {
   507         boundary_values[s] = this->manufactured_solution_function->value (pos, s);
   508         boundary_gradients[s] = this->manufactured_solution_function->gradient (pos, s);
   511     for (
int istate=0; istate<nstate; ++istate) {
   513         std::array<real,nstate> characteristic_dot_n = convective_eigenvalues(boundary_values, normal_int);
   514         const bool inflow = (characteristic_dot_n[istate] <= 0.);
   520             soln_bc[istate] = boundary_values[istate];
   521             soln_grad_bc[istate] = soln_grad_int[istate];
   528             soln_bc[istate] = soln_int[istate];
   535             soln_grad_bc[istate] = soln_grad_int[istate];
 real compute_mean_density(const std::array< real, nstate > &conservative_soln1, const std::array< real, nstate > &convervative_soln2) const
Mean density given two sets of conservative solutions. 
real max_convective_eigenvalue(const std::array< real, nstate > &soln) const
Maximum convective eigenvalue. 
real compute_density_from_pressure_temperature(const real pressure, const real temperature) const
Given pressure and temperature, returns NON-DIMENSIONALIZED density using free-stream non-dimensional...
std::array< real, nstate > convert_primitive_to_conservative(const std::array< real, nstate > &primitive_soln) const
std::array< real, nstate > compute_conservative_variables_from_entropy_variables(const std::array< real, nstate > &entropy_var) const
Computes the conservative variables from the entropy variables. 
real compute_temperature_from_density_pressure(const real density, const real pressure) const
Given density and pressure, returns NON-DIMENSIONALIZED temperature using free-stream non-dimensional...
real compute_pressure(const std::array< real, nstate > &conservative_soln) const
Evaluate pressure from conservative variables. 
Files for the baseline physics. 
real compute_magnetic_energy(const std::array< real, nstate > &conservative_soln) const
Evaluate Magnetic Energy. 
std::array< real, nstate > convert_conservative_to_primitive(const std::array< real, nstate > &conservative_soln) const
dealii::Tensor< 2, nstate, real > convective_flux_directional_jacobian(const std::array< real, nstate > &conservative_soln, const dealii::Tensor< 1, dim, real > &normal) const
Convective flux Jacobian: . 
real compute_entropy_measure(const std::array< real, nstate > &conservative_soln) const
Evaluate entropy from conservative variables. 
real compute_velocity_squared(const dealii::Tensor< 1, dim, real > &velocities) const
Given the velocity vector , returns the dot-product . 
dealii::Tensor< 1, dim, real > compute_velocities(const std::array< real, nstate > &conservative_soln) const
Evaluate velocities from conservative variables. 
dealii::Tensor< 1, dim, real > extract_velocities_from_primitive(const std::array< real, nstate > &primitive_soln) const
Given primitive variables, returns velocities. 
real compute_sound(const std::array< real, nstate > &conservative_soln) const
Evaluate speed of sound from conservative variables. 
real max_viscous_eigenvalue(const std::array< real, nstate > &soln) const
Maximum viscous eigenvalue. 
void boundary_face_values(const int, const dealii::Point< dim, real > &, const dealii::Tensor< 1, dim, real > &, const std::array< real, nstate > &, const std::array< dealii::Tensor< 1, dim, real >, nstate > &, std::array< real, nstate > &, std::array< dealii::Tensor< 1, dim, real >, nstate > &) const
Evaluates boundary values and gradients on the other side of the face. 
real compute_mean_specific_energy(const std::array< real, nstate > &conservative_soln1, const std::array< real, nstate > &convervative_soln2) const
Mean specific energy given two sets of conservative solutions. 
std::array< real, nstate > convective_eigenvalues(const std::array< real, nstate > &, const dealii::Tensor< 1, dim, real > &) const
Spectral radius of convective term Jacobian is 'c'. 
real compute_specific_enthalpy(const std::array< real, nstate > &conservative_soln, const real pressure) const
Evaluate pressure from conservative variables. 
std::array< dealii::Tensor< 1, dim, real >, nstate > dissipative_flux(const std::array< real, nstate > &conservative_soln, const std::array< dealii::Tensor< 1, dim, real >, nstate > &solution_gradient, const dealii::types::global_dof_index cell_index) const
Dissipative flux: 0. 
real compute_temperature(const std::array< real, nstate > &primitive_soln) const
Given primitive variables, returns NON-DIMENSIONALIZED temperature using free-stream non-dimensionali...
real compute_total_energy(const std::array< real, nstate > &primitive_soln) const
Given primitive variables, returns total energy. 
std::array< real, nstate > compute_entropy_variables(const std::array< real, nstate > &conservative_soln) const
Computes the entropy variables. 
real compute_dimensional_temperature(const std::array< real, nstate > &primitive_soln) const
Given primitive variables, returns DIMENSIONALIZED temperature using the equation of state...
dealii::Tensor< 1, dim, real > compute_mean_velocities(const std::array< real, nstate > &conservative_soln1, const std::array< real, nstate > &convervative_soln2) const
Mean velocities given two sets of conservative solutions. 
Magnetohydrodynamics (MHD) equations. Derived from PhysicsBase. 
std::array< real, nstate > source_term(const dealii::Point< dim, real > &pos, const std::array< real, nstate > &conservative_soln, const real current_time, const dealii::types::global_dof_index cell_index) const
Source term is zero or depends on manufactured solution. 
real compute_mean_pressure(const std::array< real, nstate > &conservative_soln1, const std::array< real, nstate > &convervative_soln2) const
Mean pressure given two sets of conservative solutions. 
real compute_mach_number(const std::array< real, nstate > &conservative_soln) const
Given conservative variables, returns Mach number. 
std::array< real, nstate > convective_normal_flux(const std::array< real, nstate > &conservative_soln, const dealii::Tensor< 1, dim, real > &normal) const
Convective flux: . 
std::array< dealii::Tensor< 1, dim, real >, nstate > convective_flux(const std::array< real, nstate > &conservative_soln) const
Convective flux: .