1 #ifndef __REYNOLDS_AVERAGED_NAVIER_STOKES__     2 #define __REYNOLDS_AVERAGED_NAVIER_STOKES__     5 #include "navier_stokes.h"    12 template <
int dim, 
int nstate, 
typename real>
    21         const double                                              ref_length,
    22         const double                                              gamma_gas,
    23         const double                                              mach_inf,
    24         const double                                              angle_of_attack,
    25         const double                                              side_slip_angle,
    26         const double                                              prandtl_number,
    27         const double                                              reynolds_number_inf,
    28         const bool                                                use_constant_viscosity,
    29         const double                                              constant_viscosity,
    31         const double                                              temperature_inf = 273.15,
    32         const double                                              isothermal_wall_temperature = 1.0,
    51         const std::array<real,nstate> &conservative_soln) 
const;
    55         const std::array<real,nstate> &conservative_soln,
    56         const std::array<dealii::Tensor<1,dim,real>,nstate> &solution_gradient,
    57         const dealii::types::global_dof_index cell_index) 
const;
    63         const std::array<real,nstate> &,
    64         const dealii::Tensor<1,dim,real> &) 
const;
    75         const std::array<real,nstate> &soln,
    76         const dealii::Tensor<1,dim,real> &normal) 
const;
    80         const dealii::Point<dim,real> &pos,
    81         const std::array<real,nstate> &conservative_solution,
    82         const real current_time,
    83         const dealii::types::global_dof_index cell_index) 
const;
    87         const dealii::Point<dim,real> &pos,
    88         const std::array<real,nstate> &conservative_solution,
    89         const dealii::types::global_dof_index cell_index) 
const;
    93         const dealii::Point<dim,real> &pos,
    94         const std::array<real,nstate> &conservative_solution,
    95         const std::array<dealii::Tensor<1,dim,real>,nstate> &solution_gradient,
    96         const dealii::types::global_dof_index cell_index) 
const override;
   100         const std::array<real,nstate_navier_stokes> &primitive_soln_rans,
   101         const std::array<dealii::Tensor<1,dim,real>,nstate_navier_stokes> &primitive_soln_gradient_rans,
   102         const std::array<real,nstate_turbulence_model> &primitive_soln_turbulence_model) 
const = 0;
   106         const std::array<real,nstate_navier_stokes> &primitive_soln_rans,
   107         const std::array<dealii::Tensor<1,dim,real>,nstate_navier_stokes> &primitive_soln_gradient_rans,
   108         const std::array<real,nstate_turbulence_model> &primitive_soln_turbulence_model) 
const = 0;
   112         const std::array<FadType,nstate_navier_stokes> &primitive_soln_rans,
   113         const std::array<dealii::Tensor<1,dim,FadType>,nstate_navier_stokes> &primitive_soln_gradient_rans,
   114         const std::array<FadType,nstate_turbulence_model> &primitive_soln_turbulence_model) 
const = 0;
   118         const std::array<FadType,nstate_navier_stokes> &primitive_soln_rans,
   119         const std::array<dealii::Tensor<1,dim,FadType>,nstate_navier_stokes> &primitive_soln_gradient_rans,
   120         const std::array<FadType,nstate_turbulence_model> &primitive_soln_turbulence_model) 
const = 0;
   124         const std::array<real,nstate_navier_stokes> &primitive_soln_rans,
   125         const std::array<real,nstate_turbulence_model> &primitive_soln_turbulence_model) 
const = 0;
   129         const std::array<FadType,nstate_navier_stokes> &primitive_soln_rans,
   130         const std::array<FadType,nstate_turbulence_model> &primitive_soln_turbulence_model) 
const = 0;
   134         const dealii::Point<dim,real> &pos,
   135         const std::array<real,nstate> &conservative_solution,
   136         const std::array<dealii::Tensor<1,dim,real>,nstate> &solution_gradient) 
const = 0;
   140     template<
typename real2> 
   144     template<
typename real2> 
   148     template<
typename real2>
   150         const std::array<real2,nstate> &conservative_soln,
   151         const std::array<dealii::Tensor<1,dim,real2>,nstate> &solution_gradient,
   152         const dealii::types::global_dof_index cell_index) 
const;
   155     template<
typename real2>
   157         const std::array<real2,nstate> &conservative_soln) 
const;
   160     template <
typename real2>
   162         const std::array<real2,nstate> &conservative_soln) 
const;
   165     template <
typename real2>
   167         const std::array<dealii::Tensor<1,dim,real2>,nstate> &solution_gradient) 
const;
   170     template <
typename real2>
   172         const std::array<real2,nstate_navier_stokes> &primitive_soln_rans,
   173         const std::array<real2,nstate_turbulence_model> &primitive_soln_turbulence_model,
   174         const std::array<dealii::Tensor<1,dim,real2>,nstate_turbulence_model> &primitive_solution_gradient_turbulence_model) 
const;
   178     template <
typename real2>
   180         const std::array<real2,nstate> &conservative_soln) 
const;
   184     template <
