VertexCFD  0.0-dev
VertexCFD_BoundaryState_ElectricPotentialFixed_impl.hpp
1 #ifndef VERTEXCFD_BOUNDARYSTATE_ELECTRICPOTENTIALFIXED_IMPL_HPP
2 #define VERTEXCFD_BOUNDARYSTATE_ELECTRICPOTENTIALFIXED_IMPL_HPP
3 
4 #include <utils/VertexCFD_Utils_VectorField.hpp>
5 
6 #include <Panzer_HierarchicParallelism.hpp>
7 
8 namespace VertexCFD
9 {
10 namespace BoundaryCondition
11 {
12 //---------------------------------------------------------------------------//
13 // This function should be used for Dirichlet boundary conditions. A ramping
14 // in time can be enabled.
15 //---------------------------------------------------------------------------//
16 template<class EvalType, class Traits>
17 ElectricPotentialFixed<EvalType, Traits>::ElectricPotentialFixed(
18  const panzer::IntegrationRule& ir, const Teuchos::ParameterList& bc_params)
19  : _boundary_electric_potential("BOUNDARY_electric_potential", ir.dl_scalar)
20  , _boundary_grad_electric_potential("BOUNDARY_GRAD_electric_potential",
21  ir.dl_vector)
22  , _num_grad_dim(ir.spatial_dimension)
23  , _grad_electric_potential("GRAD_electric_potential", ir.dl_vector)
24 {
25  // Calculate the coefficients 'a' and 'b' for the linear time ramping
26  // sc(t) = a * t + b
27  _time_init = bc_params.isType<double>("Time Initial")
28  ? bc_params.get<double>(
29  "Time "
30  "Initial")
31  : 0.0;
32  _time_final = bc_params.isType<double>("Time Final")
33  ? bc_params.get<double>("Time Final")
34  : 1.0E-06;
35  const double dt = _time_final - _time_init;
36  const auto sc_final = bc_params.get<double>("Final Value");
37  const auto sc_init = bc_params.isType<double>("Initial Value")
38  ? bc_params.get<double>("Initial Value")
39  : sc_final;
40  _a_sc = (sc_final - sc_init) / dt;
41  _b_sc = sc_init - _a_sc * _time_init;
42 
43  // Add evaluated fields
44  this->addEvaluatedField(_boundary_electric_potential);
45  this->addEvaluatedField(_boundary_grad_electric_potential);
46 
47  // Add dependent fields
48  this->addDependentField(_grad_electric_potential);
49 
50  this->setName("Boundary State Electric Potential Fixed "
51  + std::to_string(_num_grad_dim) + "D");
52 }
53 
54 //---------------------------------------------------------------------------//
55 template<class EvalType, class Traits>
56 void ElectricPotentialFixed<EvalType, Traits>::evaluateFields(
57  typename Traits::EvalData workset)
58 {
59  // Get time and make sure it only varies between '_time_init' and
60  // '_time_final'
61  _time = std::max(workset.time, _time_init);
62  _time = std::min(_time, _time_final);
63 
64  auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
65  workset.num_cells);
66  Kokkos::parallel_for(this->getName(), policy, *this);
67 }
68 
69 //---------------------------------------------------------------------------//
70 template<class EvalType, class Traits>
71 KOKKOS_INLINE_FUNCTION void
72 ElectricPotentialFixed<EvalType, Traits>::operator()(
73  const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team) const
74 {
75  const int cell = team.league_rank();
76  const int num_point = _grad_electric_potential.extent(1);
77 
78  Kokkos::parallel_for(
79  Kokkos::TeamThreadRange(team, 0, num_point), [&](const int point) {
80  // Assign time-dependent boundary values
81  _boundary_electric_potential(cell, point) = _a_sc * _time + _b_sc;
82 
83  // Assign gradient.
84  for (int d = 0; d < _num_grad_dim; ++d)
85  {
86  _boundary_grad_electric_potential(cell, point, d)
87  = _grad_electric_potential(cell, point, d);
88  }
89  });
90 }
91 
92 //---------------------------------------------------------------------------//
93 
94 } // end namespace BoundaryCondition
95 } // end namespace VertexCFD
96 
97 #endif // end VERTEXCFD_BOUNDARYSTATE_ELECTRICPOTENTIALFIXED_IMPL_HPP
VertexCFD
Definition: tstMethodManufacturedSolutionBC.cpp:23