Geant4  10.03.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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 // 20151016 M. Kelsey -- Replace forward declare of G4InuclElemPart w/include.
37 // 20170406 M. Kelsey -- Remove clusterHash and triedClusters registry.
38 
39 #ifndef G4CASCADE_COALESCENCE_HH
40 #define G4CASCADE_COALESCENCE_HH
41 
42 #include "globals.hh"
44 #include "G4InuclNuclei.hh"
45 #include "G4LorentzVector.hh"
46 #include <vector>
47 #include <set>
48 
49 class G4CollisionOutput;
50 
51 
53 public:
54  G4CascadeCoalescence(G4int verbose=0);
55  virtual ~G4CascadeCoalescence();
56 
57  // Final state particle list is modified directly
58  void FindClusters(G4CollisionOutput& finalState);
59 
60  void setVerboseLevel(G4int verbose) { verboseLevel = verbose; }
61 
62 private:
63  typedef std::vector<size_t> ClusterCandidate; // Indices of constituents
64 
65  G4int verboseLevel; // Control diagnostic messages
66 
67  std::vector<ClusterCandidate> allClusters; // List of candidates found
68  std::set<size_t> usedNucleons; // List of converted nucleons
69 
70  G4CollisionOutput* thisFinalState; // Pointers to current event
71  const std::vector<G4InuclElementaryParticle>* thisHadrons;
72 
73  ClusterCandidate thisCluster; // Reusable buffer for attempts
74  G4InuclNuclei thisLightIon; // Reusable construction buffer
75 
76  const G4double dpMaxDoublet; // Relative momenta for clusters
77  const G4double dpMaxTriplet;
78  const G4double dpMaxAlpha;
79 
80  // Processing stages -- search, construct, cleanup
81  void selectCandidates();
82  void createNuclei();
83  void removeNucleons();
84 
85  // Do combinatorics of given nucleons to make candidates
86  void tryClusters(size_t idx1, size_t idx2);
87  void tryClusters(size_t idx1, size_t idx2, size_t idx3);
88  void tryClusters(size_t idx1, size_t idx2, size_t idx3, size_t idx4);
89 
90  // Create cluster candidate with listed indices
91  void fillCluster(size_t idx1, size_t idx2);
92  void fillCluster(size_t idx1, size_t idx2, size_t idx3);
93  void fillCluster(size_t idx1, size_t idx2, size_t idx3, size_t idx4);
94 
95  // Check if indexed nucleon is already in a cluster
96  bool nucleonUsed(size_t idx) const {
97  return usedNucleons.find(idx) != usedNucleons.end();
98  }
99 
100  // Evaluate conditions for cluster to form light ion
101  bool allNucleons(const ClusterCandidate& clus) const;
102  bool goodCluster(const ClusterCandidate& clus) const;
103  G4int clusterType(const ClusterCandidate& aCluster) const;
104 
105  // Extract hadron from final state list
106  const G4InuclElementaryParticle& getHadron(size_t idx) const {
107  return (*thisHadrons)[idx];
108  }
109 
110  // Convert candidate nucleon set into output nucleus (true == success)
111  bool makeLightIon(const ClusterCandidate& aCluster);
112 
113  // Kinematics for cluster evaluations
114  G4LorentzVector getClusterMomentum(const ClusterCandidate& aCluster) const;
115  mutable G4LorentzVector pCluster; // Reusable buffer to reduce churn
116 
117  G4double maxDeltaP(const ClusterCandidate& aCluster) const;
118 
119  // Report cluster arguments for validation
120  void reportArgs(const G4String& name, const ClusterCandidate& clus) const;
121  void reportResult(const G4String& name, const G4InuclNuclei& nucl) const;
122 };
123 
124 #endif /* G4CASCADE_COALESCENCE_HH */
const XML_Char * name
Definition: expat.h:151
void setVerboseLevel(G4int verbose)
int G4int
Definition: G4Types.hh:78
G4CascadeCoalescence(G4int verbose=0)
void FindClusters(G4CollisionOutput &finalState)
double G4double
Definition: G4Types.hh:76