VertexCFD  0.0-dev
VertexCFD_Closure_IncompressibleLSVOFVortexVelocity_impl.hpp
1 #ifndef VERTEXCFD_CLOSURE_INCOMPRESSIBLELSVOFVORTEXVELOCITY_IMPL_HPP
2 #define VERTEXCFD_CLOSURE_INCOMPRESSIBLELSVOFVORTEXVELOCITY_IMPL_HPP
3 
4 #include "utils/VertexCFD_Utils_Constants.hpp"
5 
6 #include <Panzer_HierarchicParallelism.hpp>
7 #include <Panzer_Workset_Utilities.hpp>
8 
9 #include <Teuchos_StandardParameterEntryValidators.hpp>
10 
11 #include <string>
12 
13 namespace VertexCFD
14 {
15 namespace ClosureModel
16 {
17 //---------------------------------------------------------------------------//
18 template<class EvalType, class Traits, int NumSpaceDim>
19 IncompressibleLSVOFVortexVelocity<EvalType, Traits, NumSpaceDim>::
20  IncompressibleLSVOFVortexVelocity(
21  const panzer::IntegrationRule& ir,
22  const Teuchos::ParameterList& closure_params)
23  : _vel_0("velocity_0", ir.dl_scalar)
24  , _vel_1("velocity_1", ir.dl_scalar)
25  , _direction(1.0)
26  , _ir_degree(ir.cubature_degree)
27 {
28  // Check that simulation is 2D
29  if (num_space_dim != 2)
30  {
31  const std::string msg
32  = "ERROR: Incompressible LSVOF Vortex Velocity must be used for "
33  "2D problem.\n";
34 
35  throw std::runtime_error(msg);
36  }
37 
38  // Check vortex direction
39  if (closure_params.isType<std::string>("Direction"))
40  {
41  const auto type_validator = Teuchos::rcp(
42  new Teuchos::StringToIntegralParameterEntryValidator<Direction>(
43  Teuchos::tuple<std::string>("Forward", "Reverse"), "Forward"));
44 
45  const Direction direction = type_validator->getIntegralValue(
46  closure_params.get<std::string>("Direction"));
47 
48  if (direction == Direction::Reverse)
49  {
50  _direction = -1.0;
51  }
52  }
53 
54  // Evaluated fields
55  this->addEvaluatedField(_vel_0);
56  this->addEvaluatedField(_vel_1);
57 
58  this->setName("Incompressible LSVOF Vortex Velocity "
59  + std::to_string(num_space_dim) + "D");
60 }
61 
62 //---------------------------------------------------------------------------//
63 template<class EvalType, class Traits, int NumSpaceDim>
64 void IncompressibleLSVOFVortexVelocity<EvalType, Traits, NumSpaceDim>::
65  postRegistrationSetup(typename Traits::SetupData sd,
66  PHX::FieldManager<Traits>&)
67 {
68  _ir_index = panzer::getIntegrationRuleIndex(_ir_degree, (*sd.worksets_)[0]);
69 }
70 
71 //---------------------------------------------------------------------------//
72 template<class EvalType, class Traits, int NumSpaceDim>
73 void IncompressibleLSVOFVortexVelocity<EvalType, Traits, NumSpaceDim>::evaluateFields(
74  typename Traits::EvalData workset)
75 {
76  _ip_coords = workset.int_rules[_ir_index]->ip_coordinates;
77 
78  auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
79  workset.num_cells);
80  Kokkos::parallel_for(this->getName(), policy, *this);
81 }
82 
83 //---------------------------------------------------------------------------//
84 template<class EvalType, class Traits, int NumSpaceDim>
85 KOKKOS_INLINE_FUNCTION void
86 IncompressibleLSVOFVortexVelocity<EvalType, Traits, NumSpaceDim>::operator()(
87  const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team) const
88 {
89  const int cell = team.league_rank();
90  const int num_point = _vel_0.extent(1);
91 
92  using Kokkos::cos;
93  using Kokkos::sin;
94 
95  using VertexCFD::Constants::pi;
96 
97  Kokkos::parallel_for(
98  Kokkos::TeamThreadRange(team, 0, num_point), [&](const int point) {
99  const double x = _ip_coords(cell, point, 0);
100  const double y = _ip_coords(cell, point, 1);
101 
102  _vel_0(cell, point) = _direction * sin(2.0 * pi * y)
103  * pow(sin(pi * x), 2.0);
104  _vel_1(cell, point) = -_direction * sin(2.0 * pi * x)
105  * pow(sin(pi * y), 2.0);
106  });
107 }
108 
109 //---------------------------------------------------------------------------//
110 
111 } // end namespace ClosureModel
112 } // end namespace VertexCFD
113 
114 #endif // end
115  // VERTEXCFD_CLOSURE_INCOMPRESSIBLELSVOFVORTEXVELOCITY_IMPL_HPP
VertexCFD
Definition: tstMethodManufacturedSolutionBC.cpp:23