1 #ifndef VERTEXCFD_CONSTANTSPECIESPROPERTIES_HPP
2 #define VERTEXCFD_CONSTANTSPECIESPROPERTIES_HPP
4 #include <Teuchos_ParameterList.hpp>
6 #include <Kokkos_Core.hpp>
10 namespace SpeciesProperties
21 const Teuchos::ParameterList& model_params,
22 const Teuchos::ParameterList& reaction_params)
23 : _num_species(model_params.get<
int>(
"Number of Species"))
24 , _build_bateman(model_params.isType<
bool>(
"Build Reaction Bateman "
26 ? model_params.get<
bool>(
"Build Reaction Bateman "
29 , _build_advection(model_params.isType<
bool>(
"Build Advection")
30 ? model_params.get<
bool>(
"Build Advection")
32 , _build_diffusion(model_params.isType<
bool>(
"Build Diffusion")
33 ? model_params.get<
bool>(
"Build Diffusion")
35 , _build_fission_source(model_params.isType<
bool>(
"Build Fission "
37 ? model_params.get<
bool>(
"Build Fission "
40 , _build_transmutation_source(
41 model_params.isType<
bool>(
"Build Reaction Transmutation Source")
42 ? model_params.get<
bool>(
"Build Reaction Transmutation "
45 , _diffusion_coef(_build_diffusion
46 ? model_params.get<
double>(
"Diffusion "
48 : std::numeric_limits<double>::quiet_NaN())
49 , _decay(
"species_decay", _num_species, _num_species)
50 , _xs(_build_fission_source
51 ? reaction_params.get<
double>(
"Fission Cross-Section")
52 : std::numeric_limits<double>::quiet_NaN())
53 , _gamma(
"atoms_per_species", _num_species)
54 , _sigma(
"microscopic_cross_section", _num_species, _num_species)
60 = reaction_params.get<Teuchos::Array<double>>(
"Species Decay");
61 for (
int i = 0; i < _num_species; ++i)
63 for (
int j = 0; j < _num_species; ++j)
65 _decay(i, j) = decay[_num_species * i + j];
71 if (_build_transmutation_source)
73 const auto sigma = reaction_params.get<Teuchos::Array<double>>(
74 "Microscopic Cross-Section");
75 for (
int i = 0; i < _num_species; ++i)
77 for (
int j = 0; j < _num_species; ++j)
79 _sigma(i, j) = sigma[_num_species * i + j];
85 if (_build_fission_source)
87 auto gamma = reaction_params.get<Teuchos::Array<double>>(
88 "Number of atoms per species");
90 = Kokkos::View<double*, Kokkos::HostSpace, Kokkos::MemoryUnmanaged>(
91 gamma.data(), gamma.length());
92 Kokkos::deep_copy(_gamma, gamma_view);
97 KOKKOS_INLINE_FUNCTION
auto numSpecies()
const {
return _num_species; }
100 KOKKOS_INLINE_FUNCTION
auto batemanMatrix()
const {
return _decay; }
103 KOKKOS_INLINE_FUNCTION
bool buildBateman()
const {
return _build_bateman; }
106 KOKKOS_INLINE_FUNCTION
bool buildAdvection()
const
108 return _build_advection;
112 KOKKOS_INLINE_FUNCTION
bool buildDiffusion()
const
114 return _build_diffusion;
118 KOKKOS_INLINE_FUNCTION
bool buildFissionSource()
const
120 return _build_fission_source;
124 KOKKOS_INLINE_FUNCTION
bool buildTransmutationSource()
const
126 return _build_transmutation_source;
130 KOKKOS_INLINE_FUNCTION
double constantDiffusionCoefficient()
const
132 return _diffusion_coef;
136 KOKKOS_INLINE_FUNCTION
double fissionCrossSection()
const {
return _xs; }
139 KOKKOS_INLINE_FUNCTION
auto atomsPerSpecies()
const {
return _gamma; }
142 KOKKOS_INLINE_FUNCTION
auto microscopicCrossSection()
const
150 bool _build_advection;
151 bool _build_diffusion;
152 bool _build_fission_source;
153 bool _build_transmutation_source;
154 double _diffusion_coef;
155 Kokkos::View<double**, Kokkos::LayoutLeft, Kokkos::HostSpace> _decay;
157 Kokkos::View<double*, Kokkos::LayoutLeft, Kokkos::HostSpace> _gamma;
158 Kokkos::View<double**, Kokkos::LayoutLeft, Kokkos::HostSpace> _sigma;
164 #endif // VERTEXCFD_CONSTANTSPECIESPROPERTIES_HPP