source: trunk/source/processes/hadronic/models/high_energy/src/G4HEAntiSigmaPlusInelastic.cc@ 1036

Last change on this file since 1036 was 1007, checked in by garnier, 17 years ago

update to geant4.9.2

File size: 27.3 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//
27// $Id: G4HEAntiSigmaPlusInelastic.cc,v 1.15 2008/03/17 20:49:17 dennis Exp $
28// GEANT4 tag $Name: geant4-09-02 $
29//
30//
31
32#include "globals.hh"
33#include "G4ios.hh"
34
35//
36// G4 Process: Gheisha High Energy Collision model.
37// This includes the high energy cascading model, the two-body-resonance model
38// and the low energy two-body model. Not included are the low energy stuff like
39// nuclear reactions, nuclear fission without any cascading and all processes for
40// particles at rest.
41// First work done by J.L.Chuma and F.W.Jones, TRIUMF, June 96.
42// H. Fesefeldt, RWTH-Aachen, 23-October-1996
43// Last modified: 29-July-1998
44
45#include "G4HEAntiSigmaPlusInelastic.hh"
46
47G4HadFinalState * G4HEAntiSigmaPlusInelastic::
48ApplyYourself( const G4HadProjectile &aTrack, G4Nucleus &targetNucleus )
49 {
50 G4HEVector * pv = new G4HEVector[MAXPART];
51 const G4HadProjectile *aParticle = &aTrack;
52// G4DynamicParticle *originalTarget = targetNucleus.ReturnTargetParticle();
53 const G4double atomicWeight = targetNucleus.GetN();
54 const G4double atomicNumber = targetNucleus.GetZ();
55 G4HEVector incidentParticle(aParticle);
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 {
65 G4cout << "GHEAntiSigmaPlusInelastic: incident energy < 1 GeV" << G4endl;
66 }
67 if(verboseLevel > 1)
68 {
69 G4cout << "G4HEAntiSigmaPlusInelastic::ApplyYourself" << G4endl;
70 G4cout << "incident particle " << incidentParticle.getName()
71 << "mass " << incidentMass
72 << "kinetic energy " << incidentKineticEnergy
73 << G4endl;
74 G4cout << "target material with (A,Z) = ("
75 << atomicWeight << "," << atomicNumber << ")" << G4endl;
76 }
77
78 G4double inelasticity = NuclearInelasticity(incidentKineticEnergy,
79 atomicWeight, atomicNumber);
80 if(verboseLevel > 1)
81 G4cout << "nuclear inelasticity = " << inelasticity << G4endl;
82
83 incidentKineticEnergy -= inelasticity;
84
85 G4double excitationEnergyGNP = 0.;
86 G4double excitationEnergyDTA = 0.;
87
88 G4double excitation = NuclearExcitation(incidentKineticEnergy,
89 atomicWeight, atomicNumber,
90 excitationEnergyGNP,
91 excitationEnergyDTA);
92 if(verboseLevel > 1)
93 G4cout << "nuclear excitation = " << excitation << excitationEnergyGNP
94 << excitationEnergyDTA << G4endl;
95
96
97 incidentKineticEnergy -= excitation;
98 incidentTotalEnergy = incidentKineticEnergy + incidentMass;
99 incidentTotalMomentum = std::sqrt( (incidentTotalEnergy-incidentMass)
100 *(incidentTotalEnergy+incidentMass));
101
102
103 G4HEVector targetParticle;
104 if(G4UniformRand() < atomicNumber/atomicWeight)
105 {
106 targetParticle.setDefinition("Proton");
107 }
108 else
109 {
110 targetParticle.setDefinition("Neutron");
111 }
112
113 G4double targetMass = targetParticle.getMass();
114 G4double centerOfMassEnergy = std::sqrt( incidentMass*incidentMass + targetMass*targetMass
115 + 2.0*targetMass*incidentTotalEnergy);
116 G4double availableEnergy = centerOfMassEnergy - targetMass - incidentMass;
117
118 // this was the meaning of inElastic in the
119 // original Gheisha stand-alone version.
120// G4bool inElastic = InElasticCrossSectionInFirstInt
121// (availableEnergy, incidentCode, incidentTotalMomentum);
122 // by unknown reasons, it has been replaced
123 // to the following code in Geant???
124 G4bool inElastic = true;
125// if (G4UniformRand() < elasticCrossSection/totalCrossSection) inElastic = false;
126
127 vecLength = 0;
128
129 if(verboseLevel > 1)
130 G4cout << "ApplyYourself: CallFirstIntInCascade for particle "
131 << incidentCode << G4endl;
132
133 G4bool successful = false;
134
135 if(inElastic || (!inElastic && atomicWeight < 1.5))
136 {
137 FirstIntInCasAntiSigmaPlus(inElastic, availableEnergy, pv, vecLength,
138 incidentParticle, targetParticle, atomicWeight);
139
140 if(verboseLevel > 1)
141 G4cout << "ApplyYourself::StrangeParticlePairProduction" << G4endl;
142
143
144 if ((vecLength > 0) && (availableEnergy > 1.))
145 StrangeParticlePairProduction( availableEnergy, centerOfMassEnergy,
146 pv, vecLength,
147 incidentParticle, targetParticle);
148 HighEnergyCascading( successful, pv, vecLength,
149 excitationEnergyGNP, excitationEnergyDTA,
150 incidentParticle, targetParticle,
151 atomicWeight, atomicNumber);
152 if (!successful)
153 HighEnergyClusterProduction( successful, pv, vecLength,
154 excitationEnergyGNP, excitationEnergyDTA,
155 incidentParticle, targetParticle,
156 atomicWeight, atomicNumber);
157 if (!successful)
158 MediumEnergyCascading( successful, pv, vecLength,
159 excitationEnergyGNP, excitationEnergyDTA,
160 incidentParticle, targetParticle,
161 atomicWeight, atomicNumber);
162
163 if (!successful)
164 MediumEnergyClusterProduction( successful, pv, vecLength,
165 excitationEnergyGNP, excitationEnergyDTA,
166 incidentParticle, targetParticle,
167 atomicWeight, atomicNumber);
168 if (!successful)
169 QuasiElasticScattering( successful, pv, vecLength,
170 excitationEnergyGNP, excitationEnergyDTA,
171 incidentParticle, targetParticle,
172 atomicWeight, atomicNumber);
173 }
174 if (!successful)
175 {
176 ElasticScattering( successful, pv, vecLength,
177 incidentParticle,
178 atomicWeight, atomicNumber);
179 }
180
181 if (!successful)
182 {
183 G4cout << "GHEInelasticInteraction::ApplyYourself fails to produce final state particles" << G4endl;
184 }
185 FillParticleChange(pv, vecLength);
186 delete [] pv;
187 theParticleChange.SetStatusChange(stopAndKill);
188 return & theParticleChange;
189 }
190
191void
192G4HEAntiSigmaPlusInelastic::FirstIntInCasAntiSigmaPlus( G4bool &inElastic,
193 const G4double availableEnergy,
194 G4HEVector pv[],
195 G4int &vecLen,
196 G4HEVector incidentParticle,
197 G4HEVector targetParticle,
198 const G4double atomicWeight)
199
200// AntiSigma+ undergoes interaction with nucleon within a nucleus. Check if it is
201// energetically possible to produce pions/kaons. In not, assume nuclear excitation
202// occurs and input particle is degraded in energy. No other particles are produced.
203// If reaction is possible, find the correct number of pions/protons/neutrons
204// produced using an interpolation to multiplicity data. Replace some pions or
205// protons/neutrons by kaons or strange baryons according to the average
206// multiplicity per inelastic reaction.
207
208 {
209 static const G4double expxu = std::log(MAXFLOAT); // upper bound for arg. of exp
210 static const G4double expxl = -expxu; // lower bound for arg. of exp
211
212 static const G4double protb = 0.7;
213 static const G4double neutb = 0.7;
214 static const G4double c = 1.25;
215
216 static const G4int numMul = 1200;
217 static const G4int numMulAn = 400;
218 static const G4int numSec = 60;
219
220// G4int neutronCode = Neutron.getCode();
221 G4int protonCode = Proton.getCode();
222
223 G4int targetCode = targetParticle.getCode();
224// G4double incidentMass = incidentParticle.getMass();
225// G4double incidentEnergy = incidentParticle.getEnergy();
226 G4double incidentTotalMomentum = incidentParticle.getTotalMomentum();
227
228 static G4bool first = true;
229 static G4double protmul[numMul], protnorm[numSec]; // proton constants
230 static G4double protmulAn[numMulAn],protnormAn[numSec];
231 static G4double neutmul[numMul], neutnorm[numSec]; // neutron constants
232 static G4double neutmulAn[numMulAn],neutnormAn[numSec];
233
234 // misc. local variables
235 // np = number of pi+, nm = number of pi-, nz = number of pi0
236
237 G4int i, counter, nt, np, nm, nz;
238
239 if( first )
240 { // compute normalization constants, this will only be done once
241 first = false;
242 for( i=0; i<numMul ; i++ ) protmul[i] = 0.0;
243 for( i=0; i<numSec ; i++ ) protnorm[i] = 0.0;
244 counter = -1;
245 for( np=0; np<(numSec/3); np++ )
246 {
247 for( nm=std::max(0,np-1); nm<=(np+1); nm++ )
248 {
249 for( nz=0; nz<numSec/3; nz++ )
250 {
251 if( ++counter < numMul )
252 {
253 nt = np+nm+nz;
254 if( (nt>0) && (nt<=numSec) )
255 {
256 protmul[counter] = pmltpc(np,nm,nz,nt,protb,c);
257 protnorm[nt-1] += protmul[counter];
258 }
259 }
260 }
261 }
262 }
263 for( i=0; i<numMul; i++ )neutmul[i] = 0.0;
264 for( i=0; i<numSec; i++ )neutnorm[i] = 0.0;
265 counter = -1;
266 for( np=0; np<numSec/3; np++ )
267 {
268 for( nm=np; nm<=(np+2); nm++ )
269 {
270 for( nz=0; nz<numSec/3; nz++ )
271 {
272 if( ++counter < numMul )
273 {
274 nt = np+nm+nz;
275 if( (nt>0) && (nt<=numSec) )
276 {
277 neutmul[counter] = pmltpc(np,nm,nz,nt,neutb,c);
278 neutnorm[nt-1] += neutmul[counter];
279 }
280 }
281 }
282 }
283 }
284 for( i=0; i<numSec; i++ )
285 {
286 if( protnorm[i] > 0.0 )protnorm[i] = 1.0/protnorm[i];
287 if( neutnorm[i] > 0.0 )neutnorm[i] = 1.0/neutnorm[i];
288 }
289 // annihilation
290 for( i=0; i<numMulAn ; i++ ) protmulAn[i] = 0.0;
291 for( i=0; i<numSec ; i++ ) protnormAn[i] = 0.0;
292 counter = -1;
293 for( np=1; np<(numSec/3); np++ )
294 {
295 nm = np;
296 for( nz=0; nz<numSec/3; nz++ )
297 {
298 if( ++counter < numMulAn )
299 {
300 nt = np+nm+nz;
301 if( (nt>1) && (nt<=numSec) )
302 {
303 protmulAn[counter] = pmltpc(np,nm,nz,nt,protb,c);
304 protnormAn[nt-1] += protmulAn[counter];
305 }
306 }
307 }
308 }
309 for( i=0; i<numMulAn; i++ ) neutmulAn[i] = 0.0;
310 for( i=0; i<numSec; i++ ) neutnormAn[i] = 0.0;
311 counter = -1;
312 for( np=0; np<numSec/3; np++ )
313 {
314 nm = np+1;
315 for( nz=0; nz<numSec/3; nz++ )
316 {
317 if( ++counter < numMulAn )
318 {
319 nt = np+nm+nz;
320 if( (nt>1) && (nt<=numSec) )
321 {
322 neutmulAn[counter] = pmltpc(np,nm,nz,nt,neutb,c);
323 neutnormAn[nt-1] += neutmulAn[counter];
324 }
325 }
326 }
327 }
328 for( i=0; i<numSec; i++ )
329 {
330 if( protnormAn[i] > 0.0 )protnormAn[i] = 1.0/protnormAn[i];
331 if( neutnormAn[i] > 0.0 )neutnormAn[i] = 1.0/neutnormAn[i];
332 }
333 } // end of initialization
334
335
336 // initialize the first two places
337 // the same as beam and target
338 pv[0] = incidentParticle;
339 pv[1] = targetParticle;
340 vecLen = 2;
341
342 if( !inElastic )
343 { // some two-body reactions
344 G4double cech[] = {0.50, 0.45, 0.40, 0.35, 0.30, 0.25, 0.06, 0.04, 0.005, 0.};
345
346 G4int iplab = std::min(9, G4int( incidentTotalMomentum*2.5));
347 if( G4UniformRand() < cech[iplab]/std::pow(atomicWeight,0.42) )
348 {
349 G4double ran = G4UniformRand();
350
351 if ( targetCode == protonCode)
352 {
353 if(ran < 0.2)
354 {
355 pv[0] = Proton;
356 pv[1] = AntiSigmaPlus;
357 }
358 else if (ran < 0.4)
359 {
360 pv[0] = AntiLambda;
361 pv[1] = Neutron;
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] = AntiSigmaZero;
376 pv[1] = Neutron;
377 }
378 }
379 else
380 {
381 pv[0] = Neutron;
382 pv[1] = AntiSigmaPlus;
383 }
384 }
385 return;
386 }
387 else if (availableEnergy <= PionPlus.getMass())
388 return;
389
390 // inelastic scattering
391
392 np = 0; nm = 0; nz = 0;
393 G4double anhl[] = {1.00, 1.00, 1.00, 1.00, 1.00, 1.00, 1.00, 1.00, 0.97, 0.88,
394 0.85, 0.81, 0.75, 0.64, 0.64, 0.55, 0.55, 0.45, 0.47, 0.40,
395 0.39, 0.36, 0.33, 0.10, 0.01};
396 G4int iplab = G4int( incidentTotalMomentum*10.);
397 if ( iplab > 9) iplab = 10 + G4int( (incidentTotalMomentum -1.)*5. );
398 if ( iplab > 14) iplab = 15 + G4int( incidentTotalMomentum -2. );
399 if ( iplab > 22) iplab = 23 + G4int( (incidentTotalMomentum -10.)/10.);
400 iplab = std::min(24, iplab);
401
402 if ( G4UniformRand() > anhl[iplab] )
403 { // non- annihilation channels
404
405 // number of total particles vs. centre of mass Energy - 2*proton mass
406
407 G4double aleab = std::log(availableEnergy);
408 G4double n = 3.62567+aleab*(0.665843+aleab*(0.336514
409 + aleab*(0.117712+0.0136912*aleab))) - 2.0;
410
411 // normalization constant for kno-distribution.
412 // calculate first the sum of all constants, check for numerical problems.
413 G4double test, dum, anpn = 0.0;
414
415 for (nt=1; nt<=numSec; nt++) {
416 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
417 dum = pi*nt/(2.0*n*n);
418 if (std::fabs(dum) < 1.0) {
419 if( test >= 1.0e-10 )anpn += dum*test;
420 } else {
421 anpn += dum*test;
422 }
423 }
424
425 G4double ran = G4UniformRand();
426 G4double excs = 0.0;
427 if( targetCode == protonCode )
428 {
429 counter = -1;
430 for( np=0; np<numSec/3; np++ )
431 {
432 for( nm=std::max(0,np-1); nm<=(np+1); nm++ )
433 {
434 for( nz=0; nz<numSec/3; nz++ )
435 {
436 if( ++counter < numMul )
437 {
438 nt = np+nm+nz;
439 if ( (nt>0) && (nt<=numSec) ) {
440 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
441 dum = (pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n);
442 if (std::fabs(dum) < 1.0) {
443 if( test >= 1.0e-10 )excs += dum*test;
444 } else {
445 excs += dum*test;
446 }
447
448 if (ran < excs) goto outOfLoop; //----------------------->
449 }
450 }
451 }
452 }
453 }
454
455 // 3 previous loops continued to the end
456 inElastic = false; // quasi-elastic scattering
457 return;
458 }
459 else
460 { // target must be a neutron
461 counter = -1;
462 for( np=0; np<numSec/3; np++ )
463 {
464 for( nm=np; nm<=(np+2); nm++ )
465 {
466 for( nz=0; nz<numSec/3; nz++ )
467 {
468 if( ++counter < numMul )
469 {
470 nt = np+nm+nz;
471 if ( (nt>0) && (nt<=numSec) ) {
472 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
473 dum = (pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n);
474 if (std::fabs(dum) < 1.0) {
475 if( test >= 1.0e-10 )excs += dum*test;
476 } else {
477 excs += dum*test;
478 }
479
480 if (ran < excs) goto outOfLoop; // -------------------------->
481 }
482 }
483 }
484 }
485 }
486 // 3 previous loops continued to the end
487 inElastic = false; // quasi-elastic scattering.
488 return;
489 }
490
491 outOfLoop: // <------------------------------------------------------------------------
492
493 ran = G4UniformRand();
494
495 if( targetCode == protonCode)
496 {
497 if( np == nm)
498 {
499 if (ran < 0.50)
500 {
501 }
502 else if (ran < 0.75)
503 {
504 pv[0] = AntiSigmaZero;
505 pv[1] = Neutron;
506 }
507 else
508 {
509 pv[0] = AntiLambda;
510 pv[1] = Neutron;
511 }
512 }
513 else if (np == (nm-1))
514 {
515 if( ran < 0.50)
516 {
517 pv[0] = AntiLambda;
518 }
519 else
520 {
521 pv[0] = AntiSigmaZero;
522 }
523 }
524 else
525 {
526 pv[1] = Neutron;
527 }
528 }
529 else
530 {
531 if( np == nm)
532 {
533 }
534 else if ( np == (nm-1))
535 {
536 if (ran < 0.5)
537 {
538 pv[1] = Proton;
539 }
540 else if (ran < 0.75)
541 {
542 pv[0] = AntiLambda;
543 }
544 else
545 {
546 pv[0] = AntiSigmaZero;
547 }
548 }
549 else
550 {
551 if (ran < 0.5)
552 {
553 pv[0] = AntiLambda;
554 pv[1] = Proton;
555 }
556 else
557 {
558 pv[0] = AntiSigmaZero;
559 pv[1] = Proton;
560 }
561 }
562 }
563 }
564 else // annihilation
565 {
566 if ( availableEnergy > 2. * PionPlus.getMass() )
567 {
568
569 G4double aleab = std::log(availableEnergy);
570 G4double n = 3.62567+aleab*(0.665843+aleab*(0.336514
571 + aleab*(0.117712+0.0136912*aleab))) - 2.0;
572
573 // normalization constant for kno-distribution.
574 // calculate first the sum of all constants, check for numerical problems.
575 G4double test, dum, anpn = 0.0;
576
577 for (nt=2; nt<=numSec; nt++) {
578 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
579 dum = pi*nt/(2.0*n*n);
580 if (std::fabs(dum) < 1.0) {
581 if( test >= 1.0e-10 )anpn += dum*test;
582 } else {
583 anpn += dum*test;
584 }
585 }
586
587 G4double ran = G4UniformRand();
588 G4double excs = 0.0;
589 if( targetCode == protonCode )
590 {
591 counter = -1;
592 for( np=1; np<numSec/3; np++ )
593 {
594 nm = np;
595 for( nz=0; nz<numSec/3; nz++ )
596 {
597 if( ++counter < numMulAn )
598 {
599 nt = np+nm+nz;
600 if ( (nt>1) && (nt<=numSec) ) {
601 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
602 dum = (pi/anpn)*nt*protmulAn[counter]*protnormAn[nt-1]/(2.0*n*n);
603 if (std::fabs(dum) < 1.0) {
604 if( test >= 1.0e-10 )excs += dum*test;
605 } else {
606 excs += dum*test;
607 }
608
609 if (ran < excs) goto outOfLoopAn; //----------------------->
610 }
611 }
612 }
613 }
614 // 3 previous loops continued to the end
615 inElastic = false; // quasi-elastic scattering
616 return;
617 }
618 else
619 { // target must be a neutron
620 counter = -1;
621 for( np=0; np<numSec/3; np++ )
622 {
623 nm = np+1;
624 for( nz=0; nz<numSec/3; nz++ )
625 {
626 if( ++counter < numMulAn )
627 {
628 nt = np+nm+nz;
629 if ( (nt>1) && (nt<=numSec) ) {
630 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
631 dum = (pi/anpn)*nt*neutmulAn[counter]*neutnormAn[nt-1]/(2.0*n*n);
632 if (std::fabs(dum) < 1.0) {
633 if( test >= 1.0e-10 )excs += dum*test;
634 } else {
635 excs += dum*test;
636 }
637
638 if (ran < excs) goto outOfLoopAn; // -------------------------->
639 }
640 }
641 }
642 }
643 inElastic = false; // quasi-elastic scattering.
644 return;
645 }
646 outOfLoopAn: // <----------------------------------------
647 vecLen = 0;
648 }
649 }
650
651 nt = np + nm + nz;
652 while ( nt > 0)
653 {
654 G4double ran = G4UniformRand();
655 if ( ran < (G4double)np/nt)
656 {
657 if( np > 0 )
658 { pv[vecLen++] = PionPlus;
659 np--;
660 }
661 }
662 else if ( ran < (G4double)(np+nm)/nt)
663 {
664 if( nm > 0 )
665 {
666 pv[vecLen++] = PionMinus;
667 nm--;
668 }
669 }
670 else
671 {
672 if( nz > 0 )
673 {
674 pv[vecLen++] = PionZero;
675 nz--;
676 }
677 }
678 nt = np + nm + nz;
679 }
680 if (verboseLevel > 1)
681 {
682 G4cout << "Particles produced: " ;
683 G4cout << pv[0].getName() << " " ;
684 G4cout << pv[1].getName() << " " ;
685 for (i=2; i < vecLen; i++)
686 {
687 G4cout << pv[i].getName() << " " ;
688 }
689 G4cout << G4endl;
690 }
691 return;
692 }
693
694
695
696
697
698
699
700
701
702
Note: See TracBrowser for help on using the repository browser.