| [P]arallel [Hi]gh-order [Li]brary for [P]DEs
    Latest
    Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods | 
Navier-Stokes equations. Derived from Euler for the convective terms, which is derived from PhysicsBase. More...
#include <navier_stokes.h>
| Public Types | |
| using | thermal_boundary_condition_enum = Parameters::NavierStokesParam::ThermalBoundaryCondition | 
| using | two_point_num_flux_enum = Parameters::AllParameters::TwoPointNumericalFlux | 
|  Public Types inherited from PHiLiP::Physics::Euler< dim, nstate, real > | |
| using | two_point_num_flux_enum = Parameters::AllParameters::TwoPointNumericalFlux | 
|  Public Types inherited from PHiLiP::Physics::PhysicsBase< dim, nstate, real > | |
| using | NonPhysicalBehaviorEnum = Parameters::AllParameters::NonPhysicalBehaviorEnum | 
| Public Member Functions | |
| NavierStokes (const Parameters::AllParameters *const parameters_input, const double ref_length, const double gamma_gas, const double mach_inf, const double angle_of_attack, const double side_slip_angle, const double prandtl_number, const double reynolds_number_inf, const bool use_constant_viscosity, const double constant_viscosity, const double temperature_inf=273.15, const double isothermal_wall_temperature=1.0, const thermal_boundary_condition_enum thermal_boundary_condition_type=thermal_boundary_condition_enum::adiabatic, std::shared_ptr< ManufacturedSolutionFunction< dim, real > > manufactured_solution_function=nullptr, const two_point_num_flux_enum two_point_num_flux_type=two_point_num_flux_enum::KG) | |
| Constructor. | |
| template<typename real2 > | |
| std::array< dealii::Tensor< 1, dim, real2 >, nstate > | convert_conservative_gradient_to_primitive_gradient (const std::array< real2, nstate > &conservative_soln, const std::array< dealii::Tensor< 1, dim, real2 >, nstate > &conservative_soln_gradient) const | 
| template<typename real2 > | |
| dealii::Tensor< 1, dim, real2 > | compute_temperature_gradient (const std::array< real2, nstate > &primitive_soln, const std::array< dealii::Tensor< 1, dim, real2 >, nstate > &primitive_soln_gradient) const | 
| template<typename real2 > | |
| real2 | compute_viscosity_coefficient (const std::array< real2, nstate > &primitive_soln) const | 
| template<typename real2 > | |
| real2 | compute_viscosity_coefficient_sutherlands_law (const std::array< real2, nstate > &primitive_soln) const | 
| template<typename real2 > | |
| real2 | scale_viscosity_coefficient (const real2 viscosity_coefficient) const | 
| template<typename real2 > | |
| real2 | compute_scaled_viscosity_coefficient (const std::array< real2, nstate > &primitive_soln) const | 
| template<typename real2 > | |
| real2 | compute_scaled_heat_conductivity_given_scaled_viscosity_coefficient_and_prandtl_number (const real2 scaled_viscosity_coefficient, const double prandtl_number_input) const | 
| template<typename real2 > | |
| real2 | compute_scaled_heat_conductivity (const std::array< real2, nstate > &primitive_soln) const | 
| template<typename real2 > | |
| dealii::Tensor< 1, dim, real2 > | compute_heat_flux (const std::array< real2, nstate > &primitive_soln, const std::array< dealii::Tensor< 1, dim, real2 >, nstate > &primitive_soln_gradient) const | 
| template<typename real2 > | |
| dealii::Tensor< 1, dim, real2 > | compute_heat_flux_given_scaled_heat_conductivity_and_temperature_gradient (const real2 scaled_heat_conductivity, const dealii::Tensor< 1, dim, real2 > &temperature_gradient) const | 
| template<typename real2 > | |
| dealii::Tensor< 1, 3, real2 > | compute_vorticity (const std::array< real2, nstate > &conservative_soln, const std::array< dealii::Tensor< 1, dim, real2 >, nstate > &conservative_soln_gradient) const | 
| Evaluate vorticity from conservative variables and gradient of conservative variables. | |
| real | compute_vorticity_magnitude_sqr (const std::array< real, nstate > &conservative_soln, const std::array< dealii::Tensor< 1, dim, real >, nstate > &conservative_soln_gradient) const | 
| Evaluate vorticity magnitude squared from conservative variables and gradient of conservative variables. | |
| real | compute_vorticity_magnitude (const std::array< real, nstate > &conservative_soln, const std::array< dealii::Tensor< 1, dim, real >, nstate > &conservative_soln_gradient) const | 
| Evaluate vorticity magnitude from conservative variables and gradient of conservative variables. | |
| real | compute_enstrophy (const std::array< real, nstate > &conservative_soln, const std::array< dealii::Tensor< 1, dim, real >, nstate > &conservative_soln_gradient) const | 
| Evaluate enstrophy from conservative variables and gradient of conservative variables. | |
| real | compute_vorticity_based_dissipation_rate_from_integrated_enstrophy (const real integrated_enstrophy) const | 
| real | compute_pressure_dilatation (const std::array< real, nstate > &conservative_soln, const std::array< dealii::Tensor< 1, dim, real >, nstate > &conservative_soln_gradient) const | 
| dealii::Tensor< 2, dim, real > | compute_deviatoric_strain_rate_tensor (const std::array< real, nstate > &conservative_soln, const std::array< dealii::Tensor< 1, dim, real >, nstate > &conservative_soln_gradient) const | 
| real | compute_deviatoric_strain_rate_tensor_magnitude_sqr (const std::array< real, nstate > &conservative_soln, const std::array< dealii::Tensor< 1, dim, real >, nstate > &conservative_soln_gradient) const | 
| Evaluate the square of the deviatoric strain-rate tensor magnitude (i.e. double dot product) from conservative variables and gradient of conservative variables. | |
| real | compute_deviatoric_strain_rate_tensor_based_dissipation_rate_from_integrated_deviatoric_strain_rate_tensor_magnitude_sqr (const real integrated_deviatoric_strain_rate_tensor_magnitude_sqr) const | 
| template<typename real2 > | |
| dealii::Tensor< 2, dim, real2 > | extract_velocities_gradient_from_primitive_solution_gradient (const std::array< dealii::Tensor< 1, dim, real2 >, nstate > &primitive_soln_gradient) const | 
| template<typename real2 > | |
| dealii::Tensor< 2, dim, real2 > | compute_strain_rate_tensor (const dealii::Tensor< 2, dim, real2 > &vel_gradient) const | 
| real | compute_strain_rate_tensor_magnitude_sqr (const std::array< real, nstate > &conservative_soln, const std::array< dealii::Tensor< 1, dim, real >, nstate > &conservative_soln_gradient) const | 
| Evaluate the square of the strain-rate tensor magnitude (i.e. double dot product) from conservative variables and gradient of conservative variables. | |
| real | compute_strain_rate_tensor_based_dissipation_rate_from_integrated_strain_rate_tensor_magnitude_sqr (const real integrated_strain_rate_tensor_magnitude_sqr) const | 
| template<typename real2 > | |
| dealii::Tensor< 2, dim, real2 > | compute_viscous_stress_tensor_via_scaled_viscosity_and_strain_rate_tensor (const real2 scaled_viscosity_coefficient, const dealii::Tensor< 2, dim, real2 > &strain_rate_tensor) const | 
| template<typename real2 > | |
| dealii::Tensor< 2, dim, real2 > | compute_viscous_stress_tensor (const std::array< real2, nstate > &primitive_soln, const std::array< dealii::Tensor< 1, dim, real2 >, nstate > &primitive_soln_gradient) const | 
| 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 override | 
| dealii::Tensor< 1, dim, real > | compute_scaled_viscosity_gradient (const std::array< real, nstate > &primitive_soln, const dealii::Tensor< 1, dim, real > &temperature_gradient) const | 
| dealii::Tensor< 2, nstate, real > | dissipative_flux_directional_jacobian (const std::array< real, nstate > &conservative_soln, const std::array< dealii::Tensor< 1, dim, real >, nstate > &solution_gradient, const dealii::Tensor< 1, dim, real > &normal) const | 
| dealii::Tensor< 2, nstate, real > | dissipative_flux_directional_jacobian_wrt_gradient_component (const std::array< real, nstate > &conservative_soln, const std::array< dealii::Tensor< 1, dim, real >, nstate > &solution_gradient, const dealii::Tensor< 1, dim, real > &normal, const int d_gradient) const | 
| std::array< real, nstate > | dissipative_source_term (const dealii::Point< dim, real > &pos) const | 
| Dissipative flux contribution to the source term. | |
| std::array< real, nstate > | source_term (const dealii::Point< dim, real > &pos, const std::array< real, nstate > &conservative_soln, const real current_time) const override | 
| Source term is zero or depends on manufactured solution. | |
| dealii::Tensor< 2, nstate, real > | convective_flux_directional_jacobian_via_dfad (std::array< real, nstate > &conservative_soln, const dealii::Tensor< 1, dim, real > &normal) const | 
| real | compute_scaled_viscosity_coefficient_derivative_wrt_temperature_via_dfad (std::array< real, nstate > &conservative_soln) const | 
| template<typename real2 > | |
| std::array< dealii::Tensor< 1, dim, real2 >, nstate > | dissipative_flux_given_velocities_viscous_stress_tensor_and_heat_flux (const dealii::Tensor< 1, dim, real2 > &vel, const dealii::Tensor< 2, dim, real2 > &viscous_stress_tensor, const dealii::Tensor< 1, dim, real2 > &heat_flux) const | 
|  Public Member Functions inherited from PHiLiP::Physics::Euler< dim, nstate, real > | |
| Euler (const Parameters::AllParameters *const parameters_input, const double ref_length, const double gamma_gas, const double mach_inf, const double angle_of_attack, const double side_slip_angle, std::shared_ptr< ManufacturedSolutionFunction< dim, real > > manufactured_solution_function=nullptr, const two_point_num_flux_enum two_point_num_flux_type=two_point_num_flux_enum::KG, const bool has_nonzero_diffusion=false, const bool has_nonzero_physical_source=false) | |
| Constructor. | |
| std::array< dealii::Tensor< 1, dim, real >, nstate > | convective_flux (const std::array< real, nstate > &conservative_soln) const override | 
| Convective flux: \( \mathbf{F}_{conv} \). | |
| std::array< real, nstate > | convective_normal_flux (const std::array< real, nstate > &conservative_soln, const dealii::Tensor< 1, dim, real > &normal) const | 
| Convective normal flux: \( \mathbf{F}_{conv} \cdot \hat{n} \). | |
| 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: \( \frac{\partial \mathbf{F}_{conv}}{\partial w} \cdot \mathbf{n} \). | |
| std::array< real, nstate > | convective_eigenvalues (const std::array< real, nstate > &, const dealii::Tensor< 1, dim, real > &) const override | 
| Spectral radius of convective term Jacobian is 'c'. | |
| real | max_convective_eigenvalue (const std::array< real, nstate > &soln) const override | 
| Maximum convective eigenvalue. | |
| real | max_convective_normal_eigenvalue (const std::array< real, nstate > &soln, const dealii::Tensor< 1, dim, real > &normal) const override | 
| Maximum convective normal eigenvalue (used in Lax-Friedrichs)  More... | |
| real | max_viscous_eigenvalue (const std::array< real, nstate > &soln) const override | 
| Maximum viscous eigenvalue. | |
| 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. | |
| 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. | |
| std::array< real, nstate > | convective_source_term (const dealii::Point< dim, real > &pos) const | 
| Convective flux contribution to the source term. | |
| template<typename real2 > | |
| std::array< real2, nstate > | convert_conservative_to_primitive (const std::array< real2, nstate > &conservative_soln) const | 
| std::array< real, nstate > | convert_primitive_to_conservative (const std::array< real, nstate > &primitive_soln) const | 
| template<typename real2 > | |
| real2 | compute_pressure (const std::array< real2, nstate > &conservative_soln) const | 
| Evaluate pressure from conservative variables. | |
| template<typename real2 > | |
| real2 | compute_entropy (const real2 density, const real2 pressure) const | 
| Evaluate physical entropy = log(p ^{-}) from pressure and density. | |
| real | compute_specific_enthalpy (const std::array< real, nstate > &conservative_soln, const real pressure) const | 
| Evaluate pressure from conservative variables. | |
| real | compute_numerical_entropy_function (const std::array< real, nstate > &conservative_soln) const | 
| Compute numerical entropy function -rho s. | |
| real | compute_sound (const std::array< real, nstate > &conservative_soln) const | 
| Evaluate speed of sound from conservative variables. | |
| real | compute_sound (const real density, const real pressure) const | 
| Evaluate speed of sound from density and pressure. | |
| template<typename real2 > | |
| dealii::Tensor< 1, dim, real2 > | compute_velocities (const std::array< real2, nstate > &conservative_soln) const | 
| Evaluate velocities from conservative variables. | |
| template<typename real2 > | |
| real2 | compute_velocity_squared (const dealii::Tensor< 1, dim, real2 > &velocities) const | 
| Given the velocity vector \( \mathbf{u} \), returns the dot-product \( \mathbf{u} \cdot \mathbf{u} \). | |
| template<typename real2 > | |
| dealii::Tensor< 1, dim, real2 > | extract_velocities_from_primitive (const std::array< real2, nstate > &primitive_soln) const | 
| Given primitive variables, returns velocities. | |
| real | compute_total_energy (const std::array< real, nstate > &primitive_soln) const | 
| Given primitive variables, returns total energy.  More... | |
| real | compute_kinetic_energy_from_primitive_solution (const std::array< real, nstate > &primitive_soln) const | 
| Given primitive variables, returns kinetic energy. | |
| real | compute_kinetic_energy_from_conservative_solution (const std::array< real, nstate > &conservative_soln) const | 
| Given conservative variables, returns kinetic energy. | |
| real | compute_entropy_measure (const std::array< real, nstate > &conservative_soln) const | 
| Evaluate entropy from conservative variables.  More... | |
| real | compute_entropy_measure (const real density, const real pressure) const | 
| Evaluate entropy from density and pressure. | |
| real | compute_mach_number (const std::array< real, nstate > &conservative_soln) const | 
| Given conservative variables, returns Mach number. | |
| template<typename real2 > | |
| real2 | compute_temperature (const std::array< real2, nstate > &primitive_soln) const | 
| Given primitive variables, returns NON-DIMENSIONALIZED temperature using free-stream non-dimensionalization.  More... | |
| 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-dimensionalization.  More... | |
| 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-dimensionalization.  More... | |
| real | compute_pressure_from_density_temperature (const real density, const real temperature) const | 
| Given density and temperature, returns NON-DIMENSIONALIZED pressure using free-stream non-dimensionalization.  More... | |
| std::array< dealii::Tensor< 1, dim, real >, nstate > | convective_numerical_split_flux (const std::array< real, nstate > &conservative_soln1, const std::array< real, nstate > &conservative_soln2) const override | 
| Evaluates convective flux based on the chosen split form. | |
| std::array< real, nstate > | compute_entropy_variables (const std::array< real, nstate > &conservative_soln) const | 
| std::array< real, nstate > | compute_conservative_variables_from_entropy_variables (const std::array< real, nstate > &entropy_var) const | 
| std::array< real, nstate > | compute_kinetic_energy_variables (const std::array< real, nstate > &conservative_soln) const | 
| Computes the kinetic energy variables. | |
| 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.  More... | |
| 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.  More... | |
| 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.  More... | |
| real | compute_mean_specific_total_energy (const std::array< real, nstate > &conservative_soln1, const std::array< real, nstate > &convervative_soln2) const | 
| Mean specific total energy given two sets of conservative solutions.  More... | |
| 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 | 
| Boundary condition handler. | |
| virtual dealii::Vector< double > | post_compute_derived_quantities_vector (const dealii::Vector< double > &uh, const std::vector< dealii::Tensor< 1, dim > > &duh, const std::vector< dealii::Tensor< 2, dim > > &dduh, const dealii::Tensor< 1, dim > &normals, const dealii::Point< dim > &evaluation_points) const | 
| For post processing purposes (update comment later) | |
| virtual std::vector< std::string > | post_get_names () const | 
| For post processing purposes, sets the base names (with no prefix or suffix) of the computed quantities. | |
| virtual std::vector< dealii::DataComponentInterpretation::DataComponentInterpretation > | post_get_data_component_interpretation () const | 
| For post processing purposes, sets the interpretation of each computed quantity as either scalar or vector. | |
| virtual dealii::UpdateFlags | post_get_needed_update_flags () const | 
| For post processing purposes (update comment later) | |
|  Public Member Functions inherited from PHiLiP::Physics::PhysicsBase< dim, nstate, real > | |
| PhysicsBase (const Parameters::AllParameters *const parameters_input, const bool has_nonzero_diffusion_input, const bool has_nonzero_physical_source_input, const dealii::Tensor< 2, 3, double > input_diffusion_tensor=Parameters::ManufacturedSolutionParam::get_default_diffusion_tensor(), std::shared_ptr< ManufacturedSolutionFunction< dim, real > > manufactured_solution_function_input=nullptr) | |
| Default constructor that will set the constants. | |
| PhysicsBase (const Parameters::AllParameters *const parameters_input, const bool has_nonzero_diffusion_input, const bool has_nonzero_physical_source_input, std::shared_ptr< ManufacturedSolutionFunction< dim, real > > manufactured_solution_function_input=nullptr) | |
| Constructor that will call default constructor. | |
| virtual | ~PhysicsBase ()=default | 
| Virtual destructor required for abstract classes. | |
| virtual std::array< real, nstate > | physical_source_term (const dealii::Point< dim, real > &pos, const std::array< real, nstate > &solution, const std::array< dealii::Tensor< 1, dim, real >, nstate > &solution_gradient, const dealii::types::global_dof_index cell_index) const | 
| Physical source term that does require differentiation. | |
| virtual std::array< real, nstate > | artificial_source_term (const real viscosity_coefficient, const dealii::Point< dim, real > &pos, const std::array< real, nstate > &solution) const | 
| Artificial source term that does not require differentiation stemming from artificial dissipation. | |
| virtual dealii::Vector< double > | post_compute_derived_quantities_scalar (const double &uh, const dealii::Tensor< 1, dim > &, const dealii::Tensor< 2, dim > &, const dealii::Tensor< 1, dim > &, const dealii::Point< dim > &) const | 
| Returns current scalar solution to be used by PhysicsPostprocessor to output current solution.  More... | |
| template<typename real2 > | |
| real2 | handle_non_physical_result (const std::string message="") const | 
| Function to handle nonphysical results.  More... | |
| Public Attributes | |
| const double | viscosity_coefficient_inf | 
| Nondimensionalized viscosity coefficient at infinity. | |
| const bool | use_constant_viscosity | 
| Flag to use constant viscosity instead of Sutherland's law of viscosity. | |
| const double | constant_viscosity | 
| Nondimensionalized constant viscosity. | |
| const double | prandtl_number | 
| Prandtl number. | |
| const double | reynolds_number_inf | 
| Farfield (free stream) Reynolds number. | |
| const double | isothermal_wall_temperature | 
| Nondimensionalized isothermal wall temperature. | |
| const thermal_boundary_condition_enum | thermal_boundary_condition_type | 
| Thermal boundary condition type (adiabatic or isothermal) | |
|  Public Attributes inherited from PHiLiP::Physics::Euler< dim, nstate, real > | |
| const double | ref_length | 
| Reference length. | |
| const double | gam | 
| Constant heat capacity ratio of fluid. | |
| const double | gamm1 | 
| Constant heat capacity ratio (Gamma-1.0) used often. | |
| const double | density_inf | 
| const double | mach_inf | 
| Farfield Mach number. | |
| const double | mach_inf_sqr | 
| const double | angle_of_attack | 
| Angle of attack.  More... | |
| const double | side_slip_angle | 
| Sideslip angle.  More... | |
| const double | sound_inf | 
| Non-dimensionalized sound* at infinity. | |
| const double | pressure_inf | 
| Non-dimensionalized pressure* at infinity. | |
| const double | entropy_inf | 
| Entropy measure at infinity. | |
| const two_point_num_flux_enum | two_point_num_flux_type | 
| Two point numerical flux type (for split form) | |
| double | temperature_inf | 
| Non-dimensionalized temperature* at infinity. Should equal 1/density*(inf) | |
| double | dynamic_pressure_inf | 
| Non-dimensionalized dynamic pressure* at infinity. | |
| dealii::Tensor< 1, dim, double > | velocities_inf | 
| Non-dimensionalized Velocity vector at farfield.  More... | |
|  Public Attributes inherited from PHiLiP::Physics::PhysicsBase< dim, nstate, real > | |
| const bool | has_nonzero_diffusion | 
| Flag to signal that diffusion term is non-zero. | |
| const bool | has_nonzero_physical_source | 
| Flag to signal that physical source term is non-zero. | |
| const Parameters::AllParameters *const | all_parameters | 
| Pointer to parameters object. | |
| const NonPhysicalBehaviorEnum | non_physical_behavior_type | 
| Determines type of nonphysical behavior. | |
| std::shared_ptr< ManufacturedSolutionFunction< dim, real > > | manufactured_solution_function | 
| Manufactured solution function. | |
| const double | BIG_NUMBER = 1e100 | 
| BIG_NUMBER which is returned in place of NaN according to handle_non_physical_result()  More... | |
| Protected Member Functions | |
| template<typename real2 > | |
| std::array< dealii::Tensor< 1, dim, real2 >, nstate > | dissipative_flux_templated (const std::array< real2, nstate > &conservative_soln, const std::array< dealii::Tensor< 1, dim, real2 >, nstate > &solution_gradient) const | 
| void | boundary_wall (const dealii::Tensor< 1, dim, real > &normal_int, const std::array< real, nstate > &soln_int, const std::array< dealii::Tensor< 1, dim, real >, nstate > &soln_grad_int, std::array< real, nstate > &soln_bc, std::array< dealii::Tensor< 1, dim, real >, nstate > &soln_grad_bc) const override | 
| void | boundary_manufactured_solution (const dealii::Point< dim, real > &pos, const dealii::Tensor< 1, dim, real > &normal_int, const std::array< real, nstate > &soln_int, const std::array< dealii::Tensor< 1, dim, real >, nstate > &soln_grad_int, std::array< real, nstate > &soln_bc, std::array< dealii::Tensor< 1, dim, real >, nstate > &soln_grad_bc) const override | 
| Evaluate the manufactured solution boundary conditions. | |
|  Protected Member Functions inherited from PHiLiP::Physics::Euler< dim, nstate, real > | |
| template<typename real2 > | |
| bool | check_positive_quantity (real2 &quantity, const std::string qty_name) const | 
| Check positive quantity and modify it according to handle_non_physical_result()  More... | |
| void | boundary_slip_wall (const dealii::Tensor< 1, dim, real > &normal_int, const std::array< real, nstate > &soln_int, const std::array< dealii::Tensor< 1, dim, real >, nstate > &soln_grad_int, std::array< real, nstate > &soln_bc, std::array< dealii::Tensor< 1, dim, real >, nstate > &soln_grad_bc) const | 
| void | boundary_pressure_outflow (const real total_inlet_pressure, const real back_pressure, const std::array< real, nstate > &soln_int, std::array< real, nstate > &soln_bc) const | 
| void | boundary_inflow (const real total_inlet_pressure, const real total_inlet_temperature, const dealii::Tensor< 1, dim, real > &normal_int, const std::array< real, nstate > &soln_int, std::array< real, nstate > &soln_bc) const | 
| void | boundary_riemann (const dealii::Tensor< 1, dim, real > &normal_int, const std::array< real, nstate > &soln_int, std::array< real, nstate > &soln_bc) const | 
| void | boundary_farfield (std::array< real, nstate > &soln_bc) const | 
| Simple farfield boundary conditions based on freestream values. | |
| void | boundary_p0_extrapolation (const std::array< real, nstate > &soln_int, std::array< real, nstate > &soln_bc, std::array< dealii::Tensor< 1, dim, real >, nstate > &soln_grad_bc) const | 
| p0 extrapolation at the boundary | |
| void | boundary_custom (std::array< real, nstate > &soln_bc) const | 
| Custom boundary conditions for the left boundary of the astrophysical mach jet case where it is not hypersonic inflow. | |
| void | boundary_astrophysical_inflow (std::array< real, nstate > &soln_bc) const | 
| Boundary conditions based on user-defined values. | |
| std::array< real, nstate > | get_manufactured_solution_value (const dealii::Point< dim, real > &pos) const | 
| Get manufactured solution value. | |
| std::array< dealii::Tensor< 1, dim, real >, nstate > | get_manufactured_solution_gradient (const dealii::Point< dim, real > &pos) const | 
| Get manufactured solution gradient. | |
| std::array< dealii::Tensor< 1, dim, real >, nstate > | convective_numerical_split_flux_kennedy_gruber (const std::array< real, nstate > &conservative_soln1, const std::array< real, nstate > &conservative_soln2) const | 
| std::array< real, nstate > | compute_ismail_roe_parameter_vector_from_primitive (const std::array< real, nstate > &primitive_soln) const | 
| Compute Ismail-Roe parameter vector from primitive solution. | |
| real | compute_ismail_roe_logarithmic_mean (const real val1, const real val2) const | 
| Compute Ismail-Roe logarithmic mean. | |
| std::array< dealii::Tensor< 1, dim, real >, nstate > | convective_numerical_split_flux_ismail_roe (const std::array< real, nstate > &conservative_soln1, const std::array< real, nstate > &conservative_soln2) const | 
| std::array< dealii::Tensor< 1, dim, real >, nstate > | convective_numerical_split_flux_chandrashekar (const std::array< real, nstate > &conservative_soln1, const std::array< real, nstate > &conservative_soln2) const | 
| Chandrashekar entropy conserving flux. | |
| std::array< dealii::Tensor< 1, dim, real >, nstate > | convective_numerical_split_flux_ranocha (const std::array< real, nstate > &conservative_soln1, const std::array< real, nstate > &conservative_soln2) const | 
| Ranocha pressure equilibrium preserving, entropy and energy conserving flux. | |
| Protected Attributes | |
| const double | sutherlands_temperature | 
| Sutherland's temperature. Units: [K].  More... | |
| const double | freestream_temperature | 
| Freestream temperature. Units: [K]. | |
| const double | temperature_ratio | 
| Ratio of Sutherland's temperature to freestream temperature. | |
|  Protected Attributes inherited from PHiLiP::Physics::PhysicsBase< dim, nstate, real > | |
| dealii::ConditionalOStream | pcout | 
| ConditionalOStream.  More... | |
| dealii::Tensor< 2, dim, double > | diffusion_tensor | 
| Anisotropic diffusion matrix.  More... | |
| Private Member Functions | |
| real | get_tensor_magnitude_sqr (const dealii::Tensor< 2, dim, real > &tensor) const | 
| Returns the square of the magnitude of the tensor (i.e. the double dot product of a tensor with itself) | |
Navier-Stokes equations. Derived from Euler for the convective terms, which is derived from PhysicsBase.
Definition at line 12 of file navier_stokes.h.
| 
 | overrideprotectedvirtual | 
