[P]arallel [Hi]gh-order [Li]brary for [P]DEs  Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
PHiLiP::Physics::ReynoldsAveragedNavierStokes_SAneg< dim, nstate, real > Class Template Reference

Negative Spalart-Allmaras model. Derived from Reynolds Averaged Navier Stokes. More...

#include <negative_spalart_allmaras_rans_model.h>

Inheritance diagram for PHiLiP::Physics::ReynoldsAveragedNavierStokes_SAneg< dim, nstate, real >:
Collaboration diagram for PHiLiP::Physics::ReynoldsAveragedNavierStokes_SAneg< dim, nstate, real >:

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::ReynoldsAveragedNavierStokesBase< dim, nstate, real >
using thermal_boundary_condition_enum = Parameters::NavierStokesParam::ThermalBoundaryCondition
 
using two_point_num_flux_enum = Parameters::AllParameters::TwoPointNumericalFlux
 

Public Member Functions

 ReynoldsAveragedNavierStokes_SAneg (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 turbulent_prandtl_number, 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)
 
dealii::Tensor< 2, dim, real > compute_Reynolds_stress_tensor (const std::array< real, nstate_navier_stokes > &primitive_soln_rans, const std::array< dealii::Tensor< 1, dim, real >, nstate_navier_stokes > &primitive_soln_gradient_rans, const std::array< real, nstate_turbulence_model > &primitive_soln_turbulence_model) const
 Nondimensionalized Reynolds stress tensor, (tau^reynolds)*, for the negative SA model.
 
dealii::Tensor< 1, dim, real > compute_Reynolds_heat_flux (const std::array< real, nstate_navier_stokes > &primitive_soln_rans, const std::array< dealii::Tensor< 1, dim, real >, nstate_navier_stokes > &primitive_soln_gradient_rans, const std::array< real, nstate_turbulence_model > &primitive_soln_turbulence_model) const
 Nondimensionalized Reynolds heat flux, (q^reynolds)*, for the negative SA model.
 
dealii::Tensor< 2, dim, FadTypecompute_Reynolds_stress_tensor_fad (const std::array< FadType, nstate_navier_stokes > &primitive_soln_rans, const std::array< dealii::Tensor< 1, dim, FadType >, nstate_navier_stokes > &primitive_soln_gradient_rans, const std::array< FadType, nstate_turbulence_model > &primitive_soln_turbulence_model) const
 Nondimensionalized Reynolds stress tensor, (tau^reynolds)* (Automatic Differentiation Type: FadType), for the negative SA model.
 
dealii::Tensor< 1, dim, FadTypecompute_Reynolds_heat_flux_fad (const std::array< FadType, nstate_navier_stokes > &primitive_soln_rans, const std::array< dealii::Tensor< 1, dim, FadType >, nstate_navier_stokes > &primitive_soln_gradient_rans, const std::array< FadType, nstate_turbulence_model > &primitive_soln_turbulence_model) const
 Nondimensionalized Reynolds heat flux, (q^reynolds)* (Automatic Differentiation Type: FadType), for the negative SA model.
 
real compute_eddy_viscosity (const std::array< real, nstate_navier_stokes > &primitive_soln_rans, const std::array< real, nstate_turbulence_model > &primitive_soln_turbulence_model) const
 Nondimensionalized eddy viscosity for the negative SA model.
 
FadType compute_eddy_viscosity_fad (const std::array< FadType, nstate_navier_stokes > &primitive_soln_rans, const std::array< FadType, nstate_turbulence_model > &primitive_soln_turbulence_model) const
 Nondimensionalized eddy viscosity for the negative SA model (Automatic Differentiation Type: FadType)
 
std::array< real, nstate-(dim+2)> compute_effective_viscosity_turbulence_model (const std::array< real, nstate_navier_stokes > &primitive_soln_rans, const std::array< real, nstate_turbulence_model > &primitive_soln_turbulence_model) const
 Nondimensionalized effective (total) viscosities for the negative SA model.
 
