1 #ifndef VERTEXCFD_CLOSURE_CONDUCTIONEXACTSOLUTION_IMPL_HPP
2 #define VERTEXCFD_CLOSURE_CONDUCTIONEXACTSOLUTION_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>
23 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)
31 Utils::addEvaluatedVectorField(
32 *
this, ir.dl_scalar,
_velocity,
"Exact_velocity_");
36 _q = closure_params.get<
double>(
"Volumetric Heat Source Value");
37 _k = closure_params.get<
double>(
"Thermal Conductivity Coefficient");
38 _T_right = closure_params.get<
double>(
"Right Temperature Boundary Value");
41 this->setName(
"Conduction Exact Solution");
45 template<
class EvalType,
class Traits,
int NumSpaceDim>
47 typename Traits::SetupData sd, PHX::FieldManager<Traits>&)
49 _ir_index = panzer::getIntegrationRuleIndex(_ir_degree, (*sd.worksets_)[0]);
53 template<
class EvalType,
class Traits,
int NumSpaceDim>
55 typename Traits::EvalData workset)
57 auto policy = panzer::HP::inst().teamPolicy<
scalar_type, PHX::Device>(
60 _ip_coords = workset.int_rules[_ir_index]->ip_coordinates;
61 Kokkos::parallel_for(this->getName(), policy, *
this);
65 template<
class EvalType,
class Traits,
int NumSpaceDim>
66 KOKKOS_INLINE_FUNCTION
void
68 const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team)
const
70 const int cell = team.league_rank();
71 const int num_point = _temperature.extent(1);
75 Kokkos::TeamThreadRange(team, 0, num_point), [&](
const int point) {
79 _lagrange_pressure(cell, point) = 0.0;
80 _velocity[0](cell, point) = 0.0;
81 _velocity[1](cell, point) = 0.0;
82 if (num_space_dim == 3)
83 _velocity[2](cell, point) = 0.0;
86 const double x = _ip_coords(cell, point, 0);
87 _temperature(cell, point)
88 = _T_right * exp(_q / (6.0 * _k) * (1.0 - x * x * x));
97 #endif // VERTEXCFD_CLOSURE_CONDUCTIONEXACTSOLUTION_IMPL_HPP