VertexCFD  0.0-dev
VertexCFD_Closure_IncompressibleCLSEpsilon_impl.hpp
1 #ifndef VERTEXCFD_CLOSURE_INCOMPRESSIBLECLSEPSILON_IMPL_HPP
2 #define VERTEXCFD_CLOSURE_INCOMPRESSIBLECLSEPSILON_IMPL_HPP
3 
4 #include "utils/VertexCFD_Utils_SmoothMath.hpp"
5 
6 #include <Panzer_HierarchicParallelism.hpp>
7 
8 #include <Teuchos_StandardParameterEntryValidators.hpp>
9 
10 #include <string>
11 
12 namespace VertexCFD
13 {
14 namespace ClosureModel
15 {
16 //---------------------------------------------------------------------------//
17 template<class EvalType, class Traits, int NumSpaceDim>
18 IncompressibleCLSEpsilon<EvalType, Traits, NumSpaceDim>::IncompressibleCLSEpsilon(
19  const panzer::IntegrationRule& ir,
20  const Teuchos::ParameterList& lsvof_params)
21  : _epsilon("CLS_epsilon", ir.dl_scalar)
22  , _element_length("element_length", ir.dl_vector)
23  , _alpha(1.5)
24  , _method(EpsilonMeasurementMethod::Minimum)
25 {
26  // Check for non-default alpha coefficient
27  if (lsvof_params.isType<double>("Interface Thickness Coefficient"))
28  {
29  _alpha = lsvof_params.get<double>("Interface Thickness Coefficient");
30  }
31 
32  // Check for measurement type (default is minimum)
33  if (lsvof_params.isType<std::string>("Epsilon Measurement Method"))
34  {
35  const auto type_validator = Teuchos::rcp(
36  new Teuchos::StringToIntegralParameterEntryValidator<
37  EpsilonMeasurementMethod>(
38  Teuchos::tuple<std::string>(
39  "Minimum", "Maximum", "Magnitude", "Average"),
40  "Minimum"));
41 
42  _method = type_validator->getIntegralValue(
43  lsvof_params.get<std::string>("Epsilon Measurement Method"));
44  }
45 
46  // Evaluated fields
47  this->addEvaluatedField(_epsilon);
48 
49  // Dependent fields
50  this->addDependentField(_element_length);
51 
52  this->setName("Incompressible CLS Epsilon " + std::to_string(num_space_dim)
53  + "D");
54 }
55 
56 //---------------------------------------------------------------------------//
57 template<class EvalType, class Traits, int NumSpaceDim>
58 void IncompressibleCLSEpsilon<EvalType, Traits, NumSpaceDim>::evaluateFields(
59  typename Traits::EvalData workset)
60 {
61  auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
62  workset.num_cells);
63  Kokkos::parallel_for(this->getName(), policy, *this);
64 }
65 
66 //---------------------------------------------------------------------------//
67 template<class EvalType, class Traits, int NumSpaceDim>
68 KOKKOS_INLINE_FUNCTION void
69 IncompressibleCLSEpsilon<EvalType, Traits, NumSpaceDim>::operator()(
70  const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team) const
71 {
72  const int cell = team.league_rank();
73  const int num_point = _epsilon.extent(1);
74 
75  const double tol = 1e-10;
76 
77  using Kokkos::max;
78  using Kokkos::min;
79  using Kokkos::pow;
80 
81  Kokkos::parallel_for(
82  Kokkos::TeamThreadRange(team, 0, num_point), [&](const int point) {
83  // Minimum method
84  if (_method == EpsilonMeasurementMethod::Minimum)
85  {
86  _epsilon(cell, point) = 1e+10;
87 
88  for (int dim = 0; dim < num_space_dim; ++dim)
89  {
90  _epsilon(cell, point)
91  = min(_epsilon(cell, point),
92  _element_length(cell, point, dim));
93  }
94  }
95 
96  // Maximum method
97  else if (_method == EpsilonMeasurementMethod::Maximum)
98  {
99  _epsilon(cell, point) = tol;
100 
101  for (int dim = 0; dim < num_space_dim; ++dim)
102  {
103  _epsilon(cell, point)
104  = max(_epsilon(cell, point),
105  _element_length(cell, point, dim));
106  }
107  }
108 
109  // Magnitude method
110  else if (_method == EpsilonMeasurementMethod::Magnitude)
111  {
112  _epsilon(cell, point) = 0.0;
113 
114  for (int dim = 0; dim < num_space_dim; ++dim)
115  {
116  _epsilon(cell, point)
117  += pow(_element_length(cell, point, dim), 2.0);
118  }
119 
120  _epsilon(cell, point)
121  = pow(max(_epsilon(cell, point), tol), 0.5);
122  }
123 
124  // Average method
125  else
126  {
127  _epsilon(cell, point) = 0.0;
128 
129  for (int dim = 0; dim < num_space_dim; ++dim)
130  {
131  _epsilon(cell, point) += _element_length(cell, point, dim);
132  }
133 
134  _epsilon(cell, point) /= num_space_dim;
135  }
136 
137  _epsilon(cell, point) *= _alpha;
138  });
139 }
140 
141 //---------------------------------------------------------------------------//
142 
143 } // end namespace ClosureModel
144 } // end namespace VertexCFD
145 
146 #endif // end
147  // VERTEXCFD_CLOSURE_INCOMPRESSIBLECLSEPSILON_IMPL_HPP
VertexCFD
Definition: tstMethodManufacturedSolutionBC.cpp:23