Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
MCGIDI_distributions.hpp
Go to the documentation of this file.
1/*
2# <<BEGIN-copyright>>
3# Copyright 2019, Lawrence Livermore National Security, LLC.
4# This file is part of the gidiplus package (https://github.com/LLNL/gidiplus).
5# gidiplus is licensed under the MIT license (see https://opensource.org/licenses/MIT).
6# SPDX-License-Identifier: MIT
7# <<END-copyright>>
8*/
9
10#ifndef MCGIDI_distributions_hpp_included
11#define MCGIDI_distributions_hpp_included 1
12
13#include <LUPI_declareMacro.hpp>
14
15namespace MCGIDI {
16
17namespace Distributions {
18
22
23/*
24============================================================
25======================= Distribution =======================
26============================================================
27*/
29
30 private:
31 Type m_type; /**< Specifies the Type of the distribution. */
32 GIDI::Frame m_productFrame; /**< Specifies the frame the product data are given in. */
33 double m_projectileMass; /**< The mass of the projectile. */
34 double m_targetMass; /**< The mass of the target. */
35 double m_productMass; /**< The mass of the first product. */
36
37 public:
39 LUPI_HOST Distribution( Type a_type, GIDI::Distributions::Distribution const &a_distribution, SetupInfo &a_setupInfo );
40 LUPI_HOST Distribution( Type a_type, GIDI::Frame a_productFrame, SetupInfo &a_setupInfo );
42
43 LUPI_HOST_DEVICE Type type( ) const { return( m_type ); } /**< Returns the value of the **m_type**. */
44 LUPI_HOST_DEVICE GIDI::Frame productFrame( ) const { return( m_productFrame ); } /**< Returns the value of the **m_productFrame**. */
45
46 LUPI_HOST_DEVICE double projectileMass( ) const { return( m_projectileMass ); } /**< Returns the value of the **m_projectileMass**. */
47 LUPI_HOST_DEVICE double targetMass( ) const { return( m_targetMass ); } /**< Returns the value of the **m_targetMass**. */
48 LUPI_HOST_DEVICE double productMass( ) const { return( m_productMass ); } /**< Returns the value of the **m_productMass**. */
49
51
52 template <typename RNG>
53 LUPI_HOST_DEVICE MCGIDI_VIRTUAL_FUNCTION void sample( double a_X, Sampling::Input &a_input, RNG && a_rng ) const MCGIDI_TRUE_VIRTUAL;
54 template <typename RNG>
55 LUPI_HOST_DEVICE MCGIDI_VIRTUAL_FUNCTION double angleBiasing( Reaction const *a_reaction, double a_temperature, double a_energy_in, double a_mu_lab,
56 RNG && a_rng, double &a_energy_out ) const MCGIDI_TRUE_VIRTUAL;
58
59};
60
61/*
62============================================================
63====================== AngularTwoBody ======================
64============================================================
65*/
67
68 private:
69 double m_residualMass; /**< The mass of the second product (often the residual). */
70 double m_Q; /**< FIX ME. */
71 double m_twoBodyThreshold; /**< This is the T_1 value needed to do two-body kinematics (i.e., in the equation (K_{com,3_4} = m_2 * (K_1 - T_1) / (m_1 + m_2)). */
72 bool m_Upscatter; /**< Set to true if reaction is elastic which is the only reaction upscatter Model B is applied to. */
73 Probabilities::ProbabilityBase2d_d1 *m_angular; /**< The 2d angular probability. */
74 Sampling::Upscatter::ModelDBRC_data *m_modelDBRC_data; /**< The cross section and other data needed for neutron elastic upscatter model DBRC. */
75
76 template <typename RNG>
77 LUPI_HOST_DEVICE bool upscatterModelB( double a_kineticLab, Sampling::Input &a_input, RNG && a_rng ) const ;
78
79 public:
81 LUPI_HOST AngularTwoBody( GIDI::Distributions::AngularTwoBody const &a_angularTwoBody, SetupInfo &a_setupInfo );
83
84 LUPI_HOST_DEVICE double residualMass( ) const { return( m_residualMass ); } /**< Returns the value of the **m_residualMass**. */
85 LUPI_HOST_DEVICE double Q( ) const { return( m_Q ); } /**< Returns the value of the **m_Q**. */
86 LUPI_HOST_DEVICE Probabilities::ProbabilityBase2d_d1 *angular( ) const { return( m_angular ); } /**< Returns the value of the **m_angular**. */
87 template <typename RNG>
88 LUPI_HOST_DEVICE void sample( double a_X, Sampling::Input &a_input, RNG && a_rng ) const ;
89 template <typename RNG>
90 LUPI_HOST_DEVICE double angleBiasing( Reaction const *a_reaction, double a_temperature, double a_energy_in, double a_mu_lab,
91 RNG && a_rng, double &a_energy_out ) const ;
93 LUPI_HOST_DEVICE bool Upscatter( ) const { return( m_Upscatter ); } /**< Returns the value of the **m_Upscatter**. */
95};
96
97/*
98============================================================
99======================= Uncorrelated =======================
100============================================================
101*/
103
104 private:
105 Probabilities::ProbabilityBase2d_d1 *m_angular; /**< The angular probability P(mu|E). */
106 Probabilities::ProbabilityBase2d *m_energy; /**< The energy probability P(E'|E). */
107
108 public:
110 LUPI_HOST Uncorrelated( GIDI::Distributions::Uncorrelated const &a_uncorrelated, SetupInfo &a_setupInfo );
112
113 LUPI_HOST_DEVICE Probabilities::ProbabilityBase2d_d1 *angular( ) const { return( m_angular ); } /**< Returns the value of the **m_angular**. */
114 LUPI_HOST_DEVICE Probabilities::ProbabilityBase2d *energy( ) const { return( m_energy ); } /**< Returns the value of the **m_energy**. */
115 template <typename RNG>
116 LUPI_HOST_DEVICE void sample( double a_X, Sampling::Input &a_input, RNG && a_rng ) const ;
117 template <typename RNG>
118 LUPI_HOST_DEVICE double angleBiasing( Reaction const *a_reaction, double a_temperature, double a_energy_in, double a_mu_lab,
119 RNG && a_rng, double &a_energy_out ) const ;
121};
122
123/*
124============================================================
125======================== Branching3d =======================
126============================================================
127*/
128class Branching3d : public Distribution {
129
130 private:
131 int m_initialStateIndex;
132
133 public:
135 LUPI_HOST Branching3d( GIDI::Distributions::Branching3d const &a_branching3d, SetupInfo &a_setupInfo );
137
138 template <typename RNG>
139 LUPI_HOST_DEVICE void sample( double a_X, Sampling::Input &a_input, RNG && a_rng ) const ;
140 template <typename RNG>
141 LUPI_HOST_DEVICE double angleBiasing( Reaction const *a_reaction, double a_temperature, double a_energy_in, double a_mu_lab,
142 RNG && a_rng, double &a_energy_out ) const ;
144};
145
146/*
147============================================================
148====================== EnergyAngularMC =====================
149============================================================
150*/
152
153 private:
154 Probabilities::ProbabilityBase2d_d1 *m_energy; /**< The energy probability P(E'|E). */
155 Probabilities::ProbabilityBase3d *m_angularGivenEnergy; /**< The angular probability given E', P(mu|E,E'). */
156
157 public:
159 LUPI_HOST EnergyAngularMC( GIDI::Distributions::EnergyAngularMC const &a_energyAngularMC, SetupInfo &a_setupInfo );
161
162 LUPI_HOST_DEVICE Probabilities::ProbabilityBase2d_d1 *energy( ) const { return( m_energy ); } /**< Returns the value of the **m_energy**. */
163 LUPI_HOST_DEVICE Probabilities::ProbabilityBase3d *angularGivenEnergy( ) const { return( m_angularGivenEnergy ); } /**< Returns the value of the **m_angularGivenEnergy**. */
164 template <typename RNG>
165 LUPI_HOST_DEVICE void sample( double a_X, Sampling::Input &a_input, RNG && a_rng ) const ;
166 template <typename RNG>
167 LUPI_HOST_DEVICE double angleBiasing( Reaction const *a_reaction, double a_temperature, double a_energy_in, double a_mu_lab,
168 RNG && a_rng, double &a_energy_out ) const ;
170};
171
172/*
173============================================================
174====================== AngularEnergyMC =====================
175============================================================
176*/
178
179 private:
180 Probabilities::ProbabilityBase2d_d1 *m_angular; /**< The angular probability P(mu|E). */
181 Probabilities::ProbabilityBase3d *m_energyGivenAngular; /**< The energy probability P(E'|E,mu). */
182
183 public:
185 LUPI_HOST AngularEnergyMC( GIDI::Distributions::AngularEnergyMC const &a_angularEnergyMC, SetupInfo &a_setupInfo );
187
188 LUPI_HOST_DEVICE Probabilities::ProbabilityBase2d_d1 *angular( ) const { return( m_angular ); } /**< Returns the value of the **m_angular**. */
189 LUPI_HOST_DEVICE Probabilities::ProbabilityBase3d *energyGivenAngular( ) const { return( m_energyGivenAngular ); } /**< Returns the value of the **m_energyGivenAngular**. */
190 template <typename RNG>
191 LUPI_HOST_DEVICE void sample( double a_X, Sampling::Input &a_input, RNG && a_rng ) const ;
192 template <typename RNG>
193 LUPI_HOST_DEVICE double angleBiasing( Reaction const *a_reaction, double a_temperature, double a_energy_in, double a_mu_lab,
194 RNG && a_rng, double &a_energy_out ) const ;
196};
197
198/*
199============================================================
200======================== KalbachMann =======================
201============================================================
202*/
203class KalbachMann : public Distribution {
204
205 private:
206 double m_energyToMeVFactor; /**< The factor that converts energies to MeV. */
207 double m_eb_massFactor; /**< FIX ME */
208 Probabilities::ProbabilityBase2d_d1 *m_f; /**< The energy probability P(E'|E). */
209 Functions::Function2d *m_r; /**< The Kalbach-Mann r(E,E') function. */
210 Functions::Function2d *m_a; /**< The Kalbach-Mann a(E,E') function. */
211
212 public:
214 LUPI_HOST KalbachMann( GIDI::Distributions::KalbachMann const &a_KalbachMann, SetupInfo &a_setupInfo );
216
217 LUPI_HOST_DEVICE double energyToMeVFactor( ) const { return( m_energyToMeVFactor ); } /**< Returns the value of the **m_energyToMeVFactor**. */
218 LUPI_HOST_DEVICE double eb_massFactor( ) const { return( m_eb_massFactor ); } /**< Returns the value of the **m_eb_massFactor**. */
219 LUPI_HOST_DEVICE Probabilities::ProbabilityBase2d_d1 *f( ) const { return( m_f ); } /**< Returns the value of the **m_f**. */
220 LUPI_HOST_DEVICE Functions::Function2d *r( ) const { return( m_r ); } /**< Returns the value of the **m_r**. */
221 LUPI_HOST_DEVICE Functions::Function2d *a( ) const { return( m_a ); } /**< Returns the value of the **m_a**. */
222 template <typename RNG>
223 LUPI_HOST_DEVICE void sample( double a_X, Sampling::Input &a_input, RNG && a_rng ) const ;
224 template <typename RNG>
225 LUPI_HOST_DEVICE double angleBiasing( Reaction const *a_reaction, double a_temperature, double a_energy_in, double a_mu_lab,
226 RNG && a_rng, double &a_energy_out ) const ;
228
229 LUPI_HOST_DEVICE double evaluate( double E_in_lab, double E_out, double mu );
230};
231
232/*
233============================================================
234=============== CoherentPhotoAtomicScattering ==============
235============================================================
236*/
238
239 private:
240 bool m_anomalousDataPresent; /**< FIX ME */
241 Vector<double> m_energies; /**< FIX ME */
242 Vector<double> m_formFactor; /**< FIX ME */
243 Vector<double> m_a; /**< FIX ME */
244 Vector<double> m_integratedFormFactor; /**< FIX ME */
245 Vector<double> m_integratedFormFactorSquared; /**< FIX ME */
246 Vector<double> m_probabilityNorm1_1; /**< FIX ME */
247 Vector<double> m_probabilityNorm1_3; /**< FIX ME */
248 Vector<double> m_probabilityNorm1_5; /**< FIX ME */
249 Vector<double> m_probabilityNorm2_1; /**< FIX ME */
250 Vector<double> m_probabilityNorm2_3; /**< FIX ME */
251 Vector<double> m_probabilityNorm2_5; /**< FIX ME */
252 Functions::Function1d_d1 *m_realAnomalousFactor; /**< The real part of the anomalous scattering factor. */
253 Functions::Function1d_d1 *m_imaginaryAnomalousFactor; /**< The imaginary part of the anomalous scattering factor. */
254
255 LUPI_HOST_DEVICE double Z_a( double a_Z, double a_a ) const ;
256
257 public:
259 LUPI_HOST CoherentPhotoAtomicScattering( GIDI::Distributions::CoherentPhotoAtomicScattering const &a_coherentPhotoAtomicScattering, SetupInfo &a_setupInfo );
261
262 LUPI_HOST_DEVICE double evaluate( double a_energyIn, double a_mu ) const ;
263 LUPI_HOST_DEVICE double evaluateFormFactor( double a_energyIn, double a_mu ) const ;
264 template <typename RNG>
265 LUPI_HOST_DEVICE void sample( double a_X, Sampling::Input &a_input, RNG && a_rng ) const ;
266 template <typename RNG>
267 LUPI_HOST_DEVICE double angleBiasing( Reaction const *a_reaction, double a_temperature, double a_energy_in, double a_mu_lab,
268 RNG && a_rng, double &a_energy_out ) const ;
269
271};
272
273/*
274============================================================
275============== IncoherentPhotoAtomicScattering =============
276============================================================
277*/
279
280 private:
281 Vector<double> m_energies; /**< FIX ME */
282 Vector<double> m_scatteringFactor; /**< FIX ME */
283 Vector<double> m_a; /**< FIX ME */
284
285 public:
289
290 LUPI_HOST_DEVICE double energyRatio( double a_energyIn, double a_mu ) const ;
291 LUPI_HOST_DEVICE double evaluateKleinNishina( double a_energyIn, double a_mu ) const ;
292 LUPI_HOST_DEVICE double evaluateScatteringFactor( double a_X ) const ;
293 template <typename RNG>
294 LUPI_HOST_DEVICE void sample( double a_X, Sampling::Input &a_input, RNG && a_rng ) const ;
295 template <typename RNG>
296 LUPI_HOST_DEVICE double angleBiasing( Reaction const *a_reaction, double a_temperature, double a_energy_in, double a_mu_lab,
297 RNG && a_rng, double &a_energy_out ) const ;
299/*
300 LUPI_HOST_DEVICE double evaluate( double E_in_lab, double mu );
301*/
302};
303
304/*
305=======================================================================
306============== IncoherentBoundToFreePhotoAtomicScattering =============
307=======================================================================
308*/
310
311 private:
312 //Vector<double> m_energies;
313 //Vector<double> m_ComptonProfile;
314 Vector<double> m_occupationNumber;
315 //Vector<double> m_a;
316 Vector<double> m_pz;
317 double m_bindingEnergy;
318
319 public:
323 LUPI_HOST_DEVICE double energyRatio( double a_energyIn, double a_mu ) const ;
324 LUPI_HOST_DEVICE double evaluateKleinNishina( double a_energyIn, double a_mu ) const ;
325 LUPI_HOST_DEVICE double evaluateOccupationNumber( double a_X, double a_mu ) const;
326 template <typename RNG>
327 LUPI_HOST_DEVICE void sample( double a_X, Sampling::Input &a_input, RNG && a_rng ) const ;
328 template <typename RNG>
329 LUPI_HOST_DEVICE double angleBiasing( Reaction const *a_reaction, double a_temperature, double a_energy_in, double a_mu_lab,
330 RNG && a_rng, double &a_energy_out ) const ;
332};
333
334/*
335============================================================
336========== IncoherentPhotoAtomicScatteringElectron =========
337============================================================
338*/
340
341 public:
345
346 template <typename RNG>
347 LUPI_HOST_DEVICE void sample( double a_energy, Sampling::Input &a_input, RNG && a_rng ) const ;
348 template <typename RNG>
349 LUPI_HOST_DEVICE double angleBiasing( Reaction const *a_reaction, double a_temperature, double a_energy_in, double a_mu_lab,
350 RNG && a_rng, double &a_energy_out ) const ;
352};
353
354/*
355============================================================
356==================== PairProductionGamma ===================
357============================================================
358*/
360
361 private:
362 bool m_firstSampled; /**< When sampling photons for pair production, the photons must be emitted back-to-back. The flag help do this. */
363
364 public:
366 LUPI_HOST PairProductionGamma( SetupInfo &a_setupInfo, bool a_firstSampled );
368
369 template <typename RNG>
370 LUPI_HOST_DEVICE void sample( double a_X, Sampling::Input &a_input, RNG && a_rng ) const ;
371 template <typename RNG>
372 LUPI_HOST_DEVICE double angleBiasing( Reaction const *a_reaction, double a_temperature, double a_energy_in, double a_mu_lab,
373 RNG && a_rng, double &a_energy_out ) const ;
375};
376
377/*
378============================================================
379==================== CoherentElasticTNSL ===================
380============================================================
381*/
383
384 private:
385 Interpolation m_temperatureInterpolation;
386 Vector<double> m_temperatures;
387 Vector<double> m_energies;
388 Vector<double> m_S_table;
389
390 public:
393 SetupInfo &a_setupInfo );
395
396 template <typename RNG>
397 LUPI_HOST_DEVICE void sample( double a_energy, Sampling::Input &a_input, RNG && a_rng ) const ;
398 template <typename RNG>
399 LUPI_HOST_DEVICE double angleBiasing( Reaction const *a_reaction, double a_temperature, double a_energy_in, double a_mu_lab,
400 RNG && a_rng, double &a_energy_out ) const ;
402};
403
404/*
405============================================================
406==================== IncoherentElasticTNSL ===================
407============================================================
408*/
410
411 private:
412 double m_temperatureToMeV_K;
413 Functions::Function1d_d1 *m_DebyeWallerIntegral;
414
415 public:
418 SetupInfo &a_setupInfo );
420
421 template <typename RNG>
422 LUPI_HOST_DEVICE void sample( double a_energy, Sampling::Input &a_input, RNG && a_rng ) const ;
423 template <typename RNG>
424 LUPI_HOST_DEVICE double angleBiasing( Reaction const *a_reaction, double a_temperature, double a_energy_in, double a_mu_lab,
425 RNG && a_rng, double &a_energy_out ) const ;
427
428 Functions::Function1d *DebyeWallerIntegral( ) { return( m_DebyeWallerIntegral ); }
429 Functions::Function1d const *DebyeWallerIntegral( ) const { return( m_DebyeWallerIntegral ); }
430};
431
432/*
433============================================================
434======================= Unspecified ========================
435============================================================
436*/
437class Unspecified : public Distribution {
438
439 public:
441 LUPI_HOST Unspecified( GIDI::Distributions::Distribution const &a_distribution, SetupInfo &a_setupInfo );
443
444 template <typename RNG>
445 LUPI_HOST_DEVICE void sample( double a_X, Sampling::Input &a_input, RNG && a_rng ) const ;
446 template <typename RNG>
447 LUPI_HOST_DEVICE double angleBiasing( Reaction const *a_reaction, double a_temperature, double a_energy_in, double a_mu_lab,
448 RNG && a_rng, double &a_energy_out ) const ;
450};
451
452/*
453============================================================
454========================== Others ==========================
455============================================================
456*/
457LUPI_HOST Distribution *parseGIDI( GIDI::Suite const &a_distribution, SetupInfo &a_setupInfo, Transporting::MC const &a_settings );
458LUPI_HOST_DEVICE Type DistributionType( Distribution const *a_distribution );
459
460}
461
462}
463
464#endif // End of MCGIDI_distributions_hpp_included
#define LUPI_HOST_DEVICE
#define LUPI_HOST
#define MCGIDI_VIRTUAL_FUNCTION
Definition MCGIDI.hpp:33
#define MCGIDI_TRUE_VIRTUAL
Definition MCGIDI.hpp:34
LUPI_HOST_DEVICE Probabilities::ProbabilityBase3d * energyGivenAngular() const
LUPI_HOST_DEVICE double angleBiasing(Reaction const *a_reaction, double a_temperature, double a_energy_in, double a_mu_lab, RNG &&a_rng, double &a_energy_out) const
LUPI_HOST_DEVICE void serialize(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
LUPI_HOST_DEVICE void sample(double a_X, Sampling::Input &a_input, RNG &&a_rng) const
LUPI_HOST_DEVICE Probabilities::ProbabilityBase2d_d1 * angular() const
LUPI_HOST_DEVICE void sample(double a_X, Sampling::Input &a_input, RNG &&a_rng) const
LUPI_HOST_DEVICE double residualMass() const
LUPI_HOST_DEVICE bool Upscatter() const
LUPI_HOST_DEVICE void serialize(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
LUPI_HOST_DEVICE double angleBiasing(Reaction const *a_reaction, double a_temperature, double a_energy_in, double a_mu_lab, RNG &&a_rng, double &a_energy_out) const
LUPI_HOST void setModelDBRC_data2(Sampling::Upscatter::ModelDBRC_data *a_modelDBRC_data)
LUPI_HOST_DEVICE Probabilities::ProbabilityBase2d_d1 * angular() const
LUPI_HOST_DEVICE void sample(double a_X, Sampling::Input &a_input, RNG &&a_rng) const
LUPI_HOST_DEVICE double angleBiasing(Reaction const *a_reaction, double a_temperature, double a_energy_in, double a_mu_lab, RNG &&a_rng, double &a_energy_out) const
LUPI_HOST_DEVICE void serialize(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
LUPI_HOST_DEVICE void serialize(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
LUPI_HOST_DEVICE double angleBiasing(Reaction const *a_reaction, double a_temperature, double a_energy_in, double a_mu_lab, RNG &&a_rng, double &a_energy_out) const
LUPI_HOST_DEVICE void sample(double a_energy, Sampling::Input &a_input, RNG &&a_rng) const
LUPI_HOST_DEVICE double evaluate(double a_energyIn, double a_mu) const
LUPI_HOST_DEVICE void sample(double a_X, Sampling::Input &a_input, RNG &&a_rng) const
LUPI_HOST_DEVICE double angleBiasing(Reaction const *a_reaction, double a_temperature, double a_energy_in, double a_mu_lab, RNG &&a_rng, double &a_energy_out) const
LUPI_HOST_DEVICE void serialize(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
LUPI_HOST_DEVICE double evaluateFormFactor(double a_energyIn, double a_mu) const
LUPI_HOST_DEVICE MCGIDI_VIRTUAL_FUNCTION double angleBiasing(Reaction const *a_reaction, double a_temperature, double a_energy_in, double a_mu_lab, RNG &&a_rng, double &a_energy_out) const MCGIDI_TRUE_VIRTUAL
LUPI_HOST_DEVICE void serialize(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
LUPI_HOST_DEVICE MCGIDI_VIRTUAL_FUNCTION ~Distribution() MCGIDI_TRUE_VIRTUAL
LUPI_HOST_DEVICE double targetMass() const
LUPI_HOST void setModelDBRC_data(Sampling::Upscatter::ModelDBRC_data *a_modelDBRC_data)
LUPI_HOST_DEVICE MCGIDI_VIRTUAL_FUNCTION void sample(double a_X, Sampling::Input &a_input, RNG &&a_rng) const MCGIDI_TRUE_VIRTUAL
LUPI_HOST_DEVICE GIDI::Frame productFrame() const
LUPI_HOST_DEVICE double productMass() const
LUPI_HOST_DEVICE double projectileMass() const
LUPI_HOST_DEVICE Probabilities::ProbabilityBase2d_d1 * energy() const
LUPI_HOST_DEVICE Probabilities::ProbabilityBase3d * angularGivenEnergy() const
LUPI_HOST_DEVICE void sample(double a_X, Sampling::Input &a_input, RNG &&a_rng) const
LUPI_HOST_DEVICE double angleBiasing(Reaction const *a_reaction, double a_temperature, double a_energy_in, double a_mu_lab, RNG &&a_rng, double &a_energy_out) const
LUPI_HOST_DEVICE void serialize(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
LUPI_HOST_DEVICE void sample(double a_X, Sampling::Input &a_input, RNG &&a_rng) const
LUPI_HOST_DEVICE double energyRatio(double a_energyIn, double a_mu) const
LUPI_HOST_DEVICE void serialize(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
LUPI_HOST_DEVICE double evaluateKleinNishina(double a_energyIn, double a_mu) const
LUPI_HOST_DEVICE double angleBiasing(Reaction const *a_reaction, double a_temperature, double a_energy_in, double a_mu_lab, RNG &&a_rng, double &a_energy_out) const
LUPI_HOST_DEVICE double evaluateOccupationNumber(double a_X, double a_mu) const
Functions::Function1d const * DebyeWallerIntegral() const
LUPI_HOST_DEVICE void serialize(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
LUPI_HOST_DEVICE double angleBiasing(Reaction const *a_reaction, double a_temperature, double a_energy_in, double a_mu_lab, RNG &&a_rng, double &a_energy_out) const
LUPI_HOST_DEVICE void sample(double a_energy, Sampling::Input &a_input, RNG &&a_rng) const
LUPI_HOST_DEVICE double angleBiasing(Reaction const *a_reaction, double a_temperature, double a_energy_in, double a_mu_lab, RNG &&a_rng, double &a_energy_out) const
LUPI_HOST_DEVICE void serialize(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
LUPI_HOST_DEVICE void sample(double a_energy, Sampling::Input &a_input, RNG &&a_rng) const
LUPI_HOST_DEVICE void sample(double a_X, Sampling::Input &a_input, RNG &&a_rng) const
LUPI_HOST_DEVICE double evaluateKleinNishina(double a_energyIn, double a_mu) const
LUPI_HOST_DEVICE double angleBiasing(Reaction const *a_reaction, double a_temperature, double a_energy_in, double a_mu_lab, RNG &&a_rng, double &a_energy_out) const
LUPI_HOST_DEVICE void serialize(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
LUPI_HOST_DEVICE double energyRatio(double a_energyIn, double a_mu) const
LUPI_HOST_DEVICE double evaluateScatteringFactor(double a_X) const
LUPI_HOST_DEVICE Functions::Function2d * a() const
LUPI_HOST_DEVICE double eb_massFactor() const
LUPI_HOST_DEVICE double angleBiasing(Reaction const *a_reaction, double a_temperature, double a_energy_in, double a_mu_lab, RNG &&a_rng, double &a_energy_out) const
LUPI_HOST_DEVICE double energyToMeVFactor() const
LUPI_HOST_DEVICE double evaluate(double E_in_lab, double E_out, double mu)
LUPI_HOST_DEVICE Functions::Function2d * r() const
LUPI_HOST_DEVICE void serialize(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
LUPI_HOST_DEVICE Probabilities::ProbabilityBase2d_d1 * f() const
LUPI_HOST_DEVICE void sample(double a_X, Sampling::Input &a_input, RNG &&a_rng) const
LUPI_HOST_DEVICE double angleBiasing(Reaction const *a_reaction, double a_temperature, double a_energy_in, double a_mu_lab, RNG &&a_rng, double &a_energy_out) const
LUPI_HOST_DEVICE void sample(double a_X, Sampling::Input &a_input, RNG &&a_rng) const
LUPI_HOST_DEVICE void serialize(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
LUPI_HOST_DEVICE double angleBiasing(Reaction const *a_reaction, double a_temperature, double a_energy_in, double a_mu_lab, RNG &&a_rng, double &a_energy_out) const
LUPI_HOST_DEVICE Probabilities::ProbabilityBase2d_d1 * angular() const
LUPI_HOST_DEVICE void serialize(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
LUPI_HOST_DEVICE Probabilities::ProbabilityBase2d * energy() const
LUPI_HOST_DEVICE void sample(double a_X, Sampling::Input &a_input, RNG &&a_rng) const
LUPI_HOST_DEVICE double angleBiasing(Reaction const *a_reaction, double a_temperature, double a_energy_in, double a_mu_lab, RNG &&a_rng, double &a_energy_out) const
LUPI_HOST_DEVICE void serialize(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
LUPI_HOST_DEVICE void sample(double a_X, Sampling::Input &a_input, RNG &&a_rng) const
Frame
Definition GIDI.hpp:146
LUPI_HOST_DEVICE Type DistributionType(Distribution const *a_distribution)
LUPI_HOST Distribution * parseGIDI(GIDI::Suite const &a_distribution, SetupInfo &a_setupInfo, Transporting::MC const &a_settings)
Simple C++ string class, useful as replacement for std::string if this cannot be used,...
Definition MCGIDI.hpp:43