VertexCFD  0.0-dev
VertexCFD_Closure_RADFissionSource_impl.hpp
1 #ifndef VERTEXCFD_CLOSURE_RADFISSIONSOURCE_IMPL_HPP
2 #define VERTEXCFD_CLOSURE_RADFISSIONSOURCE_IMPL_HPP
3 
4 #include "utils/VertexCFD_Utils_ScalarFieldView.hpp"
5 
6 #include <Panzer_HierarchicParallelism.hpp>
7 #include <Teuchos_StandardParameterEntryValidators.hpp>
8 
9 #include <string>
10 
11 namespace VertexCFD
12 {
13 namespace ClosureModel
14 {
15 //---------------------------------------------------------------------------//
16 template<class EvalType, class Traits>
17 RADFissionSource<EvalType, Traits>::RADFissionSource(
18  const panzer::IntegrationRule& ir,
19  const SpeciesProperties::ConstantSpeciesProperties& species_prop,
20  const std::string& neutron_flux_name)
21  : _num_species(species_prop.numSpecies())
22  , _neutron_flux(neutron_flux_name, ir.dl_scalar)
23  , _gamma(Kokkos::ViewAllocateWithoutInitializing("atoms_per_species"),
24  _num_species)
25  , _xs(species_prop.fissionCrossSection())
26  , _avagadro(6.02214076e23)
27 {
28  // Copy number of atoms vector from host to device
29  auto gamma_host = Kokkos::create_mirror_view(Kokkos::HostSpace{}, _gamma);
30  gamma_host = species_prop.atomsPerSpecies();
31  Kokkos::deep_copy(_gamma, gamma_host);
32 
33  // Evaluated fields
34  Utils::addEvaluatedScalarFieldView(*this,
35  ir.dl_scalar,
36  _num_species,
37  _fission_source,
38  "FISSION_SOURCE_reaction_advection_"
39  "diffusion_equation_");
40 
41  // Dependent fields
42  this->addDependentField(_neutron_flux);
43 
44  this->setName("RAD Fission Source Term");
45 }
46 
47 //---------------------------------------------------------------------------//
48 template<class EvalType, class Traits>
49 void RADFissionSource<EvalType, Traits>::evaluateFields(
50  typename Traits::EvalData workset)
51 {
52  auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
53  workset.num_cells);
54  Kokkos::parallel_for(this->getName(), policy, *this);
55 }
56 
57 //---------------------------------------------------------------------------//
58 template<class EvalType, class Traits>
59 KOKKOS_INLINE_FUNCTION void RADFissionSource<EvalType, Traits>::operator()(
60  const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team) const
61 {
62  const int cell = team.league_rank();
63  const int num_point = _fission_source[0].extent(1);
64 
65  Kokkos::parallel_for(
66  Kokkos::TeamThreadRange(team, 0, num_point), [&](const int point) {
67  for (int num = 0; num < _num_species; ++num)
68  {
69  _fission_source[num](cell, point) = _neutron_flux(cell, point)
70  * _gamma[num] * _xs
71  / _avagadro;
72  }
73  });
74 }
75 
76 //---------------------------------------------------------------------------//
77 
78 } // end namespace ClosureModel
79 } // end namespace VertexCFD
80 
81 #endif // end VERTEXCFD_CLOSURE_RADFISSIONSOURCE_IMPL_HPP
VertexCFD
Definition: tstMethodManufacturedSolutionBC.cpp:23