source: trunk/source/processes/hadronic/models/high_energy/src/G4HEAntiXiMinusInelastic.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: G4HEAntiXiMinusInelastic.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 "G4HEAntiXiMinusInelastic.hh"
43
44G4HadFinalState*
45G4HEAntiXiMinusInelastic::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 << "GHEAntiXiMinusInelastic: incident energy < 1 GeV" << G4endl;
65
66 if (verboseLevel > 1) {
67 G4cout << "G4HEAntiXiMinusInelastic::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 FirstIntInCasAntiXiMinus(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 HighEnergyCascading(successful, pv, vecLength,
132 excitationEnergyGNP, excitationEnergyDTA,
133 incidentParticle, targetParticle,
134 atomicWeight, atomicNumber);
135 if (!successful)
136 HighEnergyClusterProduction(successful, pv, vecLength,
137 excitationEnergyGNP, excitationEnergyDTA,
138 incidentParticle, targetParticle,
139 atomicWeight, atomicNumber);
140 if (!successful)
141 MediumEnergyCascading(successful, pv, vecLength,
142 excitationEnergyGNP, excitationEnergyDTA,
143 incidentParticle, targetParticle,
144 atomicWeight, atomicNumber);
145
146 if (!successful)
147 MediumEnergyClusterProduction(successful, pv, vecLength,
148 excitationEnergyGNP, excitationEnergyDTA,
149 incidentParticle, targetParticle,
150 atomicWeight, atomicNumber);
151 if (!successful)
152 QuasiElasticScattering(successful, pv, vecLength,
153 excitationEnergyGNP, excitationEnergyDTA,
154 incidentParticle, targetParticle,
155 atomicWeight, atomicNumber);
156 if (!successful)
157 ElasticScattering(successful, pv, vecLength,
158 incidentParticle,
159 atomicWeight, atomicNumber);
160
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
173G4HEAntiXiMinusInelastic::FirstIntInCasAntiXiMinus(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// AntiXi- 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{
185 static const G4double expxu = std::log(MAXFLOAT); // upper bound for arg. of exp
186 static const G4double expxl = -expxu; // lower bound for arg. of exp
187
188 static const G4double protb = 0.7;
189 static const G4double neutb = 0.7;
190 static const G4double c = 1.25;
191
192 static const G4int numMul = 1200;
193 static const G4int numMulAn = 400;
194 static const G4int numSec = 60;
195
196 G4int protonCode = Proton.getCode();
197
198 G4int targetCode = targetParticle.getCode();
199 G4double incidentTotalMomentum = incidentParticle.getTotalMomentum();
200
201 static G4bool first = true;
202 static G4double protmul[numMul], protnorm[numSec]; // proton constants
203 static G4double protmulAn[numMulAn],protnormAn[numSec];
204 static G4double neutmul[numMul], neutnorm[numSec]; // neutron constants
205 static G4double neutmulAn[numMulAn],neutnormAn[numSec];
206
207 // misc. local variables
208 // np = number of pi+, nm = number of pi-, nz = number of pi0
209
210 G4int i, counter, nt, np, nm, nz;
211
212 if( first )
213 { // compute normalization constants, this will only be done once
214 first = false;
215 for( i=0; i<numMul ; i++ ) protmul[i] = 0.0;
216 for( i=0; i<numSec ; i++ ) protnorm[i] = 0.0;
217 counter = -1;
218 for( np=0; np<(numSec/3); np++ )
219 {
220 for( nm=std::max(0,np-2); nm<=(np+1); nm++ )
221 {
222 for( nz=0; nz<numSec/3; nz++ )
223 {
224 if( ++counter < numMul )
225 {
226 nt = np+nm+nz;
227 if( (nt>0) && (nt<=numSec) )
228 {
229 protmul[counter] = pmltpc(np,nm,nz,nt,protb,c);
230 protnorm[nt-1] += protmul[counter];
231 }
232 }
233 }
234 }
235 }
236 for( i=0; i<numMul; i++ )neutmul[i] = 0.0;
237 for( i=0; i<numSec; i++ )neutnorm[i] = 0.0;
238 counter = -1;
239 for( np=0; np<numSec/3; np++ )
240 {
241 for( nm=std::max(0,np-1); nm<=(np+2); nm++ )
242 {
243 for( nz=0; nz<numSec/3; nz++ )
244 {
245 if( ++counter < numMul )
246 {
247 nt = np+nm+nz;
248 if( (nt>0) && (nt<=numSec) )
249 {
250 neutmul[counter] = pmltpc(np,nm,nz,nt,neutb,c);
251 neutnorm[nt-1] += neutmul[counter];
252 }
253 }
254 }
255 }
256 }
257 for( i=0; i<numSec; i++ )
258 {
259 if( protnorm[i] > 0.0 )protnorm[i] = 1.0/protnorm[i];
260 if( neutnorm[i] > 0.0 )neutnorm[i] = 1.0/neutnorm[i];
261 }
262 // annihilation
263 for( i=0; i<numMulAn ; i++ ) protmulAn[i] = 0.0;
264 for( i=0; i<numSec ; i++ ) protnormAn[i] = 0.0;
265 counter = -1;
266 for( np=1; np<(numSec/3); np++ )
267 {
268 nm = std::max(0,np-1);
269 for( nz=0; nz<numSec/3; nz++ )
270 {
271 if( ++counter < numMulAn )
272 {
273 nt = np+nm+nz;
274 if( (nt>1) && (nt<=numSec) )
275 {
276 protmulAn[counter] = pmltpc(np,nm,nz,nt,protb,c);
277 protnormAn[nt-1] += protmulAn[counter];
278 }
279 }
280 }
281 }
282 for( i=0; i<numMulAn; i++ ) neutmulAn[i] = 0.0;
283 for( i=0; i<numSec; i++ ) neutnormAn[i] = 0.0;
284 counter = -1;
285 for( np=0; np<numSec/3; np++ )
286 {
287 nm = np;
288 for( nz=0; nz<numSec/3; nz++ )
289 {
290 if( ++counter < numMulAn )
291 {
292 nt = np+nm+nz;
293 if( (nt>1) && (nt<=numSec) )
294 {
295 neutmulAn[counter] = pmltpc(np,nm,nz,nt,neutb,c);
296 neutnormAn[nt-1] += neutmulAn[counter];
297 }
298 }
299 }
300 }
301 for( i=0; i<numSec; i++ )
302 {
303 if( protnormAn[i] > 0.0 )protnormAn[i] = 1.0/protnormAn[i];
304 if( neutnormAn[i] > 0.0 )neutnormAn[i] = 1.0/neutnormAn[i];
305 }
306 } // end of initialization
307
308
309 // initialize the first two places
310 // the same as beam and target
311 pv[0] = incidentParticle;
312 pv[1] = targetParticle;
313 vecLen = 2;
314
315 if( !inElastic )
316 { // some two-body reactions
317 G4double cech[] = {0.50, 0.45, 0.40, 0.35, 0.30, 0.25, 0.06, 0.04, 0.005, 0.};
318
319 G4int iplab = std::min(9, G4int( incidentTotalMomentum*2.5));
320 if( G4UniformRand() < cech[iplab]/std::pow(atomicWeight,0.42) )
321 {
322 G4double ran = G4UniformRand();
323
324 if ( targetCode == protonCode)
325 {
326 if(ran < 0.2)
327 {
328 pv[0] = AntiSigmaZero;
329 }
330 else if (ran < 0.4)
331 {
332 pv[0] = AntiSigmaMinus;
333 pv[1] = Neutron;
334 }
335 else if (ran < 0.6)
336 {
337 pv[0] = Proton;
338 pv[1] = AntiLambda;
339 }
340 else if (ran < 0.8)
341 {
342 pv[0] = Proton;
343 pv[1] = AntiSigmaZero;
344 }
345 else
346 {
347 pv[0] = Neutron;
348 pv[1] = AntiSigmaMinus;
349 }
350 }
351 else
352 {
353 if (ran < 0.2)
354 {
355 pv[0] = AntiSigmaZero;
356 }
357 else if (ran < 0.4)
358 {
359 pv[0] = AntiSigmaPlus;
360 pv[1] = Proton;
361 }
362 else if (ran < 0.6)
363 {
364 pv[0] = Neutron;
365 pv[1] = AntiLambda;
366 }
367 else if (ran < 0.8)
368 {
369 pv[0] = Neutron;
370 pv[1] = AntiSigmaZero;
371 }
372 else
373 {
374 pv[0] = Proton;
375 pv[1] = AntiSigmaPlus;
376 }
377 }
378 }
379 return;
380 }
381 else if (availableEnergy <= PionPlus.getMass())
382 return;
383
384 // inelastic scattering
385
386 np = 0; nm = 0; nz = 0;
387 G4double anhl[] = {1.00, 1.00, 1.00, 1.00, 1.00, 1.00, 1.00, 1.00, 0.97, 0.88,
388 0.85, 0.81, 0.75, 0.64, 0.64, 0.55, 0.55, 0.45, 0.47, 0.40,
389 0.39, 0.36, 0.33, 0.10, 0.01};
390 G4int iplab = G4int( incidentTotalMomentum*10.);
391 if ( iplab > 9) iplab = 10 + G4int( (incidentTotalMomentum -1.)*5. );
392 if ( iplab > 14) iplab = 15 + G4int( incidentTotalMomentum -2. );
393 if ( iplab > 22) iplab = 23 + G4int( (incidentTotalMomentum -10.)/10.);
394 iplab = std::min(24, iplab);
395
396 if ( G4UniformRand() > anhl[iplab] )
397 { // non- annihilation channels
398
399 // number of total particles vs. centre of mass Energy - 2*proton mass
400
401 G4double aleab = std::log(availableEnergy);
402 G4double n = 3.62567+aleab*(0.665843+aleab*(0.336514
403 + aleab*(0.117712+0.0136912*aleab))) - 2.0;
404
405 // normalization constant for kno-distribution.
406 // calculate first the sum of all constants, check for numerical problems.
407 G4double test, dum, anpn = 0.0;
408
409 for (nt=1; nt<=numSec; nt++) {
410 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
411 dum = pi*nt/(2.0*n*n);
412 if (std::fabs(dum) < 1.0) {
413 if( test >= 1.0e-10 )anpn += dum*test;
414 } else {
415 anpn += dum*test;
416 }
417 }
418
419 G4double ran = G4UniformRand();
420 G4double excs = 0.0;
421 if( targetCode == protonCode )
422 {
423 counter = -1;
424 for( np=0; np<numSec/3; np++ )
425 {
426 for( nm=std::max(0,np-2); nm<=(np+1); nm++ )
427 {
428 for( nz=0; nz<numSec/3; nz++ )
429 {
430 if( ++counter < numMul )
431 {
432 nt = np+nm+nz;
433 if ( (nt>0) && (nt<=numSec) ) {
434 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
435 dum = (pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n);
436 if (std::fabs(dum) < 1.0) {
437 if( test >= 1.0e-10 )excs += dum*test;
438 } else {
439 excs += dum*test;
440 }
441
442 if (ran < excs) goto outOfLoop; //----------------------->
443 }
444 }
445 }
446 }
447 }
448
449 // 3 previous loops continued to the end
450 inElastic = false; // quasi-elastic scattering
451 return;
452 }
453 else
454 { // target must be a neutron
455 counter = -1;
456 for( np=0; np<numSec/3; np++ )
457 {
458 for( nm=std::max(0,np-1); nm<=(np+2); nm++ )
459 {
460 for( nz=0; nz<numSec/3; nz++ )
461 {
462 if( ++counter < numMul )
463 {
464 nt = np+nm+nz;
465 if ( (nt>0) && (nt<=numSec) ) {
466 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
467 dum = (pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n);
468 if (std::fabs(dum) < 1.0) {
469 if( test >= 1.0e-10 )excs += dum*test;
470 } else {
471 excs += dum*test;
472 }
473
474 if (ran < excs) goto outOfLoop; // -------------------------->
475 }
476 }
477 }
478 }
479 }
480 // 3 previous loops continued to the end
481 inElastic = false; // quasi-elastic scattering.
482 return;
483 }
484
485 outOfLoop: // <------------------------------------------------------------------------
486
487 ran = G4UniformRand();
488
489 if( targetCode == protonCode)
490 {
491 if( np == nm)
492 {
493 if (ran < 0.40)
494 {
495 }
496 else if (ran < 0.8)
497 {
498 pv[0] = AntiSigmaZero;
499 }
500 else
501 {
502 pv[0] = AntiSigmaMinus;
503 pv[1] = Neutron;
504 }
505 }
506 else if (np == (nm+1))
507 {
508 if( ran < 0.25)
509 {
510 pv[1] = Neutron;
511 }
512 else if (ran < 0.5)
513 {
514 pv[0] = AntiSigmaZero;
515 pv[1] = Neutron;
516 }
517 else
518 {
519 pv[0] = AntiSigmaPlus;
520 }
521 }
522 else if (np == (nm-1))
523 {
524 pv[0] = AntiSigmaMinus;
525 }
526 else
527 {
528 pv[0] = AntiSigmaPlus;
529 pv[1] = Neutron;
530 }
531 }
532 else
533 {
534 if( np == nm)
535 {
536 if (ran < 0.4)
537 {
538 }
539 else if(ran < 0.8)
540 {
541 pv[0] = AntiSigmaZero;
542 }
543 else
544 {
545 pv[0] = AntiSigmaPlus;
546 pv[1] = Proton;
547 }
548 }
549 else if ( np == (nm-1))
550 {
551 if (ran < 0.5)
552 {
553 pv[0] = AntiSigmaMinus;
554 }
555 else if (ran < 0.75)
556 {
557 pv[1] = Proton;
558 }
559 else
560 {
561 pv[0] = AntiSigmaZero;
562 pv[1] = Proton;
563 }
564 }
565 else if (np == (nm+1))
566 {
567 pv[0] = AntiSigmaPlus;
568 }
569 else
570 {
571 pv[0] = AntiSigmaMinus;
572 pv[1] = Proton;
573 }
574 }
575
576 }
577 else // annihilation
578 {
579 if ( availableEnergy > 2. * PionPlus.getMass() )
580 {
581
582 G4double aleab = std::log(availableEnergy);
583 G4double n = 3.62567+aleab*(0.665843+aleab*(0.336514
584 + aleab*(0.117712+0.0136912*aleab))) - 2.0;
585
586 // normalization constant for kno-distribution.
587 // calculate first the sum of all constants, check for numerical problems.
588 G4double test, dum, anpn = 0.0;
589
590 for (nt=2; nt<=numSec; nt++) {
591 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
592 dum = pi*nt/(2.0*n*n);
593 if (std::fabs(dum) < 1.0) {
594 if( test >= 1.0e-10 )anpn += dum*test;
595 } else {
596 anpn += dum*test;
597 }
598 }
599
600 G4double ran = G4UniformRand();
601 G4double excs = 0.0;
602 if( targetCode == protonCode )
603 {
604 counter = -1;
605 for( np=1; np<numSec/3; np++ )
606 {
607 nm = np-1;
608 for( nz=0; nz<numSec/3; nz++ )
609 {
610 if( ++counter < numMulAn )
611 {
612 nt = np+nm+nz;
613 if ( (nt>1) && (nt<=numSec) ) {
614 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
615 dum = (pi/anpn)*nt*protmulAn[counter]*protnormAn[nt-1]/(2.0*n*n);
616 if (std::fabs(dum) < 1.0) {
617 if( test >= 1.0e-10 )excs += dum*test;
618 } else {
619 excs += dum*test;
620 }
621
622 if (ran < excs) goto outOfLoopAn; //----------------------->
623 }
624 }
625 }
626 }
627 // 3 previous loops continued to the end
628 inElastic = false; // quasi-elastic scattering
629 return;
630 }
631 else
632 { // target must be a neutron
633 counter = -1;
634 for( np=0; np<numSec/3; np++ )
635 {
636 nm = np;
637 for( nz=0; nz<numSec/3; nz++ )
638 {
639 if( ++counter < numMulAn )
640 {
641 nt = np+nm+nz;
642 if ( (nt>1) && (nt<=numSec) ) {
643 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
644 dum = (pi/anpn)*nt*neutmulAn[counter]*neutnormAn[nt-1]/(2.0*n*n);
645 if (std::fabs(dum) < 1.0) {
646 if( test >= 1.0e-10 )excs += dum*test;
647 } else {
648 excs += dum*test;
649 }
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].getName() << " " ;
697 G4cout << pv[1].getName() << " " ;
698 for (i=2; i < vecLen; i++)
699 {
700 G4cout << pv[i].getName() << " " ;
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.