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

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

update ti head

File size: 13.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// $Id: G4RPGOmegaMinusInelastic.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 "G4RPGOmegaMinusInelastic.hh"
31#include "Randomize.hh"
32
33G4HadFinalState*
34G4RPGOmegaMinusInelastic::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//    G4double targetMass = originalTarget->GetDefinition()->GetPDGMass();
50  G4ReactionProduct targetParticle( originalTarget->GetDefinition() );
51   
52  if( verboseLevel > 1 )
53  {
54    const G4Material *targetMaterial = aTrack.GetMaterial();
55    G4cout << "G4RPGOmegaMinusInelastic::ApplyYourself called" << G4endl;
56    G4cout << "kinetic energy = " << originalIncident->GetKineticEnergy() << "MeV, ";
57    G4cout << "target material = " << targetMaterial->GetName() << ", ";
58    G4cout << "target particle = " << originalTarget->GetDefinition()->GetParticleName()
59           << G4endl;
60  }
61    G4ReactionProduct currentParticle( const_cast<G4ParticleDefinition *>(originalIncident->GetDefinition() ));
62    currentParticle.SetMomentum( originalIncident->Get4Momentum().vect() );
63    currentParticle.SetKineticEnergy( originalIncident->GetKineticEnergy() );
64   
65    // Fermi motion and evaporation
66    // As of Geant3, the Fermi energy calculation had not been Done
67   
68    G4double ek = originalIncident->GetKineticEnergy();
69    G4double amas = originalIncident->GetDefinition()->GetPDGMass();
70   
71    G4double tkin = targetNucleus.Cinema( ek );
72    ek += tkin;
73    currentParticle.SetKineticEnergy( ek );
74    G4double et = ek + amas;
75    G4double p = std::sqrt( std::abs((et-amas)*(et+amas)) );
76    G4double pp = currentParticle.GetMomentum().mag();
77    if( pp > 0.0 )
78    {
79      G4ThreeVector momentum = currentParticle.GetMomentum();
80      currentParticle.SetMomentum( momentum * (p/pp) );
81    }
82   
83    // calculate black track energies
84   
85    tkin = targetNucleus.EvaporationEffects( ek );
86    ek -= tkin;
87    currentParticle.SetKineticEnergy( ek );
88    et = ek + amas;
89    p = std::sqrt( std::abs((et-amas)*(et+amas)) );
90    pp = currentParticle.GetMomentum().mag();
91    if( pp > 0.0 )
92    {
93      G4ThreeVector momentum = currentParticle.GetMomentum();
94      currentParticle.SetMomentum( momentum * (p/pp) );
95    }
96
97    G4ReactionProduct modifiedOriginal = currentParticle;
98
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*MeV;
109    if( currentParticle.GetKineticEnergy() > cutOff )
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 
128void G4RPGOmegaMinusInelastic::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 H. Fesefeldt's original FORTRAN code CASOM
139  // OmegaMinus 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  if( availableEnergy <= G4PionPlus::PionPlus()->GetPDGMass() )
156  {
157    quasiElastic = true;
158    return;
159  }
160    static G4bool first = true;
161    const G4int numMul = 1200;
162    const G4int numSec = 60;
163    static G4double protmul[numMul], protnorm[numSec]; // proton constants
164    static G4double neutmul[numMul], neutnorm[numSec]; // neutron constants
165    // np = number of pi+, nm = number of pi-, nz = number of pi0
166    G4int counter, nt=0, np=0, nm=0, nz=0;
167    G4double test;
168    const G4double c = 1.25;   
169    const G4double b[] = { 0.70, 0.70 };
170    if( first )    // compute normalization constants, this will only be Done once
171    {
172      first = false;
173      G4int i;
174      for( i=0; i<numMul; ++i )protmul[i] = 0.0;
175      for( i=0; i<numSec; ++i )protnorm[i] = 0.0;
176      counter = -1;
177      for( np=0; np<(numSec/3); ++np )
178      {
179        for( nm=std::max(0,np-1); nm<=(np+1); ++nm )
180        {
181          for( nz=0; nz<numSec/3; ++nz )
182          {
183            if( ++counter < numMul )
184            {
185              nt = np+nm+nz;
186              if( nt > 0 )
187              {
188                protmul[counter] = Pmltpc(np,nm,nz,nt,b[0],c);
189                protnorm[nt-1] += protmul[counter];
190              }
191            }
192          }
193        }
194      }
195      for( i=0; i<numMul; ++i )neutmul[i] = 0.0;
196      for( i=0; i<numSec; ++i )neutnorm[i] = 0.0;
197      counter = -1;
198      for( np=0; np<numSec/3; ++np )
199      {
200        for( nm=np; nm<=(np+2); ++nm )
201        {
202          for( nz=0; nz<numSec/3; ++nz )
203          {
204            if( ++counter < numMul )
205            {
206              nt = np+nm+nz;
207              if( (nt>0) && (nt<=numSec) )
208              {
209                neutmul[counter] = Pmltpc(np,nm,nz,nt,b[1],c);
210                neutnorm[nt-1] += neutmul[counter];
211              }
212            }
213          }
214        }
215      }
216      for( i=0; i<numSec; ++i )
217      {
218        if( protnorm[i] > 0.0 )protnorm[i] = 1.0/protnorm[i];
219        if( neutnorm[i] > 0.0 )neutnorm[i] = 1.0/neutnorm[i];
220      }
221    }   // end of initialization
222   
223    const G4double expxu = 82.;           // upper bound for arg. of exp
224    const G4double expxl = -expxu;        // lower bound for arg. of exp
225    G4ParticleDefinition *aNeutron = G4Neutron::Neutron();
226    G4ParticleDefinition *aProton = G4Proton::Proton();
227    G4ParticleDefinition *aKaonMinus = G4KaonMinus::KaonMinus();
228    G4ParticleDefinition *aPiMinus = G4PionMinus::PionMinus();
229    G4ParticleDefinition *aSigmaPlus = G4SigmaPlus::SigmaPlus();
230    G4ParticleDefinition *aXiZero = G4XiZero::XiZero();
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-1); 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 )
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    else  // target must be a neutron
273    {
274      counter = -1;
275      for( np=0; np<numSec/3 && ran>=excs; ++np )
276      {
277        for( nm=np; nm<=(np+2) && ran>=excs; ++nm )
278        {
279          for( nz=0; nz<numSec/3 && ran>=excs; ++nz )
280          {
281            if( ++counter < numMul )
282            {
283              nt = np+nm+nz;
284              if( (nt>=1) && (nt<=numSec) )
285              {
286                test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
287                dum = (pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n);
288                if( std::fabs(dum) < 1.0 )
289                {
290                  if( test >= 1.0e-10 )excs += dum*test;
291                }
292                else
293                  excs += dum*test;
294              }
295            }
296          }
297        }
298      }
299      if( ran >= excs )  // 3 previous loops continued to the end
300      {
301        quasiElastic = true;
302        return;
303      }
304      np--; nm--; nz--;
305    }
306    // number of secondary mesons determined by kno distribution
307    // check for total charge of final state mesons to determine
308    // the kind of baryons to be produced, taking into account
309    // charge and strangeness conservation
310    //
311    G4int nvefix = 0;
312    if( targetParticle.GetDefinition() == aProton )
313    {
314      if( nm > np )
315      {
316        if( nm == np+1 )
317        {
318          currentParticle.SetDefinitionAndUpdateE( aXiZero );
319          nvefix = 1;
320        }
321        else
322        {
323          currentParticle.SetDefinitionAndUpdateE( aSigmaPlus );
324          nvefix = 2;
325        }
326        incidentHasChanged = true;
327      }
328      else if( nm < np )
329      {
330        targetParticle.SetDefinitionAndUpdateE( aNeutron );
331        targetHasChanged = true;
332      }
333    }
334    else // target is a neutron
335    {
336      if( np+1 < nm )
337      {
338        if( nm == np+2 )
339        {
340          currentParticle.SetDefinitionAndUpdateE( aXiZero );
341          incidentHasChanged = true;
342          nvefix = 1;
343        }
344        else   // charge mismatch
345        {
346          currentParticle.SetDefinitionAndUpdateE( aSigmaPlus );
347          incidentHasChanged = true;
348          nvefix = 2;
349        }
350        targetParticle.SetDefinitionAndUpdateE( aProton );
351        targetHasChanged = true;
352      }
353      else if( nm == np+1 )
354      {
355        targetParticle.SetDefinitionAndUpdateE( aProton );
356        targetHasChanged = true;
357      }
358    }
359    SetUpPions( np, nm, nz, vec, vecLen );
360    for( G4int i=0; i<vecLen && nvefix>0; ++i )
361    {
362      if( vec[i]->GetDefinition() == aPiMinus )
363      {
364        if( nvefix >= 1 )vec[i]->SetDefinitionAndUpdateE( aKaonMinus );
365        --nvefix;
366      }
367    }
368    return;
369}
370
371 /* end of file */
372 
Note: See TracBrowser for help on using the repository browser.