source: trunk/source/processes/hadronic/models/high_energy/src/G4HEPionMinusInelastic.cc @ 1347

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

geant4 tag 9.4

File size: 18.4 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: G4HEPionMinusInelastic.cc,v 1.17 2010/11/29 05:44:44 dennis Exp $
27// GEANT4 tag $Name: geant4-09-04-ref-00 $
28//
29// 11-OCT-2007 F.W. Jones: fixed incorrect Imax (should be Imin) in
30//             sampling of charge exchange.
31
32
33#include "globals.hh"
34#include "G4ios.hh"
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
39// like nuclear reactions, nuclear fission without any cascading and all
40// processes for 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 "G4HEPionMinusInelastic.hh"
46
47G4HadFinalState*
48G4HEPionMinusInelastic::ApplyYourself(const G4HadProjectile& aTrack,
49                                      G4Nucleus& targetNucleus)
50{
51  G4HEVector* pv = new G4HEVector[MAXPART];
52  const G4HadProjectile* aParticle = &aTrack;
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    G4cout << "GHEPionMinusInelastic: incident energy < 1 GeV" << G4endl;
68
69  if (verboseLevel > 1) {
70    G4cout << "G4HEPionMinusInelastic::ApplyYourself" << G4endl;
71    G4cout << "incident particle " << incidentParticle.getName()
72           << "mass "              << incidentMass
73           << "kinetic energy "    << incidentKineticEnergy
74           << G4endl;
75    G4cout << "target material with (A,Z) = (" 
76           << atomicWeight << "," << atomicNumber << ")" << G4endl;
77  }
78
79  G4double inelasticity = NuclearInelasticity(incidentKineticEnergy, 
80                                              atomicWeight, atomicNumber);
81  if (verboseLevel > 1)
82    G4cout << "nuclear inelasticity = " << inelasticity << G4endl;
83
84  incidentKineticEnergy -= inelasticity;
85   
86  G4double excitationEnergyGNP = 0.;
87  G4double excitationEnergyDTA = 0.; 
88
89  G4double excitation = NuclearExcitation(incidentKineticEnergy,
90                                          atomicWeight, atomicNumber,
91                                          excitationEnergyGNP,
92                                          excitationEnergyDTA);
93  if (verboseLevel > 1)
94    G4cout << "nuclear excitation = " << excitation << excitationEnergyGNP
95           << excitationEnergyDTA << G4endl;
96
97  incidentKineticEnergy -= excitation;
98  incidentTotalEnergy = incidentKineticEnergy + incidentMass;
99  incidentTotalMomentum = std::sqrt( (incidentTotalEnergy-incidentMass)                   
100                                    *(incidentTotalEnergy+incidentMass));
101
102  G4HEVector targetParticle;
103 
104  if (G4UniformRand() < atomicNumber/atomicWeight) { 
105    targetParticle.setDefinition("Proton");
106  } else { 
107    targetParticle.setDefinition("Neutron");
108  }
109
110  G4double targetMass = targetParticle.getMass();
111  G4double centerOfMassEnergy = std::sqrt(incidentMass*incidentMass
112                                        + targetMass*targetMass
113                                        + 2.0*targetMass*incidentTotalEnergy);
114  G4double availableEnergy = centerOfMassEnergy - targetMass - incidentMass;
115
116  // The value of the inElastic flag was originally defined in the Gheisha
117  // stand-alone code as follows:
118  //   G4bool inElastic = InElasticCrossSectionInFirstInt
119  //                      (availableEnergy, incidentCode, incidentTotalMomentum); 
120  // For unknown reasons, it was replaced by the following code in Geant
121
122  G4bool inElastic = true;
123  //    if (G4UniformRand() < elasticCrossSection/totalCrossSection) inElastic = false;   
124
125  vecLength = 0;
126       
127  if (verboseLevel > 1)
128    G4cout << "ApplyYourself: CallFirstIntInCascade for particle "
129           << incidentCode << G4endl;
130
131  G4bool successful = false; 
132   
133  FirstIntInCasPionMinus(inElastic, availableEnergy, pv, vecLength,
134                         incidentParticle, targetParticle);
135
136  if (verboseLevel > 1)
137    G4cout << "ApplyYourself::StrangeParticlePairProduction" << G4endl;
138
139  if ((vecLength > 0) && (availableEnergy > 1.)) 
140    StrangeParticlePairProduction(availableEnergy, centerOfMassEnergy,
141                                  pv, vecLength,
142                                  incidentParticle, targetParticle);
143
144  HighEnergyCascading(successful, pv, vecLength,
145                      excitationEnergyGNP, excitationEnergyDTA,
146                      incidentParticle, targetParticle,
147                      atomicWeight, atomicNumber);
148  if (!successful)
149    HighEnergyClusterProduction(successful, pv, vecLength,
150                                excitationEnergyGNP, excitationEnergyDTA,
151                                incidentParticle, targetParticle,
152                                atomicWeight, atomicNumber);
153  if (!successful) 
154    MediumEnergyCascading(successful, pv, vecLength, 
155                          excitationEnergyGNP, excitationEnergyDTA, 
156                          incidentParticle, targetParticle,
157                          atomicWeight, atomicNumber);
158
159  if (!successful)
160    MediumEnergyClusterProduction(successful, pv, vecLength,
161                                  excitationEnergyGNP, excitationEnergyDTA,       
162                                  incidentParticle, targetParticle,
163                                  atomicWeight, atomicNumber);
164  if (!successful)
165    QuasiElasticScattering(successful, pv, vecLength,
166                           excitationEnergyGNP, excitationEnergyDTA,
167                           incidentParticle, targetParticle, 
168                           atomicWeight, atomicNumber);
169  if (!successful)
170    ElasticScattering(successful, pv, vecLength,
171                      incidentParticle,   
172                      atomicWeight, atomicNumber);
173
174  if (!successful)
175    G4cout << "GHEInelasticInteraction::ApplyYourself fails to produce final state particles" 
176           << G4endl;
177
178  FillParticleChange(pv,  vecLength);
179  delete [] pv;
180  theParticleChange.SetStatusChange(stopAndKill);
181  return &theParticleChange;
182}
183
184
185void
186G4HEPionMinusInelastic::FirstIntInCasPionMinus(G4bool& inElastic,
187                                               const G4double availableEnergy,
188                                               G4HEVector pv[],
189                                               G4int& vecLen,
190                                               const G4HEVector& incidentParticle,
191                                               const G4HEVector& targetParticle)
192
193// Pion- undergoes interaction with nucleon within a nucleus.  Check if it is
194// energetically possible to produce pions/kaons.  In not, assume nuclear excitation
195// occurs and input particle is degraded in energy. No other particles are produced.
196// If reaction is possible, find the correct number of pions/protons/neutrons
197// produced using an interpolation to multiplicity data.  Replace some pions or
198// protons/neutrons by kaons or strange baryons according to the average
199// multiplicity per inelastic reaction.
200{
201  static const G4double expxu = std::log(MAXFLOAT); // upper bound for arg. of exp
202  static const G4double expxl = -expxu;             // lower bound for arg. of exp
203
204  static const G4double protb = 0.7;
205  static const G4double neutb = 0.7;
206  static const G4double     c = 1.25;
207
208  static const G4int numMul = 1200;
209  static const G4int numSec = 60;
210
211  G4int protonCode  = Proton.getCode();
212
213  G4int targetCode = targetParticle.getCode();
214  G4double incidentTotalMomentum = incidentParticle.getTotalMomentum();
215
216  static G4bool first = true;
217  static G4double protmul[numMul], protnorm[numSec];  // proton constants
218  static G4double neutmul[numMul], neutnorm[numSec];  // neutron constants
219
220  // misc. local variables
221  // np = number of pi+,  nm = number of pi-,  nz = number of pi0
222
223  G4int i, counter, nt, np, nm, nz;
224
225   if( first ) 
226     {    // compute normalization constants, this will only be done once
227       first = false;
228       for( i=0; i<numMul; i++ )protmul[i]  = 0.0;
229       for( i=0; i<numSec; i++ )protnorm[i] = 0.0;
230       counter = -1;
231       for( np=0; np<(numSec/3); np++ ) 
232          {
233            for( nm=Imax(0,np-1); nm<=(np+1); nm++ ) 
234               {
235                 for( nz=0; nz<numSec/3; nz++ ) 
236                    {
237                      if( ++counter < numMul ) 
238                        {
239                          nt = np+nm+nz;
240                          if( (nt>0) && (nt<=numSec) ) 
241                            {
242                              protmul[counter] =
243                                    pmltpc(np,nm,nz,nt,protb,c) ;
244                              protnorm[nt-1] += protmul[counter];
245                            }
246                        }
247                    }
248               }
249          }
250       for( i=0; i<numMul; i++ )neutmul[i]  = 0.0;
251       for( i=0; i<numSec; i++ )neutnorm[i] = 0.0;
252       counter = -1;
253       for( np=0; np<numSec/3; np++ ) 
254          {
255            for( nm=np; nm<=(np+2); nm++ ) 
256               {
257                 for( nz=0; nz<numSec/3; nz++ ) 
258                    {
259                      if( ++counter < numMul ) 
260                        {
261                          nt = np+nm+nz;
262                          if( (nt>0) && (nt<=numSec) ) 
263                            {
264                               neutmul[counter] =
265                                      pmltpc(np,nm,nz,nt,neutb,c);
266                               neutnorm[nt-1] += neutmul[counter];
267                            }
268                        }
269                    }
270               }
271          }
272       for( i=0; i<numSec; i++ ) 
273          {
274            if( protnorm[i] > 0.0 )protnorm[i] = 1.0/protnorm[i];
275            if( neutnorm[i] > 0.0 )neutnorm[i] = 1.0/neutnorm[i];
276          }
277     }                                          // end of initialization
278
279         
280   // initialize the first two places
281   // the same as beam and target                                   
282   pv[0] = incidentParticle;
283   pv[1] = targetParticle;
284   vecLen = 2;
285
286   if (!inElastic || (availableEnergy <= PionPlus.getMass())) 
287      return;
288
289                                       
290   // inelastic scattering
291
292   np = 0, nm = 0, nz = 0;
293   G4double eab = availableEnergy;
294   G4int ieab = G4int( eab*5.0 );
295   
296   G4double supp[] = {0., 0.4, 0.55, 0.65, 0.75, 0.82, 0.86, 0.90, 0.94, 0.98};
297   if( (ieab <= 9) && (G4UniformRand() >= supp[ieab]) ) 
298     {
299       // suppress high multiplicity events at low momentum
300       // only one additional pion will be produced
301       G4double cech[] = {1., 0.95, 0.79, 0.32, 0.19, 0.16, 0.14, 0.12, 0.10, 0.08};
302       G4int iplab = Imin(9, G4int( incidentTotalMomentum*5.));
303       if( G4UniformRand() < cech[iplab] )
304         {
305           if( targetCode == protonCode )
306             {
307               pv[0] = PionZero;
308               pv[1] = Neutron;
309               return;
310             }
311           else
312             {
313               return;
314             }
315         }
316
317       G4double w0, wp, wm, wt, ran;
318       if( targetCode == protonCode )                    // target is a proton
319         {
320           w0 = - sqr(1.+protb)/(2.*c*c);
321           wp = w0 = std::exp(w0);
322           wm = - sqr(-1.+protb)/(2.*c*c);
323           wm = std::exp(wm);
324           wp *= 10.;
325           wt = w0+wp+wm;
326           wp = w0+wp;
327           ran = G4UniformRand();
328           if( ran < w0/wt ) 
329             { np = 0; nm = 0; nz = 1; }
330           else if ( ran < wp/wt )
331             { np = 1; nm = 0; nz = 0; }
332           else
333             { np = 0; nm = 1; nz = 0; }       
334         } 
335       else 
336         {                                           // target is a neutron
337           w0 = -sqr(1.+neutb)/(2.*c*c);
338           wm = -sqr(-1.+neutb)/(2.*c*c);
339           w0 = std::exp(w0);
340           wm = std::exp(wm);
341           if( G4UniformRand() < w0/(w0+wm) )
342             { np = 0; nm = 0; nz = 1; }       
343           else
344             { np = 0; nm = 1; nz = 0; }       
345         }
346     }
347   else
348     {
349       // number of total particles vs. centre of mass Energy - 2*proton mass
350   
351       G4double aleab = std::log(availableEnergy);
352       G4double n     = 3.62567+aleab*(0.665843+aleab*(0.336514
353                    + aleab*(0.117712+0.0136912*aleab))) - 2.0;
354   
355       // normalization constant for kno-distribution.
356       // calculate first the sum of all constants, check for numerical problems.   
357       G4double test, dum, anpn = 0.0;
358
359       for (nt=1; nt<=numSec; nt++) {
360         test = std::exp( Amin( expxu, Amax( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
361         dum = pi*nt/(2.0*n*n);
362         if (std::fabs(dum) < 1.0) { 
363           if( test >= 1.0e-10 )anpn += dum*test;
364         } else { 
365           anpn += dum*test;
366         }
367       }
368   
369       G4double ran = G4UniformRand();
370       G4double excs = 0.0;
371       if( targetCode == protonCode ) 
372         {
373           counter = -1;
374           for (np=0; np<numSec/3; np++) {
375             for (nm=Imax(0,np-1); nm<=(np+1); nm++) {
376               for (nz=0; nz<numSec/3; nz++) {
377                 if (++counter < numMul) {
378                   nt = np+nm+nz;
379                   if ( (nt>0) && (nt<=numSec) ) {
380                     test = std::exp( Amin( expxu, Amax( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
381                     dum = (pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n);
382                     if (std::fabs(dum) < 1.0) { 
383                       if( test >= 1.0e-10 )excs += dum*test;
384                     } else { 
385                       excs += dum*test;
386                     }
387                     if (ran < excs) goto outOfLoop;      //----------------------->
388                   }
389                 }
390               }
391             }
392           }
393       
394           // 3 previous loops continued to the end
395           inElastic = false;                 // quasi-elastic scattering   
396           return;
397         }
398       else   
399         {                                         // target must be a neutron
400           counter = -1;
401           for (np=0; np<numSec/3; np++) {
402             for (nm=np; nm<=(np+2); nm++) {
403               for (nz=0; nz<numSec/3; nz++) {
404                 if (++counter < numMul) {
405                   nt = np+nm+nz;
406                   if ( (nt>=1) && (nt<=numSec) ) {
407                     test = std::exp( Amin( expxu, Amax( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
408                     dum = (pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n);
409                     if (std::fabs(dum) < 1.0) { 
410                       if( test >= 1.0e-10 )excs += dum*test;
411                     } else { 
412                       excs += dum*test;
413                     }
414                     if (ran < excs) goto outOfLoop;     // ----------------->
415                   }
416                 }
417               }
418             }
419           }
420           // 3 previous loops continued to the end
421           inElastic = false;                     // quasi-elastic scattering.
422           return;
423         }
424     } 
425   outOfLoop:           //  <-----------------------------------------------   
426   
427   if( targetCode == protonCode)
428     {
429       if( np == (1+nm))
430         {
431           pv[1] = Neutron;
432         }
433       else if (np == nm)
434         {
435           if( G4UniformRand() < 0.75)
436             {
437             }
438           else
439             {
440               pv[0] = PionZero;
441               pv[1] = Neutron;
442             }
443         }
444       else     
445         {
446           pv[0] = PionZero;
447         } 
448     } 
449   else
450     {
451       if( np == (nm-1))
452         {
453           if( G4UniformRand() < 0.5)
454             {
455               pv[1] = Proton;
456             }
457           else
458             {
459               pv[0] = PionZero;
460             }
461         } 
462       else if ( np == nm)
463         {
464         }
465       else
466         {
467           pv[0] = PionZero;
468           pv[1] = Proton;
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].getName() << " " ;
506        G4cout << pv[1].getName() << " " ;
507        for (i=2; i < vecLen; i++)   
508            { 
509              G4cout << pv[i].getName() << " " ;
510            }
511         G4cout << G4endl;
512      }
513   return;
514 }
515
516
517
518
519
520
521
522
523
Note: See TracBrowser for help on using the repository browser.