source: trunk/source/processes/hadronic/models/low_energy/src/G4LEKaonPlusInelastic.cc @ 1007

Last change on this file since 1007 was 1007, checked in by garnier, 15 years ago

update to geant4.9.2

File size: 14.6 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: G4LEKaonPlusInelastic.cc,v 1.13 2006/06/29 20:45:01 gunter Exp $
28// GEANT4 tag $Name: geant4-09-02 $
29//
30 // Hadronic Process: Low Energy KaonPlus Inelastic Process
31 // J.L. Chuma, TRIUMF, 05-Feb-1997
32 // Last modified: 27-Mar-1997
33 // Modified by J.L.Chuma 30-Apr-97: added originalTarget for CalculateMomenta
34 
35#include "G4LEKaonPlusInelastic.hh"
36#include "Randomize.hh"
37 
38 G4HadFinalState *
39  G4LEKaonPlusInelastic::ApplyYourself( const G4HadProjectile &aTrack,
40                                        G4Nucleus &targetNucleus )
41  {
42    const G4HadProjectile *originalIncident = &aTrack;
43    if (originalIncident->GetKineticEnergy()<= 0.1*MeV) 
44    {
45      theParticleChange.SetStatusChange(isAlive);
46      theParticleChange.SetEnergyChange(aTrack.GetKineticEnergy());
47      theParticleChange.SetMomentumChange(aTrack.Get4Momentum().vect().unit()); 
48      return &theParticleChange;     
49    }
50   
51    // create the target particle
52   
53    G4DynamicParticle *originalTarget = targetNucleus.ReturnTargetParticle();
54    G4ReactionProduct targetParticle( originalTarget->GetDefinition() );
55   
56    if( verboseLevel > 1 )
57    {
58      const G4Material *targetMaterial = aTrack.GetMaterial();
59      G4cout << "G4LEKaonPlusInelastic::ApplyYourself called" << G4endl;
60      G4cout << "kinetic energy = " << originalIncident->GetKineticEnergy() << "MeV, ";
61      G4cout << "target material = " << targetMaterial->GetName() << ", ";
62      G4cout << "target particle = " << originalTarget->GetDefinition()->GetParticleName()
63           << G4endl;
64    }   
65    G4ReactionProduct currentParticle( const_cast<G4ParticleDefinition *>(originalIncident->GetDefinition()));
66    currentParticle.SetMomentum( originalIncident->Get4Momentum().vect() );
67    currentParticle.SetKineticEnergy( originalIncident->GetKineticEnergy() );
68   
69    // Fermi motion and evaporation
70    // As of Geant3, the Fermi energy calculation had not been Done
71   
72    G4double ek = originalIncident->GetKineticEnergy();
73    G4double amas = originalIncident->GetDefinition()->GetPDGMass();
74   
75    G4double tkin = targetNucleus.Cinema( ek );
76    ek += tkin;
77    currentParticle.SetKineticEnergy( ek );
78    G4double et = ek + amas;
79    G4double p = std::sqrt( std::abs((et-amas)*(et+amas)) );
80    G4double pp = currentParticle.GetMomentum().mag();
81    if( pp > 0.0 )
82    {
83      G4ThreeVector momentum = currentParticle.GetMomentum();
84      currentParticle.SetMomentum( momentum * (p/pp) );
85    }
86   
87    // calculate black track energies
88   
89    tkin = targetNucleus.EvaporationEffects( ek );
90    ek -= tkin;
91    currentParticle.SetKineticEnergy( ek );
92    et = ek + amas;
93    p = std::sqrt( std::abs((et-amas)*(et+amas)) );
94    pp = currentParticle.GetMomentum().mag();
95    if( pp > 0.0 )
96    {
97      G4ThreeVector momentum = currentParticle.GetMomentum();
98      currentParticle.SetMomentum( momentum * (p/pp) );
99    }
100   
101    G4ReactionProduct modifiedOriginal = currentParticle;
102   
103    currentParticle.SetSide( 1 ); // incident always goes in forward hemisphere
104    targetParticle.SetSide( -1 );  // target always goes in backward hemisphere
105    G4bool incidentHasChanged = false;
106    G4bool targetHasChanged = false;
107    G4bool quasiElastic = false;
108    G4FastVector<G4ReactionProduct,GHADLISTSIZE> vec;  // vec will contain the secondary particles
109    G4int vecLen = 0;
110    vec.Initialize( 0 );
111       
112    const G4double cutOff = 0.1*MeV;
113    if( currentParticle.GetKineticEnergy() > cutOff )
114      Cascade( vec, vecLen,
115               originalIncident, currentParticle, targetParticle,
116               incidentHasChanged, targetHasChanged, quasiElastic );
117   
118    CalculateMomenta( vec, vecLen,
119                      originalIncident, originalTarget, modifiedOriginal,
120                      targetNucleus, currentParticle, targetParticle,
121                      incidentHasChanged, targetHasChanged, quasiElastic );
122   
123    SetUpChange( vec, vecLen,
124                 currentParticle, targetParticle,
125                 incidentHasChanged );
126   
127    delete originalTarget;
128   
129    return &theParticleChange;   
130  }
131 
132 void
133  G4LEKaonPlusInelastic::Cascade(
134   G4FastVector<G4ReactionProduct,GHADLISTSIZE> &vec,
135   G4int &vecLen,
136   const G4HadProjectile *originalIncident,
137   G4ReactionProduct &currentParticle,
138   G4ReactionProduct &targetParticle,
139   G4bool &incidentHasChanged,
140   G4bool &targetHasChanged,
141   G4bool &quasiElastic )
142  {
143    // derived from original FORTRAN code CASKP by H. Fesefeldt (13-Sep-1987)
144    //
145    // K+ undergoes interaction with nucleon within a nucleus.  Check if it is
146    // energetically possible to produce pions/kaons.  In not, assume nuclear excitation
147    // occurs and input particle is degraded in energy. No other particles are produced.
148    // If reaction is possible, find the correct number of pions/protons/neutrons
149    // produced using an interpolation to multiplicity data.  Replace some pions or
150    // protons/neutrons by kaons or strange baryons according to the average
151    // multiplicity per Inelastic reaction.
152    //
153    const G4double mOriginal = originalIncident->GetDefinition()->GetPDGMass();
154    const G4double etOriginal = originalIncident->GetTotalEnergy();
155    const G4double targetMass = targetParticle.GetMass();
156    G4double centerofmassEnergy = std::sqrt( mOriginal*mOriginal +
157                                        targetMass*targetMass +
158                                        2.0*targetMass*etOriginal );
159    G4double availableEnergy = centerofmassEnergy-(targetMass+mOriginal);
160    if( availableEnergy < G4PionPlus::PionPlus()->GetPDGMass() )
161    {
162      quasiElastic = true;
163      return;
164    }
165    static G4bool first = true;
166    const G4int numMul = 1200;
167    const G4int numSec = 60;
168    static G4double protmul[numMul], protnorm[numSec]; // proton constants
169    static G4double neutmul[numMul], neutnorm[numSec]; // neutron constants
170    // np = number of pi+, nm = number of pi-, nz = number of pi0
171    G4int nt=0, np=0, nm=0, nz=0;
172    const G4double c = 1.25;   
173    const G4double b[] = { 0.70, 0.70 };
174    if( first )       // compute normalization constants, this will only be Done once
175    {
176      first = false;
177      G4int i;
178      for( i=0; i<numMul; ++i )protmul[i] = 0.0;
179      for( i=0; i<numSec; ++i )protnorm[i] = 0.0;
180      G4int counter = -1;
181      for( np=0; np<(numSec/3); ++np )
182      {
183        for( nm=std::max(0,np-2); nm<=np; ++nm )
184        {
185          for( nz=0; nz<numSec/3; ++nz )
186          {
187            if( ++counter < numMul )
188            {
189              nt = np+nm+nz;
190              if( nt > 0 )
191              {
192                protmul[counter] = Pmltpc(np,nm,nz,nt,b[0],c);
193                protnorm[nt-1] += protmul[counter];
194              }
195            }
196          }
197        }
198      }
199      for( i=0; i<numMul; ++i )neutmul[i] = 0.0;
200      for( i=0; i<numSec; ++i )neutnorm[i] = 0.0;
201      counter = -1;
202      for( np=0; np<numSec/3; ++np )
203      {
204        for( nm=std::max(0,np-1); nm<=(np+1); ++nm )
205        {
206          for( nz=0; nz<numSec/3; ++nz )
207          {
208            if( ++counter < numMul )
209            {
210              nt = np+nm+nz;
211              if( (nt>0) && (nt<=numSec) )
212              {
213                neutmul[counter] = Pmltpc(np,nm,nz,nt,b[1],c);
214                neutnorm[nt-1] += neutmul[counter];
215              }
216            }
217          }
218        }
219      }
220      for( i=0; i<numSec; ++i )
221      {
222        if( protnorm[i] > 0.0 )protnorm[i] = 1.0/protnorm[i];
223        if( neutnorm[i] > 0.0 )neutnorm[i] = 1.0/neutnorm[i];
224      }
225    }   // end of initialization
226       
227    const G4double expxu = 82.;           // upper bound for arg. of exp
228    const G4double expxl = -expxu;        // lower bound for arg. of exp
229    G4ParticleDefinition *aKaonZS = G4KaonZeroShort::KaonZeroShort();
230    G4ParticleDefinition *aKaonZL = G4KaonZeroLong::KaonZeroLong();
231    G4ParticleDefinition *aNeutron = G4Neutron::Neutron();
232    G4ParticleDefinition *aProton = G4Proton::Proton();
233    G4int ieab = static_cast<G4int>(availableEnergy*5.0/GeV);
234    const G4double supp[] = {0.,0.4,0.55,0.65,0.75,0.82,0.86,0.90,0.94,0.98};
235    G4double test, w0, wp, wt, wm;
236    if( (availableEnergy < 2.0*GeV) && (G4UniformRand() >= supp[ieab]) )
237    {
238      // suppress high multiplicity events at low momentum
239      // only one pion will be produced
240     
241      nm = np = nz = 0;
242      if( targetParticle.GetDefinition() == aProton )
243      {
244        test = std::exp( std::min( expxu, std::max( expxl, -sqr(1.0+b[0])/(2.0*c*c) ) ) );
245        w0 = test;
246        wp = test*2.0;       
247        if( G4UniformRand() < w0/(w0+wp) )
248          nz = 1;
249        else
250          np = 1;
251      }
252      else  // target is a neutron
253      {
254        test = std::exp( std::min( expxu, std::max( expxl, -sqr(1.0+b[1])/(2.0*c*c) ) ) );
255        w0 = test;
256        wp = test;
257        test = std::exp( std::min( expxu, std::max( expxl, -sqr(-1.0+b[1])/(2.0*c*c) ) ) );
258        wm = test;
259        wt = w0+wp+wm;
260        wp += w0;
261        G4double ran = G4UniformRand();
262        if( ran < w0/wt )
263          nz = 1;
264        else if( ran < wp/wt )
265          np = 1;
266        else
267          nm = 1;
268      }
269    }
270    else
271    {
272      G4double n, anpn;
273      GetNormalizationConstant( availableEnergy, n, anpn );
274      G4double ran = G4UniformRand();
275      G4double dum, excs = 0.0;
276      if( targetParticle.GetDefinition() == aProton )
277      {
278        G4int counter = -1;
279        for( np=0; (np<numSec/3) && (ran>=excs); ++np )
280        {
281          for( nm=std::max(0,np-2); (nm<=np) && (ran>=excs); ++nm )
282          {
283            for( nz=0; (nz<numSec/3) && (ran>=excs); ++nz )
284            {
285              if( ++counter < numMul )
286              {
287                nt = np+nm+nz;
288                if( nt > 0 )
289                {
290                  test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
291                  dum = (pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n);
292                  if( std::fabs(dum) < 1.0 )
293                  {
294                    if( test >= 1.0e-10 )excs += dum*test;
295                  }
296                  else
297                    excs += dum*test;
298                }
299              }
300            }
301          }
302        }
303        if( ran >= excs )return;  // 3 previous loops continued to the end
304        np--; nm--; nz--;
305      }
306      else  // target must be a neutron
307      {
308        G4int counter = -1;
309        for( np=0; (np<numSec/3) && (ran>=excs); ++np )
310        {
311          for( nm=std::max(0,np-1); (nm<=(np+1)) && (ran>=excs); ++nm )
312          {
313            for( nz=0; (nz<numSec/3) && (ran>=excs); ++nz )
314            {
315              if( ++counter < numMul )
316              {
317                nt = np+nm+nz;
318                if( (nt>=1) && (nt<=numSec) )
319                {
320                  test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
321                  dum = (pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n);
322                  if( std::fabs(dum) < 1.0 )
323                  {
324                    if( test >= 1.0e-10 )excs += dum*test;
325                  }
326                  else
327                    excs += dum*test;
328                }
329              }
330            }
331          }
332        }
333        if( ran >= excs )return;  // 3 previous loops continued to the end
334        np--; nm--; nz--;
335      }
336    }
337    if( targetParticle.GetDefinition() == aProton )
338    {
339      switch( np-nm )
340      {
341       case 1:
342         if( G4UniformRand() < 0.5 )
343         {
344           if( G4UniformRand() < 0.5 )
345             currentParticle.SetDefinitionAndUpdateE( aKaonZS );
346           else
347             currentParticle.SetDefinitionAndUpdateE( aKaonZL );
348           incidentHasChanged = true;
349         }
350         else
351         {
352           targetParticle.SetDefinitionAndUpdateE( aNeutron );
353           targetHasChanged = true;
354         }
355         break;
356       case 2:
357         if( G4UniformRand() < 0.5 )
358           currentParticle.SetDefinitionAndUpdateE( aKaonZS );
359         else
360           currentParticle.SetDefinitionAndUpdateE( aKaonZL );
361         incidentHasChanged = true;
362         targetParticle.SetDefinitionAndUpdateE( aNeutron );
363         incidentHasChanged = true;
364         targetHasChanged = true;
365         break;
366       default:
367         break;
368      }
369    }
370    else   // target is a neutron
371    {
372      switch( np-nm )
373      {
374       case 0:
375         if( G4UniformRand() < 0.25 )
376         {
377           if( G4UniformRand() < 0.5 )
378             currentParticle.SetDefinitionAndUpdateE( aKaonZS );
379           else
380             currentParticle.SetDefinitionAndUpdateE( aKaonZL );
381           targetParticle.SetDefinitionAndUpdateE( aProton );
382           incidentHasChanged = true;
383           targetHasChanged = true;
384         }
385         break;
386       case 1:
387         if( G4UniformRand() < 0.5 )
388           currentParticle.SetDefinitionAndUpdateE( aKaonZS );
389         else
390           currentParticle.SetDefinitionAndUpdateE( aKaonZL );
391         incidentHasChanged = true;
392         break;
393       default: // assumes nm = np+1 so charge is conserved
394         targetParticle.SetDefinitionAndUpdateE( aProton );
395         targetHasChanged = true;
396         break;
397      }
398    }
399    SetUpPions( np, nm, nz, vec, vecLen );
400    return;
401  }
402
403 /* end of file */
404 
Note: See TracBrowser for help on using the repository browser.