No-slip wall boundary conditions
Reimplemented from PHiLiP::Physics::Euler< dim, nstate, real >.
Definition at line 886 of file navier_stokes.cpp.
| dealii::Tensor< 2, dim, real > PHiLiP::Physics::NavierStokes< dim, nstate, real >::compute_deviatoric_strain_rate_tensor | ( | const std::array< real, nstate > & | conservative_soln, | 
| const std::array< dealii::Tensor< 1, dim, real >, nstate > & | conservative_soln_gradient | ||
| ) | const | 
Evaluate the deviatoric strain-rate tensor from conservative variables and gradient of conservative variables – Reference: de la Llave Plata et al. (2019). "On the performance of a high-order multiscale DG approach to LES at increasing Reynolds number."
Definition at line 347 of file navier_stokes.cpp.
| real PHiLiP::Physics::NavierStokes< dim, nstate, real >::compute_deviatoric_strain_rate_tensor_based_dissipation_rate_from_integrated_deviatoric_strain_rate_tensor_magnitude_sqr | ( | const real | integrated_deviatoric_strain_rate_tensor_magnitude_sqr | ) | const | 
Evaluate non-dimensional theoretical deviatoric strain-rate tensor based dissipation rate from integrated deviatoric strain-rate tensor magnitude squared. – Reference: Cox, Christopher, et al. "Accuracy, stability, and performance comparison between the spectral difference and flux reconstruction schemes." Computers & Fluids 221 (2021): 104922. – Equation (57a) with free-stream nondimensionalization applied
Definition at line 404 of file navier_stokes.cpp.
| dealii::Tensor< 1, dim, real2 > PHiLiP::Physics::NavierStokes< dim, nstate, real >::compute_heat_flux | ( | const std::array< real2, nstate > & | primitive_soln, | 
| const std::array< dealii::Tensor< 1, dim, real2 >, nstate > & | primitive_soln_gradient | ||
| ) | const | 
Nondimensionalized heat flux, q* Reference: Masatsuka 2018 "I do like CFD", p.148, eq.(4.14.13)
Definition at line 214 of file navier_stokes.cpp.
| template dealii::Tensor< 1, PHILIP_DIM, FadType > PHiLiP::Physics::NavierStokes< dim, nstate, real >::compute_heat_flux_given_scaled_heat_conductivity_and_temperature_gradient< FadType > | ( | const real2 | scaled_heat_conductivity, | 
| const dealii::Tensor< 1, dim, real2 > & | temperature_gradient | ||
| ) | const | 
Nondimensionalized heat flux, q*, given the scaled heat conductivity and temperature gradient Reference: Masatsuka 2018 "I do like CFD", p.148, eq.(4.14.13)
Definition at line 231 of file navier_stokes.cpp.
| real PHiLiP::Physics::NavierStokes< dim, nstate, real >::compute_pressure_dilatation | ( | const std::array< real, nstate > & | conservative_soln, | 
| const std::array< dealii::Tensor< 1, dim, real >, nstate > & | conservative_soln_gradient | ||
| ) | const | 
Evaluate pressure dilatation from conservative variables and gradient of conservative variables – Reference: Cox, Christopher, et al. "Accuracy, stability, and performance comparison between the spectral difference and flux reconstruction schemes." Computers & Fluids 221 (2021): 104922.
Definition at line 324 of file navier_stokes.cpp.
| 
 | inline | 
