1 #ifndef VERTEXCFD_INITIALCONDITION_MHDVORTEXPROBLEM_IMPL_HPP
2 #define VERTEXCFD_INITIALCONDITION_MHDVORTEXPROBLEM_IMPL_HPP
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>
10 #include "utils/VertexCFD_Utils_Constants.hpp"
17 namespace InitialCondition
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())
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)
30 _vel_0[dim] = vel_0[dim];
31 _xy_0[dim] = xy_0[dim];
34 this->addEvaluatedField(_lagrange_pressure);
35 this->addUnsharedField(_lagrange_pressure.fieldTag().clone());
37 Utils::addEvaluatedVectorField(
38 *
this, basis.functional, _velocity,
"velocity_",
true);
40 Utils::addEvaluatedVectorField(*
this,
42 _induced_magnetic_field,
43 "induced_magnetic_field_",
46 this->setName(
"MHDVortexProblem " + std::to_string(num_space_dim)
47 +
"D Initial Condition");
51 template<
class EvalType,
class Traits,
int NumSpaceDim>
52 void MHDVortexProblem<EvalType, Traits, NumSpaceDim>::postRegistrationSetup(
53 typename Traits::SetupData sd, PHX::FieldManager<Traits>&)
55 _basis_index = panzer::getPureBasisIndex(
56 _basis_name, (*sd.worksets_)[0], this->wda);
60 template<
class EvalType,
class Traits,
int NumSpaceDim>
61 void MHDVortexProblem<EvalType, Traits, NumSpaceDim>::evaluateFields(
62 typename Traits::EvalData workset)
64 _basis_coords = this->wda(workset).bases[_basis_index]->basis_coordinates;
65 auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
67 Kokkos::parallel_for(this->getName(), policy, *
this);
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
76 const int cell = team.league_rank();
77 const int num_basis = _lagrange_pressure.extent(1);
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]);
88 _lagrange_pressure(cell, basis)
89 = 1.0 + 0.5 * exp(1.) * (1. - r2 * exp(-r2));
91 _induced_magnetic_field[0](cell, basis) = exp(0.5 * (1.0 - r2))
93 _induced_magnetic_field[1](cell, basis) = exp(0.5 * (1.0 - r2))
96 _velocity[0](cell, basis) = _induced_magnetic_field[0](cell, basis)
98 _velocity[1](cell, basis) = _induced_magnetic_field[1](cell, basis)
100 if (num_space_dim > 2)
102 _induced_magnetic_field[2](cell, basis) = 0.0;
103 _velocity[2](cell, basis) = 0.0;
113 #endif // end VERTEXCFD_INITIALCONDITION_MHDVORTEXPROBLEM_IMPL_HPP