VertexCFD  0.0-dev
VertexCFD_BoundaryState_TurbulenceSymmetry_impl.hpp
1 #ifndef VERTEXCFD_BOUNDARYSTATE_TURBULENCESYMMETRY_IMPL_HPP
2 #define VERTEXCFD_BOUNDARYSTATE_TURBULENCESYMMETRY_IMPL_HPP
3 
4 #include <Panzer_HierarchicParallelism.hpp>
5 
6 namespace VertexCFD
7 {
8 namespace BoundaryCondition
9 {
10 //---------------------------------------------------------------------------//
11 // Symmetry boundary condition for turbulence quantities
12 //---------------------------------------------------------------------------//
13 template<class EvalType, class Traits>
14 TurbulenceSymmetry<EvalType, Traits>::TurbulenceSymmetry(
15  const panzer::IntegrationRule& ir, const std::string variable_name)
16  : _boundary_variable("BOUNDARY_" + variable_name, ir.dl_scalar)
17  , _boundary_grad_variable("BOUNDARY_GRAD_" + variable_name, ir.dl_vector)
18  , _variable(variable_name, ir.dl_scalar)
19  , _grad_variable("GRAD_" + variable_name, ir.dl_vector)
20  , _normals("Side Normal", ir.dl_vector)
21  , _num_grad_dim(ir.spatial_dimension)
22 {
23  // Add evaluated fields
24  this->addEvaluatedField(_boundary_variable);
25  this->addEvaluatedField(_boundary_grad_variable);
26 
27  // Add dependent fields
28  this->addDependentField(_variable);
29  this->addDependentField(_grad_variable);
30  this->addDependentField(_normals);
31 
32  this->setName(variable_name + " Boundary State Turbulence Model Symmetry "
33  + std::to_string(_num_grad_dim) + "D");
34 }
35 
36 //---------------------------------------------------------------------------//
37 template<class EvalType, class Traits>
38 void TurbulenceSymmetry<EvalType, Traits>::evaluateFields(
39  typename Traits::EvalData workset)
40 {
41  auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
42  workset.num_cells);
43  Kokkos::parallel_for(this->getName(), policy, *this);
44 }
45 
46 //---------------------------------------------------------------------------//
47 template<class EvalType, class Traits>
48 KOKKOS_INLINE_FUNCTION void TurbulenceSymmetry<EvalType, Traits>::operator()(
49  const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team) const
50 {
51  const int cell = team.league_rank();
52  const int num_point = _grad_variable.extent(1);
53 
54  Kokkos::parallel_for(
55  Kokkos::TeamThreadRange(team, 0, num_point), [&](const int point) {
56  // Assign boundary values
57  _boundary_variable(cell, point) = _variable(cell, point);
58 
59  // Compute boundary normal gradient
60  scalar_type grad_var_dot_n = 0.0;
61 
62  for (int d = 0; d < _num_grad_dim; ++d)
63  {
64  grad_var_dot_n += _grad_variable(cell, point, d)
65  * _normals(cell, point, d);
66  }
67 
68  // Assign gradient
69  for (int d = 0; d < _num_grad_dim; ++d)
70  {
71  _boundary_grad_variable(cell, point, d)
72  = _grad_variable(cell, point, d)
73  - grad_var_dot_n * _normals(cell, point, d);
74  }
75  });
76 }
77 
78 //---------------------------------------------------------------------------//
79 
80 } // end namespace BoundaryCondition
81 } // end namespace VertexCFD
82 
83 #endif // end VERTEXCFD_BOUNDARYSTATE_TURBULENCESYMMETRY_IMPL_HPP
VertexCFD
Definition: tstMethodManufacturedSolutionBC.cpp:23