VertexCFD  0.0-dev
VertexCFD_InitialCondition_MHDVortexProblem_impl.hpp
1 #ifndef VERTEXCFD_INITIALCONDITION_MHDVORTEXPROBLEM_IMPL_HPP
2 #define VERTEXCFD_INITIALCONDITION_MHDVORTEXPROBLEM_IMPL_HPP
3 
4 #include <Panzer_GlobalIndexer.hpp>
5 #include <Panzer_HierarchicParallelism.hpp>
6 #include <Panzer_PureBasis.hpp>
7 #include <Panzer_Workset_Utilities.hpp>
8 #include <utils/VertexCFD_Utils_VectorField.hpp>
9 
10 #include "utils/VertexCFD_Utils_Constants.hpp"
11 
12 #include <cmath>
13 #include <string>
14 
15 namespace VertexCFD
16 {
17 namespace InitialCondition
18 {
19 //---------------------------------------------------------------------------//
20 template<class EvalType, class Traits, int NumSpaceDim>
21 MHDVortexProblem<EvalType, Traits, NumSpaceDim>::MHDVortexProblem(
22  const Teuchos::ParameterList& params, const panzer::PureBasis& basis)
23  : _lagrange_pressure("lagrange_pressure", basis.functional)
24  , _basis_name(basis.name())
25 {
26  const auto vel_0 = params.get<Teuchos::Array<double>>("velocity_0");
27  const auto xy_0 = params.get<Teuchos::Array<double>>("center_0");
28  for (int dim = 0; dim < 2; ++dim)
29  {
30  _vel_0[dim] = vel_0[dim];
31  _xy_0[dim] = xy_0[dim];
32  }
33 
34  this->addEvaluatedField(_lagrange_pressure);
35  this->addUnsharedField(_lagrange_pressure.fieldTag().clone());
36 
37  Utils::addEvaluatedVectorField(
38  *this, basis.functional, _velocity, "velocity_", true);
39 
40  Utils::addEvaluatedVectorField(*this,
41  basis.functional,
42  _induced_magnetic_field,
43  "induced_magnetic_field_",
44  true);
45 
46  this->setName("MHDVortexProblem " + std::to_string(num_space_dim)
47  + "D Initial Condition");
48 }
49 
50 //---------------------------------------------------------------------------//
51 template<class EvalType, class Traits, int NumSpaceDim>
52 void MHDVortexProblem<EvalType, Traits, NumSpaceDim>::postRegistrationSetup(
53  typename Traits::SetupData sd, PHX::FieldManager<Traits>&)
54 {
55  _basis_index = panzer::getPureBasisIndex(
56  _basis_name, (*sd.worksets_)[0], this->wda);
57 }
58 
59 //---------------------------------------------------------------------------//
60 template<class EvalType, class Traits, int NumSpaceDim>
61 void MHDVortexProblem<EvalType, Traits, NumSpaceDim>::evaluateFields(
62  typename Traits::EvalData workset)
63 {
64  _basis_coords = this->wda(workset).bases[_basis_index]->basis_coordinates;
65  auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
66  workset.num_cells);
67  Kokkos::parallel_for(this->getName(), policy, *this);
68 }
69 
70 //---------------------------------------------------------------------------//
71 template<class EvalType, class Traits, int NumSpaceDim>
72 KOKKOS_INLINE_FUNCTION void
73 MHDVortexProblem<EvalType, Traits, NumSpaceDim>::operator()(
74  const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team) const
75 {
76  const int cell = team.league_rank();
77  const int num_basis = _lagrange_pressure.extent(1);
78  using std::exp;
79  using std::sqrt;
80 
81  Kokkos::parallel_for(
82  Kokkos::TeamThreadRange(team, 0, num_basis), [&](const int basis) {
83  const double x = _basis_coords(cell, basis, 0);
84  const double y = _basis_coords(cell, basis, 1);
85  const double r2 = (x - _xy_0[0]) * (x - _xy_0[0])
86  + (y - _xy_0[1]) * (y - _xy_0[1]);
87 
88  _lagrange_pressure(cell, basis)
89  = 1.0 + 0.5 * exp(1.) * (1. - r2 * exp(-r2));
90 
91  _induced_magnetic_field[0](cell, basis) = exp(0.5 * (1.0 - r2))
92  * (_xy_0[1] - y);
93  _induced_magnetic_field[1](cell, basis) = exp(0.5 * (1.0 - r2))
94  * (x - _xy_0[0]);
95 
96  _velocity[0](cell, basis) = _induced_magnetic_field[0](cell, basis)
97  + _vel_0[0];
98  _velocity[1](cell, basis) = _induced_magnetic_field[1](cell, basis)
99  + _vel_0[1];
100  if (num_space_dim > 2)
101  {
102  _induced_magnetic_field[2](cell, basis) = 0.0;
103  _velocity[2](cell, basis) = 0.0;
104  }
105  });
106 }
107 
108 //---------------------------------------------------------------------------//
109 
110 } // end namespace InitialCondition
111 } // end namespace VertexCFD
112 
113 #endif // end VERTEXCFD_INITIALCONDITION_MHDVORTEXPROBLEM_IMPL_HPP
VertexCFD
Definition: tstMethodManufacturedSolutionBC.cpp:23