Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
MCGIDI_protareComposite.cc
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#include <map>
11
12#include "MCGIDI.hpp"
13
14namespace MCGIDI {
15
16/*! \class ProtareComposite
17 * Class representing a **GNDS** <**reactionSuite**> node with only data needed for Monte Carlo transport. The
18 * data are also stored in a way that is better suited for Monte Carlo transport. For example, cross section data
19 * for each reaction are not stored with its reaction, but within the HeatedCrossSections member of the Protare.
20 */
21
22/* *********************************************************************************************************//**
23 * Default constructor used when broadcasting a Protare as needed by MPI or GPUs.
24 ***********************************************************************************************************/
25
28 m_numberOfReactions( 0 ),
29 m_numberOfOrphanProducts( 0 ),
30 m_minimumEnergy( 0.0 ),
31 m_maximumEnergy( 0.0 ) {
32
33}
34
35/* *********************************************************************************************************//**
36 * @param a_smr [Out] If errors are not to be thrown, then the error is reported via this instance.
37 * @param a_protare [in] The GIDI::Protare whose data is to be used to construct *this*.
38 * @param a_pops [in] A PoPs Database instance used to get particle intids and possibly other particle information.
39 * @param a_settings [in] Used to pass user options to the *this* to instruct it which data are desired.
40 * @param a_particles [in] List of transporting particles and their information (e.g., multi-group boundaries and fluxes).
41 * @param a_domainHash [in] The hash data used when looking up a cross section.
42 * @param a_temperatureInfos [in] The list of temperature data to extract from *a_protare*.
43 * @param a_reactionsToExclude [in] A list of reaction to not include in the MCGIDI::Protare. This currently does not work for ProtareComposite.
44 * @param a_reactionsToExcludeOffset [in] The starting index for the reactions in this ProtareSingle.
45 * @param a_allowFixedGrid [in] For internal (i.e., MCGIDI) use only. Users must use the default value.
46 ***********************************************************************************************************/
47
49 Transporting::MC &a_settings, GIDI::Transporting::Particles const &a_particles, DomainHash const &a_domainHash,
50 GIDI::Styles::TemperatureInfos const &a_temperatureInfos, GIDI::ExcludeReactionsSet const &a_reactionsToExclude,
51 std::size_t a_reactionsToExcludeOffset, LUPI_maybeUnused bool a_allowFixedGrid ) :
52 Protare( ProtareType::composite, a_protare, a_settings, a_pops ),
53 m_numberOfReactions( 0 ),
54 m_numberOfOrphanProducts( 0 ) {
55
56 std::vector<GIDI::Protare *> &protares = static_cast<std::vector<GIDI::Protare *> &>( const_cast<GIDI::ProtareComposite &>( a_protare ).protares( ) );
57 std::size_t length = static_cast<std::size_t>( protares.size( ) );
58
59 std::set<int> product_intids;
60 std::set<int> product_intids_transportable;
61 std::set<int> product_indices;
62 std::set<int> product_indices_transportable;
63
64 m_protares.resize( length );
65 for( std::size_t i1 = 0; i1 < length; ++i1 ) {
66 GIDI::Protare const *protare = protares[i1];
67
68 m_protares[i1] = static_cast<ProtareSingle *>( protareFromGIDIProtare( a_smr, *protare, a_pops, a_settings, a_particles, a_domainHash, a_temperatureInfos, a_reactionsToExclude, a_reactionsToExcludeOffset, false ) );
69
70 m_numberOfReactions += m_protares[i1]->numberOfReactions( );
71 m_numberOfOrphanProducts += m_protares[i1]->numberOfOrphanProducts( );
72
73 if( i1 == 0 ) {
74 m_minimumEnergy = m_protares[0]->minimumEnergy( );
75 m_maximumEnergy = m_protares[0]->maximumEnergy( );
76 }
77 if( m_protares[i1]->minimumEnergy( ) < m_minimumEnergy ) m_minimumEnergy = m_protares[i1]->minimumEnergy( );
78 if( m_protares[i1]->maximumEnergy( ) > m_maximumEnergy ) m_maximumEnergy = m_protares[i1]->maximumEnergy( );
79
80 addVectorItemsToSet( m_protares[i1]->productIntids( false ), product_intids );
81 addVectorItemsToSet( m_protares[i1]->productIntids( true ), product_intids_transportable );
82 addVectorItemsToSet( m_protares[i1]->productIndices( false ), product_indices );
83 addVectorItemsToSet( m_protares[i1]->productIndices( true ), product_indices_transportable );
84
85 a_reactionsToExcludeOffset += protare->numberOfReactions( );
86 }
87
88 productIntidsAndIndices( product_intids, product_intids_transportable, product_indices, product_indices_transportable );
89}
90
91/* *********************************************************************************************************//**
92 ***********************************************************************************************************/
93
95
96 std::size_t length = static_cast<std::size_t>( m_protares.size( ) );
97
98 for( std::size_t i1 = 0; i1 < length; ++i1 ) delete m_protares[i1];
99}
100
101/* *********************************************************************************************************//**
102 * Updates the m_userParticleIndex to *a_userParticleIndex* for all particles with PoPs index *a_particleIndex*.
103 *
104 * @param a_particleIndex [in] The PoPs index of the particle whose user index is to be set.
105 * @param a_userParticleIndex [in] The particle index specified by the user.
106 ***********************************************************************************************************/
107
108LUPI_HOST void ProtareComposite::setUserParticleIndex2( int a_particleIndex, int a_userParticleIndex ) {
109
110 for( auto iter = m_protares.begin( ); iter != m_protares.end( ); ++iter ) (*iter)->setUserParticleIndex( a_particleIndex, a_userParticleIndex );
111}
112
113/* *********************************************************************************************************//**
114 * Updates the m_userParticleIndex to *a_userParticleIndex* for all particles with PoPs intid *a_particleIntid*.
115 *
116 * @param a_particleIntid [in] The PoPs intid of the particle whose user index is to be set.
117 * @param a_userParticleIndex [in] The particle index specified by the user.
118 ***********************************************************************************************************/
119
120LUPI_HOST void ProtareComposite::setUserParticleIndexViaIntid2( int a_particleIntid, int a_userParticleIndex ) {
121
122 for( auto iter = m_protares.begin( ); iter != m_protares.end( ); ++iter ) (*iter)->setUserParticleIndexViaIntid( a_particleIntid, a_userParticleIndex );
123}
124
125/* *********************************************************************************************************//**
126 * Returns the pointer representing the (a_index - 1)th **ProtareSingle**.
127 *
128 * @param a_index [in] Index of the **ProtareSingle** to return.
129 *
130 * @return Pointer to the requested protare or nullptr if invalid *a_index*..
131 ***********************************************************************************************************/
132
133LUPI_HOST_DEVICE ProtareSingle const *ProtareComposite::protare( std::size_t a_index ) const {
134
135 for( std::size_t i1 = 0; i1 < m_protares.size( ); ++i1 ) {
136 std::size_t number = m_protares[i1]->numberOfProtares( );
137
138 if( number > a_index ) return( m_protares[i1]->protare( a_index ) );
139 a_index -= number;
140 }
141
142 return( nullptr );
143}
144
145/* *********************************************************************************************************//**
146 * Returns the pointer representing the (a_index - 1)th **ProtareSingle**.
147 *
148 * @param a_index [in] Index of the **ProtareSingle** to return.
149 *
150 * @return Pointer to the requested protare or nullptr if invalid *a_index*..
151 ***********************************************************************************************************/
152
154
155 for( std::size_t i1 = 0; i1 < m_protares.size( ); ++i1 ) {
156 std::size_t number = m_protares[i1]->numberOfProtares( );
157
158 if( number > a_index ) return( m_protares[i1]->protare( a_index ) );
159 a_index -= number;
160 }
161
162 return( nullptr );
163}
164
165/* *********************************************************************************************************//**
166 * Returns the pointer to the **ProtareSingle** that contains the (a_index - 1)th reaction.
167 *
168 * @param a_index [in] Index of the reaction.
169 *
170 * @return Pointer to the requested protare or nullptr if invalid *a_index*..
171 ***********************************************************************************************************/
172
174
175 for( std::size_t i1 = 0; i1 < m_protares.size( ); ++i1 ) {
176 std::size_t numberOfReactions = m_protares[i1]->numberOfReactions( );
177
178 if( a_index < numberOfReactions ) return( m_protares[i1] );
179 a_index -= numberOfReactions;
180 }
181
182 return( nullptr );
183}
184
185/* *********************************************************************************************************//**
186 * Returns the list of temperatures for the requested ProtareSingle.
187 *
188 * @param a_index [in] Index of the reqested ProtareSingle.
189 *
190 * @return Vector of doubles.
191 ***********************************************************************************************************/
192
194
195 for( std::size_t i1 = 0; i1 < m_protares.size( ); ++i1 ) {
196 std::size_t number = m_protares[i1]->numberOfProtares( );
197
198 if( number > a_index ) return( m_protares[i1]->temperatures( a_index ) );
199 a_index -= number;
200 }
201
202 LUPI_THROW( "ProtareSingle::temperatures: a_index not in range." );
203
204 Vector<double> temps; // Only to stop compilers from complaining.
205 return( temps );
206}
207
208/* *********************************************************************************************************//**
209 * Returns the reaction at index *a_index*. If *a_index* is negative, the reaction of the TNSL protare at index -*a_index* is
210 * returned; otherwise, the reaction from the regular protare at index *a_index* is returned.
211 *
212 * @param a_index [in] The index of the reaction to return.
213 *
214 * @return The reaction at index *a_index*.
215 ***********************************************************************************************************/
216
217LUPI_HOST_DEVICE Reaction const *ProtareComposite::reaction( std::size_t a_index ) const {
218
219 std::size_t length = static_cast<std::size_t>( m_protares.size( ) );
220
221 for( std::size_t i1 = 0; i1 < length; ++i1 ) {
222 std::size_t numberOfReactions = m_protares[i1]->numberOfReactions( );
223
224 if( a_index < numberOfReactions ) return( m_protares[i1]->reaction( a_index ) );
225 a_index -= numberOfReactions;
226 }
227
228 return( nullptr );
229}
230
231/* *********************************************************************************************************//**
232 * Returns the reaction at index *a_index*. If *a_index* is negative, the reaction of the TNSL protare at index -*a_index* is
233 * returned; otherwise, the reaction from the regular protare at index *a_index* is returned.
234 *
235 * @param a_index [in] The index of the reaction to return.
236 *
237 * @return The reaction at index *a_index*.
238 ***********************************************************************************************************/
239
240LUPI_HOST_DEVICE Reaction const *ProtareComposite::orphanProduct( std::size_t a_index ) const {
241
242 std::size_t length = static_cast<std::size_t>( m_protares.size( ) );
243
244 for( std::size_t i1 = 0; i1 < length; ++i1 ) {
245 std::size_t numberOfReactions = m_protares[i1]->numberOfOrphanProducts( );
246
247 if( a_index < numberOfReactions ) return( m_protares[i1]->orphanProduct( a_index ) );
248 a_index -= numberOfReactions;
249 }
250
251 return( nullptr );
252}
253
254/* *********************************************************************************************************//**
255 * Returns *true* is one of the protares has a fission channel and *false* otherwise.
256 *
257 * @return *true* is one of the protares has a fission channel and *false* otherwise.
258 ***********************************************************************************************************/
259
261
262 std::size_t length = static_cast<std::size_t>( m_protares.size( ) );
263
264 for( std::size_t i1 = 0; i1 < length; ++i1 ) {
265 if( m_protares[i1]->hasFission( ) ) return( true );
266 }
267
268 return( false );
269}
270
271/* *********************************************************************************************************//**
272 * Returns true if *this* has a photoatomic incoherent doppler broadened reaction and false otherwise.
273 *
274 * @return *true* is one of the protares has a fission channel and *false* otherwise.
275 ***********************************************************************************************************/
276
278
279 std::size_t length = static_cast<std::size_t>( m_protares.size( ) );
280
281 for( std::size_t i1 = 0; i1 < length; ++i1 ) {
282 if( m_protares[i1]->hasIncoherentDoppler( ) ) return( true );
283 }
284
285 return( false );
286}
287
288/* *********************************************************************************************************//**
289 * Returns true if one the the sub-protares of *this* has a unresolved resonance region (URR) data and false otherwise.
290 *
291 * @return true is if *this* has a URR data.
292 ***********************************************************************************************************/
293
295
296 std::size_t length = static_cast<std::size_t>( m_protares.size( ) );
297
298 for( std::size_t i1 = 0; i1 < length; ++i1 ) {
299 if( m_protares[i1]->hasURR_probabilityTables( ) ) return( true );
300 }
301
302 return( false );
303}
304
305/* *********************************************************************************************************//**
306 * Returns the minimum energy for the unresolved resonance region (URR) domain. If no URR data present, returns -1.
307 *
308 * @return The energy or -1 if not URR data present.
309 ***********************************************************************************************************/
310
312
313 std::size_t length = static_cast<std::size_t>( m_protares.size( ) );
314 double URR_domain_min = 1e32;
315
316 for( std::size_t i1 = 0; i1 < length; ++i1 ) {
317 if( m_protares[i1]->hasURR_probabilityTables( ) ) {
318 if( URR_domain_min > m_protares[i1]->URR_domainMin( ) ) URR_domain_min = m_protares[i1]->URR_domainMin( );
319 }
320 }
321
322 if( URR_domain_min == 1e32 ) URR_domain_min = -1.0;
323
324 return( URR_domain_min );
325}
326
327/* *********************************************************************************************************//**
328 * Returns the maximum energy for the unresolved resonance region (URR) domain. If no URR data present, returns -1.
329 *
330 * @return true is if *this* has a URR data.
331 ***********************************************************************************************************/
332
334
335 std::size_t length = static_cast<std::size_t>( m_protares.size( ) );
336 double URR_domain_max = -1.0;
337
338 for( std::size_t i1 = 0; i1 < length; ++i1 ) {
339 if( m_protares[i1]->hasURR_probabilityTables( ) ) {
340 if( URR_domain_max < m_protares[i1]->URR_domainMax( ) ) URR_domain_max = m_protares[i1]->URR_domainMax( );
341
342 }
343 }
344
345 return( URR_domain_max );
346}
347
348/* *********************************************************************************************************//**
349 * Returns *true* if the reaction at index *a_index* has URR robability tables and false otherwise.
350 *
351 * @param a_index [in] The index of the reaction.
352 *
353 * @return *true* if the reaction has URR robability tables and false otherwise.
354 ***********************************************************************************************************/
355
357
358 std::size_t length = m_protares.size( );
359
360 for( std::size_t i1 = 0; i1 < length; ++i1 ) {
361 std::size_t numberOfReactions = m_protares[i1]->numberOfReactions( );
362
363 if( a_index < numberOfReactions ) return( m_protares[i1]->reactionHasURR_probabilityTables( a_index ) );
364 a_index -= numberOfReactions;
365 }
366
367 return( false );
368}
369
370/* *********************************************************************************************************//**
371 * Returns the threshold for the reaction at index *a_index*. If *a_index* is negative, it is set to 0 before the
372 * threshold in the regular protare is returned.
373 *
374 * @param a_index [in] The index of the reaction.
375 *
376 * @return The threshold for reaction at index *a_index*.
377 ***********************************************************************************************************/
378
379LUPI_HOST_DEVICE double ProtareComposite::threshold( std::size_t a_index ) const {
380
381 std::size_t length = static_cast<std::size_t>( m_protares.size( ) );
382
383 for( std::size_t i1 = 0; i1 < length; ++i1 ) {
384 std::size_t numberOfReactions = m_protares[i1]->numberOfReactions( );
385
386 if( a_index < numberOfReactions ) return( m_protares[i1]->threshold( a_index ) );
387 a_index -= numberOfReactions;
388 }
389
390 return( 0.0 );
391}
392
393/* *********************************************************************************************************//**
394 * Returns the total cross section.
395 *
396 * @param a_URR_protareInfos [in] URR information.
397 * @param a_hashIndex [in] The cross section hash index.
398 * @param a_temperature [in] The target temperature.
399 * @param a_energy [in] The projectile energy.
400 * @param a_sampling [in] Only used for multi-group cross sections. When sampling, the cross section in the group where threshold
401 * is present the cross section is augmented.
402 *
403 * @return The total cross section.
404 ***********************************************************************************************************/
405
406LUPI_HOST_DEVICE double ProtareComposite::crossSection( URR_protareInfos const &a_URR_protareInfos, std::size_t a_hashIndex, double a_temperature, double a_energy, bool a_sampling ) const {
407
408 std::size_t length = static_cast<std::size_t>( m_protares.size( ) );
409 double cross_section = 0.0;
410
411 for( std::size_t i1 = 0; i1 < length; ++i1 ) cross_section += m_protares[i1]->crossSection( a_URR_protareInfos, a_hashIndex, a_temperature, a_energy, a_sampling );
412
413 return( cross_section );
414}
415
416/* *********************************************************************************************************//**
417 * Adds the energy dependent, total cross section corresponding to the temperature *a_temperature* multiplied by *a_userFactor* to *a_crossSectionVector*.
418 *
419 * @param a_temperature [in] Specifies the temperature of the material.
420 * @param a_userFactor [in] User factor which all cross sections are multiplied by.
421 * @param a_numberAllocated [in] The length of memory allocated for *a_crossSectionVector*.
422 * @param a_crossSectionVector [in/out] The energy dependent, total cross section to add cross section data to.
423 ***********************************************************************************************************/
424
425LUPI_HOST_DEVICE void ProtareComposite::crossSectionVector( double a_temperature, double a_userFactor, std::size_t a_numberAllocated, double *a_crossSectionVector ) const {
426
427 std::size_t length = static_cast<std::size_t>( m_protares.size( ) );
428
429 for( std::size_t i1 = 0; i1 < length; ++i1 ) m_protares[i1]->crossSectionVector( a_temperature, a_userFactor, a_numberAllocated, a_crossSectionVector );
430}
431
432/* *********************************************************************************************************//**
433 * Returns the cross section for reaction at index *a_reactionIndex*.
434 *
435 * @param a_reactionIndex [in] The index of the reaction.
436 * @param a_URR_protareInfos [in] URR information.
437 * @param a_hashIndex [in] The cross section hash index.
438 * @param a_temperature [in] The target temperature.
439 * @param a_energy [in] The projectile energy.
440 * @param a_sampling [in] Only used for multi-group cross sections. When sampling, the cross section in the group where threshold
441 * is present the cross section is augmented.
442 *
443 * @return The total cross section.
444 ***********************************************************************************************************/
445
446LUPI_HOST_DEVICE double ProtareComposite::reactionCrossSection( std::size_t a_reactionIndex, URR_protareInfos const &a_URR_protareInfos, std::size_t a_hashIndex, double a_temperature, double a_energy, bool a_sampling ) const {
447
448 std::size_t length = static_cast<std::size_t>( m_protares.size( ) );
449 double cross_section = 0.0;
450
451 for( std::size_t i1 = 0; i1 < length; ++i1 ) {
452 std::size_t numberOfReactions = m_protares[i1]->numberOfReactions( );
453
454 if( a_reactionIndex < numberOfReactions ) {
455 cross_section = m_protares[i1]->reactionCrossSection( a_reactionIndex, a_URR_protareInfos, a_hashIndex, a_temperature, a_energy, a_sampling );
456 break;
457 }
458 a_reactionIndex -= numberOfReactions;
459 }
460
461 return( cross_section );
462}
463
464/* *********************************************************************************************************//**
465 * Returns the cross section for reaction at index *a_reactionIndex*.
466 *
467 * @param a_reactionIndex [in] The index of the reaction.
468 * @param a_URR_protareInfos [in] URR information.
469 * @param a_temperature [in] The target temperature.
470 * @param a_energy [in] The projectile energy.
471 *
472 * @return The total cross section.
473 ***********************************************************************************************************/
474
475LUPI_HOST_DEVICE double ProtareComposite::reactionCrossSection( std::size_t a_reactionIndex, URR_protareInfos const &a_URR_protareInfos, double a_temperature, double a_energy ) const {
476
477 std::size_t length = static_cast<std::size_t>( m_protares.size( ) );
478 double cross_section = 0.0;
479
480 for( std::size_t i1 = 0; i1 < length; ++i1 ) {
481 std::size_t numberOfReactions = m_protares[i1]->numberOfReactions( );
482
483 if( a_reactionIndex < numberOfReactions ) {
484 cross_section = m_protares[i1]->reactionCrossSection( a_reactionIndex, a_URR_protareInfos, a_temperature, a_energy );
485 break;
486 }
487 a_reactionIndex -= numberOfReactions;
488 }
489
490 return( cross_section );
491}
492
493/* *********************************************************************************************************//**
494 * Returns the total deposition energy.
495 *
496 * @param a_hashIndex [in] The cross section hash index.
497 * @param a_temperature [in] The target temperature.
498 * @param a_energy [in] The projectile energy.
499 *
500 * @return The total deposition energy.
501 ***********************************************************************************************************/
502
503LUPI_HOST_DEVICE double ProtareComposite::depositionEnergy( std::size_t a_hashIndex, double a_temperature, double a_energy ) const {
504
505 std::size_t length = static_cast<std::size_t>( m_protares.size( ) );
506 double deposition_energy = 0.0;
507
508 for( std::size_t i1 = 0; i1 < length; ++i1 ) deposition_energy += m_protares[i1]->depositionEnergy( a_hashIndex, a_temperature, a_energy );
509
510 return( deposition_energy );
511}
512
513/* *********************************************************************************************************//**
514 * Returns the total deposition momentum.
515 *
516 * @param a_hashIndex [in] The cross section hash index.
517 * @param a_temperature [in] The target temperature.
518 * @param a_energy [in] The projectile energy.
519 *
520 * @return The total deposition momentum.
521 ***********************************************************************************************************/
522
523LUPI_HOST_DEVICE double ProtareComposite::depositionMomentum( std::size_t a_hashIndex, double a_temperature, double a_energy ) const {
524
525 std::size_t length = static_cast<std::size_t>( m_protares.size( ) );
526 double deposition_momentum = 0.0;
527
528 for( std::size_t i1 = 0; i1 < length; ++i1 ) deposition_momentum += m_protares[i1]->depositionMomentum( a_hashIndex, a_temperature, a_energy );
529
530 return( deposition_momentum );
531}
532
533/* *********************************************************************************************************//**
534 * Returns the total production energy.
535 *
536 * @param a_hashIndex [in] The cross section hash index.
537 * @param a_temperature [in] The target temperature.
538 * @param a_energy [in] The projectile energy.
539 *
540 * @return The total production energy.
541 ***********************************************************************************************************/
542
543LUPI_HOST_DEVICE double ProtareComposite::productionEnergy( std::size_t a_hashIndex, double a_temperature, double a_energy ) const {
544
545 std::size_t length = static_cast<std::size_t>( m_protares.size( ) );
546 double production_energy = 0.0;
547
548 for( std::size_t i1 = 0; i1 < length; ++i1 ) production_energy += m_protares[i1]->productionEnergy( a_hashIndex, a_temperature, a_energy );
549
550 return( production_energy );
551}
552
553/* *********************************************************************************************************//**
554 * Returns the gain for particle with index *a_particleIndex*.
555 *
556 * @param a_hashIndex [in] The continuous energy hash or multi-group index.
557 * @param a_temperature [in] The temperature of the target.
558 * @param a_energy [in] The projectile energy.
559 * @param a_particleIndex [in] The index of the particle whose gain is to be returned.
560 *
561 * @return [in] A vector of the length of the number of multi-group groups.
562 ***********************************************************************************************************/
563
564LUPI_HOST_DEVICE double ProtareComposite::gain( std::size_t a_hashIndex, double a_temperature, double a_energy, int a_particleIndex ) const {
565
566 std::size_t length = static_cast<std::size_t>( m_protares.size( ) );
567 double gain1 = m_protares[0]->gain( a_hashIndex, a_temperature, a_energy, a_particleIndex );
568
569 for( std::size_t i1 = 1; i1 < length; ++i1 ) gain1 += m_protares[i1]->gain( a_hashIndex, a_temperature, a_energy, a_particleIndex );
570
571 return( gain1 );
572}
573
574/* *********************************************************************************************************//**
575 * Returns the gain for particle with intid *a_particleIntid*.
576 *
577 * @param a_hashIndex [in] The continuous energy hash or multi-group index.
578 * @param a_temperature [in] The temperature of the target.
579 * @param a_energy [in] The projectile energy.
580 * @param a_particleIntid [in] The intid of the particle whose gain is to be returned.
581 *
582 * @return [in] A vector of the length of the number of multi-group groups.
583 ***********************************************************************************************************/
584
585LUPI_HOST_DEVICE double ProtareComposite::gainViaIntid( std::size_t a_hashIndex, double a_temperature, double a_energy, int a_particleIntid ) const {
586
587 std::size_t length = static_cast<std::size_t>( m_protares.size( ) );
588 double gain1 = m_protares[0]->gainViaIntid( a_hashIndex, a_temperature, a_energy, a_particleIntid );
589
590 for( std::size_t i1 = 1; i1 < length; ++i1 ) gain1 += m_protares[i1]->gainViaIntid( a_hashIndex, a_temperature, a_energy, a_particleIntid );
591
592 return( gain1 );
593}
594
595/* *********************************************************************************************************//**
596 * This method serializes *this* for broadcasting as needed for MPI and GPUs. The method can count the number of required
597 * bytes, pack *this* or unpack *this* depending on *a_mode*.
598 *
599 * @param a_buffer [in] The buffer to read or write data to depending on *a_mode*.
600 * @param a_mode [in] Specifies the action of this method.
601 ***********************************************************************************************************/
602
604
605 std::size_t vectorSize = m_protares.size( );
606 int vectorSizeInt = static_cast<int>( vectorSize );
607 LUPI::DataBuffer *workingBuffer = &a_buffer;
608
609 DATA_MEMBER_SIZE_T( m_numberOfReactions, a_buffer, a_mode );
610 DATA_MEMBER_SIZE_T( m_numberOfOrphanProducts, a_buffer, a_mode );
611 DATA_MEMBER_DOUBLE( m_minimumEnergy, a_buffer, a_mode );
612 DATA_MEMBER_DOUBLE( m_maximumEnergy, a_buffer, a_mode );
613
614 DATA_MEMBER_INT( vectorSizeInt, *workingBuffer, a_mode );
615 vectorSize = static_cast<std::size_t>( vectorSizeInt );
616
617 if( a_mode == LUPI::DataBuffer::Mode::Unpack ) {
618 m_protares.resize( vectorSize, &(workingBuffer->m_placement) );
619 for( std::size_t vectorIndex = 0; vectorIndex < vectorSize; ++vectorIndex ) {
620 if( workingBuffer->m_placement != nullptr ) {
621 m_protares[vectorIndex] = new(workingBuffer->m_placement) ProtareSingle;
622 workingBuffer->incrementPlacement( sizeof( ProtareSingle ) ); }
623 else {
624 m_protares[vectorIndex] = new ProtareSingle;
625 }
626 }
627 }
628 if( a_mode == LUPI::DataBuffer::Mode::Memory ) {
629 a_buffer.m_placement += m_protares.internalSize();
630 a_buffer.incrementPlacement( sizeof( ProtareSingle ) * vectorSize );
631 }
632
633 for( std::size_t i1 = 0; i1 < vectorSize; ++i1 ) {
634 m_protares[i1]->serializeCommon( a_buffer, a_mode );
635 m_protares[i1]->serialize2( a_buffer, a_mode );
636 }
637}
638
639}
#define DATA_MEMBER_SIZE_T(member, buf, mode)
#define DATA_MEMBER_DOUBLE(member, buf, mode)
#define DATA_MEMBER_INT( member, buf, mode)
#define LUPI_HOST_DEVICE
#define LUPI_THROW(arg)
#define LUPI_HOST
#define LUPI_maybeUnused
std::vector< Protare * > & protares()
Definition GIDI.hpp:4914
LUPI_HOST_DEVICE void incrementPlacement(std::size_t a_delta)
LUPI_HOST_DEVICE Reaction const * reaction(std::size_t a_index) const
LUPI_HOST_DEVICE ProtareSingle const * protareWithReaction(std::size_t a_index) const
LUPI_HOST_DEVICE Vector< double > temperatures(std::size_t a_index=0) const
LUPI_HOST_DEVICE Reaction const * orphanProduct(std::size_t a_index) const
LUPI_HOST void setUserParticleIndexViaIntid2(int a_particleIntid, int a_userParticleIndex)
LUPI_HOST_DEVICE double reactionCrossSection(std::size_t a_reactionIndex, URR_protareInfos const &a_URR_protareInfos, std::size_t a_hashIndex, double a_temperature, double a_energy, bool a_sampling=false) const
LUPI_HOST_DEVICE double depositionEnergy(std::size_t a_hashIndex, double a_temperature, double a_energy) const
LUPI_HOST_DEVICE double gain(std::size_t a_hashIndex, double a_temperature, double a_energy, int a_particleIndex) const
LUPI_HOST_DEVICE ProtareSingle const * protare(std::size_t a_index) const
LUPI_HOST_DEVICE double maximumEnergy() const
Definition MCGIDI.hpp:1778
LUPI_HOST_DEVICE double URR_domainMin() const
LUPI_HOST_DEVICE void serialize2(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
LUPI_HOST_DEVICE double crossSection(URR_protareInfos const &a_URR_protareInfos, std::size_t a_hashIndex, double a_temperature, double a_energy, bool a_sampling=false) const
LUPI_HOST_DEVICE bool hasFission() const
LUPI_HOST_DEVICE bool hasURR_probabilityTables() const
LUPI_HOST_DEVICE double threshold(std::size_t a_index) const
LUPI_HOST_DEVICE double gainViaIntid(std::size_t a_hashIndex, double a_temperature, double a_energy, int a_particleIntid) const
LUPI_HOST_DEVICE double productionEnergy(std::size_t a_hashIndex, double a_temperature, double a_energy) const
LUPI_HOST_DEVICE double URR_domainMax() const
LUPI_HOST_DEVICE bool reactionHasURR_probabilityTables(std::size_t a_index) const
LUPI_HOST_DEVICE void crossSectionVector(double a_temperature, double a_userFactor, std::size_t a_numberAllocated, double *a_crossSectionVector) const
LUPI_HOST_DEVICE double depositionMomentum(std::size_t a_hashIndex, double a_temperature, double a_energy) const
Vector< ProtareSingle * > protares() const
Definition MCGIDI.hpp:1766
LUPI_HOST void setUserParticleIndex2(int a_particleIndex, int a_userParticleIndex)
LUPI_HOST_DEVICE std::size_t numberOfReactions() const
Definition MCGIDI.hpp:1786
LUPI_HOST_DEVICE bool hasIncoherentDoppler() const
LUPI_HOST_DEVICE double minimumEnergy() const
Definition MCGIDI.hpp:1777
LUPI_HOST_DEVICE std::size_t numberOfReactions() const
Definition MCGIDI.hpp:1694
LUPI_HOST Vector< int > const & productIndices(bool a_transportablesOnly) const
friend ProtareSingle
Definition MCGIDI.hpp:1599
LUPI_HOST Vector< int > const & productIntids(bool a_transportablesOnly) const
LUPI_HOST_DEVICE Protare(ProtareType a_protareType)
std::vector< Styles::TemperatureInfo > TemperatureInfos
Definition GIDI.hpp:3440
std::set< std::size_t > ExcludeReactionsSet
Definition GIDI.hpp:47
Simple C++ string class, useful as replacement for std::string if this cannot be used,...
Definition MCGIDI.hpp:43
LUPI_HOST void addVectorItemsToSet(Vector< int > const &a_from, std::set< int > &a_to)
LUPI_HOST Protare * protareFromGIDIProtare(LUPI::StatusMessageReporting &a_smr, GIDI::Protare const &a_protare, PoPI::Database const &a_pops, Transporting::MC &a_settings, GIDI::Transporting::Particles const &a_particles, DomainHash const &a_domainHash, GIDI::Styles::TemperatureInfos const &a_temperatureInfos, GIDI::ExcludeReactionsSet const &a_reactionsToExclude, std::size_t a_reactionsToExcludeOffset=0, bool a_allowFixedGrid=true)
ProtareType
Definition MCGIDI.hpp:79