Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
MCGIDI::OutputChannel Class Reference

#include <MCGIDI.hpp>

Public Member Functions

LUPI_HOST_DEVICE OutputChannel ()
LUPI_HOST OutputChannel (GIDI::OutputChannel const *a_outputChannel, SetupInfo &a_setupInfo, Transporting::MC const &a_settings, GIDI::Transporting::Particles const &a_particles)
LUPI_HOST_DEVICE ~OutputChannel ()
LUPI_HOST_DEVICE Productoperator[] (std::size_t a_index)
LUPI_HOST_DEVICE bool isTwoBody () const
LUPI_HOST_DEVICE double finalQ (double a_x1) const
LUPI_HOST_DEVICE bool isFission () const
LUPI_HOST_DEVICE bool hasFission () const
LUPI_HOST_DEVICE Functions::Function1d_d1Q ()
LUPI_HOST void setUserParticleIndex (int a_particleIndex, int a_userParticleIndex)
LUPI_HOST void setUserParticleIndexViaIntid (int a_particleIntid, int a_userParticleIndex)
LUPI_HOST void setModelDBRC_data (Sampling::Upscatter::ModelDBRC_data *a_modelDBRC_data)
LUPI_HOST_DEVICE Vector< Product * > const & products () const
Vector< DelayedNeutron * > delayedNeutrons () const
LUPI_HOST_DEVICE DelayedNeutron const * delayedNeutron (std::size_t a_index) const
LUPI_HOST void moveProductsEtAlToReaction (std::vector< Product * > &a_products, Functions::Function1d **a_totalDelayedNeutronMultiplicity, std::vector< DelayedNeutron * > &a_delayedNeutrons, std::vector< Functions::Function1d_d1 * > &a_Qs)
LUPI_HOST_DEVICE double productAverageMultiplicity (int a_index, double a_projectileEnergy) const
LUPI_HOST_DEVICE double productAverageMultiplicityViaIntid (int a_intid, double a_projectileEnergy) const
template<typename RNG, typename PUSHBACK>
LUPI_HOST_DEVICE void sampleProducts (ProtareSingle const *a_protare, double a_projectileEnergy, Sampling::Input &a_input, RNG &&a_rng, PUSHBACK &&a_push_back, Sampling::ProductHandler &a_products) const
template<typename RNG>
LUPI_HOST_DEVICE void angleBiasing (Reaction const *a_reaction, int a_pid, double a_temperature, double a_energy_in, double a_mu_lab, double &a_probability, double &a_energy_out, RNG &&a_rng, double &a_cumulative_weight) const
template<typename RNG>
LUPI_HOST_DEVICE void angleBiasingViaIntid (Reaction const *a_reaction, int a_intid, double a_temperature, double a_energy_in, double a_mu_lab, double &a_probability, double &a_energy_out, RNG &&a_rng, double &a_cumulative_weight) const
LUPI_HOST_DEVICE void serialize (LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)

Detailed Description

Definition at line 1272 of file MCGIDI.hpp.

Constructor & Destructor Documentation

◆ OutputChannel() [1/2]

LUPI_HOST_DEVICE MCGIDI::OutputChannel::OutputChannel ( )

Default constructor used when broadcasting a Protare as needed by MPI or GPUs.

Definition at line 23 of file MCGIDI_outputChannel.cc.

23 :
24 m_channelType( ChannelType::none ),
25 m_neutronIndex( -1 ),
26 m_isFission( false ),
27 m_hasFinalStatePhotons( false ),
28 m_Q( nullptr ),
29 m_products( ),
30 m_totalDelayedNeutronMultiplicity( nullptr ) {
31
32}

◆ OutputChannel() [2/2]

LUPI_HOST MCGIDI::OutputChannel::OutputChannel ( GIDI::OutputChannel const * a_outputChannel,
SetupInfo & a_setupInfo,
Transporting::MC const & a_settings,
GIDI::Transporting::Particles const & a_particles )
Parameters
a_outputChannel[in] The GIDI::OutputChannel whose data is to be used to construct this.
a_setupInfo[in] Used internally when constructing a Protare to pass information to other constructors.
a_settings[in] Used to pass user options to the this to instruct it which data are desired.
a_particles[in] List of transporting particles and their information (e.g., multi-group boundaries and fluxes).

Definition at line 41 of file MCGIDI_outputChannel.cc.

