VertexCFD  0.0-dev
VertexCFD_Closure_IncompressibleLSVOFNthPhaseFraction_impl.hpp
1 #ifndef VERTEXCFD_CLOSURE_INCOMPRESSIBLELSVOFNTHPHASEFRACTION_IMPL_HPP
2 #define VERTEXCFD_CLOSURE_INCOMPRESSIBLELSVOFNTHPHASEFRACTION_IMPL_HPP
3 
4 #include "utils/VertexCFD_Utils_PhaseLayout.hpp"
5 #include "utils/VertexCFD_Utils_SmoothMath.hpp"
6 
7 #include <Panzer_HierarchicParallelism.hpp>
8 #include <Panzer_Workset_Utilities.hpp>
9 
10 #include <Phalanx_DataLayout_DynamicLayout.hpp>
11 
12 #include <Sacado_Traits.hpp>
13 
14 #include <string>
15 
16 namespace VertexCFD
17 {
18 namespace ClosureModel
19 {
20 //---------------------------------------------------------------------------//
21 template<class EvalType, class Traits>
22 IncompressibleLSVOFNthPhaseFraction<EvalType, Traits>::
23  IncompressibleLSVOFNthPhaseFraction(
24  const panzer::IntegrationRule& ir,
25  const std::vector<std::string>& phase_names,
26  const std::string& field_prefix)
27  : _alpha_n(field_prefix + phase_names[phase_names.size() - 1], ir.dl_scalar)
28  , _phase_layout(
29  Utils::buildPhaseLayout(ir.dl_scalar, phase_names.size() - 1))
30  , _alphas(field_prefix + "volume_fractions", _phase_layout)
31 {
32  // Evaluated fields
33  this->addEvaluatedField(_alpha_n);
34 
35  // Dependent fields
36  this->addDependentField(_alphas);
37 
38  this->setName("Incompressible LSVOF Nth Phase Fraction");
39 }
40 
41 //---------------------------------------------------------------------------//
42 template<class EvalType, class Traits>
43 void IncompressibleLSVOFNthPhaseFraction<EvalType, Traits>::evaluateFields(
44  typename Traits::EvalData workset)
45 {
46  auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
47  workset.num_cells);
48  Kokkos::parallel_for(this->getName(), policy, *this);
49 }
50 
51 //---------------------------------------------------------------------------//
52 template<class EvalType, class Traits>
53 KOKKOS_INLINE_FUNCTION void
54 IncompressibleLSVOFNthPhaseFraction<EvalType, Traits>::operator()(
55  const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team) const
56 {
57  const int cell = team.league_rank();
58  const int num_point = _alpha_n.extent(1);
59 
60  Kokkos::parallel_for(
61  Kokkos::TeamThreadRange(team, 0, num_point), [&](const int point) {
62  // Set alpha_n = 1.0 - alpha_1 - ... - alpha_{n-1}
63  _alpha_n(cell, point) = 1.0;
64 
65  for (size_t phase = 0; phase < _alphas.extent(2); ++phase)
66  {
67  _alpha_n(cell, point) -= _alphas(cell, point, phase);
68  }
69 
70  // Limit alpha_n to [0, 1]
71  _alpha_n(cell, point) = SmoothMath::max(
72  SmoothMath::min(_alpha_n(cell, point), 1.0, 0.0), 0.0, 0.0);
73  });
74 }
75 
76 //---------------------------------------------------------------------------//
77 
78 } // end namespace ClosureModel
79 } // end namespace VertexCFD
80 
81 #endif // end VERTEXCFD_CLOSURE_INCOMPRESSIBLELSVOFNTHPHASEFRACTION_IMPL_HPP
VertexCFD
Definition: tstMethodManufacturedSolutionBC.cpp:23