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

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

update ti head

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