1 #ifndef VERTEXCFD_BOUNDARYSTATE_INCOMPRESSIBLEROTATINGWALL_IMPL_HPP
2 #define VERTEXCFD_BOUNDARYSTATE_INCOMPRESSIBLEROTATINGWALL_IMPL_HPP
4 #include "utils/VertexCFD_Utils_VectorField.hpp"
6 #include <Panzer_HierarchicParallelism.hpp>
7 #include <Panzer_Workset_Utilities.hpp>
11 namespace BoundaryCondition
14 template<
class EvalType,
class Traits,
int NumSpaceDim>
15 IncompressibleRotatingWall<EvalType, Traits, NumSpaceDim>::IncompressibleRotatingWall(
16 const panzer::IntegrationRule& ir,
17 const Teuchos::ParameterList& fluid_params,
18 const Teuchos::ParameterList& bc_params,
20 : _boundary_lagrange_pressure(
"BOUNDARY_lagrange_pressure", ir.dl_scalar)
21 , _boundary_grad_lagrange_pressure(
"BOUNDARY_GRAD_lagrange_pressure",
23 , _boundary_temperature(
"BOUNDARY_temperature", ir.dl_scalar)
24 , _boundary_grad_temperature(
"BOUNDARY_GRAD_temperature", ir.dl_vector)
25 , _ir_degree(ir.cubature_degree)
26 , _set_lagrange_pressure(bc_params.isType<double>(
"Lagrange Pressure"))
27 , _lp_wall(std::numeric_limits<double>::quiet_NaN())
28 , _solve_temp(fluid_params.get<bool>(
"Build Temperature Equation"))
30 , _T_wall(std::numeric_limits<double>::quiet_NaN())
31 , _lagrange_pressure(
"lagrange_pressure", ir.dl_scalar)
32 , _grad_lagrange_pressure(
"GRAD_lagrange_pressure", ir.dl_vector)
33 , _grad_temperature(
"GRAD_temperature", ir.dl_vector)
34 , _normals(
"Side Normal", ir.dl_vector)
38 _time_init = bc_params.isType<
double>(
"Time Initial")
39 ? bc_params.get<
double>(
"Time Initial")
41 _time_final = bc_params.isType<
double>(
"Time Final")
42 ? bc_params.get<
double>(
"Time Final")
45 if (_set_lagrange_pressure)
46 _lp_wall = bc_params.get<
double>(
"Lagrange Pressure");
48 const auto angular_velocity_final
49 = bc_params.get<
double>(
"Angular Velocity");
50 const auto angular_velocity_init = bc_params.isType<
double>(
53 ? bc_params.get<
double>(
56 : angular_velocity_final;
57 const double dt = _time_final - _time_init;
58 _a_vel = (angular_velocity_final - angular_velocity_init) / dt;
59 _b_vel = angular_velocity_init - _a_vel * _time_init;
62 this->addEvaluatedField(_boundary_lagrange_pressure);
64 this->addEvaluatedField(_boundary_grad_lagrange_pressure);
66 Utils::addEvaluatedVectorField(
67 *
this, ir.dl_scalar, _boundary_velocity,
"BOUNDARY_velocity_");
69 Utils::addEvaluatedVectorField(*
this,
71 _boundary_grad_velocity,
72 "BOUNDARY_GRAD_velocity_");
75 _T_wall = bc_params.get<
double>(
"Wall Temperature");
76 this->addEvaluatedField(_boundary_temperature);
77 this->addEvaluatedField(_boundary_grad_temperature);
81 this->addDependentField(_lagrange_pressure);
84 this->addDependentField(_grad_lagrange_pressure);
85 this->addDependentField(_normals);
88 Utils::addDependentVectorField(
89 *
this, ir.dl_vector, _grad_velocity,
"GRAD_velocity_");
92 this->addDependentField(_grad_temperature);
94 this->setName(
"Boundary State Incompressible Rotating Wall "
95 + std::to_string(num_space_dim) +
"D");
99 template<
class EvalType,
class Traits,
int NumSpaceDim>
100 void IncompressibleRotatingWall<EvalType, Traits, NumSpaceDim>::postRegistrationSetup(
101 typename Traits::SetupData sd, PHX::FieldManager<Traits>&)
103 _ir_index = panzer::getIntegrationRuleIndex(_ir_degree, (*sd.worksets_)[0]);
107 template<
class EvalType,
class Traits,
int NumSpaceDim>
108 void IncompressibleRotatingWall<EvalType, Traits, NumSpaceDim>::evaluateFields(
109 typename Traits::EvalData workset)
113 const double time = workset.time < _time_init ? _time_init
114 : workset.time >= _time_final ? _time_final
116 _angular_velocity = _a_vel * time + _b_vel;
118 _ip_coords = workset.int_rules[_ir_index]->ip_coordinates;
119 auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
121 Kokkos::parallel_for(this->getName(), policy, *
this);
125 template<
class EvalType,
class Traits,
int NumSpaceDim>
126 KOKKOS_INLINE_FUNCTION
void
127 IncompressibleRotatingWall<EvalType, Traits, NumSpaceDim>::operator()(
128 const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team)
const
130 const int cell = team.league_rank();
131 const int num_point = _lagrange_pressure.extent(1);
132 const int num_grad_dim = _boundary_grad_velocity[0].extent(2);
134 Kokkos::parallel_for(
135 Kokkos::TeamThreadRange(team, 0, num_point), [&](
const int point) {
137 if (_set_lagrange_pressure)
138 _boundary_lagrange_pressure(cell, point) = _lp_wall;
141 _boundary_lagrange_pressure(cell, point)
142 = _lagrange_pressure(cell, point);
146 _boundary_temperature(cell, point) = _T_wall;
149 Kokkos::Array<double, num_space_dim> vel_bnd{};
150 vel_bnd[0] = -_angular_velocity * _ip_coords(cell, point, 1);
151 vel_bnd[1] = _angular_velocity * _ip_coords(cell, point, 0);
152 for (
int vel_dim = 0; vel_dim < num_space_dim; ++vel_dim)
153 _boundary_velocity[vel_dim](cell, point) = vel_bnd[vel_dim];
156 for (
int d = 0; d < num_grad_dim; ++d)
161 _boundary_grad_lagrange_pressure(cell, point, d)
162 = _grad_lagrange_pressure(cell, point, d);
164 for (
int grad_dim = 0; grad_dim < num_space_dim; ++grad_dim)
166 _boundary_grad_lagrange_pressure(cell, point, d)
167 -= (_grad_lagrange_pressure(cell, point, grad_dim)
168 * _normals(cell, point, grad_dim))
169 * _normals(cell, point, d);
173 for (
int vel_dim = 0; vel_dim < num_space_dim; ++vel_dim)
175 _boundary_grad_velocity[vel_dim](cell, point, d)
176 = _grad_velocity[vel_dim](cell, point, d);
181 _boundary_grad_temperature(cell, point, d)
182 = _grad_temperature(cell, point, d);
193 #endif // VERTEXCFD_BOUNDARYSTATE_INCOMPRESSIBLEROTATINGWALL_IMPL_HPP