source: trunk/source/processes/hadronic/models/high_energy/src/G4HELambdaInelastic.cc@ 1350

Last change on this file since 1350 was 1347, checked in by garnier, 15 years ago

geant4 tag 9.4

File size: 18.8 KB
Line 
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: G4HELambdaInelastic.cc,v 1.17 2010/11/29 05:44:44 dennis Exp $
27// GEANT4 tag $Name: geant4-09-04-ref-00 $
28//
29
30#include "globals.hh"
31#include "G4ios.hh"
32
33// G4 Process: Gheisha High Energy Collision model.
34// This includes the high energy cascading model, the two-body-resonance model
35// and the low energy two-body model. Not included are the low energy stuff
36// like nuclear reactions, nuclear fission without any cascading and all
37// processes for particles at rest.
38// First work done by J.L.Chuma and F.W.Jones, TRIUMF, June 96.
39// H. Fesefeldt, RWTH-Aachen, 23-October-1996
40// Last modified: 29-July-1998
41
42#include "G4HELambdaInelastic.hh"
43
44G4HadFinalState*
45G4HELambdaInelastic::ApplyYourself(const G4HadProjectile& aTrack,
46 G4Nucleus& targetNucleus)
47{
48 G4HEVector* pv = new G4HEVector[MAXPART];
49 const G4HadProjectile* aParticle = &aTrack;
50 const G4double A = targetNucleus.GetN();
51 const G4double Z = targetNucleus.GetZ();
52 G4HEVector incidentParticle(aParticle);
53
54 G4double atomicNumber = Z;
55 G4double atomicWeight = A;
56
57 G4int incidentCode = incidentParticle.getCode();
58 G4double incidentMass = incidentParticle.getMass();
59 G4double incidentTotalEnergy = incidentParticle.getEnergy();
60 G4double incidentTotalMomentum = incidentParticle.getTotalMomentum();
61 G4double incidentKineticEnergy = incidentTotalEnergy - incidentMass;
62
63 if (incidentKineticEnergy < 1.)
64 G4cout << "GHELambdaInelastic: incident energy < 1 GeV" << G4endl;
65
66 if (verboseLevel > 1) {
67 G4cout << "G4HELambdaInelastic::ApplyYourself" << G4endl;
68 G4cout << "incident particle " << incidentParticle.getName()
69 << "mass " << incidentMass
70 << "kinetic energy " << incidentKineticEnergy
71 << G4endl;
72 G4cout << "target material with (A,Z) = ("
73 << atomicWeight << "," << atomicNumber << ")" << G4endl;
74 }
75
76 G4double inelasticity = NuclearInelasticity(incidentKineticEnergy,
77 atomicWeight, atomicNumber);
78 if (verboseLevel > 1)
79 G4cout << "nuclear inelasticity = " << inelasticity << G4endl;
80
81 incidentKineticEnergy -= inelasticity;
82
83 G4double excitationEnergyGNP = 0.;
84 G4double excitationEnergyDTA = 0.;
85
86 G4double excitation = NuclearExcitation(incidentKineticEnergy,
87 atomicWeight, atomicNumber,
88 excitationEnergyGNP,
89 excitationEnergyDTA);
90 if (verboseLevel > 1)
91 G4cout << "nuclear excitation = " << excitation << excitationEnergyGNP
92 << excitationEnergyDTA << G4endl;
93
94 incidentKineticEnergy -= excitation;
95 incidentTotalEnergy = incidentKineticEnergy + incidentMass;
96 incidentTotalMomentum = std::sqrt( (incidentTotalEnergy-incidentMass)
97 *(incidentTotalEnergy+incidentMass));
98
99 G4HEVector targetParticle;
100 if (G4UniformRand() < atomicNumber/atomicWeight) {
101 targetParticle.setDefinition("Proton");
102 } else {
103 targetParticle.setDefinition("Neutron");
104 }
105
106 G4double targetMass = targetParticle.getMass();
107 G4double centerOfMassEnergy = std::sqrt(incidentMass*incidentMass
108 + targetMass*targetMass
109 + 2.0*targetMass*incidentTotalEnergy);
110 G4double availableEnergy = centerOfMassEnergy - targetMass - incidentMass;
111
112 G4bool inElastic = true;
113 vecLength = 0;
114
115 if (verboseLevel > 1)
116 G4cout << "ApplyYourself: CallFirstIntInCascade for particle "
117 << incidentCode << G4endl;
118
119 G4bool successful = false;
120
121 FirstIntInCasLambda(inElastic, availableEnergy, pv, vecLength,
122 incidentParticle, targetParticle, atomicWeight);
123
124 if (verboseLevel > 1)
125 G4cout << "ApplyYourself::StrangeParticlePairProduction" << G4endl;
126
127 if ((vecLength > 0) && (availableEnergy > 1.))
128 StrangeParticlePairProduction(availableEnergy, centerOfMassEnergy,
129 pv, vecLength,
130 incidentParticle, targetParticle);
131
132 HighEnergyCascading(successful, pv, vecLength,
133 excitationEnergyGNP, excitationEnergyDTA,
134 incidentParticle, targetParticle,
135 atomicWeight, atomicNumber);
136 if (!successful)
137 HighEnergyClusterProduction(successful, pv, vecLength,
138 excitationEnergyGNP, excitationEnergyDTA,
139 incidentParticle, targetParticle,
140 atomicWeight, atomicNumber);
141 if (!successful)
142 MediumEnergyCascading(successful, pv, vecLength,
143 excitationEnergyGNP, excitationEnergyDTA,
144 incidentParticle, targetParticle,
145 atomicWeight, atomicNumber);
146
147 if (!successful)
148 MediumEnergyClusterProduction(successful, pv, vecLength,
149 excitationEnergyGNP, excitationEnergyDTA,
150 incidentParticle, targetParticle,
151 atomicWeight, atomicNumber);
152 if (!successful)
153 QuasiElasticScattering(successful, pv, vecLength,
154 excitationEnergyGNP, excitationEnergyDTA,
155 incidentParticle, targetParticle,
156 atomicWeight, atomicNumber);
157 if (!successful)
158 ElasticScattering(successful, pv, vecLength,
159 incidentParticle,
160 atomicWeight, atomicNumber);
161
162 if (!successful)
163 G4cout << "GHEInelasticInteraction::ApplyYourself fails to produce final state particles"
164 << G4endl;
165 FillParticleChange(pv, vecLength);
166 delete [] pv;
167 theParticleChange.SetStatusChange(stopAndKill);
168 return &theParticleChange;
169}
170
171
172void
173G4HELambdaInelastic::FirstIntInCasLambda(G4bool& inElastic,
174 const G4double availableEnergy,
175 G4HEVector pv[],
176 G4int& vecLen,
177 const G4HEVector& incidentParticle,
178 const G4HEVector& targetParticle,
179 const G4double atomicWeight)
180
181// Lambda undergoes interaction with nucleon within a nucleus. Check if it is
182// energetically possible to produce pions/kaons. In not, assume nuclear
183// excitation occurs and input particle is degraded in energy. No other
184// particles are produced. If reaction is possible, find the correct number
185// of pions/protons/neutrons produced using an interpolation to multiplicity
186// data. Replace some pions or protons/neutrons by kaons or strange baryons
187// according to the average multiplicity per inelastic reaction.
188{
189 static const G4double expxu = std::log(MAXFLOAT); // upper bound for arg. of exp
190 static const G4double expxl = -expxu; // lower bound for arg. of exp
191
192 static const G4double protb = 0.7;
193 static const G4double neutb = 0.7;
194 static const G4double c = 1.25;
195
196 static const G4int numMul = 1200;
197 static const G4int numSec = 60;
198
199 G4int protonCode = Proton.getCode();
200
201 G4int targetCode = targetParticle.getCode();
202 G4double incidentTotalMomentum = incidentParticle.getTotalMomentum();
203
204 static G4bool first = true;
205 static G4double protmul[numMul], protnorm[numSec]; // proton constants
206 static G4double neutmul[numMul], neutnorm[numSec]; // neutron constants
207
208 // misc. local variables
209 // np = number of pi+, nm = number of pi-, nz = number of pi0
210
211 G4int i, counter, nt, np, nm, nz;
212
213 if( first )
214 { // compute normalization constants, this will only be done once
215 first = false;
216 for( i=0; i<numMul; i++ )protmul[i] = 0.0;
217 for( i=0; i<numSec; i++ )protnorm[i] = 0.0;
218 counter = -1;
219 for( np=0; np<(numSec/3); np++ )
220 {
221 for( nm=std::max(0,np-2); nm<=(np+1); nm++ )
222 {
223 for( nz=0; nz<numSec/3; nz++ )
224 {
225 if( ++counter < numMul )
226 {
227 nt = np+nm+nz;
228 if( (nt>0) && (nt<=numSec) )
229 {
230 protmul[counter] = pmltpc(np,nm,nz,nt,protb,c);
231 protnorm[nt-1] += protmul[counter];
232 }
233 }
234 }
235 }
236 }
237 for( i=0; i<numMul; i++ )neutmul[i] = 0.0;
238 for( i=0; i<numSec; i++ )neutnorm[i] = 0.0;
239 counter = -1;
240 for( np=0; np<numSec/3; np++ )
241 {
242 for( nm=std::max(0,np-1); nm<=(np+2); nm++ )
243 {
244 for( nz=0; nz<numSec/3; nz++ )
245 {
246 if( ++counter < numMul )
247 {
248 nt = np+nm+nz;
249 if( (nt>0) && (nt<=numSec) )
250 {
251 neutmul[counter] = pmltpc(np,nm,nz,nt,neutb,c);
252 neutnorm[nt-1] += neutmul[counter];
253 }
254 }
255 }
256 }
257 }
258 for( i=0; i<numSec; i++ )
259 {
260 if( protnorm[i] > 0.0 )protnorm[i] = 1.0/protnorm[i];
261 if( neutnorm[i] > 0.0 )neutnorm[i] = 1.0/neutnorm[i];
262 }
263 } // end of initialization
264
265
266 // initialize the first two places
267 // the same as beam and target
268 pv[0] = incidentParticle;
269 pv[1] = targetParticle;
270 vecLen = 2;
271
272 if( !inElastic )
273 { // quasi-elastic scattering, no pions produced
274 G4double cech[] = {0.50, 0.45, 0.40, 0.35, 0.30, 0.25, 0.06, 0.04, 0.005, 0.};
275 G4int iplab = G4int( std::min( 9.0, incidentTotalMomentum*2.5 ) );
276 if( G4UniformRand() < cech[iplab]/std::pow(atomicWeight,0.42) )
277 {
278 G4double ran = G4UniformRand();
279 if( targetCode == protonCode)
280 {
281 if( ran < 0.2)
282 {
283 pv[0] = SigmaPlus;
284 pv[1] = Neutron;
285 }
286 else if(ran < 0.4)
287 {
288 pv[0] = SigmaZero;
289 }
290 else if(ran < 0.6)
291 {
292 pv[0] = Proton;
293 pv[1] = Lambda;
294 }
295 else if(ran < 0.8)
296 {
297 pv[0] = Proton;
298 pv[1] = SigmaZero;
299 }
300 else
301 {
302 pv[0] = Neutron;
303 pv[1] = SigmaPlus;
304 }
305 }
306 else
307 {
308 if(ran < 0.2)
309 {
310 pv[0] = SigmaZero;
311 }
312 else if(ran < 0.4)
313 {
314 pv[0] = SigmaMinus;
315 pv[1] = Proton;
316 }
317 else if(ran < 0.6)
318 {
319 pv[0] = Neutron;
320 pv[1] = Lambda;
321 }
322 else if(ran < 0.8)
323 {
324 pv[0] = Neutron;
325 pv[1] = SigmaZero;
326 }
327 else
328 {
329 pv[0] = Proton;
330 pv[1] = SigmaMinus;
331 }
332 }
333 }
334 return;
335 }
336 else if (availableEnergy <= PionPlus.getMass())
337 return;
338
339 // inelastic scattering
340
341 np = 0; nm = 0; nz = 0;
342
343 // number of total particles vs. centre of mass Energy - 2*proton mass
344
345 G4double aleab = std::log(availableEnergy);
346 G4double n = 3.62567+aleab*(0.665843+aleab*(0.336514
347 + aleab*(0.117712+0.0136912*aleab))) - 2.0;
348
349 // normalization constant for kno-distribution.
350 // calculate first the sum of all constants, check for numerical problems.
351 G4double test, dum, anpn = 0.0;
352
353 for (nt=1; nt<=numSec; nt++) {
354 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
355 dum = pi*nt/(2.0*n*n);
356 if (std::fabs(dum) < 1.0) {
357 if( test >= 1.0e-10 )anpn += dum*test;
358 } else {
359 anpn += dum*test;
360 }
361 }
362
363 G4double ran = G4UniformRand();
364 G4double excs = 0.0;
365 if( targetCode == protonCode )
366 {
367 counter = -1;
368 for (np=0; np<numSec/3; np++) {
369 for (nm=std::max(0,np-2); nm<=(np+1); nm++) {
370 for (nz=0; nz<numSec/3; nz++) {
371 if (++counter < numMul) {
372 nt = np+nm+nz;
373 if ( (nt>0) && (nt<=numSec) ) {
374 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
375 dum = (pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n);
376 if (std::fabs(dum) < 1.0) {
377 if( test >= 1.0e-10 )excs += dum*test;
378 } else {
379 excs += dum*test;
380 }
381 if (ran < excs) goto outOfLoop; //--------------------->
382 }
383 }
384 }
385 }
386 }
387 // 3 previous loops continued to the end
388
389 inElastic = false; // quasi-elastic scattering
390 return;
391 }
392 else
393 { // target must be a neutron
394 counter = -1;
395 for (np=0; np<numSec/3; np++) {
396 for (nm=std::max(0,np-1); nm<=(np+2); nm++) {
397 for (nz=0; nz<numSec/3; nz++) {
398 if (++counter < numMul) {
399 nt = np+nm+nz;
400 if ( (nt>=1) && (nt<=numSec) ) {
401 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
402 dum = (pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n);
403 if (std::fabs(dum) < 1.0) {
404 if( test >= 1.0e-10 )excs += dum*test;
405 } else {
406 excs += dum*test;
407 }
408 if (ran < excs) goto outOfLoop; // ------------->
409 }
410 }
411 }
412 }
413 }
414 // 3 previous loops continued to the end
415 inElastic = false; // quasi-elastic scattering.
416 return;
417 }
418
419 outOfLoop: // <--------------------------------------------
420
421 ran = G4UniformRand();
422 if( targetCode == protonCode)
423 {
424 if( np == nm)
425 {
426 if (ran < 0.25)
427 {
428 }
429 else if(ran < 0.5)
430 {
431 pv[0] = SigmaZero;
432 }
433 else
434 {
435 pv[0] = SigmaPlus;
436 pv[1] = Neutron;
437 }
438 }
439 else if (np == (nm+1))
440 {
441 if( G4UniformRand() < 0.25)
442 {
443 pv[1] = Neutron;
444 }
445 else if(ran < 0.5)
446 {
447 pv[0] = SigmaZero;
448 pv[1] = Neutron;
449 }
450 else
451 {
452 pv[0] = SigmaMinus;
453 }
454 }
455 else if (np == (nm-1))
456 {
457 pv[0] = SigmaPlus;
458 }
459 else
460 {
461 pv[0] = SigmaMinus;
462 pv[1] = Neutron;
463 }
464 }
465 else
466 {
467 if (np == nm)
468 {
469 if(ran < 0.5)
470 {
471 }
472 else
473 {
474 pv[0] = SigmaMinus;
475 pv[1] = Proton;
476 }
477 }
478 else if (np == (nm-1))
479 {
480 if( ran < 0.25)
481 {
482 pv[1] = Proton;
483 }
484 else if(ran < 0.5)
485 {
486 pv[0] = SigmaZero;
487 pv[1] = Proton;
488 }
489 else
490 {
491 pv[1] = SigmaPlus;
492 }
493 }
494 else if (np == (1+nm))
495 {
496 pv[0] = SigmaMinus;
497 }
498 else
499 {
500 pv[0] = SigmaPlus;
501 pv[1] = Proton;
502 }
503 }
504
505
506 nt = np + nm + nz;
507 while ( nt > 0)
508 {
509 G4double ran = G4UniformRand();
510 if ( ran < (G4double)np/nt)
511 {
512 if( np > 0 )
513 { pv[vecLen++] = PionPlus;
514 np--;
515 }
516 }
517 else if ( ran < (G4double)(np+nm)/nt)
518 {
519 if( nm > 0 )
520 {
521 pv[vecLen++] = PionMinus;
522 nm--;
523 }
524 }
525 else
526 {
527 if( nz > 0 )
528 {
529 pv[vecLen++] = PionZero;
530 nz--;
531 }
532 }
533 nt = np + nm + nz;
534 }
535 if (verboseLevel > 1)
536 {
537 G4cout << "Particles produced: " ;
538 G4cout << pv[0].getName() << " " ;
539 G4cout << pv[1].getName() << " " ;
540 for (i=2; i < vecLen; i++)
541 {
542 G4cout << pv[i].getName() << " " ;
543 }
544 G4cout << G4endl;
545 }
546 return;
547 }
548
549
550
551
552
553
554
555
556
557
Note: See TracBrowser for help on using the repository browser.