42 :
43 m_channelType( ChannelType::none ),
44 m_neutronIndex( a_setupInfo.m_neutronIndex ),
45 m_isFission( false ),
46 m_hasFinalStatePhotons( false ),
47 m_Q( nullptr ),
48 m_products( ),
49 m_totalDelayedNeutronMultiplicity( nullptr ) {
50
51 if( a_outputChannel != nullptr ) {
52 m_channelType = a_outputChannel->twoBody( ) ? ChannelType::twoBody : ChannelType::uncorrelatedBodies;
53 m_isFission = a_outputChannel->isFission( );
54
55 m_Q = Functions::parseFunction1d_d1( a_outputChannel->Q( ).get<GIDI::Functions::Function1dForm>( 0 ) );
56 if( a_setupInfo.m_isPairProduction ) {
57 double domainMin = m_Q->domainMin( ), domainMax = m_Q->domainMax( );
58
59 delete m_Q;
60 m_Q = new Functions::Constant1d( domainMin, domainMax, 0.0 );
61 }
62 a_setupInfo.m_Q = m_Q->evaluate( 0 ); // Needed for NBodyPhaseSpace.
63
64 GIDI::Suite const &products = a_outputChannel->products( );
65 if( m_channelType == ChannelType::twoBody ) {
66 if( !a_setupInfo.m_protare.isTNSL_ProtareSingle( ) ) {
67 GIDI::Product const *product = products.get<GIDI::Product>( 1 );
68 a_setupInfo.m_product2Mass = product->particle( ).mass( "MeV/c**2" ); // Includes nuclear excitation energy.
69 }
70 }
71
72 bool electronPresent = false;
73 std::size_t size = 0;
74 std::set<std::size_t> productsToDo;
75 for( std::size_t i1 = 0; i1 < a_outputChannel->products( ).size( ); ++i1 ) {
76 GIDI::Product const *product = products.get<GIDI::Product>( i1 );
77
78 if( product->particle( ).ID( ) == PoPI::IDs::electron ) electronPresent = true;
79 if( electronPresent && !a_setupInfo.m_isPhotoAtomicIncoherentScattering ) continue;
80 if( !electronPresent && !product->isCompleteParticle( ) && ( product->outputChannel( ) == nullptr ) ) continue;
81
82 if( ( product->outputChannel( ) != nullptr ) || a_settings.sampleNonTransportingParticles( ) || a_particles.hasParticle( product->particle( ).ID( ) ) )
83 productsToDo.insert( i1 );
84 }
85 size = productsToDo.size( );
86 if( a_setupInfo.m_isPairProduction ) {
87 size += 2;
88 size = 2; // This is a kludge until the ENDL to GNDS translator is fixed.
89 }
90
91 bool addIncoherentPhotoAtomicScatteringElectron = false;
92 if( a_setupInfo.m_isPhotoAtomicIncoherentScattering && a_particles.hasParticle( PoPI::IDs::electron ) ) { // May need to add electron for legacy GNDS files.
93 if( !electronPresent ) {
94// FIXME: BRB 7/Nov/2024, Why is an electron added, this is incoherent atomic scattering which does not emit an electron?
95 addIncoherentPhotoAtomicScatteringElectron = true;
96 ++size;
97 }
98 }
99 m_products.reserve( size );
100
101 if( a_setupInfo.m_isPairProduction ) {
102 std::string ID( PoPI::IDs::photon );
103 std::string label = ID;
104
105 Product *product = new Product( a_setupInfo.m_popsUser, ID, label );
106 product->setMultiplicity( new Functions::Constant1d( a_setupInfo.m_domainMin, a_setupInfo.m_domainMax, 1.0, 0.0 ) );
107 product->distribution( new Distributions::PairProductionGamma( a_setupInfo, true ) );
108 m_products.push_back( product );
109
110 label += "__a";
111 product = new Product( a_setupInfo.m_popsUser, ID, label );
112 product->setMultiplicity( new Functions::Constant1d( a_setupInfo.m_domainMin, a_setupInfo.m_domainMax, 1.0, 0.0 ) );
113 product->distribution( new Distributions::PairProductionGamma( a_setupInfo, false ) );
114 m_products.push_back( product );
115 }
116
117 for( std::size_t i1 = 0; i1 < a_outputChannel->products( ).size( ); ++i1 ) {
118 if( productsToDo.find( i1 ) == productsToDo.end( ) ) continue;
119
120 GIDI::Product const *product = products.get<GIDI::Product>( i1 );
121
122 if( a_setupInfo.m_isPairProduction ) {
123 if( !a_settings.sampleNonTransportingParticles( ) ) continue;
124 if( a_setupInfo.m_protare.targetIntid( ) != MCGIDI_popsIntid( a_setupInfo.m_pops, product->particle( ).ID( ) ) ) continue;
125 }
126 a_setupInfo.m_twoBodyOrder = TwoBodyOrder::notApplicable;
127 if( m_channelType == ChannelType::twoBody ) a_setupInfo.m_twoBodyOrder = ( ( i1 == 0 ? TwoBodyOrder::firstParticle : TwoBodyOrder::secondParticle ) );
128 m_products.push_back( new Product( product, a_setupInfo, a_settings, a_particles, m_isFission ) );
129
130 if( addIncoherentPhotoAtomicScatteringElectron && ( product->particle( ).ID( ) == PoPI::IDs::photon ) ) {
131 addIncoherentPhotoAtomicScatteringElectron = false;
132
133 Product *product2 = new Product( a_setupInfo.m_pops, PoPI::IDs::electron, PoPI::IDs::electron );
134 product2->setMultiplicity( new Functions::Constant1d( a_setupInfo.m_domainMin, a_setupInfo.m_domainMax, 1.0, 0.0 ) );
135 product2->distribution( new Distributions::IncoherentPhotoAtomicScatteringElectron( a_setupInfo ) );
136 m_products.push_back( product2 );
137 }
138 }
139
140 if( ( a_settings.delayedNeutrons( ) == GIDI::Transporting::DelayedNeutrons::on ) && a_particles.hasParticle( PoPI::IDs::neutron ) ) {
141 GIDI::FissionFragmentData const &fissionFragmentData = a_outputChannel->fissionFragmentData( );
142 GIDI::Suite const &delayedNeutrons = fissionFragmentData.delayedNeutrons( );
143
144 if( delayedNeutrons.size( ) > 0 ) {
145 bool missingData = false;
146 GIDI::Axes axes;
147 GIDI::Functions::XYs1d totalDelayedNeutronMultiplicity( axes, ptwXY_interpolationLinLin );
148
149 m_delayedNeutrons.reserve( delayedNeutrons.size( ) );
150 for( std::size_t i1 = 0; i1 < delayedNeutrons.size( ); ++i1 ) {
151 GIDI::DelayedNeutron const *delayedNeutron = delayedNeutrons.get<GIDI::DelayedNeutron>( i1 );
152 GIDI::Product const &product = delayedNeutron->product( );
153 GIDI::Suite const &multiplicity = product.multiplicity( );
154
155 GIDI::Functions::Function1dForm const *form1d = multiplicity.get<GIDI::Functions::Function1dForm>( 0 );
156
157 if( form1d->type( ) == GIDI::FormType::unspecified1d ) {
158 missingData = true;
159 break;
160 }
161
162 if( form1d->type( ) != GIDI::FormType::XYs1d ) {
163 std::cerr << "OutputChannel::OutputChannel: GIDI::DelayedNeutron multiplicity type != GIDI::FormType::XYs1d" << std::endl;
164 missingData = true;
165 break;
166 }
167
168 GIDI::Functions::XYs1d const *multiplicityXYs1d = static_cast<GIDI::Functions::XYs1d const *>( form1d );
169 totalDelayedNeutronMultiplicity += *multiplicityXYs1d;
170
171 m_delayedNeutrons.push_back( new DelayedNeutron( static_cast<int>( i1 ), delayedNeutron, a_setupInfo, a_settings, a_particles ) );
172 }
173 if( !missingData ) m_totalDelayedNeutronMultiplicity = new Functions::XYs1d( totalDelayedNeutronMultiplicity );
174 }
175 }
176
177 m_hasFinalStatePhotons = a_setupInfo.m_hasFinalStatePhotons;
178 }
179}
FormType type() const
Definition GIDI.hpp:667
T * get(std::size_t a_Index)
Definition GIDI.hpp:2642
LUPI_HOST_DEVICE Vector< Product * > const & products() const
Definition MCGIDI.hpp:1303
Vector< DelayedNeutron * > delayedNeutrons() const
Definition MCGIDI.hpp:1305
LUPI_HOST_DEVICE DelayedNeutron const * delayedNeutron(std::size_t a_index) const
Definition MCGIDI.hpp:1306
@ fissionFragmentData
Definition GIDI.hpp:119
LUPI_HOST Function1d_d1 * parseFunction1d_d1(Transporting::MC const &a_settings, GIDI::Suite const &a_suite)
LUPI_HOST int MCGIDI_popsIntid(PoPI::Database const &a_pops, std::string const &a_ID)
@ ptwXY_interpolationLinLin
Definition ptwXY.h:37
static std::string const photon
Definition PoPI.hpp:162
static std::string const neutron
Definition PoPI.hpp:164
static std::string const electron
Definition PoPI.hpp:163

