1 #ifndef VERTEXCFD_INITIALCONDITION_INCOMPRESSIBLEVORTEXINBOX_IMPL_HPP
2 #define VERTEXCFD_INITIALCONDITION_INCOMPRESSIBLEVORTEXINBOX_IMPL_HPP
4 #include <Panzer_GlobalIndexer.hpp>
5 #include <Panzer_HierarchicParallelism.hpp>
6 #include <Panzer_PureBasis.hpp>
7 #include <Panzer_Workset_Utilities.hpp>
9 #include "utils/VertexCFD_Utils_Constants.hpp"
11 #include <Teuchos_Array.hpp>
18 namespace InitialCondition
21 template<
class EvalType,
class Traits>
22 IncompressibleVortexInBox<EvalType, Traits>::IncompressibleVortexInBox(
23 const panzer::PureBasis& basis)
24 : _velocity_0(
"velocity_0", basis.functional)
25 , _velocity_1(
"velocity_1", basis.functional)
26 , _lagrange_pressure(
"lagrange_pressure", basis.functional)
27 , _basis_name(basis.name())
29 this->addEvaluatedField(_velocity_0);
30 this->addEvaluatedField(_velocity_1);
31 this->addEvaluatedField(_lagrange_pressure);
32 this->addUnsharedField(_velocity_0.fieldTag().clone());
33 this->addUnsharedField(_velocity_1.fieldTag().clone());
34 this->addUnsharedField(_lagrange_pressure.fieldTag().clone());
35 this->setName(
"Vertex In the Box Initial Condition");
39 template<
class EvalType,
class Traits>
40 void IncompressibleVortexInBox<EvalType, Traits>::postRegistrationSetup(
41 typename Traits::SetupData sd, PHX::FieldManager<Traits>&)
43 _basis_index = panzer::getPureBasisIndex(
44 _basis_name, (*sd.worksets_)[0], this->wda);
48 template<
class EvalType,
class Traits>
49 void IncompressibleVortexInBox<EvalType, Traits>::evaluateFields(
50 typename Traits::EvalData workset)
52 _basis_coords = this->wda(workset).bases[_basis_index]->basis_coordinates;
53 auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
55 Kokkos::parallel_for(this->getName(), policy, *
this);
59 template<
class EvalType,
class Traits>
60 KOKKOS_INLINE_FUNCTION
void
61 IncompressibleVortexInBox<EvalType, Traits>::operator()(
62 const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team)
const
64 const int cell = team.league_rank();
65 const int num_basis = _velocity_0.extent(1);
72 Kokkos::TeamThreadRange(team, 0, num_basis), [&](
const int basis) {
73 const double x = _basis_coords(cell, basis, 0);
74 const double y = _basis_coords(cell, basis, 1);
75 _velocity_0(cell, basis) = -2.0 * cos(pi * y) * sin(pi * y)
76 * sin(pi * x) * sin(pi * x);
77 _velocity_1(cell, basis) = 2.0 * cos(pi * x) * sin(pi * x)
78 * sin(pi * y) * sin(pi * y);
79 _lagrange_pressure(cell, basis) = 0.0;
88 #endif // end VERTEXCFD_INITIALCONDITION_INCOMPRESSIBLEVORTEXINBOX_IMPL_HPP