source: trunk/source/processes/hadronic/models/low_energy/src/G4LEXiZeroInelastic.cc @ 1358

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

update ti head

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