source: trunk/source/processes/hadronic/models/high_energy/src/G4HEXiZeroInelastic.cc @ 1348

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

geant4 tag 9.4

File size: 18.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: G4HEXiZeroInelastic.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 "G4HEXiZeroInelastic.hh"
43
44G4HadFinalState*
45G4HEXiZeroInelastic::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 << "GHEXiZeroInelastic: incident energy < 1 GeV" << G4endl;
65
66  if (verboseLevel > 1) {
67    G4cout << "G4HEXiZeroInelastic::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  FirstIntInCasXiZero(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
174G4HEXiZeroInelastic::FirstIntInCasXiZero(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// Xi0 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               if (ran < 0.2)
283                 {
284                   pv[0] = SigmaPlus;
285                   pv[1] = SigmaZero;
286                 } 
287               else if (ran < 0.4)
288                 {
289                   pv[0] = SigmaZero;
290                   pv[1] = SigmaPlus;
291                 }
292               else if (ran < 0.6)
293                 {
294                   pv[0] = SigmaPlus;
295                   pv[1] = Lambda;
296                 }
297               else if (ran < 0.8)
298                 {
299                   pv[0] = Lambda;
300                   pv[1] = SigmaPlus;
301                 }
302               else
303                 {
304                   pv[0] = Proton;
305                   pv[1] = XiZero;
306                 }         
307             }             
308           else
309             { 
310               if (ran < 0.2)
311                 {
312                   pv[0] = Neutron;
313                   pv[1] = XiZero;
314                 } 
315               else if (ran < 0.3)
316                 {
317                   pv[0] = SigmaZero;
318                   pv[1] = SigmaZero; 
319                 }
320               else if (ran < 0.4)
321                 {
322                   pv[0] = Lambda;
323                   pv[1] = Lambda;
324                 }
325               else if (ran < 0.5)
326                 {
327                   pv[0] = SigmaZero;
328                   pv[1] = Lambda;
329                 } 
330               else if (ran < 0.6)
331                 {
332                   pv[0] = Lambda;
333                   pv[1] = SigmaZero;
334                 }
335               else if (ran < 0.7)
336                 {
337                   pv[0] = SigmaPlus;
338                   pv[1] = SigmaMinus;
339                 } 
340               else if (ran < 0.8)
341                 {
342                   pv[0] = SigmaMinus;
343                   pv[1] = SigmaPlus;
344                 }
345               else if (ran < 0.9)
346                 {
347                   pv[0] = XiMinus;
348                   pv[1] = Proton;
349                 }
350               else
351                 {
352                   pv[0] = Proton;
353                   pv[1] = XiMinus;
354                 } 
355             }
356         }
357       return;
358     }
359   else if (availableEnergy <= PionPlus.getMass())
360       return;
361
362                                                  //   inelastic scattering
363
364   np = 0; nm = 0; nz = 0;
365
366   // number of total particles vs. centre of mass Energy - 2*proton mass
367   
368   G4double aleab = std::log(availableEnergy);
369   G4double n     = 3.62567+aleab*(0.665843+aleab*(0.336514
370                    + aleab*(0.117712+0.0136912*aleab))) - 2.0;
371   
372   // normalization constant for kno-distribution.
373   // calculate first the sum of all constants, check for numerical problems.             
374   G4double test, dum, anpn = 0.0;
375
376   for (nt=1; nt<=numSec; nt++) {
377     test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
378     dum = pi*nt/(2.0*n*n);
379     if (std::fabs(dum) < 1.0) {
380       if (test >= 1.0e-10) anpn += dum*test;
381     } else { 
382       anpn += dum*test;
383     }
384   }
385   
386   G4double ran = G4UniformRand();
387   G4double excs = 0.0;
388   if( targetCode == protonCode ) 
389     {
390       counter = -1;
391       for (np=0; np<numSec/3; np++) {
392         for (nm=std::max(0,np-2); nm<=np; nm++) {
393           for (nz=0; nz<numSec/3; nz++) {
394             if (++counter < numMul) {
395               nt = np+nm+nz;
396               if ( (nt>0) && (nt<=numSec) ) {
397                 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
398                 dum = (pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n);
399                 if (std::fabs(dum) < 1.0) { 
400                   if (test >= 1.0e-10) excs += dum*test;
401                 } else { 
402                   excs += dum*test;
403                 }
404                 if (ran < excs) goto outOfLoop;      //----------------------->
405               }
406             }
407           }
408         }
409       }
410       
411       // 3 previous loops continued to the end
412
413       inElastic = false;                 // quasi-elastic scattering   
414       return;
415     }
416   else   
417     {                                         // target must be a neutron
418       counter = -1;
419       for (np=0; np<numSec/3; np++) {
420         for (nm=std::max(0,np-1); nm<=(np+1); nm++) {
421           for (nz=0; nz<numSec/3; nz++) {
422             if (++counter < numMul) {
423               nt = np+nm+nz;
424               if ( (nt>=1) && (nt<=numSec) ) {
425                 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
426                 dum = (pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n);
427                 if (std::fabs(dum) < 1.0) { 
428                   if (test >= 1.0e-10) excs += dum*test;
429                 } else { 
430                   excs += dum*test;
431                 }
432                 if (ran < excs) goto outOfLoop;       // ------------------->
433               }
434             }
435           }
436         }
437       }
438       // 3 previous loops continued to the end
439
440       inElastic = false;                         // quasi-elastic scattering.
441       return;
442     }
443 
444   outOfLoop:           //  <----------------------------------------------------   
445   
446                                                   // in the following we do not consider
447   ran = G4UniformRand();                          // strangeness transfer in high multiplicity
448   if( targetCode == protonCode)                   // events. YK combinations are added in
449     {                                             // StrangeParticlePairProduction
450       if( np == nm)
451         { 
452         }
453       else if (np == (nm+1))
454         {
455           if( ran < 0.50)
456             {
457               pv[0] = XiMinus;
458             }
459           else
460             {
461               pv[1] = Neutron;
462             }
463         }
464       else
465         {
466           pv[0] = XiMinus;
467           pv[1] = Neutron;
468         }   
469     } 
470   else
471     {
472       if (np == nm)
473         {
474           if (ran < 0.5) 
475             {
476             }
477           else
478             {
479               pv[0] = XiMinus;
480               pv[1] = Proton;
481             } 
482         } 
483       else if (np == (nm-1))
484         { 
485           pv[1] = Proton;
486         } 
487       else 
488         {
489           pv[0] = XiMinus;
490         }
491     }     
492
493
494   nt = np + nm + nz;
495   while ( nt > 0)
496       {
497         G4double ran = G4UniformRand();
498         if ( ran < (G4double)np/nt)
499            { 
500              if( np > 0 ) 
501                { pv[vecLen++] = PionPlus;
502                  np--;
503                }
504            }
505         else if ( ran < (G4double)(np+nm)/nt)
506            {   
507              if( nm > 0 )
508                { 
509                  pv[vecLen++] = PionMinus;
510                  nm--;
511                }
512            }
513         else
514            {
515              if( nz > 0 )
516                { 
517                  pv[vecLen++] = PionZero;
518                  nz--;
519                }
520            }
521         nt = np + nm + nz;
522       } 
523   if (verboseLevel > 1)
524      {
525        G4cout << "Particles produced: " ;
526        G4cout << pv[0].getCode() << " " ;
527        G4cout << pv[1].getCode() << " " ;
528        for (i=2; i < vecLen; i++)   
529            { 
530              G4cout << pv[i].getCode() << " " ;
531            }
532         G4cout << G4endl;
533      }
534   return;
535 }
536
537
538
539
540
541
542
543
544
545
Note: See TracBrowser for help on using the repository browser.