source: trunk/source/processes/hadronic/models/high_energy/src/G4HEKaonZeroInelastic.cc @ 1228

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

update geant4.9.3 tag

File size: 20.5 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: G4HEKaonZeroInelastic.cc,v 1.15 2008/03/17 20:49:17 dennis Exp $
28// GEANT4 tag $Name: geant4-09-03 $
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 "G4HEKaonZeroInelastic.hh"
46
47G4HadFinalState *  G4HEKaonZeroInelastic::
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 << "GHEKaonZeroInelastic: incident energy < 1 GeV" << G4endl;;
69      }
70    if(verboseLevel > 1)
71      {
72        G4cout << "G4HEKaonZeroInelastic::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        FirstIntInCasKaonZero(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
195G4HEKaonZeroInelastic::FirstIntInCasKaonZero( G4bool &inElastic,
196                                              const G4double availableEnergy,
197                                              G4HEVector pv[],
198                                              G4int &vecLen,
199                                              G4HEVector incidentParticle,
200                                              G4HEVector targetParticle,
201                                              const G4double atomicWeight)
202
203// Kaon0 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] =
257                                    pmltpc(np,nm,nz,nt,protb,c) ;
258                              protnorm[nt-1] += protmul[counter];
259                            }
260                        }
261                    }
262               }
263          }
264       for( i=0; i<numMul; i++ )neutmul[i]  = 0.0;
265       for( i=0; i<numSec; i++ )neutnorm[i] = 0.0;
266       counter = -1;
267       for( np=0; np<numSec/3; np++ ) 
268          {
269            for( nm=np; nm<=(np+2); nm++ ) 
270               {
271                 for( nz=0; nz<numSec/3; nz++ ) 
272                    {
273                      if( ++counter < numMul ) 
274                        {
275                          nt = np+nm+nz;
276                          if( (nt>0) && (nt<=numSec) ) 
277                            {
278                               neutmul[counter] =
279                                      pmltpc(np,nm,nz,nt,neutb,c);
280                               neutnorm[nt-1] += neutmul[counter];
281                            }
282                        }
283                    }
284               }
285          }
286       for( i=0; i<numSec; i++ ) 
287          {
288            if( protnorm[i] > 0.0 )protnorm[i] = 1.0/protnorm[i];
289            if( neutnorm[i] > 0.0 )neutnorm[i] = 1.0/neutnorm[i];
290          }
291     }                                          // end of initialization
292
293         
294                                              // initialize the first two places
295                                              // the same as beam and target                                   
296   pv[0] = incidentParticle;
297   pv[1] = targetParticle;
298   vecLen = 2;
299
300   if( !inElastic ) 
301     {                                     // quasi-elastic scattering, no pions produced
302       if( targetCode == protonCode ) 
303         {
304           G4double cech[] = {0.33,0.27,0.29,0.31,0.27,0.18,0.13,0.10,0.09,0.07};
305           G4int iplab = G4int( std::min( 9.0, incidentTotalMomentum*5. ) );
306           if( G4UniformRand() < cech[iplab]/std::pow(atomicWeight,0.42) ) 
307             {                             // charge exchange  K+ n -> K0 p
308               pv[0] = KaonPlus;
309               pv[1] = Neutron;
310             }
311         }
312       return;
313     }
314   else if (availableEnergy <= PionPlus.getMass())
315       return;
316
317                                            //  inelastic scattering
318
319   np = 0, nm = 0, nz = 0;
320   G4double eab = availableEnergy;
321   G4int ieab = G4int( eab*5.0 );
322   
323   G4double supp[] = {0., 0.4, 0.55, 0.65, 0.75, 0.82, 0.86, 0.90, 0.94, 0.98};
324   if( (ieab <= 9) && (G4UniformRand() >=  supp[ieab]) ) 
325     {
326//                                            suppress high multiplicity events at low momentum
327//                                            only one additional pion will be produced
328       G4double w0, wp, wm, wt, ran;
329       if( targetCode == neutronCode )                    // target is a neutron
330         {
331           w0 = - sqr(1.+protb)/(2.*c*c);
332           w0 = std::exp(w0);
333           wm = - sqr(-1.+protb)/(2.*c*c);
334           wm = std::exp(wm);
335           w0 = w0/2.;
336           wm = wm*1.5;
337           if( G4UniformRand() < w0/(w0+wm) )              { np = 0; nm = 0; nz = 1; }
338           else 
339             { np = 0; nm = 1; nz = 0; }       
340         } 
341       else 
342         {                                               // target is a proton
343           w0 = -sqr(1.+neutb)/(2.*c*c);
344           wp = w0 = std::exp(w0);
345           wm = -sqr(-1.+neutb)/(2.*c*c);
346           wm = std::exp(wm);
347           wt = w0+wp+wm;
348           wp = w0+wp;
349           ran = G4UniformRand();
350           if( ran < w0/wt)
351             { np = 0; nm = 0; nz = 1; }       
352           else if( ran < wp/wt)
353             { np = 1; nm = 0; nz = 0; }       
354           else
355             { np = 0; nm = 1; nz = 0; }       
356         }
357     }
358   else
359     {
360//                    number of total particles vs. centre of mass Energy - 2*proton mass
361   
362       G4double aleab = std::log(availableEnergy);
363       G4double n     = 3.62567+aleab*(0.665843+aleab*(0.336514
364                    + aleab*(0.117712+0.0136912*aleab))) - 2.0;
365   
366//                    normalization constant for kno-distribution.
367//                    calculate first the sum of all constants, check for numerical problems.   
368       G4double test, dum, anpn = 0.0;
369
370       for (nt=1; nt<=numSec; nt++) {
371         test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
372         dum = pi*nt/(2.0*n*n);
373         if (std::fabs(dum) < 1.0) { 
374           if( test >= 1.0e-10 )anpn += dum*test;
375         } else { 
376           anpn += dum*test;
377         }
378       }
379   
380       G4double ran = G4UniformRand();
381       G4double excs = 0.0;
382       if( targetCode == protonCode ) 
383         {
384           counter = -1;
385           for( np=0; np<numSec/3; np++ ) 
386              {
387                for( nm=std::max(0,np-1); nm<=(np+1); nm++ ) 
388                   {
389                     for (nz=0; nz<numSec/3; nz++) {
390                       if (++counter < numMul) {
391                         nt = np+nm+nz;
392                         if ( (nt>0) && (nt<=numSec) ) {
393                           test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
394                           dum = (pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n);
395                           if (std::fabs(dum) < 1.0) { 
396                             if( test >= 1.0e-10 )excs += dum*test;
397                           } else { 
398                             excs += dum*test;
399                           }
400                           if (ran < excs) goto outOfLoop;      //----------------------->
401                         }   
402                       }   
403                     }     
404                   }                                                                                 
405              }
406       
407                                     // 3 previous loops continued to the end
408           inElastic = false;                 // quasi-elastic scattering   
409           return;
410         }
411       else   
412         {                                         // target must be a neutron
413           counter = -1;
414           for( np=0; np<numSec/3; np++ ) 
415              {
416                for( nm=np; nm<=(np+2); nm++ ) 
417                   {
418                     for (nz=0; nz<numSec/3; nz++) {
419                       if (++counter < numMul) {
420                         nt = np+nm+nz;
421                         if ( (nt>=1) && (nt<=numSec) ) {
422                           test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
423                           dum = (pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n);
424                           if (std::fabs(dum) < 1.0) { 
425                             if( test >= 1.0e-10 )excs += dum*test;
426                           } else { 
427                             excs += dum*test;
428                           }
429                           if (ran < excs) goto outOfLoop;       // -------------------------->
430                         }
431                       }
432                     }
433                   }
434              }
435                                                  // 3 previous loops continued to the end
436           inElastic = false;                     // quasi-elastic scattering.
437           return;
438         }
439     } 
440   outOfLoop:           //  <-----------------------------------------------   
441   
442   if( targetCode == neutronCode)
443     {
444       if( np == nm)
445         {
446         }
447       else if (np == (nm-1))
448         {
449           if( G4UniformRand() < 0.5)
450             {
451               pv[0] = KaonPlus;
452             }
453           else
454             {
455               pv[1] = Proton;
456             }
457         }
458       else     
459         {
460           pv[0] = KaonPlus;
461           pv[1] = Proton;
462         } 
463     } 
464   else
465     {
466       if( np == nm )
467         {
468           if( G4UniformRand() < 0.25)
469             {
470               pv[0] = KaonPlus;
471               pv[1] = Neutron;
472             }
473           else
474             {
475             }
476         } 
477       else if ( np == (nm+1))
478         {
479           pv[1] = Neutron;
480         }
481       else
482         {
483           pv[0] = KaonPlus;
484         }
485     }     
486
487
488   nt = np + nm + nz;
489   while ( nt > 0)
490       {
491         G4double ran = G4UniformRand();
492         if ( ran < (G4double)np/nt)
493            { 
494              if( np > 0 ) 
495                { pv[vecLen++] = PionPlus;
496                  np--;
497                }
498            }
499         else if ( ran < (G4double)(np+nm)/nt)
500            {   
501              if( nm > 0 )
502                { 
503                  pv[vecLen++] = PionMinus;
504                  nm--;
505                }
506            }
507         else
508            {
509              if( nz > 0 )
510                { 
511                  pv[vecLen++] = PionZero;
512                  nz--;
513                }
514            }
515         nt = np + nm + nz;
516       } 
517   if (verboseLevel > 1)
518      {
519        G4cout << "Particles produced: " ;
520        G4cout << pv[0].getName() << " " ;
521        G4cout << pv[1].getName() << " " ;
522        for (i=2; i < vecLen; i++)   
523            { 
524              G4cout << pv[i].getName() << " " ;
525            }
526         G4cout << G4endl;
527      }
528   return; }
529
530
531
532
533
534
535
536
537
Note: See TracBrowser for help on using the repository browser.