1 #ifndef VERTEXCFD_CLOSURE_RADFISSIONSOURCE_IMPL_HPP
2 #define VERTEXCFD_CLOSURE_RADFISSIONSOURCE_IMPL_HPP
4 #include "utils/VertexCFD_Utils_ScalarFieldView.hpp"
6 #include <Panzer_HierarchicParallelism.hpp>
7 #include <Teuchos_StandardParameterEntryValidators.hpp>
13 namespace ClosureModel
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"),
25 , _xs(species_prop.fissionCrossSection())
26 , _avagadro(6.02214076e23)
29 auto gamma_host = Kokkos::create_mirror_view(Kokkos::HostSpace{}, _gamma);
30 gamma_host = species_prop.atomsPerSpecies();
31 Kokkos::deep_copy(_gamma, gamma_host);
34 Utils::addEvaluatedScalarFieldView(*
this,
38 "FISSION_SOURCE_reaction_advection_"
39 "diffusion_equation_");
42 this->addDependentField(_neutron_flux);
44 this->setName(
"RAD Fission Source Term");
48 template<
class EvalType,
class Traits>
49 void RADFissionSource<EvalType, Traits>::evaluateFields(
50 typename Traits::EvalData workset)
52 auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
54 Kokkos::parallel_for(this->getName(), policy, *
this);
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
62 const int cell = team.league_rank();
63 const int num_point = _fission_source[0].extent(1);
66 Kokkos::TeamThreadRange(team, 0, num_point), [&](
const int point) {
67 for (
int num = 0; num < _num_species; ++num)
69 _fission_source[num](cell, point) = _neutron_flux(cell, point)
81 #endif // end VERTEXCFD_CLOSURE_RADFISSIONSOURCE_IMPL_HPP