Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4CascadeCoalescence.cc
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// The relative-momentum cut offs for each cluster type may be set with
32// environment variables:
33// DPMAX_2CLUSTER 0.090 GeV/c for deuterons
34// DPMAX_3CLUSTER 0.108 GeV/c for tritons, He-3
35// DPMAX_4CLUSTER 0.115 GeV/c for alphas
36//
37// 20110917 Michael Kelsey
38// 20110920 M. Kelsey -- Use environment variables to set momentum cuts for
39// tuning, replace polymorphic argument lists with use of
40// "ClusterCandidate"
41// 20110922 M. Kelsey -- Follow G4InuclParticle::print(ostream&) migration
42// 20110927 M. Kelsey -- Bug fix; missing <iterator> header, strtof -> strtod
43// 20120822 M. Kelsey -- Move envvars to G4CascadeParameters.
44// 20130314 M. Kelsey -- Restore null initializer and if-block for _TLS_.
45// 20130326 M. Kelsey -- Replace _TLS_ with mutable data member buffer.
46// 20170406 M. Kelsey -- Remove recursive tryCluster() calls (redundant),
47// and remove use of triedClusters registry.
48
51#include "G4CollisionOutput.hh"
53#include "G4InuclNuclei.hh"
54#include "G4InuclParticle.hh"
57#include "G4ThreeVector.hh"
58#include <vector>
59#include <numeric>
60#include <algorithm>
61#include <iterator>
62
63
64// Constructor and Destructor
65
67 : verboseLevel(verbose), thisFinalState(0), thisHadrons(0),
68 dpMaxDoublet(G4CascadeParameters::dpMaxDoublet()),
69 dpMaxTriplet(G4CascadeParameters::dpMaxTriplet()),
70 dpMaxAlpha(G4CascadeParameters::dpMaxAlpha()) {}
71
73
74
75// Final state particle list is modified directly
76
78 if (verboseLevel)
79 G4cout << " >>> G4CascadeCoalescence::FindClusters()" << G4endl;
80
81 thisFinalState = &finalState; // Save pointers for use in processing
82 thisHadrons = &finalState.getOutgoingParticles();
83
84 if (verboseLevel>1) thisFinalState->printCollisionOutput(); // Before
85
86 selectCandidates();
87 createNuclei();
88 removeNucleons();
89
90 if (verboseLevel>1) thisFinalState->printCollisionOutput(); // After
91}
92
93
94// Scan list for possible nucleon clusters
95
96void G4CascadeCoalescence::selectCandidates() {
97 if (verboseLevel)
98 G4cout << " >>> G4CascadeCoalescence::selectCandidates()" << G4endl;
99
100 allClusters.clear();
101 usedNucleons.clear();
102
103 size_t nHad = thisHadrons->size();
104 for (size_t idx1=0; idx1<nHad; idx1++) {
105 if (!getHadron(idx1).nucleon()) continue;
106 for (size_t idx2=idx1+1; idx2<nHad; idx2++) {
107 if (!getHadron(idx2).nucleon()) continue;
108 for (size_t idx3=idx2+1; idx3<nHad; idx3++) {
109 if (!getHadron(idx3).nucleon()) continue;
110 for (size_t idx4=idx3+1; idx4<nHad; idx4++) {
111 if (!getHadron(idx4).nucleon()) continue;
112 tryClusters(idx1, idx2, idx3, idx4);
113 }
114 tryClusters(idx1, idx2, idx3); // If idx4 loop was empty
115 }
116 tryClusters(idx1, idx2); // If idx3 loop was empty
117 }
118 }
119
120 // All potential candidates built; report statistics
121 if (verboseLevel>1) {
122 G4cout << " Found " << allClusters.size() << " candidates, using "
123 << usedNucleons.size() << " nucleons" << G4endl;
124 }
125}
126
127
128// Do combinatorics of current set of four, skip nucleons already used
129
130void G4CascadeCoalescence::tryClusters(size_t idx1, size_t idx2,
131 size_t idx3, size_t idx4) {
132 if (nucleonUsed(idx1) || nucleonUsed(idx2) ||
133 nucleonUsed(idx3) || nucleonUsed(idx4)) return;
134
135 fillCluster(idx1,idx2,idx3,idx4);
136 if (verboseLevel>1) reportArgs("tryClusters",thisCluster);
137
138 if (goodCluster(thisCluster)) {
139 allClusters.push_back(thisCluster);
140 usedNucleons.insert(idx1);
141 usedNucleons.insert(idx2);
142 usedNucleons.insert(idx3);
143 usedNucleons.insert(idx4);
144 }
145}
146
147void
148G4CascadeCoalescence::tryClusters(size_t idx1, size_t idx2, size_t idx3) {
149 if (nucleonUsed(idx1) || nucleonUsed(idx2) || nucleonUsed(idx3)) return;
150
151 fillCluster(idx1,idx2,idx3);
152 if (verboseLevel>1) reportArgs("tryClusters",thisCluster);
153
154 if (goodCluster(thisCluster)) {
155 allClusters.push_back(thisCluster);
156 usedNucleons.insert(idx1);
157 usedNucleons.insert(idx2);
158 usedNucleons.insert(idx3);
159 }
160}
161
162void
163G4CascadeCoalescence::tryClusters(size_t idx1, size_t idx2) {
164 if (nucleonUsed(idx1) || nucleonUsed(idx2)) return;
165
166 fillCluster(idx1,idx2);
167 if (verboseLevel>1) reportArgs("tryClusters",thisCluster);
168
169 if (goodCluster(thisCluster)) {
170 allClusters.push_back(thisCluster);
171 usedNucleons.insert(idx1);
172 usedNucleons.insert(idx2);
173 }
174}
175
176
177// Process list of candidate clusters into light ions
178
179void G4CascadeCoalescence::createNuclei() {
180 if (verboseLevel)
181 G4cout << " >>> G4CascadeCoalescence::createNuclei()" << G4endl;
182
183 usedNucleons.clear();
184
185 for (size_t i=0; i<allClusters.size(); i++) {
186 if (verboseLevel>1) G4cout << " attempting candidate #" << i << G4endl;
187
188 const ClusterCandidate& cand = allClusters[i];
189 if (makeLightIon(cand)) { // Success, put ion in output
190 thisFinalState->addOutgoingNucleus(thisLightIon);
191 usedNucleons.insert(cand.begin(), cand.end());
192 }
193 }
194}
195
196
197// Remove nucleons indexed in "usedNucleons" from output
198
199void G4CascadeCoalescence::removeNucleons() {
200 if (verboseLevel>1)
201 G4cout << " >>> G4CascadeCoalescence::removeNucleons()" << G4endl;
202
203 // Remove nucleons from output from last to first (to preserve indexing)
204 for (auto usedIter = usedNucleons.crbegin(); usedIter != usedNucleons.crend(); ++usedIter)
205 thisFinalState->removeOutgoingParticle((G4int)*usedIter);
206
207 usedNucleons.clear();
208}
209
210
211// Compute momentum of whole cluster
212
214G4CascadeCoalescence::getClusterMomentum(const ClusterCandidate& aCluster) const {
215 pCluster.set(0.,0.,0.,0.);
216 for (size_t i=0; i<aCluster.size(); i++)
217 pCluster += getHadron(aCluster[i]).getMomentum();
218
219 return pCluster;
220}
221
222
223// Determine magnitude of largest momentum in CM frame
224
225G4double G4CascadeCoalescence::maxDeltaP(const ClusterCandidate& aCluster) const {
226 if (verboseLevel>1) reportArgs("maxDeltaP", aCluster);
227
228 G4ThreeVector boost = getClusterMomentum(aCluster).boostVector();
229
230 G4double dp, maxDP = -1.;
231 for (size_t i=0; i<aCluster.size(); i++) {
232 const G4InuclElementaryParticle& nucl = getHadron(aCluster[i]);
233
234 // NOTE: getMomentum() returns by value, event kinematics are not changed
235 if ((dp = nucl.getMomentum().boost(-boost).rho()) > maxDP) maxDP = dp;
236 }
237
238 if (verboseLevel>1) G4cout << " maxDP = " << maxDP << G4endl;
239
240 return maxDP;
241}
242
243
244// Compute "cluster type code" as sum of nucleon codes
245
246G4int G4CascadeCoalescence::
247clusterType(const ClusterCandidate& aCluster) const {
248 G4int type = 0;
249 for (size_t i=0; i<aCluster.size(); i++) {
250 const G4InuclElementaryParticle& had = getHadron(aCluster[i]);
251 type += had.nucleon() ? had.type() : 0;
252 }
253
254 return type;
255}
256
257
258// Create cluster candidate with listed indices
259
260void G4CascadeCoalescence::fillCluster(size_t idx1, size_t idx2) {
261 thisCluster.clear();
262 thisCluster.push_back(idx1);
263 thisCluster.push_back(idx2);
264}
265
266void G4CascadeCoalescence::fillCluster(size_t idx1, size_t idx2, size_t idx3) {
267 thisCluster.clear();
268 thisCluster.push_back(idx1);
269 thisCluster.push_back(idx2);
270 thisCluster.push_back(idx3);
271}
272
273void G4CascadeCoalescence::fillCluster(size_t idx1, size_t idx2, size_t idx3,
274 size_t idx4) {
275 thisCluster.clear();
276 thisCluster.push_back(idx1);
277 thisCluster.push_back(idx2);
278 thisCluster.push_back(idx3);
279 thisCluster.push_back(idx4);
280}
281
282
283
284// Make sure all candidates in cluster are nucleons
285
286bool G4CascadeCoalescence::allNucleons(const ClusterCandidate& clus) const {
287 bool result = true;
288 for (size_t i=0; i<clus.size(); i++)
289 result &= getHadron(clus[0]).nucleon();
290 return result;
291}
292
293
294// Determine if collection of nucleons can form a light ion
295
296bool G4CascadeCoalescence::goodCluster(const ClusterCandidate& clus) const {
297 if (verboseLevel>2) reportArgs("goodCluster?",clus);
298
299 if (!allNucleons(clus)) return false;
300
301 if (clus.size() == 2) // Deuterons (pn)
302 return (clusterType(clus) == 3 && maxDeltaP(clus) < dpMaxDoublet);
303
304 if (clus.size() == 3) // Tritons or He-3
305 return ((clusterType(clus) == 4 || clusterType(clus) == 5) // ppn OR pnn
306 && maxDeltaP(clus) < dpMaxTriplet);
307
308 if (clus.size() == 4) // Alphas (ppnn)
309 return (clusterType(clus) == 6 && maxDeltaP(clus) < dpMaxAlpha);
310
311 return false;
312}
313
314
315
316// Convert candidate nucleon set into output nucleus
317
318bool G4CascadeCoalescence::makeLightIon(const ClusterCandidate& aCluster) {
319 if (verboseLevel>1) reportArgs("makeLightIon",aCluster);
320
321 thisLightIon.clear(); // Initialize nucleus buffer before filling
322
323 if (aCluster.size()<2) return false; // Sanity check
324
325 G4int A = (G4int)aCluster.size();
326 G4int Z = -1;
327
328 G4int type = clusterType(aCluster);
329 if (A==2 && type==3) Z = 1; // Deuteron (pn)
330 if (A==3 && type==5) Z = 1; // Triton (pnn)
331 if (A==3 && type==4) Z = 2; // He-3 (ppn)
332 if (A==4 && type==6) Z = 2; // He-4/alpha (ppnn)
333
334 if (Z < 0) return false; // Invalid cluster content
335
336 // NOTE: Four-momentum will not be conserved due to binding energy
337 thisLightIon.fill(getClusterMomentum(aCluster), A, Z, 0.,
339
340 if (verboseLevel>1) reportResult("makeLightIon output",thisLightIon);
341 return true;
342}
343
344
345// Report cluster arguments for validation
346
347void G4CascadeCoalescence::reportArgs(const G4String& name,
348 const ClusterCandidate& aCluster) const {
349 G4cout << " >>> G4CascadeCoalescence::" << name << " ";
350 std::copy(aCluster.begin(), aCluster.end(),
351 std::ostream_iterator<size_t>(G4cout, " "));
352 G4cout << G4endl;
353
354 if (verboseLevel>2) {
355 for (size_t i=0; i<aCluster.size(); i++)
356 G4cout << getHadron(aCluster[i]) << G4endl;
357 }
358}
359
360void G4CascadeCoalescence::reportResult(const G4String& name,
361 const G4InuclNuclei& nucl) const {
362 G4cout << " >>> G4CascadeCoalescence::" << name << G4endl << nucl << G4endl;
363}
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4double A[17]
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
Hep3Vector boostVector() const
HepLorentzVector & boost(double, double, double)
void set(double x, double y, double z, double t)
void FindClusters(G4CollisionOutput &finalState)
G4CascadeCoalescence(G4int verbose=0)
void printCollisionOutput(std::ostream &os=G4cout) const
const std::vector< G4InuclElementaryParticle > & getOutgoingParticles() const
void addOutgoingNucleus(const G4InuclNuclei &nuclei)
void removeOutgoingParticle(G4int index)
void fill(G4int a, G4int z, G4double exc=0., Model model=DefaultModel)
G4LorentzVector getMomentum() const
const char * name(G4int ptype)