1 #ifndef VERTEXCFD_CLOSURE_INCOMPRESSIBLETAYLORGREENVORTEXEXACTSOLUTION_IMPL_HPP
2 #define VERTEXCFD_CLOSURE_INCOMPRESSIBLETAYLORGREENVORTEXEXACTSOLUTION_IMPL_HPP
4 #include <utils/VertexCFD_Utils_VectorField.hpp>
6 #include "Panzer_GlobalIndexer.hpp"
7 #include "Panzer_PureBasis.hpp"
8 #include "Panzer_Workset_Utilities.hpp"
9 #include <Panzer_HierarchicParallelism.hpp>
11 #include <Teuchos_Array.hpp>
18 namespace ClosureModel
21 template<
class EvalType,
class Traits,
int NumSpaceDim>
22 IncompressibleTaylorGreenVortexExactSolution<EvalType, Traits, NumSpaceDim>::
23 IncompressibleTaylorGreenVortexExactSolution(
24 const panzer::IntegrationRule& ir,
25 const Teuchos::ParameterList& closure_params)
26 : _lagrange_pressure(
"Exact_lagrange_pressure", ir.dl_scalar)
27 , _nu(closure_params.get<double>(
"Kinematic Viscosity"))
28 , _ir_degree(ir.cubature_degree)
30 this->addEvaluatedField(_lagrange_pressure);
31 Utils::addEvaluatedVectorField(
32 *
this, ir.dl_scalar, _velocity,
"Exact_velocity_");
34 this->setName(
"Incompressible Taylor Green Vortex Exact Solution");
38 template<
class EvalType,
class Traits,
int NumSpaceDim>
39 void IncompressibleTaylorGreenVortexExactSolution<EvalType, Traits, NumSpaceDim>::
40 postRegistrationSetup(
typename Traits::SetupData sd,
41 PHX::FieldManager<Traits>&)
43 _ir_index = panzer::getIntegrationRuleIndex(_ir_degree, (*sd.worksets_)[0]);
47 template<
class EvalType,
class Traits,
int NumSpaceDim>
48 void IncompressibleTaylorGreenVortexExactSolution<EvalType, Traits, NumSpaceDim>::
49 evaluateFields(
typename Traits::EvalData workset)
51 auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
54 _Ft = exp(-2.0 * _nu * workset.time);
56 _ip_coords = workset.int_rules[_ir_index]->ip_coordinates;
57 Kokkos::parallel_for(this->getName(), policy, *
this);
61 template<
class EvalType,
class Traits,
int NumSpaceDim>
62 KOKKOS_INLINE_FUNCTION
void
63 IncompressibleTaylorGreenVortexExactSolution<EvalType, Traits, NumSpaceDim>::
64 operator()(
const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team)
const
66 const int cell = team.league_rank();
67 const int num_point = _lagrange_pressure.extent(1);
74 Kokkos::TeamThreadRange(team, 0, num_point), [&](
const int point) {
75 const double x = _ip_coords(cell, point, 0);
76 const double y = _ip_coords(cell, point, 1);
78 _lagrange_pressure(cell, point)
79 = -0.25 * (cos(2.0 * x) + cos(2.0 * y)) * _Ft * _Ft;
80 _velocity[0](cell, point) = cos(x) * sin(y) * _Ft;
81 _velocity[1](cell, point) = -sin(x) * cos(y) * _Ft;
90 #endif // VERTEXCFD_CLOSURE_INCOMPRESSIBLETAYLORGREENVORTEXEXACTSOLUTION_IMPL_HPP