VertexCFD  0.0-dev
VertexCFD_Closure_MHDVortexProblemExact_impl.hpp
1 #ifndef VERTEXCFD_CLOSURE_MHDVORTEXPROBLEMEXACT_IMPL_HPP
2 #define VERTEXCFD_CLOSURE_MHDVORTEXPROBLEMEXACT_IMPL_HPP
3 
4 #include "utils/VertexCFD_Utils_VectorField.hpp"
5 
6 #include <Panzer_GlobalIndexer.hpp>
7 #include <Panzer_HierarchicParallelism.hpp>
8 #include <Panzer_PureBasis.hpp>
9 #include <Panzer_Workset_Utilities.hpp>
10 
11 #include <string>
12 
13 namespace VertexCFD
14 {
15 namespace ClosureModel
16 {
17 //---------------------------------------------------------------------------//
18 template<class EvalType, class Traits, int NumSpaceDim>
19 MHDVortexProblemExact<EvalType, Traits, NumSpaceDim>::MHDVortexProblemExact(
20  const panzer::IntegrationRule& ir,
21  const Teuchos::ParameterList& vortex_params)
22  : _ir_degree(ir.cubature_degree)
23  , _lagrange_pressure("Exact_lagrange_pressure", ir.dl_scalar)
24 {
25  const auto vel_0 = vortex_params.get<Teuchos::Array<double>>("velocity_0");
26  const auto xy_0 = vortex_params.get<Teuchos::Array<double>>("center_0");
27  for (int dim = 0; dim < 2; ++dim)
28  {
29  _vel_0[dim] = vel_0[dim];
30  _xy_0[dim] = xy_0[dim];
31  }
32 
33  // Add evaluated fields
34  this->addEvaluatedField(_lagrange_pressure);
35  Utils::addEvaluatedVectorField(
36  *this, ir.dl_scalar, _velocity, "Exact_velocity_");
37  Utils::addEvaluatedVectorField(*this,
38  ir.dl_scalar,
39  _induced_magnetic_field,
40  "Exact_induced_magnetic_field_");
41 
42  // Closure model name
43  this->setName("MHD Vortex Problem Exact Solution "
44  + std::to_string(num_space_dim) + "D.");
45 }
46 
47 //---------------------------------------------------------------------------//
48 template<class EvalType, class Traits, int NumSpaceDim>
49 void MHDVortexProblemExact<EvalType, Traits, NumSpaceDim>::postRegistrationSetup(
50  typename Traits::SetupData sd, PHX::FieldManager<Traits>&)
51 {
52  _ir_index = panzer::getIntegrationRuleIndex(_ir_degree, (*sd.worksets_)[0]);
53 }
54 
55 //---------------------------------------------------------------------------//
56 template<class EvalType, class Traits, int NumSpaceDim>
57 void MHDVortexProblemExact<EvalType, Traits, NumSpaceDim>::evaluateFields(
58  typename Traits::EvalData workset)
59 {
60  _time = workset.time;
61  _ip_coords = workset.int_rules[_ir_index]->ip_coordinates;
62  auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
63  workset.num_cells);
64  Kokkos::parallel_for(this->getName(), policy, *this);
65 }
66 
67 //---------------------------------------------------------------------------//
68 template<class EvalType, class Traits, int NumSpaceDim>
69 KOKKOS_INLINE_FUNCTION void
70 MHDVortexProblemExact<EvalType, Traits, NumSpaceDim>::operator()(
71  const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team) const
72 {
73  const int cell = team.league_rank();
74  const int num_point = _lagrange_pressure.extent(1);
75  using std::exp;
76  const double exp_one = exp(1.0);
77 
78  Kokkos::parallel_for(
79  Kokkos::TeamThreadRange(team, 0, num_point), [&](const int point) {
80  // Coordinnates
81  const double x = _ip_coords(cell, point, 0);
82  const double y = _ip_coords(cell, point, 1);
83  const double r2 = (x - _time * _vel_0[0] - _xy_0[0])
84  * (x - _time * _vel_0[0] - _xy_0[0])
85  + (y - _time * _vel_0[0] - _xy_0[1])
86  * (y - _time * _vel_0[0] - _xy_0[1]);
87 
88  // Exact solutions
89  _lagrange_pressure(cell, point)
90  = 1.0 + 0.5 * exp_one * (1.0 - r2 * exp(-r2));
91 
92  _induced_magnetic_field[0](cell, point) = exp(0.5 * (1.0 - r2))
93  * (_xy_0[1] - y);
94  _induced_magnetic_field[1](cell, point) = exp(0.5 * (1.0 - r2))
95  * (x - _xy_0[0]);
96 
97  _velocity[0](cell, point) = _induced_magnetic_field[0](cell, point)
98  + _vel_0[0];
99  _velocity[1](cell, point) = _induced_magnetic_field[1](cell, point);
100 
101  if (num_space_dim == 3)
102  {
103  _induced_magnetic_field[2](cell, point) = 0.0;
104  _velocity[2](cell, point) = 0.0;
105  }
106  });
107 }
108 
109 //---------------------------------------------------------------------------//
110 
111 } // end namespace ClosureModel
112 } // end namespace VertexCFD
113 
114 #endif // end
115  // VERTEXCFD_CLOSURE_MHDVORTEXPROBLEMEXACT_IMPL_HPP
VertexCFD
Definition: tstMethodManufacturedSolutionBC.cpp:23