1 #ifndef VERTEXCFD_CLOSURE_INCOMPRESSIBLECLSEPSILON_IMPL_HPP
2 #define VERTEXCFD_CLOSURE_INCOMPRESSIBLECLSEPSILON_IMPL_HPP
4 #include "utils/VertexCFD_Utils_SmoothMath.hpp"
6 #include <Panzer_HierarchicParallelism.hpp>
8 #include <Teuchos_StandardParameterEntryValidators.hpp>
14 namespace ClosureModel
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)
24 , _method(EpsilonMeasurementMethod::Minimum)
27 if (lsvof_params.isType<
double>(
"Interface Thickness Coefficient"))
29 _alpha = lsvof_params.get<
double>(
"Interface Thickness Coefficient");
33 if (lsvof_params.isType<std::string>(
"Epsilon Measurement Method"))
35 const auto type_validator = Teuchos::rcp(
36 new Teuchos::StringToIntegralParameterEntryValidator<
37 EpsilonMeasurementMethod>(
38 Teuchos::tuple<std::string>(
39 "Minimum",
"Maximum",
"Magnitude",
"Average"),
42 _method = type_validator->getIntegralValue(
43 lsvof_params.get<std::string>(
"Epsilon Measurement Method"));
47 this->addEvaluatedField(_epsilon);
50 this->addDependentField(_element_length);
52 this->setName(
"Incompressible CLS Epsilon " + std::to_string(num_space_dim)
57 template<
class EvalType,
class Traits,
int NumSpaceDim>
58 void IncompressibleCLSEpsilon<EvalType, Traits, NumSpaceDim>::evaluateFields(
59 typename Traits::EvalData workset)
61 auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
63 Kokkos::parallel_for(this->getName(), policy, *
this);
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
72 const int cell = team.league_rank();
73 const int num_point = _epsilon.extent(1);
75 const double tol = 1e-10;
82 Kokkos::TeamThreadRange(team, 0, num_point), [&](
const int point) {
84 if (_method == EpsilonMeasurementMethod::Minimum)
86 _epsilon(cell, point) = 1e+10;
88 for (
int dim = 0; dim < num_space_dim; ++dim)
91 = min(_epsilon(cell, point),
92 _element_length(cell, point, dim));
97 else if (_method == EpsilonMeasurementMethod::Maximum)
99 _epsilon(cell, point) = tol;
101 for (
int dim = 0; dim < num_space_dim; ++dim)
103 _epsilon(cell, point)
104 = max(_epsilon(cell, point),
105 _element_length(cell, point, dim));
110 else if (_method == EpsilonMeasurementMethod::Magnitude)
112 _epsilon(cell, point) = 0.0;
114 for (
int dim = 0; dim < num_space_dim; ++dim)
116 _epsilon(cell, point)
117 += pow(_element_length(cell, point, dim), 2.0);
120 _epsilon(cell, point)
121 = pow(max(_epsilon(cell, point), tol), 0.5);
127 _epsilon(cell, point) = 0.0;
129 for (
int dim = 0; dim < num_space_dim; ++dim)
131 _epsilon(cell, point) += _element_length(cell, point, dim);
134 _epsilon(cell, point) /= num_space_dim;
137 _epsilon(cell, point) *= _alpha;