VertexCFD  0.0-dev
VertexCFD_Closure_VariableTauSUPG_impl.hpp
1 #ifndef VERTEXCFD_CLOSURE_SCALARTAUSUPG_IMPL_HPP
2 #define VERTEXCFD_CLOSURE_SCALARTAUSUPG_IMPL_HPP
3 
4 #include "utils/VertexCFD_Utils_VectorField.hpp"
5 
6 #include <Panzer_HierarchicParallelism.hpp>
7 #include <Panzer_Workset_Utilities.hpp>
8 
9 #include <Teuchos_StandardParameterEntryValidators.hpp>
10 
11 #include <string>
12 
13 namespace VertexCFD
14 {
15 namespace ClosureModel
16 {
17 //---------------------------------------------------------------------------//
18 template<class EvalType, class Traits, int NumSpaceDim>
19 VariableTauSUPG<EvalType, Traits, NumSpaceDim>::VariableTauSUPG(
20  const panzer::IntegrationRule& ir,
21  const Teuchos::ParameterList& closure_params,
22  const Teuchos::ParameterList& supg_params)
23  : _field_name(closure_params.get<std::string>("Field Name"))
24  , _tau_supg("tau_supg_" + _field_name, ir.dl_scalar)
25  , _element_length("element_length", ir.dl_vector)
26  , _diffusivity("diffusivity_" + _field_name, ir.dl_scalar)
27  , _alpha(0.5)
28  , _ir_degree(ir.cubature_degree)
29 {
30  // Get tau model type for scalar equation
31  {
32  const auto type_validator = Teuchos::rcp(
33  new Teuchos::StringToIntegralParameterEntryValidator<TauModel>(
34  Teuchos::tuple<std::string>("Steady", "Transient", "NoSUPG"),
35  "Transient"));
36 
37  _tau_model = type_validator->getIntegralValue(
38  supg_params.get<std::string>("Tau model"));
39 
40  // SUPG coefficient
41  if (supg_params.isType<double>("SUPG coefficient"))
42  _alpha = supg_params.get<double>("SUPG coefficient");
43  }
44 
45  // Evaluated fields
46  this->addEvaluatedField(_tau_supg);
47 
48  // Dependent fields
49  Utils::addDependentVectorField(*this, ir.dl_scalar, _velocity, "velocity_");
50  this->addDependentField(_element_length);
51  this->addDependentField(_diffusivity);
52 
53  this->setName("Variable Tau SUPG " + std::to_string(num_space_dim) + "D");
54 }
55 
56 //---------------------------------------------------------------------------//
57 template<class EvalType, class Traits, int NumSpaceDim>
58 void VariableTauSUPG<EvalType, Traits, NumSpaceDim>::postRegistrationSetup(
59  typename Traits::SetupData sd, PHX::FieldManager<Traits>&)
60 {
61  _ir_index = panzer::getIntegrationRuleIndex(_ir_degree, (*sd.worksets_)[0]);
62 }
63 
64 //---------------------------------------------------------------------------//
65 template<class EvalType, class Traits, int NumSpaceDim>
66 void VariableTauSUPG<EvalType, Traits, NumSpaceDim>::evaluateFields(
67  typename Traits::EvalData workset)
68 {
69  // Get FEM basis order
70  _basis_order = workset.bases[_ir_index]->basis_layout->getBasis()->order();
71 
72  // Get time step
73  const double dt = workset.step_size;
74  _time_constant = 4.0 / (dt * dt);
75 
76  // Allocate space for thread-local temporaries in parallel for
77  size_t bytes;
78  if constexpr (Sacado::IsADType<scalar_type>::value)
79  {
80  const int fad_size = Kokkos::dimension_scalar(_tau_supg.get_view());
81  bytes = scratch_view::shmem_size(
82  _tau_supg.extent(1), NUM_TMPS, fad_size);
83  }
84  else
85  {
86  bytes = scratch_view::shmem_size(_tau_supg.extent(1), NUM_TMPS);
87  }
88 
89  auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
90  workset.num_cells);
91  policy.set_scratch_size(0, Kokkos::PerTeam(bytes));
92  Kokkos::parallel_for(this->getName(), policy, *this);
93 }
94 
95 //---------------------------------------------------------------------------//
96 template<class EvalType, class Traits, int NumSpaceDim>
97 KOKKOS_INLINE_FUNCTION void
98 VariableTauSUPG<EvalType, Traits, NumSpaceDim>::operator()(
99  const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team) const
100 {
101  const int cell = team.league_rank();
102  const int num_point = _tau_supg.extent(1);
103  const double tol = 1.0e-10;
104  using Kokkos::sqrt;
105 
106  scratch_view tmp_field;
107  if constexpr (Sacado::IsADType<scalar_type>::value)
108  {
109  const int fad_size = Kokkos::dimension_scalar(_tau_supg.get_view());
110  tmp_field
111  = scratch_view(team.team_shmem(), num_point, NUM_TMPS, fad_size);
112  }
113  else
114  {
115  tmp_field = scratch_view(team.team_shmem(), num_point, NUM_TMPS);
116  }
117 
118  Kokkos::parallel_for(
119  Kokkos::TeamThreadRange(team, 0, num_point), [&](const int point) {
120  // Initialize local variables
121  auto&& u2_over_h2 = tmp_field(point, U2_OVER_H2);
122  u2_over_h2 = 0.0;
123  double h2 = 0.0;
124 
125  // Compute terms `h^2` and '\frac{u^2}{h^2}`
126  for (int dim = 0; dim < num_space_dim; ++dim)
127  {
128  h2 += _element_length(cell, point, dim)
129  * _element_length(cell, point, dim);
130  u2_over_h2 += _velocity[dim](cell, point)
131  * _velocity[dim](cell, point)
132  / (_element_length(cell, point, dim)
133  * _element_length(cell, point, dim));
134  }
135 
136  // SUPG coefficient for momentum and continuity equations
137  if (_tau_model == TauModel::NoSUPG)
138  {
139  _tau_supg(cell, point) = 0.0;
140  }
141  else if (_tau_model == TauModel::Steady)
142  {
143  _tau_supg(cell, point)
144  = 1.0
145  / sqrt(4.0 * u2_over_h2
146  + 9.0 * 16.0 * _diffusivity(cell, point)
147  * _diffusivity(cell, point) / (h2 * h2)
148  + tol);
149  }
150  else // Transient
151  {
152  _tau_supg(cell, point)
153  = 1.0
154  / sqrt(_time_constant + 4.0 * u2_over_h2
155  + 9.0 * 16.0 * _diffusivity(cell, point)
156  * _diffusivity(cell, point) / (h2 * h2));
157  }
158  _tau_supg(cell, point) *= _alpha / _basis_order;
159  });
160 }
161 
162 //---------------------------------------------------------------------------//
163 
164 } // end namespace ClosureModel
165 } // end namespace VertexCFD
166 
167 #endif // end VERTEXCFD_CLOSURE_SCALARTAUSUPG_IMPL_HPP
VertexCFD
Definition: tstMethodManufacturedSolutionBC.cpp:23