VertexCFD  0.0-dev
VertexCFD_Closure_IncompressiblePlanarPoiseuilleExact_impl.hpp
1 #ifndef VERTEXCFD_CLOSURE_INCOMPRESSIBLEPLANARPOISEUILLEEXACT_IMPL_HPP
2 #define VERTEXCFD_CLOSURE_INCOMPRESSIBLEPLANARPOISEUILLEEXACT_IMPL_HPP
3 
4 #include <utils/VertexCFD_Utils_VectorField.hpp>
5 
6 #include "Panzer_GlobalIndexer.hpp"
7 #include "Panzer_PureBasis.hpp"
8 #include "Panzer_Workset_Utilities.hpp"
9 #include <Panzer_HierarchicParallelism.hpp>
10 
11 #include <string>
12 
13 namespace VertexCFD
14 {
15 namespace ClosureModel
16 {
17 //---------------------------------------------------------------------------//
18 template<class EvalType, class Traits, int NumSpaceDim>
19 IncompressiblePlanarPoiseuilleExact<EvalType, Traits, NumSpaceDim>::
20  IncompressiblePlanarPoiseuilleExact(
21  const panzer::IntegrationRule& ir,
22  const Teuchos::ParameterList& closure_params)
23  : _temperature("Exact_temperature", ir.dl_scalar)
24  , _lagrange_pressure("Exact_lagrange_pressure", ir.dl_scalar)
25  , _nu(closure_params.get<double>("Kinematic viscosity"))
26  , _rho(closure_params.get<double>("Density"))
27  , _k(closure_params.get<double>("Thermal conductivity"))
28  , _cp(closure_params.get<double>("Specific heat capacity"))
29  , _ir_degree(ir.cubature_degree)
30 {
31  // Read parameters from planar Poiseuille sublist
32  _h_min = closure_params.get<double>("H min");
33  _h_max = closure_params.get<double>("H max");
34  _T_l = closure_params.get<double>("Lower wall temperature");
35  _T_u = closure_params.get<double>("Upper wall temperature");
36  _S_u = closure_params.get<double>("Momentum source");
37 
38  // Evaluated fields
39  this->addEvaluatedField(_temperature);
40  this->addEvaluatedField(_lagrange_pressure);
41 
42  Utils::addEvaluatedVectorField(
43  *this, ir.dl_scalar, _velocity, "Exact_velocity_");
44 
45  this->setName("Exact Solution Planar Poiseuille");
46 }
47 
48 //---------------------------------------------------------------------------//
49 template<class EvalType, class Traits, int NumSpaceDim>
50 void IncompressiblePlanarPoiseuilleExact<EvalType, Traits, NumSpaceDim>::
51  postRegistrationSetup(typename Traits::SetupData sd,
52  PHX::FieldManager<Traits>&)
53 {
54  _ir_index = panzer::getIntegrationRuleIndex(_ir_degree, (*sd.worksets_)[0]);
55 }
56 
57 //---------------------------------------------------------------------------//
58 template<class EvalType, class Traits, int NumSpaceDim>
59 void IncompressiblePlanarPoiseuilleExact<EvalType, Traits, NumSpaceDim>::
60  evaluateFields(typename Traits::EvalData workset)
61 {
62  _ip_coords = workset.int_rules[_ir_index]->ip_coordinates;
63  auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
64  workset.num_cells);
65  Kokkos::parallel_for(this->getName(), policy, *this);
66 }
67 
68 //---------------------------------------------------------------------------//
69 template<class EvalType, class Traits, int NumSpaceDim>
70 KOKKOS_INLINE_FUNCTION void
71 IncompressiblePlanarPoiseuilleExact<EvalType, Traits, NumSpaceDim>::operator()(
72  const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team) const
73 {
74  const int cell = team.league_rank();
75  const int num_point = _temperature.extent(1);
76 
77  using Kokkos::pow;
78 
79  Kokkos::parallel_for(
80  Kokkos::TeamThreadRange(team, 0, num_point), [&](const int point) {
81  // Note: an analytical solution is only available for the velocity
82  // and temperature fields
83  _lagrange_pressure(cell, point) = 0.0;
84 
85  // Get geometric and thermal parameters
86  const double H = (_h_max - _h_min) / 2.0;
87  const double y = _ip_coords(cell, point, 1);
88  const double Pr = _cp * _rho * _nu / _k;
89  const double dT = _T_l - _T_u;
90  const double U = _S_u * H * H / 3.0 / _nu;
91  const double E = U * U / _cp / dT;
92 
93  // Calculate excess temperature
94  const double T_star = 0.5 * (1.0 - (y / H))
95  + 3.0 / 4.0 * Pr * E
96  * (1.0 - pow(y / H, 4.0));
97  // Back out exact temperature
98  _temperature(cell, point) = T_star * dT + _T_u;
99 
100  // Calculate exact velocity
101  _velocity[0](cell, point) = 3.0 / 2.0 * U * (1.0 - pow(y / H, 2.0));
102  _velocity[1](cell, point) = 0.0;
103 
104  if (num_space_dim == 3)
105  _velocity[2](cell, point) = 0.0;
106  });
107 }
108 
109 //---------------------------------------------------------------------------//
110 
111 } // end namespace ClosureModel
112 } // end namespace VertexCFD
113 
114 #endif // end VERTEXCFD_CLOSURE_INCOMPRESSIBLEPLANARPOISEUILLEEXACT_IMPL_HPP
VertexCFD
Definition: tstMethodManufacturedSolutionBC.cpp:23