Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4INCLClusterDecay.hh
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// INCL++ intra-nuclear cascade model
27// Pekka Kaitaniemi, CEA and Helsinki Institute of Physics
28// Davide Mancusi, CEA
29// Alain Boudard, CEA
30// Sylvie Leray, CEA
31// Joseph Cugnon, University of Liege
32//
33// INCL++ revision: v5.1.8
34//
35#define INCLXX_IN_GEANT4_MODE 1
36
37#include "globals.hh"
38
39/** \file G4INCLClusterDecay.hh
40 * \brief Static class for carrying out cluster decays
41 *
42 * \date 6th July 2011
43 * \author Davide Mancusi
44 */
45
46#ifndef G4INCLCLUSTERDECAY_HH
47#define G4INCLCLUSTERDECAY_HH
48
49#include "G4INCLCluster.hh"
51
52namespace G4INCL {
53
54 /**
55 * Pauli blocking
56 *
57 */
59 public:
60 /// \brief True if the cluster is stable.
61 static G4bool isStable(Cluster const * const c) {
62 const G4int Z = c->getZ();
63 const G4int A = c->getA();
65 }
66
67 /** \brief Carries out a cluster decay
68 *
69 * \param c cluster that should decay
70 * \return list of decay products
71 */
72 static ParticleList decay(Cluster * const c);
73
74 protected:
77
78 private:
79 /** \brief Recursively decay clusters
80 *
81 * \param c cluster that should decay
82 * \param decayProducts decay products are appended to the end of this list
83 */
84 static void recursiveDecay(Cluster * const c, ParticleList *decayProducts);
85
86 /// \brief Carries out two-body decays
87 static void twoBodyDecay(Cluster * const c, ParticleTable::ClusterDecayType theDecayMode, ParticleList *decayProducts);
88
89 /// \brief Carries out three-body decays
90 static void threeBodyDecay(Cluster * const c, ParticleTable::ClusterDecayType theDecayMode, ParticleList *decayProducts);
91
92 /** \brief Disassembles unbound nuclei using a simple phase-space model
93 *
94 * The decay products are assumed to uniformly populate the momentum space
95 * (the "phase-space" naming is a long-standing but misleading convention).
96 * The generation of the final state is done by rejection, using the
97 * Raubold-Lynch method. Parts of our implementation were "inspired" by
98 * ROOT's TGenPhaseSpace class, which in turn is a translation of CERNLIB's
99 * historical GENBOD routine [CERN report 68-15 (1968)]. The ROOT
100 * implementation is documented at the following URL:
101 *
102 * http://root.cern.ch/root/html/TGenPhaseSpace.html#TGenPhaseSpace
103 */
104 static void phaseSpaceDecay(Cluster * const c, ParticleTable::ClusterDecayType theDecayMode, ParticleList *decayProducts);
105
106 };
107
108}
109
110#endif // G4INCLCLUSTERDECAY
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
static ParticleList decay(Cluster *const c)
Carries out a cluster decay.
static G4bool isStable(Cluster const *const c)
True if the cluster is stable.
static const ClusterDecayType clusterDecayMode[clusterTableZSize][clusterTableASize]
G4int getZ() const
Returns the charge number.
G4int getA() const
Returns the baryon number.
std::list< G4INCL::Particle * > ParticleList