1 #ifndef VERTEXCFD_CLOSURE_VARIABLECONVECTIVEFLUX_IMPL_HPP
2 #define VERTEXCFD_CLOSURE_VARIABLECONVECTIVEFLUX_IMPL_HPP
4 #include <utils/VertexCFD_Utils_VectorField.hpp>
6 #include <Panzer_HierarchicParallelism.hpp>
12 namespace ClosureModel
15 template<
class EvalType,
class Traits,
int NumSpaceDim>
16 VariableConvectiveFlux<EvalType, Traits, NumSpaceDim>::VariableConvectiveFlux(
17 const panzer::IntegrationRule& ir,
18 const Teuchos::ParameterList& closure_params,
19 const std::string& flux_prefix,
20 const std::string& field_prefix)
21 : _variable_name(closure_params.get<std::string>(
"Field Name"))
22 , _equation_name(closure_params.get<std::string>(
"Equation Name"))
23 , _var_flux(flux_prefix +
"CONVECTIVE_FLUX_" + _equation_name, ir.dl_vector)
24 , _var(field_prefix + _variable_name, ir.dl_scalar)
27 this->addEvaluatedField(_var_flux);
30 this->addDependentField(_var);
32 Utils::addDependentVectorField(
33 *
this, ir.dl_scalar, _velocity, field_prefix +
"velocity_");
35 this->setName(_equation_name +
" Convective Flux "
36 + std::to_string(num_space_dim) +
"D");
40 template<
class EvalType,
class Traits,
int NumSpaceDim>
41 void VariableConvectiveFlux<EvalType, Traits, NumSpaceDim>::evaluateFields(
42 typename Traits::EvalData workset)
44 auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
46 Kokkos::parallel_for(this->getName(), policy, *
this);
50 template<
class EvalType,
class Traits,
int NumSpaceDim>
51 KOKKOS_INLINE_FUNCTION
void
52 VariableConvectiveFlux<EvalType, Traits, NumSpaceDim>::operator()(
53 const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team)
const
55 const int cell = team.league_rank();
56 const int num_point = _var_flux.extent(1);
57 const int num_grad_dim = _var_flux.extent(2);
60 Kokkos::TeamThreadRange(team, 0, num_point), [&](
const int point) {
61 for (
int dim = 0; dim < num_grad_dim; ++dim)
63 _var_flux(cell, point, dim) = _var(cell, point)
64 * _velocity[dim](cell, point);
74 #endif // end VERTEXCFD_CLOSURE_VARIABLECONVECTIVEFLUX_IMPL_HPP