VertexCFD  0.0-dev
VertexCFD_Closure_RADTimeDerivative_impl.hpp
1 #ifndef VERTEXCFD_CLOSURE_RADTIMEDERIVATIVE_IMPL_HPP
2 #define VERTEXCFD_CLOSURE_RADTIMEDERIVATIVE_IMPL_HPP
3 
4 #include "utils/VertexCFD_Utils_ScalarFieldView.hpp"
5 
6 #include <Panzer_HierarchicParallelism.hpp>
7 
8 namespace VertexCFD
9 {
10 namespace ClosureModel
11 {
12 //---------------------------------------------------------------------------//
13 template<class EvalType, class Traits>
14 RADTimeDerivative<EvalType, Traits>::RADTimeDerivative(
15  const panzer::IntegrationRule& ir,
16  const SpeciesProperties::ConstantSpeciesProperties& species_prop)
17  : _num_species(species_prop.numSpecies())
18 {
19  // Evaluated variable
20  Utils::addEvaluatedScalarFieldView(*this,
21  ir.dl_scalar,
22  _num_species,
23  _dqdt_rad,
24  "DQDT_reaction_advection_diffusion_"
25  "equation_");
26 
27  // Dependent variable
28  Utils::addDependentScalarFieldView(
29  *this, ir.dl_scalar, _num_species, _dxdt_species, "DXDT_species_");
30 
31  // Closure model name
32  this->setName("RAD Time Derivative");
33 }
34 
35 //---------------------------------------------------------------------------//
36 template<class EvalType, class Traits>
37 void RADTimeDerivative<EvalType, Traits>::evaluateFields(
38  typename Traits::EvalData workset)
39 {
40  auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
41  workset.num_cells);
42  Kokkos::parallel_for(this->getName(), policy, *this);
43 }
44 
45 //---------------------------------------------------------------------------//
46 template<class EvalType, class Traits>
47 KOKKOS_INLINE_FUNCTION void RADTimeDerivative<EvalType, Traits>::operator()(
48  const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team) const
49 {
50  const int cell = team.league_rank();
51  const int num_point = _dqdt_rad[0].extent(1);
52 
53  Kokkos::parallel_for(
54  Kokkos::TeamThreadRange(team, 0, num_point), [&](const int point) {
55  for (int num = 0; num < _num_species; ++num)
56  _dqdt_rad[num](cell, point) = _dxdt_species[num](cell, point);
57  });
58 }
59 
60 //---------------------------------------------------------------------------//
61 
62 } // end namespace ClosureModel
63 } // namespace VertexCFD
64 
65 #endif // end VERTEXCFD_CLOSURE_RADTIMEDERIVATIVE_IMPL_HPP
VertexCFD
Definition: tstMethodManufacturedSolutionBC.cpp:23