VertexCFD  0.0-dev
VertexCFD_Closure_SolidElectricPotentialExactSolution_impl.hpp
1 #ifndef VERTEXCFD_CLOSURE_SOLIDELECTRICPOTENTIALEXACTSOLUTION_IMPL_HPP
2 #define VERTEXCFD_CLOSURE_SOLIDELECTRICPOTENTIALEXACTSOLUTION_IMPL_HPP
3 
4 #include "utils/VertexCFD_Utils_VectorField.hpp"
5 
6 #include <Panzer_GlobalIndexer.hpp>
7 #include <Panzer_HierarchicParallelism.hpp>
8 #include <Panzer_PureBasis.hpp>
9 #include <Panzer_Workset_Utilities.hpp>
10 
11 #include <Teuchos_Array.hpp>
12 
13 #include <cmath>
14 #include <string>
15 
16 namespace VertexCFD
17 {
18 namespace ClosureModel
19 {
20 //---------------------------------------------------------------------------//
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)
29 {
30  // Add evaluated field
31  this->addEvaluatedField(_lagrange_pressure);
32  Utils::addEvaluatedVectorField(
33  *this, ir.dl_scalar, _velocity, "Exact_velocity_");
34  this->addEvaluatedField(_electric_potential);
35 
36  // Get heat source value and coefficient for exact solution
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");
40 
41  // Set name
42  this->setName("Conduction Exact Solution");
43 }
44 
45 //---------------------------------------------------------------------------//
46 template<class EvalType, class Traits, int NumSpaceDim>
47 void SolidElectricPotentialExactSolution<EvalType, Traits, NumSpaceDim>::
48  postRegistrationSetup(typename Traits::SetupData sd,
49  PHX::FieldManager<Traits>&)
50 {
51  _ir_index = panzer::getIntegrationRuleIndex(_ir_degree, (*sd.worksets_)[0]);
52 }
53 
54 //---------------------------------------------------------------------------//
55 template<class EvalType, class Traits, int NumSpaceDim>
56 void SolidElectricPotentialExactSolution<EvalType, Traits, NumSpaceDim>::
57  evaluateFields(typename Traits::EvalData workset)
58 {
59  auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
60  workset.num_cells);
61 
62  _ip_coords = workset.int_rules[_ir_index]->ip_coordinates;
63  Kokkos::parallel_for(this->getName(), policy, *this);
64 }
65 
66 //---------------------------------------------------------------------------//
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
71 {
72  const int cell = team.league_rank();
73  const int num_point = _electric_potential.extent(1);
74  using Kokkos::exp;
75  using Kokkos::log;
76 
77  Kokkos::parallel_for(
78  Kokkos::TeamThreadRange(team, 0, num_point), [&](const int point) {
79  // NS variables are still defined on a solid
80  // region because they are expected by the Tempus
81  // observer that compute the error norms.
82  _lagrange_pressure(cell, point) = 0.0;
83  for (int dim = 0; dim < num_space_dim; ++dim)
84  _velocity[dim](cell, point) = 0.0;
85 
86  // With electric_potential-dependent sigma
87  const double x = _ip_coords(cell, point, 0);
88  _electric_potential(cell, point)
89  = _phi_left * exp(log(_phi_right / _phi_left) * x);
90  });
91 }
92 
93 //---------------------------------------------------------------------------//
94 
95 } // end namespace ClosureModel
96 } // end namespace VertexCFD
97 
98 #endif // VERTEXCFD_CLOSURE_SOLIDELECTRICPOTENTIALEXACTSOLUTION_IMPL_HPP
VertexCFD
Definition: tstMethodManufacturedSolutionBC.cpp:23