VertexCFD  0.0-dev
VertexCFD_Closure_IncompressibleLSVOFTimeDerivative_impl.hpp
1 #ifndef VERTEXCFD_CLOSURE_INCOMPRESSIBLELSVOFTIMEDERIVATIVE_IMPL_HPP
2 #define VERTEXCFD_CLOSURE_INCOMPRESSIBLELSVOFTIMEDERIVATIVE_IMPL_HPP
3 
4 #include "utils/VertexCFD_Utils_VectorField.hpp"
5 
6 #include <Panzer_HierarchicParallelism.hpp>
7 
8 #include <string>
9 
10 namespace VertexCFD
11 {
12 namespace ClosureModel
13 {
14 //---------------------------------------------------------------------------//
15 template<class EvalType, class Traits, int NumSpaceDim>
16 IncompressibleLSVOFTimeDerivative<EvalType, Traits, NumSpaceDim>::
17  IncompressibleLSVOFTimeDerivative(const panzer::IntegrationRule& ir,
18  const Teuchos::ParameterList& lsvof_params)
19  : _dqdt_continuity("DQDT_continuity", ir.dl_scalar)
20  , _rho("density", ir.dl_scalar)
21  , _dxdt_pressure("DXDT_lagrange_pressure", ir.dl_scalar)
22  , _dxdt_density("DXDT_density", ir.dl_scalar)
23  , _beta(lsvof_params.get<double>("Mixture Artificial Compressibility"))
24 {
25  // Evaluated fields
26  this->addEvaluatedField(_dqdt_continuity);
27 
28  Utils::addEvaluatedVectorField(
29  *this, ir.dl_scalar, _dqdt_momentum, "DQDT_momentum_");
30 
31  // Dependent fields
32  this->addDependentField(_rho);
33  this->addDependentField(_dxdt_pressure);
34  this->addDependentField(_dxdt_density);
35  Utils::addDependentVectorField(*this, ir.dl_scalar, _velocity, "velocity_");
36  Utils::addDependentVectorField(
37  *this, ir.dl_scalar, _dxdt_velocity, "DXDT_velocity_");
38 
39  this->setName("Incompressible LSVOF Time Derivative "
40  + std::to_string(num_space_dim) + "D");
41 }
42 
43 //---------------------------------------------------------------------------//
44 template<class EvalType, class Traits, int NumSpaceDim>
45 void IncompressibleLSVOFTimeDerivative<EvalType, Traits, NumSpaceDim>::evaluateFields(
46  typename Traits::EvalData workset)
47 {
48  auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
49  workset.num_cells);
50  Kokkos::parallel_for(this->getName(), policy, *this);
51 }
52 
53 //---------------------------------------------------------------------------//
54 template<class EvalType, class Traits, int NumSpaceDim>
55 KOKKOS_INLINE_FUNCTION void
56 IncompressibleLSVOFTimeDerivative<EvalType, Traits, NumSpaceDim>::operator()(
57  const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team) const
58 {
59  const int cell = team.league_rank();
60  const int num_point = _dqdt_continuity.extent(1);
61 
62  Kokkos::parallel_for(
63  Kokkos::TeamThreadRange(team, 0, num_point), [&](const int point) {
64  _dqdt_continuity(cell, point) = _dxdt_pressure(cell, point) / _beta;
65 
66  for (int dim = 0; dim < num_space_dim; ++dim)
67  {
68  _dqdt_momentum[dim](cell, point)
69  = _rho(cell, point) * _dxdt_velocity[dim](cell, point)
70  + _dxdt_density(cell, point)
71  * _velocity[dim](cell, point);
72  }
73  });
74 }
75 
76 //---------------------------------------------------------------------------//
77 
78 } // end namespace ClosureModel
79 } // end namespace VertexCFD
80 
81 #endif // end VERTEXCFD_CLOSURE_INCOMPRESSIBLELSVOFTIMEDERIVATIVE_IMPL_HPP
VertexCFD
Definition: tstMethodManufacturedSolutionBC.cpp:23