| [P]arallel [Hi]gh-order [Li]brary for [P]DEs
    Latest
    Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods | 
Vreman eddy viscosity model. Derived from LargeEddySimulation_Smagorinsky for only modifying compute_eddy_viscosity. More...
#include <large_eddy_simulation.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::LargeEddySimulation_Smagorinsky< dim, nstate, real > | |
| using | thermal_boundary_condition_enum = Parameters::NavierStokesParam::ThermalBoundaryCondition | 
| using | two_point_num_flux_enum = Parameters::AllParameters::TwoPointNumericalFlux | 
|  Public Types inherited from PHiLiP::Physics::LargeEddySimulationBase< dim, nstate, real > | |
| using | thermal_boundary_condition_enum = Parameters::NavierStokesParam::ThermalBoundaryCondition | 
| using | two_point_num_flux_enum = Parameters::AllParameters::TwoPointNumericalFlux | 
| Public Member Functions | |
| LargeEddySimulation_Vreman (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, const double turbulent_prandtl_number, const double ratio_of_filter_width_to_cell_size, const double model_constant, 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) | |
| real | compute_eddy_viscosity (const std::array< real, nstate > &primitive_soln, const std::array< dealii::Tensor< 1, dim, real >, nstate > &primitive_soln_gradient, const dealii::types::global_dof_index cell_index) const override | 
| FadType | compute_eddy_viscosity_fad (const std::array< FadType, nstate > &primitive_soln, const std::array< dealii::Tensor< 1, dim, FadType >, nstate > &primitive_soln_gradient, const dealii::types::global_dof_index cell_index) const override | 
|  Public Member Functions inherited from PHiLiP::Physics::LargeEddySimulation_Smagorinsky< dim, nstate, real > | |
| LargeEddySimulation_Smagorinsky (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, const double turbulent_prandtl_number, const double ratio_of_filter_width_to_cell_size, const double model_constant, 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) | |
| double | get_model_constant_times_filter_width (const dealii::types::global_dof_index cell_index) const | 
| Returns the product of the eddy viscosity model constant and the filter width. | |
| dealii::Tensor< 2, dim, real > | compute_SGS_stress_tensor (const std::array< real, nstate > &primitive_soln, const std::array< dealii::Tensor< 1, dim, real >, nstate > &primitive_soln_gradient, const dealii::types::global_dof_index cell_index) const | 
| Nondimensionalized sub-grid scale (SGS) stress tensor, (tau^sgs)*. | |
| dealii::Tensor< 1, dim, real > | compute_SGS_heat_flux (const std::array< real, nstate > &primitive_soln, const std::array< dealii::Tensor< 1, dim, real >, nstate > &primitive_soln_gradient, const dealii::types::global_dof_index cell_index) const | 
| Nondimensionalized sub-grid scale (SGS) heat flux, (q^sgs)*. | |
| dealii::Tensor< 2, dim, FadType > | compute_SGS_stress_tensor_fad (const std::array< FadType, nstate > &primitive_soln, const std::array< dealii::Tensor< 1, dim, FadType >, nstate > &primitive_soln_gradient, const dealii::types::global_dof_index cell_index) const | 
| Nondimensionalized sub-grid scale (SGS) stress tensor, (tau^sgs)* (Automatic Differentiation Type: FadType) | |
| dealii::Tensor< 1, dim, FadType > | compute_SGS_heat_flux_fad (const std::array< FadType, nstate > &primitive_soln, const std::array< dealii::Tensor< 1, dim, FadType >, nstate > &primitive_soln_gradient, const dealii::types::global_dof_index cell_index) const | 
| Nondimensionalized sub-grid scale (SGS) heat flux, (q^sgs)* (Automatic Differentiation Type: FadType) | |
|  Public Member Functions inherited from PHiLiP::Physics::LargeEddySimulationBase< dim, nstate, real > | |
| LargeEddySimulationBase (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, const double turbulent_prandtl_number, const double ratio_of_filter_width_to_cell_size, 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 | 
| Convective flux: \( \mathbf{F}_{conv} \). | |
| 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 (i.e. viscous) flux: \( \mathbf{F}_{diss} \). | |
| std::array< real, nstate > | convective_eigenvalues (const std::array< real, nstate > &, const dealii::Tensor< 1, dim, real > &) const override | 
| 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 > &solution, const real current_time, const dealii::types::global_dof_index cell_index) const | 
| Source term for manufactured solution functions. | |
| double | get_filter_width (const dealii::types::global_dof_index cell_index) const | 
| Compute the nondimensionalized filter width used by the SGS model given a cell index.  More... | |
|  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. | |
| 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 terms additional to the baseline physics (including physical source terms in additional PDEs of model) | |
| 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. | |
| virtual 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 | 
| Returns current vector solution to be used by PhysicsPostprocessor to output current solution.  More... | |
| virtual std::vector< dealii::DataComponentInterpretation::DataComponentInterpretation > | post_get_data_component_interpretation () const | 
| Returns DataComponentInterpretation of the solution to be used by PhysicsPostprocessor to output current solution.  More... | |
| virtual std::vector< std::string > | post_get_names () const | 
| Returns names of the solution to be used by PhysicsPostprocessor to output current solution.  More... | |
| Private Member Functions | |
| template<typename real2 > | |
| real2 | compute_eddy_viscosity_templated (const std::array< real2, nstate > &primitive_soln, const std::array< dealii::Tensor< 1, dim, real2 >, nstate > &primitive_soln_gradient, const dealii::types::global_dof_index cell_index) const | 
| Additional Inherited Members | |
|  Public Attributes inherited from PHiLiP::Physics::LargeEddySimulation_Smagorinsky< dim, nstate, real > | |
| const double | model_constant | 
| SGS model constant. | |
|  Public Attributes inherited from PHiLiP::Physics::LargeEddySimulationBase< dim, nstate, real > | |
| const double | turbulent_prandtl_number | 
| Turbulent Prandtl number. | |
| const double | ratio_of_filter_width_to_cell_size | 
| Ratio of filter width to cell size. | |
| std::unique_ptr< NavierStokes< dim, nstate, 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 Member Functions inherited from PHiLiP::Physics::LargeEddySimulation_Smagorinsky< dim, nstate, real > | |
| template<typename real2 > | |
| dealii::Tensor< 2, dim, real2 > | compute_SGS_stress_tensor_templated (const std::array< real2, nstate > &primitive_soln, const std::array< dealii::Tensor< 1, dim, real2 >, nstate > &primitive_soln_gradient, const dealii::types::global_dof_index cell_index) const | 
| Templated nondimensionalized sub-grid scale (SGS) stress tensor, (tau^sgs)*. | |
| template<typename real2 > | |
| dealii::Tensor< 1, dim, real2 > | compute_SGS_heat_flux_templated (const std::array< real2, nstate > &primitive_soln, const std::array< dealii::Tensor< 1, dim, real2 >, nstate > &primitive_soln_gradient, const dealii::types::global_dof_index cell_index) const | 
| Templated nondimensionalized sub-grid scale (SGS) heat flux, (q^sgs)*. | |
| template<typename real2 > | |
| real2 | scale_eddy_viscosity_templated (const std::array< real2, nstate > &primitive_soln, const real2 eddy_viscosity) const | 
| Templated scale nondimensionalized eddy viscosity for Smagorinsky model. | |
|  Protected Member Functions inherited from PHiLiP::Physics::LargeEddySimulationBase< dim, nstate, real > | |
| 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 dissipative (i.e. viscous) flux: \( \mathbf{F}_{diss} \). | |
| 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 (repeated from Euler) | |
| std::array< dealii::Tensor< 1, dim, real >, nstate > | get_manufactured_solution_gradient (const dealii::Point< dim, real > &pos) const | 
| Get manufactured solution value (repeated from Euler) | |
| std::array< real, nstate > | dissipative_source_term (const dealii::Point< dim, real > &pos, const dealii::types::global_dof_index cell_index) const | 
| Dissipative flux contribution to the source term (repeated from NavierStokes)  More... | |
|  Protected Member Functions inherited from PHiLiP::Physics::ModelBase< dim, nstate, real > | |
| virtual 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 | 
| Evaluate the manufactured solution boundary conditions. | |
| virtual void | boundary_wall (std::array< real, nstate > &soln_bc, std::array< dealii::Tensor< 1, dim, real >, nstate > &soln_grad_bc) const | 
| Wall boundary condition. | |
| virtual 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 | 
| Outflow Boundary Condition. | |
| virtual 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 | 
| Inflow boundary conditions. | |
| virtual void | boundary_farfield (std::array< real, nstate > &soln_bc) const | 
| Farfield boundary conditions based on freestream values. | |
| 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 | 
|  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. | |
Vreman eddy viscosity model. Derived from LargeEddySimulation_Smagorinsky for only modifying compute_eddy_viscosity.
Definition at line 312 of file large_eddy_simulation.h.
| PHiLiP::Physics::LargeEddySimulation_Vreman< dim, nstate, real >::LargeEddySimulation_Vreman | ( | 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, | ||
| const double | turbulent_prandtl_number, | ||
| const double | ratio_of_filter_width_to_cell_size, | ||
| const double | model_constant, | ||
| 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 sub-grid scale (SGS) model: Vreman Reference: Vreman, A. W. (2004) "An eddy-viscosity subgrid-scale model for turbulent shear flow: Algebraic theory and applications."
Definition at line 767 of file large_eddy_simulation.cpp.
| 
 | overridevirtual | 
Nondimensionalized eddy viscosity for the Vreman model. Reference: Vreman, A. W. (2004) "An eddy-viscosity subgrid-scale model for turbulent shear flow: Algebraic theory and applications."
Reimplemented from PHiLiP::Physics::LargeEddySimulation_Smagorinsky< dim, nstate, real >.
Definition at line 808 of file large_eddy_simulation.cpp.
| 
 | overridevirtual | 
Nondimensionalized eddy viscosity for the Vreman model. (Automatic Differentiation Type: FadType) Reference: Vreman, A. W. (2004) "An eddy-viscosity subgrid-scale model for turbulent shear flow: Algebraic theory and applications."
Reimplemented from PHiLiP::Physics::LargeEddySimulation_Smagorinsky< dim, nstate, real >.
Definition at line 818 of file large_eddy_simulation.cpp.
| 
 | private | 
Templated nondimensionalized eddy viscosity for the Vreman model. Reference: Vreman, A. W. (2004) "An eddy-viscosity subgrid-scale model for turbulent shear flow: Algebraic theory and applications."
Reference: Vreman (2004) - Equation (7)
Eddy viscosity is zero in the absence of turbulent fluctuations, i.e. zero strain rate and zero rotation rate, also B is positive-semidefinite (i.e. beta_tensor_determinant>=0). See Vreman (2004). Since the denominator in this eddy viscosity model will go to zero, we must explicitly set the eddy viscosity to zero to avoid a division by zero. Or equivalently, update it from its zero initialization only if there is turbulence.
Definition at line 829 of file large_eddy_simulation.cpp.