std::array< FadType, nstate-(dim+2)> compute_effective_viscosity_turbulence_model_fad (const std::array< FadType, nstate_navier_stokes > &primitive_soln_rans, const std::array< FadType, nstate_turbulence_model > &primitive_soln_turbulence_model) const
 Nondimensionalized effective (total) viscosities for the negative SA model (Automatic Differentiation Type: FadType)
 
std::array< real, nstate > compute_production_dissipation_cross_term (const dealii::Point< dim, real > &pos, const std::array< real, nstate > &conservative_solution, const std::array< dealii::Tensor< 1, dim, real >, nstate > &solution_gradient) const
 Physical source term (production, dissipation source terms and source term with cross derivatives) for the negative SA model.
 
dealii::Vector< double > post_compute_derived_quantities_vector (const dealii::Vector< double > &uh, const std::vector< dealii::Tensor< 1, dim > > &, const std::vector< dealii::Tensor< 2, dim > > &, const dealii::Tensor< 1, dim > &, const dealii::Point< dim > &) const override
 For post processing purposes, returns conservative and primitive solution variables for the negative SA model.
 
std::vector< dealii::DataComponentInterpretation::DataComponentInterpretation > post_get_data_component_interpretation () const override
 For post processing purposes, sets the interpretation of each computed quantity as either scalar or vector for the negative SA model.
 
std::vector< std::string > post_get_names () const override
 For post processing purposes, sets the base names (with no prefix or suffix) of the computed quantities for the negative SA model.
 
