VertexCFD  0.0-dev
VertexCFD_ConstantSpeciesProperties.hpp
1 #ifndef VERTEXCFD_CONSTANTSPECIESPROPERTIES_HPP
2 #define VERTEXCFD_CONSTANTSPECIESPROPERTIES_HPP
3 
4 #include <Teuchos_ParameterList.hpp>
5 
6 #include <Kokkos_Core.hpp>
7 
8 namespace VertexCFD
9 {
10 namespace SpeciesProperties
11 {
12 //---------------------------------------------------------------------------//
13 // Constant species properties
14 //---------------------------------------------------------------------------//
15 
17 {
18  public:
19  ConstantSpeciesProperties() = default;
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 "
25  "Source")
26  ? model_params.get<bool>("Build Reaction Bateman "
27  "Source")
28  : false)
29  , _build_advection(model_params.isType<bool>("Build Advection")
30  ? model_params.get<bool>("Build Advection")
31  : false)
32  , _build_diffusion(model_params.isType<bool>("Build Diffusion")
33  ? model_params.get<bool>("Build Diffusion")
34  : false)
35  , _build_fission_source(model_params.isType<bool>("Build Fission "
36  "Source")
37  ? model_params.get<bool>("Build Fission "
38  "Source")
39  : false)
40  , _build_transmutation_source(
41  model_params.isType<bool>("Build Reaction Transmutation Source")
42  ? model_params.get<bool>("Build Reaction Transmutation "
43  "Source")
44  : false)
45  , _diffusion_coef(_build_diffusion
46  ? model_params.get<double>("Diffusion "
47  "Coefficient")
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)
55  {
56  // Bateman matrix
57  if (_build_bateman)
58  {
59  const auto decay
60  = reaction_params.get<Teuchos::Array<double>>("Species Decay");
61  for (int i = 0; i < _num_species; ++i)
62  {
63  for (int j = 0; j < _num_species; ++j)
64  {
65  _decay(i, j) = decay[_num_species * i + j];
66  }
67  }
68  }
69 
70  // Microscopic Cross-Section
71  if (_build_transmutation_source)
72  {
73  const auto sigma = reaction_params.get<Teuchos::Array<double>>(
74  "Microscopic Cross-Section");
75  for (int i = 0; i < _num_species; ++i)
76  {
77  for (int j = 0; j < _num_species; ++j)
78  {
79  _sigma(i, j) = sigma[_num_species * i + j];
80  }
81  }
82  }
83 
84  // Atoms per species
85  if (_build_fission_source)
86  {
87  auto gamma = reaction_params.get<Teuchos::Array<double>>(
88  "Number of atoms per species");
89  const auto gamma_view
90  = Kokkos::View<double*, Kokkos::HostSpace, Kokkos::MemoryUnmanaged>(
91  gamma.data(), gamma.length());
92  Kokkos::deep_copy(_gamma, gamma_view);
93  }
94  }
95 
96  // Bateman matrix
97  KOKKOS_INLINE_FUNCTION auto numSpecies() const { return _num_species; }
98 
99  // Bateman matrix
100  KOKKOS_INLINE_FUNCTION auto batemanMatrix() const { return _decay; }
101 
102  // Include reaction term
103  KOKKOS_INLINE_FUNCTION bool buildBateman() const { return _build_bateman; }
104 
105  // Include advection term
106  KOKKOS_INLINE_FUNCTION bool buildAdvection() const
107  {
108  return _build_advection;
109  }
110 
111  // Include diffusion term
112  KOKKOS_INLINE_FUNCTION bool buildDiffusion() const
113  {
114  return _build_diffusion;
115  }
116 
117  // Include fission source term
118  KOKKOS_INLINE_FUNCTION bool buildFissionSource() const
119  {
120  return _build_fission_source;
121  }
122 
123  // Include transmutation term
124  KOKKOS_INLINE_FUNCTION bool buildTransmutationSource() const
125  {
126  return _build_transmutation_source;
127  }
128 
129  // Constant diffusion coefficient
130  KOKKOS_INLINE_FUNCTION double constantDiffusionCoefficient() const
131  {
132  return _diffusion_coef;
133  }
134 
135  // Fission cross section
136  KOKKOS_INLINE_FUNCTION double fissionCrossSection() const { return _xs; }
137 
138  // Number of atoms per species
139  KOKKOS_INLINE_FUNCTION auto atomsPerSpecies() const { return _gamma; }
140 
141  // Microscopic cross section
142  KOKKOS_INLINE_FUNCTION auto microscopicCrossSection() const
143  {
144  return _sigma;
145  }
146 
147  private:
148  int _num_species;
149  bool _build_bateman;
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;
156  double _xs;
157  Kokkos::View<double*, Kokkos::LayoutLeft, Kokkos::HostSpace> _gamma;
158  Kokkos::View<double**, Kokkos::LayoutLeft, Kokkos::HostSpace> _sigma;
159 };
160 
161 } // namespace SpeciesProperties
162 } // namespace VertexCFD
163 
164 #endif // VERTEXCFD_CONSTANTSPECIESPROPERTIES_HPP
VertexCFD
Definition: tstMethodManufacturedSolutionBC.cpp:23
VertexCFD::SpeciesProperties::ConstantSpeciesProperties
Definition: VertexCFD_ConstantSpeciesProperties.hpp:17