VertexCFD  0.0-dev
VertexCFD_BoundaryState_TurbulenceKEpsilonWallFunction.hpp
1 #ifndef VERTEXCFD_BOUNDARYSTATE_TURBULENCEKEPSILONWALLFUNCTION_HPP
2 #define VERTEXCFD_BOUNDARYSTATE_TURBULENCEKEPSILONWALLFUNCTION_HPP
3 
4 #include <Panzer_Dimension.hpp>
5 #include <Panzer_Evaluator_WithBaseImpl.hpp>
6 
7 #include <Phalanx_Evaluator_Derived.hpp>
8 #include <Phalanx_Evaluator_WithBaseImpl.hpp>
9 #include <Phalanx_FieldManager.hpp>
10 #include <Phalanx_config.hpp>
11 
12 #include <string>
13 
14 namespace VertexCFD
15 {
16 namespace BoundaryCondition
17 {
18 //---------------------------------------------------------------------------//
19 // Wall function boundary conditions for K-Epsilon family of turbulence
20 // models as outlined by Kuzmin et al. (2007)
21 //---------------------------------------------------------------------------//
22 template<class EvalType, class Traits, int NumSpaceDim>
24  : public panzer::EvaluatorWithBaseImpl<Traits>,
25  public PHX::EvaluatorDerived<EvalType, Traits>
26 {
27  public:
28  using scalar_type = typename EvalType::ScalarT;
29  static constexpr int num_space_dim = NumSpaceDim;
30 
31  TurbulenceKEpsilonWallFunction(const panzer::IntegrationRule& ir,
32  const Teuchos::ParameterList& bc_params);
33 
34  void evaluateFields(typename Traits::EvalData workset) override;
35 
36  KOKKOS_INLINE_FUNCTION
37  void operator()(
38  const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team) const;
39 
40  public:
41  PHX::MDField<scalar_type, panzer::Cell, panzer::Point> _boundary_k;
42  PHX::MDField<scalar_type, panzer::Cell, panzer::Point> _boundary_e;
43 
44  PHX::MDField<scalar_type, panzer::Cell, panzer::Point, panzer::Dim>
45  _boundary_grad_k;
46 
47  PHX::MDField<scalar_type, panzer::Cell, panzer::Point, panzer::Dim>
48  _boundary_grad_e;
49 
50  PHX::MDField<scalar_type, panzer::Cell, panzer::Point> _boundary_u_tau;
51  PHX::MDField<scalar_type, panzer::Cell, panzer::Point> _boundary_y_plus;
52  PHX::MDField<scalar_type, panzer::Cell, panzer::Point> _wall_func_nu_t;
53 
54  private:
55  PHX::MDField<const scalar_type, panzer::Cell, panzer::Point> _k;
56  PHX::MDField<const scalar_type, panzer::Cell, panzer::Point> _e;
57  Kokkos::Array<PHX::MDField<const scalar_type, panzer::Cell, panzer::Point>,
58  num_space_dim>
59  _velocity;
60  PHX::MDField<const scalar_type, panzer::Cell, panzer::Point, panzer::Dim>
61  _grad_k;
62  PHX::MDField<const scalar_type, panzer::Cell, panzer::Point, panzer::Dim>
63  _grad_e;
64  PHX::MDField<const scalar_type, panzer::Cell, panzer::Point, panzer::Dim>
65  _normals;
66  PHX::MDField<const scalar_type, panzer::Cell, panzer::Point> _nu;
67 
68  int _num_grad_dim;
69  double _C_mu;
70  double _kappa;
71  double _yp_tr;
72  bool _neumann;
73 };
74 
75 //---------------------------------------------------------------------------//
76 
77 } // end namespace BoundaryCondition
78 } // end namespace VertexCFD
79 
80 #endif // VERTEXCFD_BOUNDARYSTATE_TURBULENCEKEPSILONWALLFUNCTION_HPP
VertexCFD
Definition: tstMethodManufacturedSolutionBC.cpp:23
VertexCFD::BoundaryCondition::TurbulenceKEpsilonWallFunction
Definition: VertexCFD_BoundaryState_TurbulenceKEpsilonWallFunction.hpp:26