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