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 | // $Id: G4RPGTwoCluster.cc,v 1.2 2007/08/15 20:38:48 dennis Exp $ |
---|
27 | // GEANT4 tag $Name: geant4-09-01-patch-02 $ |
---|
28 | // |
---|
29 | |
---|
30 | #include "G4RPGTwoCluster.hh" |
---|
31 | #include "Randomize.hh" |
---|
32 | #include "G4Poisson.hh" |
---|
33 | #include <iostream> |
---|
34 | #include "G4HadReentrentException.hh" |
---|
35 | #include <signal.h> |
---|
36 | |
---|
37 | |
---|
38 | G4RPGTwoCluster::G4RPGTwoCluster() |
---|
39 | : G4RPGReaction() {} |
---|
40 | |
---|
41 | |
---|
42 | G4bool G4RPGTwoCluster:: |
---|
43 | ReactionStage(const G4HadProjectile* originalIncident, |
---|
44 | G4ReactionProduct& modifiedOriginal, |
---|
45 | G4bool& incidentHasChanged, |
---|
46 | const G4DynamicParticle* originalTarget, |
---|
47 | G4ReactionProduct& targetParticle, |
---|
48 | G4bool& targetHasChanged, |
---|
49 | const G4Nucleus& targetNucleus, |
---|
50 | G4ReactionProduct& currentParticle, |
---|
51 | G4FastVector<G4ReactionProduct,256>& vec, |
---|
52 | G4int& vecLen, |
---|
53 | G4bool leadFlag, |
---|
54 | G4ReactionProduct& leadingStrangeParticle) |
---|
55 | { |
---|
56 | // Derived from H. Fesefeldt's FORTRAN code TWOCLU |
---|
57 | // |
---|
58 | // A simple two cluster model is used to generate x- and pt- values for |
---|
59 | // incident, target, and all secondary particles. |
---|
60 | // This should be sufficient for low energy interactions. |
---|
61 | // |
---|
62 | |
---|
63 | G4int i; |
---|
64 | G4ParticleDefinition *aProton = G4Proton::Proton(); |
---|
65 | G4ParticleDefinition *aNeutron = G4Neutron::Neutron(); |
---|
66 | G4ParticleDefinition *aPiPlus = G4PionPlus::PionPlus(); |
---|
67 | G4ParticleDefinition *aPiMinus = G4PionMinus::PionMinus(); |
---|
68 | G4ParticleDefinition *aPiZero = G4PionZero::PionZero(); |
---|
69 | G4bool veryForward = false; |
---|
70 | |
---|
71 | const G4double protonMass = aProton->GetPDGMass()/MeV; |
---|
72 | const G4double ekOriginal = modifiedOriginal.GetKineticEnergy()/GeV; |
---|
73 | const G4double etOriginal = modifiedOriginal.GetTotalEnergy()/GeV; |
---|
74 | const G4double mOriginal = modifiedOriginal.GetMass()/GeV; |
---|
75 | const G4double pOriginal = modifiedOriginal.GetMomentum().mag()/GeV; |
---|
76 | G4double targetMass = targetParticle.GetDefinition()->GetPDGMass()/GeV; |
---|
77 | G4double centerofmassEnergy = std::sqrt( mOriginal*mOriginal + |
---|
78 | targetMass*targetMass + |
---|
79 | 2.0*targetMass*etOriginal ); // GeV |
---|
80 | G4double currentMass = currentParticle.GetMass()/GeV; |
---|
81 | targetMass = targetParticle.GetMass()/GeV; |
---|
82 | |
---|
83 | if( currentMass == 0.0 && targetMass == 0.0 ) |
---|
84 | { |
---|
85 | G4double ek = currentParticle.GetKineticEnergy(); |
---|
86 | G4ThreeVector m = currentParticle.GetMomentum(); |
---|
87 | currentParticle = *vec[0]; |
---|
88 | targetParticle = *vec[1]; |
---|
89 | for( i=0; i<(vecLen-2); ++i )*vec[i] = *vec[i+2]; |
---|
90 | if(vecLen<2) |
---|
91 | { |
---|
92 | for(G4int i=0; i<vecLen; i++) delete vec[i]; |
---|
93 | vecLen = 0; |
---|
94 | throw G4HadReentrentException(__FILE__, __LINE__, |
---|
95 | "G4RPGTwoCluster::ReactionStage : Negative number of particles"); |
---|
96 | } |
---|
97 | delete vec[vecLen-1]; |
---|
98 | delete vec[vecLen-2]; |
---|
99 | vecLen -= 2; |
---|
100 | currentMass = currentParticle.GetMass()/GeV; |
---|
101 | targetMass = targetParticle.GetMass()/GeV; |
---|
102 | incidentHasChanged = true; |
---|
103 | targetHasChanged = true; |
---|
104 | currentParticle.SetKineticEnergy( ek ); |
---|
105 | currentParticle.SetMomentum( m ); |
---|
106 | veryForward = true; |
---|
107 | } |
---|
108 | |
---|
109 | const G4double atomicWeight = targetNucleus.GetN(); |
---|
110 | const G4double atomicNumber = targetNucleus.GetZ(); |
---|
111 | // |
---|
112 | // particles have been distributed in forward and backward hemispheres |
---|
113 | // in center of mass system of the hadron nucleon interaction |
---|
114 | // |
---|
115 | |
---|
116 | // Incident particle always in forward hemisphere |
---|
117 | |
---|
118 | G4int forwardCount = 1; // number of particles in forward hemisphere |
---|
119 | currentParticle.SetSide( 1 ); |
---|
120 | G4double forwardMass = currentParticle.GetMass()/GeV; |
---|
121 | G4double cMass = forwardMass; |
---|
122 | |
---|
123 | // Target particle always in backward hemisphere |
---|
124 | |
---|
125 | G4int backwardCount = 1; // number of particles in backward hemisphere |
---|
126 | targetParticle.SetSide( -1 ); |
---|
127 | G4double backwardMass = targetParticle.GetMass()/GeV; |
---|
128 | G4double bMass = backwardMass; |
---|
129 | |
---|
130 | G4int backwardNucleonCount = 1; // number of nucleons in backward hemisphere |
---|
131 | |
---|
132 | for( i=0; i<vecLen; ++i ) |
---|
133 | { |
---|
134 | if( vec[i]->GetSide() < 0 )vec[i]->SetSide( -1 ); // to take care of |
---|
135 | // case where vec has been preprocessed by GenerateXandPt |
---|
136 | // and some of them have been set to -2 or -3 |
---|
137 | if( vec[i]->GetSide() == -1 ) |
---|
138 | { |
---|
139 | ++backwardCount; |
---|
140 | backwardMass += vec[i]->GetMass()/GeV; |
---|
141 | } |
---|
142 | else |
---|
143 | { |
---|
144 | ++forwardCount; |
---|
145 | forwardMass += vec[i]->GetMass()/GeV; |
---|
146 | } |
---|
147 | } |
---|
148 | |
---|
149 | // |
---|
150 | // Add nucleons and some pions from intra-nuclear cascade |
---|
151 | // |
---|
152 | |
---|
153 | G4double term1 = std::log(centerofmassEnergy*centerofmassEnergy); |
---|
154 | if(term1 < 0) term1 = 0.0001; // making sure xtarg<0; |
---|
155 | const G4double afc = 0.312 + 0.2 * std::log(term1); |
---|
156 | G4double xtarg; |
---|
157 | if( centerofmassEnergy < 2.0+G4UniformRand() ) // added +2 below, JLC 4Jul97 |
---|
158 | xtarg = afc * (std::pow(atomicWeight,0.33)-1.0) * (2*backwardCount+vecLen+2)/2.0; |
---|
159 | else |
---|
160 | xtarg = afc * (std::pow(atomicWeight,0.33)-1.0) * (2*backwardCount); |
---|
161 | if( xtarg <= 0.0 )xtarg = 0.01; |
---|
162 | G4int nuclearExcitationCount = G4Poisson( xtarg ); |
---|
163 | |
---|
164 | if(atomicWeight<1.0001) nuclearExcitationCount = 0; |
---|
165 | G4int extraNucleonCount = 0; |
---|
166 | G4double extraMass = 0.0; |
---|
167 | G4double extraNucleonMass = 0.0; |
---|
168 | if( nuclearExcitationCount > 0 ) |
---|
169 | { |
---|
170 | G4int momentumBin = std::min( 4, G4int(pOriginal/3.0) ); |
---|
171 | const G4double nucsup[] = { 1.0, 0.8, 0.6, 0.5, 0.4 }; |
---|
172 | // |
---|
173 | // NOTE: in TWOCLU, these new particles were given negative codes |
---|
174 | // here we use NewlyAdded = true instead |
---|
175 | // |
---|
176 | for( i=0; i<nuclearExcitationCount; ++i ) |
---|
177 | { |
---|
178 | G4ReactionProduct* pVec = new G4ReactionProduct(); |
---|
179 | if( G4UniformRand() < nucsup[momentumBin] ) // add proton or neutron |
---|
180 | { |
---|
181 | if( G4UniformRand() > 1.0-atomicNumber/atomicWeight ) |
---|
182 | pVec->SetDefinition( aProton ); |
---|
183 | else |
---|
184 | pVec->SetDefinition( aNeutron ); |
---|
185 | ++backwardNucleonCount; |
---|
186 | ++extraNucleonCount; |
---|
187 | extraNucleonMass += pVec->GetMass()/GeV; |
---|
188 | } |
---|
189 | else |
---|
190 | { // add a pion |
---|
191 | G4double ran = G4UniformRand(); |
---|
192 | if( ran < 0.3181 ) |
---|
193 | pVec->SetDefinition( aPiPlus ); |
---|
194 | else if( ran < 0.6819 ) |
---|
195 | pVec->SetDefinition( aPiZero ); |
---|
196 | else |
---|
197 | pVec->SetDefinition( aPiMinus ); |
---|
198 | } |
---|
199 | pVec->SetSide( -2 ); // backside particle |
---|
200 | extraMass += pVec->GetMass()/GeV; |
---|
201 | pVec->SetNewlyAdded( true ); |
---|
202 | vec.SetElement( vecLen++, pVec ); |
---|
203 | } |
---|
204 | } |
---|
205 | |
---|
206 | // Masses of particles added from cascade not included in energy balance |
---|
207 | G4double forwardEnergy = (centerofmassEnergy-cMass-bMass)/2.0 +cMass - forwardMass; |
---|
208 | G4double backwardEnergy = (centerofmassEnergy-cMass-bMass)/2.0 +bMass - backwardMass; |
---|
209 | G4double eAvailable = centerofmassEnergy - (forwardMass+backwardMass); |
---|
210 | G4bool secondaryDeleted; |
---|
211 | G4double pMass; |
---|
212 | |
---|
213 | while( eAvailable <= 0.0 ) // must eliminate a particle |
---|
214 | { |
---|
215 | secondaryDeleted = false; |
---|
216 | for( i=(vecLen-1); i>=0; --i ) |
---|
217 | { |
---|
218 | if( vec[i]->GetSide() == 1 && vec[i]->GetMayBeKilled()) |
---|
219 | { |
---|
220 | pMass = vec[i]->GetMass()/GeV; |
---|
221 | for( G4int j=i; j<(vecLen-1); ++j )*vec[j] = *vec[j+1]; // shift up |
---|
222 | --forwardCount; |
---|
223 | forwardEnergy += pMass; |
---|
224 | forwardMass -= pMass; |
---|
225 | secondaryDeleted = true; |
---|
226 | break; |
---|
227 | } |
---|
228 | else if( vec[i]->GetSide() == -1 && vec[i]->GetMayBeKilled()) |
---|
229 | { |
---|
230 | pMass = vec[i]->GetMass()/GeV; |
---|
231 | for( G4int j=i; j<(vecLen-1); ++j )*vec[j] = *vec[j+1]; // shift up |
---|
232 | --backwardCount; |
---|
233 | backwardEnergy += pMass; |
---|
234 | backwardMass -= pMass; |
---|
235 | secondaryDeleted = true; |
---|
236 | break; |
---|
237 | } |
---|
238 | } // breaks go down to here |
---|
239 | |
---|
240 | if( secondaryDeleted ) |
---|
241 | { |
---|
242 | delete vec[vecLen-1]; |
---|
243 | --vecLen; |
---|
244 | // DEBUGGING --> DumpFrames::DumpFrame(vec, vecLen); |
---|
245 | } |
---|
246 | else |
---|
247 | { |
---|
248 | if( vecLen == 0 ) return false; // all secondaries have been eliminated |
---|
249 | if( targetParticle.GetSide() == -1 ) |
---|
250 | { |
---|
251 | pMass = targetParticle.GetMass()/GeV; |
---|
252 | targetParticle = *vec[0]; |
---|
253 | for( G4int j=0; j<(vecLen-1); ++j )*vec[j] = *vec[j+1]; // shift up |
---|
254 | --backwardCount; |
---|
255 | backwardEnergy += pMass; |
---|
256 | backwardMass -= pMass; |
---|
257 | secondaryDeleted = true; |
---|
258 | } |
---|
259 | else if( targetParticle.GetSide() == 1 ) |
---|
260 | { |
---|
261 | pMass = targetParticle.GetMass()/GeV; |
---|
262 | targetParticle = *vec[0]; |
---|
263 | for( G4int j=0; j<(vecLen-1); ++j )*vec[j] = *vec[j+1]; // shift up |
---|
264 | --forwardCount; |
---|
265 | forwardEnergy += pMass; |
---|
266 | forwardMass -= pMass; |
---|
267 | secondaryDeleted = true; |
---|
268 | } |
---|
269 | |
---|
270 | if( secondaryDeleted ) |
---|
271 | { |
---|
272 | delete vec[vecLen-1]; |
---|
273 | --vecLen; |
---|
274 | } |
---|
275 | else |
---|
276 | { |
---|
277 | if( currentParticle.GetSide() == -1 ) |
---|
278 | { |
---|
279 | pMass = currentParticle.GetMass()/GeV; |
---|
280 | currentParticle = *vec[0]; |
---|
281 | for( G4int j=0; j<(vecLen-1); ++j )*vec[j] = *vec[j+1]; // shift up |
---|
282 | --backwardCount; |
---|
283 | backwardEnergy += pMass; |
---|
284 | backwardMass -= pMass; |
---|
285 | secondaryDeleted = true; |
---|
286 | } |
---|
287 | else if( currentParticle.GetSide() == 1 ) |
---|
288 | { |
---|
289 | pMass = currentParticle.GetMass()/GeV; |
---|
290 | currentParticle = *vec[0]; |
---|
291 | for( G4int j=0; j<(vecLen-1); ++j )*vec[j] = *vec[j+1]; // shift up |
---|
292 | --forwardCount; |
---|
293 | forwardEnergy += pMass; |
---|
294 | forwardMass -= pMass; |
---|
295 | secondaryDeleted = true; |
---|
296 | } |
---|
297 | if( secondaryDeleted ) |
---|
298 | { |
---|
299 | delete vec[vecLen-1]; |
---|
300 | --vecLen; |
---|
301 | } |
---|
302 | else break; |
---|
303 | |
---|
304 | } // secondary not deleted |
---|
305 | } // secondary not deleted |
---|
306 | |
---|
307 | eAvailable = centerofmassEnergy - (forwardMass+backwardMass); |
---|
308 | } // while |
---|
309 | |
---|
310 | // |
---|
311 | // This is the start of the TwoCluster function |
---|
312 | // Choose multi-particle resonance masses by sampling |
---|
313 | // P(M) = gc[g(M-M0)]**(c-1) *exp[-(g(M-M0))**c] |
---|
314 | // for M > M0 |
---|
315 | // |
---|
316 | // Use for the forward and backward clusters, but not |
---|
317 | // the cascade cluster |
---|
318 | |
---|
319 | const G4double cpar[] = { 1.60, 1.35, 1.15, 1.10 }; |
---|
320 | const G4double gpar[] = { 2.60, 1.80, 1.30, 1.20 }; |
---|
321 | G4int ntc = 0; |
---|
322 | |
---|
323 | if (forwardCount < 1 || backwardCount < 1) return false; // array bounds protection |
---|
324 | |
---|
325 | G4double rmc = forwardMass; |
---|
326 | if (forwardCount > 1) { |
---|
327 | ntc = std::min(3,forwardCount-2); |
---|
328 | rmc += std::pow(-std::log(1.0-G4UniformRand()),1./cpar[ntc])/gpar[ntc]; |
---|
329 | } |
---|
330 | G4double rmd = backwardMass; |
---|
331 | if( backwardCount > 1 ) { |
---|
332 | ntc = std::min(3,backwardCount-2); |
---|
333 | rmd += std::pow(-std::log(1.0-G4UniformRand()),1./cpar[ntc])/gpar[ntc]; |
---|
334 | } |
---|
335 | |
---|
336 | while( rmc+rmd > centerofmassEnergy ) |
---|
337 | { |
---|
338 | if( (rmc <= forwardMass) && (rmd <= backwardMass) ) |
---|
339 | { |
---|
340 | G4double temp = 0.999*centerofmassEnergy/(rmc+rmd); |
---|
341 | rmc *= temp; |
---|
342 | rmd *= temp; |
---|
343 | } |
---|
344 | else |
---|
345 | { |
---|
346 | rmc = 0.1*forwardMass + 0.9*rmc; |
---|
347 | rmd = 0.1*backwardMass + 0.9*rmd; |
---|
348 | } |
---|
349 | } |
---|
350 | |
---|
351 | G4ReactionProduct pseudoParticle[8]; |
---|
352 | for( i=0; i<8; ++i )pseudoParticle[i].SetZero(); |
---|
353 | |
---|
354 | pseudoParticle[1].SetMass( mOriginal*GeV ); |
---|
355 | pseudoParticle[1].SetTotalEnergy( etOriginal*GeV ); |
---|
356 | pseudoParticle[1].SetMomentum( 0.0, 0.0, pOriginal*GeV ); |
---|
357 | |
---|
358 | pseudoParticle[2].SetMass( protonMass*MeV ); |
---|
359 | pseudoParticle[2].SetTotalEnergy( protonMass*MeV ); |
---|
360 | pseudoParticle[2].SetMomentum( 0.0, 0.0, 0.0 ); |
---|
361 | // |
---|
362 | // transform into center of mass system |
---|
363 | // |
---|
364 | pseudoParticle[0] = pseudoParticle[1] + pseudoParticle[2]; |
---|
365 | pseudoParticle[1].Lorentz( pseudoParticle[1], pseudoParticle[0] ); |
---|
366 | pseudoParticle[2].Lorentz( pseudoParticle[2], pseudoParticle[0] ); |
---|
367 | |
---|
368 | // Calculate cm momentum for forward and backward masses |
---|
369 | // W = sqrt(pf*pf + rmc*rmc) + sqrt(pf*pf + rmd*rmd) |
---|
370 | // Solve for pf |
---|
371 | |
---|
372 | const G4double pfMin = 0.0001; |
---|
373 | G4double pf = (centerofmassEnergy*centerofmassEnergy+rmd*rmd-rmc*rmc); |
---|
374 | pf *= pf; |
---|
375 | pf -= 4*centerofmassEnergy*centerofmassEnergy*rmd*rmd; |
---|
376 | pf = std::sqrt( std::max(pf,pfMin) )/(2.0*centerofmassEnergy); |
---|
377 | // |
---|
378 | // set final state masses and energies in centre of mass system |
---|
379 | // |
---|
380 | pseudoParticle[3].SetMass( rmc*GeV ); |
---|
381 | pseudoParticle[3].SetTotalEnergy( std::sqrt(pf*pf+rmc*rmc)*GeV ); |
---|
382 | |
---|
383 | pseudoParticle[4].SetMass( rmd*GeV ); |
---|
384 | pseudoParticle[4].SetTotalEnergy( std::sqrt(pf*pf+rmd*rmd)*GeV ); |
---|
385 | |
---|
386 | // |
---|
387 | // Get cm scattering angle by sampling t from tmin to tmax |
---|
388 | // |
---|
389 | const G4double bMin = 0.01; |
---|
390 | const G4double b1 = 4.0; |
---|
391 | const G4double b2 = 1.6; |
---|
392 | G4double pin = pseudoParticle[1].GetMomentum().mag()/GeV; |
---|
393 | G4double dtb = 4.0*pin*pf*std::max( bMin, b1+b2*std::log(pOriginal) ); |
---|
394 | G4double factor = 1.0 - std::exp(-dtb); |
---|
395 | G4double costheta = 1.0 + 2.0*std::log(1.0 - G4UniformRand()*factor) / dtb; |
---|
396 | |
---|
397 | costheta = std::max(-1.0, std::min(1.0, costheta)); |
---|
398 | G4double sintheta = std::sqrt((1.0-costheta)*(1.0+costheta)); |
---|
399 | G4double phi = G4UniformRand() * twopi; |
---|
400 | // |
---|
401 | // calculate final state momenta in centre of mass system |
---|
402 | // |
---|
403 | pseudoParticle[3].SetMomentum( pf*sintheta*std::cos(phi)*GeV, |
---|
404 | pf*sintheta*std::sin(phi)*GeV, |
---|
405 | pf*costheta*GeV ); |
---|
406 | pseudoParticle[4].SetMomentum( -pseudoParticle[3].GetMomentum()); |
---|
407 | |
---|
408 | // Backward cluster of nucleons and pions from intra-nuclear cascade |
---|
409 | // Set up in lab system and transform to cms |
---|
410 | |
---|
411 | G4double pp, pp1; |
---|
412 | if( nuclearExcitationCount > 0 ) |
---|
413 | { |
---|
414 | const G4double ga = 1.2; |
---|
415 | G4double ekit1 = 0.04; |
---|
416 | G4double ekit2 = 0.6; // Max KE of cascade particle |
---|
417 | if( ekOriginal <= 5.0 ) |
---|
418 | { |
---|
419 | ekit1 *= ekOriginal*ekOriginal/25.0; |
---|
420 | ekit2 *= ekOriginal*ekOriginal/25.0; |
---|
421 | } |
---|
422 | G4double scale = std::pow(ekit2/ekit1, 1.0-ga) - 1.0; |
---|
423 | for( i=0; i<vecLen; ++i ) |
---|
424 | { |
---|
425 | if( vec[i]->GetSide() == -2 ) |
---|
426 | { |
---|
427 | G4double kineticE = ekit1*std::pow((1.0 + G4UniformRand()*scale), 1.0/(1.0-ga) ); |
---|
428 | vec[i]->SetKineticEnergy( kineticE*GeV ); |
---|
429 | G4double vMass = vec[i]->GetMass()/MeV; |
---|
430 | G4double totalE = kineticE*GeV + vMass; |
---|
431 | pp = std::sqrt( std::abs(totalE*totalE-vMass*vMass) ); |
---|
432 | G4double cost = std::min( 1.0, std::max( -1.0, std::log(2.23*G4UniformRand()+0.383)/0.96 ) ); |
---|
433 | G4double sint = std::sqrt(1.0-cost*cost); |
---|
434 | phi = twopi*G4UniformRand(); |
---|
435 | vec[i]->SetMomentum( pp*sint*std::cos(phi)*MeV, |
---|
436 | pp*sint*std::sin(phi)*MeV, |
---|
437 | pp*cost*MeV ); |
---|
438 | vec[i]->Lorentz( *vec[i], pseudoParticle[0] ); |
---|
439 | } |
---|
440 | } |
---|
441 | } |
---|
442 | |
---|
443 | // |
---|
444 | // Fragmentation of forward and backward clusters |
---|
445 | // |
---|
446 | |
---|
447 | currentParticle.SetMomentum( pseudoParticle[3].GetMomentum() ); |
---|
448 | currentParticle.SetTotalEnergy( pseudoParticle[3].GetTotalEnergy() ); |
---|
449 | |
---|
450 | targetParticle.SetMomentum( pseudoParticle[4].GetMomentum() ); |
---|
451 | targetParticle.SetTotalEnergy( pseudoParticle[4].GetTotalEnergy() ); |
---|
452 | |
---|
453 | pseudoParticle[5].SetMomentum( pseudoParticle[3].GetMomentum() * (-1.0) ); |
---|
454 | pseudoParticle[5].SetMass( pseudoParticle[3].GetMass() ); |
---|
455 | pseudoParticle[5].SetTotalEnergy( pseudoParticle[3].GetTotalEnergy() ); |
---|
456 | |
---|
457 | pseudoParticle[6].SetMomentum( pseudoParticle[4].GetMomentum() * (-1.0) ); |
---|
458 | pseudoParticle[6].SetMass( pseudoParticle[4].GetMass() ); |
---|
459 | pseudoParticle[6].SetTotalEnergy( pseudoParticle[4].GetTotalEnergy() ); |
---|
460 | |
---|
461 | G4double wgt; |
---|
462 | // DEBUGGING --> DumpFrames::DumpFrame(vec, vecLen); |
---|
463 | if( forwardCount > 1 ) // tempV will contain the forward particles |
---|
464 | { |
---|
465 | G4FastVector<G4ReactionProduct,256> tempV; |
---|
466 | tempV.Initialize( forwardCount ); |
---|
467 | G4bool constantCrossSection = true; |
---|
468 | G4int tempLen = 0; |
---|
469 | if( currentParticle.GetSide() == 1 ) |
---|
470 | tempV.SetElement( tempLen++, ¤tParticle ); |
---|
471 | if( targetParticle.GetSide() == 1 ) |
---|
472 | tempV.SetElement( tempLen++, &targetParticle ); |
---|
473 | for( i=0; i<vecLen; ++i ) |
---|
474 | { |
---|
475 | if( vec[i]->GetSide() == 1 ) |
---|
476 | { |
---|
477 | if( tempLen < 18 ) |
---|
478 | tempV.SetElement( tempLen++, vec[i] ); |
---|
479 | else |
---|
480 | { |
---|
481 | vec[i]->SetSide( -1 ); |
---|
482 | continue; |
---|
483 | } |
---|
484 | } |
---|
485 | } |
---|
486 | if( tempLen >= 2 ) |
---|
487 | { |
---|
488 | wgt = GenerateNBodyEvent( pseudoParticle[3].GetMass()/MeV, |
---|
489 | constantCrossSection, tempV, tempLen ); |
---|
490 | if( currentParticle.GetSide() == 1 ) |
---|
491 | currentParticle.Lorentz( currentParticle, pseudoParticle[5] ); |
---|
492 | if( targetParticle.GetSide() == 1 ) |
---|
493 | targetParticle.Lorentz( targetParticle, pseudoParticle[5] ); |
---|
494 | for( i=0; i<vecLen; ++i ) |
---|
495 | { |
---|
496 | if( vec[i]->GetSide() == 1 )vec[i]->Lorentz( *vec[i], pseudoParticle[5] ); |
---|
497 | } |
---|
498 | } |
---|
499 | } |
---|
500 | // DEBUGGING --> DumpFrames::DumpFrame(vec, vecLen); |
---|
501 | if( backwardCount > 1 ) // tempV will contain the backward particles, |
---|
502 | { // but not those created from the intranuclear cascade |
---|
503 | G4FastVector<G4ReactionProduct,256> tempV; |
---|
504 | tempV.Initialize( backwardCount ); |
---|
505 | G4bool constantCrossSection = true; |
---|
506 | G4int tempLen = 0; |
---|
507 | if( currentParticle.GetSide() == -1 ) |
---|
508 | tempV.SetElement( tempLen++, ¤tParticle ); |
---|
509 | if( targetParticle.GetSide() == -1 ) |
---|
510 | tempV.SetElement( tempLen++, &targetParticle ); |
---|
511 | for( i=0; i<vecLen; ++i ) |
---|
512 | { |
---|
513 | if( vec[i]->GetSide() == -1 ) |
---|
514 | { |
---|
515 | if( tempLen < 18 ) |
---|
516 | tempV.SetElement( tempLen++, vec[i] ); |
---|
517 | else |
---|
518 | { |
---|
519 | vec[i]->SetSide( -2 ); |
---|
520 | vec[i]->SetKineticEnergy( 0.0 ); |
---|
521 | vec[i]->SetMomentum( 0.0, 0.0, 0.0 ); |
---|
522 | continue; |
---|
523 | } |
---|
524 | } |
---|
525 | } |
---|
526 | if( tempLen >= 2 ) |
---|
527 | { |
---|
528 | wgt = GenerateNBodyEvent( pseudoParticle[4].GetMass()/MeV, |
---|
529 | constantCrossSection, tempV, tempLen ); |
---|
530 | if( currentParticle.GetSide() == -1 ) |
---|
531 | currentParticle.Lorentz( currentParticle, pseudoParticle[6] ); |
---|
532 | if( targetParticle.GetSide() == -1 ) |
---|
533 | targetParticle.Lorentz( targetParticle, pseudoParticle[6] ); |
---|
534 | for( i=0; i<vecLen; ++i ) |
---|
535 | { |
---|
536 | if( vec[i]->GetSide() == -1 )vec[i]->Lorentz( *vec[i], pseudoParticle[6] ); |
---|
537 | } |
---|
538 | } |
---|
539 | } |
---|
540 | |
---|
541 | // DEBUGGING --> DumpFrames::DumpFrame(vec, vecLen); |
---|
542 | // |
---|
543 | // Lorentz transformation in lab system |
---|
544 | // |
---|
545 | currentParticle.Lorentz( currentParticle, pseudoParticle[2] ); |
---|
546 | targetParticle.Lorentz( targetParticle, pseudoParticle[2] ); |
---|
547 | for( i=0; i<vecLen; ++i ) vec[i]->Lorentz( *vec[i], pseudoParticle[2] ); |
---|
548 | |
---|
549 | // DEBUGGING --> DumpFrames::DumpFrame(vec, vecLen); |
---|
550 | // |
---|
551 | // sometimes the leading strange particle is lost, set it back |
---|
552 | // |
---|
553 | G4bool dum = true; |
---|
554 | if( leadFlag ) |
---|
555 | { |
---|
556 | // leadFlag will be true |
---|
557 | // iff original particle is strange AND if incident particle is strange |
---|
558 | // leadFlag is set to the incident particle |
---|
559 | // or |
---|
560 | // target particle is strange leadFlag is set to the target particle |
---|
561 | |
---|
562 | if( currentParticle.GetDefinition() == leadingStrangeParticle.GetDefinition() ) |
---|
563 | dum = false; |
---|
564 | else if( targetParticle.GetDefinition() == leadingStrangeParticle.GetDefinition() ) |
---|
565 | dum = false; |
---|
566 | else |
---|
567 | { |
---|
568 | for( i=0; i<vecLen; ++i ) |
---|
569 | { |
---|
570 | if( vec[i]->GetDefinition() == leadingStrangeParticle.GetDefinition() ) |
---|
571 | { |
---|
572 | dum = false; |
---|
573 | break; |
---|
574 | } |
---|
575 | } |
---|
576 | } |
---|
577 | if( dum ) |
---|
578 | { |
---|
579 | G4double leadMass = leadingStrangeParticle.GetMass()/MeV; |
---|
580 | G4double ekin; |
---|
581 | if( ((leadMass < protonMass) && (targetParticle.GetMass()/MeV < protonMass)) || |
---|
582 | ((leadMass >= protonMass) && (targetParticle.GetMass()/MeV >= protonMass)) ) |
---|
583 | { |
---|
584 | ekin = targetParticle.GetKineticEnergy()/GeV; |
---|
585 | pp1 = targetParticle.GetMomentum().mag()/MeV; // old momentum |
---|
586 | targetParticle.SetDefinition( leadingStrangeParticle.GetDefinition() ); |
---|
587 | targetParticle.SetKineticEnergy( ekin*GeV ); |
---|
588 | pp = targetParticle.GetTotalMomentum()/MeV; // new momentum |
---|
589 | if( pp1 < 1.0e-3 ) { |
---|
590 | G4ThreeVector iso = Isotropic(pp); |
---|
591 | targetParticle.SetMomentum( iso.x(), iso.y(), iso.z() ); |
---|
592 | } else { |
---|
593 | targetParticle.SetMomentum( targetParticle.GetMomentum() * (pp/pp1) ); |
---|
594 | } |
---|
595 | targetHasChanged = true; |
---|
596 | } |
---|
597 | else |
---|
598 | { |
---|
599 | ekin = currentParticle.GetKineticEnergy()/GeV; |
---|
600 | pp1 = currentParticle.GetMomentum().mag()/MeV; |
---|
601 | currentParticle.SetDefinition( leadingStrangeParticle.GetDefinition() ); |
---|
602 | currentParticle.SetKineticEnergy( ekin*GeV ); |
---|
603 | pp = currentParticle.GetTotalMomentum()/MeV; |
---|
604 | if( pp1 < 1.0e-3 ) { |
---|
605 | G4ThreeVector iso = Isotropic(pp); |
---|
606 | currentParticle.SetMomentum( iso.x(), iso.y(), iso.z() ); |
---|
607 | } else { |
---|
608 | currentParticle.SetMomentum( currentParticle.GetMomentum() * (pp/pp1) ); |
---|
609 | } |
---|
610 | incidentHasChanged = true; |
---|
611 | } |
---|
612 | } |
---|
613 | } // end of if( leadFlag ) |
---|
614 | |
---|
615 | // Get number of final state nucleons and nucleons remaining in |
---|
616 | // target nucleus |
---|
617 | |
---|
618 | std::pair<G4int, G4int> finalStateNucleons = |
---|
619 | GetFinalStateNucleons(originalTarget, vec, vecLen); |
---|
620 | |
---|
621 | G4int protonsInFinalState = finalStateNucleons.first; |
---|
622 | G4int neutronsInFinalState = finalStateNucleons.second; |
---|
623 | |
---|
624 | G4int numberofFinalStateNucleons = |
---|
625 | protonsInFinalState + neutronsInFinalState; |
---|
626 | |
---|
627 | if (currentParticle.GetDefinition()->GetBaryonNumber() == 1 && |
---|
628 | targetParticle.GetDefinition()->GetBaryonNumber() == 1 && |
---|
629 | originalIncident->GetDefinition()->GetPDGMass() < |
---|
630 | G4Lambda::Lambda()->GetPDGMass()) |
---|
631 | numberofFinalStateNucleons++; |
---|
632 | |
---|
633 | numberofFinalStateNucleons = std::max(1, numberofFinalStateNucleons); |
---|
634 | |
---|
635 | G4int PinNucleus = std::max(0, |
---|
636 | G4int(targetNucleus.GetZ()) - protonsInFinalState); |
---|
637 | G4int NinNucleus = std::max(0, |
---|
638 | G4int(targetNucleus.GetN()-targetNucleus.GetZ()) - neutronsInFinalState); |
---|
639 | // |
---|
640 | // for various reasons, the energy balance is not sufficient, |
---|
641 | // check that, energy balance, angle of final system, etc. |
---|
642 | // |
---|
643 | pseudoParticle[4].SetMass( mOriginal*GeV ); |
---|
644 | pseudoParticle[4].SetTotalEnergy( etOriginal*GeV ); |
---|
645 | pseudoParticle[4].SetMomentum( 0.0, 0.0, pOriginal*GeV ); |
---|
646 | |
---|
647 | G4ParticleDefinition * aOrgDef = modifiedOriginal.GetDefinition(); |
---|
648 | G4int diff = 0; |
---|
649 | if(aOrgDef == G4Proton::Proton() || aOrgDef == G4Neutron::Neutron() ) diff = 1; |
---|
650 | if(numberofFinalStateNucleons == 1) diff = 0; |
---|
651 | pseudoParticle[5].SetMomentum( 0.0, 0.0, 0.0 ); |
---|
652 | pseudoParticle[5].SetMass( protonMass*(numberofFinalStateNucleons-diff)*MeV); |
---|
653 | pseudoParticle[5].SetTotalEnergy( protonMass*(numberofFinalStateNucleons-diff)*MeV); |
---|
654 | |
---|
655 | G4double theoreticalKinetic = |
---|
656 | pseudoParticle[4].GetTotalEnergy()/GeV + pseudoParticle[5].GetTotalEnergy()/GeV; |
---|
657 | |
---|
658 | pseudoParticle[6] = pseudoParticle[4] + pseudoParticle[5]; |
---|
659 | pseudoParticle[4].Lorentz( pseudoParticle[4], pseudoParticle[6] ); |
---|
660 | pseudoParticle[5].Lorentz( pseudoParticle[5], pseudoParticle[6] ); |
---|
661 | |
---|
662 | if( vecLen < 16 ) |
---|
663 | { |
---|
664 | G4ReactionProduct tempR[130]; |
---|
665 | tempR[0] = currentParticle; |
---|
666 | tempR[1] = targetParticle; |
---|
667 | for( i=0; i<vecLen; ++i )tempR[i+2] = *vec[i]; |
---|
668 | |
---|
669 | G4FastVector<G4ReactionProduct,256> tempV; |
---|
670 | tempV.Initialize( vecLen+2 ); |
---|
671 | G4bool constantCrossSection = true; |
---|
672 | G4int tempLen = 0; |
---|
673 | for( i=0; i<vecLen+2; ++i )tempV.SetElement( tempLen++, &tempR[i] ); |
---|
674 | |
---|
675 | if( tempLen >= 2 ) |
---|
676 | { |
---|
677 | // DEBUGGING --> DumpFrames::DumpFrame(vec, vecLen); |
---|
678 | wgt = GenerateNBodyEvent( pseudoParticle[4].GetTotalEnergy()/MeV + |
---|
679 | pseudoParticle[5].GetTotalEnergy()/MeV, |
---|
680 | constantCrossSection, tempV, tempLen ); |
---|
681 | if (wgt == -1) { |
---|
682 | G4double Qvalue = 0; |
---|
683 | for (i = 0; i < tempLen; i++) Qvalue += tempV[i]->GetMass(); |
---|
684 | wgt = GenerateNBodyEvent( Qvalue/MeV, |
---|
685 | constantCrossSection, tempV, tempLen ); |
---|
686 | } |
---|
687 | theoreticalKinetic = 0.0; |
---|
688 | for( i=0; i<vecLen+2; ++i ) |
---|
689 | { |
---|
690 | pseudoParticle[7].SetMomentum( tempV[i]->GetMomentum() ); |
---|
691 | pseudoParticle[7].SetMass( tempV[i]->GetMass() ); |
---|
692 | pseudoParticle[7].SetTotalEnergy( tempV[i]->GetTotalEnergy() ); |
---|
693 | pseudoParticle[7].Lorentz( pseudoParticle[7], pseudoParticle[5] ); |
---|
694 | theoreticalKinetic += pseudoParticle[7].GetKineticEnergy()/GeV; |
---|
695 | } |
---|
696 | } |
---|
697 | // DEBUGGING --> DumpFrames::DumpFrame(vec, vecLen); |
---|
698 | } |
---|
699 | else |
---|
700 | { |
---|
701 | theoreticalKinetic -= |
---|
702 | ( currentParticle.GetMass()/GeV + targetParticle.GetMass()/GeV ); |
---|
703 | for( i=0; i<vecLen; ++i )theoreticalKinetic -= vec[i]->GetMass()/GeV; |
---|
704 | } |
---|
705 | G4double simulatedKinetic = |
---|
706 | currentParticle.GetKineticEnergy()/GeV + targetParticle.GetKineticEnergy()/GeV; |
---|
707 | for( i=0; i<vecLen; ++i )simulatedKinetic += vec[i]->GetKineticEnergy()/GeV; |
---|
708 | |
---|
709 | // make sure that kinetic energies are correct |
---|
710 | // the backward nucleon cluster is not produced within proper kinematics!!! |
---|
711 | |
---|
712 | if( simulatedKinetic != 0.0 ) |
---|
713 | { |
---|
714 | wgt = (theoreticalKinetic)/simulatedKinetic; |
---|
715 | currentParticle.SetKineticEnergy( wgt*currentParticle.GetKineticEnergy() ); |
---|
716 | pp = currentParticle.GetTotalMomentum()/MeV; |
---|
717 | pp1 = currentParticle.GetMomentum().mag()/MeV; |
---|
718 | if( pp1 < 0.001*MeV ) { |
---|
719 | G4ThreeVector iso = Isotropic(pp); |
---|
720 | currentParticle.SetMomentum( iso.x(), iso.y(), iso.z() ); |
---|
721 | } else { |
---|
722 | currentParticle.SetMomentum( currentParticle.GetMomentum() * (pp/pp1) ); |
---|
723 | } |
---|
724 | |
---|
725 | targetParticle.SetKineticEnergy( wgt*targetParticle.GetKineticEnergy() ); |
---|
726 | pp = targetParticle.GetTotalMomentum()/MeV; |
---|
727 | pp1 = targetParticle.GetMomentum().mag()/MeV; |
---|
728 | if( pp1 < 0.001*MeV ) { |
---|
729 | G4ThreeVector iso = Isotropic(pp); |
---|
730 | targetParticle.SetMomentum( iso.x(), iso.y(), iso.z() ); |
---|
731 | } else { |
---|
732 | targetParticle.SetMomentum( targetParticle.GetMomentum() * (pp/pp1) ); |
---|
733 | } |
---|
734 | |
---|
735 | for( i=0; i<vecLen; ++i ) |
---|
736 | { |
---|
737 | vec[i]->SetKineticEnergy( wgt*vec[i]->GetKineticEnergy() ); |
---|
738 | pp = vec[i]->GetTotalMomentum()/MeV; |
---|
739 | pp1 = vec[i]->GetMomentum().mag()/MeV; |
---|
740 | if( pp1 < 0.001 ) { |
---|
741 | G4ThreeVector iso = Isotropic(pp); |
---|
742 | vec[i]->SetMomentum( iso.x(), iso.y(), iso.z() ); |
---|
743 | } else { |
---|
744 | vec[i]->SetMomentum( vec[i]->GetMomentum() * (pp/pp1) ); |
---|
745 | } |
---|
746 | } |
---|
747 | } |
---|
748 | // DEBUGGING --> DumpFrames::DumpFrame(vec, vecLen); |
---|
749 | |
---|
750 | Rotate( numberofFinalStateNucleons, pseudoParticle[4].GetMomentum(), |
---|
751 | modifiedOriginal, originalIncident, targetNucleus, |
---|
752 | currentParticle, targetParticle, vec, vecLen ); |
---|
753 | |
---|
754 | // Add black track particles |
---|
755 | // the total number of particles produced is restricted to 198 |
---|
756 | // this may have influence on very high energies |
---|
757 | |
---|
758 | if( atomicWeight >= 1.5 ) |
---|
759 | { |
---|
760 | // npnb is number of proton/neutron black track particles |
---|
761 | // ndta is the number of deuterons, tritons, and alphas produced |
---|
762 | // epnb is the kinetic energy available for proton/neutron black track |
---|
763 | // particles |
---|
764 | // edta is the kinetic energy available for deuteron/triton/alpha |
---|
765 | // particles |
---|
766 | |
---|
767 | G4int npnb = 0; |
---|
768 | G4int ndta = 0; |
---|
769 | |
---|
770 | G4double epnb, edta; |
---|
771 | if (veryForward) { |
---|
772 | epnb = targetNucleus.GetAnnihilationPNBlackTrackEnergy(); |
---|
773 | edta = targetNucleus.GetAnnihilationDTABlackTrackEnergy(); |
---|
774 | } else { |
---|
775 | epnb = targetNucleus.GetPNBlackTrackEnergy(); |
---|
776 | edta = targetNucleus.GetDTABlackTrackEnergy(); |
---|
777 | } |
---|
778 | |
---|
779 | const G4double pnCutOff = 0.001; // GeV |
---|
780 | const G4double dtaCutOff = 0.001; // GeV |
---|
781 | const G4double kineticMinimum = 1.e-6; |
---|
782 | const G4double kineticFactor = -0.005; |
---|
783 | |
---|
784 | G4double sprob = 0.0; // sprob = probability of self-absorption in |
---|
785 | // heavy molecules |
---|
786 | const G4double ekIncident = originalIncident->GetKineticEnergy()/GeV; |
---|
787 | if( ekIncident >= 5.0 )sprob = std::min( 1.0, 0.6*std::log(ekIncident-4.0) ); |
---|
788 | |
---|
789 | if( epnb >= pnCutOff ) |
---|
790 | { |
---|
791 | npnb = G4Poisson((1.5+1.25*numberofFinalStateNucleons)*epnb/(epnb+edta)); |
---|
792 | if( numberofFinalStateNucleons + npnb > atomicWeight ) |
---|
793 | npnb = G4int(atomicWeight - numberofFinalStateNucleons); |
---|
794 | npnb = std::min( npnb, 127-vecLen ); |
---|
795 | } |
---|
796 | if( edta >= dtaCutOff ) |
---|
797 | { |
---|
798 | ndta = G4Poisson( (1.5+1.25*numberofFinalStateNucleons)*edta/(epnb+edta) ); |
---|
799 | ndta = std::min( ndta, 127-vecLen ); |
---|
800 | } |
---|
801 | if (npnb == 0 && ndta == 0) npnb = 1; |
---|
802 | |
---|
803 | // DEBUGGING --> DumpFrames::DumpFrame(vec, vecLen); |
---|
804 | |
---|
805 | AddBlackTrackParticles(epnb, npnb, edta, ndta, sprob, kineticMinimum, |
---|
806 | kineticFactor, modifiedOriginal, |
---|
807 | PinNucleus, NinNucleus, targetNucleus, |
---|
808 | vec, vecLen ); |
---|
809 | // DEBUGGING --> DumpFrames::DumpFrame(vec, vecLen); |
---|
810 | } |
---|
811 | |
---|
812 | //if( centerofmassEnergy <= (4.0+G4UniformRand()) ) |
---|
813 | // MomentumCheck( modifiedOriginal, currentParticle, targetParticle, vec, vecLen ); |
---|
814 | // |
---|
815 | // calculate time delay for nuclear reactions |
---|
816 | // |
---|
817 | if( (atomicWeight >= 1.5) && (atomicWeight <= 230.0) && (ekOriginal <= 0.2) ) |
---|
818 | currentParticle.SetTOF( 1.0-500.0*std::exp(-ekOriginal/0.04)*std::log(G4UniformRand()) ); |
---|
819 | else |
---|
820 | currentParticle.SetTOF( 1.0 ); |
---|
821 | |
---|
822 | return true; |
---|
823 | } |
---|
824 | |
---|
825 | /* end of file */ |
---|