Scaled nondimensionalized heat conductivity, hat{kappa*} Reference: Masatsuka 2018 "I do like CFD", p.148, eq.(4.14.13)
Definition at line 199 of file navier_stokes.cpp.
| 
 | inline | 
Scaled nondimensionalized heat conductivity, hat{kappa*}, given scaled nondimensionalized viscosity coefficient and Prandtl number Reference: Masatsuka 2018 "I do like CFD", p.148, eq.(4.14.13)
Definition at line 186 of file navier_stokes.cpp.
| 
 | inline | 
Scaled nondimensionalized viscosity coefficient, hat{mu*} Reference: Masatsuka 2018 "I do like CFD", p.148, eq.(4.14.14)
Definition at line 172 of file navier_stokes.cpp.
| 
 | inline | 
Scaled viscosity coefficient derivative wrt temperature via dfad (automatic differentiation) – Only used for verifying the basic dfad procedure that is extended in convective_flux_directional_jacobian_via_dfad()
Definition at line 800 of file navier_stokes.cpp.
| dealii::Tensor< 1, dim, real > PHiLiP::Physics::NavierStokes< dim, nstate, real >::compute_scaled_viscosity_gradient | ( | const std::array< real, nstate > & | primitive_soln, | 
| const dealii::Tensor< 1, dim, real > & | temperature_gradient | ||
| ) | const | 
Gradient of the scaled nondimensionalized viscosity coefficient Reference: Masatsuka 2018 "I do like CFD", p.148, eq.(4.14.14 and 4.14.17)
Definition at line 540 of file navier_stokes.cpp.
| template dealii::Tensor< 2, PHILIP_DIM, FadType > PHiLiP::Physics::NavierStokes< dim, nstate, real >::compute_strain_rate_tensor< FadType > | ( | const dealii::Tensor< 2, dim, real2 > & | vel_gradient | ) | const | 
Nondimensionalized strain rate tensor, S* Reference: Masatsuka 2018 "I do like CFD", p.148, extracted from eq.(4.14.12)
Definition at line 429 of file navier_stokes.cpp.
| real PHiLiP::Physics::NavierStokes< dim, nstate, real >::compute_strain_rate_tensor_based_dissipation_rate_from_integrated_strain_rate_tensor_magnitude_sqr | ( | const real | integrated_strain_rate_tensor_magnitude_sqr | ) | const | 
Evaluate non-dimensional theoretical strain-rate tensor based dissipation rate from integrated strain-rate tensor magnitude squared. – Reference: Navah, Farshad, et al. "A High-Order Variational Multiscale Approach to Turbulence for Compact Nodal Schemes." – Equation (E.9) with free-stream nondimensionalization applied
Definition at line 463 of file navier_stokes.cpp.
| template dealii::Tensor< 1, PHILIP_DIM, FadType > PHiLiP::Physics::NavierStokes< dim, nstate, real >::compute_temperature_gradient< FadType > | ( | const std::array< real2, nstate > & | primitive_soln, | 
| const std::array< dealii::Tensor< 1, dim, real2 >, nstate > & | primitive_soln_gradient | ||
| ) | const | 
Nondimensionalized temperature gradient
Definition at line 107 of file navier_stokes.cpp.
| 
 | inline | 