◆ ~OutputChannel()

LUPI_HOST_DEVICE MCGIDI::OutputChannel::~OutputChannel ( )

Definition at line 184 of file MCGIDI_outputChannel.cc.

184 {
185
186 delete m_Q;
187 for( std::size_t i1 = 0; i1 < m_products.size( ); ++i1 ) delete m_products[i1];
188
189 delete m_totalDelayedNeutronMultiplicity;
190 for( std::size_t i1 = 0; i1 < m_delayedNeutrons.size( ); ++i1 ) delete m_delayedNeutrons[i1];
191}

Member Function Documentation

◆ angleBiasing()

template<typename RNG>
LUPI_HOST_DEVICE void MCGIDI::OutputChannel::angleBiasing ( Reaction const * a_reaction,
int a_index,
double a_temperature,
double a_energy_in,
double a_mu_lab,
double & a_weight,
double & a_energy_out,
RNG && a_rng,
double & a_cumulative_weight ) const
inline

Returns the probability for a project with energy a_energy_in to cause this channel to emitted a particle of index a_index at angle a_mu_lab as seen in the lab frame. If a particle is emitted, a_energy_out is its sampled outgoing energy.

Parameters
a_reaction[in] The reaction containing the particle which this distribution describes.
a_index[in] The index of the particle to emit.
a_temperature[in] Specifies the temperature of the material.
a_energy_in[in] The energy of the incident particle.
a_mu_lab[in] The desired mu in the lab frame for the emitted particle.
a_weight[in] The probability of emitting outgoing particle into lab angle a_mu_lab.
a_energy_out[in] The energy of the emitted outgoing particle.
a_rng[in] The random number generator function that returns a double in the range [0, 1.0).
a_cumulative_weight[in] The sum of the multiplicity for other outgoing particles with index a_index.

