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

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

geant4 tag 9.4

File size: 26.5 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: G4HEAntiXiZeroInelastic.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 "G4HEAntiXiZeroInelastic.hh"
43
44G4HadFinalState*
45G4HEAntiXiZeroInelastic::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 << "GHEAntiXiZeroInelastic: incident energy < 1 GeV" << G4endl;
65
66 if (verboseLevel > 1) {
67 G4cout << "G4HEAntiXiZeroInelastic::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 FirstIntInCasAntiXiZero(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 if (!successful)
162 G4cout << "GHEInelasticInteraction::ApplyYourself fails to produce final state particles"
163 << G4endl;
164
165 FillParticleChange(pv, vecLength);
166 delete [] pv;
167 theParticleChange.SetStatusChange(stopAndKill);
168 return &theParticleChange;
169}
170
171
172void
173G4HEAntiXiZeroInelastic::FirstIntInCasAntiXiZero(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// AntiXi0 undergoes interaction with nucleon within a nucleus.
182// As in Geant3, we think that this routine has absolutely no influence
183// on the whole performance of the program. Take AntiLambda instaed.
184// ( decay Xi0 -> L Pi > 99 % )
185{
186 static const G4double expxu = std::log(MAXFLOAT); // upper bound for arg. of exp
187 static const G4double expxl = -expxu; // lower bound for arg. of exp
188
189 static const G4double protb = 0.7;
190 static const G4double neutb = 0.7;
191 static const G4double c = 1.25;
192
193 static const G4int numMul = 1200;
194 static const G4int numMulAn = 400;
195 static const G4int numSec = 60;
196
197 G4int protonCode = Proton.getCode();
198
199 G4int targetCode = targetParticle.getCode();
200 G4double incidentTotalMomentum = incidentParticle.getTotalMomentum();
201
202 static G4bool first = true;
203 static G4double protmul[numMul], protnorm[numSec]; // proton constants
204 static G4double protmulAn[numMulAn],protnormAn[numSec];
205 static G4double neutmul[numMul], neutnorm[numSec]; // neutron constants
206 static G4double neutmulAn[numMulAn],neutnormAn[numSec];
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 // annihilation
264 for( i=0; i<numMulAn ; i++ ) protmulAn[i] = 0.0;
265 for( i=0; i<numSec ; i++ ) protnormAn[i] = 0.0;
266 counter = -1;
267 for( np=1; np<(numSec/3); np++ )
268 {
269 nm = std::max(0,np-1);
270 for( nz=0; nz<numSec/3; nz++ )
271 {
272 if( ++counter < numMulAn )
273 {
274 nt = np+nm+nz;
275 if( (nt>1) && (nt<=numSec) )
276 {
277 protmulAn[counter] = pmltpc(np,nm,nz,nt,protb,c);
278 protnormAn[nt-1] += protmulAn[counter];
279 }
280 }
281 }
282 }
283 for( i=0; i<numMulAn; i++ ) neutmulAn[i] = 0.0;
284 for( i=0; i<numSec; i++ ) neutnormAn[i] = 0.0;
285 counter = -1;
286 for( np=0; np<numSec/3; np++ )
287 {
288 nm = np;
289 for( nz=0; nz<numSec/3; nz++ )
290 {
291 if( ++counter < numMulAn )
292 {
293 nt = np+nm+nz;
294 if( (nt>1) && (nt<=numSec) )
295 {
296 neutmulAn[counter] = pmltpc(np,nm,nz,nt,neutb,c);
297 neutnormAn[nt-1] += neutmulAn[counter];
298 }
299 }
300 }
301 }
302 for( i=0; i<numSec; i++ )
303 {
304 if( protnormAn[i] > 0.0 )protnormAn[i] = 1.0/protnormAn[i];
305 if( neutnormAn[i] > 0.0 )neutnormAn[i] = 1.0/neutnormAn[i];
306 }
307 } // end of initialization
308
309
310 // initialize the first two places
311 // the same as beam and target
312 pv[0] = incidentParticle;
313 pv[1] = targetParticle;
314 vecLen = 2;
315
316 if( !inElastic )
317 { // some two-body reactions
318 G4double cech[] = {0.50, 0.45, 0.40, 0.35, 0.30, 0.25, 0.06, 0.04, 0.005, 0.};
319
320 G4int iplab = std::min(9, G4int( incidentTotalMomentum*2.5 ));
321 if( G4UniformRand() < cech[iplab]/std::pow(atomicWeight,0.42) )
322 {
323 G4double ran = G4UniformRand();
324
325 if ( targetCode == protonCode)
326 {
327 if(ran < 0.2)
328 {
329 pv[0] = AntiSigmaZero;
330 }
331 else if (ran < 0.4)
332 {
333 pv[0] = AntiSigmaMinus;
334 pv[1] = Neutron;
335 }
336 else if (ran < 0.6)
337 {
338 pv[0] = Proton;
339 pv[1] = AntiLambda;
340 }
341 else if (ran < 0.8)
342 {
343 pv[0] = Proton;
344 pv[1] = AntiSigmaZero;
345 }
346 else
347 {
348 pv[0] = Neutron;
349 pv[1] = AntiSigmaMinus;
350 }
351 }
352 else
353 {
354 if (ran < 0.2)
355 {
356 pv[0] = AntiSigmaZero;
357 }
358 else if (ran < 0.4)
359 {
360 pv[0] = AntiSigmaPlus;
361 pv[1] = Proton;
362 }
363 else if (ran < 0.6)
364 {
365 pv[0] = Neutron;
366 pv[1] = AntiLambda;
367 }
368 else if (ran < 0.8)
369 {
370 pv[0] = Neutron;
371 pv[1] = AntiSigmaZero;
372 }
373 else
374 {
375 pv[0] = Proton;
376 pv[1] = AntiSigmaPlus;
377 }
378 }
379 }
380 return;
381 }
382 else if (availableEnergy <= PionPlus.getMass())
383 return;
384
385 // inelastic scattering
386
387 np = 0; nm = 0; nz = 0;
388 G4double anhl[] = {1.00, 1.00, 1.00, 1.00, 1.00, 1.00, 1.00, 1.00, 0.97, 0.88,
389 0.85, 0.81, 0.75, 0.64, 0.64, 0.55, 0.55, 0.45, 0.47, 0.40,
390 0.39, 0.36, 0.33, 0.10, 0.01};
391 G4int iplab = G4int( incidentTotalMomentum*10.);
392 if ( iplab > 9) iplab = 10 + G4int( (incidentTotalMomentum -1.)*5. );
393 if ( iplab > 14) iplab = 15 + G4int( incidentTotalMomentum -2. );
394 if ( iplab > 22) iplab = 23 + G4int( (incidentTotalMomentum -10.)/10.);
395 iplab = std::min(24, iplab);
396
397 if ( G4UniformRand() > anhl[iplab] )
398 { // non- annihilation channels
399
400 // number of total particles vs. centre of mass Energy - 2*proton mass
401
402 G4double aleab = std::log(availableEnergy);
403 G4double n = 3.62567+aleab*(0.665843+aleab*(0.336514
404 + aleab*(0.117712+0.0136912*aleab))) - 2.0;
405
406 // normalization constant for kno-distribution.
407 // calculate first the sum of all constants, check for numerical problems.
408 G4double test, dum, anpn = 0.0;
409
410 for (nt=1; nt<=numSec; nt++) {
411 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
412 dum = pi*nt/(2.0*n*n);
413 if (std::fabs(dum) < 1.0) {
414 if( test >= 1.0e-10 )anpn += dum*test;
415 } else {
416 anpn += dum*test;
417 }
418 }
419
420 G4double ran = G4UniformRand();
421 G4double excs = 0.0;
422 if( targetCode == protonCode )
423 {
424 counter = -1;
425 for( np=0; np<numSec/3; np++ )
426 {
427 for( nm=std::max(0,np-2); nm<=(np+1); nm++ )
428 {
429 for( nz=0; nz<numSec/3; nz++ )
430 {
431 if( ++counter < numMul )
432 {
433 nt = np+nm+nz;
434 if ( (nt>0) && (nt<=numSec) ) {
435 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
436 dum = (pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n);
437 if (std::fabs(dum) < 1.0) {
438 if( test >= 1.0e-10 )excs += dum*test;
439 } else {
440 excs += dum*test;
441 }
442
443 if (ran < excs) goto outOfLoop; //----------------------->
444 }
445 }
446 }
447 }
448 }
449
450 // 3 previous loops continued to the end
451 inElastic = false; // quasi-elastic scattering
452 return;
453 }
454 else
455 { // target must be a neutron
456 counter = -1;
457 for( np=0; np<numSec/3; np++ )
458 {
459 for( nm=std::max(0,np-1); nm<=(np+2); nm++ )
460 {
461 for( nz=0; nz<numSec/3; nz++ )
462 {
463 if( ++counter < numMul )
464 {
465 nt = np+nm+nz;
466 if ( (nt>0) && (nt<=numSec) ) {
467 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
468 dum = (pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n);
469 if (std::fabs(dum) < 1.0) {
470 if( test >= 1.0e-10 )excs += dum*test;
471 } else {
472 excs += dum*test;
473 }
474
475 if (ran < excs) goto outOfLoop; // -------------------------->
476 }
477 }
478 }
479 }
480 }
481 // 3 previous loops continued to the end
482 inElastic = false; // quasi-elastic scattering.
483 return;
484 }
485
486 outOfLoop: // <------------------------------------------------------------------------
487
488 ran = G4UniformRand();
489
490 if( targetCode == protonCode)
491 {
492 if( np == nm)
493 {
494 if (ran < 0.40)
495 {
496 }
497 else if (ran < 0.8)
498 {
499 pv[0] = AntiSigmaZero;
500 }
501 else
502 {
503 pv[0] = AntiSigmaMinus;
504 pv[1] = Neutron;
505 }
506 }
507 else if (np == (nm+1))
508 {
509 if( ran < 0.25)
510 {
511 pv[1] = Neutron;
512 }
513 else if (ran < 0.5)
514 {
515 pv[0] = AntiSigmaZero;
516 pv[1] = Neutron;
517 }
518 else
519 {
520 pv[0] = AntiSigmaPlus;
521 }
522 }
523 else if (np == (nm-1))
524 {
525 pv[0] = AntiSigmaMinus;
526 }
527 else
528 {
529 pv[0] = AntiSigmaPlus;
530 pv[1] = Neutron;
531 }
532 }
533 else
534 {
535 if( np == nm)
536 {
537 if (ran < 0.4)
538 {
539 }
540 else if(ran < 0.8)
541 {
542 pv[0] = AntiSigmaZero;
543 }
544 else
545 {
546 pv[0] = AntiSigmaPlus;
547 pv[1] = Proton;
548 }
549 }
550 else if ( np == (nm-1))
551 {
552 if (ran < 0.5)
553 {
554 pv[0] = AntiSigmaMinus;
555 }
556 else if (ran < 0.75)
557 {
558 pv[1] = Proton;
559 }
560 else
561 {
562 pv[0] = AntiSigmaZero;
563 pv[1] = Proton;
564 }
565 }
566 else if (np == (nm+1))
567 {
568 pv[0] = AntiSigmaPlus;
569 }
570 else
571 {
572 pv[0] = AntiSigmaMinus;
573 pv[1] = Proton;
574 }
575 }
576
577 }
578 else // annihilation
579 {
580 if ( availableEnergy > 2. * PionPlus.getMass() )
581 {
582
583 G4double aleab = std::log(availableEnergy);
584 G4double n = 3.62567+aleab*(0.665843+aleab*(0.336514
585 + aleab*(0.117712+0.0136912*aleab))) - 2.0;
586
587 // normalization constant for kno-distribution.
588 // calculate first the sum of all constants, check for numerical problems.
589 G4double test, dum, anpn = 0.0;
590
591 for (nt=2; nt<=numSec; nt++) {
592 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
593 dum = pi*nt/(2.0*n*n);
594 if (std::fabs(dum) < 1.0) {
595 if( test >= 1.0e-10 )anpn += dum*test;
596 } else {
597 anpn += dum*test;
598 }
599 }
600
601 G4double ran = G4UniformRand();
602 G4double excs = 0.0;
603 if( targetCode == protonCode )
604 {
605 counter = -1;
606 for( np=1; np<numSec/3; np++ )
607 {
608 nm = np-1;
609 for( nz=0; nz<numSec/3; nz++ )
610 {
611 if( ++counter < numMulAn )
612 {
613 nt = np+nm+nz;
614 if ( (nt>1) && (nt<=numSec) ) {
615 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
616 dum = (pi/anpn)*nt*protmulAn[counter]*protnormAn[nt-1]/(2.0*n*n);
617 if (std::fabs(dum) < 1.0) {
618 if( test >= 1.0e-10 )excs += dum*test;
619 } else {
620 excs += dum*test;
621 }
622
623 if (ran < excs) goto outOfLoopAn; //----------------------->
624 }
625 }
626 }
627 }
628 // 3 previous loops continued to the end
629 inElastic = false; // quasi-elastic scattering
630 return;
631 }
632 else
633 { // target must be a neutron
634 counter = -1;
635 for( np=0; np<numSec/3; np++ )
636 {
637 nm = np;
638 for( nz=0; nz<numSec/3; nz++ )
639 {
640 if( ++counter < numMulAn )
641 {
642 nt = np+nm+nz;
643 if ( (nt>1) && (nt<=numSec) ) {
644 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
645 dum = (pi/anpn)*nt*neutmulAn[counter]*neutnormAn[nt-1]/(2.0*n*n);
646 if (std::fabs(dum) < 1.0) {
647 if( test >= 1.0e-10 )excs += dum*test;
648 } else {
649 excs += dum*test;
650 }
651 if (ran < excs) goto outOfLoopAn; // -------------------------->
652 }
653 }
654 }
655 }
656 inElastic = false; // quasi-elastic scattering.
657 return;
658 }
659 outOfLoopAn: // <------------------------------------------------------------------
660 vecLen = 0;
661 }
662 }
663
664 nt = np + nm + nz;
665 while ( nt > 0)
666 {
667 G4double ran = G4UniformRand();
668 if ( ran < (G4double)np/nt)
669 {
670 if( np > 0 )
671 { pv[vecLen++] = PionPlus;
672 np--;
673 }
674 }
675 else if ( ran < (G4double)(np+nm)/nt)
676 {
677 if( nm > 0 )
678 {
679 pv[vecLen++] = PionMinus;
680 nm--;
681 }
682 }
683 else
684 {
685 if( nz > 0 )
686 {
687 pv[vecLen++] = PionZero;
688 nz--;
689 }
690 }
691 nt = np + nm + nz;
692 }
693 if (verboseLevel > 1)
694 {
695 G4cout << "Particles produced: " ;
696 G4cout << pv[0].getCode() << " " ;
697 G4cout << pv[1].getCode() << " " ;
698 for (i=2; i < vecLen; i++)
699 {
700 G4cout << pv[i].getCode() << " " ;
701 }
702 G4cout << G4endl;
703 }
704 return;
705 }
706
707
708
709
710
711
712
713
714
715
Note: See TracBrowser for help on using the repository browser.