Nondimensionalized viscosity coefficient, mu* Based on the use_constant_viscosity flag, it returns a value based on either: (1) Sutherland's viscosity law, or (2) Constant nondimensionalized viscosity value
Definition at line 124 of file navier_stokes.cpp.
| 
 | inline | 
Nondimensionalized viscosity coefficient, mu* Reference: Masatsuka 2018 "I do like CFD", p.148, eq.(4.14.16)
Based on Sutherland's law for viscosity
Definition at line 140 of file navier_stokes.cpp.
| dealii::Tensor< 2, dim, real2 > PHiLiP::Physics::NavierStokes< dim, nstate, real >::compute_viscous_stress_tensor | ( | const std::array< real2, nstate > & | primitive_soln, | 
| const std::array< dealii::Tensor< 1, dim, real2 >, nstate > & | primitive_soln_gradient | ||
| ) | const | 
Nondimensionalized viscous stress tensor, tau* Reference: Masatsuka 2018 "I do like CFD", p.148, eq.(4.14.12)
Definition at line 507 of file navier_stokes.cpp.
| template dealii::Tensor< 2, PHILIP_DIM, FadType > PHiLiP::Physics::NavierStokes< dim, nstate, real >::compute_viscous_stress_tensor_via_scaled_viscosity_and_strain_rate_tensor< FadType > | ( | const real2 | scaled_viscosity_coefficient, | 
| const dealii::Tensor< 2, dim, real2 > & | strain_rate_tensor | ||
| ) | const | 
Nondimensionalized viscous stress tensor, tau* Reference: Masatsuka 2018 "I do like CFD", p.148, eq.(4.14.12)
Definition at line 473 of file navier_stokes.cpp.
| real PHiLiP::Physics::NavierStokes< dim, nstate, real >::compute_vorticity_based_dissipation_rate_from_integrated_enstrophy | ( | const real | integrated_enstrophy | ) | const | 
Evaluate non-dimensional theoretical vorticity-based dissipation rate integrated enstrophy. Note: For incompressible flows or when dilatation effects are negligible – Reference: Cox, Christopher, et al. "Accuracy, stability, and performance comparison between the spectral difference and flux reconstruction schemes." Computers & Fluids 221 (2021): 104922. – Equation (56) with free-stream nondimensionalization applied
Definition at line 315 of file navier_stokes.cpp.
| dealii::Tensor< 2, nstate, real > PHiLiP::Physics::NavierStokes< dim, nstate, real >::convective_flux_directional_jacobian_via_dfad | ( | std::array< real, nstate > & | conservative_soln, | 
| const dealii::Tensor< 1, dim, real > & | normal | ||
| ) | const | 
Convective flux Jacobian computed via dfad (automatic differentiation) – Only used for verifying the dfad procedure used in dissipative flux jacobian
Definition at line 749 of file navier_stokes.cpp.
| template std::array< dealii::Tensor< 1, PHILIP_DIM, FadType >, PHILIP_DIM+2 > PHiLiP::Physics::NavierStokes< dim, nstate, real >::convert_conservative_gradient_to_primitive_gradient< FadType > | ( | const std::array< real2, nstate > & | conservative_soln, | 
| const std::array< dealii::Tensor< 1, dim, real2 >, nstate > & | conservative_soln_gradient | ||
| ) | const | 
Obtain gradient of primitive variables from gradient of conservative variables
Definition at line 59 of file navier_stokes.cpp.
| 
 | overridevirtual | 
