1 #ifndef VERTEXCFD_CLOSURE_RADFISSIONSOURCEEXACTSOLUTION_IMPL_HPP
2 #define VERTEXCFD_CLOSURE_RADFISSIONSOURCEEXACTSOLUTION_IMPL_HPP
4 #include "utils/VertexCFD_Utils_ScalarFieldView.hpp"
6 #include <Panzer_GlobalIndexer.hpp>
7 #include <Panzer_HierarchicParallelism.hpp>
8 #include <Panzer_PureBasis.hpp>
9 #include <Panzer_Workset_Utilities.hpp>
11 #include <Teuchos_Array.hpp>
18 namespace ClosureModel
21 template<
class EvalType,
class Traits>
22 RADFissionSourceExactSolution<EvalType, Traits>::RADFissionSourceExactSolution(
23 const panzer::IntegrationRule& ir,
24 const SpeciesProperties::ConstantSpeciesProperties& species_prop,
25 const Teuchos::ParameterList& closure_params)
26 : _num_species(species_prop.numSpecies())
27 , _neutron_flux(
"neutron_flux", ir.dl_scalar)
28 , _xs(species_prop.fissionCrossSection())
29 , _avagadro(6.02214076e23)
30 , _gamma(Kokkos::ViewAllocateWithoutInitializing(
"atoms_per_species"),
34 auto gamma_host = Kokkos::create_mirror_view(Kokkos::HostSpace{}, _gamma);
35 gamma_host = species_prop.atomsPerSpecies();
36 Kokkos::deep_copy(_gamma, gamma_host);
38 _a = closure_params.get<
double>(
"Time Shape Function Ramp Scale");
39 _b = closure_params.get<
double>(
"Time Shape Function Ramp Offset");
40 _kappa = closure_params.get<
double>(
"Time Shape Function Decline Scale");
41 _beta = closure_params.get<
double>(
"Time Shape Function Decline Offset");
42 _flux_amp = closure_params.get<
double>(
"Neutron Flux Amplitude");
43 _calc_flux = closure_params.get<
bool>(
"Calculate Analytical Neutron Flux");
44 const auto initial_amp = closure_params.get<Teuchos::Array<double>>(
45 "Initial Species Concentration");
46 for (
int i = 0; i < _num_species; ++i)
48 _initial_amp[i] = initial_amp[i];
51 Utils::addEvaluatedScalarFieldView(
52 *
this, ir.dl_scalar, _num_species, _exact_species,
"Exact_species_");
55 this->addEvaluatedField(_neutron_flux);
57 this->setName(
"RAD Fission Source Exact Solution");
61 template<
class EvalType,
class Traits>
62 void RADFissionSourceExactSolution<EvalType, Traits>::evaluateFields(
63 typename Traits::EvalData workset)
67 if constexpr (Sacado::IsADType<scalar_type>::value)
70 = Kokkos::dimension_scalar(_exact_species[0].get_view());
71 bytes = scratch_view::shmem_size(
72 _exact_species[0].extent(1), NUM_TMPS, fad_size);
76 bytes = scratch_view::shmem_size(_exact_species[0].extent(1), NUM_TMPS);
79 auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
81 policy.set_scratch_size(0, Kokkos::PerTeam(bytes));
83 Kokkos::parallel_for(this->getName(), policy, *
this);
87 template<
class EvalType,
class Traits>
88 KOKKOS_INLINE_FUNCTION
void
89 RADFissionSourceExactSolution<EvalType, Traits>::operator()(
90 const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team)
const
92 const int cell = team.league_rank();
93 const int num_point = _exact_species[0].extent(1);
99 scratch_view tmp_field;
100 if constexpr (Sacado::IsADType<scalar_type>::value)
103 = Kokkos::dimension_scalar(_exact_species[0].get_view());
105 = scratch_view(team.team_shmem(), num_point, NUM_TMPS, fad_size);
109 tmp_field = scratch_view(team.team_shmem(), num_point, NUM_TMPS);
112 Kokkos::parallel_for(
113 Kokkos::TeamThreadRange(team, 0, num_point), [&](
const int point) {
114 auto&& const_s = tmp_field(point, CONSTS);
115 auto&& const_c = tmp_field(point, CONSTC);
117 _neutron_flux(cell, point)
119 * ((tanh(_a * _time - _b) + 1.0)
120 - (tanh(_kappa * _time - _beta) + 1.0));
121 for (
int num = 0; num < _num_species; ++num)
123 const_s = _flux_amp * _gamma[num] * _xs / _avagadro;
124 const_c = -const_s * log(cosh(_b)) / (2 * _a)
125 + const_s * log(cosh(-_beta)) / (2 * _kappa)
127 _exact_species[num](cell, point)
128 = const_s * log(cosh(_b - _a * _time)) / (2 * _a)
129 - const_s * log(cosh(_kappa * _time - _beta))
141 #endif // VERTEXCFD_CLOSURE_RADFISSIONSOURCEEXACTSOLUTION_IMPL_HPP