source: trunk/source/processes/hadronic/models/high_energy/src/G4HEXiMinusInelastic.cc@ 1201

Last change on this file since 1201 was 1196, checked in by garnier, 16 years ago

update CVS release candidate geant4.9.3.01

File size: 20.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: G4HEXiMinusInelastic.cc,v 1.15 2008/03/17 20:49:17 dennis Exp $
28// GEANT4 tag $Name: geant4-09-03-cand-01 $
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 "G4HEXiMinusInelastic.hh"
46
47G4HadFinalState * G4HEXiMinusInelastic::
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 A = targetNucleus.GetN();
54 const G4double Z = targetNucleus.GetZ();
55 G4HEVector incidentParticle(aParticle);
56
57 G4double atomicNumber = Z;
58 G4double atomicWeight = A;
59
60 G4int incidentCode = incidentParticle.getCode();
61 G4double incidentMass = incidentParticle.getMass();
62 G4double incidentTotalEnergy = incidentParticle.getEnergy();
63 G4double incidentTotalMomentum = incidentParticle.getTotalMomentum();
64 G4double incidentKineticEnergy = incidentTotalEnergy - incidentMass;
65
66 if(incidentKineticEnergy < 1.)
67 {
68 G4cout << "GHEXiMinusInelastic: incident energy < 1 GeV" << G4endl;
69 }
70 if(verboseLevel > 1)
71 {
72 G4cout << "G4HEXiMinusInelastic::ApplyYourself" << G4endl;
73 G4cout << "incident particle " << incidentParticle.getName()
74 << "mass " << incidentMass
75 << "kinetic energy " << incidentKineticEnergy
76 << G4endl;
77 G4cout << "target material with (A,Z) = ("
78 << atomicWeight << "," << atomicNumber << ")" << G4endl;
79 }
80
81 G4double inelasticity = NuclearInelasticity(incidentKineticEnergy,
82 atomicWeight, atomicNumber);
83 if(verboseLevel > 1)
84 G4cout << "nuclear inelasticity = " << inelasticity << G4endl;
85
86 incidentKineticEnergy -= inelasticity;
87
88 G4double excitationEnergyGNP = 0.;
89 G4double excitationEnergyDTA = 0.;
90
91 G4double excitation = NuclearExcitation(incidentKineticEnergy,
92 atomicWeight, atomicNumber,
93 excitationEnergyGNP,
94 excitationEnergyDTA);
95 if(verboseLevel > 1)
96 G4cout << "nuclear excitation = " << excitation << excitationEnergyGNP
97 << excitationEnergyDTA << G4endl;
98
99
100 incidentKineticEnergy -= excitation;
101 incidentTotalEnergy = incidentKineticEnergy + incidentMass;
102 incidentTotalMomentum = std::sqrt( (incidentTotalEnergy-incidentMass)
103 *(incidentTotalEnergy+incidentMass));
104
105
106 G4HEVector targetParticle;
107 if(G4UniformRand() < atomicNumber/atomicWeight)
108 {
109 targetParticle.setDefinition("Proton");
110 }
111 else
112 {
113 targetParticle.setDefinition("Neutron");
114 }
115
116 G4double targetMass = targetParticle.getMass();
117 G4double centerOfMassEnergy = std::sqrt( incidentMass*incidentMass + targetMass*targetMass
118 + 2.0*targetMass*incidentTotalEnergy);
119 G4double availableEnergy = centerOfMassEnergy - targetMass - incidentMass;
120
121 // this was the meaning of inElastic in the
122 // original Gheisha stand-alone version.
123// G4bool inElastic = InElasticCrossSectionInFirstInt
124// (availableEnergy, incidentCode, incidentTotalMomentum);
125 // by unknown reasons, it has been replaced
126 // to the following code in Geant???
127 G4bool inElastic = true;
128// if (G4UniformRand() < elasticCrossSection/totalCrossSection) inElastic = false;
129
130 vecLength = 0;
131
132 if(verboseLevel > 1)
133 G4cout << "ApplyYourself: CallFirstIntInCascade for particle "
134 << incidentCode << G4endl;
135
136 G4bool successful = false;
137
138 if(inElastic || (!inElastic && atomicWeight < 1.5))
139 {
140 FirstIntInCasXiMinus(inElastic, availableEnergy, pv, vecLength,
141 incidentParticle, targetParticle, atomicWeight);
142
143 if(verboseLevel > 1)
144 G4cout << "ApplyYourself::StrangeParticlePairProduction" << G4endl;
145
146
147 if ((vecLength > 0) && (availableEnergy > 1.))
148 StrangeParticlePairProduction( availableEnergy, centerOfMassEnergy,
149 pv, vecLength,
150 incidentParticle, targetParticle);
151 HighEnergyCascading( successful, pv, vecLength,
152 excitationEnergyGNP, excitationEnergyDTA,
153 incidentParticle, targetParticle,
154 atomicWeight, atomicNumber);
155 if (!successful)
156 HighEnergyClusterProduction( successful, pv, vecLength,
157 excitationEnergyGNP, excitationEnergyDTA,
158 incidentParticle, targetParticle,
159 atomicWeight, atomicNumber);
160 if (!successful)
161 MediumEnergyCascading( successful, pv, vecLength,
162 excitationEnergyGNP, excitationEnergyDTA,
163 incidentParticle, targetParticle,
164 atomicWeight, atomicNumber);
165
166 if (!successful)
167 MediumEnergyClusterProduction( successful, pv, vecLength,
168 excitationEnergyGNP, excitationEnergyDTA,
169 incidentParticle, targetParticle,
170 atomicWeight, atomicNumber);
171 if (!successful)
172 QuasiElasticScattering( successful, pv, vecLength,
173 excitationEnergyGNP, excitationEnergyDTA,
174 incidentParticle, targetParticle,
175 atomicWeight, atomicNumber);
176 }
177 if (!successful)
178 {
179 ElasticScattering( successful, pv, vecLength,
180 incidentParticle,
181 atomicWeight, atomicNumber);
182 }
183
184 if (!successful)
185 {
186 G4cout << "GHEInelasticInteraction::ApplyYourself fails to produce final state particles" << G4endl;
187 }
188 FillParticleChange(pv, vecLength);
189 delete [] pv;
190 theParticleChange.SetStatusChange(stopAndKill);
191 return & theParticleChange;
192 }
193
194void
195G4HEXiMinusInelastic::FirstIntInCasXiMinus( G4bool &inElastic,
196 const G4double availableEnergy,
197 G4HEVector pv[],
198 G4int &vecLen,
199 G4HEVector incidentParticle,
200 G4HEVector targetParticle,
201 const G4double atomicWeight)
202
203// Xi0 undergoes interaction with nucleon within a nucleus. Check if it is
204// energetically possible to produce pions/kaons. In not, assume nuclear excitation
205// occurs and input particle is degraded in energy. No other particles are produced.
206// If reaction is possible, find the correct number of pions/protons/neutrons
207// produced using an interpolation to multiplicity data. Replace some pions or
208// protons/neutrons by kaons or strange baryons according to the average
209// multiplicity per inelastic reaction.
210
211 {
212 static const G4double expxu = std::log(MAXFLOAT); // upper bound for arg. of exp
213 static const G4double expxl = -expxu; // lower bound for arg. of exp
214
215 static const G4double protb = 0.7;
216 static const G4double neutb = 0.7;
217 static const G4double c = 1.25;
218
219 static const G4int numMul = 1200;
220 static const G4int numSec = 60;
221
222 G4int neutronCode = Neutron.getCode();
223 G4int protonCode = Proton.getCode();
224
225 G4int targetCode = targetParticle.getCode();
226// G4double incidentMass = incidentParticle.getMass();
227// G4double incidentEnergy = incidentParticle.getEnergy();
228 G4double incidentTotalMomentum = incidentParticle.getTotalMomentum();
229
230 static G4bool first = true;
231 static G4double protmul[numMul], protnorm[numSec]; // proton constants
232 static G4double neutmul[numMul], neutnorm[numSec]; // neutron constants
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 } // end of initialization
290
291
292 // initialize the first two places
293 // the same as beam and target
294 pv[0] = incidentParticle;
295 pv[1] = targetParticle;
296 vecLen = 2;
297
298 if( !inElastic )
299 { // quasi-elastic scattering, no pions produced
300 G4double cech[] = {0.50, 0.45, 0.40, 0.35, 0.30, 0.25, 0.06, 0.04, 0.005, 0.};
301 G4int iplab = G4int( std::min( 9.0, incidentTotalMomentum*2.5 ) );
302 if( G4UniformRand() < cech[iplab]/std::pow(atomicWeight,0.42) )
303 {
304 G4double ran = G4UniformRand();
305 if( targetCode == neutronCode)
306 {
307 if (ran < 0.2)
308 {
309 pv[0] = SigmaMinus;
310 pv[1] = SigmaZero;
311 }
312 else if (ran < 0.4)
313 {
314 pv[0] = SigmaZero;
315 pv[1] = SigmaMinus;
316 }
317 else if (ran < 0.6)
318 {
319 pv[0] = SigmaMinus;
320 pv[1] = Lambda;
321 }
322 else if (ran < 0.8)
323 {
324 pv[0] = Lambda;
325 pv[1] = SigmaMinus;
326 }
327 else
328 {
329 pv[0] = Neutron;
330 pv[1] = XiMinus;
331 }
332 }
333 else
334 {
335 if (ran < 0.2)
336 {
337 pv[0] = SigmaZero;
338 pv[1] = SigmaZero;
339 }
340 else if (ran < 0.3)
341 {
342 pv[0] = Lambda;
343 pv[1] = Lambda;
344 }
345 else if (ran < 0.4)
346 {
347 pv[0] = SigmaZero;
348 pv[1] = Lambda;
349 }
350 else if (ran < 0.5)
351 {
352 pv[0] = Lambda;
353 pv[1] = SigmaZero;
354 }
355 else if (ran < 0.7)
356 {
357 pv[0] = SigmaZero;
358 pv[1] = Neutron;
359 }
360 else if (ran < 0.8)
361 {
362 pv[0] = Neutron;
363 pv[1] = SigmaZero;
364 }
365 else
366 {
367 pv[0] = Proton;
368 pv[1] = XiMinus;
369 }
370 }
371 }
372 return;
373 }
374 else if (availableEnergy <= PionPlus.getMass())
375 return;
376
377 // inelastic scattering
378
379 np = 0; nm = 0; nz = 0;
380
381// number of total particles vs. centre of mass Energy - 2*proton mass
382
383 G4double aleab = std::log(availableEnergy);
384 G4double n = 3.62567+aleab*(0.665843+aleab*(0.336514
385 + aleab*(0.117712+0.0136912*aleab))) - 2.0;
386
387// normalization constant for kno-distribution.
388// calculate first the sum of all constants, check for numerical problems.
389 G4double test, dum, anpn = 0.0;
390
391 for (nt=1; nt<=numSec; nt++) {
392 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
393 dum = pi*nt/(2.0*n*n);
394 if (std::fabs(dum) < 1.0) {
395 if (test >= 1.0e-10) anpn += dum*test;
396 } else {
397 anpn += dum*test;
398 }
399 }
400
401 G4double ran = G4UniformRand();
402 G4double excs = 0.0;
403 if( targetCode == protonCode )
404 {
405 counter = -1;
406 for (np=0; np<numSec/3; np++) {
407 for (nm=std::max(0,np-1); nm<=(np+1); nm++) {
408 for (nz=0; nz<numSec/3; nz++) {
409 if (++counter < numMul) {
410 nt = np+nm+nz;
411 if ( (nt>0) && (nt<=numSec) ) {
412 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
413 dum = (pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n);
414 if (std::fabs(dum) < 1.0) {
415 if (test >= 1.0e-10) excs += dum*test;
416 } else {
417 excs += dum*test;
418 }
419 if (ran < excs) goto outOfLoop; //----------------------->
420 }
421 }
422 }
423 }
424 }
425
426 // 3 previous loops continued to the end
427 inElastic = false; // quasi-elastic scattering
428 return;
429 }
430 else
431 { // target must be a neutron
432 counter = -1;
433 for (np=0; np<numSec/3; np++) {
434 for (nm=np; nm<=(np+2); nm++) {
435 for( nz=0; nz<numSec/3; nz++) {
436 if (++counter < numMul) {
437 nt = np+nm+nz;
438 if ( (nt>=1) && (nt<=numSec) ) {
439 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
440 dum = (pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n);
441 if (std::fabs(dum) < 1.0) {
442 if( test >= 1.0e-10 )excs += dum*test;
443 } else {
444 excs += dum*test;
445 }
446 if (ran < excs) goto outOfLoop; // --------------------->
447 }
448 }
449 }
450 }
451 }
452 // 3 previous loops continued to the end
453 inElastic = false; // quasi-elastic scattering.
454 return;
455 }
456
457 outOfLoop: // <---------------------------------------------------
458
459 // in the following we do not consider
460 ran = G4UniformRand(); // strangeness transfer in high multiplicity
461 if( targetCode == neutronCode) // events. YK combinations are added in
462 { // StrangeParticlePairProduction
463 if( np == nm)
464 {
465 }
466 else if (np == (nm-1))
467 {
468 if( ran < 0.50)
469 {
470 pv[0] = XiZero;
471 }
472 else
473 {
474 pv[1] = Proton;
475 }
476 }
477 else
478 {
479 pv[0] = XiZero;
480 pv[1] = Proton;
481 }
482 }
483 else
484 {
485 if (np == nm)
486 {
487 if (ran < 0.5)
488 {
489 }
490 else
491 {
492 pv[0] = XiZero;
493 pv[1] = Neutron;
494 }
495 }
496 else if (np == (nm+1))
497 {
498 pv[1] = Neutron;
499 }
500 else
501 {
502 pv[0] = XiZero;
503 }
504 }
505
506
507 nt = np + nm + nz;
508 while ( nt > 0)
509 {
510 G4double ran = G4UniformRand();
511 if ( ran < (G4double)np/nt)
512 {
513 if( np > 0 )
514 { pv[vecLen++] = PionPlus;
515 np--;
516 }
517 }
518 else if ( ran < (G4double)(np+nm)/nt)
519 {
520 if( nm > 0 )
521 {
522 pv[vecLen++] = PionMinus;
523 nm--;
524 }
525 }
526 else
527 {
528 if( nz > 0 )
529 {
530 pv[vecLen++] = PionZero;
531 nz--;
532 }
533 }
534 nt = np + nm + nz;
535 }
536 if (verboseLevel > 1)
537 {
538 G4cout << "Particles produced: " ;
539 G4cout << pv[0].getName() << " " ;
540 G4cout << pv[1].getName() << " " ;
541 for (i=2; i < vecLen; i++)
542 {
543 G4cout << pv[i].getName() << " " ;
544 }
545 G4cout << G4endl;
546 }
547 return;
548 }
549
550
551
552
553
554
555
556
557
558
Note: See TracBrowser for help on using the repository browser.