- Public Member Functions inherited from PHiLiP::Physics::ReynoldsAveragedNavierStokesBase< dim, nstate, real >
 ReynoldsAveragedNavierStokesBase (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 turbulent_prandtl_number, 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.
 
std::array< dealii::Tensor< 1, dim, real >, nstate > convective_flux (const std::array< real, nstate > &conservative_soln) const
 Additional convective flux of RANS + convective flux of turbulence model.
 
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
 Additional viscous flux of RANS + viscous flux of turbulence model.
 
std::array< real, nstate > convective_eigenvalues (const std::array< real, nstate > &, const dealii::Tensor< 1, dim, real > &) const
 Convective eigenvalues of the additional models' PDEs. More...
 
real max_convective_eigenvalue (const std::array< real, nstate > &soln) const
 Maximum convective eigenvalue of the additional models' PDEs. More...
 
real max_convective_normal_eigenvalue (const std::array< real, nstate > &soln, const dealii::Tensor< 1, dim, real > &normal) const
 Maximum convective normal eigenvalue (used in Lax-Friedrichs) of the additional models' PDEs. More...
 
std::array< real, nstate > source_term (const dealii::Point< dim, real > &pos, const std::array< real, nstate > &conservative_solution, const real current_time, const dealii::types::global_dof_index cell_index) const
 Source term for manufactured solution functions.
 
std::array< real, nstate > convective_dissipative_source_term (const dealii::Point< dim, real > &pos, const std::array< real, nstate > &conservative_solution, const dealii::types::global_dof_index cell_index) const
 Convective and dissipative source term for manufactured solution functions.
 
std::array< real, nstate > physical_source_term (const dealii::Point< dim, real > &pos, const std::array< real, nstate > &conservative_solution, const std::array< dealii::Tensor< 1, dim, real >, nstate > &solution_gradient, const dealii::types::global_dof_index cell_index) const override
 Physical source term.
 
- Public Member Functions inherited from PHiLiP::Physics::ModelBase< dim, nstate, real >
 ModelBase (std::shared_ptr< ManufacturedSolutionFunction< dim, real > > manufactured_solution_function_input=nullptr)
 Constructor.
 
virtual ~ModelBase ()=default
 Virtual destructor required for abstract classes.
 
void boundary_face_values (const int boundary_type, const dealii::Point< dim, real > &pos, const dealii::Tensor< 1, dim, real > &normal, 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
 Boundary condition handler.
 

Static Public Attributes

static const int nstate_navier_stokes = dim+2
 Number of PDEs for RANS equations.
 
static const int nstate_turbulence_model = nstate-(dim+2)
 Number of PDEs for RANS turbulence model.
 
- Static Public Attributes inherited from PHiLiP::Physics::ReynoldsAveragedNavierStokesBase< dim, nstate, real >
static const int nstate_navier_stokes = dim+2
 Number of PDEs for RANS equations.
 
static const int nstate_turbulence_model = nstate-(dim+2)
 Number of PDEs for RANS turbulence model.
 

Protected Member Functions

template<typename real2 >
dealii::Tensor< 2, dim, real2 > compute_Reynolds_stress_tensor_templated (const std::array< real2, nstate_navier_stokes > &primitive_soln_rans, const std::array< dealii::Tensor< 1, dim, real2 >, nstate_navier_stokes > &primitive_soln_gradient_rans, const std::array< real2, nstate_turbulence_model > &primitive_soln_turbulence_model) const
 Templated nondimensionalized Reynolds stress tensor, (tau^reynolds)* for the negative SA model.
 
template<typename real2 >
dealii::Tensor< 1, dim, real2 > compute_Reynolds_heat_flux_templated (const std::array< real2, nstate_navier_stokes > &primitive_soln_rans, const std::array< dealii::Tensor< 1, dim, real2 >, nstate_navier_stokes > &primitive_soln_gradient_rans, const std::array< real2, nstate_turbulence_model > &primitive_soln_turbulence_model) const
 Templated nondimensionalized Reynolds heat flux, (q^reynolds)* for the negative SA model.
 
template<typename real2 >
real2 scale_coefficient (const real2 coefficient) const
 Templated nondimensionalized variables scaled by reynolds_number_inf for the negative SA model.
 
void boundary_wall (std::array< real, nstate > &soln_bc, std::array< dealii::Tensor< 1, dim, real >, nstate > &soln_grad_bc) const override
 Wall boundary condition. More...
 
void boundary_outflow (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
 Outflow boundary condition.
 
void boundary_inflow (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
 Inflow boundary condition. More...
 
void boundary_farfield (std::array< real, nstate > &soln_bc) const override
 Farfield boundary conditions based on freestream values. More...
 
- Protected Member Functions inherited from PHiLiP::Physics::ReynoldsAveragedNavierStokesBase< dim, nstate, real >
template<typename real2 >
real2 get_vector_magnitude_sqr (const dealii::Tensor< 1, 3, real2 > &vector) const
 Returns the square of the magnitude of the vector.
 
template<typename real2 >
real2 get_tensor_magnitude_sqr (const dealii::Tensor< 2, dim, real2 > &tensor) const
 Returns the square of the magnitude of the tensor (i.e. the double dot product of a tensor with itself)
 
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 dealii::types::global_dof_index cell_index) const
 Templated additional dissipative (i.e. viscous) flux.
 
template<typename real2 >
std::array< dealii::Tensor< 1, dim, real2 >, nstate > convective_flux_templated (const std::array< real2, nstate > &conservative_soln) const
 Templated additional convective flux.
 
template<typename real2 >
std::array< real2, dim+2 > extract_rans_conservative_solution (const std::array< real2, nstate > &conservative_soln) const
 Returns the conservative solutions of Reynolds-averaged Navier-Stokes equations (without additional RANS turbulence model)
 
template<typename real2 >
std::array< dealii::Tensor< 1, dim, real2 >, dim+2 > extract_rans_solution_gradient (const std::array< dealii::Tensor< 1, dim, real2 >, nstate > &solution_gradient) const
 Returns the conservative solutions gradient of Reynolds-averaged Navier-Stokes equations (without additional RANS turbulence model)
 
template<typename real2 >
std::array< dealii::Tensor< 1, dim, real2 >, nstate-(dim+2)> dissipative_flux_turbulence_model (const std::array< real2, nstate_navier_stokes > &primitive_soln_rans, const std::array< real2, nstate_turbulence_model > &primitive_soln_turbulence_model, const std::array< dealii::Tensor< 1, dim, real2 >, nstate_turbulence_model > &primitive_solution_gradient_turbulence_model) const
 Templated Additional viscous flux of RANS + viscous flux of turbulence model.
 
template<typename real2 >
std::array< real2, nstate-(dim+2)> convert_conservative_to_primitive_turbulence_model (const std::array< real2, nstate > &conservative_soln) const
 
template<typename real2 >
std::array< dealii::Tensor< 1, dim, real2 >, nstate-(dim+2)> convert_conservative_gradient_to_primitive_gradient_turbulence_model (const std::array< real2, nstate > &conservative_soln, const std::array< dealii::Tensor< 1, dim, real2 >, nstate > &solution_gradient) const
 
std::array< real, nstate-(dim+2)> compute_mean_turbulence_property (const std::array< real, nstate > &conservative_soln1, const std::array< real, nstate > &conservative_soln2) const
 Mean turbulence properties given two sets of conservative solutions. More...
 
dealii::Tensor< 2, nstate, real > convective_flux_directional_jacobian (const std::array< real, nstate > &conservative_soln, const dealii::Tensor< 1, dim, real > &normal) 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::types::global_dof_index cell_index) 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 dealii::types::global_dof_index cell_index) const
 
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 value.
 
std::array< real, nstate > convective_source_term_computed_from_manufactured_solution (const dealii::Point< dim, real > &pos) const
 
std::array< real, nstate > dissipative_source_term_computed_from_manufactured_solution (const dealii::Point< dim, real > &pos, const dealii::types::global_dof_index cell_index) const
 
std::array< real, nstate > physical_source_term_computed_from_manufactured_solution (const dealii::Point< dim, real > &pos, const dealii::types::global_dof_index cell_index) const
 
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::ModelBase< dim, nstate, real >
virtual 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
 
virtual 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
 

Private Member Functions

template<typename real2 >
real2 compute_eddy_viscosity_templated (const std::array< real2, nstate_navier_stokes > &primitive_soln_rans, const std::array< real2, nstate_turbulence_model > &primitive_soln_turbulence_model) const
 Templated nondimensionalized eddy viscosity for the negative SA model. More...
 
template<typename real2 >
std::array< real2, nstate-(dim+2)> compute_effective_viscosity_turbulence_model_templated (const std::array< real2, nstate_navier_stokes > &primitive_soln_rans, const std::array< real2, nstate_turbulence_model > &primitive_soln_turbulence_model) const
 Templated nondimensionalized effective (total) viscosities for the negative SA model. More...
 
template<typename real2 >
real2 compute_coefficient_Chi (const real2 nu_tilde, const real2 laminar_kinematic_viscosity) const
 Templated coefficient Chi for the negative SA model. More...
 
template<typename real2 >
real2 compute_coefficient_f_v1 (const real2 coefficient_Chi) const
 Templated coefficient f_v1 for the negative SA model. More...
 
real compute_coefficient_f_v2 (const real coefficient_Chi) const
 Coefficient f_v2 for the negative SA model. More...
 
template<typename real2 >
real2 compute_coefficient_f_n (const real2 nu_tilde, const real2 laminar_kinematic_viscosity) const
 Templated coefficient f_n for the negative SA model. More...
 
real compute_coefficient_f_t2 (const real coefficient_Chi) const
 Coefficient f_t2 for the negative SA model. More...
 
real compute_coefficient_f_w (const real coefficient_g) const
 Coefficient f_t2 for the negative SA model. More...
 
real compute_coefficient_r (const real nu_tilde, const real d_wall, const real s_tilde) const
 Coefficient r for the negative SA model. More...
 
real compute_coefficient_g (const real coefficient_r) const
 Coefficient g for the negative SA model. More...
 
real compute_s (const std::array< real, nstate_navier_stokes > &conservative_soln_rans, const std::array< dealii::Tensor< 1, dim, real >, nstate_navier_stokes > &conservative_soln_gradient_rans) const
 Vorticity magnitude for the negative SA model, sqrt(vorticity_x^2+vorticity_y^2+vorticity_z^2)
 
real compute_s_bar (const real coefficient_Chi, const real nu_tilde, const real d_wall) const
 Correction of vorticity magnitude for the negative SA model. More...
 
real compute_s_tilde (const real coefficient_Chi, const real nu_tilde, const real d_wall, const real s) const
 Modified vorticity magnitude for the negative SA model. More...
 
std::array< real, nstate > compute_production_source (const real coefficient_f_t2, const real density, const real nu_tilde, const real s, const real s_tilde) const
 Production source term for the negative SA model. More...
 
std::array< real, nstate > compute_dissipation_source (const real coefficient_f_t2, const real density, const real nu_tilde, const real d_wall, const real s_tilde) const
 Dissipation source term for the negative SA model. More...
 
std::array< real, nstate > compute_cross_source (const real density, const real nu_tilde, const real laminar_kinematic_viscosity, const std::array< dealii::Tensor< 1, dim, real >, nstate_navier_stokes > &primitive_soln_gradient_rans, const std::array< dealii::Tensor< 1, dim, real >, nstate_turbulence_model > &primitive_solution_gradient_turbulence_model) const
 Source term with cross derivatives for the negative SA model. More...
 

Private Attributes

const real c_b1 = 0.1355
 Constant coefficients for the negative SA model. More...
 
const real sigma = 2.0/3.0
 
const real c_b2 = 0.622
 
const real kappa = 0.41
 
const real c_w1 = c_b1/(kappa*kappa)+(1+c_b2)/sigma
 
const real c_w2 = 0.3
 
const real c_w3 = 2.0
 
const real c_v1 = 7.1
 
const real c_v2 = 0.7
 
const real c_v3 = 0.9
 
const real c_t3 = 1.2
 
const real c_t4 = 0.5
 
const real c_n1 = 16.0
 
const real r_lim = 10.0
 
const FadType sigma_fad = 2.0/3.0
 
const FadType c_v1_fad = 7.1
 
const FadType c_n1_fad = 16.0
 

Additional Inherited Members

- Public Attributes inherited from PHiLiP::Physics::ReynoldsAveragedNavierStokesBase< dim, nstate, real >
const double turbulent_prandtl_number
 Turbulent Prandtl number.
 
std::unique_ptr< NavierStokes< dim, nstate_navier_stokes, real > > navier_stokes_physics
 Pointer to Navier-Stokes physics object.
 
- Public Attributes inherited from PHiLiP::Physics::ModelBase< dim, nstate, real >
std::shared_ptr< ManufacturedSolutionFunction< dim, real > > manufactured_solution_function
 Manufactured solution function.
 
dealii::LinearAlgebra::distributed::Vector< int > cellwise_poly_degree
 Cellwise polynomial degree.
 
dealii::LinearAlgebra::distributed::Vector< double > cellwise_volume
 
- Protected Attributes inherited from PHiLiP::Physics::ModelBase< dim, nstate, real >
const MPI_Comm mpi_communicator
 MPI communicator.
 
dealii::ConditionalOStream pcout
 Parallel std::cout that only outputs on mpi_rank==0.
 

Detailed Description

template<int dim, int nstate, typename real>
class PHiLiP::Physics::ReynoldsAveragedNavierStokes_SAneg< dim, nstate, real >

Negative Spalart-Allmaras model. Derived from Reynolds Averaged Navier Stokes.

Definition at line 14 of file negative_spalart_allmaras_rans_model.h.

Constructor & Destructor Documentation

◆ ReynoldsAveragedNavierStokes_SAneg()

template<int dim, int nstate, typename real >
PHiLiP::Physics::ReynoldsAveragedNavierStokes_SAneg< dim, nstate, real >::ReynoldsAveragedNavierStokes_SAneg ( 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  turbulent_prandtl_number,
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 for the Reynolds-averaged Navier-Stokes model: negative SA Reference: Steven R. Allmaras. (2012). "Modifications and Clarifications for the Implementation of the Spalart-Allmaras Turbulence Model."

Definition at line 18 of file negative_spalart_allmaras_rans_model.cpp.

Member Function Documentation

◆ boundary_farfield()

template<int dim, int nstate, typename real >
void PHiLiP::Physics::ReynoldsAveragedNavierStokes_SAneg< dim, nstate, real >::boundary_farfield ( std::array< real, nstate > &  soln_bc) const
overrideprotectedvirtual

Farfield boundary conditions based on freestream values.

Reference: Steven R. Allmaras. (2012). "Modifications and Clarifications for the Implementation of the Spalart-Allmaras Turbulence Model." eq.(8)

Reimplemented from PHiLiP::Physics::ModelBase< dim, nstate, real >.

Definition at line 646 of file negative_spalart_allmaras_rans_model.cpp.

◆ boundary_inflow()

template<int dim, int nstate, typename real >
void PHiLiP::Physics::ReynoldsAveragedNavierStokes_SAneg< dim, nstate, real >::boundary_inflow ( 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
overrideprotectedvirtual

Inflow boundary condition.

Reference: Philippe R. Spalart. (2007). "Effective Inflow Conditions for Turbulence Models in Aerodynamc Calculations." Section III. Ambient Values: A. Turbulent Reynolds Numbers

Reimplemented from PHiLiP::Physics::ModelBase< dim, nstate, real >.

Definition at line 629 of file negative_spalart_allmaras_rans_model.cpp.

◆ boundary_wall()

template<int dim, int nstate, typename real >
void PHiLiP::Physics::ReynoldsAveragedNavierStokes_SAneg< dim, nstate, real >::boundary_wall ( std::array< real, nstate > &  soln_bc,
std::array< dealii::Tensor< 1, dim, real >, nstate > &  soln_grad_bc 
) const
overrideprotectedvirtual

Wall boundary condition.

Reference: Steven R. Allmaras. (2012). "Modifications and Clarifications for the Implementation of the Spalart-Allmaras Turbulence Model." eq.(7)

Reimplemented from PHiLiP::Physics::ModelBase< dim, nstate, real >.

Definition at line 601 of file negative_spalart_allmaras_rans_model.cpp.

◆ compute_coefficient_Chi()

template<int dim, int nstate, typename real >
template<typename real2 >
real2 PHiLiP::Physics::ReynoldsAveragedNavierStokes_SAneg< dim, nstate, real >::compute_coefficient_Chi ( const real2  nu_tilde,
const real2  laminar_kinematic_viscosity 
) const
private

Templated coefficient Chi for the negative SA model.

Reference: Steven R. Allmaras. (2012). "Modifications and Clarifications for the Implementation of the Spalart-Allmaras Turbulence Model." eq.(1)

Definition at line 168 of file negative_spalart_allmaras_rans_model.cpp.

◆ compute_coefficient_f_n()

template<int dim, int nstate, typename real >
template<typename real2 >
real2 PHiLiP::Physics::ReynoldsAveragedNavierStokes_SAneg< dim, nstate, real >::compute_coefficient_f_n ( const real2  nu_tilde,
const real2  laminar_kinematic_viscosity 
) const
private

Templated coefficient f_n for the negative SA model.

Reference: Steven R. Allmaras. (2012). "Modifications and Clarifications for the Implementation of the Spalart-Allmaras Turbulence Model." eq.(21)

Definition at line 219 of file negative_spalart_allmaras_rans_model.cpp.

◆ compute_coefficient_f_t2()

template<int dim, int nstate, typename real >
real PHiLiP::Physics::ReynoldsAveragedNavierStokes_SAneg< dim, nstate, real >::compute_coefficient_f_t2 ( const real  coefficient_Chi) const
private

Coefficient f_t2 for the negative SA model.

Reference: Steven R. Allmaras. (2012). "Modifications and Clarifications for the Implementation of the Spalart-Allmaras Turbulence Model." eq.(6)

Definition at line 250 of file negative_spalart_allmaras_rans_model.cpp.

◆ compute_coefficient_f_v1()

template<int dim, int nstate, typename real >
template<typename real2 >
real2 PHiLiP::Physics::ReynoldsAveragedNavierStokes_SAneg< dim, nstate, real >::compute_coefficient_f_v1 ( const real2  coefficient_Chi) const
private

Templated coefficient f_v1 for the negative SA model.

Reference: Steven R. Allmaras. (2012). "Modifications and Clarifications for the Implementation of the Spalart-Allmaras Turbulence Model." eq.(1)

Definition at line 181 of file negative_spalart_allmaras_rans_model.cpp.

◆ compute_coefficient_f_v2()

template<int dim, int nstate, typename real >
real PHiLiP::Physics::ReynoldsAveragedNavierStokes_SAneg< dim, nstate, real >::compute_coefficient_f_v2 ( const real  coefficient_Chi) const
private

Coefficient f_v2 for the negative SA model.

Reference: Steven R. Allmaras. (2012). "Modifications and Clarifications for the Implementation of the Spalart-Allmaras Turbulence Model." eq.(4)

Definition at line 204 of file negative_spalart_allmaras_rans_model.cpp.

◆ compute_coefficient_f_w()

template<int dim, int nstate, typename real >
real PHiLiP::Physics::ReynoldsAveragedNavierStokes_SAneg< dim, nstate, real >::compute_coefficient_f_w ( const real  coefficient_g) const
private

Coefficient f_t2 for the negative SA model.

Reference: Steven R. Allmaras. (2012). "Modifications and Clarifications for the Implementation of the Spalart-Allmaras Turbulence Model." eq.(5)

Definition at line 263 of file negative_spalart_allmaras_rans_model.cpp.

◆ compute_coefficient_g()

template<int dim, int nstate, typename real >
real PHiLiP::Physics::ReynoldsAveragedNavierStokes_SAneg< dim, nstate, real >::compute_coefficient_g ( const real  coefficient_r) const
private

Coefficient g for the negative SA model.

Reference: Steven R. Allmaras. (2012). "Modifications and Clarifications for the Implementation of the Spalart-Allmaras Turbulence Model." eq.(5)

Definition at line 293 of file negative_spalart_allmaras_rans_model.cpp.

◆ compute_coefficient_r()

template<int dim, int nstate, typename real >
real PHiLiP::Physics::ReynoldsAveragedNavierStokes_SAneg< dim, nstate, real >::compute_coefficient_r ( const real  nu_tilde,
const real  d_wall,
const real  s_tilde 
) const
private

Coefficient r for the negative SA model.

Reference: Steven R. Allmaras. (2012). "Modifications and Clarifications for the Implementation of the Spalart-Allmaras Turbulence Model." eq.(5)

Definition at line 276 of file negative_spalart_allmaras_rans_model.cpp.

◆ compute_cross_source()

template<int dim, int nstate, typename real >
std::array< real, nstate > PHiLiP::Physics::ReynoldsAveragedNavierStokes_SAneg< dim, nstate, real >::compute_cross_source ( const real  density,
const real  nu_tilde,
const real  laminar_kinematic_viscosity,
const std::array< dealii::Tensor< 1, dim, real >, nstate_navier_stokes > &  primitive_soln_gradient_rans,
const std::array< dealii::Tensor< 1, dim, real >, nstate_turbulence_model > &  primitive_solution_gradient_turbulence_model 
) const
private

Source term with cross derivatives for the negative SA model.

Reference: Steven R. Allmaras. (2012). "Modifications and Clarifications for the Implementation of the Spalart-Allmaras Turbulence Model." eq.(17) for positive nu_tilde

Definition at line 418 of file negative_spalart_allmaras_rans_model.cpp.

◆ compute_dissipation_source()

template<int dim, int nstate, typename real >
std::array< real, nstate > PHiLiP::Physics::ReynoldsAveragedNavierStokes_SAneg< dim, nstate, real >::compute_dissipation_source ( const real  coefficient_f_t2,
const real  density,
const real  nu_tilde,
const real  d_wall,
const real  s_tilde 
) const
private

Dissipation source term for the negative SA model.

Reference: Steven R. Allmaras. (2012). "Modifications and Clarifications for the Implementation of the Spalart-Allmaras Turbulence Model." eq.(3) for positive nu_tilde eq.(22) for negative nu_tilde

Definition at line 389 of file negative_spalart_allmaras_rans_model.cpp.

◆ compute_eddy_viscosity_templated()

template<int dim, int nstate, typename real >
template<typename real2 >
real2 PHiLiP::Physics::ReynoldsAveragedNavierStokes_SAneg< dim, nstate, real >::compute_eddy_viscosity_templated ( const std::array< real2, nstate_navier_stokes > &  primitive_soln_rans,
const std::array< real2, nstate_turbulence_model > &  primitive_soln_turbulence_model 
) const
private

Templated nondimensionalized eddy viscosity for the negative SA model.

Reference: Steven R. Allmaras. (2012). "Modifications and Clarifications for the Implementation of the Spalart-Allmaras Turbulence Model." eq.(1)

Definition at line 74 of file negative_spalart_allmaras_rans_model.cpp.

◆ compute_effective_viscosity_turbulence_model_templated()

template<int dim, int nstate, typename real >
template<typename real2 >
std::array< real2, nstate-(dim+2)> PHiLiP::Physics::ReynoldsAveragedNavierStokes_SAneg< dim, nstate, real >::compute_effective_viscosity_turbulence_model_templated ( const std::array< real2, nstate_navier_stokes > &  primitive_soln_rans,
const std::array< real2, nstate_turbulence_model > &  primitive_soln_turbulence_model 
) const
private

Templated nondimensionalized effective (total) viscosities for the negative SA model.

Reference: Steven R. Allmaras. (2012). "Modifications and Clarifications for the Implementation of the Spalart-Allmaras Turbulence Model." eq.(14)

Definition at line 137 of file negative_spalart_allmaras_rans_model.cpp.

◆ compute_production_source()

template<int dim, int nstate, typename real >
std::array< real, nstate > PHiLiP::Physics::ReynoldsAveragedNavierStokes_SAneg< dim, nstate, real >::compute_production_source ( const real  coefficient_f_t2,
const real  density,
const real  nu_tilde,
const real  s,
const real  s_tilde 
) const
private

Production source term for the negative SA model.

Reference: Steven R. Allmaras. (2012). "Modifications and Clarifications for the Implementation of the Spalart-Allmaras Turbulence Model." eq.(3) for positive nu_tilde eq.(22) for negative nu_tilde

Definition at line 364 of file negative_spalart_allmaras_rans_model.cpp.

◆ compute_s_bar()

template<int dim, int nstate, typename real >
real PHiLiP::Physics::ReynoldsAveragedNavierStokes_SAneg< dim, nstate, real >::compute_s_bar ( const real  coefficient_Chi,
const real  nu_tilde,
const real  d_wall 
) const
private

Correction of vorticity magnitude for the negative SA model.

Reference: Steven R. Allmaras. (2012). "Modifications and Clarifications for the Implementation of the Spalart-Allmaras Turbulence Model." eq.(11)

Definition at line 324 of file negative_spalart_allmaras_rans_model.cpp.

◆ compute_s_tilde()

template<int dim, int nstate, typename real >
real PHiLiP::Physics::ReynoldsAveragedNavierStokes_SAneg< dim, nstate, real >::compute_s_tilde ( const real  coefficient_Chi,
const real  nu_tilde,
const real  d_wall,
const real  s 
) const
private

Modified vorticity magnitude for the negative SA model.

Reference: Steven R. Allmaras. (2012). "Modifications and Clarifications for the Implementation of the Spalart-Allmaras Turbulence Model." eq.(12)

Definition at line 340 of file negative_spalart_allmaras_rans_model.cpp.

Member Data Documentation

◆ c_b1

template<int dim, int nstate, typename real>
const real PHiLiP::Physics::ReynoldsAveragedNavierStokes_SAneg< dim, nstate, real >::c_b1 = 0.1355
private

Constant coefficients for the negative SA model.

Reference: Steven R. Allmaras. (2012). "Modifications and Clarifications for the Implementation of the Spalart-Allmaras Turbulence Model."

Definition at line 298 of file negative_spalart_allmaras_rans_model.h.


The documentation for this class was generated from the following files: