source: trunk/source/processes/hadronic/models/low_energy/src/G4LEAntiOmegaMinusInelastic.cc @ 1340

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

update ti head

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