typename real2>
   186         const std::array<real2,nstate> &conservative_soln,
   187         const std::array<dealii::Tensor<1,dim,real2>,nstate> &solution_gradient) 
const;
   193         const std::array<real,nstate> &conservative_soln1,
   194         const std::array<real,nstate> &conservative_soln2) 
const;
   201         const std::array<real,nstate> &conservative_soln,
   202         const dealii::Tensor<1,dim,real> &normal) 
const;
   209         const std::array<real,nstate> &conservative_soln,
   210         const std::array<dealii::Tensor<1,dim,real>,nstate> &solution_gradient,
   211         const dealii::Tensor<1,dim,real> &normal,
   212         const dealii::types::global_dof_index cell_index) 
const;
   219         const std::array<real,nstate> &conservative_soln,
   220         const std::array<dealii::Tensor<1,dim,real>,nstate> &solution_gradient,
   221         const dealii::Tensor<1,dim,real> &normal,
   222         const int d_gradient,
   223         const dealii::types::global_dof_index cell_index) 
const;
   227         const dealii::Point<dim,real> &pos) 
const;
   231         const dealii::Point<dim,real> &pos) 
const;
   238         const dealii::Point<dim,real> &pos) 
const;
   245         const dealii::Point<dim,real> &pos,
   246         const dealii::types::global_dof_index cell_index) 
const;
   252         const dealii::Point<dim,real> &pos,
   253         const dealii::types::global_dof_index cell_index) 
const;
   257         const dealii::Point<dim, real> &pos,
   258         const dealii::Tensor<1,dim,real> &normal_int,
   259         const std::array<real,nstate> &soln_int,
   260         const std::array<dealii::Tensor<1,dim,real>,nstate> &soln_grad_int,
   261         std::array<real,nstate> &soln_bc,
   262         std::array<dealii::Tensor<1,dim,real>,nstate> &soln_grad_bc) 
const override;
 virtual 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 =0
Physical source term (production, dissipation source terms and source term with cross derivatives) in...
std::shared_ptr< ManufacturedSolutionFunction< dim, real > > manufactured_solution_function
Manufactured solution function. 
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. 
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. 
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< real2, nstate-(dim+2)> convert_conservative_to_primitive_turbulence_model(const std::array< real2, nstate > &conservative_soln) 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
virtual std::array< real, nstate_turbulence_model > 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 =0
Nondimensionalized effective (total) viscosities for the turbulence model. 
Manufactured solution used for grid studies to check convergence orders. 
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, 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. 
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 > get_manufactured_solution_value(const dealii::Point< dim, real > &pos) const
Get manufactured solution value. 
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
Files for the baseline physics. 
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. 
real max_convective_eigenvalue(const std::array< real, nstate > &soln) const
Maximum convective eigenvalue of the additional models' PDEs. 
virtual std::array< FadType, nstate_turbulence_model > 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 =0
Nondimensionalized effective (total) viscosities for the turbulence model (Automatic Differentiation ...
virtual dealii::Tensor< 2, dim, FadType > compute_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 =0
Nondimensionalized Reynolds stress tensor, (tau^reynolds)* (Automatic Differentiation Type: FadType) ...
Physics model additional terms and equations to the baseline physics. 
Main parameter class that contains the various other sub-parameter classes. 
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< 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 add...
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. 
virtual dealii::Tensor< 1, dim, FadType > compute_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 =0
Nondimensionalized Reynolds heat flux, (q^reynolds)* (Automatic Differentiation Type: FadType) ...
TwoPointNumericalFlux
Two point numerical flux type for split form. 
std::array< dealii::Tensor< 1, dim, real2 >, nstate > convective_flux_templated(const std::array< real2, nstate > &conservative_soln) const
Templated additional convective flux. 
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. 
virtual 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 =0
Nondimensionalized Reynolds stress tensor, (tau^reynolds)*. 
std::unique_ptr< NavierStokes< dim, nstate_navier_stokes, real > > navier_stokes_physics
Pointer to Navier-Stokes physics object. 
static const int nstate_navier_stokes
Number of PDEs for RANS equations. 
Reynolds-Averaged Navier-Stokes (RANS) equations. Derived from Navier-Stokes for modifying the stress...
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. 
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...
dealii::Tensor< 2, nstate, real > convective_flux_directional_jacobian(const std::array< real, nstate > &conservative_soln, const dealii::Tensor< 1, dim, real > &normal) const
const double turbulent_prandtl_number
Turbulent Prandtl number. 
real2 get_vector_magnitude_sqr(const dealii::Tensor< 1, 3, real2 > &vector) const
Returns the square of the magnitude of the vector. 
virtual 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 =0
Nondimensionalized Reynolds heat flux, (q^reynolds)*. 
std::array< real, nstate > convective_source_term_computed_from_manufactured_solution(const dealii::Point< dim, real > &pos) const
static const int nstate_turbulence_model
Number of PDEs for RANS turbulence model. 
ThermalBoundaryCondition
Types of thermal boundary conditions available. 
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
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 itsel...
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 R...
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 > 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. 
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.