Geant4  10.01.p01
G4CascadeCoalescence.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 // G4CascadeCoalescence: Factory model for final-state interactions to
27 // produce light ions from cascade nucleons. The algorithm implemented
28 // here is descirbed in Section 2.3 of the LAQGSM documentation (p. 11-12)
29 // [http://lib-www.lanl.gov/la-pubs/00818645.pdf].
30 //
31 // 20110917 Michael Kelsey
32 // 20110920 M. Kelsey -- Use environment variables to set momentum cuts for tuning,
33 // replace polymorphic argument lists with use of "ClusterCandidate"
34 // 20140116 M. Kelsey -- Move statics to const data members to avoid weird
35 // interactions with MT.
36 
37 #ifndef G4CASCADE_COALESCENCE_HH
38 #define G4CASCADE_COALESCENCE_HH
39 
40 #include "globals.hh"
41 #include "G4InuclNuclei.hh"
42 #include "G4LorentzVector.hh"
43 #include <vector>
44 #include <set>
45 
46 class G4CollisionOutput;
48 
49 
51 public:
52  G4CascadeCoalescence(G4int verbose=0);
53  virtual ~G4CascadeCoalescence();
54 
55  // Final state particle list is modified directly
56  void FindClusters(G4CollisionOutput& finalState);
57 
58  void setVerboseLevel(G4int verbose) { verboseLevel = verbose; }
59 
60 private:
61  typedef std::vector<size_t> ClusterCandidate; // Indices of constituents
62 
63  G4int verboseLevel; // Control diagnostic messages
64 
65  std::vector<ClusterCandidate> allClusters; // List of candidates found
66  std::set<size_t> triedClusters; // Hashes of combinatorics
67  std::set<size_t> usedNucleons; // List of converted nucleons
68 
69  G4CollisionOutput* thisFinalState; // Pointers to current event
70  const std::vector<G4InuclElementaryParticle>* thisHadrons;
71 
72  ClusterCandidate thisCluster; // Reusable buffer for attempts
73  G4InuclNuclei thisLightIon; // Reusable construction buffer
74 
75  const G4double dpMaxDoublet; // Relative momenta for clusters
78 
79  // Processing stages -- search, construct, cleanup
80  void selectCandidates();
81  void createNuclei();
82  void removeNucleons();
83 
84  // Do combinatorics of given nucleons to make candidates
85  void tryClusters(size_t idx1, size_t idx2);
86  void tryClusters(size_t idx1, size_t idx2, size_t idx3);
87  void tryClusters(size_t idx1, size_t idx2, size_t idx3, size_t idx4);
88 
89  // Create cluster candidate with listed indices
90  void fillCluster(size_t idx1, size_t idx2);
91  void fillCluster(size_t idx1, size_t idx2, size_t idx3);
92  void fillCluster(size_t idx1, size_t idx2, size_t idx3, size_t idx4);
93 
94  // Convert cluster to hash index (for combinatoric reduction)
95  size_t clusterHash(const ClusterCandidate& clus) const;
96 
97  // Check if candidate cluster has already been evaluated
98  bool clusterTried(const ClusterCandidate& clus) const {
99  return triedClusters.find(clusterHash(clus)) != triedClusters.end();
100  }
101 
102  // Check if indexed nucleon is already in a cluster
103  bool nucleonUsed(size_t idx) const {
104  return usedNucleons.find(idx) != usedNucleons.end();
105  }
106 
107  // Evaluate conditions for cluster to form light ion
108  bool allNucleons(const ClusterCandidate& clus) const;
109  bool goodCluster(const ClusterCandidate& clus) const;
110  G4int clusterType(const ClusterCandidate& aCluster) const;
111 
112  // Extract hadron from final state list
113  const G4InuclElementaryParticle& getHadron(size_t idx) const {
114  return (*thisHadrons)[idx];
115  }
116 
117  // Convert candidate nucleon set into output nucleus (true == success)
118  bool makeLightIon(const ClusterCandidate& aCluster);
119 
120  // Kinematics for cluster evaluations
121  G4LorentzVector getClusterMomentum(const ClusterCandidate& aCluster) const;
122  mutable G4LorentzVector pCluster; // Reusable buffer to reduce churn
123 
124  G4double maxDeltaP(const ClusterCandidate& aCluster) const;
125 
126  // Report cluster arguments for validation
127  void reportArgs(const G4String& name, const ClusterCandidate& clus) const;
128  void reportResult(const G4String& name, const G4InuclNuclei& nucl) const;
129 };
130 
131 #endif /* G4CASCADE_COALESCENCE_HH */
bool allNucleons(const ClusterCandidate &clus) const
void setVerboseLevel(G4int verbose)
size_t clusterHash(const ClusterCandidate &clus) const
G4String name
Definition: TRTMaterials.hh:40
std::vector< size_t > ClusterCandidate
std::set< size_t > triedClusters
G4double maxDeltaP(const ClusterCandidate &aCluster) const
bool makeLightIon(const ClusterCandidate &aCluster)
std::set< size_t > usedNucleons
int G4int
Definition: G4Types.hh:78
void reportArgs(const G4String &name, const ClusterCandidate &clus) const
bool goodCluster(const ClusterCandidate &clus) const
std::vector< ClusterCandidate > allClusters
const std::vector< G4InuclElementaryParticle > * thisHadrons
G4CollisionOutput * thisFinalState
const G4InuclElementaryParticle & getHadron(size_t idx) const
bool nucleonUsed(size_t idx) const
G4int clusterType(const ClusterCandidate &aCluster) const
void reportResult(const G4String &name, const G4InuclNuclei &nucl) const
bool clusterTried(const ClusterCandidate &clus) const
G4CascadeCoalescence(G4int verbose=0)
void tryClusters(size_t idx1, size_t idx2)
ClusterCandidate thisCluster
void FindClusters(G4CollisionOutput &finalState)
double G4double
Definition: G4Types.hh:76
G4LorentzVector getClusterMomentum(const ClusterCandidate &aCluster) const
CLHEP::HepLorentzVector G4LorentzVector
void fillCluster(size_t idx1, size_t idx2)