Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
MCGIDI_protare.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/* *********************************************************************************************************//**
17 * Returns the proper **MCGIDI** protare base on the type of **GIDI** protare.
18 *
19 * @param a_smr [Out] If errors are not to be thrown, then the error is reported via this instance.
20 * @param a_protare [in] The GIDI::Protare whose data is to be used to construct *this*.
21 * @param a_pops [in] A PoPs Database instance used to get particle intids and possibly other particle information.
22 * @param a_settings [in] Used to pass user options to the *this* to instruct it which data are desired.
23 * @param a_particles [in] List of transporting particles and their information (e.g., multi-group boundaries and fluxes).
24 * @param a_domainHash [in] The hash data used when looking up a cross section.
25 * @param a_temperatureInfos [in] The list of temperature data to extract from *a_protare*.
26 * @param a_reactionsToExclude [in] A list of reaction to not include in the MCGIDI::Protare.
27 * @param a_reactionsToExcludeOffset [in] The starting index for the reactions in this ProtareSingle.
28 * @param a_allowFixedGrid [in] For internal (i.e., MCGIDI) use only. Users must use the default value.
29 ***********************************************************************************************************/
30
32 GIDI::Transporting::Particles const &a_particles, DomainHash const &a_domainHash, GIDI::Styles::TemperatureInfos const &a_temperatureInfos,
33 GIDI::ExcludeReactionsSet const &a_reactionsToExclude, std::size_t a_reactionsToExcludeOffset, bool a_allowFixedGrid ) {
34
35 Protare *protare( nullptr );
36
37 if( a_protare.protareType( ) == GIDI::ProtareType::single ) {
38 protare = new ProtareSingle( a_smr, static_cast<GIDI::ProtareSingle const &>( a_protare ), a_pops, a_settings, a_particles, a_domainHash,
39 a_temperatureInfos, a_reactionsToExclude, a_reactionsToExcludeOffset, a_allowFixedGrid ); }
40 else if( a_protare.protareType( ) == GIDI::ProtareType::composite ) {
41 protare = new ProtareComposite( a_smr, static_cast<GIDI::ProtareComposite const &>( a_protare ), a_pops, a_settings, a_particles, a_domainHash,
42 a_temperatureInfos, a_reactionsToExclude, a_reactionsToExcludeOffset, false ); }
43 else if( a_protare.protareType( ) == GIDI::ProtareType::TNSL ) {
44 protare = new ProtareTNSL( a_smr, static_cast<GIDI::ProtareTNSL const &>( a_protare ), a_pops, a_settings, a_particles, a_domainHash,
45 a_temperatureInfos, a_reactionsToExclude, a_reactionsToExcludeOffset, false );
46 }
47
48 return( protare );
49}
50
51/*! \class Protare
52 * Base class for the *MCGIDI* protare classes.
53 */
54
55/* *********************************************************************************************************//**
56 * @param a_protareType [in] The enum for the type of Protare (i.e., single, composite or TNSL).
57 *
58 * Default constructor used when broadcasting a Protare as needed by MPI or GPUs.
59 ***********************************************************************************************************/
60
62 m_protareType( a_protareType ),
63 m_projectileID( ),
64 m_projectileIntid( -1 ),
65 m_projectileIndex( -1 ),
66 m_projectileUserIndex( -1 ),
67 m_projectileMass( 0.0 ),
68 m_projectileExcitationEnergy( 0.0 ),
69
70 m_targetID( ),
71 m_targetIntid( -1 ),
72 m_targetIndex( -1 ),
73 m_targetUserIndex( -1 ),
74 m_targetMass( 0.0 ),
75 m_targetExcitationEnergy( 0.0 ),
76
77 m_neutronIndex( -1 ),
78 m_userNeutronIndex( -1 ),
79 m_photonIndex( -1 ),
80 m_userPhotonIndex( -1 ),
81
82 m_evaluation( ),
83 m_projectileFrame( GIDI::Frame::lab ),
84
85 m_isTNSL_ProtareSingle( false ) {
86
87}
88
89/* *********************************************************************************************************//**
90 * Default base Protare constructor.
91 *
92 * @param a_protareType [in] The enum for the type of Protare (i.e., single, composite or TNSL).
93 * @param a_protare [in] The GIDI::Protare whose data is to be used to construct *this*.
94 * @param a_settings [in] Used to pass user options to the *this* to instruct it which data are desired.
95 ***********************************************************************************************************/
96
97LUPI_HOST Protare::Protare( ProtareType a_protareType, GIDI::Protare const &a_protare, LUPI_maybeUnused Transporting::MC const &a_settings, PoPI::Database const &a_pops ) :
98 m_protareType( a_protareType ),
99 m_projectileID( a_protare.projectile( ).ID( ).c_str( ) ),
100 m_projectileIntid( -1 ),
101 m_projectileIndex( MCGIDI_popsIndex( a_pops, m_projectileID.c_str( ) ) ),
102 m_projectileUserIndex( -1 ),
103 m_projectileMass( a_protare.projectile( ).mass( "MeV/c**2" ) ), // Includes nuclear excitation energy.
104 m_projectileExcitationEnergy( a_protare.projectile( ).excitationEnergy( ).value( ) ),
105
106 m_targetID( a_protare.target( ).ID( ).c_str( ) ),
107 m_targetIntid( -1 ),
108 m_targetIndex( -1 ),
109 m_targetUserIndex( -1 ),
110 m_targetMass( a_protare.target( ).mass( "MeV/c**2" ) ), // Includes nuclear excitation energy.
111 m_targetExcitationEnergy( a_protare.target( ).excitationEnergy( ).value( ) ),
112
113 m_neutronIndex( MCGIDI_popsIndex( a_pops, PoPI::IDs::neutron ) ),
114 m_userNeutronIndex( -1 ),
115 m_photonIndex( MCGIDI_popsIndex( a_pops, PoPI::IDs::photon ) ),
116 m_userPhotonIndex( -1 ),
117
118 m_evaluation( a_protare.evaluation( ).c_str( ) ),
119 m_projectileFrame( a_protare.projectileFrame( ) ),
120 m_productIntids( 0 ),
121 m_userProductIndices( 0 ),
122 m_productIntidsTransportable( 0 ),
123 m_userProductIndicesTransportable( 0 ),
124 m_isTNSL_ProtareSingle( a_protare.isTNSL_ProtareSingle( ) ) {
125
126 PoPI::Database const &pops = a_protare.protare( 0 )->internalPoPs( );
127 m_projectileIntid = MCGIDI_popsIntid( pops, m_projectileID.c_str( ) );
128
129 if( a_protare.protare( 0 )->interaction( ) != GIDI_MapInteractionTNSLChars ) {
130 m_targetIntid = MCGIDI_popsIntid( pops, m_targetID.c_str( ) );
131 m_targetIndex = MCGIDI_popsIndex( a_pops, m_targetID.c_str( ) );
132 }
133}
134
135/* *********************************************************************************************************//**
136 ***********************************************************************************************************/
137
141
142/* *********************************************************************************************************//**
143 * Sets *this* members *m_productIntids* and *m_productIntidsTransportable* to *a_intids* and *a_transportableIntids* respectively.
144 * And, sets *this* members *m_productIndices* and *m_productIndicesTransportable* to *a_indices* and *a_transportableIndices* respectively.
145 *
146 * @param a_intids [out] The list of intids for the outgoing particles (i.e., products).
147 * @param a_transportableIntids [in] The list of transportable intids for the outgoing particles (i.e., products).
148 * @param a_indices [in] The list of indices for the outgoing particles (i.e., products).
149 * @param a_transportableIndices [in] The list of transportable indices for the outgoing particles (i.e., products).
150 ***********************************************************************************************************/
151
152LUPI_HOST void Protare::productIntidsAndIndices( std::set<int> const &a_intids, std::set<int> const &a_transportableIntids,
153 std::set<int> const &a_indices, std::set<int> const &a_transportableIndices ) {
154
155 m_productIntids.reserve( a_intids.size( ) );
156 m_userProductIndices.reserve( a_intids.size( ) );
157 for( std::set<int>::const_iterator iter = a_intids.begin( ); iter != a_intids.end( ); ++iter ) {
158 m_productIntids.push_back( *iter );
159 m_userProductIndices.push_back( -1 );
160 }
161
162 m_productIndices.reserve( a_indices.size( ) );
163 for( auto iter = a_indices.begin( ); iter != a_indices.end( ); ++iter ) m_productIndices.push_back( *iter );
164
165 m_productIntidsTransportable.reserve( a_transportableIntids.size( ) );
166 m_userProductIndicesTransportable.reserve( a_transportableIntids.size( ) );
167 for( std::set<int>::const_iterator iter = a_transportableIntids.begin( ); iter != a_transportableIntids.end( ); ++iter ) {
168 m_productIntidsTransportable.push_back( *iter );
169 m_userProductIndicesTransportable.push_back( -1 );
170 }
171
172 m_productIndicesTransportable.reserve( a_transportableIndices.size( ) );
173 for( auto iter = a_transportableIndices.begin( ); iter != a_transportableIndices.end( ); ++iter ) m_productIndicesTransportable.push_back( *iter );
174}
175
176/* *********************************************************************************************************//**
177 * Returns the list product intids. If *a_transportablesOnly* is true, the list only includes transportable particle.
178 *
179 * @param a_transportablesOnly [in] If **true**, a reference to *m_productIntidsTransportable* is returned; otherwise a reference to *m_productIntids* is returned.
180 ***********************************************************************************************************/
181
182LUPI_HOST Vector<int> const &Protare::productIntids( bool a_transportablesOnly ) const {
183
184 if( a_transportablesOnly ) return( m_productIntidsTransportable );
185 return( m_productIntids );
186}
187
188/* *********************************************************************************************************//**
189 * Returns the list product indices. If *a_transportablesOnly* is true, the list only includes transportable particle.
190 *
191 * @param a_transportablesOnly [in] If **true**, a reference to *m_productIndicesTransportable* is returned; otherwise a reference to *m_productIndices* is returned.
192 ***********************************************************************************************************/
193
194LUPI_HOST Vector<int> const &Protare::productIndices( bool a_transportablesOnly ) const {
195
196 if( a_transportablesOnly ) return( m_productIndicesTransportable );
197 return( m_productIndices );
198}
199
200/* *********************************************************************************************************//**
201 * Returns the list user product indices. If *a_transportablesOnly* is true, the list only includes transportable particle.
202 *
203 * @param a_transportablesOnly [in] If **true**, a reference to *m_userProductIndicesTransportable* is returned; otherwise a reference to *m_userProductIndices* is returned.
204 ***********************************************************************************************************/
205
206LUPI_HOST Vector<int> const &Protare::userProductIndices( bool a_transportablesOnly ) const {
207
208 if( a_transportablesOnly ) return( m_userProductIndicesTransportable );
209 return( m_userProductIndices );
210}
211
212/* *********************************************************************************************************//**
213 * Updates the m_userParticleIndex to *a_userParticleIndex* for all particles with PoPs index *a_particleIndex*.
214 *
215 * @param a_particleIndex [in] The PoPs index of the particle whose user index is to be set.
216 * @param a_userParticleIndex [in] The particle index specified by the user.
217 ***********************************************************************************************************/
218
219LUPI_HOST void Protare::setUserParticleIndex( int a_particleIndex, int a_userParticleIndex ) {
220
221 if( m_projectileIndex == a_particleIndex ) m_projectileUserIndex = a_userParticleIndex;
222 if( m_targetIndex == a_particleIndex ) m_targetUserIndex = a_userParticleIndex;
223
224 if( m_photonIndex == a_particleIndex ) m_userPhotonIndex = a_userParticleIndex;
225
226 for( std::size_t i1 = 0; i1 < m_productIndices.size( ); ++i1 ) {
227 if( m_productIndices[i1] == a_particleIndex ) m_userProductIndices[i1] = a_userParticleIndex;
228 }
229
230 for( std::size_t i1 = 0; i1 < m_productIndicesTransportable.size( ); ++i1 ) {
231 if( m_productIndicesTransportable[i1] == a_particleIndex ) m_userProductIndicesTransportable[i1] = a_userParticleIndex;
232 }
233
234 switch( m_protareType ) {
236 static_cast<ProtareSingle *>( this )->setUserParticleIndex2( a_particleIndex, a_userParticleIndex );
237 break;
239 static_cast<ProtareComposite *>( this )->setUserParticleIndex2( a_particleIndex, a_userParticleIndex );
240 break;
242 static_cast<ProtareTNSL *>( this )->setUserParticleIndex2( a_particleIndex, a_userParticleIndex );
243 break;
244 }
245}
246
247/* *********************************************************************************************************//**
248 * Updates the m_userParticleIndex to *a_userParticleIndex* for all particles with PoPs intid *a_particleIntid*.
249 *
250 * @param a_particleIntid [in] The PoPs intid of the particle whose user index is to be set.
251 * @param a_userParticleIndex [in] The particle index specified by the user.
252 ***********************************************************************************************************/
253
254LUPI_HOST void Protare::setUserParticleIndexViaIntid( int a_particleIntid, int a_userParticleIndex ) {
255
256 if( m_projectileIntid == a_particleIntid ) m_projectileUserIndex = a_userParticleIndex;
257 if( m_targetIntid == a_particleIntid ) m_targetUserIndex = a_userParticleIndex;
258
259 if( PoPI::Intids::photon == a_particleIntid ) m_userPhotonIndex = a_userParticleIndex;
260
261 for( std::size_t i1 = 0; i1 < m_productIntids.size( ); ++i1 ) {
262 if( m_productIntids[i1] == a_particleIntid ) m_userProductIndices[i1] = a_userParticleIndex;
263 }
264
265 for( std::size_t i1 = 0; i1 < m_productIntidsTransportable.size( ); ++i1 ) {
266 if( m_productIntidsTransportable[i1] == a_particleIntid ) m_userProductIndicesTransportable[i1] = a_userParticleIndex;
267 }
268
269 switch( m_protareType ) {
271 static_cast<ProtareSingle *>( this )->setUserParticleIndexViaIntid2( a_particleIntid, a_userParticleIndex );
272 break;
274 static_cast<ProtareComposite *>( this )->setUserParticleIndexViaIntid2( a_particleIntid, a_userParticleIndex );
275 break;
277 static_cast<ProtareTNSL *>( this )->setUserParticleIndexViaIntid2( a_particleIntid, a_userParticleIndex );
278 break;
279 }
280}
281
282/* *********************************************************************************************************//**
283 * This method serializes *this* for broadcasting as needed for MPI and GPUs. The method can count the number of required
284 * bytes, pack *this* or unpack *this* depending on *a_mode*.
285 *
286 * @param a_buffer [in] The buffer to read or write data to depending on *a_mode*.
287 * @param a_mode [in] Specifies the action of this method.
288 ***********************************************************************************************************/
289
291
292 serializeCommon( a_buffer, a_mode );
293
294 switch( m_protareType ) {
296 static_cast<ProtareSingle *>( this )->serialize2( a_buffer, a_mode );
297 break;
299 static_cast<ProtareComposite *>( this )->serialize2( a_buffer, a_mode );
300 break;
302 static_cast<ProtareTNSL *>( this )->serialize2( a_buffer, a_mode );
303 break;
304 }
305}
306
307/* *********************************************************************************************************//**
308 * This method serializes *this* for broadcasting as needed for MPI and GPUs. The method can count the number of required
309 * bytes, pack *this* or unpack *this* depending on *a_mode*.
310 *
311 * @param a_buffer [in] The buffer to read or write data to depending on *a_mode*.
312 * @param a_mode [in] Specifies the action of this method.
313 ***********************************************************************************************************/
314
316
317 serializeCommon( a_buffer, a_mode );
318
319 switch( m_protareType ) {
321 static_cast<ProtareSingle *>( this )->serialize2( a_buffer, a_mode );
322 break;
324 LUPI_THROW( "Protare::serialize2:: Oops1, this should not happend." );
325 break;
327 LUPI_THROW( "Protare::serialize2:: Oops2, this should not happend." );
328 break;
329 }
330}
331
332/* *********************************************************************************************************//**
333 * This method serializes *this* for broadcasting as needed for MPI and GPUs. The method can count the number of required
334 * bytes, pack *this* or unpack *this* depending on *a_mode*.
335 *
336 * @param a_buffer [in] The buffer to read or write data to depending on *a_mode*.
337 * @param a_mode [in] Specifies the action of this method.
338 ***********************************************************************************************************/
339
341
342 int protareType1 = 0;
343 if( a_mode != LUPI::DataBuffer::Mode::Unpack ) {
344 switch( m_protareType ) {
346 break;
348 protareType1 = 1;
349 break;
350 case ProtareType::TNSL :
351 protareType1 = 2;
352 break;
353 }
354 }
355 DATA_MEMBER_INT( protareType1, a_buffer, a_mode );
356 if( a_mode == LUPI::DataBuffer::Mode::Unpack ) {
357 switch( protareType1 ) {
358 case 0 :
359 m_protareType = ProtareType::single;
360 break;
361 case 1 :
362 m_protareType = ProtareType::composite;
363 break;
364 case 2 :
365 m_protareType = ProtareType::TNSL;
366 break;
367 }
368 }
369
370 DATA_MEMBER_STRING( m_projectileID, a_buffer, a_mode );
371 DATA_MEMBER_INT( m_projectileIntid, a_buffer, a_mode );
372 DATA_MEMBER_INT( m_projectileIndex, a_buffer, a_mode );
373 DATA_MEMBER_INT( m_projectileUserIndex, a_buffer, a_mode );
374 DATA_MEMBER_DOUBLE( m_projectileMass, a_buffer, a_mode );
375 DATA_MEMBER_DOUBLE( m_projectileExcitationEnergy, a_buffer, a_mode );
376
377 DATA_MEMBER_STRING( m_targetID, a_buffer, a_mode );
378 DATA_MEMBER_INT( m_targetIntid, a_buffer, a_mode );
379 DATA_MEMBER_INT( m_targetIndex, a_buffer, a_mode );
380 DATA_MEMBER_INT( m_targetUserIndex, a_buffer, a_mode );
381 DATA_MEMBER_DOUBLE( m_targetMass, a_buffer, a_mode );
382 DATA_MEMBER_DOUBLE( m_targetExcitationEnergy, a_buffer, a_mode );
383
384 DATA_MEMBER_INT( m_photonIndex, a_buffer, a_mode );
385 DATA_MEMBER_INT( m_userPhotonIndex, a_buffer, a_mode );
386
387 DATA_MEMBER_STRING( m_evaluation, a_buffer, a_mode );
388
389 int frame = 0;
390 if( m_projectileFrame == GIDI::Frame::centerOfMass ) frame = 1;
391 DATA_MEMBER_INT( frame, a_buffer, a_mode );
392 m_projectileFrame = GIDI::Frame::lab;
393 if( frame == 1 ) m_projectileFrame = GIDI::Frame::centerOfMass;
394
395 DATA_MEMBER_VECTOR_INT( m_productIntids, a_buffer, a_mode );
396 DATA_MEMBER_VECTOR_INT( m_productIndices, a_buffer, a_mode );
397 DATA_MEMBER_VECTOR_INT( m_userProductIndices, a_buffer, a_mode );
398
399 DATA_MEMBER_VECTOR_INT( m_productIntidsTransportable, a_buffer, a_mode );
400 DATA_MEMBER_VECTOR_INT( m_productIndicesTransportable, a_buffer, a_mode );
401 DATA_MEMBER_VECTOR_INT( m_userProductIndicesTransportable, a_buffer, a_mode );
402
403 DATA_MEMBER_CAST( m_isTNSL_ProtareSingle, a_buffer, a_mode, bool );
404}
405
406/* *********************************************************************************************************//**
407 * Returns the number of memory bytes used by *this*.
408 *
409 * @return The number of bytes used by *this*.
410 ***********************************************************************************************************/
411
413
414 long sizeOf1 = 0;
415
416 switch( m_protareType ) {
418 sizeOf1 = static_cast<ProtareSingle const *>( this )->sizeOf2( );
419 break;
421 sizeOf1 = static_cast<ProtareComposite const *>( this )->sizeOf2( );
422 break;
423 case ProtareType::TNSL:
424 sizeOf1 = static_cast<ProtareTNSL const *>( this )->sizeOf2( );
425 break;
426 }
427
428 return( sizeOf1 );
429}
430
431/* *********************************************************************************************************//**
432 * This method counts the number of bytes of memory allocated by *this*.
433 * This is an improvement to the internalSize() method of getting memory size.
434 ***********************************************************************************************************/
436
438 // Written this way for debugger to modify buf.m_placementStart here for easier double checking.
439 buf.m_placement = buf.m_placementStart + sizeOf();
441 return( ( buf.m_placement - buf.m_placementStart ) + ( buf.m_sharedPlacement - buf.m_sharedPlacementStart ) );
442}
443
444/* *********************************************************************************************************//**
445 * This method counts the number of bytes of memory allocated by *this* and puts it into a_totalMemory.
446 * If shared memory is used, the size of shared memory is a_sharedMemory. If using shared memory,
447 * the host code only needs to allocate (a_totalMemory - a_sharedMemory) in main memory.
448 ***********************************************************************************************************/
449LUPI_HOST_DEVICE void Protare::incrementMemorySize( long &a_totalMemory, long &a_sharedMemory ) {
450
451 LUPI::DataBuffer buf; // Written this way for debugger to modify buf.m_placementStart here for easier double checking.
452
453 buf.m_placement = buf.m_placementStart + sizeOf( );
455 a_totalMemory += buf.m_placement - buf.m_placementStart;
456 a_sharedMemory += buf.m_sharedPlacement - buf.m_sharedPlacementStart;
457}
458
459/* *********************************************************************************************************//**
460 * Returns the number of protares contained in *this*.
461 *
462 * @return Integer number of protares.
463 ***********************************************************************************************************/
464
466
467 std::size_t numberOfProtares2 = 0;
468
469 switch( protareType( ) ) {
471 numberOfProtares2 = static_cast<ProtareSingle const *>( this )->numberOfProtares( );
472 break;
474 numberOfProtares2 = static_cast<ProtareComposite const *>( this )->numberOfProtares( );
475 break;
477 numberOfProtares2 = static_cast<ProtareTNSL const *>( this )->numberOfProtares( );
478 break;
479 }
480
481 return( numberOfProtares2 );
482}
483
484/* *********************************************************************************************************//**
485 * Returns the const pointer representing the protare at index *a_index*.
486 *
487 * @param a_index [in] Index of protare in *this*.
488 *
489 * @return Returns the const pointer representing the protare.
490 ***********************************************************************************************************/
491
492LUPI_HOST_DEVICE ProtareSingle const *Protare::protare( std::size_t a_index ) const {
493
494 ProtareSingle const *protare1 = nullptr;
495
496 switch( protareType( ) ) {
498 protare1 = static_cast<ProtareSingle const *>( this )->protare( a_index );
499 break;
501 protare1 = static_cast<ProtareComposite const *>( this )->protare( a_index );
502 break;
504 protare1 = static_cast<ProtareTNSL const *>( this )->protare( a_index );
505 break;
506 }
507
508 return( protare1 );
509}
510
511/* *********************************************************************************************************//**
512 * Returns the pointer representing the protare at index *a_index*.
513 *
514 * @param a_index [in] Index of protare in *this*.
515 *
516 * @return Returns the pointer representing the protare.
517 ***********************************************************************************************************/
518
520
521 ProtareSingle *protare1 = nullptr;
522
523 switch( protareType( ) ) {
525 protare1 = static_cast<ProtareSingle *>( this )->protare( a_index );
526 break;
528 protare1 = static_cast<ProtareComposite *>( this )->protare( a_index );
529 break;
531 protare1 = static_cast<ProtareTNSL *>( this )->protare( a_index );
532 break;
533 }
534
535 return( protare1 );
536}
537
538/* *********************************************************************************************************//**
539 * Returns the pointer to the **ProtareSingle** that contains the (a_index - 1)th reaction.
540 *
541 * @param a_index [in] Index of the reaction.
542 *
543 * @return Pointer to the requested protare or nullptr if invalid *a_index*..
544 ***********************************************************************************************************/
545
547
548 ProtareSingle const *protare1 = nullptr;
549
550 switch( protareType( ) ) {
552 protare1 = static_cast<ProtareSingle const *>( this )->protareWithReaction( a_index );
553 break;
555 protare1 = static_cast<ProtareComposite const *>( this )->protareWithReaction( a_index );
556 break;
558 protare1 = static_cast<ProtareTNSL const *>( this )->protareWithReaction( a_index );
559 break;
560 }
561
562 return( protare1 );
563}
564
565/* *********************************************************************************************************//**
566 * Returns the minimum cross section domain for all reaction..
567 *
568 * @return Returns the minimum cross section domain.
569 ***********************************************************************************************************/
570
572
573 double minimumEnergy1 = 0.0;
574
575 switch( protareType( ) ) {
577 minimumEnergy1 = static_cast<ProtareSingle const *>( this )->minimumEnergy( );
578 break;
580 minimumEnergy1 = static_cast<ProtareComposite const *>( this )->minimumEnergy( );
581 break;
583 minimumEnergy1 = static_cast<ProtareTNSL const *>( this )->minimumEnergy( );
584 break;
585 }
586
587 return( minimumEnergy1 );
588}
589
590/* *********************************************************************************************************//**
591 * Returns the maximum cross section domain for all reaction..
592 *
593 * @return Returns the maximum cross section domain.
594 ***********************************************************************************************************/
595
597
598 double maximumEnergy1 = 0.0;
599
600 switch( protareType( ) ) {
602 maximumEnergy1 = static_cast<ProtareSingle const *>( this )->maximumEnergy( );
603 break;
605 maximumEnergy1 = static_cast<ProtareComposite const *>( this )->maximumEnergy( );
606 break;
608 maximumEnergy1 = static_cast<ProtareTNSL const *>( this )->maximumEnergy( );
609 break;
610 }
611
612 return( maximumEnergy1 );
613}
614
615/* *********************************************************************************************************//**
616 * Returns the list of temperatures for the requested ProtareSingle.
617 *
618 * @param a_index [in] Index of the reqested ProtareSingle.
619 *
620 * @return Vector of doubles.
621 ***********************************************************************************************************/
622
624
625 ProtareSingle const *protareSingle = protare( a_index );
626 return( protareSingle->temperatures( ) );
627}
628
629/* *********************************************************************************************************//**
630 * Returns the projectile's multi-group boundaries that were read from the file (i.e., pre-collapse).
631 *
632 * @param a_index [in] Index of the reqested ProtareSingle.
633 *
634 * @return Vector of doubles.
635 ***********************************************************************************************************/
636
638
639 ProtareSingle const *protareSingle = protare( 0 );
640 return( protareSingle->projectileMultiGroupBoundaries( ) );
641}
642
643
644/* *********************************************************************************************************//**
645 * Returns the projectile's collapsed multi-group boundaries (i.e., those used for transport).
646 *
647 * @param a_index [in] Index of the reqested ProtareSingle.
648 *
649 * @return Vector of doubles.
650 ***********************************************************************************************************/
651
653
654 ProtareSingle const *protareSingle = protare( 0 );
655 return( protareSingle->projectileMultiGroupBoundariesCollapsed( ) );
656}
657
658/* *********************************************************************************************************//**
659 * Returns the number of reactions of *this*.
660 *
661 * @return Number of reactions of *this*.
662 ***********************************************************************************************************/
663
665
666 std::size_t numberOfReactions1 = 0;
667
668 switch( protareType( ) ) {
670 numberOfReactions1 = static_cast<ProtareSingle const *>( this )->numberOfReactions( );
671 break;
673 numberOfReactions1 = static_cast<ProtareComposite const *>( this )->numberOfReactions( );
674 break;
676 numberOfReactions1 = static_cast<ProtareTNSL const *>( this )->numberOfReactions( );
677 break;
678 }
679
680 return( numberOfReactions1 );
681}
682
683/* *********************************************************************************************************//**
684 * Returns the reaction at index *a_index*.
685 *
686 * @param a_index [in] The index of the reaction to return.
687 *
688 * @return The reaction at index *a_index*.
689 ***********************************************************************************************************/
690
691LUPI_HOST_DEVICE Reaction const *Protare::reaction( std::size_t a_index ) const {
692
693 Reaction const *reaction1 = nullptr;
694
695 switch( protareType( ) ) {
697 reaction1 = static_cast<ProtareSingle const *>( this )->reaction( a_index );
698 break;
700 reaction1 = static_cast<ProtareComposite const *>( this )->reaction( a_index );
701 break;
703 reaction1 = static_cast<ProtareTNSL const *>( this )->reaction( a_index );
704 break;
705 }
706
707 return( reaction1 );
708}
709
710/* *********************************************************************************************************//**
711 * Returns the number of orphanProducts of *this*.
712 *
713 * @return Number of orphanProducts of *this*.
714 ***********************************************************************************************************/
715
717
718 std::size_t numberOfReactions1 = 0;
719
720 switch( protareType( ) ) {
722 numberOfReactions1 = static_cast<ProtareSingle const *>( this )->numberOfOrphanProducts( );
723 break;
725 numberOfReactions1 = static_cast<ProtareComposite const *>( this )->numberOfOrphanProducts( );
726 break;
728 numberOfReactions1 = static_cast<ProtareTNSL const *>( this )->numberOfOrphanProducts( );
729 break;
730 }
731
732 return( numberOfReactions1 );
733}
734
735/* *********************************************************************************************************//**
736 * Returns the orphanProduct at index *a_index*.
737 *
738 * @param a_index [in] The index of the orphanProduct to return.
739 *
740 * @return The orphanProduct at index *a_index*.
741 ***********************************************************************************************************/
742
743LUPI_HOST_DEVICE Reaction const *Protare::orphanProduct( std::size_t a_index ) const {
744
745 Reaction const *orphanProduct1 = nullptr;
746
747 switch( protareType( ) ) {
749 orphanProduct1 = static_cast<ProtareSingle const *>( this )->orphanProduct( a_index );
750 break;
752 orphanProduct1 = static_cast<ProtareComposite const *>( this )->orphanProduct( a_index );
753 break;
755 orphanProduct1 = static_cast<ProtareTNSL const *>( this )->orphanProduct( a_index );
756 break;
757 }
758
759 return( orphanProduct1 );
760}
761
762/* *********************************************************************************************************//**
763 * Returns true if *this* has a fission reaction and false otherwise.
764 *
765 * @return true is if *this* has a fission reaction and false otherwise.
766 ***********************************************************************************************************/
767
769
770 bool hasFission1 = false;
771
772 switch( protareType( ) ) {
774 hasFission1 = static_cast<ProtareSingle const *>( this )->hasFission( );
775 break;
777 hasFission1 = static_cast<ProtareComposite const *>( this )->hasFission( );
778 break;
780 hasFission1 = static_cast<ProtareTNSL const *>( this )->hasFission( );
781 break;
782 }
783
784 return( hasFission1 );
785}
786
787/* *********************************************************************************************************//**
788 * Returns true if *this* has a photoatomic incoherent doppler broadened reaction and false otherwise.
789 *
790 * @return true is if *this* has a specified reaction and false otherwise.
791 ***********************************************************************************************************/
792
794
795 bool hasIncoherentDoppler1 = false;
796
797 switch( protareType( ) ) {
799 hasIncoherentDoppler1 = static_cast<ProtareSingle const *>( this )->hasIncoherentDoppler( );
800 break;
802 hasIncoherentDoppler1 = static_cast<ProtareComposite const *>( this )->hasIncoherentDoppler( );
803 break;
805 hasIncoherentDoppler1 = static_cast<ProtareTNSL const *>( this )->hasIncoherentDoppler( );
806 break;
807 }
808
809 return( hasIncoherentDoppler1 );
810}
811
812/* *********************************************************************************************************//**
813 * Returns URR index of *this*.
814 *
815 * @return Integer URR index of *this*.
816 ***********************************************************************************************************/
817
819
820 int URR_index1 = 0;
821
822 switch( protareType( ) ) {
824 URR_index1 = static_cast<ProtareSingle const *>( this )->URR_index( );
825 break;
827 URR_index1 = static_cast<ProtareComposite const *>( this )->URR_index( );
828 break;
830 URR_index1 = static_cast<ProtareTNSL const *>( this )->URR_index( );
831 break;
832 }
833
834 return( URR_index1 );
835}
836
837/* *********************************************************************************************************//**
838 * Returns **true** if *this* has URR probability tables and **false** otherwise.
839 *
840 * @return boolean.
841 ***********************************************************************************************************/
842
844
845 bool hasURR_probabilityTables1 = false;
846
847 switch( protareType( ) ) {
849 hasURR_probabilityTables1 = static_cast<ProtareSingle const *>( this )->hasURR_probabilityTables( );
850 break;
852 hasURR_probabilityTables1 = static_cast<ProtareComposite const *>( this )->hasURR_probabilityTables( );
853 break;
855 hasURR_probabilityTables1 = static_cast<ProtareTNSL const *>( this )->hasURR_probabilityTables( );
856 break;
857 }
858
859 return( hasURR_probabilityTables1 );
860}
861
862/* *********************************************************************************************************//**
863 * Returns the minimum energy for the unresolved resonance region (URR) domain. If no URR data present, returns -1.
864 *
865 * @return Minimum energy for the unresolved resonance region (URR) domain.
866 ***********************************************************************************************************/
867
869
870 double URR_domainMin1 = 0.0;
871
872 switch( protareType( ) ) {
874 URR_domainMin1 = static_cast<ProtareSingle const *>( this )->URR_domainMin( );
875 break;
877 URR_domainMin1 = static_cast<ProtareComposite const *>( this )->URR_domainMin( );
878 break;
880 URR_domainMin1 = static_cast<ProtareTNSL const *>( this )->URR_domainMin( );
881 break;
882 }
883
884 return( URR_domainMin1 );
885}
886
887/* *********************************************************************************************************//**
888 * Returns the maximum energy for the unresolved resonance region (URR) domain. If no URR data present, returns -1.
889 *
890 * @return Maximum energy for the unresolved resonance region (URR) domain.
891 ***********************************************************************************************************/
892
894
895 double URR_domainMax1 = 0.0;
896
897 switch( protareType( ) ) {
899 URR_domainMax1 = static_cast<ProtareSingle const *>( this )->URR_domainMax( );
900 break;
902 URR_domainMax1 = static_cast<ProtareComposite const *>( this )->URR_domainMax( );
903 break;
905 URR_domainMax1 = static_cast<ProtareTNSL const *>( this )->URR_domainMax( );
906 break;
907 }
908
909 return( URR_domainMax1 );
910}
911
912/* *********************************************************************************************************//**
913 * Returns *true* if the reaction at index *a_index* has URR robability tables and *false* otherwise.
914 *
915 * @param a_index [in] The index of the reaction.
916 *
917 * @return *true* if the reaction has URR robability tables and *false* otherwise.
918 ***********************************************************************************************************/
919
921
922 bool reactionHasURR_probabilityTables1 = false;
923
924 switch( protareType( ) ) {
926 reactionHasURR_probabilityTables1 = static_cast<ProtareSingle const *>( this )->reactionHasURR_probabilityTables( a_index );
927 break;
929 reactionHasURR_probabilityTables1 = static_cast<ProtareComposite const *>( this )->reactionHasURR_probabilityTables( a_index );
930 break;
932 reactionHasURR_probabilityTables1 = static_cast<ProtareTNSL const *>( this )->reactionHasURR_probabilityTables( a_index );
933 break;
934 }
935
936 return( reactionHasURR_probabilityTables1 );
937}
938
939/* *********************************************************************************************************//**
940 * Returns the threshold for the reaction at index *a_index*. If *a_index* is negative, it is set to 0 before the
941 * threshold in the regular protare is returned.
942 *
943 * @param a_index [in] The index of the reaction.
944 *
945 * @return The threshold for reaction at index *a_index*.
946 ***********************************************************************************************************/
947
948LUPI_HOST_DEVICE double Protare::threshold( std::size_t a_index ) const {
949
950 double threshold1 = 0.0;
951
952 switch( protareType( ) ) {
954 threshold1 = static_cast<ProtareSingle const *>( this )->threshold( a_index );
955 break;
957 threshold1 = static_cast<ProtareComposite const *>( this )->threshold( a_index );
958 break;
960 threshold1 = static_cast<ProtareTNSL const *>( this )->threshold( a_index );
961 break;
962 }
963
964 return( threshold1 );
965}
966
967/* *********************************************************************************************************//**
968 * Returns the total cross section.
969 *
970 * @param a_URR_protareInfos [in] URR information.
971 * @param a_hashIndex [in] The cross section hash index.
972 * @param a_temperature [in] The target temperature.
973 * @param a_energy [in] The projectile energy.
974 * @param a_sampling [in] Only used for multi-group cross sections. When sampling, the cross section in the group where threshold
975 * is present the cross section is augmented.
976 *
977 * @return The total cross section.
978 ***********************************************************************************************************/
979
980LUPI_HOST_DEVICE double Protare::crossSection( URR_protareInfos const &a_URR_protareInfos, std::size_t a_hashIndex, double a_temperature, double a_energy, bool a_sampling ) const {
981
982 double crossSection1 = 0.0;
983
984 switch( protareType( ) ) {
986 crossSection1 = static_cast<ProtareSingle const *>( this )->crossSection( a_URR_protareInfos, a_hashIndex, a_temperature, a_energy, a_sampling );
987 break;
989 crossSection1 = static_cast<ProtareComposite const *>( this )->crossSection( a_URR_protareInfos, a_hashIndex, a_temperature, a_energy, a_sampling );
990 break;
992 crossSection1 = static_cast<ProtareTNSL const *>( this )->crossSection( a_URR_protareInfos, a_hashIndex, a_temperature, a_energy, a_sampling );
993 break;
994 }
995
996 return( crossSection1 );
997}
998
999/* *********************************************************************************************************//**
1000 * Adds the energy dependent, total cross section corresponding to the temperature *a_temperature* multiplied by *a_userFactor* to *a_crossSectionVector*.
1001 *
1002 * @param a_temperature [in] Specifies the temperature of the material.
1003 * @param a_userFactor [in] User factor which all cross sections are multiplied by.
1004 * @param a_numberAllocated [in] The length of memory allocated for *a_crossSectionVector*.
1005 * @param a_crossSectionVector [in/out] The energy dependent, total cross section to add cross section data to.
1006 ***********************************************************************************************************/
1007
1008LUPI_HOST_DEVICE void Protare::crossSectionVector( double a_temperature, double a_userFactor, std::size_t a_numberAllocated,
1009 double *a_crossSectionVector ) const {
1010
1011 switch( protareType( ) ) {
1012 case ProtareType::single:
1013 static_cast<ProtareSingle const *>( this )->crossSectionVector( a_temperature, a_userFactor, a_numberAllocated, a_crossSectionVector );
1014 break;
1016 static_cast<ProtareComposite const *>( this )->crossSectionVector( a_temperature, a_userFactor, a_numberAllocated, a_crossSectionVector );
1017 break;
1018 case ProtareType::TNSL:
1019 static_cast<ProtareTNSL const *>( this )->crossSectionVector( a_temperature, a_userFactor, a_numberAllocated, a_crossSectionVector );
1020 break;
1021 }
1022}
1023
1024/* *********************************************************************************************************//**
1025 * Returns the cross section for reaction at index *a_reactionIndex*.
1026 *
1027 * @param a_reactionIndex [in] The index of the reaction.
1028 * @param a_URR_protareInfos [in] URR information.
1029 * @param a_hashIndex [in] The cross section hash index.
1030 * @param a_temperature [in] The target temperature.
1031 * @param a_energy [in] The projectile energy.
1032 * @param a_sampling [in] Only used for multi-group cross sections. When sampling, the cross section in the group where threshold
1033 * is present the cross section is augmented.
1034 *
1035 * @return The total cross section.
1036 ***********************************************************************************************************/
1037
1038LUPI_HOST_DEVICE double Protare::reactionCrossSection( std::size_t a_reactionIndex, URR_protareInfos const &a_URR_protareInfos, std::size_t a_hashIndex,
1039 double a_temperature, double a_energy, bool a_sampling ) const {
1040
1041 double reactionCrossSection1 = 0.0;
1042
1043 switch( protareType( ) ) {
1044 case ProtareType::single:
1045 reactionCrossSection1 = static_cast<ProtareSingle const *>( this )->reactionCrossSection( a_reactionIndex, a_URR_protareInfos, a_hashIndex,
1046 a_temperature, a_energy, a_sampling );
1047 break;
1049 reactionCrossSection1 = static_cast<ProtareComposite const *>( this )->reactionCrossSection( a_reactionIndex, a_URR_protareInfos, a_hashIndex,
1050 a_temperature, a_energy, a_sampling );
1051 break;
1052 case ProtareType::TNSL:
1053 reactionCrossSection1 = static_cast<ProtareTNSL const *>( this )->reactionCrossSection( a_reactionIndex, a_URR_protareInfos, a_hashIndex,
1054 a_temperature, a_energy, a_sampling );
1055 break;
1056 }
1057
1058 return( reactionCrossSection1 );
1059}
1060
1061/* *********************************************************************************************************//**
1062 * Returns the cross section for reaction at index *a_reactionIndex*.
1063 *
1064 * @param a_reactionIndex [in] The index of the reaction.
1065 * @param a_URR_protareInfos [in] URR information.
1066 * @param a_temperature [in] The target temperature.
1067 * @param a_energy [in] The projectile energy.
1068 *
1069 * @return The total cross section.
1070 ***********************************************************************************************************/
1071
1072LUPI_HOST_DEVICE double Protare::reactionCrossSection( std::size_t a_reactionIndex, URR_protareInfos const &a_URR_protareInfos, double a_temperature, double a_energy ) const {
1073
1074 double reactionCrossSection1 = 0.0;
1075
1076 switch( protareType( ) ) {
1077 case ProtareType::single:
1078 reactionCrossSection1 = static_cast<ProtareSingle const *>( this )->reactionCrossSection( a_reactionIndex, a_URR_protareInfos, a_temperature, a_energy );
1079 break;
1081 reactionCrossSection1 = static_cast<ProtareComposite const *>( this )->reactionCrossSection( a_reactionIndex, a_URR_protareInfos, a_temperature, a_energy );
1082 break;
1083 case ProtareType::TNSL:
1084 reactionCrossSection1 = static_cast<ProtareTNSL const *>( this )->reactionCrossSection( a_reactionIndex, a_URR_protareInfos, a_temperature, a_energy );
1085 break;
1086 }
1087
1088 return( reactionCrossSection1 );
1089}
1090
1091/* *********************************************************************************************************//**
1092 * Returns the total deposition energy.
1093 *
1094 * @param a_hashIndex [in] The cross section hash index.
1095 * @param a_temperature [in] The target temperature.
1096 * @param a_energy [in] The projectile energy.
1097 *
1098 * @return The total deposition energy.
1099 ***********************************************************************************************************/
1100
1101LUPI_HOST_DEVICE double Protare::depositionEnergy( std::size_t a_hashIndex, double a_temperature, double a_energy ) const {
1102
1103 double depositionEnergy1 = 0.0;
1104
1105 switch( protareType( ) ) {
1106 case ProtareType::single:
1107 depositionEnergy1 = static_cast<ProtareSingle const *>( this )->depositionEnergy( a_hashIndex, a_temperature, a_energy );
1108 break;
1110 depositionEnergy1 = static_cast<ProtareComposite const *>( this )->depositionEnergy( a_hashIndex, a_temperature, a_energy );
1111 break;
1112 case ProtareType::TNSL:
1113 depositionEnergy1 = static_cast<ProtareTNSL const *>( this )->depositionEnergy( a_hashIndex, a_temperature, a_energy );
1114 break;
1115 }
1116
1117 return( depositionEnergy1 );
1118}
1119
1120/* *********************************************************************************************************//**
1121 * Returns the total deposition momentum.
1122 *
1123 * @param a_hashIndex [in] The cross section hash index.
1124 * @param a_temperature [in] The target temperature.
1125 * @param a_energy [in] The projectile energy.
1126 *
1127 * @return The total deposition momentum.
1128 ***********************************************************************************************************/
1129
1130LUPI_HOST_DEVICE double Protare::depositionMomentum( std::size_t a_hashIndex, double a_temperature, double a_energy ) const {
1131
1132 double depositionMomentum1 = 0.0;
1133
1134 switch( protareType( ) ) {
1135 case ProtareType::single:
1136 depositionMomentum1 = static_cast<ProtareSingle const *>( this )->depositionMomentum( a_hashIndex, a_temperature, a_energy );
1137 break;
1139 depositionMomentum1 = static_cast<ProtareComposite const *>( this )->depositionMomentum( a_hashIndex, a_temperature, a_energy );
1140 break;
1141 case ProtareType::TNSL:
1142 depositionMomentum1 = static_cast<ProtareTNSL const *>( this )->depositionMomentum( a_hashIndex, a_temperature, a_energy );
1143 break;
1144 }
1145
1146 return( depositionMomentum1 );
1147}
1148
1149/* *********************************************************************************************************//**
1150 * Returns the total production energy.
1151 *
1152 * @param a_hashIndex [in] The cross section hash index.
1153 * @param a_temperature [in] The target temperature.
1154 * @param a_energy [in] The projectile energy.
1155 *
1156 * @return The total production energy.
1157 ***********************************************************************************************************/
1158
1159LUPI_HOST_DEVICE double Protare::productionEnergy( std::size_t a_hashIndex, double a_temperature, double a_energy ) const {
1160
1161 double productionEnergy1 = 0.0;
1162
1163 switch( protareType( ) ) {
1164 case ProtareType::single:
1165 productionEnergy1 = static_cast<ProtareSingle const *>( this )->productionEnergy( a_hashIndex, a_temperature, a_energy );
1166 break;
1168 productionEnergy1 = static_cast<ProtareComposite const *>( this )->productionEnergy( a_hashIndex, a_temperature, a_energy );
1169 break;
1170 case ProtareType::TNSL:
1171 productionEnergy1 = static_cast<ProtareTNSL const *>( this )->productionEnergy( a_hashIndex, a_temperature, a_energy );
1172 break;
1173 }
1174
1175 return( productionEnergy1 );
1176}
1177
1178/* *********************************************************************************************************//**
1179 * Returns the gain for particle with index *a_particleIndex*.
1180 *
1181 * @param a_hashIndex [in] The continuous energy hash or multi-group index.
1182 * @param a_temperature [in] The temperature of the target.
1183 * @param a_energy [in] The projectile energy.
1184 * @param a_particleIndex [in] The index of the particle whose gain is to be returned.
1185 *
1186 * @return A vector of the length of the number of multi-group groups.
1187 ***********************************************************************************************************/
1188
1189LUPI_HOST_DEVICE double Protare::gain( std::size_t a_hashIndex, double a_temperature, double a_energy, int a_particleIndex ) const {
1190
1191 double gain1 = 0.0;
1192
1193 switch( protareType( ) ) {
1194 case ProtareType::single:
1195 gain1 = static_cast<ProtareSingle const *>( this )->gain( a_hashIndex, a_temperature, a_energy, a_particleIndex );
1196 break;
1198 gain1 = static_cast<ProtareComposite const *>( this )->gain( a_hashIndex, a_temperature, a_energy, a_particleIndex );
1199 break;
1200 case ProtareType::TNSL:
1201 gain1 = static_cast<ProtareTNSL const *>( this )->gain( a_hashIndex, a_temperature, a_energy, a_particleIndex );
1202 break;
1203 }
1204
1205 return( gain1 );
1206}
1207
1208/* *********************************************************************************************************//**
1209 * Returns the gain for particle with intid *a_particleIntid*.
1210 *
1211 * @param a_hashIndex [in] The continuous energy hash or multi-group index.
1212 * @param a_temperature [in] The temperature of the target.
1213 * @param a_energy [in] The projectile energy.
1214 * @param a_particleIntid [in] The intid of the particle whose gain is to be returned.
1215 *
1216 * @return A vector of the length of the number of multi-group groups.
1217 ***********************************************************************************************************/
1218
1219LUPI_HOST_DEVICE double Protare::gainViaIntid( std::size_t a_hashIndex, double a_temperature, double a_energy, int a_particleIntid ) const {
1220
1221 double gain1 = 0.0;
1222
1223 switch( protareType( ) ) {
1225 gain1 = static_cast<ProtareSingle const *>( this )->gainViaIntid( a_hashIndex, a_temperature, a_energy, a_particleIntid );
1226 break;
1228 gain1 = static_cast<ProtareComposite const *>( this )->gainViaIntid( a_hashIndex, a_temperature, a_energy, a_particleIntid );
1229 break;
1230 case ProtareType::TNSL:
1231 gain1 = static_cast<ProtareTNSL const *>( this )->gainViaIntid( a_hashIndex, a_temperature, a_energy, a_particleIntid );
1232 break;
1233 }
1234
1235 return( gain1 );
1236}
1237
1238/* *********************************************************************************************************//**
1239 * Returns a reference to the **m_upscatterModelAGroupVelocities**.
1240 *
1241 * @return Reference to the upscatter model A group velocities.
1242 ***********************************************************************************************************/
1243
1245
1246 Vector<double> const *upscatterModelAGroupVelocities1 = nullptr;
1247 switch( protareType( ) ) {
1248 case ProtareType::single:
1249 upscatterModelAGroupVelocities1 = &static_cast<ProtareSingle const *>( this )->upscatterModelAGroupVelocities( );
1250 break;
1252 upscatterModelAGroupVelocities1 = &static_cast<ProtareComposite const *>( this )->upscatterModelAGroupVelocities( );
1253 break;
1254 case ProtareType::TNSL:
1255 upscatterModelAGroupVelocities1 = &static_cast<ProtareTNSL const *>( this )->upscatterModelAGroupVelocities( );
1256 break;
1257 }
1258
1259 return( *upscatterModelAGroupVelocities1 );
1260}
1261
1262/*! \class ProtareSingle
1263 * Class representing a **GNDS** <**reactionSuite**> node with only data needed for Monte Carlo transport. The
1264 * data are also stored in a way that is better suited for Monte Carlo transport. For example, cross section data
1265 * for each reaction are not stored with its reaction, but within the HeatedCrossSections member of the Protare.
1266 */
1267
1268/* *********************************************************************************************************//**
1269 * Default constructor used when broadcasting a Protare as needed by MPI or GPUs.
1270 ***********************************************************************************************************/
1271
1274 m_URR_index( -1 ),
1275 m_hasURR_probabilityTables( false ),
1276 m_URR_domainMin( -1.0 ),
1277 m_URR_domainMax( -1.0 ),
1278 m_upscatterModelASupported( false ),
1279 m_projectileMultiGroupBoundaries( 0 ),
1280 m_projectileMultiGroupBoundariesCollapsed( 0 ),
1281 m_reactions( 0 ),
1282 m_orphanProducts( 0 ) {
1283
1284}
1285
1286/* *********************************************************************************************************//**
1287 * @param a_smr [Out] If errors are not to be thrown, then the error is reported via this instance.
1288 * @param a_protare [in] The GIDI::Protare whose data is to be used to construct *this*.
1289 * @param a_pops [in] A PoPs Database instance used to get particle intids and possibly other particle information.
1290 * @param a_settings [in] Used to pass user options to the *this* to instruct it which data are desired.
1291 * @param a_particles [in] List of transporting particles and their information (e.g., multi-group boundaries and fluxes).
1292 * @param a_domainHash [in] The hash data used when looking up a cross section.
1293 * @param a_temperatureInfos [in] The list of temperature data to extract from *a_protare*.
1294 * @param a_reactionsToExclude [in] A list of reaction to not include in the MCGIDI::Protare.
1295 * @param a_reactionsToExcludeOffset [in] The starting index for the reactions in this ProtareSingle.
1296 * @param a_allowFixedGrid [in] For internal (i.e., MCGIDI) use only. Users must use the default value.
1297 ***********************************************************************************************************/
1298
1300 Transporting::MC &a_settings, GIDI::Transporting::Particles const &a_particles, DomainHash const &a_domainHash,
1301 GIDI::Styles::TemperatureInfos const &a_temperatureInfos, GIDI::ExcludeReactionsSet const &a_reactionsToExclude,
1302 std::size_t a_reactionsToExcludeOffset, bool a_allowFixedGrid ) :
1303 Protare( ProtareType::single, a_protare, a_settings, a_pops ),
1304 m_interaction( a_protare.interaction( ).c_str( ) ),
1305 m_URR_index( -1 ),
1306 m_hasURR_probabilityTables( false ),
1307 m_URR_domainMin( -1.0 ),
1308 m_URR_domainMax( -1.0 ),
1309 m_domainHash( a_domainHash ),
1310 m_upscatterModelASupported( ( projectileIntid( ) != PoPI::Intids::photon ) &&
1311 ( projectileIntid( ) != PoPI::Intids::electron ) &&
1313 m_projectileMultiGroupBoundaries( 0 ),
1314 m_projectileMultiGroupBoundariesCollapsed( 0 ),
1315 m_reactions( 0 ),
1316 m_orphanProducts( 0 ),
1317 m_isPhotoAtomic( a_protare.isPhotoAtomic( ) ),
1318 m_heatedCrossSections( ),
1319 m_heatedMultigroupCrossSections( ) {
1320
1321 a_protare.updateReactionIndices( 0 ); // This is not correct as the offset should be passed as an arguent.
1322 PoPI::Database const &pops = a_protare.protare( 0 )->internalPoPs( );
1323
1324 if( !a_protare.isPhotoAtomic( ) ) {
1325 std::set<std::string> incompleteParticles;
1326 a_protare.incompleteParticles( a_settings, incompleteParticles );
1327 for( auto particle = a_particles.particles( ).begin( ); particle != a_particles.particles( ).end( ); ++particle ) {
1328 if( incompleteParticles.count( particle->first ) != 0 ) {
1329 std::string message = "Requested particle '" + particle->first + "' is incomplete in '" + a_protare.realFileName( ) + "'.";
1330 if( a_settings.throwOnError( ) ) {
1331 throw std::runtime_error( message.c_str( ) ); }
1332 else {
1333 smr_setReportError2p( a_smr.smr( ), 0, 0, message.c_str( ) );
1334 }
1335 }
1336 }
1337 }
1338
1339 SetupInfo setupInfo( *this, a_protare, a_pops, pops );
1340 setupInfo.m_formatVersion = a_protare.formatVersion( );
1341 setupInfo.m_GRIN_continuumGammas = a_protare.GRIN_continuumGammas2( );
1342
1344 for( std::map<std::string, GIDI::Transporting::Particle>::const_iterator particle = a_particles.particles( ).begin( ); particle != a_particles.particles( ).end( ); ++particle ) {
1345 setupInfo.m_particleIntids[particle->first] = MCGIDI_popsIntid( pops, particle->first );
1346 setupInfo.m_particleIndices[particle->first] = MCGIDI_popsIndex( a_pops, particle->first );
1347
1348 if( ( m_interaction == GIDI_MapInteractionAtomicChars ) &&
1349 !( ( particle->first == PoPI::IDs::photon ) || ( particle->first == PoPI::IDs::electron ) ) ) continue;
1350 particles.add( particle->second );
1351 }
1352
1353 GIDI::Transporting::MG multiGroupSettings( a_settings.projectileID( ), GIDI::Transporting::Mode::MonteCarloContinuousEnergy, a_settings.delayedNeutrons( ) );
1354 multiGroupSettings.setThrowOnError( a_settings.throwOnError( ) );
1355
1356 setupInfo.m_distributionLabel = a_temperatureInfos[0].griddedCrossSection( );
1357
1358 a_settings.styles( &a_protare.styles( ) );
1359
1360 switch( a_settings.crossSectionLookupMode( ) ) {
1362 m_continuousEnergy = true;
1363 break;
1365 m_continuousEnergy = false;
1366 if( a_settings.upscatterModel( ) == Sampling::Upscatter::Model::B ) {
1367 multiGroupSettings.setMode( GIDI::Transporting::Mode::multiGroupWithSnElasticUpScatter ); }
1368 else {
1369 multiGroupSettings.setMode( GIDI::Transporting::Mode::multiGroup );
1370 }
1371 break;
1372 default :
1373 throw std::runtime_error( "ProtareSingle::ProtareSingle: invalid lookupMode" );
1374 }
1375 m_fixedGrid = a_allowFixedGrid && ( a_protare.projectile( ).ID( ) == PoPI::IDs::photon ) && ( a_settings.fixedGridPoints( ).size( ) > 0 );
1376
1377 setupNuclideGammaBranchStateInfos( setupInfo, a_protare, a_settings.makePhotonEmissionProbabilitiesOne( ),
1378 a_settings.zeroNuclearLevelEnergyWidth( ) );
1379 convertACE_URR_probabilityTablesFromGIDI( a_protare, a_settings, setupInfo );
1380
1381 if( ( a_settings.crossSectionLookupMode( ) == Transporting::LookupMode::Data1d::multiGroup ) ||
1382 ( a_settings.other1dDataLookupMode( ) == Transporting::LookupMode::Data1d::multiGroup ) ) {
1383
1384 GIDI::Suite const *transportables = nullptr;
1385 if( setupInfo.m_formatVersion.major( ) > 1 ) {
1386 GIDI::Styles::HeatedMultiGroup const &heatedMultiGroup = *a_protare.styles( ).get<GIDI::Styles::HeatedMultiGroup>( a_settings.label( ) );
1387 transportables = &heatedMultiGroup.transportables( ); }
1388 else {
1389 std::vector<GIDI::Suite::const_iterator> tags = a_protare.styles( ).findAllOfMoniker( GIDI_multiGroupStyleChars );
1390
1391 if( tags.size( ) != 1 ) throw std::runtime_error( "MCGIDI::ProtareSingle::ProtareSingle: What is going on here?" );
1392 GIDI::Styles::MultiGroup const &multiGroup = static_cast<GIDI::Styles::MultiGroup const &>( **tags[0] );
1393 transportables = &multiGroup.transportables( );
1394 }
1395
1396 GIDI::Transportable const transportable = *transportables->get<GIDI::Transportable>( a_protare.projectile( ).ID( ) );
1397 m_projectileMultiGroupBoundaries = transportable.groupBoundaries( );
1398 GIDI::Transporting::Particle const *particle = a_particles.particle( a_protare.projectile( ).ID( ) );
1399 m_projectileMultiGroupBoundariesCollapsed = particle->multiGroup( ).boundaries( );
1400 }
1401
1402 std::vector<GIDI::Reaction const *> GIDI_reactions;
1403 std::set<std::string> product_ids;
1404 std::set<std::string> product_ids_transportable;
1405 GIDI::Reaction const *nuclearPlusCoulombInterferenceReaction = nullptr;
1406 if( a_settings.nuclearPlusCoulombInterferenceOnly( ) ) nuclearPlusCoulombInterferenceReaction = a_protare.nuclearPlusCoulombInterferenceOnlyReaction( );
1407
1408 for( std::size_t reactionIndex = 0; reactionIndex < a_protare.reactions( ).size( ); ++reactionIndex ) {
1409 if( a_reactionsToExclude.find( reactionIndex + a_reactionsToExcludeOffset ) != a_reactionsToExclude.end( ) ) continue;
1410
1411 GIDI::Reaction const *GIDI_reaction = a_protare.reaction( reactionIndex );
1412
1413 if( !GIDI_reaction->active( ) ) continue;
1414
1415 if( m_continuousEnergy ) {
1416 if( GIDI_reaction->crossSectionThreshold( ) >= a_settings.energyDomainMax( ) ) continue; }
1417 else {
1418 GIDI::Vector multi_group_cross_section = GIDI_reaction->multiGroupCrossSection( a_smr, multiGroupSettings, a_temperatureInfos[0] );
1419 GIDI::Vector vector = GIDI::collapse( multi_group_cross_section, a_settings, a_particles, 0.0 );
1420
1421 std::size_t i1 = 0;
1422 for( ; i1 < vector.size( ); ++i1 ) if( vector[i1] != 0.0 ) break;
1423 if( i1 == vector.size( ) ) continue;
1424 }
1425 if( a_settings.ignoreENDF_MT5( ) && ( GIDI_reaction->ENDF_MT( ) == 5 ) && ( a_reactionsToExclude.size( ) == 0 ) ) continue;
1426
1427 GIDI_reaction->productIDs( product_ids, particles, false );
1428 GIDI_reaction->productIDs( product_ids_transportable, particles, true );
1429
1430 if( a_settings.nuclearPlusCoulombInterferenceOnly( ) && GIDI_reaction->RutherfordScatteringPresent( ) ) {
1431 if( nuclearPlusCoulombInterferenceReaction != nullptr ) GIDI_reactions.push_back( nuclearPlusCoulombInterferenceReaction ); }
1432 else {
1433 GIDI_reactions.push_back( GIDI_reaction );
1434 }
1435 }
1436
1437 bool zeroReactions = GIDI_reactions.size( ) == 0; // Happens when all reactions are skipped in the prior loop.
1438 if( zeroReactions ) GIDI_reactions.push_back( a_protare.reaction( 0 ) ); // Special case where no reaction in the protare is wanted so the first one is used but its cross section is set to 0.0 at all energies.
1439
1440 setupInfo.m_reactionType = Transporting::Reaction::Type::Reactions;
1441 m_reactions.reserve( GIDI_reactions.size( ) );
1442 for( auto GIDI_reaction = GIDI_reactions.begin( ); GIDI_reaction != GIDI_reactions.end( ); ++GIDI_reaction ) {
1443 setupInfo.m_reaction = *GIDI_reaction;
1444 setupInfo.m_isPairProduction = (*GIDI_reaction)->isPairProduction( );
1445 setupInfo.m_isPhotoAtomicIncoherentScattering = (*GIDI_reaction)->isPhotoAtomicIncoherentScattering( );
1446 setupInfo.m_initialStateIndex = -1;
1447 Reaction *reaction = new Reaction( **GIDI_reaction, setupInfo, a_settings, particles, a_temperatureInfos );
1448 setupInfo.m_initialStateIndices[(*GIDI_reaction)->label( )] = setupInfo.m_initialStateIndex;
1449 reaction->updateProtareSingleInfo( this, m_reactions.size( ) );
1450 m_reactions.push_back( reaction );
1451 }
1452
1453 std::set<int> product_intids;
1454 std::set<int> product_indices;
1455 for( std::set<std::string>::iterator iter = product_ids.begin( ); iter != product_ids.end( ); ++iter ) {
1456 product_intids.insert( MCGIDI_popsIntid( pops, *iter ) );
1457 product_indices.insert( MCGIDI_popsIndex( a_pops, *iter ) );
1458 }
1459 std::set<int> product_intids_transportable;
1460 std::set<int> product_indices_transportable;
1461 for( std::set<std::string>::iterator iter = product_ids_transportable.begin( ); iter != product_ids_transportable.end( ); ++iter ) {
1462 product_intids_transportable.insert( MCGIDI_popsIntid( pops, *iter ) );
1463 product_indices_transportable.insert( MCGIDI_popsIndex( a_pops, *iter ) );
1464 }
1465 productIntidsAndIndices( product_intids, product_intids_transportable, product_indices, product_indices_transportable );
1466
1467 if( a_settings.sampleNonTransportingParticles( ) || particles.hasParticle( PoPI::IDs::photon ) ) {
1468 setupInfo.m_reactionType = Transporting::Reaction::Type::OrphanProducts;
1469 m_orphanProducts.reserve( a_protare.orphanProducts( ).size( ) );
1470 std::vector< std::vector<std::size_t> > associatedOrphanProductIndices( m_reactions.size( ) );
1471
1472 for( std::size_t orphanProductIndex = 0; orphanProductIndex < a_protare.orphanProducts( ).size( ); ++orphanProductIndex ) {
1473 GIDI::Reaction const *GIDI_reaction = a_protare.orphanProduct( orphanProductIndex );
1474
1475 if( GIDI_reaction->crossSectionThreshold( ) >= a_settings.energyDomainMax( ) ) continue;
1476
1477 setupInfo.m_reaction = GIDI_reaction;
1478 Reaction *orphanProductReaction = new Reaction( *GIDI_reaction, setupInfo, a_settings, particles, a_temperatureInfos );
1479 orphanProductReaction->updateProtareSingleInfo( this, m_orphanProducts.size( ) );
1480 m_orphanProducts.push_back( orphanProductReaction );
1481
1482 GIDI::Functions::Reference1d const *reference( GIDI_reaction->crossSection( ).get<GIDI::Functions::Reference1d>( 0 ) );
1483 std::string xlink = reference->xlink( );
1484 GUPI::Ancestry const *ancestry = a_protare.findInAncestry( xlink );
1485 if( ancestry == nullptr ) throw std::runtime_error( "Could not find xlink for orphan product - 1." );
1486 ancestry = ancestry->ancestor( );
1487 if( ancestry == nullptr ) throw std::runtime_error( "Could not find xlink for orphan product - 2." );
1488 if( ancestry->moniker( ) != GIDI_crossSectionSumChars ) {
1489 ancestry = ancestry->ancestor( );
1490 if( ancestry == nullptr ) throw std::runtime_error( "Could not find xlink for orphan product - 3." );
1491 }
1492 GIDI::Sums::CrossSectionSum const *crossSectionSum = static_cast<GIDI::Sums::CrossSectionSum const *>( ancestry );
1493 GIDI::Sums::Summands const &summands = crossSectionSum->summands( );
1494 for( std::size_t i1 = 0; i1 < summands.size( ); ++i1 ) {
1495 GIDI::Sums::Summand::Base const *summand = summands[i1];
1496
1497 ancestry = a_protare.findInAncestry( summand->href( ) );
1498 if( ancestry == nullptr ) throw std::runtime_error( "Could not find href for summand - 1." );
1499 ancestry = ancestry->ancestor( );
1500 if( ancestry == nullptr ) throw std::runtime_error( "Could not find href for summand - 2." );
1501
1502 GIDI::Reaction const *GIDI_reaction2 = static_cast<GIDI::Reaction const *>( ancestry );
1503 for( std::size_t reactionIndex = 0; reactionIndex < m_reactions.size( ); ++reactionIndex ) {
1504 std::string label( m_reactions[reactionIndex]->label( ).c_str( ) );
1505
1506 if( label == GIDI_reaction2->label( ) ) {
1507 associatedOrphanProductIndices[reactionIndex].push_back( m_orphanProducts.size( ) - 1 );
1508 break;
1509 }
1510 }
1511 }
1512 }
1513
1514 for( std::size_t reactionIndex = 0; reactionIndex < m_reactions.size( ); ++reactionIndex ) {
1515 Reaction *reaction = m_reactions[reactionIndex];
1516 std::size_t size = associatedOrphanProductIndices[reactionIndex].size( );
1517 if( size > 0 ) {
1518 std::vector<Product *> associatedOrphanProducts;
1519 for( std::size_t index1 = 0; index1 < size; ++index1 ) {
1520 std::size_t associatedOrphanProductIndex = associatedOrphanProductIndices[reactionIndex][index1];
1521 m_orphanProducts[associatedOrphanProductIndex]->addOrphanProductToProductList( associatedOrphanProducts );
1522 }
1523 reaction->setOrphanProductData( associatedOrphanProductIndices[reactionIndex], associatedOrphanProducts );
1524 }
1525 }
1526 }
1527
1528 std::vector<GIDI::Reaction const *> GIDI_orphanProducts;
1529 for( std::size_t reactionIndex = 0; reactionIndex < a_protare.orphanProducts( ).size( ); ++reactionIndex ) {
1530 GIDI::Reaction const *GIDI_reaction = a_protare.orphanProduct( reactionIndex );
1531
1532 if( GIDI_reaction->crossSectionThreshold( ) >= a_settings.energyDomainMax( ) ) continue;
1533 GIDI_orphanProducts.push_back( GIDI_reaction );
1534 }
1535
1536 bool removeContinuousEnergyData = false;
1537 if( m_continuousEnergy ) {
1538 m_heatedCrossSections.update( a_smr, setupInfo, a_settings, particles, a_domainHash, a_temperatureInfos, GIDI_reactions, GIDI_orphanProducts,
1539 m_fixedGrid, zeroReactions );
1540 m_hasURR_probabilityTables = m_heatedCrossSections.hasURR_probabilityTables( );
1541 m_URR_domainMin = m_heatedCrossSections.URR_domainMin( );
1542 m_URR_domainMax = m_heatedCrossSections.URR_domainMax( ); }
1543 else {
1544 m_heatedMultigroupCrossSections.update( a_smr, a_protare, setupInfo, a_settings, particles, a_temperatureInfos, GIDI_reactions,
1545 GIDI_orphanProducts, zeroReactions, a_reactionsToExclude );
1546
1547 if( ( a_settings.upscatterModelAGroupBoundaries().size( ) > 0 ) ||
1548 ( a_settings.upscatterModel( ) == Sampling::Upscatter::Model::DBRC ) ) { // Load pointwise data to recompute Model A cross sections on user-defined grid or for DBRC.
1549 removeContinuousEnergyData = true;
1550 m_heatedCrossSections.update( a_smr, setupInfo, a_settings, particles, a_domainHash, a_temperatureInfos, GIDI_reactions, GIDI_orphanProducts,
1551 m_fixedGrid, zeroReactions );
1552 }
1553 }
1554
1555 if( m_upscatterModelASupported && ( a_settings.upscatterModel( ) == Sampling::Upscatter::Model::A ) ) {
1556 std::vector<double> const &upscatterModelAGroupBoundaries = a_settings.upscatterModelAGroupBoundaries( );
1557 if( upscatterModelAGroupBoundaries.size( ) == 0 ) {
1558 GIDI::Styles::Base const *style = a_protare.styles( ).get<GIDI::Styles::Base>( a_temperatureInfos[0].heatedMultiGroup( ) );
1559
1560 if( style->moniker( ) == GIDI_SnElasticUpScatterStyleChars ) style = a_protare.styles( ).get<GIDI::Styles::Base>( style->derivedStyle( ) );
1561 if( style->moniker( ) != GIDI_heatedMultiGroupStyleChars ) throw GIDI::Exception( "Label does not yield a heatedMultiGroup style." );
1562
1563 GIDI::Styles::HeatedMultiGroup const &heatedMultiGroup = *static_cast<GIDI::Styles::HeatedMultiGroup const *>( style );
1564 std::vector<double> const &boundaries = heatedMultiGroup.groupBoundaries( a_protare.projectile( ).ID( ) );
1565
1566 m_upscatterModelAGroupEnergies.resize( boundaries.size( ) );
1567 m_upscatterModelAGroupVelocities.resize( boundaries.size( ) );
1568 for( std::size_t i1 = 0; i1 < boundaries.size( ); ++i1 ) {
1569 m_upscatterModelAGroupEnergies[i1] = boundaries[i1];
1570 m_upscatterModelAGroupVelocities[i1] = MCGIDI_particleBeta( projectileMass( ), boundaries[i1] );
1571 }
1572
1573 GIDI::ExcludeReactionsSet reactionsToExclude;
1574 auto upscatterModelACrossSectionForm = a_protare.multiGroupCrossSection( a_smr, multiGroupSettings, a_temperatureInfos[0],
1575 reactionsToExclude, a_temperatureInfos[0].heatedMultiGroup( ) );
1576 m_upscatterModelACrossSection.resize( upscatterModelACrossSectionForm.size( ) );
1577 for( std::size_t i1 = 0; i1 < upscatterModelACrossSectionForm.size( ); ++i1 )
1578 m_upscatterModelACrossSection[i1] = upscatterModelACrossSectionForm[i1]; }
1579 else {
1580
1581 m_upscatterModelAGroupEnergies.reserve( upscatterModelAGroupBoundaries.size( ) );
1582 m_upscatterModelAGroupVelocities.reserve( upscatterModelAGroupBoundaries.size( ) );
1583 for( auto iter = upscatterModelAGroupBoundaries.begin( ); iter != upscatterModelAGroupBoundaries.end( ); ++iter ) {
1584 m_upscatterModelAGroupEnergies.push_back( *iter );
1585 m_upscatterModelAGroupVelocities.push_back( MCGIDI_particleBeta( projectileMass( ), *iter ) );
1586 }
1587
1588 GIDI::Transporting::MultiGroup boundaries( "Model A", upscatterModelAGroupBoundaries );
1589
1590 GIDI::Functions::XYs1d crossSectionXYs1d = m_heatedCrossSections.crossSectionAsGIDI_XYs1d( 0.0 );
1591
1592 GIDI::Transporting::Flux flux( "Model A", 0.0 );
1593 std::vector<double> energies, fluxes;
1594 energies.push_back( upscatterModelAGroupBoundaries[0] );
1595 energies.push_back( upscatterModelAGroupBoundaries.back( ) );
1596 fluxes.push_back( 1.0 );
1597 fluxes.push_back( 1.0 );
1598 flux.addFluxOrder( GIDI::Transporting::Flux_order( 0, energies, fluxes ) );
1599
1600 GIDI::Vector crossSectionVector = multiGroupXYs1d( boundaries, crossSectionXYs1d, flux );
1601 m_upscatterModelACrossSection.resize( crossSectionVector.size( ) );
1602 for( std::size_t index = 0; index < crossSectionVector.size( ); ++index )
1603 m_upscatterModelACrossSection[index] = crossSectionVector[index];
1604 }
1605 if( !m_continuousEnergy ) m_multiGroupHash = MultiGroupHash( m_projectileMultiGroupBoundariesCollapsed );
1606 }
1607
1608 if( ( PoPI::Intids::neutron == projectileIntid( ) ) && ( a_settings.upscatterModel( ) == Sampling::Upscatter::Model::DBRC ) ) {
1609 std::size_t reactionIndex = 0;
1610 for( auto reactionIter = m_reactions.begin( ); reactionIter != m_reactions.end( ); ++reactionIter, ++reactionIndex ) {
1611 if( (*reactionIter)->ENDF_MT( ) == 2 ) {
1612 Reaction *reaction = *reactionIter;
1613
1614 HeatedCrossSectionContinuousEnergy const *heatedCrossSectionContinuousEnergy = m_heatedCrossSections.heatedCrossSections( )[0];
1615 HeatedReactionCrossSectionContinuousEnergy const *heatedReactionCrossSectionContinuousEnergy
1616 = heatedCrossSectionContinuousEnergy->reactionCrossSection( reactionIndex );
1617
1618 Vector<double> const &energies = heatedCrossSectionContinuousEnergy->energies( );
1619 Vector<MCGIDI_FLOAT> const &crossSectionsFloat = heatedReactionCrossSectionContinuousEnergy->crossSections( );
1620 Vector<double> crossSections( crossSectionsFloat.size( ) );
1621 std::size_t index = 0;
1622 for( auto iter = crossSectionsFloat.begin( ); iter != crossSectionsFloat.end( ); ++iter, ++index )
1623 crossSections[index] = *iter;
1624
1625 Sampling::Upscatter::ModelDBRC_data *modelDBRC_data =
1626 new Sampling::Upscatter::ModelDBRC_data( projectileMass( ), targetMass( ), energies, crossSections, a_domainHash );
1627 reaction->setModelDBRC_data( modelDBRC_data );
1628 break;
1629 }
1630 }
1631 }
1632 if( removeContinuousEnergyData ) {
1633 m_heatedCrossSections.clear( );
1634 }
1635}
1636
1637/* *********************************************************************************************************//**
1638 ***********************************************************************************************************/
1639
1641
1642 for( Vector<NuclideGammaBranchInfo *>::const_iterator iter = m_branches.begin( ); iter < m_branches.end( ); ++iter ) delete *iter;
1643 for( Vector<NuclideGammaBranchStateInfo *>::const_iterator iter = m_nuclideGammaBranchStateInfos.begin( ); iter < m_nuclideGammaBranchStateInfos.end( ); ++iter ) delete *iter;
1644 for( Vector<Reaction *>::const_iterator iter = m_reactions.begin( ); iter < m_reactions.end( ); ++iter ) delete *iter;
1645 for( Vector<Reaction *>::const_iterator iter = m_orphanProducts.begin( ); iter < m_orphanProducts.end( ); ++iter ) delete *iter;
1646}
1647
1648/* *********************************************************************************************************//**
1649 * Updates the m_userParticleIndex to *a_userParticleIndex* for all particles with PoPs index *a_particleIndex*.
1650 *
1651 * @param a_particleIndex [in] The PoPs index of the particle whose user index is to be set.
1652 * @param a_userParticleIndex [in] The particle index specified by the user.
1653 ***********************************************************************************************************/
1654
1655LUPI_HOST void ProtareSingle::setUserParticleIndex2( int a_particleIndex, int a_userParticleIndex ) {
1656
1657 m_heatedCrossSections.setUserParticleIndex( a_particleIndex, a_userParticleIndex );
1658 m_heatedMultigroupCrossSections.setUserParticleIndex( a_particleIndex, a_userParticleIndex );
1659 for( auto iter = m_reactions.begin( ); iter < m_reactions.end( ); ++iter ) (*iter)->setUserParticleIndex( a_particleIndex, a_userParticleIndex );
1660 for( auto iter = m_orphanProducts.begin( ); iter < m_orphanProducts.end( ); ++iter ) (*iter)->setUserParticleIndex( a_particleIndex, a_userParticleIndex );
1661}
1662
1663/* *********************************************************************************************************//**
1664 * Updates the m_userParticleIndex to *a_userParticleIndex* for all particles with PoPs intid *a_particleIntid*.
1665 *
1666 * @param a_particleIndex [in] The PoPs index of the particle whose user index is to be set.
1667 * @param a_userParticleIndex [in] The particle index specified by the user.
1668 ***********************************************************************************************************/
1669
1670LUPI_HOST void ProtareSingle::setUserParticleIndexViaIntid2( int a_particleIntid, int a_userParticleIndex ) {
1671
1672 m_heatedCrossSections.setUserParticleIndexViaIntid( a_particleIntid, a_userParticleIndex );
1673 m_heatedMultigroupCrossSections.setUserParticleIndexViaIntid( a_particleIntid, a_userParticleIndex );
1674 for( auto iter = m_reactions.begin( ); iter < m_reactions.end( ); ++iter ) (*iter)->setUserParticleIndexViaIntid( a_particleIntid, a_userParticleIndex );
1675 for( auto iter = m_orphanProducts.begin( ); iter < m_orphanProducts.end( ); ++iter ) (*iter)->setUserParticleIndexViaIntid( a_particleIntid, a_userParticleIndex );
1676}
1677
1678/* *********************************************************************************************************//**
1679 * Returns the pointer representing the protare (i.e., *this*) if *a_index* is 0 and nullptr otherwise.
1680 *
1681 * @param a_index [in] Must always be 0.
1682 *
1683 * @return Returns the pointer representing *this*.
1684 ***********************************************************************************************************/
1685
1686LUPI_HOST_DEVICE ProtareSingle const *ProtareSingle::protare( std::size_t a_index ) const {
1687
1688 if( a_index != 0 ) return( nullptr );
1689 return( this );
1690}
1691
1692/* *********************************************************************************************************//**
1693 * Returns the pointer representing the protare (i.e., *this*) if *a_index* is 0 and nullptr otherwise.
1694 *
1695 * @param a_index [in] Must always be 0.
1696 *
1697 * @return Returns the pointer representing *this*.
1698 ***********************************************************************************************************/
1699
1701
1702 if( a_index != 0 ) return( nullptr );
1703 return( this );
1704}
1705
1706/* *********************************************************************************************************//**
1707 * Returns the pointer to the **this** if (*a_index* - 1)th is a value reaction index and nullptr otherwise.
1708 *
1709 * @param a_index [in] Index of the reaction.
1710 *
1711 * @return Pointer to the requested protare or nullptr if invalid *a_index*..
1712 ***********************************************************************************************************/
1713
1715
1716 if( static_cast<std::size_t>( a_index ) < numberOfReactions( ) ) return( this );
1717 return( nullptr );
1718}
1719
1720/* *********************************************************************************************************//**
1721 * Returns the list of temperatures for *this*.
1722 *
1723 * @param a_index [in] Index of the reqested ProtareSingle. Must be 0.
1724 *
1725 * @return Vector of doubles.
1726 ***********************************************************************************************************/
1727
1729
1730 if( a_index != 0 ) LUPI_THROW( "ProtareSingle::temperatures: a_index not 0." );
1731 if( m_continuousEnergy ) return( m_heatedCrossSections.temperatures( ) );
1732 return( m_heatedMultigroupCrossSections.temperatures( ) );
1733}
1734
1735/* *********************************************************************************************************//**
1736 * Sets up the nuclear gamma branching data needed to sample gamma decays.
1737 *
1738 * @param a_setupInfo [in] Used internally when constructing a Protare to pass information to other constructors.
1739 * @param a_protare [in] The **GIDI::Protare** whose data are to be used to construct gamma branching data.
1740 * @param a_makePhotonEmissionProbabilitiesOne [in] If true, all photon emission probabilities are set to 1.0 (i.e., all ICCs are set to 0.0).
1741 ***********************************************************************************************************/
1742
1743LUPI_HOST void ProtareSingle::setupNuclideGammaBranchStateInfos( SetupInfo &a_setupInfo, GIDI::ProtareSingle const &a_protare,
1744 bool a_makePhotonEmissionProbabilitiesOne, bool a_zeroNuclearLevelEnergyWidth ) {
1745
1746 PoPI::NuclideGammaBranchStateInfos const &nuclideGammaBranchStateInfos = a_protare.nuclideGammaBranchStateInfos( );
1747 std::vector<NuclideGammaBranchInfo *> nuclideGammaBranchInfos;
1748
1749 for( std::size_t i1 = 0; i1 < nuclideGammaBranchStateInfos.size( ); ++i1 ) {
1750 a_setupInfo.m_stateNamesToIndices[nuclideGammaBranchStateInfos[i1]->state( )] = (int) i1;
1751 PoPI::NuclideGammaBranchStateInfo const *nuclideGammaBranchStateInfo = nuclideGammaBranchStateInfos.find( nuclideGammaBranchStateInfos[i1]->state( ) );
1752 a_setupInfo.m_nuclearLevelEnergies[nuclideGammaBranchStateInfos[i1]->state( )] = nuclideGammaBranchStateInfo->nuclearLevelEnergy( );
1753 }
1754
1755 m_nuclideGammaBranchStateInfos.reserve( nuclideGammaBranchStateInfos.size( ) );
1756 for( std::size_t i1 = 0; i1 < nuclideGammaBranchStateInfos.size( ); ++i1 ) {
1757 m_nuclideGammaBranchStateInfos.push_back( new NuclideGammaBranchStateInfo( *nuclideGammaBranchStateInfos[i1], nuclideGammaBranchInfos,
1758 a_setupInfo.m_stateNamesToIndices, a_makePhotonEmissionProbabilitiesOne, a_zeroNuclearLevelEnergyWidth ) );
1759 }
1760
1761 m_branches.reserve( nuclideGammaBranchInfos.size( ) );
1762 for( std::size_t i1 = 0; i1 < nuclideGammaBranchInfos.size( ); ++i1 ) m_branches.push_back( nuclideGammaBranchInfos[i1] );
1763}
1764
1765/* *********************************************************************************************************//**
1766 * Returns true if *this* has a fission reaction and false otherwise.
1767 *
1768 * @return true is if *this* has a fission reaction and false otherwise.
1769 ***********************************************************************************************************/
1770
1772
1773 for( Vector<Reaction *>::const_iterator iter = m_reactions.begin( ); iter < m_reactions.end( ); ++iter ) {
1774 if( (*iter)->hasFission( ) ) return( true );
1775 }
1776 return( false );
1777}
1778
1779/* *********************************************************************************************************//**
1780 * Returns true if *this* has an incoherent photoatomic doppler broadened reaction and false otherwise.
1781 *
1782 * @return true is if *this* has a specified reaction and false otherwise.
1783 ***********************************************************************************************************/
1784
1786
1787 for( Vector<Reaction *>::const_iterator iter = m_reactions.begin( ); iter < m_reactions.end( ); ++iter ) {
1788 if( (*iter)->ENDF_MT( ) == 1534 ) return( true );
1789 }
1790 return( false );
1791}
1792
1793/* *********************************************************************************************************//**
1794 * Returns true if *a_energy* with unresolved resonance region (URR) of *this* and false otherwise.
1795 *
1796 * @return true is if *this* has a URR data.
1797 ***********************************************************************************************************/
1798
1799LUPI_HOST_DEVICE bool ProtareSingle::inURR( double a_energy ) const {
1800
1801 if( a_energy < m_URR_domainMin ) return( false );
1802 if( a_energy > m_URR_domainMax ) return( false );
1803
1804 return( true );
1805}
1806
1807/* *********************************************************************************************************//**
1808 * Returns the total cross section for target temperature *a_temperature* and projectile energy *a_energy*.
1809 * *a_sampling* is only used for multi-group cross section look up.
1810 *
1811 * @param a_URR_protareInfos [in] URR information.
1812 * @param a_hashIndex [in] Specifies the continuous energy or multi-group index.
1813 * @param a_temperature [in] The temperature of the target.
1814 * @param a_energy [in] The energy of the projectile.
1815 * @param a_sampling [in] Used for multi-group look up. If *true*, use augmented cross sections.
1816 ***********************************************************************************************************/
1817
1818LUPI_HOST_DEVICE double ProtareSingle::crossSection( URR_protareInfos const &a_URR_protareInfos, std::size_t a_hashIndex, double a_temperature, double a_energy, bool a_sampling ) const {
1819
1820 if( m_continuousEnergy ) return( m_heatedCrossSections.crossSection( a_URR_protareInfos, m_URR_index, a_hashIndex, a_temperature, a_energy ) );
1821
1822 return( m_heatedMultigroupCrossSections.crossSection( a_hashIndex, a_temperature, a_sampling ) );
1823}
1824
1825/* *********************************************************************************************************//**
1826 * Adds the energy dependent, total cross section corresponding to the temperature *a_temperature* multiplied by *a_userFactor* to *a_crossSectionVector*.
1827 *
1828 * @param a_temperature [in] Specifies the temperature of the material.
1829 * @param a_userFactor [in] User factor which all cross sections are multiplied by.
1830 * @param a_numberAllocated [in] The length of memory allocated for *a_crossSectionVector*.
1831 * @param a_crossSectionVector [in/out] The energy dependent, total cross section to add cross section data to.
1832 ***********************************************************************************************************/
1833
1834LUPI_HOST_DEVICE void ProtareSingle::crossSectionVector( double a_temperature, double a_userFactor, std::size_t a_numberAllocated,
1835 double *a_crossSectionVector ) const {
1836
1837 if( m_continuousEnergy ) {
1838 if( !m_fixedGrid ) LUPI_THROW( "ProtareSingle::crossSectionVector: continuous energy cannot be supported." );
1839 m_heatedCrossSections.crossSectionVector( a_temperature, a_userFactor, a_numberAllocated, a_crossSectionVector ); }
1840 else {
1841 m_heatedMultigroupCrossSections.crossSectionVector( a_temperature, a_userFactor, a_numberAllocated, a_crossSectionVector );
1842 }
1843}
1844
1845/* *********************************************************************************************************//**
1846 * Returns the reaction's cross section for the reaction at index *a_reactionIndex*, for target temperature *a_temperature* and projectile energy *a_energy*.
1847 * *a_sampling* is only used for multi-group cross section look up.
1848 *
1849 * @param a_reactionIndex [in] The index of the reaction.
1850 * @param a_URR_protareInfos [in] URR information.
1851 * @param a_hashIndex [in] Specifies the continuous energy or multi-group index.
1852 * @param a_temperature [in] The temperature of the target.
1853 * @param a_energy [in] The energy of the projectile.
1854 * @param a_sampling [in] Used for multi-group look up. If *true*, use augmented cross sections.
1855 ***********************************************************************************************************/
1856
1857LUPI_HOST_DEVICE double ProtareSingle::reactionCrossSection( std::size_t a_reactionIndex, URR_protareInfos const &a_URR_protareInfos, std::size_t a_hashIndex,
1858 double a_temperature, double a_energy, bool a_sampling ) const {
1859
1860 if( m_continuousEnergy ) return( m_heatedCrossSections.reactionCrossSection( a_reactionIndex, a_URR_protareInfos, m_URR_index, a_hashIndex,
1861 a_temperature, a_energy ) );
1862
1863 return( m_heatedMultigroupCrossSections.reactionCrossSection( a_reactionIndex, a_hashIndex, a_temperature, a_sampling ) );
1864}
1865
1866/* *********************************************************************************************************//**
1867 * Returns the reaction's cross section for the reaction at index *a_reactionIndex*, for target temperature *a_temperature* and projectile energy *a_energy*.
1868 *
1869 * @param a_reactionIndex [in] The index of the reaction.
1870 * @param a_URR_protareInfos [in] URR information.
1871 * @param a_temperature [in] The temperature of the target.
1872 * @param a_energy [in] The energy of the projectile.
1873 ***********************************************************************************************************/
1874
1875LUPI_HOST_DEVICE double ProtareSingle::reactionCrossSection( std::size_t a_reactionIndex, URR_protareInfos const &a_URR_protareInfos, double a_temperature, double a_energy ) const {
1876
1877 if( m_continuousEnergy ) return( m_heatedCrossSections.reactionCrossSection( a_reactionIndex, a_URR_protareInfos, m_URR_index, a_temperature, a_energy ) );
1878
1879 return( m_heatedMultigroupCrossSections.reactionCrossSection( a_reactionIndex, a_temperature, a_energy ) );
1880}
1881
1882/* *********************************************************************************************************//**
1883 * Returns the index of a sampled reaction for a target with termpature *a_temperature*, a projectile with energy *a_energy* and total cross section
1884 * *a_crossSection*. Random numbers are obtained via *a_userrng* and *a_rngState*.
1885 *
1886 * @param a_hashIndex [in] Specifies the continuous energy or multi-group index.
1887 * @param a_temperature [in] The temperature of the target.
1888 * @param a_energy [in] The energy of the projectile.
1889 ***********************************************************************************************************/
1890
1891LUPI_HOST_DEVICE double ProtareSingle::depositionEnergy( std::size_t a_hashIndex, double a_temperature, double a_energy ) const {
1892
1893 if( m_continuousEnergy ) return( m_heatedCrossSections.depositionEnergy( a_hashIndex, a_temperature, a_energy ) );
1894
1895 return( m_heatedMultigroupCrossSections.depositionEnergy( a_hashIndex, a_temperature ) );
1896}
1897
1898/* *********************************************************************************************************//**
1899 * Returns the index of a sampled reaction for a target with termpature *a_temperature*, a projectile with energy *a_energy* and total cross section
1900 * *a_crossSection*. Random numbers are obtained via *a_userrng* and *a_rngState*.
1901 *
1902 * @param a_hashIndex [in] Specifies the continuous energy or multi-group index.
1903 * @param a_temperature [in] The temperature of the target.
1904 * @param a_energy [in] The energy of the projectile.
1905 ***********************************************************************************************************/
1906
1907LUPI_HOST_DEVICE double ProtareSingle::depositionMomentum( std::size_t a_hashIndex, double a_temperature, double a_energy ) const {
1908
1909 if( m_continuousEnergy ) return( m_heatedCrossSections.depositionMomentum( a_hashIndex, a_temperature, a_energy ) );
1910
1911 return( m_heatedMultigroupCrossSections.depositionMomentum( a_hashIndex, a_temperature ) );
1912}
1913
1914/* *********************************************************************************************************//**
1915 * Returns the index of a sampled reaction for a target with termpature *a_temperature*, a projectile with energy *a_energy* and total cross section
1916 * *a_crossSection*. Random numbers are obtained via *a_userrng* and *a_rngState*.
1917 *
1918 * @param a_hashIndex [in] Specifies the continuous energy or multi-group index.
1919 * @param a_temperature [in] The temperature of the target.
1920 * @param a_energy [in] The energy of the projectile.
1921 ***********************************************************************************************************/
1922
1923LUPI_HOST_DEVICE double ProtareSingle::productionEnergy( std::size_t a_hashIndex, double a_temperature, double a_energy ) const {
1924
1925 if( m_continuousEnergy ) return( m_heatedCrossSections.productionEnergy( a_hashIndex, a_temperature, a_energy ) );
1926
1927 return( m_heatedMultigroupCrossSections.productionEnergy( a_hashIndex, a_temperature ) );
1928}
1929
1930/* *********************************************************************************************************//**
1931 * Returns the index of a sampled reaction for a target with termpature *a_temperature*, a projectile with energy *a_energy* and total cross section
1932 * *a_crossSection*. Random numbers are obtained via *a_userrng* and *a_rngState*.
1933 *
1934 * @param a_hashIndex [in] Specifies the continuous energy or multi-group index.
1935 * @param a_temperature [in] The temperature of the target.
1936 * @param a_energy [in] The energy of the projectile.
1937 * @param a_particleIndex [in] The index of the particle whose gain is to be returned.
1938 ***********************************************************************************************************/
1939
1940LUPI_HOST_DEVICE double ProtareSingle::gain( std::size_t a_hashIndex, double a_temperature, double a_energy, int a_particleIndex ) const {
1941
1942 if( m_continuousEnergy ) return( m_heatedCrossSections.gain( a_hashIndex, a_temperature, a_energy, a_particleIndex ) );
1943
1944 return( m_heatedMultigroupCrossSections.gain( a_hashIndex, a_temperature, a_particleIndex ) );
1945}
1946
1947/* *********************************************************************************************************//**
1948 * Returns the intid of a sampled reaction for a target with termpature *a_temperature*, a projectile with energy *a_energy* and total cross section
1949 * *a_crossSection*. Random numbers are obtained via *a_userrng* and *a_rngState*.
1950 *
1951 * @param a_hashIndex [in] Specifies the continuous energy or multi-group index.
1952 * @param a_temperature [in] The temperature of the target.
1953 * @param a_energy [in] The energy of the projectile.
1954 * @param a_particleIntid [in] The intid of the particle whose gain is to be returned.
1955 ***********************************************************************************************************/
1956
1957LUPI_HOST_DEVICE double ProtareSingle::gainViaIntid( std::size_t a_hashIndex, double a_temperature, double a_energy, int a_particleIntid ) const {
1958
1959 if( m_continuousEnergy ) return( m_heatedCrossSections.gainViaIntid( a_hashIndex, a_temperature, a_energy, a_particleIntid ) );
1960
1961 return( m_heatedMultigroupCrossSections.gainViaIntid( a_hashIndex, a_temperature, a_particleIntid ) );
1962}
1963
1964/* *********************************************************************************************************//**
1965 * This method serializes *this* for broadcasting as needed for MPI and GPUs. The method can count the number of required
1966 * bytes, pack *this* or unpack *this* depending on *a_mode*.
1967 *
1968 * @param a_buffer [in] The buffer to read or write data to depending on *a_mode*.
1969 * @param a_mode [in] Specifies the action of this method.
1970 ***********************************************************************************************************/
1971
1973
1974 std::size_t vectorSize;
1975 LUPI::DataBuffer *workingBuffer = &a_buffer;
1976
1977 DATA_MEMBER_STRING( m_interaction, a_buffer, a_mode );
1978 DATA_MEMBER_INT( m_URR_index, a_buffer, a_mode );
1979 DATA_MEMBER_CAST( m_hasURR_probabilityTables, a_buffer, a_mode, bool );
1980 DATA_MEMBER_DOUBLE( m_URR_domainMin, a_buffer, a_mode );
1981 DATA_MEMBER_DOUBLE( m_URR_domainMax, a_buffer, a_mode );
1982 m_domainHash.serialize( a_buffer, a_mode );
1983 DATA_MEMBER_VECTOR_DOUBLE( m_projectileMultiGroupBoundaries, a_buffer, a_mode );
1984 DATA_MEMBER_VECTOR_DOUBLE( m_projectileMultiGroupBoundariesCollapsed, a_buffer, a_mode );
1985 DATA_MEMBER_CAST( m_upscatterModelASupported, a_buffer, a_mode, bool );
1986 DATA_MEMBER_VECTOR_DOUBLE( m_upscatterModelAGroupEnergies, a_buffer, a_mode );
1987 DATA_MEMBER_VECTOR_DOUBLE( m_upscatterModelAGroupVelocities, a_buffer, a_mode );
1988 DATA_MEMBER_VECTOR_DOUBLE( m_upscatterModelACrossSection, a_buffer, a_mode );
1989 m_multiGroupHash.serialize( a_buffer, a_mode );
1990
1991 vectorSize = m_nuclideGammaBranchStateInfos.size( );
1992 int vectorSizeInt = (int) vectorSize;
1993 DATA_MEMBER_INT( vectorSizeInt, *workingBuffer, a_mode );
1994 vectorSize = (std::size_t) vectorSizeInt;
1995
1996 if( a_mode == LUPI::DataBuffer::Mode::Unpack ) {
1997 m_nuclideGammaBranchStateInfos.resize( vectorSize, &(workingBuffer->m_placement) );
1998 for( std::size_t vectorIndex = 0; vectorIndex < vectorSize; ++vectorIndex ) {
1999 if (workingBuffer->m_placement != nullptr) {
2000 m_nuclideGammaBranchStateInfos[vectorIndex] = new(workingBuffer->m_placement) NuclideGammaBranchStateInfo;
2001 workingBuffer->incrementPlacement( sizeof( NuclideGammaBranchStateInfo ) );
2002 }
2003 else {
2004 m_nuclideGammaBranchStateInfos[vectorIndex] = new NuclideGammaBranchStateInfo;
2005 }
2006 }
2007 }
2008 if( a_mode == LUPI::DataBuffer::Mode::Memory ) {
2009 a_buffer.m_placement += m_nuclideGammaBranchStateInfos.internalSize();
2010 a_buffer.incrementPlacement( sizeof( NuclideGammaBranchStateInfo ) * vectorSize );
2011 }
2012 for( std::size_t vectorIndex = 0; vectorIndex < vectorSize; ++vectorIndex ) {
2013 m_nuclideGammaBranchStateInfos[vectorIndex]->serialize( *workingBuffer, a_mode );
2014 }
2015
2016 vectorSize = m_branches.size( );
2017 vectorSizeInt = (int) vectorSize;
2018 DATA_MEMBER_INT( vectorSizeInt, *workingBuffer, a_mode );
2019 vectorSize = (std::size_t) vectorSizeInt;
2020
2021 if( a_mode == LUPI::DataBuffer::Mode::Unpack ) {
2022 m_branches.resize( vectorSize, &(workingBuffer->m_placement) );
2023 for( std::size_t vectorIndex = 0; vectorIndex < vectorSize; ++vectorIndex ) {
2024 if (workingBuffer->m_placement != nullptr) {
2025 m_branches[vectorIndex] = new(workingBuffer->m_placement) NuclideGammaBranchInfo;
2026 workingBuffer->incrementPlacement( sizeof( NuclideGammaBranchInfo ) );
2027 }
2028 else {
2029 m_branches[vectorIndex] = new NuclideGammaBranchInfo;
2030 }
2031 }
2032 }
2033 if( a_mode == LUPI::DataBuffer::Mode::Memory ) {
2034 a_buffer.m_placement += m_branches.internalSize();
2035 workingBuffer->incrementPlacement( sizeof( NuclideGammaBranchInfo ) * vectorSize );
2036 }
2037 for( std::size_t vectorIndex = 0; vectorIndex < vectorSize; ++vectorIndex ) {
2038 m_branches[vectorIndex]->serialize( *workingBuffer, a_mode );
2039 }
2040
2041 vectorSize = m_reactions.size( );
2042 vectorSizeInt = (int) vectorSize;
2043 DATA_MEMBER_INT( vectorSizeInt, *workingBuffer, a_mode );
2044 vectorSize = (std::size_t) vectorSizeInt;
2045
2046 if( a_mode == LUPI::DataBuffer::Mode::Unpack ) {
2047 m_reactions.resize( vectorSize, &(workingBuffer->m_placement) );
2048 for( std::size_t vectorIndex = 0; vectorIndex < vectorSize; ++vectorIndex ) {
2049 if (workingBuffer->m_placement != nullptr) {
2050 m_reactions[vectorIndex] = new(workingBuffer->m_placement) Reaction;
2051 workingBuffer->incrementPlacement( sizeof(Reaction));
2052 }
2053 else {
2054 m_reactions[vectorIndex] = new Reaction;
2055 }
2056 }
2057 }
2058 if( a_mode == LUPI::DataBuffer::Mode::Memory ) {
2059 a_buffer.m_placement += m_reactions.internalSize();
2060 a_buffer.incrementPlacement( sizeof(Reaction) * vectorSize);
2061 }
2062 for( std::size_t vectorIndex = 0; vectorIndex < vectorSize; ++vectorIndex ) {
2063 m_reactions[vectorIndex]->serialize( *workingBuffer, a_mode );
2064 m_reactions[vectorIndex]->updateProtareSingleInfo( this, vectorIndex );
2065 }
2066
2067 vectorSize = m_orphanProducts.size( );
2068 vectorSizeInt = (int) vectorSize;
2069 DATA_MEMBER_INT( vectorSizeInt, *workingBuffer, a_mode );
2070 vectorSize = (std::size_t) vectorSizeInt;
2071
2072 if( a_mode == LUPI::DataBuffer::Mode::Unpack ) {
2073 m_orphanProducts.resize( vectorSize, &(workingBuffer->m_placement) );
2074 for( std::size_t vectorIndex = 0; vectorIndex < vectorSize; ++vectorIndex ) {
2075 if (workingBuffer->m_placement != nullptr) {
2076 m_orphanProducts[vectorIndex] = new(workingBuffer->m_placement) Reaction;
2077 workingBuffer->incrementPlacement( sizeof(Reaction));
2078 }
2079 else {
2080 m_orphanProducts[vectorIndex] = new Reaction;
2081 }
2082 }
2083 }
2084
2085 if( a_mode == LUPI::DataBuffer::Mode::Memory ) {
2086 a_buffer.m_placement += m_orphanProducts.internalSize( );
2087 a_buffer.incrementPlacement( sizeof( Reaction ) * vectorSize );
2088 }
2089
2090 for( std::size_t vectorIndex = 0; vectorIndex < vectorSize; ++vectorIndex ) {
2091 m_orphanProducts[vectorIndex]->serialize( *workingBuffer, a_mode );
2092 m_orphanProducts[vectorIndex]->updateProtareSingleInfo( this, vectorIndex );
2093 }
2094
2095 if( a_mode == LUPI::DataBuffer::Mode::Unpack ) {
2096 for( auto reactionIter = m_reactions.begin( ); reactionIter != m_reactions.end( ); ++reactionIter ) {
2097 (*reactionIter)->addOrphanProductToProductList( m_orphanProducts );
2098 }
2099 }
2100
2101 DATA_MEMBER_CAST( m_isPhotoAtomic, *workingBuffer, a_mode, bool );
2102 DATA_MEMBER_CAST( m_continuousEnergy, *workingBuffer, a_mode, bool );
2103 DATA_MEMBER_CAST( m_fixedGrid, *workingBuffer, a_mode, bool );
2104 m_heatedCrossSections.serialize( *workingBuffer, a_mode );
2105 m_heatedMultigroupCrossSections.serialize( *workingBuffer, a_mode );
2106}
2107
2108}
#define GIDI_multiGroupStyleChars
Definition GIDI.hpp:246
#define GIDI_MapInteractionTNSLChars
Definition GIDI.hpp:5165
#define GIDI_MapInteractionAtomicChars
Definition GIDI.hpp:5164
#define GIDI_heatedMultiGroupStyleChars
Definition GIDI.hpp:253
#define GIDI_crossSectionSumChars
Definition GIDI.hpp:210
#define GIDI_SnElasticUpScatterStyleChars
Definition GIDI.hpp:254
#define DATA_MEMBER_STRING(member, buf, mode)
#define DATA_MEMBER_VECTOR_INT(member, buf, mode)
#define DATA_MEMBER_CAST(member, buf, mode, someType)
#define DATA_MEMBER_VECTOR_DOUBLE(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
#define MCGIDI_particleBeta(a_mass_unitOfEnergy, a_kineticEnergy)
Definition MCGIDI.hpp:71
std::string const & label() const
Definition GIDI.hpp:658
std::string const & ID() const
Definition GIDI.hpp:760
void updateReactionIndices(std::size_t a_offset) const
std::string const & interaction(LUPI_maybeUnused std::size_t a_index=0) const
Definition GIDI.hpp:4787
bool isPhotoAtomic() const
Definition GIDI.hpp:4746
ProtareSingle * protare(std::size_t a_index)
PoPI::NuclideGammaBranchStateInfos const & nuclideGammaBranchStateInfos() const
Definition GIDI.hpp:4731
virtual Vector multiGroupCrossSection(LUPI::StatusMessageReporting &a_smr, Transporting::MG const &a_settings, Styles::TemperatureInfo const &a_temperatureInfo, ExcludeReactionsSet const &a_reactionsToExclude=ExcludeReactionsSet {}, std::string const &a_label="") const =0
virtual Reaction * orphanProduct(std::size_t a_index)=0
virtual Styles::Suite & styles()=0
virtual ProtareType protareType() const =0
virtual Reaction * reaction(std::size_t a_index)=0
ParticleInfo const & projectile() const
Definition GIDI.hpp:4541
virtual LUPI::FormatVersion const & formatVersion(std::size_t a_index=0) const =0
virtual ProtareSingle * protare(std::size_t a_index)=0
bool active() const
Definition GIDI.hpp:4280
bool RutherfordScatteringPresent() const
Definition GIDI.hpp:4290
Vector multiGroupCrossSection(LUPI::StatusMessageReporting &a_smr, Transporting::MG const &a_settings, Styles::TemperatureInfo const &a_temperatureInfo, std::string const &a_label="") const
int ENDF_MT() const
Definition GIDI.hpp:4284
Component & crossSection()
Definition GIDI.hpp:4299
void productIDs(std::set< std::string > &a_ids, Transporting::Particles const &a_particles, bool a_transportablesOnly) const
double crossSectionThreshold() const
Definition GIDI.hpp:4324
std::vector< double > groupBoundaries(std::string const &a_ID) const
GIDI::Suite const & transportables() const
Definition GIDI.hpp:3323
std::vector< iterator > findAllOfMoniker(std::string const &a_moniker)
T * get(std::size_t a_Index)
Definition GIDI.hpp:2642
std::string const & href() const
Definition GIDI.hpp:4393
std::vector< double > const & boundaries() const
Definition GIDI.hpp:3501
MultiGroup multiGroup() const
Definition GIDI.hpp:3668
bool add(Particle const &a_particle)
bool hasParticle(std::string const &a_id) const
std::size_t size() const
Definition GIDI_data.hpp:79
Ancestry * findInAncestry(std::string const &a_href)
Ancestry * ancestor()
Definition GUPI.hpp:104
std::string const & moniker() const
Definition GUPI.hpp:102
LUPI_HOST_DEVICE void incrementPlacement(std::size_t a_delta)
LUPI_HOST_DEVICE ProtareSingle()
LUPI_HOST_DEVICE bool inURR(double a_energy) const
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 double productionEnergy(std::size_t a_hashIndex, double a_temperature, double a_energy) const
LUPI_HOST Vector< double > const & projectileMultiGroupBoundaries() const
Definition MCGIDI.hpp:1689
LUPI_HOST_DEVICE Vector< double > temperatures(std::size_t a_index=0) const
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 void setUserParticleIndexViaIntid2(int a_particleIntid, int a_userParticleIndex)
LUPI_HOST_DEVICE double depositionMomentum(std::size_t a_hashIndex, double a_temperature, double a_energy) const
LUPI_HOST_DEVICE void serialize2(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
LUPI_HOST_DEVICE double gain(std::size_t a_hashIndex, double a_temperature, double a_energy, int a_particleIndex) const
LUPI_HOST_DEVICE bool isPhotoAtomic() const
Definition MCGIDI.hpp:1650
LUPI_HOST_DEVICE void crossSectionVector(double a_temperature, double a_userFactor, std::size_t a_numberAllocated, double *a_crossSectionVector) const
LUPI_HOST void setUserParticleIndex2(int a_particleIndex, int a_userParticleIndex)
LUPI_HOST_DEVICE ~ProtareSingle()
LUPI_HOST_DEVICE bool hasIncoherentDoppler() const
LUPI_HOST Vector< double > const & projectileMultiGroupBoundariesCollapsed() const
Definition MCGIDI.hpp:1691
LUPI_HOST_DEVICE double gainViaIntid(std::size_t a_hashIndex, double a_temperature, double a_energy, int a_particleIntid) const
LUPI_HOST_DEVICE double depositionEnergy(std::size_t a_hashIndex, double a_temperature, double a_energy) const
LUPI_HOST_DEVICE String interaction() const
Definition MCGIDI.hpp:1700
LUPI_HOST_DEVICE std::size_t numberOfReactions() const
Definition MCGIDI.hpp:1694
LUPI_HOST_DEVICE ProtareSingle const * protare(std::size_t a_index) const
LUPI_HOST_DEVICE ProtareSingle const * protareWithReaction(std::size_t a_index) const
LUPI_HOST_DEVICE bool hasFission() const
LUPI_HOST_DEVICE int projectileIntid() const
Definition MCGIDI.hpp:1517
LUPI_HOST void setUserParticleIndex(int a_particleIndex, int a_userParticleIndex)
LUPI_HOST_DEVICE ProtareType protareType() const
Definition MCGIDI.hpp:1514
LUPI_HOST Vector< int > const & productIndices(bool a_transportablesOnly) const
MCGIDI_VIRTUAL_FUNCTION LUPI_HOST_DEVICE double depositionMomentum(std::size_t a_hashIndex, double a_temperature, double a_energy) const MCGIDI_TRUE_VIRTUAL
LUPI_HOST_DEVICE void serialize(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
MCGIDI_VIRTUAL_FUNCTION LUPI_HOST_DEVICE bool hasIncoherentDoppler() const MCGIDI_TRUE_VIRTUAL
LUPI_HOST void setUserParticleIndexViaIntid(int a_particleIntid, int a_userParticleIndex)
MCGIDI_VIRTUAL_FUNCTION LUPI_HOST_DEVICE ProtareSingle const * protareWithReaction(std::size_t a_index) const MCGIDI_TRUE_VIRTUAL
LUPI_HOST_DEVICE long memorySize()
friend ProtareTNSL
Definition MCGIDI.hpp:1601
MCGIDI_VIRTUAL_FUNCTION LUPI_HOST Vector< double > const & projectileMultiGroupBoundariesCollapsed() const MCGIDI_TRUE_VIRTUAL
MCGIDI_VIRTUAL_FUNCTION LUPI_HOST_DEVICE Vector< double > temperatures(std::size_t a_index=0) const MCGIDI_TRUE_VIRTUAL
MCGIDI_VIRTUAL_FUNCTION LUPI_HOST_DEVICE double minimumEnergy() const MCGIDI_TRUE_VIRTUAL
LUPI_HOST_DEVICE long sizeOf() const
LUPI_HOST Vector< int > const & userProductIndices(bool a_transportablesOnly) const
MCGIDI_VIRTUAL_FUNCTION LUPI_HOST_DEVICE double depositionEnergy(std::size_t a_hashIndex, double a_temperature, double a_energy) const MCGIDI_TRUE_VIRTUAL
friend ProtareComposite
Definition MCGIDI.hpp:1600
MCGIDI_VIRTUAL_FUNCTION 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 MCGIDI_TRUE_VIRTUAL
MCGIDI_VIRTUAL_FUNCTION LUPI_HOST_DEVICE std::size_t numberOfReactions() const MCGIDI_TRUE_VIRTUAL
MCGIDI_VIRTUAL_FUNCTION LUPI_HOST_DEVICE Reaction const * reaction(std::size_t a_index) const MCGIDI_TRUE_VIRTUAL
MCGIDI_VIRTUAL_FUNCTION LUPI_HOST_DEVICE double gain(std::size_t a_hashIndex, double a_temperature, double a_energy, int a_particleIndex) const MCGIDI_TRUE_VIRTUAL
MCGIDI_VIRTUAL_FUNCTION LUPI_HOST_DEVICE double URR_domainMin() const MCGIDI_TRUE_VIRTUAL
MCGIDI_VIRTUAL_FUNCTION LUPI_HOST_DEVICE std::size_t numberOfOrphanProducts() const MCGIDI_TRUE_VIRTUAL
friend ProtareSingle
Definition MCGIDI.hpp:1599
MCGIDI_VIRTUAL_FUNCTION 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 MCGIDI_TRUE_VIRTUAL
MCGIDI_VIRTUAL_FUNCTION LUPI_HOST_DEVICE Reaction const * orphanProduct(std::size_t a_index) const MCGIDI_TRUE_VIRTUAL
MCGIDI_VIRTUAL_FUNCTION LUPI_HOST_DEVICE double URR_domainMax() const MCGIDI_TRUE_VIRTUAL
LUPI_HOST_DEVICE void serializeCommon(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
LUPI_HOST_DEVICE void serialize2(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
MCGIDI_VIRTUAL_FUNCTION LUPI_HOST_DEVICE std::size_t numberOfProtares() const MCGIDI_TRUE_VIRTUAL
MCGIDI_VIRTUAL_FUNCTION LUPI_HOST_DEVICE Vector< double > const & upscatterModelAGroupVelocities() const MCGIDI_TRUE_VIRTUAL
MCGIDI_VIRTUAL_FUNCTION LUPI_HOST_DEVICE bool hasURR_probabilityTables() const MCGIDI_TRUE_VIRTUAL
MCGIDI_VIRTUAL_FUNCTION LUPI_HOST Vector< double > const & projectileMultiGroupBoundaries() const MCGIDI_TRUE_VIRTUAL
LUPI_HOST Vector< int > const & productIntids(bool a_transportablesOnly) const
MCGIDI_VIRTUAL_FUNCTION LUPI_HOST_DEVICE bool hasFission() const MCGIDI_TRUE_VIRTUAL
LUPI_HOST_DEVICE Protare(ProtareType a_protareType)
MCGIDI_VIRTUAL_FUNCTION LUPI_HOST_DEVICE int URR_index() const MCGIDI_TRUE_VIRTUAL
MCGIDI_VIRTUAL_FUNCTION LUPI_HOST_DEVICE double threshold(std::size_t a_index) const MCGIDI_TRUE_VIRTUAL
LUPI_HOST_DEVICE void incrementMemorySize(long &a_totalMemory, long &a_sharedMemory)
virtual LUPI_HOST_DEVICE ~Protare()
LUPI_HOST_DEVICE bool isTNSL_ProtareSingle() const
Definition MCGIDI.hpp:1541
MCGIDI_VIRTUAL_FUNCTION LUPI_HOST_DEVICE bool reactionHasURR_probabilityTables(std::size_t a_index) const MCGIDI_TRUE_VIRTUAL
MCGIDI_VIRTUAL_FUNCTION LUPI_HOST_DEVICE double maximumEnergy() const MCGIDI_TRUE_VIRTUAL
MCGIDI_VIRTUAL_FUNCTION LUPI_HOST_DEVICE double gainViaIntid(std::size_t a_hashIndex, double a_temperature, double a_energy, int a_particleIntid) const MCGIDI_TRUE_VIRTUAL
MCGIDI_VIRTUAL_FUNCTION LUPI_HOST_DEVICE double productionEnergy(std::size_t a_hashIndex, double a_temperature, double a_energy) const MCGIDI_TRUE_VIRTUAL
MCGIDI_VIRTUAL_FUNCTION LUPI_HOST_DEVICE ProtareSingle const * protare(std::size_t a_index) const MCGIDI_TRUE_VIRTUAL
MCGIDI_VIRTUAL_FUNCTION LUPI_HOST_DEVICE void crossSectionVector(double a_temperature, double a_userFactor, std::size_t a_numberAllocated, double *a_crossSectionVector) const MCGIDI_TRUE_VIRTUAL
std::map< std::string, int > m_stateNamesToIndices
Definition MCGIDI.hpp:277
std::map< std::string, double > m_nuclearLevelEnergies
Definition MCGIDI.hpp:278
LUPI_HOST_DEVICE const char * c_str() const
raw data
LUPI_HOST_DEVICE void push_back(const T &dataElem)
LUPI_HOST_DEVICE void resize(std::size_t s, char **address=nullptr, bool mem_flag=CPU_MEM)
LUPI_HOST_DEVICE void reserve(std::size_t s, char **address=nullptr, bool mem_flag=CPU_MEM)
double nuclearLevelEnergy() const
Definition PoPI.hpp:592
std::size_t size() const
Definition PoPI.hpp:622
NuclideGammaBranchStateInfo * find(std::string const &a_state)
std::vector< Styles::TemperatureInfo > TemperatureInfos
Definition GIDI.hpp:3440
Definition GIDI.hpp:32
Vector multiGroupXYs1d(Transporting::MultiGroup const &a_boundaries, Functions::XYs1d const &a_function, Transporting::Flux const &a_flux)
@ centerOfMass
Definition GIDI.hpp:146
Vector collapse(Vector const &a_vector, Transporting::Settings const &a_settings, Transporting::Particles const &a_particles, double a_temperature)
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 int MCGIDI_popsIntid(PoPI::Database const &a_pops, std::string const &a_ID)
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
LUPI_HOST int MCGIDI_popsIndex(PoPI::Database const &a_pops, std::string const &a_ID)
LUPI_HOST void convertACE_URR_probabilityTablesFromGIDI(GIDI::ProtareSingle const &a_protare, Transporting::MC &a_settings, SetupInfo &a_setupInfo)
Definition PoPI.hpp:28
static std::string const photon
Definition PoPI.hpp:162
static std::string const electron
Definition PoPI.hpp:163
static int constexpr photon
Definition PoPI.hpp:178
static int constexpr neutron
Definition PoPI.hpp:177