1 #ifndef VERTEXCFD_CLOSURE_INCOMPRESSIBLELSVOFCONVECTIVEFLUX_IMPL_HPP
2 #define VERTEXCFD_CLOSURE_INCOMPRESSIBLELSVOFCONVECTIVEFLUX_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 IncompressibleLSVOFConvectiveFlux<EvalType, Traits, NumSpaceDim>::
17 IncompressibleLSVOFConvectiveFlux(
const panzer::IntegrationRule& ir,
18 const std::string& flux_prefix,
19 const std::string& field_prefix)
20 : _continuity_flux(flux_prefix +
"CONVECTIVE_FLUX_continuity", ir.dl_vector)
21 , _rho(field_prefix +
"density", ir.dl_scalar)
22 , _pressure(field_prefix +
"lagrange_pressure", ir.dl_scalar)
25 this->addEvaluatedField(_continuity_flux);
27 Utils::addEvaluatedVectorField(*
this, ir.dl_vector, _momentum_flux,
28 flux_prefix +
"CONVECTIVE_FLUX_"
32 this->addDependentField(_rho);
33 this->addDependentField(_pressure);
34 Utils::addDependentVectorField(
35 *
this, ir.dl_scalar, _velocity, field_prefix +
"velocity_");
37 this->setName(
"Incompressible LSVOF Convective Flux "
38 + std::to_string(num_space_dim) +
"D");
42 template<
class EvalType,
class Traits,
int NumSpaceDim>
43 void IncompressibleLSVOFConvectiveFlux<EvalType, Traits, NumSpaceDim>::evaluateFields(
44 typename Traits::EvalData workset)
46 auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
48 Kokkos::parallel_for(this->getName(), policy, *
this);
52 template<
class EvalType,
class Traits,
int NumSpaceDim>
53 KOKKOS_INLINE_FUNCTION
void
54 IncompressibleLSVOFConvectiveFlux<EvalType, Traits, NumSpaceDim>::operator()(
55 const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team)
const
57 const int cell = team.league_rank();
58 const int num_point = _continuity_flux.extent(1);
61 Kokkos::TeamThreadRange(team, 0, num_point), [&](
const int point) {
62 for (
int dim = 0; dim < num_space_dim; ++dim)
64 _continuity_flux(cell, point, dim)
65 = _rho(cell, point) * _velocity[dim](cell, point);
67 for (
int mom_dim = 0; mom_dim < num_space_dim; ++mom_dim)
69 _momentum_flux[mom_dim](cell, point, dim)
70 = _rho(cell, point) * _velocity[dim](cell, point)
71 * _velocity[mom_dim](cell, point);
74 _momentum_flux[mom_dim](cell, point, dim)
75 += _pressure(cell, point);
87 #endif // end VERTEXCFD_CLOSURE_INCOMPRESSIBLELSVOFCONVECTIVEFLUX_IMPL_HPP