1 #ifndef VERTEXCFD_CLOSURE_SOLIDELECTRICPOTENTIALEXACTSOLUTION_IMPL_HPP
2 #define VERTEXCFD_CLOSURE_SOLIDELECTRICPOTENTIALEXACTSOLUTION_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 SolidElectricPotentialExactSolution<EvalType, Traits, NumSpaceDim>::
23 SolidElectricPotentialExactSolution(
24 const panzer::IntegrationRule& ir,
25 const Teuchos::ParameterList& closure_params)
26 : _lagrange_pressure(
"Exact_lagrange_pressure", ir.dl_scalar)
27 , _electric_potential(
"Exact_electric_potential", ir.dl_scalar)
28 , _ir_degree(ir.cubature_degree)
31 this->addEvaluatedField(_lagrange_pressure);
32 Utils::addEvaluatedVectorField(
33 *
this, ir.dl_scalar, _velocity,
"Exact_velocity_");
34 this->addEvaluatedField(_electric_potential);
37 _sigma = closure_params.get<
double>(
"Electrical Conductivity Coefficient");
38 _phi_right = closure_params.get<
double>(
"Right Electric Boundary Value");
39 _phi_left = closure_params.get<
double>(
"Left Electric Boundary Value");
42 this->setName(
"Conduction Exact Solution");
46 template<
class EvalType,
class Traits,
int NumSpaceDim>
47 void SolidElectricPotentialExactSolution<EvalType, Traits, NumSpaceDim>::
48 postRegistrationSetup(
typename Traits::SetupData sd,
49 PHX::FieldManager<Traits>&)
51 _ir_index = panzer::getIntegrationRuleIndex(_ir_degree, (*sd.worksets_)[0]);
55 template<
class EvalType,
class Traits,
int NumSpaceDim>
56 void SolidElectricPotentialExactSolution<EvalType, Traits, NumSpaceDim>::
57 evaluateFields(
typename Traits::EvalData workset)
59 auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
62 _ip_coords = workset.int_rules[_ir_index]->ip_coordinates;
63 Kokkos::parallel_for(this->getName(), policy, *
this);
67 template<
class EvalType,
class Traits,
int NumSpaceDim>
68 KOKKOS_INLINE_FUNCTION
void
69 SolidElectricPotentialExactSolution<EvalType, Traits, NumSpaceDim>::operator()(
70 const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team)
const
72 const int cell = team.league_rank();
73 const int num_point = _electric_potential.extent(1);
78 Kokkos::TeamThreadRange(team, 0, num_point), [&](
const int point) {
82 _lagrange_pressure(cell, point) = 0.0;
83 for (
int dim = 0; dim < num_space_dim; ++dim)
84 _velocity[dim](cell, point) = 0.0;
87 const double x = _ip_coords(cell, point, 0);
88 _electric_potential(cell, point)
89 = _phi_left * exp(log(_phi_right / _phi_left) * x);
98 #endif // VERTEXCFD_CLOSURE_SOLIDELECTRICPOTENTIALEXACTSOLUTION_IMPL_HPP