1 #ifndef VERTEXCFD_CLOSURE_INCOMPRESSIBLESUPGEXACTSOLUTION_IMPL_HPP
2 #define VERTEXCFD_CLOSURE_INCOMPRESSIBLESUPGEXACTSOLUTION_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 <Teuchos_Array.hpp>
18 namespace ClosureModel
21 template<
class EvalType,
class Traits,
int NumSpaceDim>
22 IncompressibleSUPGExactSolution<EvalType, Traits, NumSpaceDim>::
23 IncompressibleSUPGExactSolution(
const panzer::IntegrationRule& ir,
24 const Teuchos::ParameterList& closure_params)
25 : _lagrange_pressure(
"Exact_lagrange_pressure", ir.dl_scalar)
26 , _temperature(
"Exact_temperature", ir.dl_scalar)
27 , _ir_degree(ir.cubature_degree)
28 , _D(closure_params.get<double>(
"Thermal diffusivity"))
29 , _u(closure_params.get<double>(
"Uniform advection velocity"))
31 this->addEvaluatedField(_lagrange_pressure);
32 this->addEvaluatedField(_temperature);
33 Utils::addEvaluatedVectorField(
34 *
this, ir.dl_scalar, _velocity,
"Exact_velocity_");
36 this->setName(
"Incompressible SUPG Exact Solution");
40 template<
class EvalType,
class Traits,
int NumSpaceDim>
41 void IncompressibleSUPGExactSolution<EvalType, Traits, NumSpaceDim>::
42 postRegistrationSetup(
typename Traits::SetupData sd,
43 PHX::FieldManager<Traits>&)
45 _ir_index = panzer::getIntegrationRuleIndex(_ir_degree, (*sd.worksets_)[0]);
49 template<
class EvalType,
class Traits,
int NumSpaceDim>
50 void IncompressibleSUPGExactSolution<EvalType, Traits, NumSpaceDim>::evaluateFields(
51 typename Traits::EvalData workset)
53 auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
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 IncompressibleSUPGExactSolution<EvalType, Traits, NumSpaceDim>::operator()(
64 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);
72 Kokkos::TeamThreadRange(team, 0, num_point), [&](
const int point) {
73 const double x = _ip_coords(cell, point, 0);
76 _lagrange_pressure(cell, point) = 0.0;
77 _velocity[0](cell, point) = 0.0;
78 _velocity[1](cell, point) = 0.0;
79 _temperature(cell, point)
80 = (x - (1.0 - exp(_u * x / _D)) / (1.0 - exp(_u / _D))) / _u;
89 #endif // VERTEXCFD_CLOSURE_INCOMPRESSIBLESUPGEXACTSOLUTION_IMPL_HPP