Definition at line 2469 of file MCGIDI_headerSource.hpp.

2470 {
2471
2472 for( Vector<Product *>::const_iterator iter = m_products.begin( ); iter != m_products.end( ); ++iter )
2473 (*iter)->angleBiasing( a_reaction, a_index, a_temperature, a_energy_in, a_mu_lab, a_weight, a_energy_out, a_rng, a_cumulative_weight );
2474
2475 if( ( m_totalDelayedNeutronMultiplicity != nullptr ) && ( a_index == m_neutronIndex ) ) {
2476 for( std::size_t i1 = 0; i1 < (std::size_t) m_delayedNeutrons.size( ); ++i1 ) {
2477 DelayedNeutron const *delayedNeutron1( delayedNeutron( i1 ) );
2478 Product const &product = delayedNeutron1->product( );
2479
2480 product.angleBiasing( a_reaction, a_index, a_temperature, a_energy_in, a_mu_lab, a_weight, a_energy_out, a_rng, a_cumulative_weight );
2481 }
2482 }
2483}

◆ angleBiasingViaIntid()

template<typename RNG>
LUPI_HOST_DEVICE void MCGIDI::OutputChannel::angleBiasingViaIntid ( Reaction const * a_reaction,
int a_intid,
double a_temperature,
double a_energy_in,
double a_mu_lab,
double & a_weight,
double & a_energy_out,
RNG && a_rng,
double & a_cumulative_weight ) const
inline

Returns the probability for a project with energy a_energy_in to cause this channel to emitted a particle of intid a_intid at angle a_mu_lab as seen in the lab frame. If a particle is emitted, a_energy_out is its sampled outgoing energy.

Parameters
a_reaction[in] The reaction containing the particle which this distribution describes.
a_intid[in] The intid of the particle to emit.
a_temperature[in] Specifies the temperature of the material.
a_energy_in[in] The energy of the incident particle.
a_mu_lab[in] The desired mu in the lab frame for the emitted particle.
a_weight[in] The probability of emitting outgoing particle into lab angle a_mu_lab.
a_energy_out[in] The energy of the emitted outgoing particle.
a_rng[in] The random number generator function that returns a double in the range [0, 1.0).
a_cumulative_weight[in] The sum of the multiplicity for other outgoing particles with intid a_intid.

