source: trunk/source/processes/hadronic/models/rpg/src/G4RPGAntiOmegaMinusInelastic.cc @ 1340

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

update ti head

File size: 13.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: G4RPGAntiOmegaMinusInelastic.cc,v 1.1 2007/07/18 21:04:20 dennis Exp $
27// GEANT4 tag $Name: geant4-09-03-ref-09 $
28//
29//
30// NOTE:  The FORTRAN version of the cascade, CASAOM, simply called the
31//        routine for the OmegaMinus particle.  Hence, the Cascade function
32//        below is just a copy of the Cascade from the OmegaMinus particle.
33 
34#include "G4RPGAntiOmegaMinusInelastic.hh"
35#include "Randomize.hh"
36
37G4HadFinalState*
38G4RPGAntiOmegaMinusInelastic::ApplyYourself( const G4HadProjectile &aTrack,
39                                             G4Nucleus &targetNucleus )
40{ 
41  const G4HadProjectile *originalIncident = &aTrack;
42  if (originalIncident->GetKineticEnergy()<= 0.1*MeV) 
43  {
44    theParticleChange.SetStatusChange(isAlive);
45    theParticleChange.SetEnergyChange(aTrack.GetKineticEnergy());
46    theParticleChange.SetMomentumChange(aTrack.Get4Momentum().vect().unit()); 
47    return &theParticleChange;     
48  }
49
50  // create the target particle
51
52  G4DynamicParticle *originalTarget = targetNucleus.ReturnTargetParticle();
53   
54    if( verboseLevel > 1 )
55    {
56      const G4Material *targetMaterial = aTrack.GetMaterial();
57      G4cout << "kinetic energy = " << originalIncident->GetKineticEnergy()/MeV << "MeV, ";
58      G4cout << "target material = " << targetMaterial->GetName() << ", ";
59      G4cout << "target particle = " << originalTarget->GetDefinition()->GetParticleName()
60           << G4endl;
61    }
62    //
63    // Fermi motion and evaporation
64    // As of Geant3, the Fermi energy calculation had not been Done
65    //
66    G4double ek = originalIncident->GetKineticEnergy()/MeV;
67    G4double amas = originalIncident->GetDefinition()->GetPDGMass()/MeV;
68    G4ReactionProduct modifiedOriginal;
69    modifiedOriginal = *originalIncident;
70   
71    G4double tkin = targetNucleus.Cinema( ek );
72    ek += tkin;
73    modifiedOriginal.SetKineticEnergy( ek*MeV );
74    G4double et = ek + amas;
75    G4double p = std::sqrt( std::abs((et-amas)*(et+amas)) );
76    G4double pp = modifiedOriginal.GetMomentum().mag()/MeV;
77    if( pp > 0.0 )
78    {
79      G4ThreeVector momentum = modifiedOriginal.GetMomentum();
80      modifiedOriginal.SetMomentum( momentum * (p/pp) );
81    }
82    //
83    // calculate black track energies
84    //
85    tkin = targetNucleus.EvaporationEffects( ek );
86    ek -= tkin;
87    modifiedOriginal.SetKineticEnergy( ek*MeV );
88    et = ek + amas;
89    p = std::sqrt( std::abs((et-amas)*(et+amas)) );
90    pp = modifiedOriginal.GetMomentum().mag()/MeV;
91    if( pp > 0.0 )
92    {
93      G4ThreeVector momentum = modifiedOriginal.GetMomentum();
94      modifiedOriginal.SetMomentum( momentum * (p/pp) );
95    }
96    G4ReactionProduct currentParticle = modifiedOriginal;
97    G4ReactionProduct targetParticle;
98    targetParticle = *originalTarget;
99    currentParticle.SetSide( 1 ); // incident always goes in forward hemisphere
100    targetParticle.SetSide( -1 );  // target always goes in backward hemisphere
101    G4bool incidentHasChanged = false;
102    G4bool targetHasChanged = false;
103    G4bool quasiElastic = false;
104    G4FastVector<G4ReactionProduct,GHADLISTSIZE> vec;  // vec will contain the secondary particles
105    G4int vecLen = 0;
106    vec.Initialize( 0 );
107       
108    const G4double cutOff = 0.1;
109    const G4double anni = std::min( 1.3*currentParticle.GetTotalMomentum()/GeV, 0.4 );
110   
111    if( (currentParticle.GetKineticEnergy()/MeV > cutOff) || (G4UniformRand() > anni) )
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
130void G4RPGAntiOmegaMinusInelastic::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 H. Fesefeldt's original FORTRAN code CASOM
141  // AntiOmegaMinus 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 pOriginal = originalIncident->GetTotalMomentum()/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  {  // not energetically possible to produce pion(s)
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-1); 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=np; 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 *aXiZero = G4XiZero::XiZero();
232    G4double n, anpn;
233    GetNormalizationConstant( availableEnergy, n, anpn );
234    G4double ran = G4UniformRand();
235    G4double dum, excs = 0.0;
236    G4int nvefix = 0;
237    if( targetParticle.GetDefinition() == aProton )
238    {
239      counter = -1;
240      for( np=0; np<numSec/3 && ran>=excs; ++np )
241      {
242        for( nm=std::max(0,np-1); nm<=(np+1) && ran>=excs; ++nm )
243        {
244          for( nz=0; nz<numSec/3 && ran>=excs; ++nz )
245          {
246            if( ++counter < numMul )
247            {
248              nt = np+nm+nz;
249              if( nt>0 && nt<=numSec )
250              {
251                test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
252                dum = (pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n);
253                if( std::fabs(dum) < 1.0 )
254                {
255                  if( test >= 1.0e-10 )excs += dum*test;
256                }
257                else
258                  excs += dum*test;
259              }
260            }
261          }
262        }
263      }
264      if( ran >= excs )  // 3 previous loops continued to the end
265      {
266        quasiElastic = true;
267        return;
268      }
269      np--; nm--; nz--;
270      //
271      // number of secondary mesons determined by kno distribution
272      // check for total charge of final state mesons to determine
273      // the kind of baryons to be produced, taking into account
274      // charge and strangeness conservation
275      //
276      if( np < nm )
277      {
278        if( np+1 == nm )
279        {
280          currentParticle.SetDefinitionAndUpdateE( aXiZero );
281          incidentHasChanged = true;
282          nvefix = 1;
283        }
284        else   // charge mismatch
285        {
286          currentParticle.SetDefinitionAndUpdateE( aSigmaPlus );
287          incidentHasChanged = true;
288          nvefix = 2;
289        }
290      }
291      else if( np > nm )
292      {
293        targetParticle.SetDefinitionAndUpdateE( aNeutron );
294        targetHasChanged = true;
295      }
296    }
297    else  // target must be a neutron
298    {
299      counter = -1;
300      for( np=0; np<numSec/3 && ran>=excs; ++np )
301      {
302        for( nm=np; nm<=(np+2) && ran>=excs; ++nm )
303        {
304          for( nz=0; nz<numSec/3 && ran>=excs; ++nz )
305          {
306            if( ++counter < numMul )
307            {
308              nt = np+nm+nz;
309              if( nt>0 && nt<=numSec )
310              {
311                test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
312                dum = (pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n);
313                if( std::fabs(dum) < 1.0 )
314                {
315                  if( test >= 1.0e-10 )excs += dum*test;
316                }
317                else
318                  excs += dum*test;
319              }
320            }
321          }
322        }
323      }
324      if( ran >= excs )  // 3 previous loops continued to the end
325      {
326        quasiElastic = true;
327        return;
328      }
329      np--; nm--; nz--;
330      if( np+1 < nm )
331      {
332        if( np+2 == nm )
333        {
334          currentParticle.SetDefinitionAndUpdateE( aXiZero );
335          incidentHasChanged = true;
336          nvefix = 1;
337        }
338        else   // charge mismatch
339        {
340          currentParticle.SetDefinitionAndUpdateE( aSigmaPlus );
341          incidentHasChanged = true;
342          nvefix = 2;
343        }
344        targetParticle.SetDefinitionAndUpdateE( aProton );
345        targetHasChanged = true;
346      }
347      else if( np+1 == nm )
348      {
349        targetParticle.SetDefinitionAndUpdateE( aProton );
350        targetHasChanged = true;
351      }
352    }
353    SetUpPions( np, nm, nz, vec, vecLen );
354    for( G4int i=0; i<vecLen && nvefix>0; ++i )
355    {
356      if( vec[i]->GetDefinition() == G4PionMinus::PionMinus() )
357      {
358        //
359        // correct the strangeness by replacing a pi- by a kaon-
360        //
361        if( nvefix >= 1 )vec[i]->SetDefinitionAndUpdateE( aKaonMinus );
362        --nvefix;
363      }
364    }
365  return;
366}
367
368 /* end of file */
369 
Note: See TracBrowser for help on using the repository browser.