VertexCFD  0.0-dev
VertexCFD_Closure_HartmannProblemExact_impl.hpp
1 #ifndef VERTEXCFD_CLOSURE_HARTMANNPROBLEMEXACT_IMPL_HPP
2 #define VERTEXCFD_CLOSURE_HARTMANNPROBLEMEXACT_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 HartmannProblemExact<EvalType, Traits, NumSpaceDim>::HartmannProblemExact(
20  const panzer::IntegrationRule& ir,
21  const Teuchos::ParameterList& closure_params,
22  const Teuchos::ParameterList& user_params)
23  : _exact_lagrange_pressure("Exact_lagrange_pressure", ir.dl_scalar)
24  , _exact_elec_pot("Exact_electric_potential", ir.dl_scalar)
25  , _nu(closure_params.get<double>("Kinematic Viscosity"))
26  , _sigma(closure_params.get<double>("Electrical Conductivity"))
27  , _rho(closure_params.get<double>("Density"))
28  , _L(closure_params.get<double>("Hartmann Solution Half-Width Channel"))
29  , _B(0.0)
30  , _ir_degree(ir.cubature_degree)
31 {
32  // Evaluated fields
33  this->addEvaluatedField(_exact_lagrange_pressure);
34  Utils::addEvaluatedVectorField(
35  *this, ir.dl_scalar, _exact_velocity, "Exact_velocity_");
36  this->addEvaluatedField(_exact_elec_pot);
37 
38  // Get external magnetic vector
39  const auto ext_magn_vct
40  = user_params.get<Teuchos::Array<double>>("External Magnetic Field");
41  for (int dim = 0; dim < 3; ++dim)
42  {
43  _B += ext_magn_vct[dim] * ext_magn_vct[dim];
44  }
45  _B = std::sqrt(_B);
46 
47  // Hartmann Number
48  _Ha = _B * _L * sqrt(_sigma / (_rho * _nu));
49 
50  // Closure model name
51  this->setName("Exact Solution Hartmann Problem");
52 }
53 
54 //---------------------------------------------------------------------------//
55 template<class EvalType, class Traits, int NumSpaceDim>
56 void HartmannProblemExact<EvalType, Traits, NumSpaceDim>::postRegistrationSetup(
57  typename Traits::SetupData sd, PHX::FieldManager<Traits>&)
58 {
59  _ir_index = panzer::getIntegrationRuleIndex(_ir_degree, (*sd.worksets_)[0]);
60 }
61 
62 //---------------------------------------------------------------------------//
63 template<class EvalType, class Traits, int NumSpaceDim>
64 void HartmannProblemExact<EvalType, Traits, NumSpaceDim>::evaluateFields(
65  typename Traits::EvalData workset)
66 {
67  _ip_coords = workset.int_rules[_ir_index]->ip_coordinates;
68  auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
69  workset.num_cells);
70  Kokkos::parallel_for(this->getName(), policy, *this);
71 }
72 
73 //---------------------------------------------------------------------------//
74 template<class EvalType, class Traits, int NumSpaceDim>
75 KOKKOS_INLINE_FUNCTION void
76 HartmannProblemExact<EvalType, Traits, NumSpaceDim>::operator()(
77  const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team) const
78 {
79  const int cell = team.league_rank();
80  const int num_point = _exact_elec_pot.extent(1);
81  using Kokkos::cosh;
82  using Kokkos::sqrt;
83 
84  Kokkos::parallel_for(
85  Kokkos::TeamThreadRange(team, 0, num_point), [&](const int point) {
86  // Note: an analytical solution is only available for
87  // the velocity and the electric potential.
88  _exact_lagrange_pressure(cell, point) = 0.0;
89  const double y = _ip_coords(cell, point, 1);
90 
91  _exact_velocity[0](cell, point) = (cosh(_Ha) - cosh(_Ha * y / _L))
92  / (cosh(_Ha) - 1.0);
93  _exact_velocity[1](cell, point) = 0.0;
94  if (num_space_dim == 3)
95  _exact_velocity[2](cell, point) = 0.0;
96  _exact_elec_pot(cell, point) = 0.0;
97  });
98 }
99 
100 //---------------------------------------------------------------------------//
101 
102 } // end namespace ClosureModel
103 } // end namespace VertexCFD
104 
105 #endif // end
106  // VERTEXCFD_CLOSURE_HARTMANNPROBLEMEXACT_IMPL_HPP
VertexCFD
Definition: tstMethodManufacturedSolutionBC.cpp:23