VertexCFD  0.0-dev
VertexCFD_Closure_IncompressibleCLSNonReconstructedNormal_impl.hpp
1 #ifndef VERTEXCFD_CLOSURE_INCOMPRESSIBLECLSNONRECONSTRUCTEDNORMAL_IMPL_HPP
2 #define VERTEXCFD_CLOSURE_INCOMPRESSIBLECLSNONRECONSTRUCTEDNORMAL_IMPL_HPP
3 
4 #include <Panzer_HierarchicParallelism.hpp>
5 
6 #include <string>
7 
8 namespace VertexCFD
9 {
10 namespace ClosureModel
11 {
12 //---------------------------------------------------------------------------//
13 template<class EvalType, class Traits, int NumSpaceDim>
14 IncompressibleCLSNonReconstructedNormal<EvalType, Traits, NumSpaceDim>::
15  IncompressibleCLSNonReconstructedNormal(const panzer::IntegrationRule& ir)
16  : _q("CLS_interface_normal", ir.dl_vector)
17  , _grad_phi_0("STAR_GRAD_level_set", ir.dl_vector)
18 {
19  // Evaluated fields
20  this->addEvaluatedField(_q);
21 
22  // Dependent fields
23  this->addDependentField(_grad_phi_0);
24 
25  this->setName("Incompressible CLS Non-Reconstructed Interface Normal "
26  + std::to_string(num_space_dim) + "D");
27 }
28 
29 //---------------------------------------------------------------------------//
30 template<class EvalType, class Traits, int NumSpaceDim>
31 void IncompressibleCLSNonReconstructedNormal<EvalType, Traits, NumSpaceDim>::
32  evaluateFields(typename Traits::EvalData workset)
33 {
34  auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
35  workset.num_cells);
36  Kokkos::parallel_for(this->getName(), policy, *this);
37 }
38 
39 //---------------------------------------------------------------------------//
40 template<class EvalType, class Traits, int NumSpaceDim>
41 KOKKOS_INLINE_FUNCTION void
42 IncompressibleCLSNonReconstructedNormal<EvalType, Traits, NumSpaceDim>::operator()(
43  const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team) const
44 {
45  const int cell = team.league_rank();
46  const int num_point = _q.extent(1);
47 
48  const double tol = 1e-15;
49  const int N = num_space_dim - 1;
50 
51  using Kokkos::pow;
52 
53  Kokkos::parallel_for(
54  Kokkos::TeamThreadRange(team, 0, num_point), [&](const int point) {
55  // Store mag(grad(phi_0)) in _q[N]
56  _q(cell, point, N) = 0.0;
57 
58  for (int dim = 0; dim < num_space_dim; ++dim)
59  {
60  _q(cell, point, N) += pow(_grad_phi_0(cell, point, dim), 2.0);
61  }
62 
63  _q(cell, point, N) = pow(_q(cell, point, N) + tol, 0.5);
64 
65  // Evaluate interface normal
66  for (int dim = 0; dim < num_space_dim; ++dim)
67  {
68  _q(cell, point, dim) = _grad_phi_0(cell, point, dim)
69  / _q(cell, point, N);
70  }
71  });
72 }
73 
74 //---------------------------------------------------------------------------//
75 
76 } // end namespace ClosureModel
77 } // end namespace VertexCFD
78 
79 #endif // end
80  // VERTEXCFD_CLOSURE_INCOMPRESSIBLECLSNONRECONSTRUCTEDNORMAL_IMPL_HPP
VertexCFD
Definition: tstMethodManufacturedSolutionBC.cpp:23