Definition at line 2501 of file MCGIDI_headerSource.hpp.

2502 {
2503
2504 for( Vector<Product *>::const_iterator iter = m_products.begin( ); iter != m_products.end( ); ++iter )
2505 (*iter)->angleBiasingViaIntid( a_reaction, a_intid, a_temperature, a_energy_in, a_mu_lab, a_weight, a_energy_out, a_rng, a_cumulative_weight );
2506
2507 if( ( m_totalDelayedNeutronMultiplicity != nullptr ) && ( a_intid == PoPI::Intids::neutron ) ) {
2508 for( std::size_t i1 = 0; i1 < (std::size_t) m_delayedNeutrons.size( ); ++i1 ) {
2509 DelayedNeutron const *delayedNeutron1( delayedNeutron( i1 ) );
2510 Product const &product = delayedNeutron1->product( );
2511
2512 product.angleBiasingViaIntid( a_reaction, a_intid, a_temperature, a_energy_in, a_mu_lab, a_weight, a_energy_out, a_rng, a_cumulative_weight );
2513 }
2514 }
2515}
static int constexpr neutron
Definition PoPI.hpp:177

◆ delayedNeutron()

LUPI_HOST_DEVICE DelayedNeutron const * MCGIDI::OutputChannel::delayedNeutron ( std::size_t a_index) const
inline

Definition at line 1306 of file MCGIDI.hpp.

1306{ return( m_delayedNeutrons[a_index] ); }

Referenced by angleBiasing(), angleBiasingViaIntid(), OutputChannel(), and sampleProducts().

◆ delayedNeutrons()

Vector< DelayedNeutron * > MCGIDI::OutputChannel::delayedNeutrons ( ) const
inline

Definition at line 1305 of file MCGIDI.hpp.

1305{ return( m_delayedNeutrons ); }

Referenced by OutputChannel().

◆ finalQ()

LUPI_HOST_DEVICE double MCGIDI::OutputChannel::finalQ ( double a_x1) const

This method returns the final Q for this by getting its final Q plus any sub-output channel's finalQ.

Parameters
a_x1[in] The energy of the projectile.
Returns
The Q-value at product energy a_x1.

Definition at line 201 of file MCGIDI_outputChannel.cc.

201 {
202
203 double final_Q = m_Q->evaluate( a_x1 );
204
205 for( std::size_t i1 = 0; i1 < m_products.size( ); ++i1 ) final_Q += m_products[i1]->finalQ( a_x1 );
206 return( final_Q );
207}
LUPI_HOST_DEVICE double finalQ(double a_x1) const

Referenced by finalQ(), and MCGIDI::Reaction::finalQ().

◆ hasFission()

LUPI_HOST_DEVICE bool MCGIDI::OutputChannel::hasFission ( ) const

This method returns true if the output channel or any of its sub-output channels is a fission channel and false otherwise.

Returns
true if this or any sub-output channel is a fission channel and false otherwise.

Definition at line 215 of file MCGIDI_outputChannel.cc.

215 {
216
217 if( m_isFission ) return( true );
218 for( std::size_t i1 = 0; i1 < m_products.size( ); ++i1 ) {
219 if( m_products[i1]->hasFission( ) ) return( true );
220 }
221 return( false );
222}
LUPI_HOST_DEVICE bool hasFission() const

Referenced by hasFission().

◆ isFission()

LUPI_HOST_DEVICE bool MCGIDI::OutputChannel::isFission ( ) const
inline

Returns the value of the m_isFission.

Definition at line 1294 of file MCGIDI.hpp.

◆ isTwoBody()

LUPI_HOST_DEVICE bool MCGIDI::OutputChannel::isTwoBody ( ) const
inline

Returns true if output channel is two-body and false otherwise.

Definition at line 1292 of file MCGIDI.hpp.

◆ moveProductsEtAlToReaction()

LUPI_HOST void MCGIDI::OutputChannel::moveProductsEtAlToReaction ( std::vector< Product * > & a_products,
Functions::Function1d ** a_totalDelayedNeutronMultiplicity,
std::vector< DelayedNeutron * > & a_delayedNeutrons,
std::vector< Functions::Function1d_d1 * > & a_Qs )

Returns the energy dependent multiplicity for outgoing particle with pops id a_id. The returned value may not be an integer. Energy dependent multiplicity mainly occurs for photons and fission neutrons.

