source: trunk/source/processes/hadronic/models/high_energy/src/G4HEKaonMinusInelastic.cc

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

geant4 tag 9.4

File size: 20.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: G4HEKaonMinusInelastic.cc,v 1.16 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 is 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 "G4HEKaonMinusInelastic.hh"
43
44G4HadFinalState*
45G4HEKaonMinusInelastic::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 << "GHEKaonMinusInelastic: incident energy < 1 GeV" << G4endl;
65
66  if (verboseLevel > 1) {
67    G4cout << "G4HEKaonMinusInelastic::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
114  vecLength = 0;
115
116  if (verboseLevel > 1)
117    G4cout << "ApplyYourself: CallFirstIntInCascade for particle "
118           << incidentCode << G4endl;
119
120  G4bool successful = false; 
121   
122  FirstIntInCasKaonMinus(inElastic, availableEnergy, pv, vecLength,
123                         incidentParticle, targetParticle);
124
125  if (verboseLevel > 1)
126    G4cout << "ApplyYourself::StrangeParticlePairProduction" << G4endl;
127
128  if ((vecLength > 0) && (availableEnergy > 1.)) 
129    StrangeParticlePairProduction(availableEnergy, centerOfMassEnergy,
130                                  pv, vecLength,
131                                  incidentParticle, targetParticle);
132
133  HighEnergyCascading(successful, pv, vecLength,
134                      excitationEnergyGNP, excitationEnergyDTA,
135                      incidentParticle, targetParticle,
136                      atomicWeight, atomicNumber);
137  if (!successful)
138    HighEnergyClusterProduction(successful, pv, vecLength,
139                                excitationEnergyGNP, excitationEnergyDTA,
140                                incidentParticle, targetParticle,
141                                atomicWeight, atomicNumber);
142  if (!successful) 
143    MediumEnergyCascading(successful, pv, vecLength, 
144                          excitationEnergyGNP, excitationEnergyDTA, 
145                          incidentParticle, targetParticle,
146                          atomicWeight, atomicNumber);
147
148  if (!successful)
149    MediumEnergyClusterProduction(successful, pv, vecLength,
150                                  excitationEnergyGNP, excitationEnergyDTA,       
151                                  incidentParticle, targetParticle,
152                                  atomicWeight, atomicNumber);
153  if (!successful)
154    QuasiElasticScattering(successful, pv, vecLength,
155                           excitationEnergyGNP, excitationEnergyDTA,
156                           incidentParticle, targetParticle, 
157                           atomicWeight, atomicNumber);
158  if (!successful)
159    ElasticScattering(successful, pv, vecLength,
160                      incidentParticle,   
161                      atomicWeight, atomicNumber);
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
174G4HEKaonMinusInelastic::FirstIntInCasKaonMinus(G4bool& inElastic,
175                                               const G4double availableEnergy,
176                                               G4HEVector pv[],
177                                               G4int& vecLen,
178                                               const G4HEVector& incidentParticle,
179                                               const G4HEVector& targetParticle)
180
181// Kaon- undergoes interaction with nucleon within a nucleus.  Check if it is
182// energetically possible to produce pions/kaons.  In not, assume nuclear excitation
183// occurs and input particle is degraded in energy. No other particles are produced.
184// If reaction is possible, find the correct number of pions/protons/neutrons
185// produced using an interpolation to multiplicity data.  Replace some pions or
186// protons/neutrons by kaons or strange baryons according to the average
187// multiplicity per inelastic reaction.
188{
189  static const G4double expxu = std::log(MAXFLOAT); // upper bound for arg. of exp
190  static const G4double expxl = -expxu;         // lower bound for arg. of exp
191
192  static const G4double protb = 0.7;
193  static const G4double neutb = 0.7;
194  static const G4double     c = 1.25;
195
196  static const G4int numMul = 1200;
197  static const G4int numSec = 60;
198
199  G4int neutronCode = Neutron.getCode();
200  G4int protonCode = Proton.getCode();
201  G4int kaonMinusCode = KaonMinus.getCode();
202  G4int kaonZeroCode = KaonZero.getCode();
203  G4int antiKaonZeroCode = AntiKaonZero.getCode(); 
204
205  G4int targetCode = targetParticle.getCode();
206  G4double incidentTotalMomentum = incidentParticle.getTotalMomentum();
207
208  static G4bool first = true;
209  static G4double protmul[numMul], protnorm[numSec];  // proton constants
210  static G4double neutmul[numMul], neutnorm[numSec];  // neutron constants
211
212  // misc. local variables
213  // np = number of pi+,  nm = number of pi-,  nz = number of pi0
214
215  G4int i, counter, nt, np, nm, nz;
216
217   if( first ) 
218     {       // compute normalization constants, this will only be done once
219       first = false;
220       for( i=0; i<numMul; i++ )protmul[i]  = 0.0;
221       for( i=0; i<numSec; i++ )protnorm[i] = 0.0;
222       counter = -1;
223       for( np=0; np<(numSec/3); np++ ) 
224          {
225            for( nm=Imax(0,np-1); nm<=np+1; nm++ ) 
226               {
227                 for( nz=0; nz<numSec/3; nz++ ) 
228                    {
229                      if( ++counter < numMul ) 
230                        {
231                          nt = np+nm+nz;
232                          if( (nt>0) && (nt<=numSec) ) 
233                            {
234                              protmul[counter] =
235                                    pmltpc(np,nm,nz,nt,protb,c) ;
236                              protnorm[nt-1] += protmul[counter];
237                            }
238                        }
239                    }
240               }
241          }
242       for( i=0; i<numMul; i++ )neutmul[i]  = 0.0;
243       for( i=0; i<numSec; i++ )neutnorm[i] = 0.0;
244       counter = -1;
245       for( np=0; np<numSec/3; np++ ) 
246          {
247            for( nm=np; nm<=(np+2); 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                               neutmul[counter] =
257                                      pmltpc(np,nm,nz,nt,neutb,c);
258                               neutnorm[nt-1] += neutmul[counter];
259                            }
260                        }
261                    }
262               }
263          }
264       for( i=0; i<numSec; i++ ) 
265          {
266            if( protnorm[i] > 0.0 )protnorm[i] = 1.0/protnorm[i];
267            if( neutnorm[i] > 0.0 )neutnorm[i] = 1.0/neutnorm[i];
268          }
269     }                                          // end of initialization
270
271         
272                                              // initialize the first two places
273                                              // the same as beam and target                                   
274   pv[0] = incidentParticle;
275   pv[1] = targetParticle;
276   vecLen = 2;
277
278   if (!inElastic || (availableEnergy <= PionPlus.getMass())) 
279      return;
280
281                                       
282//                                            inelastic scattering
283
284   np = 0, nm = 0, nz = 0;
285   G4double cech[] = { 1., 1., 1., 0.70, 0.60, 0.55, 0.35, 0.25, 0.18, 0.15};
286   G4int iplab = G4int( incidentTotalMomentum*5.);
287   if( (iplab < 10) && (G4UniformRand() < cech[iplab])) 
288     {
289       G4int     iplab = Imin(19, G4int( incidentTotalMomentum*5.));
290       G4double cnk0[] = {0.17, 0.18, 0.17, 0.24, 0.26, 0.20, 0.22, 0.21, 0.34, 0.45,
291                          0.58, 0.55, 0.36, 0.29, 0.29, 0.32, 0.32, 0.33, 0.33, 0.33};
292       if( G4UniformRand() < cnk0[iplab] )
293         {
294           if( targetCode == protonCode )
295             {
296               pv[0] = AntiKaonZero;
297               pv[1] = Neutron;
298               return;
299             }
300           else
301             {
302               return;
303             }
304         }
305       G4double ran = G4UniformRand();
306       if( targetCode == protonCode )                    // target is a proton
307         {
308           if( ran < 0.25 )
309             { 
310               pv[0] = PionMinus;
311               pv[1] = SigmaPlus;
312             } 
313           else if (ran < 0.50)
314             {
315               pv[0] = PionZero;
316               pv[1] = SigmaZero;
317             }
318           else if (ran < 0.75)
319             { 
320               pv[0] = PionPlus;
321               pv[1] = SigmaMinus;
322             }
323           else
324             {
325               pv[0] = PionZero;
326               pv[1] = Lambda;
327             }
328         } 
329       else 
330         {                                               // target is a neutron
331           if( ran < 0.25 )
332             { 
333             } 
334           else if (ran < 0.50)
335             {
336               pv[0] = PionMinus;
337               pv[1] = SigmaZero;
338             }
339           else if (ran < 0.75)
340             { 
341               pv[0] = PionZero;
342               pv[1] = SigmaMinus;
343             }
344           else
345             {
346               pv[0] = PionMinus;
347               pv[1] = Lambda;
348             }
349         }
350       return;
351     }
352   else
353     {
354//                    number of total particles vs. centre of mass Energy - 2*proton mass
355   
356       G4double aleab = std::log(availableEnergy);
357       G4double n     = 3.62567+aleab*(0.665843+aleab*(0.336514
358                    + aleab*(0.117712+0.0136912*aleab))) - 2.0;
359   
360//                    normalization constant for kno-distribution.
361//                    calculate first the sum of all constants, check for numerical problems.   
362       G4double test, dum, anpn = 0.0;
363
364       for (nt=1; nt<=numSec; nt++) {
365         test = std::exp( Amin( expxu, Amax( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
366         dum = pi*nt/(2.0*n*n);
367         if (std::fabs(dum) < 1.0) { 
368           if( test >= 1.0e-10 )anpn += dum*test;
369         } else { 
370           anpn += dum*test;
371         }
372       }
373   
374       G4double ran = G4UniformRand();
375       G4double excs = 0.0;
376       if( targetCode == protonCode ) 
377         {
378           counter = -1;
379           for( np=0; np<numSec/3; np++ ) 
380              {
381                for( nm=Imax(0,np-1); nm<=np+1; nm++ ) 
382                   {
383                     for (nz=0; nz<numSec/3; nz++) {
384                       if (++counter < numMul) {
385                         nt = np+nm+nz;
386                         if ( (nt>0) && (nt<=numSec) ) {
387                           test = std::exp( Amin( expxu, Amax( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
388                           dum = (pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n);
389                           if (std::fabs(dum) < 1.0) { 
390                             if( test >= 1.0e-10 )excs += dum*test;
391                           } else { 
392                             excs += dum*test;
393                           }
394
395                           if (ran < excs) goto outOfLoop;      //----------------------->
396                         }   
397                       }   
398                     }     
399                   }                                                                                 
400              }
401       
402                                  // 3 previous loops continued to the end
403           inElastic = false;     // quasi-elastic scattering   
404           return;
405         }
406       else   
407         {                                         // target must be a neutron
408           counter = -1;
409           for( np=0; np<numSec/3; np++ ) 
410              {
411                for( nm=np; nm<=(np+2); nm++ ) 
412                   {
413                     for (nz=0; nz<numSec/3; nz++) {
414                       if (++counter < numMul) {
415                         nt = np+nm+nz;
416                         if ( (nt>=1) && (nt<=numSec) ) {
417                           test = std::exp( Amin( expxu, Amax( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
418                           dum = (pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n);
419                           if (std::fabs(dum) < 1.0) { 
420                             if( test >= 1.0e-10 )excs += dum*test;
421                           } else { 
422                             excs += dum*test;
423                           }
424                           if (ran < excs) goto outOfLoop;       // -------------------------->
425                         }
426                       }
427                     }
428                   }
429              }
430                                 // 3 previous loops continued to the end
431           inElastic = false;    // quasi-elastic scattering.
432           return;
433         }
434     } 
435   outOfLoop:           //  <---------------------------------------------   
436   
437   if( targetCode == protonCode)
438     {
439       if( np == (1+nm))
440         {
441           pv[1] = Neutron;
442         }
443       else if (np == nm)
444         {
445           if( G4UniformRand() < 0.75)
446             {
447             }
448           else
449             {
450               pv[0] = AntiKaonZero;
451               pv[1] = Neutron;
452             }
453         }
454       else     
455         {
456           pv[0] = AntiKaonZero;
457         } 
458     } 
459   else
460     {
461       if( np == (nm-1))
462         {
463           if( G4UniformRand() < 0.5)
464             {
465               pv[1] = Proton;
466             }
467           else
468             {
469               pv[0] = AntiKaonZero;
470             }
471         } 
472       else if ( np == nm)
473         {
474         }
475       else
476         {
477           pv[0] = AntiKaonZero;
478           pv[1] = Proton;
479         }
480     }     
481
482   if( G4UniformRand() < 0.5 )
483     {
484       if(    (    (pv[0].getCode() == kaonMinusCode)
485                && (pv[1].getCode() == neutronCode)  )
486           || (    (pv[0].getCode() == kaonZeroCode)
487                && (pv[1].getCode() == protonCode)   )
488           || (    (pv[0].getCode() == antiKaonZeroCode)
489                && (pv[1].getCode() == protonCode)   )   )
490         {
491           G4double ran = G4UniformRand();
492           if( pv[1].getCode() == protonCode)
493             { 
494               if(ran < 0.68)
495                 {
496                   pv[0] = PionPlus;
497                   pv[1] = Lambda;
498                 }
499               else if (ran < 0.84)
500                 {
501                   pv[0] = PionZero;
502                   pv[1] = SigmaPlus;
503                 }
504               else
505                 {
506                   pv[0] = PionPlus;
507                   pv[1] = SigmaZero;
508                 }
509             }
510           else
511             {
512               if(ran < 0.68)
513                 {
514                   pv[0] = PionMinus;
515                   pv[1] = Lambda;
516                 }
517               else if (ran < 0.84)
518                 {
519                   pv[0] = PionMinus;
520                   pv[1] = SigmaZero;
521                 }
522               else
523                 {
524                   pv[0] = PionZero;
525                   pv[1] = SigmaMinus;
526                 }
527             }
528         } 
529       else
530         {
531           G4double ran = G4UniformRand();
532           if (ran < 0.67)
533              {
534                pv[0] = PionZero;
535                pv[1] = Lambda;
536              }
537           else if (ran < 0.78)
538              {
539                pv[0] = PionMinus;
540                pv[1] = SigmaPlus;
541              }
542           else if (ran < 0.89)
543              {
544                pv[0] = PionZero;
545                pv[1] = SigmaZero;
546              }
547           else
548              {
549                pv[0] = PionPlus;
550                pv[1] = SigmaMinus;
551              }
552         }
553     }
554               
555
556   nt = np + nm + nz;
557   while ( nt > 0)
558       {
559         G4double ran = G4UniformRand();
560         if ( ran < (G4double)np/nt)
561            { 
562              if( np > 0 ) 
563                { pv[vecLen++] = PionPlus;
564                  np--;
565                }
566            }
567         else if ( ran < (G4double)(np+nm)/nt)
568            {   
569              if( nm > 0 )
570                { 
571                  pv[vecLen++] = PionMinus;
572                  nm--;
573                }
574            }
575         else
576            {
577              if( nz > 0 )
578                { 
579                  pv[vecLen++] = PionZero;
580                  nz--;
581                }
582            }
583         nt = np + nm + nz;
584       } 
585   if (verboseLevel > 1)
586      {
587        G4cout << "Particles produced: " ;
588        G4cout << pv[0].getName() << " " ;
589        G4cout << pv[1].getName() << " " ;
590        for (i=2; i < vecLen; i++)   
591            { 
592              G4cout << pv[i].getName() << " " ;
593            }
594         G4cout << G4endl;
595      }
596   return;
597 }
598
599
600
601
602
603
604
605
606
Note: See TracBrowser for help on using the repository browser.