VertexCFD  0.0-dev
VertexCFD_BoundaryState_IncompressibleDirichlet_impl.hpp
1 #ifndef VERTEXCFD_BOUNDARYSTATE_INCOMPRESSIBLEDIRICHLET_IMPL_HPP
2 #define VERTEXCFD_BOUNDARYSTATE_INCOMPRESSIBLEDIRICHLET_IMPL_HPP
3 
4 #include <utils/VertexCFD_Utils_VectorField.hpp>
5 
6 #include <Panzer_HierarchicParallelism.hpp>
7 
8 #include <string>
9 
10 namespace VertexCFD
11 {
12 namespace BoundaryCondition
13 {
14 //---------------------------------------------------------------------------//
15 // This function should be used for Dirichlet boundary conditions. A ramping
16 // in time can be enabled for all variables or only a few. Be aware that the
17 // ramping in time does not include any logic with the characteristics and thus
18 // should only be used for non-transient runs.
19 //---------------------------------------------------------------------------//
20 template<class EvalType, class Traits, int NumSpaceDim>
21 IncompressibleDirichlet<EvalType, Traits, NumSpaceDim>::IncompressibleDirichlet(
22  const panzer::IntegrationRule& ir,
23  const Teuchos::ParameterList& fluid_params,
24  const Teuchos::ParameterList& bc_params,
25  const bool is_edac)
26  : _boundary_lagrange_pressure("BOUNDARY_lagrange_pressure", ir.dl_scalar)
27  , _boundary_grad_lagrange_pressure("BOUNDARY_GRAD_lagrange_pressure",
28  ir.dl_vector)
29  , _boundary_temperature("BOUNDARY_temperature", ir.dl_scalar)
30  , _boundary_grad_temperature("BOUNDARY_GRAD_temperature", ir.dl_vector)
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  , _solve_temp(fluid_params.get<bool>("Build Temperature Equation"))
35  , _is_edac(is_edac)
36 {
37  // Calculate the coefficients 'a' and 'b' for the linear time ramping
38  // f(t) = a * t + b for each variable
39  _time_init = bc_params.isType<double>("Time Initial")
40  ? bc_params.get<double>(
41  "Time "
42  "Initial")
43  : 0.0;
44  _time_final = bc_params.isType<double>("Time Final")
45  ? bc_params.get<double>("Time Final")
46  : 1.0E-06;
47  const double dt = _time_final - _time_init;
48 
49  for (int vel_dim = 0; vel_dim < num_space_dim; ++vel_dim)
50  {
51  const std::string vel_string = "velocity_" + std::to_string(vel_dim);
52  const auto vel_final = bc_params.get<double>(vel_string);
53  const auto vel_init = bc_params.isType<double>(vel_string + "_init")
54  ? bc_params.get<double>(vel_string + "_init")
55  : vel_final;
56  _a_vel[vel_dim] = (vel_final - vel_init) / dt;
57  _b_vel[vel_dim] = vel_init - _a_vel[vel_dim] * _time_init;
58  }
59 
60  // Add evaluated fields
61  this->addEvaluatedField(_boundary_lagrange_pressure);
62  if (_is_edac)
63  this->addEvaluatedField(_boundary_grad_lagrange_pressure);
64  if (_solve_temp)
65  {
66  _T_dirichlet = bc_params.get<double>("temperature");
67  this->addEvaluatedField(_boundary_temperature);
68  this->addEvaluatedField(_boundary_grad_temperature);
69  }
70  Utils::addEvaluatedVectorField(
71  *this, ir.dl_scalar, _boundary_velocity, "BOUNDARY_velocity_");
72 
73  Utils::addEvaluatedVectorField(*this,
74  ir.dl_vector,
75  _boundary_grad_velocity,
76  "BOUNDARY_GRAD_velocity_");
77 
78  // Add dependent fields
79  Utils::addDependentVectorField(
80  *this, ir.dl_vector, _grad_velocity, "GRAD_velocity_");
81  this->addDependentField(_lagrange_pressure);
82  if (_is_edac)
83  this->addDependentField(_grad_lagrange_pressure);
84  if (_solve_temp)
85  this->addDependentField(_grad_temperature);
86 
87  this->setName("Boundary State Incompressible Dirichlet "
88  + std::to_string(num_space_dim) + "D");
89 }
90 
91 //---------------------------------------------------------------------------//
92 template<class EvalType, class Traits, int NumSpaceDim>
93 void IncompressibleDirichlet<EvalType, Traits, NumSpaceDim>::evaluateFields(
94  typename Traits::EvalData workset)
95 {
96  // Get time and make sure it only varies between '_time_init' and
97  // '_time_final'
98  _time = std::max(workset.time, _time_init);
99  _time = std::min(_time, _time_final);
100 
101  auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
102  workset.num_cells);
103  Kokkos::parallel_for(this->getName(), policy, *this);
104 }
105 
106 //---------------------------------------------------------------------------//
107 template<class EvalType, class Traits, int NumSpaceDim>
108 KOKKOS_INLINE_FUNCTION void
109 IncompressibleDirichlet<EvalType, Traits, NumSpaceDim>::operator()(
110  const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team) const
111 {
112  const int cell = team.league_rank();
113  const int num_point = _lagrange_pressure.extent(1);
114  const int num_grad_dim = _boundary_grad_velocity[0].extent(2);
115 
116  Kokkos::parallel_for(
117  Kokkos::TeamThreadRange(team, 0, num_point), [&](const int point) {
118  // Assign time-dependent boundary values
119  _boundary_lagrange_pressure(cell, point)
120  = _lagrange_pressure(cell, point);
121  for (int vel_dim = 0; vel_dim < num_space_dim; ++vel_dim)
122  {
123  _boundary_velocity[vel_dim](cell, point)
124  = _a_vel[vel_dim] * _time + _b_vel[vel_dim];
125  }
126  if (_solve_temp)
127  _boundary_temperature(cell, point) = _T_dirichlet;
128 
129  // Set gradients
130  for (int d = 0; d < num_grad_dim; ++d)
131  {
132  if (_is_edac)
133  {
134  _boundary_grad_lagrange_pressure(cell, point, d)
135  = _grad_lagrange_pressure(cell, point, d);
136  }
137 
138  for (int vel_dim = 0; vel_dim < num_space_dim; ++vel_dim)
139  {
140  _boundary_grad_velocity[vel_dim](cell, point, d)
141  = _grad_velocity[vel_dim](cell, point, d);
142  }
143 
144  if (_solve_temp)
145  {
146  _boundary_grad_temperature(cell, point, d)
147  = _grad_temperature(cell, point, d);
148  }
149  }
150  });
151 }
152 
153 //---------------------------------------------------------------------------//
154 
155 } // end namespace BoundaryCondition
156 } // end namespace VertexCFD
157 
158 #endif // end VERTEXCFD_BOUNDARYSTATE_INCOMPRESSIBLEDIRICHLET_IMPL_HPP
VertexCFD
Definition: tstMethodManufacturedSolutionBC.cpp:23