1 #ifndef VERTEXCFD_BOUNDARYSTATE_VARIABLEFIXED_IMPL_HPP
2 #define VERTEXCFD_BOUNDARYSTATE_VARIABLEFIXED_IMPL_HPP
4 #include <Panzer_HierarchicParallelism.hpp>
8 namespace BoundaryCondition
13 template<
class EvalType,
class Traits>
14 VariableFixed<EvalType, Traits>::VariableFixed(
15 const panzer::IntegrationRule& ir,
16 const Teuchos::ParameterList& bc_params,
17 const std::string variable_name)
18 : _boundary_variable(
"BOUNDARY_" + variable_name, ir.dl_scalar)
19 , _boundary_grad_variable(
"BOUNDARY_GRAD_" + variable_name, ir.dl_vector)
20 , _grad_variable(
"GRAD_" + variable_name, ir.dl_vector)
21 , _fixed_value(bc_params.get<double>(variable_name +
" Value"))
22 , _num_grad_dim(ir.spatial_dimension)
26 _time_init = bc_params.isType<
double>(
"Time Initial")
27 ? bc_params.get<
double>(
"Time Initial")
29 _time_final = bc_params.isType<
double>(
"Time Final")
30 ? bc_params.get<
double>(
"Time Final")
33 const auto variable_init = bc_params.isType<
double>(
34 variable_name +
" Value "
36 ? bc_params.get<
double>(
37 variable_name +
" Value "
40 const double dt = _time_final - _time_init;
41 _a = (_fixed_value - variable_init) / dt;
42 _b = variable_init - _a * _time_init;
45 this->addEvaluatedField(_boundary_variable);
46 this->addEvaluatedField(_boundary_grad_variable);
49 this->addDependentField(_grad_variable);
51 this->setName(variable_name +
" Boundary State Variable Fixed "
52 + std::to_string(_num_grad_dim) +
"D");
56 template<
class EvalType,
class Traits>
57 void VariableFixed<EvalType, Traits>::evaluateFields(
58 typename Traits::EvalData workset)
62 const double time = workset.time < _time_init ? _time_init
63 : workset.time >= _time_final ? _time_final
65 _variable_value = _a * time + _b;
67 auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
69 Kokkos::parallel_for(this->getName(), policy, *
this);
73 template<
class EvalType,
class Traits>
74 KOKKOS_INLINE_FUNCTION
void VariableFixed<EvalType, Traits>::operator()(
75 const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team)
const
77 const int cell = team.league_rank();
78 const int num_point = _grad_variable.extent(1);
80 Kokkos::parallel_for(Kokkos::TeamThreadRange(team, 0, num_point),
81 [&](
const int point) {
83 _boundary_variable(cell, point) = _variable_value;
86 for (
int d = 0; d < _num_grad_dim; ++d)
88 _boundary_grad_variable(cell, point, d)
89 = _grad_variable(cell, point, d);
99 #endif // end VERTEXCFD_BOUNDARYSTATE_VARIABLEFIXED_IMPL_HPP