Parameters
a_products[in] The std::vector instance to add products to.
a_totalDelayedNeutronMultiplicity[in] The std::vector instance to add delayed neutron multiplicities to.
a_delayedNeutrons[in] The std::vector instance to add delayed neutrons to.
a_Qs[in] The std::vector instance to add Q functions to.

Definition at line 271 of file MCGIDI_outputChannel.cc.

272 {
273
274 if( a_totalDelayedNeutronMultiplicity != nullptr ) { /* This will not work if fission is a nested channel. Ergo "n + (R -> fission)". */
275 *a_totalDelayedNeutronMultiplicity = m_totalDelayedNeutronMultiplicity;
276 m_totalDelayedNeutronMultiplicity = nullptr;
277 for( std::size_t index = 0; index < m_delayedNeutrons.size( ); ++index ) {
278 a_delayedNeutrons.push_back( m_delayedNeutrons[index] );
279 m_delayedNeutrons[index] = nullptr;
280 }
281 }
282
283 a_Qs.push_back( m_Q );
284 m_Q = nullptr;
285 for( std::size_t productIndex = 0; productIndex < m_products.size( ); ++productIndex ) {
286 Product *product = m_products[productIndex];
287
288 if( product->outputChannel( ) != nullptr ) {
289 product->outputChannel( )->moveProductsEtAlToReaction( a_products, nullptr, a_delayedNeutrons, a_Qs );
290 delete product; }
291 else {
292 a_products.push_back( product );
293 }
294 m_products[productIndex] = nullptr;
295 }
296}

◆ operator[]()

LUPI_HOST_DEVICE Product * MCGIDI::OutputChannel::operator[] ( std::size_t a_index)
inline

Returns a pointer to the product at index a_index.

Definition at line 1290 of file MCGIDI.hpp.

◆ productAverageMultiplicity()

LUPI_HOST_DEVICE double MCGIDI::OutputChannel::productAverageMultiplicity ( int a_index,
double a_projectileEnergy ) const

Returns the energy dependent multiplicity for outgoing particle with pops index a_index. The returned value may not be an integer. Energy dependent multiplicity mainly occurs for photons and fission neutrons.

Parameters
a_index[in] The index of the requested particle.
a_projectileEnergy[in] The energy of the projectile.
Returns
The multiplicity value for the requested particle.

Definition at line 376 of file MCGIDI_outputChannel.cc.

376 {
377
378 double multiplicity = 0.0;
379
380 for( Vector<Product *>::const_iterator iter = m_products.begin( ); iter != m_products.end( ); ++iter ) {
381 multiplicity += (*iter)->productAverageMultiplicity( a_index, a_projectileEnergy );
382 }
383
384 if( m_totalDelayedNeutronMultiplicity != nullptr ) {
385 if( a_index == m_neutronIndex ) multiplicity += m_totalDelayedNeutronMultiplicity->evaluate( a_projectileEnergy );
386 }
387
388 return( multiplicity );
389}

Referenced by MCGIDI::Reaction::productAverageMultiplicity().

◆ productAverageMultiplicityViaIntid()

LUPI_HOST_DEVICE double MCGIDI::OutputChannel::productAverageMultiplicityViaIntid ( int a_intid,
double a_projectileEnergy ) const

Returns the energy dependent multiplicity for outgoing particle with pops intid a_intid. The returned value may not be an integer. Energy dependent multiplicity mainly occurs for photons and fission neutrons.

Parameters
a_intid[in] The intid of the requested particle.
a_projectileEnergy[in] The energy of the projectile.
Returns
The multiplicity value for the requested particle.

Definition at line 401 of file MCGIDI_outputChannel.cc.

401 {
402
403 double multiplicity = 0.0;
404
405 for( Vector<Product *>::const_iterator iter = m_products.begin( ); iter != m_products.end( ); ++iter ) {
406 multiplicity += (*iter)->productAverageMultiplicityViaIntid( a_intid, a_projectileEnergy );
407 }
408
409 if( m_totalDelayedNeutronMultiplicity != nullptr ) {
410 if( a_intid == PoPI::Intids::neutron ) multiplicity += m_totalDelayedNeutronMultiplicity->evaluate( a_projectileEnergy );
411 }
412
413 return( multiplicity );
414}

Referenced by MCGIDI::Reaction::productAverageMultiplicityViaIntid().

◆ products()

LUPI_HOST_DEVICE Vector< Product * > const & MCGIDI::OutputChannel::products ( ) const
inline

Returns the value of the m_products.

Definition at line 1303 of file MCGIDI.hpp.

