source: trunk/source/processes/hadronic/models/rpg/src/G4RPGAntiKZeroInelastic.cc @ 1337

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

tag geant4.9.4 beta 1 + modifs locales

File size: 18.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: G4RPGAntiKZeroInelastic.cc,v 1.1 2007/07/18 21:04:20 dennis Exp $
27// GEANT4 tag $Name: geant4-09-04-beta-01 $
28//
29 
30#include "G4RPGAntiKZeroInelastic.hh"
31#include "Randomize.hh"
32#include "G4HadReentrentException.hh"
33
34 G4HadFinalState *
35  G4RPGAntiKZeroInelastic::ApplyYourself( const G4HadProjectile &aTrack,
36                                         G4Nucleus &targetNucleus )
37  {   
38    const G4HadProjectile *originalIncident = &aTrack;
39    //
40    // create the target particle
41    //
42    G4DynamicParticle *originalTarget = targetNucleus.ReturnTargetParticle();
43   
44    if( verboseLevel > 1 )
45    {
46      const G4Material *targetMaterial = aTrack.GetMaterial();
47      G4cout << "G4RPGAntiKZeroInelastic::ApplyYourself called" << G4endl;   
48      G4cout << "kinetic energy = " << originalIncident->GetKineticEnergy()/MeV << "MeV, ";
49      G4cout << "target material = " << targetMaterial->GetName() << ", ";
50      G4cout << "target particle = " << originalTarget->GetDefinition()->GetParticleName()
51           << G4endl;
52    }
53    //
54    // Fermi motion and evaporation
55    // As of Geant3, the Fermi energy calculation had not been Done
56    //
57    G4double ek = originalIncident->GetKineticEnergy()/MeV;
58    G4double amas = originalIncident->GetDefinition()->GetPDGMass()/MeV;
59    G4ReactionProduct modifiedOriginal;
60    modifiedOriginal = *originalIncident;
61   
62    G4double tkin = targetNucleus.Cinema( ek );
63    ek += tkin;
64    modifiedOriginal.SetKineticEnergy( ek*MeV );
65    G4double et = ek + amas;
66    G4double p = std::sqrt( std::abs((et-amas)*(et+amas)) );
67    G4double pp = modifiedOriginal.GetMomentum().mag()/MeV;
68    if( pp > 0.0 )
69    {
70      G4ThreeVector momentum = modifiedOriginal.GetMomentum();
71      modifiedOriginal.SetMomentum( momentum * (p/pp) );
72    }
73    //
74    // calculate black track energies
75    //
76    tkin = targetNucleus.EvaporationEffects( ek );
77    ek -= tkin;
78    modifiedOriginal.SetKineticEnergy( ek*MeV );
79    et = ek + amas;
80    p = std::sqrt( std::abs((et-amas)*(et+amas)) );
81    pp = modifiedOriginal.GetMomentum().mag()/MeV;
82    if( pp > 0.0 )
83    {
84      G4ThreeVector momentum = modifiedOriginal.GetMomentum();
85      modifiedOriginal.SetMomentum( momentum * (p/pp) );
86    }
87    G4ReactionProduct currentParticle = modifiedOriginal;
88    G4ReactionProduct targetParticle;
89    targetParticle = *originalTarget;
90    currentParticle.SetSide( 1 ); // incident always goes in forward hemisphere
91    targetParticle.SetSide( -1 );  // target always goes in backward hemisphere
92    G4bool incidentHasChanged = false;
93    G4bool targetHasChanged = false;
94    G4bool quasiElastic = false;
95    G4FastVector<G4ReactionProduct,GHADLISTSIZE> vec;  // vec will contain the secondary particles
96    G4int vecLen = 0;
97    vec.Initialize( 0 );
98   
99    const G4double cutOff = 0.1;
100    if( currentParticle.GetKineticEnergy()/MeV > cutOff )
101      Cascade( vec, vecLen,
102               originalIncident, currentParticle, targetParticle,
103               incidentHasChanged, targetHasChanged, quasiElastic );
104   
105    try
106    {
107        CalculateMomenta( vec, vecLen,
108                      originalIncident, originalTarget, modifiedOriginal,
109                      targetNucleus, currentParticle, targetParticle,
110                      incidentHasChanged, targetHasChanged, quasiElastic );
111    }
112    catch(G4HadReentrentException aR)
113    {
114      aR.Report(G4cout);
115      throw G4HadReentrentException(__FILE__, __LINE__, "Bailing out");
116    }
117    SetUpChange( vec, vecLen,
118                 currentParticle, targetParticle,
119                 incidentHasChanged );
120   
121    delete originalTarget;
122    return &theParticleChange;
123  }
124 
125void G4RPGAntiKZeroInelastic::Cascade(
126   G4FastVector<G4ReactionProduct,GHADLISTSIZE> &vec,
127   G4int& vecLen,
128   const G4HadProjectile *originalIncident,
129   G4ReactionProduct &currentParticle,
130   G4ReactionProduct &targetParticle,
131   G4bool &incidentHasChanged,
132   G4bool &targetHasChanged,
133   G4bool &quasiElastic )
134{
135  // Derived from H. Fesefeldt's original FORTRAN code CASK0B
136  //
137  // K0Long undergoes interaction with nucleon within a nucleus.  Check if it is
138  // energetically possible to produce pions/kaons.  In not, assume nuclear excitation
139  // occurs and input particle is degraded in energy. No other particles are produced.
140  // If reaction is possible, find the correct number of pions/protons/neutrons
141  // produced using an interpolation to multiplicity data.  Replace some pions or
142  // protons/neutrons by kaons or strange baryons according to the average
143  // multiplicity per Inelastic reaction.
144
145  const G4double mOriginal = originalIncident->GetDefinition()->GetPDGMass()/MeV;
146  const G4double etOriginal = originalIncident->Get4Momentum().e()/MeV;
147  const G4double pOriginal = originalIncident->Get4Momentum().vect().mag()/MeV;
148  const G4double targetMass = targetParticle.GetMass()/MeV;
149  G4double centerofmassEnergy = std::sqrt( mOriginal*mOriginal +
150                                      targetMass*targetMass +
151                                      2.0*targetMass*etOriginal );
152  G4double availableEnergy = centerofmassEnergy-(targetMass+mOriginal);
153
154  static G4bool first = true;
155  const G4int numMul = 1200;
156  const G4int numSec = 60;
157  static G4double protmul[numMul], protnorm[numSec]; // proton constants
158  static G4double neutmul[numMul], neutnorm[numSec]; // neutron constants
159
160  // np = number of pi+, nm = number of pi-, nz = number of pi0
161
162  G4int counter, nt=0, np=0, nm=0, nz=0;
163  const G4double c = 1.25;   
164  const G4double b[] = { 0.7, 0.7 };
165  if( first )       // compute normalization constants, this will only be Done once
166  {
167      first = false;
168      G4int i;
169      for( i=0; i<numMul; ++i )protmul[i] = 0.0;
170      for( i=0; i<numSec; ++i )protnorm[i] = 0.0;
171      counter = -1;
172      for( np=0; np<numSec/3; ++np )
173      {
174        for( nm=std::max(0,np-2); nm<=np; ++nm )
175        {
176          for( nz=0; nz<numSec/3; ++nz )
177          {
178            if( ++counter < numMul )
179            {
180              nt = np+nm+nz;
181              if( nt>0 && nt<=numSec )
182              {
183                protmul[counter] = Pmltpc(np,nm,nz,nt,b[0],c);
184                protnorm[nt-1] += protmul[counter];
185              }
186            }
187          }
188        }
189      }
190      for( i=0; i<numMul; ++i )neutmul[i] = 0.0;
191      for( i=0; i<numSec; ++i )neutnorm[i] = 0.0;
192      counter = -1;
193      for( np=0; np<(numSec/3); ++np )
194      {
195        for( nm=std::max(0,np-1); nm<=(np+1); ++nm )
196        {
197          for( nz=0; nz<numSec/3; ++nz )
198          {
199            if( ++counter < numMul )
200            {
201              nt = np+nm+nz;
202              if( nt>0 && nt<=numSec )
203              {
204                neutmul[counter] = Pmltpc(np,nm,nz,nt,b[1],c);
205                neutnorm[nt-1] += neutmul[counter];
206              }
207            }
208          }
209        }
210      }
211    for( i=0; i<numSec; ++i )
212    {
213      if( protnorm[i] > 0.0 )protnorm[i] = 1.0/protnorm[i];
214      if( neutnorm[i] > 0.0 )neutnorm[i] = 1.0/neutnorm[i];
215    }
216  }   // end of initialization
217   
218  const G4double expxu = 82.;           // upper bound for arg. of exp
219  const G4double expxl = -expxu;        // lower bound for arg. of exp
220  G4ParticleDefinition *aKaonMinus = G4KaonMinus::KaonMinus();
221  G4ParticleDefinition *aKaonZS = G4KaonZeroShort::KaonZeroShort();
222  G4ParticleDefinition *aKaonZL = G4KaonZeroLong::KaonZeroLong();
223  G4ParticleDefinition *aNeutron = G4Neutron::Neutron();
224  G4ParticleDefinition *aProton = G4Proton::Proton();
225  G4ParticleDefinition *aPiPlus = G4PionPlus::PionPlus();
226  G4ParticleDefinition *aPiMinus = G4PionMinus::PionMinus();
227  G4ParticleDefinition *aPiZero = G4PionZero::PionZero();
228  G4ParticleDefinition *aLambda = G4Lambda::Lambda();
229  G4ParticleDefinition *aSigmaPlus = G4SigmaPlus::SigmaPlus();
230  G4ParticleDefinition *aSigmaMinus = G4SigmaMinus::SigmaMinus();
231  G4ParticleDefinition *aSigmaZero = G4SigmaZero::SigmaZero();
232  const G4double cech[] = {1.,1.,1.,0.70,0.60,0.55,0.35,0.25,0.18,0.15};
233  G4int iplab = G4int(std::min( 9.0, 5.0*pOriginal*MeV/GeV ));
234
235    if( (pOriginal*MeV/GeV <= 2.0) && (G4UniformRand() < cech[iplab]) )
236    {
237      np = nm = nz = nt = 0;
238      iplab = G4int(std::min( 19.0, pOriginal*MeV/GeV*10.0 ));
239      const G4double cnk0[] = {0.17,0.18,0.17,0.24,0.26,0.20,0.22,0.21,0.34,0.45,
240                               0.58,0.55,0.36,0.29,0.29,0.32,0.32,0.33,0.33,0.33};
241      if( G4UniformRand() > cnk0[iplab] )
242      {
243        G4double ran = G4UniformRand();
244        if( ran < 0.25 )         // k0Long n --> pi- s+
245        {
246          if( targetParticle.GetDefinition() == aNeutron )
247          {
248            currentParticle.SetDefinitionAndUpdateE( aPiMinus );
249            targetParticle.SetDefinitionAndUpdateE( aSigmaPlus );
250            incidentHasChanged = true;
251            targetHasChanged = true;
252          }
253        }
254        else if( ran < 0.50 )  // k0Long p --> pi+ s0  or  k0Long n --> pi0 s0
255        {
256          if( targetParticle.GetDefinition() == aNeutron )
257            currentParticle.SetDefinitionAndUpdateE( aPiZero );
258          else
259            currentParticle.SetDefinitionAndUpdateE( aPiPlus );
260          targetParticle.SetDefinitionAndUpdateE( aSigmaZero );
261          incidentHasChanged = true;
262          targetHasChanged = true;
263        }
264        else if( ran < 0.75 )  // k0Long n --> pi+ s-
265        {
266          if( targetParticle.GetDefinition() == aNeutron )
267          {
268            currentParticle.SetDefinitionAndUpdateE( aPiPlus );
269            targetParticle.SetDefinitionAndUpdateE( aSigmaMinus );
270            incidentHasChanged = true;
271            targetHasChanged = true;
272          }
273        }
274        else                   // k0Long p --> pi+ L  or  k0Long n --> pi0 L
275        {
276          if( targetParticle.GetDefinition() == aNeutron )
277            currentParticle.SetDefinitionAndUpdateE( aPiZero );
278          else
279            currentParticle.SetDefinitionAndUpdateE( aPiPlus );
280          targetParticle.SetDefinitionAndUpdateE( aLambda );
281          incidentHasChanged = true;
282          targetHasChanged = true;
283        }
284      }
285      else   // ran <= cnk0
286      {
287        quasiElastic = true;
288        if( targetParticle.GetDefinition() == aNeutron )
289        {
290          currentParticle.SetDefinitionAndUpdateE( aKaonMinus );
291          targetParticle.SetDefinitionAndUpdateE( aProton );
292          incidentHasChanged = true;
293          targetHasChanged = true;
294        }
295      }
296    }
297    else  // (pOriginal > 2.0*GeV) || (random number >= cech[iplab])
298    {
299      if( availableEnergy < aPiPlus->GetPDGMass()/MeV )
300      {
301        quasiElastic = true;
302        return;
303      }
304      G4double n, anpn;
305      GetNormalizationConstant( availableEnergy, n, anpn );
306      G4double ran = G4UniformRand();
307      G4double dum, test, excs = 0.0;
308      if( targetParticle.GetDefinition() == aProton )
309      {
310        counter = -1;
311        for( np=0; (np<numSec/3) && (ran>=excs); ++np )
312        {
313          for( nm=std::max(0,np-2); nm<=np && ran>=excs; ++nm )
314          {
315            for( nz=0; nz<numSec/3 && ran>=excs; ++nz )
316            {
317              if( ++counter < numMul )
318              {
319                nt = np+nm+nz;
320                if( nt>0 && nt<=numSec )
321                {
322                  test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
323                  dum = (pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n);
324                  if( std::fabs(dum) < 1.0 )
325                  {
326                    if( test >= 1.0e-10 )excs += dum*test;
327                  }
328                  else
329                    excs += dum*test;
330                }
331              }
332            }
333          }
334        }
335        if( ran >= excs )  // 3 previous loops continued to the end
336        {
337          quasiElastic = true;
338          return;
339        }
340        np--; nm--; nz--;
341        switch( np-nm )
342        {
343         case 1:
344           if( G4UniformRand() < 0.5 )
345           {
346             currentParticle.SetDefinitionAndUpdateE( aKaonMinus );
347             incidentHasChanged = true;
348           }
349           else
350           {
351             targetParticle.SetDefinitionAndUpdateE( aNeutron );
352             targetHasChanged = true;
353           }
354         case 0:
355           break;
356         default:
357           currentParticle.SetDefinitionAndUpdateE( aKaonMinus );
358           targetParticle.SetDefinitionAndUpdateE( aNeutron );
359           incidentHasChanged = true;
360           targetHasChanged = true;
361           break;
362        }
363      }
364      else  // target must be a neutron
365      {
366        counter = -1;
367        for( np=0; np<numSec/3 && ran>=excs; ++np )
368        {
369          for( nm=std::max(0,np-1); nm<=(np+1) && ran>=excs; ++nm )
370          {
371            for( nz=0; nz<numSec/3 && ran>=excs; ++nz )
372            {
373              if( ++counter < numMul )
374              {
375                nt = np+nm+nz;
376                if( nt>0 && nt<=numSec )
377                {
378                  test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
379                  dum = (pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n);
380                  if( std::fabs(dum) < 1.0 )
381                  {
382                    if( test >= 1.0e-10 )excs += dum*test;
383                  }
384                  else
385                    excs += dum*test;
386                }
387              }
388            }
389          }
390        }
391        if( ran >= excs )  // 3 previous loops continued to the end
392        {
393          quasiElastic = true;
394          return;
395        }
396        np--; nm--; nz--;
397        switch( np-nm )
398        {
399         case 0:
400           currentParticle.SetDefinitionAndUpdateE( aKaonMinus );
401           targetParticle.SetDefinitionAndUpdateE( aProton );
402           incidentHasChanged = true;
403           targetHasChanged = true;
404           break;
405         case 1:
406           currentParticle.SetDefinitionAndUpdateE( aKaonMinus );
407           incidentHasChanged = true;
408           break;
409         default:
410           targetParticle.SetDefinitionAndUpdateE( aProton );
411           targetHasChanged = true;
412           break;
413        }
414      }
415      if( G4UniformRand() >= 0.5 )
416      {
417        if( currentParticle.GetDefinition() == aKaonMinus &&
418            targetParticle.GetDefinition() == aNeutron )
419        {
420          ran = G4UniformRand();
421          if( ran < 0.68 )
422          {
423            currentParticle.SetDefinitionAndUpdateE( aPiMinus );
424            targetParticle.SetDefinitionAndUpdateE( aLambda );
425          }
426          else if( ran < 0.84 )
427          {
428            currentParticle.SetDefinitionAndUpdateE( aPiMinus );
429            targetParticle.SetDefinitionAndUpdateE( aSigmaZero );
430          }
431          else
432          {
433            currentParticle.SetDefinitionAndUpdateE( aPiZero );
434            targetParticle.SetDefinitionAndUpdateE( aSigmaMinus );
435          }
436        }
437        else if( (currentParticle.GetDefinition() == aKaonZS ||
438                  currentParticle.GetDefinition() == aKaonZL ) &&
439                 targetParticle.GetDefinition() == aProton )
440        {
441          ran = G4UniformRand();
442          if( ran < 0.68 )
443          {
444            currentParticle.SetDefinitionAndUpdateE( aPiPlus );
445            targetParticle.SetDefinitionAndUpdateE( aLambda );
446          }
447          else if( ran < 0.84 )
448          {
449            currentParticle.SetDefinitionAndUpdateE( aPiZero );
450            targetParticle.SetDefinitionAndUpdateE( aSigmaPlus );
451          }
452          else
453          {
454            currentParticle.SetDefinitionAndUpdateE( aPiPlus );
455            targetParticle.SetDefinitionAndUpdateE( aSigmaZero );
456          }
457        }
458        else
459        {
460          ran = G4UniformRand();
461          if( ran < 0.67 )
462          {
463            currentParticle.SetDefinitionAndUpdateE( aPiZero );
464            targetParticle.SetDefinitionAndUpdateE( aLambda );
465          }
466          else if( ran < 0.78 )
467          {
468            currentParticle.SetDefinitionAndUpdateE( aPiMinus );
469            targetParticle.SetDefinitionAndUpdateE( aSigmaPlus );
470          }
471          else if( ran < 0.89 )
472          {
473            currentParticle.SetDefinitionAndUpdateE( aPiZero );
474            targetParticle.SetDefinitionAndUpdateE( aSigmaZero );
475          }
476          else
477          {
478            currentParticle.SetDefinitionAndUpdateE( aPiPlus );
479            targetParticle.SetDefinitionAndUpdateE( aSigmaMinus );
480          }
481        }
482        incidentHasChanged = true;
483        targetHasChanged = true;
484      }
485    }
486    if( currentParticle.GetDefinition() == aKaonZL )
487    {
488      if( G4UniformRand() >= 0.5 )
489      {
490        currentParticle.SetDefinitionAndUpdateE( aKaonZS );
491        incidentHasChanged = true;
492      }
493    }
494    if( targetParticle.GetDefinition() == aKaonZL )
495    {
496      if( G4UniformRand() >= 0.5 )
497      {
498        targetParticle.SetDefinitionAndUpdateE( aKaonZS );
499        targetHasChanged = true;
500      }
501    }
502    SetUpPions( np, nm, nz, vec, vecLen );
503    return;
504}
505
506 /* end of file */
507 
Note: See TracBrowser for help on using the repository browser.