VertexCFD  0.0-dev
VertexCFD_Closure_IncompressibleLSVOFSurfaceTensionForce_impl.hpp
1 #ifndef VERTEXCFD_CLOSURE_INCOMPRESSIBLELSVOFSURFACETENSIONFORCE_IMPL_HPP
2 #define VERTEXCFD_CLOSURE_INCOMPRESSIBLELSVOFSURFACETENSIONFORCE_IMPL_HPP
3 
4 #include "utils/VertexCFD_Utils_VectorFieldView.hpp"
5 
6 #include "utils/VertexCFD_Utils_SmoothMath.hpp"
7 #include "utils/VertexCFD_Utils_VectorField.hpp"
8 
9 #include <Panzer_HierarchicParallelism.hpp>
10 
11 #include <string>
12 
13 namespace VertexCFD
14 {
15 namespace ClosureModel
16 {
17 //---------------------------------------------------------------------------//
18 template<class EvalType, class Traits, int NumSpaceDim>
19 IncompressibleLSVOFSurfaceTensionForce<EvalType, Traits, NumSpaceDim>::
20  IncompressibleLSVOFSurfaceTensionForce(
21  const panzer::IntegrationRule& ir,
22  const Teuchos::ParameterList& closure_params,
23  const std::vector<std::string>& phase_names,
24  const std::string& flux_prefix,
25  const std::string& gradient_prefix)
26  : _sigma_alpha(closure_params.get<double>("Surface Tension Coefficient"))
27  , num_phases(closure_params.get<int>("Number of Phases") - 1)
28 {
29  // Evaluated fields
30  Utils::addEvaluatedVectorField(*this, ir.dl_vector, _surface_tension_flux,
31  flux_prefix + "SURFACE_TENSION_FLUX_"
32  "momentum_");
33 
34  // Dependent fields
35  this->setName("Incompressible LSVOF Surface Tension Force "
36  + std::to_string(num_space_dim) + "D");
37 
38  Utils::addDependentVectorFieldView(*this,
39  ir.dl_vector,
40  num_phases,
41  _phase_grad,
42  gradient_prefix + "GRAD_",
43  phase_names);
44 }
45 
46 //---------------------------------------------------------------------------//
47 template<class EvalType, class Traits, int NumSpaceDim>
48 void IncompressibleLSVOFSurfaceTensionForce<EvalType, Traits, NumSpaceDim>::
49  evaluateFields(typename Traits::EvalData workset)
50 {
51  // Allocate space for thread-local temporaries
52  size_t bytes;
53  if constexpr (Sacado::IsADType<scalar_type>::value)
54  {
55  const int fad_size
56  = Kokkos::dimension_scalar(_surface_tension_flux[0].get_view());
57  bytes = scratch_view::shmem_size(_surface_tension_flux[0].extent(1),
58  fad_size);
59 
60  bytes += scratch_view_Fs::shmem_size(
61  _surface_tension_flux[0].extent(1), fad_size);
62 
63  bytes += scratch_view_normal::shmem_size(
64  _surface_tension_flux[0].extent(1), num_space_dim, fad_size);
65  }
66  else
67  {
68  bytes = scratch_view::shmem_size(_surface_tension_flux[0].extent(1));
69 
70  bytes
71  += scratch_view_Fs::shmem_size(_surface_tension_flux[0].extent(1));
72 
73  bytes += scratch_view_normal::shmem_size(
74  _surface_tension_flux[0].extent(1), num_space_dim);
75  }
76  auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
77  workset.num_cells);
78  policy.set_scratch_size(0, Kokkos::PerTeam(bytes));
79  Kokkos::parallel_for(this->getName(), policy, *this);
80 }
81 
82 //---------------------------------------------------------------------------//
83 template<class EvalType, class Traits, int NumSpaceDim>
84 KOKKOS_INLINE_FUNCTION void
85 IncompressibleLSVOFSurfaceTensionForce<EvalType, Traits, NumSpaceDim>::operator()(
86  const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team) const
87 {
88  using Kokkos::sqrt;
89 
90  const int cell = team.league_rank();
91  const int num_point = _surface_tension_flux[0].extent(1);
92  const double max_tol = 1e-10;
93 
94  scratch_view scratch_data;
95  scratch_view_Fs scratch_data_force;
96  scratch_view_normal scratch_data_normal;
97 
98  if constexpr (Sacado::IsADType<scalar_type>::value)
99  {
100  const int fad_size
101  = Kokkos::dimension_scalar(_surface_tension_flux[0].get_view());
102 
103  scratch_data = scratch_view(team.team_shmem(), num_point, fad_size);
104 
105  scratch_data_force
106  = scratch_view_Fs(team.team_shmem(), num_point, fad_size);
107 
108  scratch_data_normal = scratch_view_normal(
109  team.team_shmem(), num_point, num_space_dim, fad_size);
110  }
111  else
112  {
113  scratch_data = scratch_view(team.team_shmem(), num_point);
114 
115  scratch_data_force = scratch_view_Fs(team.team_shmem(), num_point);
116 
117  scratch_data_normal
118  = scratch_view_normal(team.team_shmem(), num_point, num_space_dim);
119  }
120 
121  Kokkos::parallel_for(
122  Kokkos::TeamThreadRange(team, 0, num_point), [&](const int point) {
123  auto&& tmp_surf_force = scratch_data_force(point);
124  auto&& mag_grad_scalar = scratch_data(point);
125  auto&& unit_normal
126  = Kokkos::subview(scratch_data_normal, point, Kokkos::ALL);
127 
128  for (int dim = 0; dim < num_space_dim; ++dim)
129  {
130  for (int d = 0; d < num_space_dim; ++d)
131  {
132  _surface_tension_flux[d](cell, point, dim) = 0.0;
133  }
134  }
135 
136  for (int n = 0; n < num_phases; ++n)
137  {
138  mag_grad_scalar = 0.0;
139 
140  for (int dim = 0; dim < num_space_dim; ++dim)
141  {
142  unit_normal[dim] = 0.0;
143  mag_grad_scalar
144  += pow(_phase_grad[n](cell, point, dim), 2.0);
145  }
146  mag_grad_scalar
147  = sqrt(SmoothMath::max(mag_grad_scalar, max_tol, 0));
148 
149  for (int dim = 0; dim < num_space_dim; ++dim)
150  {
151  unit_normal[dim] = _phase_grad[n](cell, point, dim)
152  / mag_grad_scalar;
153  }
154 
155  for (int dim = 0; dim < num_space_dim; ++dim)
156  {
157  for (int d = 0; d < num_space_dim; ++d)
158  {
159  tmp_surf_force = -unit_normal[dim] * unit_normal[d];
160 
161  if (d == dim)
162  {
163  tmp_surf_force += 1.0;
164  }
165 
166  tmp_surf_force *= _sigma_alpha * mag_grad_scalar;
167 
168  _surface_tension_flux[d](cell, point, dim)
169  += tmp_surf_force;
170  }
171  }
172  }
173  });
174 }
175 
176 //---------------------------------------------------------------------------//
177 
178 } // end namespace ClosureModel
179 } // end namespace VertexCFD
180 
181 #endif // end
182  // VERTEXCFD_CLOSURE_INCOMPRESSIBLELSVOFSURFACETENSIONFORCE_IMPL_HPP
VertexCFD
Definition: tstMethodManufacturedSolutionBC.cpp:23