62typedef std::vector<G4InuclElementaryParticle>
hadronList;
77 : verboseLevel(verbose), thisFinalState(0), thisHadrons(0) {}
86 G4cout <<
" >>> G4CascadeCoalescence::FindClusters()" <<
G4endl;
88 thisFinalState = &finalState;
103void G4CascadeCoalescence::selectCandidates() {
105 G4cout <<
" >>> G4CascadeCoalescence::selectCandidates()" <<
G4endl;
107 triedClusters.clear();
109 usedNucleons.clear();
111 size_t nHad = thisHadrons->size();
112 for (
size_t idx1=0; idx1<nHad; idx1++) {
113 if (!getHadron(idx1).nucleon())
continue;
114 for (
size_t idx2=idx1+1; idx2<nHad; idx2++) {
115 if (!getHadron(idx2).nucleon())
continue;
116 for (
size_t idx3=idx2+1; idx3<nHad; idx3++) {
117 if (!getHadron(idx3).nucleon())
continue;
118 for (
size_t idx4=idx3+1; idx4<nHad; idx4++) {
119 if (!getHadron(idx4).nucleon())
continue;
120 tryClusters(idx1, idx2, idx3, idx4);
122 tryClusters(idx1, idx2, idx3);
124 tryClusters(idx1, idx2);
129 if (verboseLevel>1) {
130 G4cout <<
" Found " << allClusters.size() <<
" candidates, using "
131 << usedNucleons.size() <<
" nucleons" <<
G4endl;
138void G4CascadeCoalescence::tryClusters(
size_t idx1,
size_t idx2,
139 size_t idx3,
size_t idx4) {
140 fillCluster(idx1,idx2,idx3,idx4);
141 if (clusterTried(thisCluster))
return;
144 G4cout <<
" >>> G4CascadeCoalescence::tryClusters " << idx1 <<
" " << idx2
145 <<
" " << idx3 <<
" " << idx4 <<
G4endl;
147 triedClusters.insert(clusterHash(thisCluster));
149 if (!nucleonUsed(idx1) && !nucleonUsed(idx2) &&
150 !nucleonUsed(idx3) && !nucleonUsed(idx4)) {
151 if (goodCluster(thisCluster)) {
152 allClusters.push_back(thisCluster);
153 usedNucleons.insert(idx1);
154 usedNucleons.insert(idx2);
155 usedNucleons.insert(idx3);
156 usedNucleons.insert(idx4);
161 tryClusters(idx2,idx3,idx4);
162 tryClusters(idx1,idx3,idx4);
163 tryClusters(idx1,idx2,idx4);
164 tryClusters(idx1,idx2,idx3);
168G4CascadeCoalescence::tryClusters(
size_t idx1,
size_t idx2,
size_t idx3) {
169 fillCluster(idx1,idx2,idx3);
170 if (clusterTried(thisCluster))
return;
173 G4cout <<
" >>> G4CascadeCoalescence::tryClusters " << idx1 <<
" " << idx2
176 triedClusters.insert(clusterHash(thisCluster));
178 if (!nucleonUsed(idx1) && !nucleonUsed(idx2) && !nucleonUsed(idx3)) {
179 fillCluster(idx1,idx2,idx3);
180 if (goodCluster(thisCluster)) {
181 allClusters.push_back(thisCluster);
182 usedNucleons.insert(idx1);
183 usedNucleons.insert(idx2);
184 usedNucleons.insert(idx3);
189 tryClusters(idx2,idx3);
190 tryClusters(idx1,idx3);
191 tryClusters(idx1,idx2);
195G4CascadeCoalescence::tryClusters(
size_t idx1,
size_t idx2) {
196 if (nucleonUsed(idx1) || nucleonUsed(idx2))
return;
198 fillCluster(idx1,idx2);
199 if (clusterTried(thisCluster))
return;
202 G4cout <<
" >>> G4CascadeCoalescence::tryClusters " << idx1 <<
" " << idx2
205 triedClusters.insert(clusterHash(thisCluster));
207 fillCluster(idx1,idx2);
208 if (goodCluster(thisCluster)) {
209 allClusters.push_back(thisCluster);
210 usedNucleons.insert(idx1);
211 usedNucleons.insert(idx2);
218void G4CascadeCoalescence::createNuclei() {
220 G4cout <<
" >>> G4CascadeCoalescence::createNuclei()" <<
G4endl;
222 usedNucleons.clear();
224 for (
size_t i=0; i<allClusters.size(); i++) {
225 if (verboseLevel>1)
G4cout <<
" attempting candidate #" << i <<
G4endl;
227 const ClusterCandidate& cand = allClusters[i];
228 if (makeLightIon(cand)) {
230 usedNucleons.insert(cand.begin(), cand.end());
238void G4CascadeCoalescence::removeNucleons() {
240 G4cout <<
" >>> G4CascadeCoalescence::removeNucleons()" <<
G4endl;
243 std::set<size_t>::reverse_iterator usedIter;
244 for (usedIter = usedNucleons.rbegin(); usedIter != usedNucleons.rend(); ++usedIter)
247 usedNucleons.clear();
254G4CascadeCoalescence::getClusterMomentum(
const ClusterCandidate& aCluster)
const {
256 ptot.
set(0.,0.,0.,0.);
257 for (
size_t i=0; i<aCluster.size(); i++)
266G4double G4CascadeCoalescence::maxDeltaP(
const ClusterCandidate& aCluster)
const {
267 if (verboseLevel>1) reportArgs(
"maxDeltaP", aCluster);
273 for (
size_t i=0; i<aCluster.size(); i++) {
280 if (verboseLevel>1)
G4cout <<
" maxDP = " << maxDP <<
G4endl;
288G4int G4CascadeCoalescence::
289clusterType(
const ClusterCandidate& aCluster)
const {
291 for (
size_t i=0; i<aCluster.size(); i++) {
301size_t G4CascadeCoalescence::
302clusterHash(
const ClusterCandidate& aCluster)
const {
304 for (
size_t i=0; i<aCluster.size(); i++) {
305 hash = hash*1000 + aCluster[i];
314void G4CascadeCoalescence::fillCluster(
size_t idx1,
size_t idx2) {
316 thisCluster.push_back(idx1);
317 thisCluster.push_back(idx2);
320void G4CascadeCoalescence::fillCluster(
size_t idx1,
size_t idx2,
size_t idx3) {
322 thisCluster.push_back(idx1);
323 thisCluster.push_back(idx2);
324 thisCluster.push_back(idx3);
327void G4CascadeCoalescence::fillCluster(
size_t idx1,
size_t idx2,
size_t idx3,
330 thisCluster.push_back(idx1);
331 thisCluster.push_back(idx2);
332 thisCluster.push_back(idx3);
333 thisCluster.push_back(idx4);
340bool G4CascadeCoalescence::allNucleons(
const ClusterCandidate& clus)
const {
342 for (
size_t i=0; i<clus.size(); i++)
343 result &= getHadron(clus[0]).
nucleon();
350bool G4CascadeCoalescence::goodCluster(
const ClusterCandidate& clus)
const {
351 if (verboseLevel>2) reportArgs(
"goodCluster?",clus);
353 if (!allNucleons(clus))
return false;
355 if (clus.size() == 2)
356 return (clusterType(clus) == 3 && maxDeltaP(clus) < dpMaxDoublet);
358 if (clus.size() == 3)
359 return ((clusterType(clus) == 4 || clusterType(clus) == 5)
360 && maxDeltaP(clus) < dpMaxTriplet);
362 if (clus.size() == 4)
363 return (clusterType(clus) == 6 && maxDeltaP(clus) < dpMaxAlpha);
372bool G4CascadeCoalescence::makeLightIon(
const ClusterCandidate& aCluster) {
373 if (verboseLevel>1) reportArgs(
"makeLightIon",aCluster);
375 thisLightIon.
clear();
377 if (aCluster.size()<2)
return false;
379 G4int A = aCluster.size();
382 G4int type = clusterType(aCluster);
383 if (A==2 && type==3) Z = 1;
384 if (A==3 && type==5) Z = 1;
385 if (A==3 && type==4) Z = 2;
386 if (A==4 && type==6) Z = 2;
388 if (Z < 0)
return false;
391 thisLightIon.
fill(getClusterMomentum(aCluster), A, Z, 0.,
394 if (verboseLevel>1) reportResult(
"makeLightIon output",thisLightIon);
401void G4CascadeCoalescence::reportArgs(
const G4String& name,
402 const ClusterCandidate& aCluster)
const {
403 G4cout <<
" >>> G4CascadeCoalescence::" << name <<
" ";
404 std::copy(aCluster.begin(), aCluster.end(),
405 std::ostream_iterator<size_t>(
G4cout,
" "));
408 if (verboseLevel>2) {
409 for (
size_t i=0; i<aCluster.size(); i++)
414void G4CascadeCoalescence::reportResult(
const G4String& name,
std::vector< G4InuclElementaryParticle > hadronList
hadronList::const_iterator hadronIter
G4DLLIMPORT 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)
virtual ~G4CascadeCoalescence()
static G4double dpMaxDoublet()
static G4double dpMaxTriplet()
static G4double dpMaxAlpha()
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