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

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

update ti head

File size: 20.3 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: G4HELambdaInelastic.cc,v 1.15 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 "G4HELambdaInelastic.hh"
46
47G4HadFinalState *  G4HELambdaInelastic::
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 << "GHELambdaInelastic: incident energy < 1 GeV" << G4endl;
69      }
70    if(verboseLevel > 1)
71      {
72        G4cout << "G4HELambdaInelastic::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        FirstIntInCasLambda(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
195G4HELambdaInelastic::FirstIntInCasLambda( G4bool &inElastic,
196                                          const G4double availableEnergy,
197                                          G4HEVector pv[],
198                                          G4int &vecLen,
199                                          G4HEVector incidentParticle,
200                                          G4HEVector targetParticle,
201                                          const G4double atomicWeight)
202
203// Lambda undergoes interaction with nucleon within a nucleus.  Check if it is
204// energetically possible to produce pions/kaons.  In not, assume nuclear
205// excitation occurs and input particle is degraded in energy. No other
206// particles are produced.  If reaction is possible, find the correct number
207// of pions/protons/neutrons produced using an interpolation to multiplicity
208// data.  Replace some pions or protons/neutrons by kaons or strange baryons
209// according to the average 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-2); 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] = pmltpc(np,nm,nz,nt,protb,c);
257                              protnorm[nt-1] += protmul[counter];
258                            }
259                        }
260                    }
261               }
262          }
263       for( i=0; i<numMul; i++ )neutmul[i]  = 0.0;
264       for( i=0; i<numSec; i++ )neutnorm[i] = 0.0;
265       counter = -1;
266       for( np=0; np<numSec/3; np++ ) 
267          {
268            for( nm=std::max(0,np-1); nm<=(np+2); nm++ ) 
269               {
270                 for( nz=0; nz<numSec/3; nz++ ) 
271                    {
272                      if( ++counter < numMul ) 
273                        {
274                          nt = np+nm+nz;
275                          if( (nt>0) && (nt<=numSec) ) 
276                            {
277                               neutmul[counter] = pmltpc(np,nm,nz,nt,neutb,c);
278                               neutnorm[nt-1] += neutmul[counter];
279                            }
280                        }
281                    }
282               }
283          }
284       for( i=0; i<numSec; i++ ) 
285          {
286            if( protnorm[i] > 0.0 )protnorm[i] = 1.0/protnorm[i];
287            if( neutnorm[i] > 0.0 )neutnorm[i] = 1.0/neutnorm[i];
288          }
289     }                                      // end of initialization
290
291         
292                                            // initialize the first two places
293                                            // the same as beam and target
294   pv[0] = incidentParticle;
295   pv[1] = targetParticle;
296   vecLen = 2;
297
298   if( !inElastic ) 
299     {                       // quasi-elastic scattering, no pions produced
300       G4double cech[] = {0.50, 0.45, 0.40, 0.35, 0.30, 0.25, 0.06, 0.04, 0.005, 0.};
301       G4int iplab = G4int( std::min( 9.0, incidentTotalMomentum*2.5 ) );
302       if( G4UniformRand() < cech[iplab]/std::pow(atomicWeight,0.42) ) 
303         {                                           
304           G4double ran = G4UniformRand();
305           if( targetCode == protonCode)
306             { 
307               if( ran < 0.2)
308                 {
309                   pv[0] = SigmaPlus;
310                   pv[1] = Neutron;
311                 }
312               else if(ran < 0.4)
313                 {
314                   pv[0] = SigmaZero;
315                 }
316               else if(ran < 0.6)
317                 {
318                   pv[0] = Proton;
319                   pv[1] = Lambda;
320                 }
321               else if(ran < 0.8)
322                 {
323                   pv[0] = Proton;
324                   pv[1] = SigmaZero;
325                 }
326               else
327                 {
328                   pv[0] = Neutron;
329                   pv[1] = SigmaPlus;
330                 }
331             }             
332           else
333             { 
334               if(ran < 0.2)
335                 {
336                   pv[0] = SigmaZero;
337                 }
338               else if(ran < 0.4)
339                 {
340                   pv[0] = SigmaMinus;
341                   pv[1] = Proton;
342                 }
343               else if(ran < 0.6)
344                 {
345                   pv[0] = Neutron;
346                   pv[1] = Lambda;
347                 } 
348               else if(ran < 0.8)
349                 {
350                   pv[0] = Neutron;
351                   pv[1] = SigmaZero;
352                 }
353               else
354                 {
355                   pv[0] = Proton;
356                   pv[1] = SigmaMinus;
357                 } 
358             }
359         }
360       return;
361     }
362   else if (availableEnergy <= PionPlus.getMass())
363       return;
364
365   //   inelastic scattering
366
367   np = 0; nm = 0; nz = 0;
368
369   // number of total particles vs. centre of mass Energy - 2*proton mass
370   
371   G4double aleab = std::log(availableEnergy);
372   G4double n     = 3.62567+aleab*(0.665843+aleab*(0.336514
373                    + aleab*(0.117712+0.0136912*aleab))) - 2.0;
374   
375   // normalization constant for kno-distribution.
376   // calculate first the sum of all constants, check for numerical problems. 
377   G4double test, dum, anpn = 0.0;
378
379   for (nt=1; nt<=numSec; nt++) {
380     test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
381     dum = pi*nt/(2.0*n*n);
382     if (std::fabs(dum) < 1.0) {
383       if( test >= 1.0e-10 )anpn += dum*test;
384     } else { 
385       anpn += dum*test;
386     }
387   }
388   
389   G4double ran = G4UniformRand();
390   G4double excs = 0.0;
391   if( targetCode == protonCode ) 
392     {
393       counter = -1;
394       for (np=0; np<numSec/3; np++) {
395         for (nm=std::max(0,np-2); nm<=(np+1); nm++) {
396           for (nz=0; nz<numSec/3; nz++) {
397             if (++counter < numMul) {
398               nt = np+nm+nz;
399               if ( (nt>0) && (nt<=numSec) ) {
400                 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
401                 dum = (pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n);
402                 if (std::fabs(dum) < 1.0) { 
403                   if( test >= 1.0e-10 )excs += dum*test;
404                 } else { 
405                   excs += dum*test;
406                 }
407                 if (ran < excs) goto outOfLoop;      //--------------------->
408               }   
409             }
410           }
411         }
412       }
413       // 3 previous loops continued to the end
414
415       inElastic = false;                 // quasi-elastic scattering   
416       return;
417     }
418   else   
419     {                                         // target must be a neutron
420       counter = -1;
421       for (np=0; np<numSec/3; np++) {
422         for (nm=std::max(0,np-1); nm<=(np+2); nm++) {
423           for (nz=0; nz<numSec/3; nz++) {
424             if (++counter < numMul) {
425               nt = np+nm+nz;
426               if ( (nt>=1) && (nt<=numSec) ) {
427                 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
428                 dum = (pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n);
429                 if (std::fabs(dum) < 1.0) { 
430                   if( test >= 1.0e-10 )excs += dum*test;
431                 } else { 
432                   excs += dum*test;
433                 }
434                 if (ran < excs) goto outOfLoop;       // ------------->
435               }
436             }
437           }
438         }
439       }
440                                  // 3 previous loops continued to the end
441       inElastic = false;         // quasi-elastic scattering.
442       return;
443     }
444 
445   outOfLoop:           //  <--------------------------------------------   
446   
447   ran = G4UniformRand();
448   if( targetCode == protonCode)
449     {
450       if( np == nm)
451         { 
452           if (ran < 0.25)
453             { 
454             }
455           else if(ran < 0.5)
456             {
457               pv[0] = SigmaZero;
458             } 
459           else 
460             {
461               pv[0] = SigmaPlus;   
462               pv[1] = Neutron;
463             }
464         }
465       else if (np == (nm+1))
466         {
467           if( G4UniformRand() < 0.25)
468             {
469               pv[1] = Neutron;
470             }
471           else if(ran < 0.5)
472             {
473               pv[0] = SigmaZero;
474               pv[1] = Neutron; 
475             }
476           else
477             {
478               pv[0] = SigmaMinus;
479             }
480         }
481       else if (np == (nm-1))     
482         {
483           pv[0] = SigmaPlus;
484         } 
485       else
486         {
487           pv[0] = SigmaMinus;
488           pv[1] = Neutron;
489         } 
490     } 
491   else
492     {
493       if (np == nm)
494         {
495           if(ran < 0.5) 
496             {
497             }
498           else
499             {
500               pv[0] = SigmaMinus;
501               pv[1] = Proton;
502             }
503         } 
504       else if (np == (nm-1))
505         { 
506           if( ran < 0.25)
507             {
508               pv[1] = Proton;
509             }
510           else if(ran < 0.5)               
511             {
512               pv[0] = SigmaZero;
513               pv[1] = Proton;
514             }
515           else
516             {
517               pv[1] = SigmaPlus;
518             }
519         } 
520       else if (np == (1+nm))
521         {
522           pv[0] = SigmaMinus;
523         }
524       else
525         {
526           pv[0] = SigmaPlus; 
527           pv[1] = Proton;
528         }
529     }     
530
531
532   nt = np + nm + nz;
533   while ( nt > 0)
534       {
535         G4double ran = G4UniformRand();
536         if ( ran < (G4double)np/nt)
537            { 
538              if( np > 0 ) 
539                { pv[vecLen++] = PionPlus;
540                  np--;
541                }
542            }
543         else if ( ran < (G4double)(np+nm)/nt)
544            {   
545              if( nm > 0 )
546                { 
547                  pv[vecLen++] = PionMinus;
548                  nm--;
549                }
550            }
551         else
552            {
553              if( nz > 0 )
554                { 
555                  pv[vecLen++] = PionZero;
556                  nz--;
557                }
558            }
559         nt = np + nm + nz;
560       } 
561   if (verboseLevel > 1)
562      {
563        G4cout << "Particles produced: " ;
564        G4cout << pv[0].getName() << " " ;
565        G4cout << pv[1].getName() << " " ;
566        for (i=2; i < vecLen; i++)   
567            { 
568              G4cout << pv[i].getName() << " " ;
569            }
570         G4cout << G4endl;
571      }
572   return;
573 }
574
575
576
577
578
579
580
581
582
583
Note: See TracBrowser for help on using the repository browser.