VertexCFD  0.0-dev
VertexCFD_Closure_RADAdvectionFlux_impl.hpp
1 #ifndef VERTEXCFD_CLOSURE_RADADVECTIONFLUX_IMPL_HPP
2 #define VERTEXCFD_CLOSURE_RADADVECTIONFLUX_IMPL_HPP
3 
4 #include "utils/VertexCFD_Utils_ScalarFieldView.hpp"
5 #include "utils/VertexCFD_Utils_VectorField.hpp"
6 
7 #include <Panzer_HierarchicParallelism.hpp>
8 
9 #include <string>
10 
11 namespace VertexCFD
12 {
13 namespace ClosureModel
14 {
15 //---------------------------------------------------------------------------//
16 template<class EvalType, class Traits, int NumSpaceDim>
17 RADAdvectionFlux<EvalType, Traits, NumSpaceDim>::RADAdvectionFlux(
18  const panzer::IntegrationRule& ir,
19  const SpeciesProperties::ConstantSpeciesProperties& species_prop,
20  const std::string& flux_prefix,
21  const std::string& field_prefix)
22  : _num_species(species_prop.numSpecies())
23 {
24  // Evaluated fields
25  Utils::addEvaluatedScalarFieldView(*this, ir.dl_vector, _num_species, _species_flux,
26  flux_prefix + "ADVECTION_FLUX_"
27  "reaction_advection_diffusion_equation_");
28 
29  // Dependent fields
30  Utils::addDependentScalarFieldView(
31  *this, ir.dl_scalar, _num_species, _species, field_prefix + "species_");
32 
33  Utils::addDependentVectorField(
34  *this, ir.dl_scalar, _velocity, field_prefix + "velocity_");
35 
36  this->setName("RAD Equation Advection Flux "
37  + std::to_string(num_space_dim) + "D");
38 }
39 
40 //---------------------------------------------------------------------------//
41 template<class EvalType, class Traits, int NumSpaceDim>
42 void RADAdvectionFlux<EvalType, Traits, NumSpaceDim>::evaluateFields(
43  typename Traits::EvalData workset)
44 {
45  auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
46  workset.num_cells);
47  Kokkos::parallel_for(this->getName(), policy, *this);
48 }
49 
50 //---------------------------------------------------------------------------//
51 template<class EvalType, class Traits, int NumSpaceDim>
52 KOKKOS_INLINE_FUNCTION void
53 RADAdvectionFlux<EvalType, Traits, NumSpaceDim>::operator()(
54  const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team) const
55 {
56  const int cell = team.league_rank();
57  const int num_point = _species[0].extent(1);
58 
59  Kokkos::parallel_for(Kokkos::TeamThreadRange(team, 0, num_point),
60  [&](const int point) {
61  for (int dim = 0; dim < num_space_dim; ++dim)
62  {
63  for (int sp = 0; sp < _num_species; ++sp)
64  {
65  _species_flux[sp](cell, point, dim)
66  = _species[sp](cell, point)
67  * _velocity[dim](cell, point);
68  }
69  }
70  });
71 }
72 
73 //---------------------------------------------------------------------------//
74 
75 } // end namespace ClosureModel
76 } // end namespace VertexCFD
77 
78 #endif // end VERTEXCFD_CLOSURE_RADADVECTIONFLUX_IMPL_HPP
VertexCFD
Definition: tstMethodManufacturedSolutionBC.cpp:23