Nondimensionalized viscous flux (i.e. dissipative flux) Reference: Masatsuka 2018 "I do like CFD", p.142, eq.(4.12.1-4.12.4)
Reimplemented from PHiLiP::Physics::Euler< dim, nstate, real >.
Definition at line 527 of file navier_stokes.cpp.
| dealii::Tensor< 2, nstate, real > PHiLiP::Physics::NavierStokes< dim, nstate, real >::dissipative_flux_directional_jacobian | ( | const std::array< real, nstate > & | conservative_soln, | 
| const std::array< dealii::Tensor< 1, dim, real >, nstate > & | solution_gradient, | ||
| const dealii::Tensor< 1, dim, real > & | normal | ||
| ) | const | 
Dissipative flux Jacobian Note: Only used for computing the manufactured solution source term; computed using automatic differentiation
Definition at line 584 of file navier_stokes.cpp.
| dealii::Tensor< 2, nstate, real > PHiLiP::Physics::NavierStokes< dim, nstate, real >::dissipative_flux_directional_jacobian_wrt_gradient_component | ( | const std::array< real, nstate > & | conservative_soln, | 
| const std::array< dealii::Tensor< 1, dim, real >, nstate > & | solution_gradient, | ||
| const dealii::Tensor< 1, dim, real > & | normal, | ||
| const int | d_gradient | ||
| ) | const | 
Dissipative flux Jacobian wrt gradient component Note: Only used for computing the manufactured solution source term; computed using automatic differentiation
Definition at line 622 of file navier_stokes.cpp.
| template std::array< dealii::Tensor< 1, PHILIP_DIM, FadType >, PHILIP_DIM+2 > PHiLiP::Physics::NavierStokes< dim, nstate, real >::dissipative_flux_given_velocities_viscous_stress_tensor_and_heat_flux< FadType > | ( | const dealii::Tensor< 1, dim, real2 > & | vel, | 
| const dealii::Tensor< 2, dim, real2 > & | viscous_stress_tensor, | ||
| const dealii::Tensor< 1, dim, real2 > & | heat_flux | ||
| ) | const | 
Nondimensionalized viscous flux (i.e. dissipative flux) computed via given velocities, viscous stress tensor, and heat flux. Reference: Masatsuka 2018 "I do like CFD", p.142, eq.(4.12.1-4.12.4)
Definition at line 854 of file navier_stokes.cpp.
| 
 | protected | 
