1 #ifndef VERTEXCFD_INITIALCONDITION_INCOMPRESSIBLETAYLORGREENVORTEX_IMPL_HPP
2 #define VERTEXCFD_INITIALCONDITION_INCOMPRESSIBLETAYLORGREENVORTEX_IMPL_HPP
4 #include <utils/VertexCFD_Utils_VectorField.hpp>
6 #include <Panzer_GlobalIndexer.hpp>
7 #include <Panzer_HierarchicParallelism.hpp>
8 #include <Panzer_PureBasis.hpp>
9 #include <Panzer_Workset_Utilities.hpp>
11 #include "utils/VertexCFD_Utils_Constants.hpp"
13 #include <Teuchos_Array.hpp>
20 namespace InitialCondition
23 template<
class EvalType,
class Traits,
int NumSpaceDim>
24 IncompressibleTaylorGreenVortex<EvalType, Traits, NumSpaceDim>::
25 IncompressibleTaylorGreenVortex(
const panzer::PureBasis& basis)
26 : _lagrange_pressure(
"lagrange_pressure", basis.functional)
27 , _basis_name(basis.name())
29 this->addEvaluatedField(_lagrange_pressure);
30 this->addUnsharedField(_lagrange_pressure.fieldTag().clone());
32 Utils::addEvaluatedVectorField(
33 *
this, basis.functional, _velocity,
"velocity_",
true);
35 this->setName(
"Taylor-Green Vortex Initial Condition");
39 template<
class EvalType,
class Traits,
int NumSpaceDim>
40 void IncompressibleTaylorGreenVortex<EvalType, Traits, NumSpaceDim>::
41 postRegistrationSetup(
typename Traits::SetupData sd,
42 PHX::FieldManager<Traits>&)
44 _basis_index = panzer::getPureBasisIndex(
45 _basis_name, (*sd.worksets_)[0], this->wda);
49 template<
class EvalType,
class Traits,
int NumSpaceDim>
50 void IncompressibleTaylorGreenVortex<EvalType, Traits, NumSpaceDim>::evaluateFields(
51 typename Traits::EvalData workset)
53 _basis_coords = this->wda(workset).bases[_basis_index]->basis_coordinates;
54 auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
56 Kokkos::parallel_for(this->getName(), policy, *
this);
60 template<
class EvalType,
class Traits,
int NumSpaceDim>
61 KOKKOS_INLINE_FUNCTION
void
62 IncompressibleTaylorGreenVortex<EvalType, Traits, NumSpaceDim>::operator()(
63 const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team)
const
65 const int cell = team.league_rank();
66 const int num_basis = _lagrange_pressure.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 _lagrange_pressure(cell, basis) = -0.25
76 * (cos(2.0 * x) + cos(2.0 * y));
77 _velocity[0](cell, basis) = cos(x) * sin(y);
78 _velocity[1](cell, basis) = -sin(x) * cos(y);
80 if (num_space_dim == 3)
82 const double z = _basis_coords(cell, basis, 2);
83 _lagrange_pressure(cell, basis) *= (0.25 * cos(2.0 * z + 2.0));
84 _velocity[0](cell, basis) *= cos(z);
85 _velocity[1](cell, basis) *= cos(z);