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 | // |
---|
27 | // $Id: G4LEKaonZeroInelastic.cc,v 1.9 2006/06/29 20:45:03 gunter Exp $ |
---|
28 | // GEANT4 tag $Name: geant4-09-03-ref-09 $ |
---|
29 | // |
---|
30 | // Hadronic Process: Low Energy KaonZeroShort Inelastic Process |
---|
31 | // J.L. Chuma, TRIUMF, 11-Feb-1997 |
---|
32 | // Last modified: 27-Mar-1997 |
---|
33 | // Modified by J.L.Chuma 30-Apr-97: added originalTarget for CalculateMomenta |
---|
34 | |
---|
35 | #include "G4LEKaonZeroInelastic.hh" |
---|
36 | #include "Randomize.hh" |
---|
37 | |
---|
38 | G4HadFinalState * |
---|
39 | G4LEKaonZeroInelastic::ApplyYourself( const G4HadProjectile &aTrack, |
---|
40 | G4Nucleus &targetNucleus ) |
---|
41 | { |
---|
42 | const G4HadProjectile *originalIncident = &aTrack; |
---|
43 | // |
---|
44 | // create the target particle |
---|
45 | // |
---|
46 | G4DynamicParticle *originalTarget = targetNucleus.ReturnTargetParticle(); |
---|
47 | |
---|
48 | if( verboseLevel > 1 ) |
---|
49 | { |
---|
50 | const G4Material *targetMaterial = aTrack.GetMaterial(); |
---|
51 | G4cout << "G4LEKaonZeroInelastic::ApplyYourself called" << G4endl; |
---|
52 | G4cout << "kinetic energy = " << originalIncident->GetKineticEnergy()/MeV << "MeV, "; |
---|
53 | G4cout << "target material = " << targetMaterial->GetName() << ", "; |
---|
54 | G4cout << "target particle = " << originalTarget->GetDefinition()->GetParticleName() |
---|
55 | << G4endl; |
---|
56 | } |
---|
57 | // |
---|
58 | // Fermi motion and evaporation |
---|
59 | // As of Geant3, the Fermi energy calculation had not been Done |
---|
60 | // |
---|
61 | G4double ek = originalIncident->GetKineticEnergy()/MeV; |
---|
62 | G4double amas = originalIncident->GetDefinition()->GetPDGMass()/MeV; |
---|
63 | G4ReactionProduct modifiedOriginal; |
---|
64 | modifiedOriginal = *originalIncident; |
---|
65 | |
---|
66 | G4double tkin = targetNucleus.Cinema( ek ); |
---|
67 | ek += tkin; |
---|
68 | modifiedOriginal.SetKineticEnergy( ek*MeV ); |
---|
69 | G4double et = ek + amas; |
---|
70 | G4double p = std::sqrt( std::abs((et-amas)*(et+amas)) ); |
---|
71 | G4double pp = modifiedOriginal.GetMomentum().mag()/MeV; |
---|
72 | if( pp > 0.0 ) |
---|
73 | { |
---|
74 | G4ThreeVector momentum = modifiedOriginal.GetMomentum(); |
---|
75 | modifiedOriginal.SetMomentum( momentum * (p/pp) ); |
---|
76 | } |
---|
77 | // |
---|
78 | // calculate black track energies |
---|
79 | // |
---|
80 | tkin = targetNucleus.EvaporationEffects( ek ); |
---|
81 | ek -= tkin; |
---|
82 | modifiedOriginal.SetKineticEnergy( ek*MeV ); |
---|
83 | et = ek + amas; |
---|
84 | p = std::sqrt( std::abs((et-amas)*(et+amas)) ); |
---|
85 | pp = modifiedOriginal.GetMomentum().mag()/MeV; |
---|
86 | if( pp > 0.0 ) |
---|
87 | { |
---|
88 | G4ThreeVector momentum = modifiedOriginal.GetMomentum(); |
---|
89 | modifiedOriginal.SetMomentum( momentum * (p/pp) ); |
---|
90 | } |
---|
91 | G4ReactionProduct currentParticle = modifiedOriginal; |
---|
92 | G4ReactionProduct targetParticle; |
---|
93 | targetParticle = *originalTarget; |
---|
94 | currentParticle.SetSide( 1 ); // incident always goes in forward hemisphere |
---|
95 | targetParticle.SetSide( -1 ); // target always goes in backward hemisphere |
---|
96 | G4bool incidentHasChanged = false; |
---|
97 | G4bool targetHasChanged = false; |
---|
98 | G4bool quasiElastic = false; |
---|
99 | G4FastVector<G4ReactionProduct,GHADLISTSIZE> vec; // vec will contain the secondary particles |
---|
100 | G4int vecLen = 0; |
---|
101 | vec.Initialize( 0 ); |
---|
102 | |
---|
103 | const G4double cutOff = 0.1; |
---|
104 | if( currentParticle.GetKineticEnergy()/MeV > cutOff ) |
---|
105 | Cascade( vec, vecLen, |
---|
106 | originalIncident, currentParticle, targetParticle, |
---|
107 | incidentHasChanged, targetHasChanged, quasiElastic ); |
---|
108 | |
---|
109 | CalculateMomenta( vec, vecLen, |
---|
110 | originalIncident, originalTarget, modifiedOriginal, |
---|
111 | targetNucleus, currentParticle, targetParticle, |
---|
112 | incidentHasChanged, targetHasChanged, quasiElastic ); |
---|
113 | |
---|
114 | SetUpChange( vec, vecLen, |
---|
115 | currentParticle, targetParticle, |
---|
116 | incidentHasChanged ); |
---|
117 | |
---|
118 | delete originalTarget; |
---|
119 | return &theParticleChange; |
---|
120 | } |
---|
121 | |
---|
122 | void |
---|
123 | G4LEKaonZeroInelastic::Cascade( |
---|
124 | G4FastVector<G4ReactionProduct,GHADLISTSIZE> &vec, |
---|
125 | G4int& vecLen, |
---|
126 | const G4HadProjectile *originalIncident, |
---|
127 | G4ReactionProduct ¤tParticle, |
---|
128 | G4ReactionProduct &targetParticle, |
---|
129 | G4bool &incidentHasChanged, |
---|
130 | G4bool &targetHasChanged, |
---|
131 | G4bool &quasiElastic ) |
---|
132 | { |
---|
133 | // derived from original FORTRAN code CASK0 by H. Fesefeldt (13-Sep-1987) |
---|
134 | // |
---|
135 | // K0Short undergoes interaction with nucleon within a nucleus. Check if it is |
---|
136 | // energetically possible to produce pions/kaons. In not, assume nuclear excitation |
---|
137 | // occurs and input particle is degraded in energy. No other particles are produced. |
---|
138 | // If reaction is possible, find the correct number of pions/protons/neutrons |
---|
139 | // produced using an interpolation to multiplicity data. Replace some pions or |
---|
140 | // protons/neutrons by kaons or strange baryons according to the average |
---|
141 | // multiplicity per Inelastic reaction. |
---|
142 | // |
---|
143 | const G4double mOriginal = originalIncident->GetDefinition()->GetPDGMass()/MeV; |
---|
144 | const G4double etOriginal = originalIncident->GetTotalEnergy()/MeV; |
---|
145 | const G4double targetMass = targetParticle.GetMass()/MeV; |
---|
146 | G4double centerofmassEnergy = std::sqrt( mOriginal*mOriginal + |
---|
147 | targetMass*targetMass + |
---|
148 | 2.0*targetMass*etOriginal ); |
---|
149 | G4double availableEnergy = centerofmassEnergy-(targetMass+mOriginal); |
---|
150 | if( availableEnergy <= G4PionPlus::PionPlus()->GetPDGMass()/MeV ) |
---|
151 | { |
---|
152 | quasiElastic = true; |
---|
153 | return; |
---|
154 | } |
---|
155 | static G4bool first = true; |
---|
156 | const G4int numMul = 1200; |
---|
157 | const G4int numSec = 60; |
---|
158 | static G4double protmul[numMul], protnorm[numSec]; // proton constants |
---|
159 | static G4double neutmul[numMul], neutnorm[numSec]; // neutron constants |
---|
160 | // np = number of pi+, nm = number of pi-, nz = number of pi0 |
---|
161 | G4int counter, nt=0, np=0, nm=0, nz=0; |
---|
162 | const G4double c = 1.25; |
---|
163 | const G4double b[] = { 0.7, 0.7 }; |
---|
164 | if( first ) // compute normalization constants, this will only be Done once |
---|
165 | { |
---|
166 | first = false; |
---|
167 | G4int i; |
---|
168 | for( i=0; i<numMul; ++i )protmul[i] = 0.0; |
---|
169 | for( i=0; i<numSec; ++i )protnorm[i] = 0.0; |
---|
170 | counter = -1; |
---|
171 | for( np=0; np<(numSec/3); ++np ) |
---|
172 | { |
---|
173 | for( nm=std::max(0,np-1); nm<=(np+1); ++nm ) |
---|
174 | { |
---|
175 | for( nz=0; nz<numSec/3; ++nz ) |
---|
176 | { |
---|
177 | if( ++counter < numMul ) |
---|
178 | { |
---|
179 | nt = np+nm+nz; |
---|
180 | if( nt>0 && nt<=numSec ) |
---|
181 | { |
---|
182 | protmul[counter] = Pmltpc(np,nm,nz,nt,b[0],c); |
---|
183 | protnorm[nt-1] += protmul[counter]; |
---|
184 | } |
---|
185 | } |
---|
186 | } |
---|
187 | } |
---|
188 | } |
---|
189 | for( i=0; i<numMul; ++i )neutmul[i] = 0.0; |
---|
190 | for( i=0; i<numSec; ++i )neutnorm[i] = 0.0; |
---|
191 | counter = -1; |
---|
192 | for( np=0; np<numSec/3; ++np ) |
---|
193 | { |
---|
194 | for( nm=std::max(0,np-2); nm<=np; ++nm ) |
---|
195 | { |
---|
196 | for( nz=0; nz<numSec/3; ++nz ) |
---|
197 | { |
---|
198 | if( ++counter < numMul ) |
---|
199 | { |
---|
200 | nt = np+nm+nz; |
---|
201 | if( nt>0 && nt<=numSec ) |
---|
202 | { |
---|
203 | neutmul[counter] = Pmltpc(np,nm,nz,nt,b[1],c); |
---|
204 | neutnorm[nt-1] += neutmul[counter]; |
---|
205 | } |
---|
206 | } |
---|
207 | } |
---|
208 | } |
---|
209 | } |
---|
210 | for( i=0; i<numSec; ++i ) |
---|
211 | { |
---|
212 | if( protnorm[i] > 0.0 )protnorm[i] = 1.0/protnorm[i]; |
---|
213 | if( neutnorm[i] > 0.0 )neutnorm[i] = 1.0/neutnorm[i]; |
---|
214 | } |
---|
215 | } // end of initialization |
---|
216 | |
---|
217 | const G4double expxu = 82.; // upper bound for arg. of exp |
---|
218 | const G4double expxl = -expxu; // lower bound for arg. of exp |
---|
219 | G4ParticleDefinition *aKaonPlus = G4KaonPlus::KaonPlus(); |
---|
220 | G4ParticleDefinition *aKaonZL = G4KaonZeroLong::KaonZeroLong(); |
---|
221 | G4ParticleDefinition *aKaonZS = G4KaonZeroShort::KaonZeroShort(); |
---|
222 | G4ParticleDefinition *aNeutron = G4Neutron::Neutron(); |
---|
223 | G4ParticleDefinition *aProton = G4Proton::Proton(); |
---|
224 | G4int ieab = static_cast<G4int>(5.0*availableEnergy*MeV/GeV); |
---|
225 | const G4double supp[] = {0.,0.4,0.55,0.65,0.75,0.82,0.86,0.90,0.94,0.98}; |
---|
226 | G4double test, w0, wp, wt, wm; |
---|
227 | if( (availableEnergy*MeV/GeV < 2.0) && (G4UniformRand() >= supp[ieab]) ) |
---|
228 | { |
---|
229 | // |
---|
230 | // suppress high multiplicity events at low momentum |
---|
231 | // only one pion will be produced |
---|
232 | // |
---|
233 | nm = np = nz = 0; |
---|
234 | if( targetParticle.GetDefinition() == aNeutron ) |
---|
235 | { |
---|
236 | test = std::exp( std::min( expxu, std::max( expxl, -(1.0+b[0])*(1.0+b[0])/(2.0*c*c) ) ) ); |
---|
237 | w0 = test/2.0; |
---|
238 | test = std::exp( std::min( expxu, std::max( expxl, -(-1.0+b[0])*(1.0+b[0])/(2.0*c*c) ) ) ); |
---|
239 | wm = test*1.5; |
---|
240 | if( G4UniformRand() < w0/(w0+wm) ) |
---|
241 | nz = 1; |
---|
242 | else |
---|
243 | nm = 1; |
---|
244 | } |
---|
245 | else // target is a proton |
---|
246 | { |
---|
247 | test = std::exp( std::min( expxu, std::max( expxl, -(1.0+b[1])*(1.0+b[1])/(2.0*c*c) ) ) ); |
---|
248 | w0 = test; |
---|
249 | wp = test; |
---|
250 | test = std::exp( std::min( expxu, std::max( expxl, -(-1.0+b[1])*(-1.0+b[1])/(2.0*c*c) ) ) ); |
---|
251 | wm = test; |
---|
252 | wt = w0+wp+wm; |
---|
253 | wp += w0; |
---|
254 | G4double ran = G4UniformRand(); |
---|
255 | if( ran < w0/wt ) |
---|
256 | nz = 1; |
---|
257 | else if( ran < wp/wt ) |
---|
258 | np = 1; |
---|
259 | else |
---|
260 | nm = 1; |
---|
261 | } |
---|
262 | } |
---|
263 | else // (availableEnergy*MeV/GeV >= 2.0) || (G4UniformRand() < supp[ieab]) |
---|
264 | { |
---|
265 | G4double n, anpn; |
---|
266 | GetNormalizationConstant( availableEnergy, n, anpn ); |
---|
267 | G4double ran = G4UniformRand(); |
---|
268 | G4double dum, excs = 0.0; |
---|
269 | if( targetParticle.GetDefinition() == aProton ) |
---|
270 | { |
---|
271 | counter = -1; |
---|
272 | for( np=0; np<numSec/3 && ran>=excs; ++np ) |
---|
273 | { |
---|
274 | for( nm=std::max(0,np-1); nm<=(np+1) && ran>=excs; ++nm ) |
---|
275 | { |
---|
276 | for( nz=0; nz<numSec/3 && ran>=excs; ++nz ) |
---|
277 | { |
---|
278 | if( ++counter < numMul ) |
---|
279 | { |
---|
280 | nt = np+nm+nz; |
---|
281 | if( nt>0 && nt<=numSec ) |
---|
282 | { |
---|
283 | test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) ); |
---|
284 | dum = (pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n); |
---|
285 | if( std::fabs(dum) < 1.0 ) |
---|
286 | { |
---|
287 | if( test >= 1.0e-10 )excs += dum*test; |
---|
288 | } |
---|
289 | else |
---|
290 | excs += dum*test; |
---|
291 | } |
---|
292 | } |
---|
293 | } |
---|
294 | } |
---|
295 | } |
---|
296 | if( ran >= excs ) // 3 previous loops continued to the end |
---|
297 | { |
---|
298 | quasiElastic = true; |
---|
299 | return; |
---|
300 | } |
---|
301 | np--; nm--; nz--; |
---|
302 | } |
---|
303 | else // target must be a neutron |
---|
304 | { |
---|
305 | counter = -1; |
---|
306 | for( np=0; np<numSec/3 && ran>=excs; ++np ) |
---|
307 | { |
---|
308 | for( nm=std::max(0,np-2); nm<=np && ran>=excs; ++nm ) |
---|
309 | { |
---|
310 | for( nz=0; nz<numSec/3 && ran>=excs; ++nz ) |
---|
311 | { |
---|
312 | if( ++counter < numMul ) |
---|
313 | { |
---|
314 | nt = np+nm+nz; |
---|
315 | if( nt>0 && nt<=numSec ) |
---|
316 | { |
---|
317 | test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) ); |
---|
318 | dum = (pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n); |
---|
319 | if( std::fabs(dum) < 1.0 ) |
---|
320 | { |
---|
321 | if( test >= 1.0e-10 )excs += dum*test; |
---|
322 | } |
---|
323 | else |
---|
324 | excs += dum*test; |
---|
325 | } |
---|
326 | } |
---|
327 | } |
---|
328 | } |
---|
329 | } |
---|
330 | if( ran >= excs ) // 3 previous loops continued to the end |
---|
331 | { |
---|
332 | quasiElastic = true; |
---|
333 | return; |
---|
334 | } |
---|
335 | np--; nm--; nz--; |
---|
336 | } |
---|
337 | } |
---|
338 | if( targetParticle.GetDefinition() == aProton ) |
---|
339 | { |
---|
340 | switch( np-nm ) |
---|
341 | { |
---|
342 | case 0: |
---|
343 | if( G4UniformRand() < 0.25 ) |
---|
344 | { |
---|
345 | currentParticle.SetDefinitionAndUpdateE( aKaonPlus ); |
---|
346 | targetParticle.SetDefinitionAndUpdateE( aNeutron ); |
---|
347 | incidentHasChanged = true; |
---|
348 | targetHasChanged = true; |
---|
349 | } |
---|
350 | break; |
---|
351 | case 1: |
---|
352 | targetParticle.SetDefinitionAndUpdateE( aNeutron ); |
---|
353 | targetHasChanged = true; |
---|
354 | break; |
---|
355 | default: |
---|
356 | targetParticle.SetDefinitionAndUpdateE( aNeutron ); |
---|
357 | targetHasChanged = true; |
---|
358 | break; |
---|
359 | } |
---|
360 | } |
---|
361 | else // targetParticle is a neutron |
---|
362 | { |
---|
363 | switch( np-nm ) // seems wrong, charge not conserved |
---|
364 | { |
---|
365 | case 1: |
---|
366 | if( G4UniformRand() < 0.5 ) |
---|
367 | { |
---|
368 | currentParticle.SetDefinitionAndUpdateE( aKaonPlus ); |
---|
369 | incidentHasChanged = true; |
---|
370 | } |
---|
371 | else |
---|
372 | { |
---|
373 | targetParticle.SetDefinitionAndUpdateE( aProton ); |
---|
374 | targetHasChanged = true; |
---|
375 | } |
---|
376 | break; |
---|
377 | case 2: |
---|
378 | currentParticle.SetDefinitionAndUpdateE( aKaonPlus ); |
---|
379 | incidentHasChanged = true; |
---|
380 | targetParticle.SetDefinitionAndUpdateE( aProton ); |
---|
381 | targetHasChanged = true; |
---|
382 | break; |
---|
383 | default: |
---|
384 | break; |
---|
385 | } |
---|
386 | } |
---|
387 | if( currentParticle.GetDefinition() == aKaonZS ) |
---|
388 | { |
---|
389 | if( G4UniformRand() >= 0.5 ) |
---|
390 | { |
---|
391 | currentParticle.SetDefinitionAndUpdateE( aKaonZL); |
---|
392 | incidentHasChanged = true; |
---|
393 | } |
---|
394 | } |
---|
395 | if( targetParticle.GetDefinition() == aKaonZS ) |
---|
396 | { |
---|
397 | if( G4UniformRand() >= 0.5 ) |
---|
398 | { |
---|
399 | targetParticle.SetDefinitionAndUpdateE( aKaonZL ); |
---|
400 | targetHasChanged = true; |
---|
401 | } |
---|
402 | } |
---|
403 | SetUpPions( np, nm, nz, vec, vecLen ); |
---|
404 | return; |
---|
405 | } |
---|
406 | |
---|
407 | /* end of file */ |
---|
408 | |
---|