Nondimensionalized viscous flux (i.e. dissipative flux) Reference: Masatsuka 2018 "I do like CFD", p.142, eq.(4.12.1-4.12.4)
Definition at line 827 of file navier_stokes.cpp.
| template dealii::Tensor< 2, PHILIP_DIM, FadType > PHiLiP::Physics::NavierStokes< dim, nstate, real >::extract_velocities_gradient_from_primitive_solution_gradient< FadType > | ( | const std::array< dealii::Tensor< 1, dim, real2 >, nstate > & | primitive_soln_gradient | ) | const | 
Extract gradient of velocities
Definition at line 414 of file navier_stokes.cpp.
| 
 | inline | 
Scaled nondimensionalized viscosity coefficient, hat{mu*}, given nondimensionalized viscosity coefficient Reference: Masatsuka 2018 "I do like CFD", p.148, eq.(4.14.14)
Definition at line 159 of file navier_stokes.cpp.
| 
 | protected | 
Sutherland's temperature. Units: [K].
Constants for Sutherland's law for viscosity Reference: Sutherland, W. (1893), "The viscosity of gases and molecular force", Philosophical Magazine, S. 5, 36, pp. 507-531 (1893) Values: https://www.cfd-online.com/Wiki/Sutherland%27s_law
Definition at line 66 of file navier_stokes.h.