67 : verboseLevel(verbose), thisFinalState(0), thisHadrons(0),
79 G4cout <<
" >>> G4CascadeCoalescence::FindClusters()" <<
G4endl;
81 thisFinalState = &finalState;
96 void G4CascadeCoalescence::selectCandidates() {
98 G4cout <<
" >>> G4CascadeCoalescence::selectCandidates()" <<
G4endl;
101 usedNucleons.clear();
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);
114 tryClusters(idx1, idx2, idx3);
116 tryClusters(idx1, idx2);
121 if (verboseLevel>1) {
122 G4cout <<
" Found " << allClusters.size() <<
" candidates, using "
123 << usedNucleons.size() <<
" nucleons" <<
G4endl;
130 void 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;
135 fillCluster(idx1,idx2,idx3,idx4);
136 if (verboseLevel>1) reportArgs(
"tryClusters",thisCluster);
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);
148 G4CascadeCoalescence::tryClusters(
size_t idx1,
size_t idx2,
size_t idx3) {
149 if (nucleonUsed(idx1) || nucleonUsed(idx2) || nucleonUsed(idx3))
return;
151 fillCluster(idx1,idx2,idx3);
152 if (verboseLevel>1) reportArgs(
"tryClusters",thisCluster);
154 if (goodCluster(thisCluster)) {
155 allClusters.push_back(thisCluster);
156 usedNucleons.insert(idx1);
157 usedNucleons.insert(idx2);
158 usedNucleons.insert(idx3);
163 G4CascadeCoalescence::tryClusters(
size_t idx1,
size_t idx2) {
164 if (nucleonUsed(idx1) || nucleonUsed(idx2))
return;
166 fillCluster(idx1,idx2);
167 if (verboseLevel>1) reportArgs(
"tryClusters",thisCluster);
169 if (goodCluster(thisCluster)) {
170 allClusters.push_back(thisCluster);
171 usedNucleons.insert(idx1);
172 usedNucleons.insert(idx2);
179 void G4CascadeCoalescence::createNuclei() {
181 G4cout <<
" >>> G4CascadeCoalescence::createNuclei()" <<
G4endl;
183 usedNucleons.clear();
185 for (
size_t i=0; i<allClusters.size(); i++) {
186 if (verboseLevel>1)
G4cout <<
" attempting candidate #" << i <<
G4endl;
188 const ClusterCandidate& cand = allClusters[i];
189 if (makeLightIon(cand)) {
191 usedNucleons.insert(cand.begin(), cand.end());
199 void G4CascadeCoalescence::removeNucleons() {
201 G4cout <<
" >>> G4CascadeCoalescence::removeNucleons()" <<
G4endl;
204 std::set<size_t>::reverse_iterator usedIter;
205 for (usedIter = usedNucleons.rbegin(); usedIter != usedNucleons.rend(); ++usedIter)
208 usedNucleons.clear();
215 G4CascadeCoalescence::getClusterMomentum(
const ClusterCandidate& aCluster)
const {
216 pCluster.
set(0.,0.,0.,0.);
217 for (
size_t i=0; i<aCluster.size(); i++)
226 G4double G4CascadeCoalescence::maxDeltaP(
const ClusterCandidate& aCluster)
const {
227 if (verboseLevel>1) reportArgs(
"maxDeltaP", aCluster);
232 for (
size_t i=0; i<aCluster.size(); i++) {
239 if (verboseLevel>1)
G4cout <<
" maxDP = " << maxDP <<
G4endl;
247 G4int G4CascadeCoalescence::
248 clusterType(
const ClusterCandidate& aCluster)
const {
250 for (
size_t i=0; i<aCluster.size(); i++) {
261 void G4CascadeCoalescence::fillCluster(
size_t idx1,
size_t idx2) {
263 thisCluster.push_back(idx1);
264 thisCluster.push_back(idx2);
267 void G4CascadeCoalescence::fillCluster(
size_t idx1,
size_t idx2,
size_t idx3) {
269 thisCluster.push_back(idx1);
270 thisCluster.push_back(idx2);
271 thisCluster.push_back(idx3);
274 void G4CascadeCoalescence::fillCluster(
size_t idx1,
size_t idx2,
size_t idx3,
277 thisCluster.push_back(idx1);
278 thisCluster.push_back(idx2);
279 thisCluster.push_back(idx3);
280 thisCluster.push_back(idx4);
287 bool G4CascadeCoalescence::allNucleons(
const ClusterCandidate& clus)
const {
289 for (
size_t i=0; i<clus.size(); i++)
290 result &= getHadron(clus[0]).
nucleon();
297 bool G4CascadeCoalescence::goodCluster(
const ClusterCandidate& clus)
const {
298 if (verboseLevel>2) reportArgs(
"goodCluster?",clus);
300 if (!allNucleons(clus))
return false;
302 if (clus.size() == 2)
303 return (clusterType(clus) == 3 && maxDeltaP(clus) < dpMaxDoublet);
305 if (clus.size() == 3)
306 return ((clusterType(clus) == 4 || clusterType(clus) == 5)
307 && maxDeltaP(clus) < dpMaxTriplet);
309 if (clus.size() == 4)
310 return (clusterType(clus) == 6 && maxDeltaP(clus) < dpMaxAlpha);
319 bool G4CascadeCoalescence::makeLightIon(
const ClusterCandidate& aCluster) {
320 if (verboseLevel>1) reportArgs(
"makeLightIon",aCluster);
322 thisLightIon.
clear();
324 if (aCluster.size()<2)
return false;
326 G4int A = aCluster.size();
329 G4int type = clusterType(aCluster);
330 if (A==2 && type==3) Z = 1;
331 if (A==3 && type==5) Z = 1;
332 if (A==3 && type==4) Z = 2;
333 if (A==4 && type==6) Z = 2;
335 if (Z < 0)
return false;
338 thisLightIon.
fill(getClusterMomentum(aCluster), A, Z, 0.,
341 if (verboseLevel>1) reportResult(
"makeLightIon output",thisLightIon);
348 void G4CascadeCoalescence::reportArgs(
const G4String&
name,
349 const ClusterCandidate& aCluster)
const {
350 G4cout <<
" >>> G4CascadeCoalescence::" << name <<
" ";
351 std::copy(aCluster.begin(), aCluster.end(),
352 std::ostream_iterator<size_t>(
G4cout,
" "));
355 if (verboseLevel>2) {
356 for (
size_t i=0; i<aCluster.size(); i++)
361 void G4CascadeCoalescence::reportResult(
const G4String& name,
363 G4cout <<
" >>> G4CascadeCoalescence::" << name << G4endl << nucl <<
G4endl;
G4double G4ParticleHPJENDLHEData::G4double result
void fill(G4int a, G4int z, G4double exc=0., Model model=DefaultModel)
Hep3Vector boostVector() const
G4LorentzVector getMomentum() const
void printCollisionOutput(std::ostream &os=G4cout) const
G4GLOB_DLL std::ostream G4cout
double A(double temperature)
G4bool nucleon(G4int ityp)
HepLorentzVector & boost(double, double, double)
void removeOutgoingParticle(G4int index)
G4CascadeCoalescence(G4int verbose=0)
void set(double x, double y, double z, double t)
const std::vector< G4InuclElementaryParticle > & getOutgoingParticles() const
void addOutgoingNucleus(const G4InuclNuclei &nuclei)
void FindClusters(G4CollisionOutput &finalState)
virtual ~G4CascadeCoalescence()