Referenced by OutputChannel().

◆ Q()

LUPI_HOST_DEVICE Functions::Function1d_d1 * MCGIDI::OutputChannel::Q ( )
inline

Returns the pointer of the m_Q member.

Definition at line 1297 of file MCGIDI.hpp.

◆ sampleProducts()

template<typename RNG, typename PUSHBACK>
LUPI_HOST_DEVICE void MCGIDI::OutputChannel::sampleProducts ( ProtareSingle const * a_protare,
double a_projectileEnergy,
Sampling::Input & a_input,
RNG && a_rng,
PUSHBACK && a_push_back,
Sampling::ProductHandler & a_products ) const
inline

This method adds sampled products to a_products.

Parameters
a_protare[in] The Protare this Reaction belongs to.
a_projectileEnergy[in] The energy of the projectile.
a_input[in] Sample options requested by user.
a_rng[in] The random number generator function that returns a double in the range [0, 1.0).
a_products[in] The object to add all sampled products to.

Definition at line 2406 of file MCGIDI_headerSource.hpp.

2407 {
2408
2409 if( m_hasFinalStatePhotons ) {
2410 double random = a_rng( );
2411 double cumulative = 0.0;
2412 bool sampled = false;
2413 for( auto productIter = m_products.begin( ); productIter != m_products.end( ); ++productIter ) {
2414 cumulative += (*productIter)->multiplicity( )->evaluate( a_projectileEnergy );
2415 if( cumulative >= random ) {
2416 (*productIter)->sampleFinalState( a_protare, a_projectileEnergy, a_input, a_rng, a_push_back, a_products );
2417 sampled = true;
2418 break;
2419 }
2420 }
2421 if( !sampled ) { // BRB: FIXME: still need to code for continuum photon.
2422 } }
2423 else {
2424 for( Vector<Product *>::const_iterator iter = m_products.begin( ); iter != m_products.end( ); ++iter )
2425 (*iter)->sampleProducts( a_protare, a_projectileEnergy, a_input, a_rng, a_push_back, a_products );
2426 }
2427
2428 if( m_totalDelayedNeutronMultiplicity != nullptr ) {
2429 double totalDelayedNeutronMultiplicity = m_totalDelayedNeutronMultiplicity->evaluate( a_projectileEnergy );
2430
2431 if( a_rng( ) < totalDelayedNeutronMultiplicity ) { // Assumes that totalDelayedNeutronMultiplicity < 1.0, which it is.
2432 double sum = 0.0;
2433
2434 totalDelayedNeutronMultiplicity *= a_rng( );
2435 for( std::size_t i1 = 0; i1 < (std::size_t) m_delayedNeutrons.size( ); ++i1 ) {
2436 DelayedNeutron const *delayedNeutron1( delayedNeutron( i1 ) );
2437 Product const &product = delayedNeutron1->product( );
2438
2439 sum += product.multiplicity( )->evaluate( a_projectileEnergy );
2440 if( sum >= totalDelayedNeutronMultiplicity ) {
2441 product.distribution( )->sample( a_projectileEnergy, a_input, a_rng );
2442 a_input.m_delayedNeutronIndex = delayedNeutron1->delayedNeutronIndex( );
2443 a_input.m_delayedNeutronDecayRate = delayedNeutron1->rate( );
2444 a_products.add( a_projectileEnergy, product.intid( ), product.index( ), product.userParticleIndex( ), product.mass( ),
2445 a_input, a_rng, a_push_back, false );
2446 break;
2447 }
2448 }
2449 }
2450 }
2451}

◆ serialize()

LUPI_HOST_DEVICE void MCGIDI::OutputChannel::serialize ( LUPI::DataBuffer & a_buffer,
LUPI::DataBuffer::Mode a_mode )

This method serializes this for broadcasting as needed for MPI and GPUs. The method can count the number of required bytes, pack this or unpack this depending on a_mode.

Parameters
a_buffer[in] The buffer to read or write data to depending on a_mode.
a_mode[in] Specifies the action of this method.

Definition at line 424 of file MCGIDI_outputChannel.cc.

