source: trunk/source/processes/hadronic/models/rpg/src/G4RPGKMinusInelastic.cc

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

update ti head

File size: 19.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// $Id: G4RPGKMinusInelastic.cc,v 1.1 2007/07/18 21:04:20 dennis Exp $
27// GEANT4 tag $Name: geant4-09-03-ref-09 $
28//
29 
30#include "G4RPGKMinusInelastic.hh"
31#include "Randomize.hh"
32
33G4HadFinalState*
34G4RPGKMinusInelastic::ApplyYourself( const G4HadProjectile &aTrack,
35                                     G4Nucleus &targetNucleus )
36{
37  const G4HadProjectile *originalIncident = &aTrack;
38  if (originalIncident->GetKineticEnergy()<= 0.1*MeV) 
39  {
40    theParticleChange.SetStatusChange(isAlive);
41    theParticleChange.SetEnergyChange(aTrack.GetKineticEnergy());
42    theParticleChange.SetMomentumChange(aTrack.Get4Momentum().vect().unit()); 
43    return &theParticleChange;     
44  }
45
46  // create the target particle
47   
48  G4DynamicParticle *originalTarget = targetNucleus.ReturnTargetParticle();
49  G4ReactionProduct targetParticle( originalTarget->GetDefinition() );
50   
51  if( verboseLevel > 1 )
52  {
53    const G4Material *targetMaterial = aTrack.GetMaterial();
54    G4cout << "G4RPGKMinusInelastic::ApplyYourself called" << G4endl;
55    G4cout << "kinetic energy = " << originalIncident->GetKineticEnergy() << "MeV, ";
56    G4cout << "target material = " << targetMaterial->GetName() << ", ";
57    G4cout << "target particle = " << originalTarget->GetDefinition()->GetParticleName()
58           << G4endl;
59  }
60  G4ReactionProduct currentParticle( const_cast<G4ParticleDefinition *>(originalIncident->GetDefinition()) );
61  currentParticle.SetMomentum( originalIncident->Get4Momentum().vect() );
62  currentParticle.SetKineticEnergy( originalIncident->GetKineticEnergy() );
63   
64  // Fermi motion and evaporation
65  // As of Geant3, the Fermi energy calculation had not been Done
66   
67  G4double ek = originalIncident->GetKineticEnergy();
68  G4double amas = originalIncident->GetDefinition()->GetPDGMass();
69   
70  G4double tkin = targetNucleus.Cinema( ek );
71  ek += tkin;
72  currentParticle.SetKineticEnergy( ek );
73  G4double et = ek + amas;
74  G4double p = std::sqrt( std::abs((et-amas)*(et+amas)) );
75  G4double pp = currentParticle.GetMomentum().mag();
76  if( pp > 0.0 )
77  {
78    G4ThreeVector momentum = currentParticle.GetMomentum();
79    currentParticle.SetMomentum( momentum * (p/pp) );
80  }
81   
82  // calculate black track energies
83   
84  tkin = targetNucleus.EvaporationEffects( ek );
85  ek -= tkin;
86  currentParticle.SetKineticEnergy( ek );
87  et = ek + amas;
88  p = std::sqrt( std::abs((et-amas)*(et+amas)) );
89  pp = currentParticle.GetMomentum().mag();
90  if( pp > 0.0 )
91  {
92    G4ThreeVector momentum = currentParticle.GetMomentum();
93    currentParticle.SetMomentum( momentum * (p/pp) );
94  }
95
96  G4ReactionProduct modifiedOriginal = currentParticle;
97
98  currentParticle.SetSide( 1 ); // incident always goes in forward hemisphere
99  targetParticle.SetSide( -1 );  // target always goes in backward hemisphere
100  G4bool incidentHasChanged = false;
101  G4bool targetHasChanged = false;
102  G4bool quasiElastic = false;
103  G4FastVector<G4ReactionProduct,GHADLISTSIZE> vec;  // vec will contain the secondary particles
104  G4int vecLen = 0;
105  vec.Initialize( 0 );
106   
107  const G4double cutOff = 0.1*MeV;
108  if( currentParticle.GetKineticEnergy() > cutOff )
109    Cascade( vec, vecLen,
110             originalIncident, currentParticle, targetParticle,
111             incidentHasChanged, targetHasChanged, quasiElastic );
112   
113  CalculateMomenta( vec, vecLen,
114                    originalIncident, originalTarget, modifiedOriginal,
115                    targetNucleus, currentParticle, targetParticle,
116                    incidentHasChanged, targetHasChanged, quasiElastic );
117   
118  SetUpChange( vec, vecLen,
119               currentParticle, targetParticle,
120               incidentHasChanged );
121   
122  delete originalTarget;
123  return &theParticleChange;   
124}
125
126 
127void G4RPGKMinusInelastic::Cascade(
128   G4FastVector<G4ReactionProduct,GHADLISTSIZE> &vec,
129   G4int& vecLen,
130   const G4HadProjectile *originalIncident,
131   G4ReactionProduct &currentParticle,
132   G4ReactionProduct &targetParticle,
133   G4bool &incidentHasChanged,
134   G4bool &targetHasChanged,
135   G4bool &quasiElastic )
136{
137    // Derived from H. Fesefeldt's original FORTRAN code CASKM
138    //
139    // K- undergoes interaction with nucleon within a nucleus.  Check if it is
140    // energetically possible to produce pions/kaons.  In not, assume nuclear excitation
141    // occurs and input particle is degraded in energy. No other particles are produced.
142    // If reaction is possible, find the correct number of pions/protons/neutrons
143    // produced using an interpolation to multiplicity data.  Replace some pions or
144    // protons/neutrons by kaons or strange baryons according to the average
145    // multiplicity per Inelastic reaction.
146    //
147    const G4double mOriginal = originalIncident->GetDefinition()->GetPDGMass();
148    const G4double etOriginal = originalIncident->GetTotalEnergy();
149    const G4double pOriginal = originalIncident->GetTotalMomentum();
150    const G4double targetMass = targetParticle.GetMass();
151    G4double centerofmassEnergy = std::sqrt( mOriginal*mOriginal +
152                                        targetMass*targetMass +
153                                        2.0*targetMass*etOriginal );
154    G4double availableEnergy = centerofmassEnergy-(targetMass+mOriginal);
155   
156    static G4bool first = true;
157    const G4int numMul = 1200;
158    const G4int numSec = 60;
159    static G4double protmul[numMul], protnorm[numSec]; // proton constants
160    static G4double neutmul[numMul], neutnorm[numSec]; // neutron constants
161    // np = number of pi+, nm = number of pi-, nz = number of pi0
162    G4int nt(0), np(0), nm(0), nz(0);
163    const G4double c = 1.25;   
164    const G4double b[] = { 0.70, 0.70 };
165    if( first )       // compute normalization constants, this will only be Done once
166    {
167      first = false;
168      G4int i;
169      for( i=0; i<numMul; ++i )protmul[i] = 0.0;
170      for( i=0; i<numSec; ++i )protnorm[i] = 0.0;
171      G4int counter = -1;
172      for( np=0; np<(numSec/3); ++np )
173      {
174        for( nm=std::max(0,np-1); nm<=(np+1); ++nm )
175        {
176          for( nz=0; nz<numSec/3; ++nz )
177          {
178            if( ++counter < numMul )
179            {
180              nt = np+nm+nz;
181              if( (nt>0) && (nt<=numSec) )
182              {
183                protmul[counter] = Pmltpc(np,nm,nz,nt,b[0],c);
184                protnorm[nt-1] += protmul[counter];
185              }
186            }
187          }
188        }
189      }
190      for( i=0; i<numMul; ++i )neutmul[i] = 0.0;
191      for( i=0; i<numSec; ++i )neutnorm[i] = 0.0;
192      counter = -1;
193      for( np=0; np<numSec/3; ++np )
194      {
195        for( nm=np; nm<=(np+2); ++nm )
196        {
197          for( nz=0; nz<numSec/3; ++nz )
198          {
199            if( ++counter < numMul )
200            {
201              nt = np+nm+nz;
202              if( (nt>0) && (nt<=numSec) )
203              {
204                neutmul[counter] = Pmltpc(np,nm,nz,nt,b[1],c);
205                neutnorm[nt-1] += neutmul[counter];
206              }
207            }
208          }
209        }
210      }
211      for( i=0; i<numSec; ++i )
212      {
213        if( protnorm[i] > 0.0 )protnorm[i] = 1.0/protnorm[i];
214        if( neutnorm[i] > 0.0 )neutnorm[i] = 1.0/neutnorm[i];
215      }
216    }   // end of initialization
217   
218    const G4double expxu = 82.;           // upper bound for arg. of exp
219    const G4double expxl = -expxu;        // lower bound for arg. of exp
220    G4ParticleDefinition *aKaonMinus = G4KaonMinus::KaonMinus();
221    G4ParticleDefinition *aKaonZS = G4KaonZeroShort::KaonZeroShort();
222    G4ParticleDefinition *aKaonZL = G4KaonZeroLong::KaonZeroLong();
223    G4ParticleDefinition *aNeutron = G4Neutron::Neutron();
224    G4ParticleDefinition *aProton = G4Proton::Proton();
225    G4ParticleDefinition *aPiPlus = G4PionPlus::PionPlus();
226    G4ParticleDefinition *aPiMinus = G4PionMinus::PionMinus();
227    G4ParticleDefinition *aPiZero = G4PionZero::PionZero();
228    G4ParticleDefinition *aLambda = G4Lambda::Lambda();
229    G4ParticleDefinition *aSigmaPlus = G4SigmaPlus::SigmaPlus();
230    G4ParticleDefinition *aSigmaMinus = G4SigmaMinus::SigmaMinus();
231    G4ParticleDefinition *aSigmaZero = G4SigmaZero::SigmaZero();
232    const G4double cech[] = {1.,1.,1.,0.70,0.60,0.55,0.35,0.25,0.18,0.15};
233    G4int iplab = G4int(std::min( 9.0, pOriginal/GeV*5.0 ));
234    if( (pOriginal <= 2.0*GeV) && (G4UniformRand() < cech[iplab]) )
235    {
236      np = nm = nz = nt = 0;
237      iplab = G4int(std::min( 19.0, pOriginal/GeV*10.0 ));
238      const G4double cnk0[] = {0.17,0.18,0.17,0.24,0.26,0.20,0.22,0.21,0.34,0.45,
239                               0.58,0.55,0.36,0.29,0.29,0.32,0.32,0.33,0.33,0.33};
240      if( G4UniformRand() <= cnk0[iplab] )
241      {
242        quasiElastic = true;
243        if( targetParticle.GetDefinition() == aProton )
244        {
245          currentParticle.SetDefinitionAndUpdateE( aKaonZL );
246          incidentHasChanged = true;
247          targetParticle.SetDefinitionAndUpdateE( aNeutron );
248          targetHasChanged = true;
249        }
250      }
251      else  // random number > cnk0[iplab]
252      {
253        G4double ran = G4UniformRand();
254        if( ran < 0.25 )         // k- p --> pi- s+
255        {
256          if( targetParticle.GetDefinition() == aProton )
257          {
258            currentParticle.SetDefinitionAndUpdateE( aPiMinus );
259            targetParticle.SetDefinitionAndUpdateE( aSigmaPlus );
260            incidentHasChanged = true;
261            targetHasChanged = true;
262          }
263        }
264        else if( ran < 0.50 )  // k- p --> pi0 s0  or  k- n --> pi- s0
265        {
266          if( targetParticle.GetDefinition() == aNeutron )
267            currentParticle.SetDefinitionAndUpdateE( aPiMinus );
268          else
269            currentParticle.SetDefinitionAndUpdateE( aPiZero );
270          targetParticle.SetDefinitionAndUpdateE( aSigmaZero );
271          incidentHasChanged = true;
272          targetHasChanged = true;
273        }
274        else if( ran < 0.75 )  // k- p --> pi+ s-  or  k- n --> pi0 s-
275        {
276          if( targetParticle.GetDefinition() == aNeutron )
277            currentParticle.SetDefinitionAndUpdateE( aPiZero );
278          else
279            currentParticle.SetDefinitionAndUpdateE( aPiPlus );
280          targetParticle.SetDefinitionAndUpdateE( aSigmaMinus );
281          incidentHasChanged = true;
282          targetHasChanged = true;
283        }
284        else                   // k- p --> pi0 L  or  k- n --> pi- L
285        {
286          if( targetParticle.GetDefinition() == aNeutron )
287            currentParticle.SetDefinitionAndUpdateE( aPiMinus );
288          else
289            currentParticle.SetDefinitionAndUpdateE( aPiZero );
290          targetParticle.SetDefinitionAndUpdateE( aLambda );
291          incidentHasChanged = true;
292          targetHasChanged = true;
293        }
294      }
295    }
296    else  // (pOriginal > 2.0*GeV) || (random number >= cech[iplab])
297    {
298      if( availableEnergy < aPiPlus->GetPDGMass() )
299      {               // not energetically possible to produce pion(s)
300        quasiElastic = true;
301        return;
302      }
303      G4double n, anpn;
304      GetNormalizationConstant( availableEnergy, n, anpn );
305      G4double ran = G4UniformRand();
306      G4double dum, test, excs = 0.0;
307      if( targetParticle.GetDefinition() == aProton )
308      {
309        G4int counter = -1;
310        for( np=0; np<numSec/3 && ran>=excs; ++np )
311        {
312          for( nm=std::max(0,np-1); nm<=(np+1) && ran>=excs; ++nm )
313          {
314            for( nz=0; nz<numSec/3 && ran>=excs; ++nz )
315            {
316              if( ++counter < numMul )
317              {
318                nt = np+nm+nz;
319                if( nt > 0 )
320                {
321                  test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
322                  dum = (pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n);
323                  if( std::fabs(dum) < 1.0 )
324                  {
325                    if( test >= 1.0e-10 )excs += dum*test;
326                  }
327                  else
328                    excs += dum*test;
329                }
330              }
331            }
332          }
333        }
334        if( ran >= excs )  // 3 previous loops continued to the end
335        {
336          quasiElastic = true;
337          return;
338        }
339        np--; nm--; nz--;
340        if( np == nm )
341        {
342          if( G4UniformRand() >= 0.75 )
343          {
344            currentParticle.SetDefinitionAndUpdateE( aKaonZL );
345            targetParticle.SetDefinitionAndUpdateE( aNeutron );
346            incidentHasChanged = true;
347            targetHasChanged = true;
348          }
349        }
350        else if( np == nm+1 )
351        {
352          targetParticle.SetDefinitionAndUpdateE( aNeutron );
353          targetHasChanged = true;
354        }
355        else
356        {
357          currentParticle.SetDefinitionAndUpdateE( aKaonZL );
358          incidentHasChanged = true;
359        }
360      }
361      else   // target must be a neutron
362      {
363        G4int counter = -1;
364        for( np=0; np<numSec/3 && ran>=excs; ++np )
365        {
366          for( nm=np; nm<=(np+2) && ran>=excs; ++nm )
367          {
368            for( nz=0; nz<numSec/3 && ran>=excs; ++nz )
369            {
370              if( ++counter < numMul )
371              {
372                nt = np+nm+nz;
373                if( (nt>=1) && (nt<=numSec) )
374                {
375                  test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
376                  dum = (pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n);
377                  if( std::fabs(dum) < 1.0 )
378                  {
379                    if( test >= 1.0e-10 )excs += dum*test;
380                  }
381                  else
382                    excs += dum*test;
383                }
384              }
385            }
386          }
387        }
388        if( ran >= excs )  // 3 previous loops continued to the end
389        {
390          quasiElastic = true;
391          return;
392        }
393        np--; nm--; nz--;
394        if( np == nm-1 )
395        {
396          if( G4UniformRand() < 0.5 )
397          {
398            targetParticle.SetDefinitionAndUpdateE( aProton );
399            targetHasChanged = true;
400          }
401          else
402          {
403            currentParticle.SetDefinitionAndUpdateE( aKaonZL );
404            incidentHasChanged = true;
405          }
406        }
407        else if( np != nm )
408        {
409          currentParticle.SetDefinitionAndUpdateE( aKaonZL );
410          incidentHasChanged = true;
411        }
412      }
413      if( G4UniformRand() >= 0.5 )
414      {
415        if( (currentParticle.GetDefinition() == aKaonMinus &&
416             targetParticle.GetDefinition() == aNeutron )     ||
417            (currentParticle.GetDefinition() == aKaonZL &&
418             targetParticle.GetDefinition() == aProton ) )
419        {
420          ran = G4UniformRand();
421          if( ran < 0.68 )
422          {
423            if( targetParticle.GetDefinition() == aProton )
424            {
425              currentParticle.SetDefinitionAndUpdateE( aPiPlus );
426              targetParticle.SetDefinitionAndUpdateE( aLambda );
427              incidentHasChanged = true;
428              targetHasChanged = true;
429            }
430            else
431            {
432              currentParticle.SetDefinitionAndUpdateE( aPiMinus );
433              targetParticle.SetDefinitionAndUpdateE( aLambda );
434              incidentHasChanged = true;
435              targetHasChanged = true;
436            }
437          }
438          else if( ran < 0.84 )
439          {
440            if( targetParticle.GetDefinition() == aProton )
441            {
442              currentParticle.SetDefinitionAndUpdateE( aPiZero );
443              targetParticle.SetDefinitionAndUpdateE( aSigmaPlus );
444              incidentHasChanged = true;
445              targetHasChanged = true;
446            }
447            else
448            {
449              currentParticle.SetDefinitionAndUpdateE( aPiMinus );
450              targetParticle.SetDefinitionAndUpdateE( aSigmaZero );
451              incidentHasChanged = true;
452              targetHasChanged = true;
453            }
454          }
455          else
456          {
457            if( targetParticle.GetDefinition() == aProton )
458            {
459              currentParticle.SetDefinitionAndUpdateE( aPiPlus );
460              targetParticle.SetDefinitionAndUpdateE( aSigmaZero );
461              incidentHasChanged = true;
462              targetHasChanged = true;
463            }
464            else
465            {
466              currentParticle.SetDefinitionAndUpdateE( aPiZero );
467              targetParticle.SetDefinitionAndUpdateE( aSigmaMinus );
468              incidentHasChanged = true;
469              targetHasChanged = true;
470            }
471          }
472        }
473        else  // ( current != aKaonMinus || target != aNeutron ) &&
474              // ( current != aKaonZL    || target != aProton  )
475        {
476          ran = G4UniformRand();
477          if( ran < 0.67 )
478          {
479            currentParticle.SetDefinitionAndUpdateE( aPiZero );
480            targetParticle.SetDefinitionAndUpdateE( aLambda );
481            incidentHasChanged = true;
482            targetHasChanged = true;
483          }
484          else if( ran < 0.78 )
485          {
486            currentParticle.SetDefinitionAndUpdateE( aPiMinus );
487            targetParticle.SetDefinitionAndUpdateE( aSigmaPlus );
488            incidentHasChanged = true;
489            targetHasChanged = true;
490          }
491          else if( ran < 0.89 )
492          {
493            currentParticle.SetDefinitionAndUpdateE( aPiZero );
494            targetParticle.SetDefinitionAndUpdateE( aSigmaZero );
495            incidentHasChanged = true;
496            targetHasChanged = true;
497          }
498          else
499          {
500            currentParticle.SetDefinitionAndUpdateE( aPiPlus );
501            targetParticle.SetDefinitionAndUpdateE( aSigmaMinus );
502            incidentHasChanged = true;
503            targetHasChanged = true;
504          }
505        }
506      }
507    }
508    if( currentParticle.GetDefinition() == aKaonZL )
509    {
510      if( G4UniformRand() >= 0.5 )
511      {
512        currentParticle.SetDefinitionAndUpdateE( aKaonZS );
513        incidentHasChanged = true;
514      }
515    }
516    if( targetParticle.GetDefinition() == aKaonZL )
517    {
518      if( G4UniformRand() >= 0.5 )
519      {
520        targetParticle.SetDefinitionAndUpdateE( aKaonZS );
521        targetHasChanged = true;
522      }
523    }
524    SetUpPions( np, nm, nz, vec, vecLen );
525    return;
526}
527
528 /* end of file */
529 
530
531
Note: See TracBrowser for help on using the repository browser.