VertexCFD  0.0-dev
VertexCFD_Closure_VariableDiffusionFlux_impl.hpp
1 #ifndef VERTEXCFD_CLOSURE_VARIABLEDIFFUSIONFLUX_IMPL_HPP
2 #define VERTEXCFD_CLOSURE_VARIABLEDIFFUSIONFLUX_IMPL_HPP
3 
4 #include <utils/VertexCFD_Utils_VectorField.hpp>
5 
6 #include <Panzer_HierarchicParallelism.hpp>
7 
8 namespace VertexCFD
9 {
10 namespace ClosureModel
11 {
12 //---------------------------------------------------------------------------//
13 template<class EvalType, class Traits>
14 VariableDiffusionFlux<EvalType, Traits>::VariableDiffusionFlux(
15  const panzer::IntegrationRule& ir,
16  const Teuchos::ParameterList& closure_params,
17  const std::string& flux_prefix,
18  const std::string& gradient_prefix)
19  : _variable_name(closure_params.get<std::string>("Field Name"))
20  , _equation_name(closure_params.get<std::string>("Equation Name"))
21  , _var_diff_flux(flux_prefix + "DIFFUSION_FLUX_" + _equation_name,
22  ir.dl_vector)
23  , _diffusivity_var("diffusivity_" + _variable_name, ir.dl_scalar)
24  , _grad_var(gradient_prefix + "GRAD_" + _variable_name, ir.dl_vector)
25 {
26  // Add evaludated fields
27  this->addEvaluatedField(_var_diff_flux);
28 
29  // Add dependent fields
30  this->addDependentField(_diffusivity_var);
31  this->addDependentField(_grad_var);
32 
33  // Closure model name
34  this->setName(_equation_name + " Diffusion Flux");
35 }
36 
37 //---------------------------------------------------------------------------//
38 template<class EvalType, class Traits>
39 void VariableDiffusionFlux<EvalType, Traits>::evaluateFields(
40  typename Traits::EvalData workset)
41 {
42  auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
43  workset.num_cells);
44  Kokkos::parallel_for(this->getName(), policy, *this);
45 }
46 
47 //---------------------------------------------------------------------------//
48 template<class EvalType, class Traits>
49 KOKKOS_INLINE_FUNCTION void VariableDiffusionFlux<EvalType, Traits>::operator()(
50  const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team) const
51 {
52  const int cell = team.league_rank();
53  const int num_point = _var_diff_flux.extent(1);
54  const int num_grad_dim = _var_diff_flux.extent(2);
55 
56  Kokkos::parallel_for(
57  Kokkos::TeamThreadRange(team, 0, num_point), [&](const int point) {
58  // Loop over spatial dimension
59  for (int i = 0; i < num_grad_dim; ++i)
60  {
61  _var_diff_flux(cell, point, i) = _diffusivity_var(cell, point)
62  * _grad_var(cell, point, i);
63  }
64  });
65 }
66 
67 //---------------------------------------------------------------------------//
68 
69 } // end namespace ClosureModel
70 } // end namespace VertexCFD
71 
72 #endif // end VERTEXCFD_CLOSURE_VARIABLEDIFFUSIONFLUX_IMPL_HPP
VertexCFD
Definition: tstMethodManufacturedSolutionBC.cpp:23