1 #ifndef VERTEXCFD_CLOSURE_INCOMPRESSIBLECLSREGULARIZATION_IMPL_HPP
2 #define VERTEXCFD_CLOSURE_INCOMPRESSIBLECLSREGULARIZATION_IMPL_HPP
4 #include <Panzer_HierarchicParallelism.hpp>
10 namespace ClosureModel
13 template<
class EvalType,
class Traits,
int NumSpaceDim>
14 IncompressibleCLSRegularization<EvalType, Traits, NumSpaceDim>::
15 IncompressibleCLSRegularization(
const panzer::IntegrationRule& ir)
16 : _cls_reg(
"CLS_regularization", ir.dl_vector)
17 , _lambda(
"CLS_lambda", ir.dl_scalar)
18 , _grad_phi(
"GRAD_level_set", ir.dl_vector)
19 , _q(
"CLS_interface_normal", ir.dl_vector)
22 this->addEvaluatedField(_cls_reg);
25 this->addDependentField(_lambda);
26 this->addDependentField(_grad_phi);
27 this->addDependentField(_q);
29 this->setName(
"Incompressible CLS Regularization "
30 + std::to_string(num_space_dim) +
"D");
34 template<
class EvalType,
class Traits,
int NumSpaceDim>
35 void IncompressibleCLSRegularization<EvalType, Traits, NumSpaceDim>::evaluateFields(
36 typename Traits::EvalData workset)
38 auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
40 Kokkos::parallel_for(this->getName(), policy, *
this);
44 template<
class EvalType,
class Traits,
int NumSpaceDim>
45 KOKKOS_INLINE_FUNCTION
void
46 IncompressibleCLSRegularization<EvalType, Traits, NumSpaceDim>::operator()(
47 const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team)
const
49 const int cell = team.league_rank();
50 const int num_point = _cls_reg.extent(1);
53 Kokkos::TeamThreadRange(team, 0, num_point), [&](
const int point) {
54 for (
int dim = 0; dim < num_space_dim; ++dim)
56 _cls_reg(cell, point, dim)
57 = -1.0 * _lambda(cell, point)
58 * (_grad_phi(cell, point, dim) - _q(cell, point, dim));