source: trunk/source/processes/hadronic/models/high_energy/src/G4HESigmaMinusInelastic.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.0 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: G4HESigmaMinusInelastic.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 "G4HESigmaMinusInelastic.hh"
43
44G4HadFinalState*
45G4HESigmaMinusInelastic::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 << "GHESigmaMinusInelastic: incident energy < 1 GeV" << G4endl;
65
66 if (verboseLevel > 1) {
67 G4cout << "G4HESigmaMinusInelastic::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 FirstIntInCasSigmaMinus(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
166 FillParticleChange(pv, vecLength);
167 delete [] pv;
168 theParticleChange.SetStatusChange(stopAndKill);
169 return &theParticleChange;
170}
171
172
173void
174G4HESigmaMinusInelastic::FirstIntInCasSigmaMinus(G4bool& inElastic,
175 const G4double availableEnergy,
176 G4HEVector pv[],
177 G4int& vecLen,
178 const G4HEVector& incidentParticle,
179 const G4HEVector& targetParticle,
180 const G4double atomicWeight)
181
182// Sigma- undergoes interaction with nucleon within a nucleus. Check if it is
183// energetically possible to produce pions/kaons. In not, assume nuclear excitation
184// occurs and input particle is degraded in energy. No other particles are produced.
185// If reaction is possible, find the correct number of pions/protons/neutrons
186// produced using an interpolation to multiplicity data. Replace some pions or
187// protons/neutrons by kaons or strange baryons according to the average
188// multiplicity per inelastic reaction.
189{
190 static const G4double expxu = std::log(MAXFLOAT); // upper bound for arg. of exp
191 static const G4double expxl = -expxu; // lower bound for arg. of exp
192
193 static const G4double protb = 0.7;
194 static const G4double neutb = 0.7;
195 static const G4double c = 1.25;
196
197 static const G4int numMul = 1200;
198 static const G4int numSec = 60;
199
200 G4int neutronCode = Neutron.getCode();
201 G4int protonCode = Proton.getCode();
202
203 G4int targetCode = targetParticle.getCode();
204 G4double incidentTotalMomentum = incidentParticle.getTotalMomentum();
205
206 static G4bool first = true;
207 static G4double protmul[numMul], protnorm[numSec]; // proton constants
208 static G4double neutmul[numMul], neutnorm[numSec]; // neutron constants
209
210 // misc. local variables
211 // np = number of pi+, nm = number of pi-, nz = number of pi0
212
213 G4int i, counter, nt, np, nm, nz;
214
215 if( first )
216 { // compute normalization constants, this will only be done once
217 first = false;
218 for( i=0; i<numMul; i++ )protmul[i] = 0.0;
219 for( i=0; i<numSec; i++ )protnorm[i] = 0.0;
220 counter = -1;
221 for( np=0; np<(numSec/3); np++ )
222 {
223 for( nm=std::max(0,np-1); nm<=(np+1); nm++ )
224 {
225 for( nz=0; nz<numSec/3; nz++ )
226 {
227 if( ++counter < numMul )
228 {
229 nt = np+nm+nz;
230 if( (nt>0) && (nt<=numSec) )
231 {
232 protmul[counter] = pmltpc(np,nm,nz,nt,protb,c);
233 protnorm[nt-1] += protmul[counter];
234 }
235 }
236 }
237 }
238 }
239 for( i=0; i<numMul; i++ )neutmul[i] = 0.0;
240 for( i=0; i<numSec; i++ )neutnorm[i] = 0.0;
241 counter = -1;
242 for( np=0; np<numSec/3; np++ )
243 {
244 for( nm=np; nm<=(np+2); nm++ )
245 {
246 for( nz=0; nz<numSec/3; nz++ )
247 {
248 if( ++counter < numMul )
249 {
250 nt = np+nm+nz;
251 if( (nt>0) && (nt<=numSec) )
252 {
253 neutmul[counter] = pmltpc(np,nm,nz,nt,neutb,c);
254 neutnorm[nt-1] += neutmul[counter];
255 }
256 }
257 }
258 }
259 }
260 for( i=0; i<numSec; i++ )
261 {
262 if( protnorm[i] > 0.0 )protnorm[i] = 1.0/protnorm[i];
263 if( neutnorm[i] > 0.0 )neutnorm[i] = 1.0/neutnorm[i];
264 }
265 } // end of initialization
266
267
268 // initialize the first two places
269 // the same as beam and target
270 pv[0] = incidentParticle;
271 pv[1] = targetParticle;
272 vecLen = 2;
273
274 if( !inElastic )
275 { // quasi-elastic scattering, no pions produced
276 G4double cech[] = {0.50, 0.45, 0.40, 0.35, 0.30, 0.25, 0.06, 0.04, 0.005, 0.};
277 G4int iplab = G4int( std::min( 9.0, incidentTotalMomentum*2.5 ) );
278 if( G4UniformRand() < cech[iplab]/std::pow(atomicWeight,0.42) )
279 {
280 G4double ran = G4UniformRand();
281 if( targetCode == neutronCode)
282 {
283 pv[0] = Neutron;
284 pv[1] = SigmaMinus;
285 }
286 else
287 {
288 if(ran < 0.2)
289 {
290 pv[0] = SigmaZero;
291 pv[1] = Neutron;
292 }
293 else if(ran < 0.4)
294 {
295 pv[0] = Lambda;
296 pv[1] = Neutron;
297 }
298 else if(ran < 0.6)
299 {
300 pv[0] = Proton;
301 pv[1] = SigmaMinus;
302 }
303 else if(ran < 0.8)
304 {
305 pv[0] = Neutron;
306 pv[1] = SigmaZero;
307 }
308 else
309 {
310 pv[0] = Neutron;
311 pv[1] = Lambda;
312 }
313 }
314 }
315 return;
316 }
317 else if (availableEnergy <= PionPlus.getMass())
318 return;
319
320 // inelastic scattering
321
322 np = 0, nm = 0, nz = 0;
323
324 // number of total particles vs. centre of mass Energy - 2*proton mass
325
326 G4double aleab = std::log(availableEnergy);
327 G4double n = 3.62567+aleab*(0.665843+aleab*(0.336514
328 + aleab*(0.117712+0.0136912*aleab))) - 2.0;
329
330 // normalization constant for kno-distribution.
331 // calculate first the sum of all constants, check for numerical problems.
332 G4double test, dum, anpn = 0.0;
333
334 for (nt=1; nt<=numSec; nt++) {
335 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
336 dum = pi*nt/(2.0*n*n);
337 if (std::fabs(dum) < 1.0) {
338 if (test >= 1.0e-10) anpn += dum*test;
339 } else {
340 anpn += dum*test;
341 }
342 }
343
344 G4double ran = G4UniformRand();
345 G4double excs = 0.0;
346 if( targetCode == protonCode )
347 {
348 counter = -1;
349 for (np=0; np<numSec/3; np++) {
350 for (nm=std::max(0,np-1); nm<=(np+1); nm++) {
351 for (nz=0; nz<numSec/3; nz++) {
352 if (++counter < numMul) {
353 nt = np+nm+nz;
354 if ( (nt>0) && (nt<=numSec) ) {
355 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
356 dum = (pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n);
357 if (std::fabs(dum) < 1.0) {
358 if( test >= 1.0e-10 )excs += dum*test;
359 } else {
360 excs += dum*test;
361 }
362 if (ran < excs) goto outOfLoop; //----------------------->
363 }
364 }
365 }
366 }
367 }
368
369 // 3 previous loops continued to the end
370 inElastic = false; // quasi-elastic scattering
371 return;
372 }
373 else
374 { // target must be a neutron
375 counter = -1;
376 for (np=0; np<numSec/3; np++) {
377 for (nm=np; nm<=(np+2); nm++) {
378 for (nz=0; nz<numSec/3; nz++) {
379 if (++counter < numMul) {
380 nt = np+nm+nz;
381 if ( (nt>=1) && (nt<=numSec) ) {
382 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
383 dum = (pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n);
384 if (std::fabs(dum) < 1.0) {
385 if( test >= 1.0e-10 )excs += dum*test;
386 } else {
387 excs += dum*test;
388 }
389 if (ran < excs) goto outOfLoop; // -------------------->
390 }
391 }
392 }
393 }
394 }
395 // 3 previous loops continued to the end
396 inElastic = false; // quasi-elastic scattering.
397 return;
398 }
399
400 outOfLoop: // <-------------------------------------------
401
402 ran = G4UniformRand();
403 if( targetCode == neutronCode)
404 {
405 if( np == nm)
406 {
407 }
408 else if (np == (nm-1))
409 {
410 if( ran < 0.25)
411 {
412 pv[0] = SigmaZero;
413 }
414 else if (ran < 0.5)
415 {
416 pv[0] = Lambda;
417 }
418 else
419 {
420 pv[1] = Proton;
421 }
422 }
423 else
424 {
425 if(ran < 0.5)
426 {
427 pv[0] = SigmaZero;
428 pv[1] = Proton;
429 }
430 else
431 {
432 pv[0] = Lambda;
433 pv[1] = Proton;
434 }
435 }
436 }
437 else
438 {
439 if (np == nm)
440 {
441 if (ran < 0.5)
442 {
443 }
444 else if (ran < 0.75)
445 {
446 pv[0] = SigmaZero;
447 pv[1] = Neutron;
448 }
449 else
450 {
451 pv[0] = Lambda;
452 pv[1] = Neutron;
453 }
454 }
455 else if (np == (nm+1))
456 {
457 pv[1] = Neutron;
458 }
459 else
460 {
461 if (ran < 0.5)
462 {
463 pv[0] = SigmaZero;
464 }
465 else
466 {
467 pv[0] = Lambda;
468 }
469 }
470 }
471
472
473 nt = np + nm + nz;
474 while ( nt > 0)
475 {
476 G4double ran = G4UniformRand();
477 if ( ran < (G4double)np/nt)
478 {
479 if( np > 0 )
480 { pv[vecLen++] = PionPlus;
481 np--;
482 }
483 }
484 else if ( ran < (G4double)(np+nm)/nt)
485 {
486 if( nm > 0 )
487 {
488 pv[vecLen++] = PionMinus;
489 nm--;
490 }
491 }
492 else
493 {
494 if( nz > 0 )
495 {
496 pv[vecLen++] = PionZero;
497 nz--;
498 }
499 }
500 nt = np + nm + nz;
501 }
502 if (verboseLevel > 1)
503 {
504 G4cout << "Particles produced: " ;
505 G4cout << pv[0].getCode() << " " ;
506 G4cout << pv[1].getCode() << " " ;
507 for (i=2; i < vecLen; i++)
508 {
509 G4cout << pv[i].getCode() << " " ;
510 }
511 G4cout << G4endl;
512 }
513 return;
514 }
515
516
517
518
519
520
521
522
523
524
Note: See TracBrowser for help on using the repository browser.