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

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

geant4 tag 9.4

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