source: trunk/source/processes/hadronic/models/high_energy/src/G4HEProtonInelastic.cc @ 1340

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

update ti head

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