Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
PoPI.hpp
Go to the documentation of this file.
1/*
2# <<BEGIN-copyright>>
3# Copyright 2019, Lawrence Livermore National Security, LLC.
4# This file is part of the gidiplus package (https://github.com/LLNL/gidiplus).
5# gidiplus is licensed under the MIT license (see https://opensource.org/licenses/MIT).
6# SPDX-License-Identifier: MIT
7# <<END-copyright>>
8*/
9
10#ifndef PoPI_hpp_included
11#define PoPI_hpp_included 1
12
13#include <string>
14#include <map>
15#include <vector>
16#include <list>
17#include <iostream>
18#include <stdexcept>
19#include <typeinfo>
20#include <fstream>
21#include <exception>
22#include <utility>
23#include <stddef.h>
24
25#include <LUPI.hpp>
26#include <HAPI.hpp>
27
28namespace PoPI {
29
30#define PoPI_AMU2MeV_c2 931.494028
31#define PoPI_electronMass_MeV_c2 0.5109989461
32
33#define PoPI_formatVersion_0_1_Chars "0.1"
34#define PoPI_formatVersion_1_10_Chars "1.10"
35#define PoPI_formatVersion_2_0_LLNL_3_Chars "2.0.LLNL_3"
36
37#define PoPI_PoPsChars "PoPs"
38
39#define PoPI_idChars "id"
40#define PoPI_symbolChars "symbol"
41#define PoPI_chemicalElementsChars "chemicalElements"
42#define PoPI_chemicalElementChars "chemicalElement"
43#define PoPI_isotopesChars "isotopes"
44#define PoPI_isotopeChars "isotope"
45#define PoPI_gaugeBosonChars "gaugeBoson"
46#define PoPI_leptonChars "lepton"
47#define PoPI_baryonChars "baryon"
48#define PoPI_nuclidesChars "nuclides"
49#define PoPI_nuclideChars "nuclide"
50#define PoPI_nucleusChars "nucleus"
51#define PoPI_unorthodoxChars "unorthodox"
52#define PoPI_aliasesChars "aliases"
53
54/*! \enum Particle_class
55 * This enum represents the various type of allowed particle types.
56 */
57
58enum class Particle_class { nuclide, /**< Specifies that the particle is a nuclide. */
59 nucleus, /**< Specifies that the particle is a nucleus. */
60 gaugeBoson, /**< Specifies that the particle is a gauge boson. */
61 lepton, /**< Specifies that the particle is a lepton. */
62 baryon, /**< Specifies that the particle is a baryon. */
63 nuclideMetaStable, /**< Specifies that the particle is a nuclide meta-stable alias. */
64 nucleusMetaStable, /**< Specifies that the particle is a nucleus meta-stable alias. */
65 TNSL, /**< Specifies that the particle is a TNSL target. Currently not used. */
66 ENDL_fissionProduct, /**< Specifies that the particle is an ENDL fissiont product (e.g., 99120, 99125). */
67 unorthodox, /**< Specifies that the particle is an unorthodox. */
68 alias, /**< Specifies that the particle is a alias. */
69 chemicalElement, /**< Specifies that the particle is a chemicalElement. */
70 isotope, /**< Specifies that the particle is a isotope. */
71 unknown /**< Specifies that the particle is a unknown. */ };
72
73#define PoPI_massChars "mass"
74#define PoPI_spinChars "spin"
75#define PoPI_parityChars "parity"
76#define PoPI_chargeChars "charge"
77#define PoPI_halflifeChars "halflife"
78
79#define PoPI_doubleChars "double"
80#define PoPI_integerChars "integer"
81#define PoPI_fractionChars "fraction"
82#define PoPI_stringChars "string"
83#define PoPI_shellChars "shell"
84#define PoPI_decayDataChars "decayData"
85#define PoPI_gammaDecayDataChars "gammaDecayData"
86
87#define PoPI_decayModeElectroMagnetic "electroMagnetic"
88
89#define PoPI_formatChars "format"
90#define PoPI_labelChars "label"
91#define PoPI_indexChars "index"
92#define PoPI_pidChars "pid"
93#define PoPI_nameChars "name"
94#define PoPI_versionChars "version"
95#define PoPI_aliasChars "alias"
96#define PoPI_metaStableChars "metaStable"
97#define PoPI_particleChars "particle"
98#define PoPI_discreteChars "discrete"
99#define PoPI_continuumChars "continuum"
100
101/*! \enum PQ_class
102 * This enum represents the various type of allowed physcial quantity types.
103 */
104
105enum class PQ_class { Double, /**< Specifies that the physcial quantity is a double. */
106 integer, /**< Specifies that the physcial quantity is an integer. */
107 fraction, /**< Specifies that the physcial quantity is a fraction. */
108 string, /**< Specifies that the physcial quantity is a string. */
109 shell /**< Specifies that the physcial quantity is a shell. */ };
110
111/*! \enum SpecialParticleID_mode
112 * This enum specifies how the light charged particle ids are handled. The light charged particles ids are familiarly known as
113 * p, d, t, h and a.
114 */
115
116enum class SpecialParticleID_mode { familiar, /**< Treat ids as the familiar p, d, t, h and a. */
117 nuclide, /**< Treat ids as the familiar p, d, t, h and a as h1, h2, h3, he3 and he4, respectively. */
118 nucleus /**< Treat ids as the familiar p, d, t, h and a as H1, H2, H3, He3 and He4, respectively. */ };
119
120class NuclideGammaBranchStateInfos;
121class Base;
122class SymbolBase;
123class Decay;
124class DecayMode;
125class DecayData;
126class Particle;
127class MetaStable;
128class Alias;
129class Baryon;
130class GaugeBoson;
131class Lepton;
132class Nuclide;
133class Nucleus;
134class Unorthodox;
135
136class Isotope;
137class ChemicalElement;
138class Database;
139
140void appendXMLEnd( std::vector<std::string> &a_XMLList, std::string const &a_label );
141
142int particleZ( Base const &a_particle, bool a_isNeutronProtonANucleon = false );
143int particleZ( Database const &a_pops, std::size_t a_index, bool a_isNeutronProtonANucleon = false );
144int particleZ( Database const &a_pops, std::string const &a_id, bool a_isNeutronProtonANucleon = false );
145
146int particleA( Base const &a_particle, bool a_isNeutronProtonANucleon = false );
147int particleA( Database const &a_pops, std::size_t a_index, bool a_isNeutronProtonANucleon = false );
148int particleA( Database const &a_pops, std::string const &a_id, bool a_isNeutronProtonANucleon = false );
149
150int particleZA( Base const &a_particle, bool a_isNeutronProtonANucleon = false );
151int particleZA( Database const &a_pops, std::size_t a_index, bool a_isNeutronProtonANucleon = false );
152int particleZA( Database const &a_pops, std::string const &a_id, bool a_isNeutronProtonANucleon = false );
153
154int particleMetaStableIndex( Base const &a_particle );
155int particleMetaStableIndex( Database const &a_pops, std::size_t a_index );
156int particleMetaStableIndex( Database const &a_pops, std::string const &a_id );
157
158std::string specialParticleID( SpecialParticleID_mode a_mode, std::string const &a_id );
159bool compareSpecialParticleIDs( std::string const &a_id1, std::string const &a_id2 );
160
161struct IDs {
162 static std::string const photon;
163 static std::string const electron;
164 static std::string const neutron;
165 static std::string const proton;
166 static std::string const familiarPhoton;
167 static std::string const familiarDeuteron;
168 static std::string const familiarTriton;
169 static std::string const familiarHelion;
170 static std::string const familiarAlpha;
171 static std::string const FissionProductENDL99120;
172 static std::string const FissionProductENDL99125;
173 static std::string const anti;
174};
175
176struct Intids {
177 static int constexpr neutron = 1020000000;
178 static int constexpr photon = 1000000000;
179 static int constexpr electron = 1010000000;
180 static int constexpr FissionProductENDL99120 = 1990099120;
181 static int constexpr FissionProductENDL99125 = 1990099125;
182};
183
184extern std::map<std::string, std::string> supportedNucleusAliases;
185
186typedef std::vector<Base *> ParticleList;
187typedef std::vector<SymbolBase *> SymbolList;
188
189/*
190============================================================
191======================== Exception =========================
192============================================================
193*/
194
195class Exception : public std::runtime_error {
196
197 public :
198 explicit Exception( std::string const &a_message );
199
200};
201
202/*
203============================================================
204====================== ParseIntidInfo ======================
205============================================================
206*/
207
209
210 private:
211 int m_intid; /**< The intid for the rest of the data. */
212 Particle_class m_family; /**< The family of the particle. */
213 bool m_isAnti; /**< **true** if particle is an anti-particle and **false** otherwise. */
214// The following are for nuclear like particles.
215 bool m_isNuclear; /**< *true* if the particle is a nuclear particle and *false* otherwise. */
216 int m_AAA; /**< For a nuclear particle, its AAA value (i.e., atomic mass number). */
217 int m_ZZZ; /**< For a nuclear particle, its ZZZ value (i.e., atomic number). */
218 int m_III; /**< For a nuclear particle, its III value (nuclear excitation level index or meta-stable index. */
219 int m_nuclearLevelIndex; /**< For a nuclear particle, nuclear excitation level index. */
220 int m_metaStableIndex; /**< For a nuclear meta-stable particle, its meta-stable index. */
221// The following are for leptons.
222 int m_generation; /**< For a lepton, its generation. */
223 bool m_isNeutrino; /**< For a lepton, **true** if leption is a neutrino and **false** otherwise. */
224// The follow are for baryons.
225 int m_baryonGroup; /**< For a baryon, its baryon group. */
226 int m_baryonId; /**< For a baryon, its id within a baryon group. */
227// Other data.
228 int m_familyId; /**< For non-nuclear particles, the particle's indentifier within its family. */
229
230 public:
231 ParseIntidInfo( int a_intid, bool a_GRIN_mode = false );
232
233 int intid( ) { return( m_intid ); } /**< Returns the value of the *m_intid* member. */
234 Particle_class family( ) { return( m_family ); } /**< Returns the value of the *m_family* member. */
235 bool isAnti( ) { return( m_isAnti ); } /**< Returns the value of the *m_isAnti * member. */
236
237 bool isNuclear( ) { return( m_isNuclear ); } /**< Returns the value of the *m_isNuclear* member. */
238 int AAA( ) { return( m_AAA ); } /**< Returns the value of the *m_AAA* member. */
239 int ZZZ( ) { return( m_ZZZ ); } /**< Returns the value of the *m_ZZZ* member. */
240 int III( ) { return( m_III ); } /**< Returns the value of the *m_III* member. */
241 bool isNuclearMetaStable( ) { return( ( m_family == Particle_class::nuclideMetaStable ) || ( m_family == Particle_class::nucleusMetaStable ) ); }
242 /**< Returns **true** if particle is a nuclear meta-stable alias. */
243 int metaStableIndex( ) { return( m_metaStableIndex ); } /**< Returns the value of the *m_metaStableIndex* member. */
244
245 int generation( ) { return( m_generation ); } /**< Returns the value of the *m_generation* member. */
246 bool isNeutrino( ) { return( m_isNeutrino ); } /**< Returns the value of the *m_isNeutrino* member. */
247
248 int baryonGroup( ) { return( m_baryonGroup ); } /**< Returns the value of the *m_baryonGroup* member. */
249 int baryonId( ) { return( m_baryonId ); } /**< Returns the value of the *m_baryonId * member. */
250
251 int familyId( ) { return( m_familyId ); } /**< Returns the value of the *m_familyId* member. */
252
253 std::string id( );
254};
255
256/*
257============================================================
258======================= ParseIdInfo ========================
259============================================================
260*/
261
263
264 private:
265 bool m_isSupported; /**< If **true** the particle's id was parsed and the other member of *this* are valid. Otherwise, parsing of the particle's id is currently not supported. */
266 std::string m_id; /**< The id for the particles. */
267 bool m_isNuclear; /**< **true** if particle is a valid nuclear name (i.e., nuclide, nucleus of meta-stable) and **false** otherwise. */
268 bool m_isNucleus; /**< **true** if particle is a nucleus and **false** otherwise. */
269 bool m_isChemicalElement; /**< **true** if id is only a chemical element symbol. */
270 bool m_isAnti; /**< **true** if particle is an anti-particle and **false** otherwise. */
271 bool m_isMetaStable; /**< **true** if particle is a meta-stable alias and **false** otherwise. */
272 std::string m_symbol; /**< The chemical element symbol part of *m_id*. This will always be the nuclide symbol even if *m_id* is for a nucleus. */
273 int m_Z; /**< The atomic number of the *m_id*. */
274 int m_A; /**< The atomic mass number of the *m_id*. */
275 int m_index; /**< The nuclear level index of *m_id*. */
276 std::string m_qualifier;
277
278 std::string boolToString( bool a_value, std::string const &a_prefix ) const;
279
280 public:
281 ParseIdInfo( std::string const &a_id );
282
283 bool isSupported( ) { return( m_isSupported ); } /**< Returns the value of the *m_isSupported. */
284 std::string const &Id( ) { return( m_id ); } /**< Returns a reference to the *m_id* member. */
285 bool isNuclear( ) { return( m_isNuclear ); } /**< Returns the value of the *m_isNuclear* member. */
286 bool isNucleus( ) { return( m_isNucleus ); } /**< Returns the value of the *m_isNucleus* member. */
287 bool isChemicalElement( ) { return( m_isChemicalElement ); } /**< Returns the value of the *m_isChemicalElement* member. */
288 bool isAnti( ) { return( m_isAnti ); } /**< Returns the value of the *m_isAnti* member. */
289 std::string const &symbol( ) { return( m_symbol); } /**< Returns a reference to the *m_symbol* member. */
290 int Z( ) { return( m_Z ); } /**< Returns the value of the *m_Z* member. */
291 int A( ) { return( m_A ); } /**< Returns the value of the *m_A* member. */
292 int index( ) { return( m_index ); } /**< Returns the value of the *m_index* member. */
293 std::string const &qualifier( ) { return( m_qualifier ); } /**< Returns a reference to the *m_qualifier* member. */
294
295 void print( bool a_terse, std::string const &a_indent = "" ) const ;
296};
297
298/*! \class Suite
299 * This is the base class for all suite like members.
300 */
301
302/*
303============================================================
304========================== Suite ===========================
305============================================================
306*/
307
308template <class T, class T2>
309class Suite {
310
311 private:
312 std::string m_moniker; /**< The moniker (i.e., name) of the suite. */
313 std::vector<T *> m_items; /**< The list of all items in the suite. */
314
315 public:
316 Suite( std::string const &a_moniker ) : m_moniker( a_moniker ) { }
318 void appendFromParentNode( HAPI::Node const &a_node, Database *a_DB, T2 *a_parent );
319 void appendFromParentNode2( HAPI::Node const &a_node, T2 *a_parent );
320
321 std::string::size_type size( void ) const { return( m_items.size( ) ); } /**< Returns the number of items in the suite. */
322 T &operator[]( std::size_t a_index ) const { return( *m_items[a_index] ); } /**< Returns the item at index *a_index*. */
323 std::string const &moniker( void ) { return( m_moniker ); } /**< Returns the value of the *m_moniker* member. */
324
325 void toXMLList( std::vector<std::string> &a_XMLList, std::string const &a_indent1 ) const ;
326};
327
328/* *********************************************************************************************************//**
329 ***********************************************************************************************************/
330
331template <class T, class T2>
333
334 std::string::size_type i1, _size = m_items.size( );
335
336 for( i1 = 0; i1 < _size; ++i1 ) delete m_items[i1];
337
338// Ask Adam why next line does not work.
339// for( std::vector<T *>::iterator iter = m_items.begin( ); iter != m_items.end( ); ++iter ) delete *iter;
340}
341
342/* *********************************************************************************************************//**
343 * Adds the children of *a_node* to the suite and to *a_DB*.
344 *
345 * @param a_node [in] The **HAPI::Node** to be parsed.
346 * @param a_DB [in] The **PoPI::Database** instance to add the constructed items to.
347 * @param a_parent [in] The parent suite that will contain *this*.
348 ***********************************************************************************************************/
349
350template <class T, class T2>
351void Suite<T, T2>::appendFromParentNode( HAPI::Node const &a_node, Database *a_DB, T2 *a_parent ) {
352
353 for( HAPI::Node child = a_node.first_child( ); !child.empty( ); child.to_next_sibling( ) ) {
354 T *item = new T( child, a_DB, a_parent );
355 m_items.push_back( item );
356 }
357}
358
359/* *********************************************************************************************************//**
360 * Adds the children of *a_node* to the suite.
361 *
362 * @param a_node [in] The **HAPI::Node** to be parsed.
363 * @param a_parent [in] The parent suite that will contain *this*.
364 ***********************************************************************************************************/
365
366template <class T, class T2>
367void Suite<T, T2>::appendFromParentNode2( HAPI::Node const &a_node, T2 *a_parent ) {
368
369 for( HAPI::Node child = a_node.first_child( ); !child.empty( ); child.to_next_sibling( ) ) {
370 T *item = new T( child, a_parent );
371 m_items.push_back( item );
372 }
373}
374
375/* *********************************************************************************************************//**
376 * Creates an XML representation of the suite.
377 *
378 * @param a_XMLList [in] The list the XML lines are added to.
379 * @param a_indent1 [in] The amount to indent the XML text.
380 ***********************************************************************************************************/
381
382template <class T, class T2>
383void Suite<T, T2>::toXMLList( std::vector<std::string> &a_XMLList, std::string const &a_indent1 ) const {
384
385 std::string::size_type _size = m_items.size( );
386 std::string indent2 = a_indent1 + " ";
387
388 if( _size == 0 ) return;
389
390 std::string header = a_indent1 + "<" + m_moniker + ">";
391 a_XMLList.push_back( std::move( header ) );
392 for( std::string::size_type i1 = 0; i1 < _size; ++i1 ) m_items[i1]->toXMLList( a_XMLList, indent2 );
393
394 appendXMLEnd( a_XMLList, m_moniker );
395}
396
397/*
398============================================================
399===================== PhysicalQuantity =====================
400============================================================
401*/
402
404
405 private:
406 PQ_class m_class; /**< The class for the physical quanity. */
407 std::string m_tag; /**< The name of the physical quanity. */
408 std::string m_label; /**< The label for the physical quanity. */
409 std::string m_valueString; /**< The string value of the physical quanity. */
410 std::string m_unit; /**< The unit of the physical quanity. */
411
412 public:
413 PhysicalQuantity( HAPI::Node const &a_node, PQ_class a_class );
414 virtual ~PhysicalQuantity( );
415
416 PQ_class Class( void ) const { return( m_class ); } /**< Returns the value of the *m_class* member. */
417 std::string const &tag( void ) const { return( m_tag ); } /**< Returns the value of the *m_tag* member. */
418 std::string const &label( void ) const { return( m_label ); } /**< Returns the value of the *m_label* member. */
419 std::string const &valueString( void ) const { return( m_valueString ); } /**< Returns the value of the *valueString* member. */
420 std::string const &unit( void ) const { return( m_unit ); } /**< Returns the value of the *m_unit* member. */
421
422 void toXMLList( std::vector<std::string> &a_XMLList, std::string const &a_indent1 ) const ;
423 virtual std::string valueToString( void ) const = 0;
424};
425
426/*
427============================================================
428========================= PQ_double ========================
429============================================================
430*/
431
433
434 private:
435 double m_value; /**< The double value of the physical quanity. */
436 void initialize( );
437
438 public:
439 PQ_double( HAPI::Node const &a_node );
440 PQ_double( HAPI::Node const &a_node, PQ_class a_class );
441 virtual ~PQ_double( );
442
443 double value( void ) const { return( m_value ); } /**< Returns the value of the *m_value* member. */
444 double value( char const *a_unit ) const ;
445 double value( std::string const &a_unit ) const { return( value( a_unit.c_str( ) ) ); }
446 /**< Returns the value of the *m_value* member in units of *a_unit*. */
447 virtual std::string valueToString( void ) const ;
448};
449
450/*
451============================================================
452========================= PQ_integer =======================
453============================================================
454*/
455
457
458 private:
459 int m_value; /**< The integer value of the physical quanity. */
460
461 public:
462 PQ_integer ( HAPI::Node const &a_node );
463 virtual ~PQ_integer( );
464
465 int value( void ) const { return( m_value ); } /**< Returns the value of the *m_value* member. */
466 int value( char const *a_unit ) const ;
467 int value( std::string const &a_unit ) const { return( value( a_unit.c_str( ) ) ); }
468 virtual std::string valueToString( void ) const ;
469};
470
471/*
472==============================================================
473========================= PQ_fraction ========================
474==============================================================
475*/
476
478
479 public:
480 PQ_fraction( HAPI::Node const &a_node );
481 virtual ~PQ_fraction( );
482
483 std::string value( void ) const ;
484 std::string value( char const *a_unit ) const ;
485 std::string value( std::string const &a_unit ) const { return( value( a_unit.c_str( ) ) ); }
486 /**< Returns the value of the *m_value* member in units of *a_unit*. */
487 virtual std::string valueToString( void ) const ;
488};
489
490/*
491============================================================
492========================= PQ_string ========================
493============================================================
494*/
495
497
498 public:
499 PQ_string( HAPI::Node const &a_node );
500 virtual ~PQ_string( );
501
502 std::string value( void ) const { return( valueString( ) ); } /**< Returns the value returned by calling the **valueString** methods. */
503 std::string value( char const *a_unit ) const ;
504 std::string value( std::string const &a_unit ) const { return( value( a_unit.c_str( ) ) ); }
505 /**< Returns the value of the *m_value* member in units of *a_unit*. */
506 virtual std::string valueToString( void ) const ;
507};
508
509/*
510============================================================
511========================= PQ_shell =========================
512============================================================
513*/
514
515class PQ_shell : public PQ_double {
516
517 public:
518 PQ_shell( HAPI::Node const &a_node );
519 ~PQ_shell( );
520};
521
522/*
523============================================================
524========================= PQ_suite =========================
525============================================================
526*/
527
528class PQ_suite : public std::vector<PhysicalQuantity *> {
529
530 private:
531 std::string m_label;
532
533 public:
534 PQ_suite( HAPI::Node const &a_node );
535 ~PQ_suite( );
536
537 std::string &label( void ) { return( m_label ); }
538
539 void toXMLList( std::vector<std::string> &a_XMLList, std::string const &a_indent1 ) const ;
540};
541
542/*
543============================================================
544================== NuclideGammaBranchInfo ==================
545============================================================
546*/
547
549
550 private:
551 double m_probability; /**< The probability that the level decays to state *m_residualState*. */
552 double m_photonEmissionProbability; /**< The conditional probability the the decay emitted a photon. */
553 double m_gammaEnergy; /**< The energy of the emitted photon. */
554 std::string m_residualState; /**< The state the residual is left in after photon decay. */
555
556 public:
557 NuclideGammaBranchInfo( double a_probability, double a_photonEmissionProbability, double a_gammaEnergy, std::string const &a_residualState );
558 NuclideGammaBranchInfo( NuclideGammaBranchInfo const &a_nuclideGammaBranchInfo );
559
560 double probability( ) const { return( m_probability ); } /**< Returns the value of the *m_probability* member. */
561 double photonEmissionProbability( ) const { return( m_photonEmissionProbability ); }
562 /**< Returns the value of the *m_photonEmissionProbability* member. */
563 double gammaEnergy( ) const { return( m_gammaEnergy ); } /**< Returns the value of the *m_gammaEnergy* member. */
564 std::string const &residualState( ) const { return( m_residualState ); } /**< Returns the value of the *m_residualState* member. */
565};
566
567/*
568==============================================================
569================= NuclideGammaBranchStateInfo ================
570==============================================================
571*/
572
574
575 private:
576 std::string m_state; /**< The inital state the decay starts from. */
577 int m_intid; /**< The intid for the inital state. */
578 std::string m_kind; /**< The kind of the particle. Currently can be 'discrete' or 'continuum'. */
579 double m_nuclearLevelEnergy; /**< The nuclear level excitation energy of the level (state). */
580 double m_nuclearLevelEnergyWidth; /**< This is 0.0 except for GRIN realized continuum levels where this is the energy width from this level to the next higher level. */
581 bool m_derivedCalculated; /**< For internal use to determine if other members have been set or not. */
582 double m_multiplicity; /**< The average number of photons emitted when transitioning from the initial to the final state. Data derived from m_branches data. */
583 double m_averageGammaEnergy; /**< The average energy per decay from the initial to the final state. Data derived from m_branches data. */
584 std::vector<NuclideGammaBranchInfo> m_branches;
585
586 public:
587 NuclideGammaBranchStateInfo( std::string const &a_state, int a_intid, std::string const &a_kind, double a_nuclearLevelEnergy );
588
589 std::string const &state( ) const { return( m_state ); } /**< Returns the value of the *m_state* member. */
590 int intid( ) const { return( m_intid ); } /**< Returns the value of the *m_intid* member. */
591 std::string const &kind( ) const { return( m_kind ); } /**< Returns the value of the *m_kind* member. */
592 double nuclearLevelEnergy( ) const { return( m_nuclearLevelEnergy ); } /**< Returns the value of the *m_nuclearLevelEnergy* member. */
593 double nuclearLevelEnergyWidth( ) const { return( m_nuclearLevelEnergyWidth ); }
594 /**< Returns the value of the *m_nuclearLevelEnergyWidth* member. */
595 void setNuclearLevelEnergyWidth( double a_nuclearLevelEnergyWidth ) { m_nuclearLevelEnergyWidth = a_nuclearLevelEnergyWidth; }
596 /**< Set the value of the *m_nuclearLevelEnergyWidth* member to *a_nuclearLevelEnergyWidth*. */
597 bool derivedCalculated( ) const { return( m_derivedCalculated ); } /**< Returns the value of the *m_derivedCalculated* member. */
598 double multiplicity( ) const { return( m_multiplicity ); } /**< Returns the value of the *m_multiplicity* member. */
599 double averageGammaEnergy( ) const { return( m_averageGammaEnergy ); } /**< Returns the value of the *m_averageGammaEnergy* member. */
600 std::vector<NuclideGammaBranchInfo> const &branches( ) const { return( m_branches ); }
601 /**< Returns a reference to the *m_branches* member. */
602
603 void add( NuclideGammaBranchInfo const &a_nuclideGammaBranchInfo );
604 void calculateDerivedData( NuclideGammaBranchStateInfos &a_nuclideGammaBranchStateInfos );
605};
606
607/*
608==============================================================
609================ NuclideGammaBranchStateInfos ================
610==============================================================
611*/
612
614
615 private:
616 std::vector<NuclideGammaBranchStateInfo *> m_nuclideGammaBranchStateInfos;
617
618 public:
621
622 std::size_t size( ) const { return( m_nuclideGammaBranchStateInfos.size( ) ); }
623 NuclideGammaBranchStateInfo *operator[]( std::size_t a_index ) { return( m_nuclideGammaBranchStateInfos[a_index] ); }
624 NuclideGammaBranchStateInfo const *operator[]( std::size_t a_index ) const { return( m_nuclideGammaBranchStateInfos[a_index] ); }
625 std::vector<NuclideGammaBranchStateInfo *> &nuclideGammaBranchStateInfos( ) { return( m_nuclideGammaBranchStateInfos ); }
626 void add( NuclideGammaBranchStateInfo *a_nuclideGammaBranchStateInfo );
627 NuclideGammaBranchStateInfo *find( std::string const &a_state );
628 NuclideGammaBranchStateInfo const *find( std::string const &a_state ) const ;
629};
630
631/*
632============================================================
633=========================== Base ===========================
634============================================================
635*/
636
637class Base {
638
639 private:
640 std::string m_id; /**< The **PoPs** id for the particle or **PoPs** symbol for a chemicalElement or isotope. */
641 Particle_class m_class; /**< The **Particle_class** for the particle, chemicalElement or isotope. */
642 std::size_t m_index; /**< The for the particle, chemicalElement or isotope. */
643 int m_intid; /**< The unique integer id for a particle or a meta-stable alias. For a non meta-stable alias, an isotope or chemical element, this is -1. */
644
645 void setIntid( int a_intid ) { m_intid = a_intid; } /**< Sets the value of the *m_intid* member to *a_intid*. */
646
647 public:
648 Base( std::string const &a_id, Particle_class a_class );
649 Base( HAPI::Node const &a_node, std::string const &a_label, Particle_class a_class );
650 virtual ~Base( );
651
652 std::string const &ID( void ) const { return( m_id ); } /**< Returns a *const* reference to the *m_id* member of *this*. */
653 std::size_t index( void ) const { return( m_index ); } /**< Returns the value of the *m_index* member of *this*. */
654 void setIndex( std::size_t a_index ) { m_index = a_index; } /**< Sets the value of the *m_index* member of *this* to *a_index*. */
655 int intid( ) const { return( m_intid ); } /**< Returns the value of the *m_intid* member. */
656 Particle_class Class( void ) const { return( m_class ); } /**< Returns the value of the *m_class* member of *this*. */
657 virtual bool isParticle( ) const { return( true ); } /**< Returns **true** if *this* is a **Particle** and **false** it *this* is a **ChemicalElement** or **Isotope** instance. */
658 bool isAlias( void ) const { return( ( m_class == Particle_class::alias ) || isMetaStableAlias( ) ); }
659 /**< Returns **true** if *this* is an **Alias** or **MetaStable** instance and **false** otherwise. */
660 bool isMetaStableAlias( void ) const { return( ( m_class == Particle_class::nuclideMetaStable ) || ( m_class == Particle_class::nucleusMetaStable ) ); }
661 /**< Returns **true** if *this* is a **MetaStable** instance and **false** otherwise. */
662
663 bool isGaugeBoson( ) const { return( m_class == Particle_class::gaugeBoson ); } /**< Returns **true** if *this* is a **GaugeBoson** instance and **false** otherwise. */
664 bool isLepton( ) const { return( m_class == Particle_class::lepton ); } /**< Returns **true** if *this* is a **Lepton** instance and **false** otherwise. */
665 bool isBaryon( ) const { return( m_class == Particle_class::baryon ); } /**< Returns **true** if *this* is a **Baryon** instance and **false** otherwise. */
666 bool isUnorthodox( ) const { return( m_class == Particle_class::unorthodox ); } /**< Returns **true** if *this* is a **Unorthodox** instance and **false** otherwise. */
667 bool isNucleus( ) const { return( m_class == Particle_class::nucleus ); } /**< Returns **true** if *this* is a **Nucleus** instance and **false** otherwise. */
668 bool isNuclide( ) const { return( m_class == Particle_class::nuclide ); } /**< Returns **true** if *this* is a **Nuclide** instance and **false** otherwise. */
669 bool isIsotope( ) const { return( m_class == Particle_class::isotope ); } /**< Returns **true** if *this* is a **Isotope** instance and **false** otherwise. */
670 bool isChemicalElement( ) const { return( m_class == Particle_class::chemicalElement ); }
671 /**< Returns **true** if *this* is a **ChemicalElement** instance and **false** otherwise. */
672
674 friend Alias;
675 friend Baryon;
677 friend Lepton;
678 friend Nucleus;
679 friend Nuclide;
681};
682
683/*
684============================================================
685========================== IDBase ==========================
686============================================================
687*/
688
689class IDBase : public Base {
690
691 public:
692 IDBase( std::string const &a_id, Particle_class a_class );
693 IDBase( HAPI::Node const &a_node, Particle_class a_class );
694 virtual ~IDBase( ); // BRB This should be virtual but I cannot get it to work without crashing.
695
696 std::size_t addToDatabase( Database *a_DB );
697 double massValue2( Database const &a_DB, std::string const &a_unit ) const ;
698};
699
700/*
701============================================================
702======================== SymbolBase ========================
703============================================================
704*/
705
706class SymbolBase : public Base {
707
708 public:
709 SymbolBase( HAPI::Node const &a_node, Particle_class a_class );
710 ~SymbolBase( );
711
712 std::string const &symbol( ) const { return( ID( ) ); } /**< Returns the value of the symbol. */
713
714 std::size_t addToSymbols( Database *a_DB );
715 bool isParticle( ) const { return( false ); }
716};
717
718/*
719============================================================
720========================= Product ==========================
721============================================================
722*/
723
724class Product {
725
726 private:
727 int m_id;
728 std::string m_pid;
729 std::string m_label;
730
731 public:
732 Product( HAPI::Node const &a_node, Decay *a_DB );
733 ~Product( );
734
735 int ID( ) const { return( m_id ); }
736 std::string const &pid( ) const { return( m_pid ); }
737 std::string const &label( ) const { return( m_label ); }
738
739 void toXMLList( std::vector<std::string> &a_XMLList, std::string const &a_indent1 ) const ;
740};
741
742/*
743============================================================
744========================== Decay ===========================
745============================================================
746*/
747
748class Decay {
749
750 private:
751 int m_index;
752 std::string m_mode;
753 bool m_complete;
754 Suite<Product, Decay> m_products;
755
756 public:
757 Decay( HAPI::Node const &a_node, DecayMode const *a_decayMode );
758 ~Decay( );
759
760 int index( void ) const { return( m_index ); }
761 std::string const &mode( ) const { return( m_mode ); }
762 bool complete( ) const { return( m_complete ); }
763 Suite<Product, Decay> const &products( void ) const { return( m_products ); }
764 void toXMLList( std::vector<std::string> &a_XMLList, std::string const &a_indent1 ) const ;
765};
766
767/*
768============================================================
769======================== DecayMode =========================
770============================================================
771*/
772
774
775 private:
776 std::string m_label;
777 std::string m_mode;
778 PQ_suite m_probability;
779 PQ_suite m_photonEmissionProbabilities;
780 Suite<Decay, DecayMode> m_decayPath;
781
782 public:
783 DecayMode( HAPI::Node const &a_node, DecayData const *a_decayData );
784 ~DecayMode( );
785
786 std::string const &label( ) const { return( m_label ); }
787 std::string const &mode( ) const { return( m_mode ); }
788 PQ_suite const &probability( ) const { return( m_probability ); }
789 PQ_suite const &photonEmissionProbabilities( ) const { return( m_photonEmissionProbabilities ); }
790 Suite<Decay, DecayMode> const &decayPath( ) const { return( m_decayPath ); }
791
792 void calculateNuclideGammaBranchStateInfo( PoPI::Database const &a_pops, NuclideGammaBranchStateInfo &a_nuclideGammaBranchStateInfo ) const ;
793 void toXMLList( std::vector<std::string> &a_XMLList, std::string const &a_indent1 ) const ;
794};
795
796/*
797============================================================
798======================== DecayData =========================
799============================================================
800*/
801
803
804 private:
805 Suite<DecayMode, DecayData> m_decayModes;
806
807 public:
808 DecayData( HAPI::Node const &a_node );
809 ~DecayData( );
810
811 Suite<DecayMode, DecayData> const &decayModes( void ) const { return( m_decayModes ); }
812
813 void calculateNuclideGammaBranchStateInfo( PoPI::Database const &a_pops, NuclideGammaBranchStateInfo &a_nuclideGammaBranchStateInfo ) const ;
814 void toXMLList( std::vector<std::string> &a_XMLList, std::string const &a_indent1 ) const ;
815};
816
817/*
818============================================================
819====================== GammaDecayData ======================
820============================================================
821*/
822
824
825 private:
826 std::string m_kind; /**< The kind of the particle. Currently can be 'discrete' or 'continuum' but may be 'experimental', 'evaluated' or 'modelled' in the future. */
827 int m_rows; /**< The number of nuclides listed. */
828 int m_columns; /**< The number of items in a row. */
829 std::vector<std::string> m_ids; /**< The list of nulcides. */
830 std::vector<double> m_probabilities; /**< The list of probabilities for each nuclide. This must sum to 1. */
831 std::vector<double> m_photonEmissionProbabilities; /**< The list of photon emission probabilities for each nuclide. */
832
833 public:
834 GammaDecayData( HAPI::Node const &a_node );
836
837 std::string const &kind( ) const { return( m_kind ); } /**< Returns a const reference to the *m_kind* member. */
838 int rows( ) const { return( m_rows ); }
839 int colunms( ) const { return( m_columns ); }
840 std::vector<std::string> const &ids( ) const { return( m_ids ); }
841 std::vector<double> const &probabilities( ) const { return( m_probabilities ); }
842 std::vector<double> const &photonEmissionProbabilities( ) const { return( m_photonEmissionProbabilities ); }
843
844 void calculateNuclideGammaBranchStateInfo( PoPI::Database const &a_pops, NuclideGammaBranchStateInfo &nuclideGammaBranchStateInfo ) const ;
845};
846
847/*
848============================================================
849========================= Particle =========================
850============================================================
851*/
852
853class Particle : public IDBase {
854
855 private:
856 std::string m_baseId; /**< The base part of the id (i.e., without the anti and quailifier). */
857 std::string m_family; /**< The family of the particle. */
858 std::string m_anti; /**< The string "_anti" if particle is an anti-particle and an empty string otherwise. */
859 int m_hasNucleus; /**< Indicates if the particle is or contains a nucleus. 0 = no, -1 = yes and 1 = is nucleus. */
860 PQ_suite m_mass; /**< A suite storing the mass physical quantities for the particle. */
861 PQ_suite m_spin; /**< A suite storing the spin physical quantities for the particle. */
862 PQ_suite m_parity; /**< A suite storing the parity physical quantities for the particle. */
863 PQ_suite m_charge; /**< A suite storing the charge physical quantities for the particle. */
864 PQ_suite m_halflife; /**< A suite storing the halflife physical quantities for the particle. */
865 DecayData m_decayData; /**< Stores the decay data for the particle. */
866
867 void setHasNucleus( bool a_hasNucleus ) { m_hasNucleus = a_hasNucleus; }
868
869 public:
870 Particle( HAPI::Node const &a_node, Particle_class a_class, std::string const &a_family, int a_hasNucleus = 0 );
871 virtual ~Particle( );
872
873 std::string const &baseId( void ) const { return( m_baseId ); } /**< Returns a *const* reference to the *m_baseId* member. */
874 std::string const &family( void ) const { return( m_family ); } /**< Returns a *const* reference to the *m_family* member. */
875 bool isAnti( ) const { return( m_anti == IDs::anti ); } /**< Returns the value of the *m_anti* member. */
876 int hasNucleus( void ) const { return( m_hasNucleus ); } /**< Returns the value of the *m_hasNucleus* member. */
877
878 virtual PQ_suite const &mass( void ) const { return( m_mass ); } /**< Returns a *const* reference to the *m_mass* member. */
879 virtual double massValue( char const *a_unit ) const ;
880 double massValue( std::string const &a_unit ) const { return( massValue( a_unit.c_str( ) ) ); }
881 /**< Returns the value of massValue( a_unit.c_str( ) ). */
882
883 PQ_suite const &spin( ) const { return( m_spin ); } /**< Returns a *const* reference to the *m_spin* member. */
884 PQ_suite const &parity( ) const { return( m_parity ); } /**< Returns a *const* reference to the *m_parity* member. */
885 PQ_suite const &charge( ) const { return( m_charge ); } /**< Returns a *const* reference to the *m_charge* member. */
886 PQ_suite const &halflife( ) const { return( m_halflife ); } /**< Returns a *const* reference to the *m_halflife* member. */
887 DecayData const &decayData( ) const { return( m_decayData ); } /**< Returns a *const* reference to the *m_decayData* member. */
888
889 void toXMLList( std::vector<std::string> &a_XMLList, std::string const &a_indent1 ) const ;
890 virtual std::string toXMLListExtraAttributes( void ) const ;
891 virtual void toXMLListExtraElements( std::vector<std::string> &a_XMLList, std::string const &a_indent1 ) const ;
892
893 friend class Unorthodox;
894};
895
896/*
897============================================================
898======================== GaugeBoson ========================
899============================================================
900*/
901
902class GaugeBoson : public Particle {
903
904 public:
905 GaugeBoson( HAPI::Node const &a_node, Database *a_DB, Database *a_parent );
906 virtual ~GaugeBoson( );
907};
908
909/*
910============================================================
911========================== Lepton ==========================
912============================================================
913*/
914
915class Lepton : public Particle {
916
917 private:
918 std::string m_generation; /**< The generation of the lepton (i.e., electronic, muonic or tauonic). */
919
920 public:
921 Lepton( HAPI::Node const &a_node, Database *a_DB, Database *a_parent );
922 virtual ~Lepton( );
923
924 std::string const &generation( void ) const { return( m_generation ); }
925 virtual std::string toXMLListExtraAttributes( void ) const ;
926};
927
928/*
929============================================================
930========================== Baryon ==========================
931============================================================
932*/
933
934class Baryon : public Particle {
935
936 public:
937 Baryon( HAPI::Node const &a_node, Database *a_DB, Database *a_parent );
938 virtual ~Baryon( );
939};
940
941/*
942============================================================
943======================== Unorthodox ========================
944============================================================
945*/
946
947class Unorthodox : public Particle {
948
949 public:
950 Unorthodox( HAPI::Node const &a_node, Database *a_DB, Database *a_parent );
951 virtual ~Unorthodox( );
952};
953
954/*
955============================================================
956========================== Nucleus =========================
957============================================================
958*/
959
960class Nucleus : public Particle {
961
962 private:
963 Nuclide *m_nuclide; /**< The parent nuclide of *this*. */
964 int m_Z; /**< The atomic number of the parent nuclide. */
965 int m_A; /**< The atomic mass number of the parent nuclide. */
966 std::string m_levelName; /**< The string representationn of *m_levelIndex*. */
967 int m_levelIndex; /**< The index of the excited nucleus state. */
968 PQ_suite m_energy; /**< A suite storing the physical quantities representing the nucleus excited energy for the particle.. */
969
970 public:
971 Nucleus( HAPI::Node const &node, Database *a_DB, Nuclide *a_parent );
972 virtual ~Nucleus( );
973
974 Nuclide const *nuclide( ) const { return( m_nuclide ); } /**< Returns a *const* reference to the *m_nuclide* member. */
975 int Z( void ) const { return( m_Z ); } /**< Returns a *const* reference to the *m_Z* member. */
976 int A( void ) const { return( m_A ); } /**< Returns a *const* reference to the *m_A* member. */
977 std::string const &levelName( ) const { return( m_levelName ); } /**< Returns a *const* reference to the *m_levelName* member of *this*. */
978 int levelIndex( void ) const { return( m_levelIndex ); } /**< Returns a *const* reference to the *m_levelIndex* member. */
979 std::string const &atomsID( void ) const ;
980
981 double massValue( char const *a_unit ) const ;
982 PQ_suite const &energy( void ) const { return( m_energy ); } /**< Returns a *const* reference to the *m_energy* member. */
983 double energy( std::string const &a_unit ) const ;
984 virtual std::string toXMLListExtraAttributes( void ) const ;
985 virtual void toXMLListExtraElements( std::vector<std::string> &a_XMLList, std::string const &a_indent1 ) const ;
986};
987
988/*
989============================================================
990========================== Nuclide =========================
991============================================================
992*/
993
994class Nuclide : public Particle {
995
996 private:
997 Isotope *m_isotope; /**< A pointer to the parent isotope. */
998 Nucleus m_nucleus; /**< The nucleus for *this* nuclide. */
999 GammaDecayData m_gammaDecayData; /**< . */
1000
1001 public:
1002 Nuclide( HAPI::Node const &a_node, Database *a_DB, Isotope *a_parent );
1003 virtual ~Nuclide( );
1004
1005 int Z( void ) const;
1006 int A( void ) const;
1007 std::string const &levelName( void ) const { return( m_nucleus.levelName( ) ); }
1008 /**< Returns the result of calling m_nucleus.levelName( ). */
1009 int levelIndex( void ) const { return( m_nucleus.levelIndex( ) ); } /**< Returns the result of calling m_nucleus.levelIndex( ). */
1010 std::string const &atomsID( ) const ;
1011 std::string const &kind( ) const { return( m_gammaDecayData.kind( ) ); }
1012 GammaDecayData const &gammaDecayData( ) const { return( m_gammaDecayData ); }
1013
1014 Isotope const *isotope( ) const { return( m_isotope ); } /**< Returns a *const* reference to the *m_isotope* member. */
1015 Nucleus const &nucleus( ) const { return( m_nucleus ); } /**< Returns a *const* reference to the *m_nucleus* member. */
1016
1017 PQ_suite const &baseMass( void ) const ;
1018 double massValue( char const *a_unit ) const ;
1019 double levelEnergy( std::string const &a_unit ) const { return( m_nucleus.energy( a_unit ) ); }
1020 /**< Returns the result of calling m_nucleus.energy( a_unit ). */
1021
1022 void calculateNuclideGammaBranchStateInfos( PoPI::Database const &a_pops, NuclideGammaBranchStateInfos &a_nuclideGammaBranchStateInfos,
1023 bool a_alwaysAdd = false ) const ;
1024 virtual void toXMLListExtraElements( std::vector<std::string> &a_XMLList, std::string const &a_indent1 ) const ;
1025};
1026
1027/*
1028============================================================
1029========================= Isotope ==========================
1030============================================================
1031*/
1032
1033class Isotope : public SymbolBase {
1034
1035 private:
1036 ChemicalElement *m_chemicalElement; /**< A pointer to the parent chemical element. */
1037 int m_Z; /**< A atomic number for the isotope. */
1038 int m_A; /**< The atomic mass number for the isotope. */
1039 Suite<Nuclide, Isotope> m_nuclides; /**< The suite of nuclides for this isotope. */
1040
1041 public:
1042 Isotope( HAPI::Node const &a_node, Database *a_DB, ChemicalElement *a_parent );
1043 virtual ~Isotope( );
1044
1045 ChemicalElement const *chemicalElement( ) const { return( m_chemicalElement ); } /**< Returns a *const* reference to the *m_isotope* member. */
1046 int Z( void ) const { return( m_Z ); } /**< Returns the value of the *m_Z* member. */
1047 int A( void ) const { return( m_A ); } /**< Returns the value of the *m_A* member. */
1048 Suite<Nuclide, Isotope> const &nuclides( ) const { return( m_nuclides ); } /**< Returns a *const* reference to the *m_nuclides* member. */
1049
1050 void calculateNuclideGammaBranchStateInfos( PoPI::Database const &a_pops, NuclideGammaBranchStateInfos &a_nuclideGammaBranchStateInfos ) const ;
1051 void toXMLList( std::vector<std::string> &a_XMLList, std::string const &a_indent1 ) const ;
1052};
1053
1054/*
1055============================================================
1056===================== ChemicalElement ======================
1057============================================================
1058*/
1059
1061
1062 private:
1063 int m_Z; /**< A atomic number for all isotopes in *thie* chemical element. */
1064 std::string m_name; /**< The name of the chemical element. */
1065 Suite<Isotope, ChemicalElement> m_isotopes; /**< The suite of isotopes for this chemical element. */
1066
1067 public:
1068 ChemicalElement( HAPI::Node const &a_node, Database *a_DB, Database *a_parent );
1069 virtual ~ChemicalElement( );
1070
1071 int Z( void ) const { return( m_Z ); } /**< Returns the value of the *m_Z* member. */
1072 std::string const &name( void ) const { return( m_name ); } /**< Returns the value of the *m_name* member. */
1073
1074 Suite<Isotope, ChemicalElement> const &isotopes( ) const { return( m_isotopes ); } /**< Returns a *const* reference to the *m_isotopes* member. */
1075
1076 void calculateNuclideGammaBranchStateInfos( PoPI::Database const &a_pops, NuclideGammaBranchStateInfos &a_nuclideGammaBranchStateInfos ) const ;
1077 void toXMLList( std::vector<std::string> &a_XMLList, std::string const &a_indent1 ) const ;
1078};
1079
1080/*
1081============================================================
1082=========================== Alias ========================== // FIXME, there should be an alias base class that Alias and MetaStable inherit from.
1083============================================================
1084*/
1085
1086class Alias : public IDBase {
1087
1088 private:
1089 std::string m_pid; /**< The id of the particle *this* is an alias for. */
1090 std::size_t m_pidIndex; /**< The index of the particle with id *m_pid*. */
1091
1092 public:
1093 Alias( HAPI::Node const &a_node, Database *a_DB, Particle_class a_class = Particle_class::alias );
1094 virtual ~Alias( );
1095
1096 std::string const &pid( void ) const { return( m_pid ); } /**< Returns a *const* reference to the *m_pid* member of *this*. */
1097 std::size_t pidIndex( void ) const { return( m_pidIndex ); } /**< Returns a *const* reference to the *m_pidIndex* member of *this*. */
1098 void setPidIndex( std::size_t a_index ) { m_pidIndex = a_index; } /**< Set the member *m_pidIndex* to *a_index*. */
1099
1100 void toXMLList( std::vector<std::string> &a_XMLList, std::string const &a_indent1 ) const ;
1101};
1102
1103/*
1104============================================================
1105======================== MetaStable ========================
1106============================================================
1107*/
1108
1109class MetaStable : public Alias {
1110
1111 private:
1112 int m_metaStableIndex; /**< The meta-stable index for *this*. */
1113
1114 public:
1115 MetaStable( HAPI::Node const &a_node, Database *a_DB );
1116 virtual ~MetaStable( );
1117
1118 int metaStableIndex( void ) const { return( m_metaStableIndex ); } /**< Returns the value of the *m_metaStableIndex* member. */
1119 void toXMLList( std::vector<std::string> &a_XMLList, std::string const &a_indent1 ) const ;
1120};
1121
1122/*
1123============================================================
1124========================= Database =========================
1125============================================================
1126*/
1127
1129
1130 private:
1131 LUPI::FormatVersion m_formatVersion; /**< The **GNDS** **format** attribute of the first file read. */
1132 std::string m_name; /**< The **GNDS** **name** of the first file read in. */
1133 std::string m_version; /**< The **GNDS** **version** of the first file read in. */
1134 ParticleList m_list; /**< The internal list of the particles. */
1135 std::map<std::string, std::size_t> m_idsMap; // Be careful with this as a map[key] will add key if it is not in the map.
1136 std::map<int, std::size_t> m_intidsMap; // Be careful with this as a map[key] will add key if it is not in the map.
1137 /**< This maps each particle id to a unique index. */
1138 SymbolList m_symbolList; /**< The internal list of the symbols. */
1139 std::map<std::string, std::size_t> m_symbolMap; // Be careful with this as a map[key] will add key if it is not in the map.
1140 /**< This maps each symbol to a unique index. */
1141
1142 std::vector<Alias *> m_unresolvedAliases; /**< This is used internally to store aliases when a **PoPs** node is being parsed as the aliases onde is parsed before the particles are parsed. */
1143 std::vector<Alias *> m_aliases; /**< Represents the **PoPs** *aliases* node which contains a list of the **PoPs** **alias** and **metaStable** nodes. */
1144
1145 Suite<GaugeBoson, Database> m_gaugeBosons; /**< Represents the **PoPs** **gaugeBosons** node which contains a list of the **PoPs** **gaugeBoson** nodes. */
1146 Suite<Lepton, Database> m_leptons; /**< Represents the **PoPs** **leptons** node which contains a list of the **PoPs** **lepton** nodes. */
1147 Suite<Baryon, Database> m_baryons; /**< Represents the **PoPs** **baryons** node which contains a list of the **PoPs** **baryon** nodes. */
1148 Suite<ChemicalElement, Database> m_chemicalElements; /**< Represents the **PoPs** **chemicalElements** node which contains a list of the **PoPs** **chemicalElement** nodes. */
1149 Suite<Unorthodox, Database> m_unorthodoxes; /**< Represents the **PoPs** **unorthodoxes** node which contains a list of the **PoPs** **unorthodox** nodes. */
1150
1151 public:
1152 Database( );
1153 Database( std::string const &a_fileName );
1154 Database( HAPI::Node const &a_database );
1155 ~Database( );
1156
1157 LUPI::FormatVersion const &formatVersion( void ) const { return( m_formatVersion ); } /**< Returns a *const* *reference* to the *m_formatVersion* variable of *this*. */
1158 std::string const &name( void ) const { return( m_name ); } /**< Returns a *const* *reference* to the *m_name* variable of *this*. */
1159 std::string const &version( void ) const { return( m_version ); } /**< Returns a *const* *reference* to the *m_version* variable of *this*. */
1160
1161 std::vector<Alias *> unresolvedAliases( ) { return( m_unresolvedAliases ); } /**< Returns a *reference* to the *m_unresolvedAliases* member of *this*. */
1162 std::size_t numberOfUnresolvedAliases( ) { return( m_unresolvedAliases.size( ) ); } /**< Returns the number of unresolved aliases. */
1163 std::vector<std::string> unresolvedAliasIds( ) const ;
1164 std::vector<Alias *> &aliases( ) { return( m_aliases ); } /**< Returns a *const* *reference* to the *m_aliases* variable of *this*. */
1165
1166 void addFile( char const *a_fileName, bool a_warnIfDuplicate );
1167 void addFile( std::string const &a_fileName, bool a_warnIfDuplicate );
1168 void addDatabase( std::string const &a_string, bool a_warnIfDuplicate );
1169 void addDatabase( HAPI::Node const &a_database, bool a_warnIfDuplicate );
1170 void addAlias( Alias *a_alias ) { m_aliases.push_back( a_alias ); } /**< Added the **Alias** *a_alias* to *this*. */
1171
1172 std::string::size_type size( void ) const { return( m_list.size( ) ); } /**< Returns the number of particle in *this*. */
1173 ParticleList const &list( ) { return( m_list ); } /**< Returns a *const* *reference* to the *m_list* member. */
1174 SymbolList const &symbolList( ) { return( m_symbolList ); } /**< Returns a *const* *reference* to the *m_symbolList* member. */
1175 std::size_t operator[]( std::string const &a_id ) const ;
1176 template<typename T> T const &get( std::string const &a_id ) const ;
1177 template<typename T> T const &get( std::size_t a_index ) const ;
1178 Particle const &particle( std::string const &a_id ) const { return( get<Particle>( a_id ) ); } /**< Returns a *const* *reference* to the particle with id *a_id*. */
1179 Particle const &particle( std::size_t a_index ) const { return( get<Particle>( a_index ) ); } /**< Returns a *const* *reference* to the particle with index *a_index*. */
1180 IDBase const &idBase( std::string const &a_id ) const { return( get<IDBase>( a_id ) ); } /**< Returns a *const* *reference* to a **IDBase** instance with id *a_id*. */
1181 IDBase const &idBase( std::size_t &a_index ) const { return( get<IDBase>( a_index ) ); } /**< Returns a *const* *reference* to a **IDBase** instance with id *a_index*. */
1182 ParticleList const &particleList( ) const { return( m_list ); } /**< Returns a *const* *reference* to the *m_list* variable of *this*. */
1183 SymbolList symbolList( ) const { return( m_symbolList ); } /**< Returns a *const* *reference* to the *m_symbolList* variable of *this*. */
1184
1185 bool exists( std::string const &a_id ) const ;
1186 bool exists( std::size_t a_index ) const ;
1187 bool existsIntid( int a_intid ) const ;
1188
1189 Suite<ChemicalElement, Database> const &chemicalElements( ) const { return( m_chemicalElements ); }
1190 /**< Returns a *const* *reference* to the *m_chemicalElements* variable of *this*. */
1191
1192 bool isParticle( std::string const &a_id ) const { return( get<Base>( a_id ).isParticle( ) ); } /**< Returns **true** if *a_id* is a particle and **false** otherwise. */
1193 bool isParticle( std::size_t a_index ) const { return( m_list[a_index]->isParticle( ) ); } /**< Returns **true** if *a_index* is a particle and **false** otherwise. */
1194 bool isAlias( std::string const &a_id ) const { return( get<Base>( a_id ).isAlias( ) ); } /**< Returns **true** if *a_id* is an alias and **false** otherwise. */
1195 bool isAlias( std::size_t a_index ) const { return( m_list[a_index]->isAlias( ) ); } /**< Returns **true** if *a_index* is an alias and **false** otherwise. */
1196 bool isMetaStableAlias( std::string const &a_id ) const { return( get<Base>( a_id ).isMetaStableAlias( ) ); }
1197 /**< Returns **true** if *a_id* is a meta-stable and **false** otherwise. */
1198 bool isMetaStableAlias( std::size_t a_index ) const { return( m_list[a_index]->isMetaStableAlias( ) ); }
1199 /**< Returns **true** if *a_index* is a meta-stable and **false** otherwise. */
1200 std::vector<std::string> aliasReferences( std::string const &a_id );
1201
1202 std::string final( std::string const &a_id, bool a_returnAtMetaStableAlias = false ) const ;
1203 std::size_t final( std::size_t a_index, bool a_returnAtMetaStableAlias = false ) const ;
1204
1205 std::string chemicalElementSymbol( std::string const &a_id ) const ;
1206 std::string isotopeSymbol( std::string const &a_id ) const ;
1207 int intid( std::string const &a_id ) const ;
1208 int intid( std::size_t a_index ) const ;
1209 std::size_t indexFromIntid( int a_intid ) const ;
1210
1211 std::size_t add( Base *a_item );
1212 std::size_t addSymbol( SymbolBase *a_item );
1213
1214 void calculateNuclideGammaBranchStateInfos( NuclideGammaBranchStateInfos &a_nuclideGammaBranchStateInfos, Database const *a_pops2,
1215 std::vector<std::string> &a_extraGammaBranchStates ) const ;
1216 void calculateNuclideGammaBranchStateInfos2( NuclideGammaBranchStateInfos &a_nuclideGammaBranchStateInfos ) const ;
1217
1218 double massValue( std::string const &a_id, std::string const &a_unit ) const ;
1219
1220 void saveAs( std::string const &a_fileName ) const ;
1221 void toXMLList( std::vector<std::string> &a_XMLList, std::string const &a_indent1 ) const ;
1222 void print( bool a_printIndices );
1223};
1224
1225/* *********************************************************************************************************//**
1226 * Returns the partile in *this* that has index *a_index*.
1227 *
1228 * @param a_index [in] The index of the particle to return.
1229 *
1230 * @return A *const* reference to the particle at index *a_index*.
1231 ***********************************************************************************************************/
1232
1233template<typename T> T const &Database::get( std::size_t a_index ) const {
1234
1235 Base *particle = m_list[a_index];
1236 if( particle == nullptr ) throw std::range_error( std::string( "particle not in database" ) );
1237 T const *object = dynamic_cast<T const *>( particle );
1238 if( object == nullptr ) throw std::bad_cast( );
1239
1240 return( *object );
1241}
1242
1243/* *********************************************************************************************************//**
1244 * Returns the partile in *this* that has index *a_index*.
1245 *
1246 * @param a_id [in] The **PoPs** id of the particle to return.
1247 *
1248 * @return A *const* reference to the particle with id *a_id*.
1249 ***********************************************************************************************************/
1250
1251template<typename T> T const &Database::get( std::string const &a_id ) const {
1252
1253 auto index = (*this)[a_id];
1254 Base *particle = m_list[index];
1255 T const *object = dynamic_cast<T const *>( particle );
1256 if( object == nullptr ) throw std::bad_cast( );
1257
1258 return( *object );
1259}
1260
1261double getPhysicalQuantityAsDouble( PhysicalQuantity const &a_physicalQuantity );
1262double getPhysicalQuantityOfSuiteAsDouble( PQ_suite const &a_suite, bool a_allowEmpty = false, double a_emptyValue = 0.0 );
1263bool supportedFormat( LUPI::FormatVersion const &a_formatVersion );
1264std::string baseAntiQualifierFromID( std::string const &a_id, std::string &a_anti, std::string *a_qualifier = nullptr );
1265
1267std::string chemicalElementInfoFromZ( int a_Z, bool a_wantSymbol, bool a_asNucleus = false );
1268std::string const &chemicalElementSymbolFromZ( int a_Z );
1269int Z_FromChemicalElementSymbol( std::string const &a_symbol );
1270
1271int family2Integer( Particle_class a_family );
1272int intidHelper( bool a_isAnti, Particle_class a_family, int a_SSSSSSS );
1273
1274}
1275
1276#endif // End of PoPI_hpp_included
bool empty() const
Definition HAPI_Node.cc:150
Node first_child() const
Definition HAPI_Node.cc:82
Alias(HAPI::Node const &a_node, Database *a_DB, Particle_class a_class=Particle_class::alias)
Definition PoPI_alias.cc:28
void toXMLList(std::vector< std::string > &a_XMLList, std::string const &a_indent1) const
Definition PoPI_alias.cc:55
void setPidIndex(std::size_t a_index)
Definition PoPI.hpp:1098
std::size_t pidIndex(void) const
Definition PoPI.hpp:1097
std::string const & pid(void) const
Definition PoPI.hpp:1096
virtual ~Alias()
Definition PoPI_alias.cc:44
Baryon(HAPI::Node const &a_node, Database *a_DB, Database *a_parent)
virtual ~Baryon()
bool isAlias(void) const
Definition PoPI.hpp:658
Base(std::string const &a_id, Particle_class a_class)
Definition PoPI_base.cc:23
bool isLepton() const
Definition PoPI.hpp:664
void setIndex(std::size_t a_index)
Definition PoPI.hpp:654
friend MetaStable
Definition PoPI.hpp:673
std::size_t index(void) const
Definition PoPI.hpp:653
bool isUnorthodox() const
Definition PoPI.hpp:666
Particle_class Class(void) const
Definition PoPI.hpp:656
friend Nuclide
Definition PoPI.hpp:679
bool isBaryon() const
Definition PoPI.hpp:665
friend Nucleus
Definition PoPI.hpp:678
friend Alias
Definition PoPI.hpp:674
bool isGaugeBoson() const
Definition PoPI.hpp:663
friend Baryon
Definition PoPI.hpp:675
bool isChemicalElement() const
Definition PoPI.hpp:670
int intid() const
Definition PoPI.hpp:655
bool isNucleus() const
Definition PoPI.hpp:667
bool isNuclide() const
Definition PoPI.hpp:668
friend GaugeBoson
Definition PoPI.hpp:676
bool isIsotope() const
Definition PoPI.hpp:669
friend Unorthodox
Definition PoPI.hpp:680
virtual bool isParticle() const
Definition PoPI.hpp:657
friend Lepton
Definition PoPI.hpp:677
virtual ~Base()
Definition PoPI_base.cc:50
std::string const & ID(void) const
Definition PoPI.hpp:652
bool isMetaStableAlias(void) const
Definition PoPI.hpp:660
std::string const & name(void) const
Definition PoPI.hpp:1072
Suite< Isotope, ChemicalElement > const & isotopes() const
Definition PoPI.hpp:1074
void calculateNuclideGammaBranchStateInfos(PoPI::Database const &a_pops, NuclideGammaBranchStateInfos &a_nuclideGammaBranchStateInfos) const
void toXMLList(std::vector< std::string > &a_XMLList, std::string const &a_indent1) const
ChemicalElement(HAPI::Node const &a_node, Database *a_DB, Database *a_parent)
int Z(void) const
Definition PoPI.hpp:1071
std::vector< std::string > aliasReferences(std::string const &a_id)
std::string isotopeSymbol(std::string const &a_id) const
bool isMetaStableAlias(std::size_t a_index) const
Definition PoPI.hpp:1198
bool existsIntid(int a_intid) const
std::size_t addSymbol(SymbolBase *a_item)
LUPI::FormatVersion const & formatVersion(void) const
Definition PoPI.hpp:1157
void toXMLList(std::vector< std::string > &a_XMLList, std::string const &a_indent1) const
void addFile(char const *a_fileName, bool a_warnIfDuplicate)
std::string const & version(void) const
Definition PoPI.hpp:1159
double massValue(std::string const &a_id, std::string const &a_unit) const
void print(bool a_printIndices)
ParticleList const & particleList() const
Definition PoPI.hpp:1182
std::size_t operator[](std::string const &a_id) const
bool isParticle(std::size_t a_index) const
Definition PoPI.hpp:1193
T const & get(std::string const &a_id) const
Particle const & particle(std::size_t a_index) const
Definition PoPI.hpp:1179
bool exists(std::string const &a_id) const
void addAlias(Alias *a_alias)
Definition PoPI.hpp:1170
bool isParticle(std::string const &a_id) const
Definition PoPI.hpp:1192
void addDatabase(std::string const &a_string, bool a_warnIfDuplicate)
bool isMetaStableAlias(std::string const &a_id) const
Definition PoPI.hpp:1196
void calculateNuclideGammaBranchStateInfos2(NuclideGammaBranchStateInfos &a_nuclideGammaBranchStateInfos) const
Particle const & particle(std::string const &a_id) const
Definition PoPI.hpp:1178
std::size_t numberOfUnresolvedAliases()
Definition PoPI.hpp:1162
void saveAs(std::string const &a_fileName) const
IDBase const & idBase(std::string const &a_id) const
Definition PoPI.hpp:1180
std::string chemicalElementSymbol(std::string const &a_id) const
SymbolList const & symbolList()
Definition PoPI.hpp:1174
Suite< ChemicalElement, Database > const & chemicalElements() const
Definition PoPI.hpp:1189
std::size_t indexFromIntid(int a_intid) const
std::vector< std::string > unresolvedAliasIds() const
SymbolList symbolList() const
Definition PoPI.hpp:1183
std::string::size_type size(void) const
Definition PoPI.hpp:1172
bool isAlias(std::size_t a_index) const
Definition PoPI.hpp:1195
std::string const & name(void) const
Definition PoPI.hpp:1158
IDBase const & idBase(std::size_t &a_index) const
Definition PoPI.hpp:1181
int intid(std::string const &a_id) const
void addDatabase(HAPI::Node const &a_database, bool a_warnIfDuplicate)
ParticleList const & list()
Definition PoPI.hpp:1173
bool isAlias(std::string const &a_id) const
Definition PoPI.hpp:1194
std::vector< Alias * > unresolvedAliases()
Definition PoPI.hpp:1161
std::size_t add(Base *a_item)
void calculateNuclideGammaBranchStateInfos(NuclideGammaBranchStateInfos &a_nuclideGammaBranchStateInfos, Database const *a_pops2, std::vector< std::string > &a_extraGammaBranchStates) const
std::vector< Alias * > & aliases()
Definition PoPI.hpp:1164
T const & get(std::size_t a_index) const
DecayData(HAPI::Node const &a_node)
Suite< DecayMode, DecayData > const & decayModes(void) const
Definition PoPI.hpp:811
void calculateNuclideGammaBranchStateInfo(PoPI::Database const &a_pops, NuclideGammaBranchStateInfo &a_nuclideGammaBranchStateInfo) const
void toXMLList(std::vector< std::string > &a_XMLList, std::string const &a_indent1) const
void calculateNuclideGammaBranchStateInfo(PoPI::Database const &a_pops, NuclideGammaBranchStateInfo &a_nuclideGammaBranchStateInfo) const
DecayMode(HAPI::Node const &a_node, DecayData const *a_decayData)
void toXMLList(std::vector< std::string > &a_XMLList, std::string const &a_indent1) const
std::string const & mode() const
Definition PoPI.hpp:787
PQ_suite const & photonEmissionProbabilities() const
Definition PoPI.hpp:789
std::string const & label() const
Definition PoPI.hpp:786
Suite< Decay, DecayMode > const & decayPath() const
Definition PoPI.hpp:790
PQ_suite const & probability() const
Definition PoPI.hpp:788
std::string const & mode() const
Definition PoPI.hpp:761
Suite< Product, Decay > const & products(void) const
Definition PoPI.hpp:763
void toXMLList(std::vector< std::string > &a_XMLList, std::string const &a_indent1) const
Decay(HAPI::Node const &a_node, DecayMode const *a_decayMode)
bool complete() const
Definition PoPI.hpp:762
int index(void) const
Definition PoPI.hpp:760
Exception(std::string const &a_message)
Definition PoPI_misc.cc:494
int rows() const
Definition PoPI.hpp:838
std::vector< double > const & probabilities() const
Definition PoPI.hpp:841
void calculateNuclideGammaBranchStateInfo(PoPI::Database const &a_pops, NuclideGammaBranchStateInfo &nuclideGammaBranchStateInfo) const
std::vector< std::string > const & ids() const
Definition PoPI.hpp:840
GammaDecayData(HAPI::Node const &a_node)
std::string const & kind() const
Definition PoPI.hpp:837
int colunms() const
Definition PoPI.hpp:839
std::vector< double > const & photonEmissionProbabilities() const
Definition PoPI.hpp:842
GaugeBoson(HAPI::Node const &a_node, Database *a_DB, Database *a_parent)
double massValue2(Database const &a_DB, std::string const &a_unit) const
virtual ~IDBase()
Definition PoPI_base.cc:83
std::size_t addToDatabase(Database *a_DB)
Definition PoPI_base.cc:95
IDBase(std::string const &a_id, Particle_class a_class)
Definition PoPI_base.cc:63
virtual ~Isotope()
void calculateNuclideGammaBranchStateInfos(PoPI::Database const &a_pops, NuclideGammaBranchStateInfos &a_nuclideGammaBranchStateInfos) const
int Z(void) const
Definition PoPI.hpp:1046
Suite< Nuclide, Isotope > const & nuclides() const
Definition PoPI.hpp:1048
ChemicalElement const * chemicalElement() const
Definition PoPI.hpp:1045
Isotope(HAPI::Node const &a_node, Database *a_DB, ChemicalElement *a_parent)
void toXMLList(std::vector< std::string > &a_XMLList, std::string const &a_indent1) const
int A(void) const
Definition PoPI.hpp:1047
std::string const & generation(void) const
Definition PoPI.hpp:924
Lepton(HAPI::Node const &a_node, Database *a_DB, Database *a_parent)
virtual ~Lepton()
virtual std::string toXMLListExtraAttributes(void) const
virtual ~MetaStable()
Definition PoPI_alias.cc:89
MetaStable(HAPI::Node const &a_node, Database *a_DB)
Definition PoPI_alias.cc:72
void toXMLList(std::vector< std::string > &a_XMLList, std::string const &a_indent1) const
int metaStableIndex(void) const
Definition PoPI.hpp:1118
std::string const & atomsID(void) const
int levelIndex(void) const
Definition PoPI.hpp:978
std::string const & levelName() const
Definition PoPI.hpp:977
int A(void) const
Definition PoPI.hpp:976
Nucleus(HAPI::Node const &node, Database *a_DB, Nuclide *a_parent)
int Z(void) const
Definition PoPI.hpp:975
virtual ~Nucleus()
virtual void toXMLListExtraElements(std::vector< std::string > &a_XMLList, std::string const &a_indent1) const
PQ_suite const & energy(void) const
Definition PoPI.hpp:982
double massValue(char const *a_unit) const
virtual std::string toXMLListExtraAttributes(void) const
Nuclide const * nuclide() const
Definition PoPI.hpp:974
double probability() const
Definition PoPI.hpp:560
NuclideGammaBranchInfo(double a_probability, double a_photonEmissionProbability, double a_gammaEnergy, std::string const &a_residualState)
double photonEmissionProbability() const
Definition PoPI.hpp:561
std::string const & residualState() const
Definition PoPI.hpp:564
double gammaEnergy() const
Definition PoPI.hpp:563
void setNuclearLevelEnergyWidth(double a_nuclearLevelEnergyWidth)
Definition PoPI.hpp:595
double averageGammaEnergy() const
Definition PoPI.hpp:599
void calculateDerivedData(NuclideGammaBranchStateInfos &a_nuclideGammaBranchStateInfos)
double nuclearLevelEnergyWidth() const
Definition PoPI.hpp:593
std::string const & state() const
Definition PoPI.hpp:589
std::string const & kind() const
Definition PoPI.hpp:591
double nuclearLevelEnergy() const
Definition PoPI.hpp:592
void add(NuclideGammaBranchInfo const &a_nuclideGammaBranchInfo)
NuclideGammaBranchStateInfo(std::string const &a_state, int a_intid, std::string const &a_kind, double a_nuclearLevelEnergy)
std::vector< NuclideGammaBranchInfo > const & branches() const
Definition PoPI.hpp:600
std::vector< NuclideGammaBranchStateInfo * > & nuclideGammaBranchStateInfos()
Definition PoPI.hpp:625
std::size_t size() const
Definition PoPI.hpp:622
NuclideGammaBranchStateInfo * operator[](std::size_t a_index)
Definition PoPI.hpp:623
NuclideGammaBranchStateInfo const * operator[](std::size_t a_index) const
Definition PoPI.hpp:624
NuclideGammaBranchStateInfo * find(std::string const &a_state)
void add(NuclideGammaBranchStateInfo *a_nuclideGammaBranchStateInfo)
virtual ~Nuclide()
GammaDecayData const & gammaDecayData() const
Definition PoPI.hpp:1012
void calculateNuclideGammaBranchStateInfos(PoPI::Database const &a_pops, NuclideGammaBranchStateInfos &a_nuclideGammaBranchStateInfos, bool a_alwaysAdd=false) const
Isotope const * isotope() const
Definition PoPI.hpp:1014
int levelIndex(void) const
Definition PoPI.hpp:1009
PQ_suite const & baseMass(void) const
virtual void toXMLListExtraElements(std::vector< std::string > &a_XMLList, std::string const &a_indent1) const
std::string const & levelName(void) const
Definition PoPI.hpp:1007
double levelEnergy(std::string const &a_unit) const
Definition PoPI.hpp:1019
std::string const & kind() const
Definition PoPI.hpp:1011
std::string const & atomsID() const
int A(void) const
int Z(void) const
Nucleus const & nucleus() const
Definition PoPI.hpp:1015
double massValue(char const *a_unit) const
Nuclide(HAPI::Node const &a_node, Database *a_DB, Isotope *a_parent)
PQ_double(HAPI::Node const &a_node)
double value(void) const
Definition PoPI.hpp:443
virtual std::string valueToString(void) const
double value(std::string const &a_unit) const
Definition PoPI.hpp:445
double value(char const *a_unit) const
std::string value(char const *a_unit) const
PQ_fraction(HAPI::Node const &a_node)
virtual std::string valueToString(void) const
std::string value(void) const
std::string value(std::string const &a_unit) const
Definition PoPI.hpp:485
PQ_integer(HAPI::Node const &a_node)
int value(char const *a_unit) const
int value(std::string const &a_unit) const
Definition PoPI.hpp:467
virtual std::string valueToString(void) const
int value(void) const
Definition PoPI.hpp:465
PQ_shell(HAPI::Node const &a_node)
std::string value(std::string const &a_unit) const
Definition PoPI.hpp:504
std::string value(void) const
Definition PoPI.hpp:502
virtual std::string valueToString(void) const
std::string value(char const *a_unit) const
PQ_string(HAPI::Node const &a_node)
std::string & label(void)
Definition PoPI.hpp:537
PQ_suite(HAPI::Node const &a_node)
void toXMLList(std::vector< std::string > &a_XMLList, std::string const &a_indent1) const
std::string const & qualifier()
Definition PoPI.hpp:293
bool isNucleus()
Definition PoPI.hpp:286
void print(bool a_terse, std::string const &a_indent="") const
std::string const & symbol()
Definition PoPI.hpp:289
std::string const & Id()
Definition PoPI.hpp:284
ParseIdInfo(std::string const &a_id)
bool isSupported()
Definition PoPI.hpp:283
bool isChemicalElement()
Definition PoPI.hpp:287
bool isNuclear()
Definition PoPI.hpp:285
ParseIntidInfo(int a_intid, bool a_GRIN_mode=false)
Definition PoPI_intId.cc:76
bool isNuclearMetaStable()
Definition PoPI.hpp:241
Particle_class family()
Definition PoPI.hpp:234
friend class Unorthodox
Definition PoPI.hpp:893
Particle(HAPI::Node const &a_node, Particle_class a_class, std::string const &a_family, int a_hasNucleus=0)
virtual std::string toXMLListExtraAttributes(void) const
virtual ~Particle()
int hasNucleus(void) const
Definition PoPI.hpp:876
virtual double massValue(char const *a_unit) const
virtual PQ_suite const & mass(void) const
Definition PoPI.hpp:878
virtual void toXMLListExtraElements(std::vector< std::string > &a_XMLList, std::string const &a_indent1) const
bool isAnti() const
Definition PoPI.hpp:875
PQ_suite const & charge() const
Definition PoPI.hpp:885
void toXMLList(std::vector< std::string > &a_XMLList, std::string const &a_indent1) const
std::string const & family(void) const
Definition PoPI.hpp:874
double massValue(std::string const &a_unit) const
Definition PoPI.hpp:880
PQ_suite const & halflife() const
Definition PoPI.hpp:886
PQ_suite const & spin() const
Definition PoPI.hpp:883
DecayData const & decayData() const
Definition PoPI.hpp:887
PQ_suite const & parity() const
Definition PoPI.hpp:884
std::string const & baseId(void) const
Definition PoPI.hpp:873
void toXMLList(std::vector< std::string > &a_XMLList, std::string const &a_indent1) const
std::string const & valueString(void) const
Definition PoPI.hpp:419
virtual std::string valueToString(void) const =0
PhysicalQuantity(HAPI::Node const &a_node, PQ_class a_class)
std::string const & label(void) const
Definition PoPI.hpp:418
std::string const & tag(void) const
Definition PoPI.hpp:417
PQ_class Class(void) const
Definition PoPI.hpp:416
std::string const & unit(void) const
Definition PoPI.hpp:420
void toXMLList(std::vector< std::string > &a_XMLList, std::string const &a_indent1) const
std::string const & label() const
Definition PoPI.hpp:737
std::string const & pid() const
Definition PoPI.hpp:736
int ID() const
Definition PoPI.hpp:735
Product(HAPI::Node const &a_node, Decay *a_DB)
std::string::size_type size(void) const
Definition PoPI.hpp:321
T & operator[](std::size_t a_index) const
Definition PoPI.hpp:322
Suite(std::string const &a_moniker)
Definition PoPI.hpp:316
void appendFromParentNode2(HAPI::Node const &a_node, T2 *a_parent)
Definition PoPI.hpp:367
void toXMLList(std::vector< std::string > &a_XMLList, std::string const &a_indent1) const
Definition PoPI.hpp:383
std::string const & moniker(void)
Definition PoPI.hpp:323
void appendFromParentNode(HAPI::Node const &a_node, Database *a_DB, T2 *a_parent)
Definition PoPI.hpp:351
SymbolBase(HAPI::Node const &a_node, Particle_class a_class)
Definition PoPI_base.cc:110
std::string const & symbol() const
Definition PoPI.hpp:712
bool isParticle() const
Definition PoPI.hpp:715
std::size_t addToSymbols(Database *a_DB)
Definition PoPI_base.cc:130
Unorthodox(HAPI::Node const &a_node, Database *a_DB, Database *a_parent)
Definition PoPI.hpp:28
double getPhysicalQuantityOfSuiteAsDouble(PQ_suite const &a_suite, bool a_allowEmpty=false, double a_emptyValue=0.0)
Definition PoPI_misc.cc:435
int particleA(Base const &a_particle, bool a_isNeutronProtonANucleon=false)
Definition PoPI_misc.cc:227
std::vector< Base * > ParticleList
Definition PoPI.hpp:186
std::string specialParticleID(SpecialParticleID_mode a_mode, std::string const &a_id)
Definition PoPI_misc.cc:79
int particleMetaStableIndex(Base const &a_particle)
Definition PoPI_misc.cc:349
PQ_class
Definition PoPI.hpp:105
int particleZ(Base const &a_particle, bool a_isNeutronProtonANucleon=false)
Definition PoPI_misc.cc:150
int maximumChemicalElementZ()
double getPhysicalQuantityAsDouble(PhysicalQuantity const &a_physicalQuantity)
Definition PoPI_misc.cc:403
bool supportedFormat(LUPI::FormatVersion const &a_formatVersion)
Definition PoPI_misc.cc:39
std::map< std::string, std::string > supportedNucleusAliases
std::string const & chemicalElementSymbolFromZ(int a_Z)
int intidHelper(bool a_isAnti, Particle_class a_family, int a_SSSSSSS)
Definition PoPI_intId.cc:47
int particleZA(Base const &a_particle, bool a_isNeutronProtonANucleon=false)
Definition PoPI_misc.cc:294
std::vector< SymbolBase * > SymbolList
Definition PoPI.hpp:187
SpecialParticleID_mode
Definition PoPI.hpp:116
std::string chemicalElementInfoFromZ(int a_Z, bool a_wantSymbol, bool a_asNucleus=false)
bool compareSpecialParticleIDs(std::string const &a_id1, std::string const &a_id2)
Definition PoPI_misc.cc:133
void appendXMLEnd(std::vector< std::string > &a_XMLList, std::string const &a_label)
Definition PoPI_misc.cc:53
int family2Integer(Particle_class a_family)
Definition PoPI_intId.cc:21
std::string baseAntiQualifierFromID(std::string const &a_id, std::string &a_anti, std::string *a_qualifier=nullptr)
Definition PoPI_misc.cc:458
Particle_class
Definition PoPI.hpp:58
int Z_FromChemicalElementSymbol(std::string const &a_symbol)
static std::string const anti
Definition PoPI.hpp:173
static std::string const photon
Definition PoPI.hpp:162
static std::string const familiarPhoton
Definition PoPI.hpp:166
static std::string const neutron
Definition PoPI.hpp:164
static std::string const FissionProductENDL99125
Definition PoPI.hpp:172
static std::string const familiarTriton
Definition PoPI.hpp:168
static std::string const FissionProductENDL99120
Definition PoPI.hpp:171
static std::string const familiarHelion
Definition PoPI.hpp:169
static std::string const proton
Definition PoPI.hpp:165
static std::string const familiarAlpha
Definition PoPI.hpp:170
static std::string const familiarDeuteron
Definition PoPI.hpp:167
static std::string const electron
Definition PoPI.hpp:163
static int constexpr FissionProductENDL99120
Definition PoPI.hpp:180
static int constexpr photon
Definition PoPI.hpp:178
static int constexpr neutron
Definition PoPI.hpp:177
static int constexpr electron
Definition PoPI.hpp:179
static int constexpr FissionProductENDL99125
Definition PoPI.hpp:181