VertexCFD  0.0-dev
VertexCFD_Closure_IncompressibleLSVOFScalarTimeDerivative_impl.hpp
1 #ifndef VERTEXCFD_CLOSURE_INCOMPRESSIBLELSVOFSCALARTIMEDERIVATIVE_IMPL_HPP
2 #define VERTEXCFD_CLOSURE_INCOMPRESSIBLELSVOFSCALARTIMEDERIVATIVE_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>
16 IncompressibleLSVOFScalarTimeDerivative<EvalType, Traits>::
17  IncompressibleLSVOFScalarTimeDerivative(const panzer::IntegrationRule& ir,
18  const std::string& dof_name,
19  const std::string& eqn_name,
20  const bool& mass_weighted)
21  : _dqdt_dof("DQDT_" + eqn_name, ir.dl_scalar)
22  , _dof(dof_name, ir.dl_scalar)
23  , _rho("density", ir.dl_scalar)
24  , _dxdt_dof("DXDT_" + dof_name, ir.dl_scalar)
25  , _dxdt_rho("DXDT_density", ir.dl_scalar)
26  , _mass_weighted(mass_weighted)
27 {
28  // Evaluated fields
29  this->addEvaluatedField(_dqdt_dof);
30 
31  // Dependent fields
32  this->addDependentField(_dxdt_dof);
33 
34  if (_mass_weighted)
35  {
36  this->addDependentField(_dof);
37  this->addDependentField(_rho);
38  this->addDependentField(_dxdt_rho);
39  }
40 
41  this->setName("Incompressible LSVOF " + dof_name + " Time Derivative");
42 }
43 
44 //---------------------------------------------------------------------------//
45 template<class EvalType, class Traits>
46 void IncompressibleLSVOFScalarTimeDerivative<EvalType, Traits>::evaluateFields(
47  typename Traits::EvalData workset)
48 {
49  auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
50  workset.num_cells);
51  Kokkos::parallel_for(this->getName(), policy, *this);
52 }
53 
54 //---------------------------------------------------------------------------//
55 template<class EvalType, class Traits>
56 KOKKOS_INLINE_FUNCTION void
57 IncompressibleLSVOFScalarTimeDerivative<EvalType, Traits>::operator()(
58  const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team) const
59 {
60  const int cell = team.league_rank();
61  const int num_point = _dqdt_dof.extent(1);
62 
63  Kokkos::parallel_for(
64  Kokkos::TeamThreadRange(team, 0, num_point), [&](const int point) {
65  if (_mass_weighted)
66  {
67  _dqdt_dof(cell, point)
68  = _rho(cell, point) * _dxdt_dof(cell, point)
69  + _dxdt_rho(cell, point) * _dof(cell, point);
70  }
71  else
72  {
73  _dqdt_dof(cell, point) = _dxdt_dof(cell, point);
74  }
75  });
76 }
77 
78 //---------------------------------------------------------------------------//
79 
80 } // end namespace ClosureModel
81 } // end namespace VertexCFD
82 
83 #endif // end
84  // VERTEXCFD_CLOSURE_INCOMPRESSIBLELSVOFSCALARTIMEDERIVATIVE_IMPL_HPP
VertexCFD
Definition: tstMethodManufacturedSolutionBC.cpp:23