424 {
425
426 int channelType = 0;
427 switch( m_channelType ) {
428 case ChannelType::none :
429 break;
431 channelType = 1;
432 break;
434 channelType = 2;
435 break;
436 }
437 DATA_MEMBER_INT( channelType, a_buffer, a_mode );
438 if( a_mode == LUPI::DataBuffer::Mode::Unpack ) {
439 switch( channelType ) {
440 case 0 :
441 m_channelType = ChannelType::none;
442 break;
443 case 1 :
444 m_channelType = ChannelType::twoBody;
445 break;
446 case 2 :
447 m_channelType = ChannelType::uncorrelatedBodies;
448 break;
449 }
450 }
451
452 DATA_MEMBER_INT( m_neutronIndex, a_buffer, a_mode );
453 DATA_MEMBER_CAST( m_isFission, a_buffer, a_mode, bool );
454 DATA_MEMBER_CAST( m_hasFinalStatePhotons, a_buffer, a_mode, bool );
455
456 m_Q = serializeFunction1d_d1( a_buffer, a_mode, m_Q );
457 serializeProducts( a_buffer, a_mode, m_products );
458 m_totalDelayedNeutronMultiplicity = serializeFunction1d( a_buffer, a_mode, m_totalDelayedNeutronMultiplicity );
459 serializeDelayedNeutrons( a_buffer, a_mode, m_delayedNeutrons );
460}
G4ConcreteNNToDeltaDelta channelType
#define DATA_MEMBER_CAST(member, buf, mode, someType)
#define DATA_MEMBER_INT( member, buf, mode)
LUPI_HOST_DEVICE void serializeDelayedNeutrons(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode, Vector< DelayedNeutron * > &a_delayedNeutrons)
LUPI_HOST_DEVICE Functions::Function1d * serializeFunction1d(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode, Functions::Function1d *a_function1d)
LUPI_HOST_DEVICE void serializeProducts(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode, Vector< Product * > &a_products)
LUPI_HOST_DEVICE Functions::Function1d_d1 * serializeFunction1d_d1(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode, Functions::Function1d_d1 *a_function1d)

Referenced by MCGIDI::Reaction::serialize().

◆ setModelDBRC_data()

LUPI_HOST void MCGIDI::OutputChannel::setModelDBRC_data ( Sampling::Upscatter::ModelDBRC_data * a_modelDBRC_data)

This method calls the **setModelDBRC_data* method on the first product of this with a_modelDBRC_data.

Parameters
a_modelDBRC_data[in] The instance storing data needed to treat the DRRC upscatter mode.

Definition at line 256 of file MCGIDI_outputChannel.cc.

256 {
257
258 m_products[0]->setModelDBRC_data( a_modelDBRC_data );
259}

Referenced by MCGIDI::Reaction::setModelDBRC_data().

◆ setUserParticleIndex()

LUPI_HOST void MCGIDI::OutputChannel::setUserParticleIndex ( int a_particleIndex,
int a_userParticleIndex )

Updates the m_userParticleIndex to a_userParticleIndex for all particles with PoPs index a_particleIndex.

Parameters
a_particleIndex[in] The PoPs index of the particle whose user index is to be set.
a_userParticleIndex[in] The particle index specified by the user.

Definition at line 231 of file MCGIDI_outputChannel.cc.

231 {
232
233 for( auto iter = m_products.begin( ); iter != m_products.end( ); ++iter ) (*iter)->setUserParticleIndex( a_particleIndex, a_userParticleIndex );
234 for( auto iter = m_delayedNeutrons.begin( ); iter != m_delayedNeutrons.end( ); ++iter ) (*iter)->setUserParticleIndex( a_particleIndex, a_userParticleIndex );
235}

Referenced by MCGIDI::Reaction::setUserParticleIndex().

◆ setUserParticleIndexViaIntid()

LUPI_HOST void MCGIDI::OutputChannel::setUserParticleIndexViaIntid ( int a_particleIntid,
int a_userParticleIndex )

Updates the m_userParticleIndex to a_userParticleIndex for all particles with PoPs intid a_particleIntid.

Parameters
a_particleIntid[in] The PoPs intid of the particle whose user index is to be set.
a_userParticleIndex[in] The particle index specified by the user.

Definition at line 244 of file MCGIDI_outputChannel.cc.

244 {
245
246 for( auto iter = m_products.begin( ); iter != m_products.end( ); ++iter ) (*iter)->setUserParticleIndexViaIntid( a_particleIntid, a_userParticleIndex );
247 for( auto iter = m_delayedNeutrons.begin( ); iter != m_delayedNeutrons.end( ); ++iter ) (*iter)->setUserParticleIndexViaIntid( a_particleIntid, a_userParticleIndex );
248}

Referenced by MCGIDI::Reaction::setUserParticleIndexViaIntid().


The documentation for this class was generated from the following files: