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

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

update ti head

File size: 19.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// $Id: G4RPGAntiSigmaPlusInelastic.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 "G4RPGAntiSigmaPlusInelastic.hh"
31#include "Randomize.hh"
32
33G4HadFinalState*
34G4RPGAntiSigmaPlusInelastic::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  // Choose the target particle
47
48  G4DynamicParticle *originalTarget = targetNucleus.ReturnTargetParticle();
49   
50  if( verboseLevel > 1 )
51  {
52    const G4Material *targetMaterial = aTrack.GetMaterial();
53    G4cout << "G4RPGAntiSigmaPlusInelastic::ApplyYourself called" << G4endl;
54    G4cout << "kinetic energy = " << originalIncident->GetKineticEnergy()/MeV << "MeV, ";
55    G4cout << "target material = " << targetMaterial->GetName() << ", ";
56    G4cout << "target particle = " << originalTarget->GetDefinition()->GetParticleName()
57           << G4endl;
58  }
59
60  // Fermi motion and evaporation
61  // As of Geant3, the Fermi energy calculation had not been Done
62
63    G4double ek = originalIncident->GetKineticEnergy()/MeV;
64    G4double amas = originalIncident->GetDefinition()->GetPDGMass()/MeV;
65    G4ReactionProduct modifiedOriginal;
66    modifiedOriginal = *originalIncident;
67   
68    G4double tkin = targetNucleus.Cinema( ek );
69    ek += tkin;
70    modifiedOriginal.SetKineticEnergy( ek*MeV );
71    G4double et = ek + amas;
72    G4double p = std::sqrt( std::abs((et-amas)*(et+amas)) );
73    G4double pp = modifiedOriginal.GetMomentum().mag()/MeV;
74    if( pp > 0.0 )
75    {
76      G4ThreeVector momentum = modifiedOriginal.GetMomentum();
77      modifiedOriginal.SetMomentum( momentum * (p/pp) );
78    }
79    //
80    // calculate black track energies
81    //
82    tkin = targetNucleus.EvaporationEffects( ek );
83    ek -= tkin;
84    modifiedOriginal.SetKineticEnergy( ek*MeV );
85    et = ek + amas;
86    p = std::sqrt( std::abs((et-amas)*(et+amas)) );
87    pp = modifiedOriginal.GetMomentum().mag()/MeV;
88    if( pp > 0.0 )
89    {
90      G4ThreeVector momentum = modifiedOriginal.GetMomentum();
91      modifiedOriginal.SetMomentum( momentum * (p/pp) );
92    }
93    G4ReactionProduct currentParticle = modifiedOriginal;
94    G4ReactionProduct targetParticle;
95    targetParticle = *originalTarget;
96    currentParticle.SetSide( 1 ); // incident always goes in forward hemisphere
97    targetParticle.SetSide( -1 );  // target always goes in backward hemisphere
98    G4bool incidentHasChanged = false;
99    G4bool targetHasChanged = false;
100    G4bool quasiElastic = false;
101    G4FastVector<G4ReactionProduct,GHADLISTSIZE> vec;  // vec will contain the secondary particles
102    G4int vecLen = 0;
103    vec.Initialize( 0 );
104   
105    const G4double cutOff = 0.1;
106    const G4double anni = std::min( 1.3*currentParticle.GetTotalMomentum()/GeV, 0.4 );
107    if( (currentParticle.GetKineticEnergy()/MeV > cutOff) || (G4UniformRand() > anni) )
108      Cascade( vec, vecLen,
109               originalIncident, currentParticle, targetParticle,
110               incidentHasChanged, targetHasChanged, quasiElastic );
111   
112    CalculateMomenta( vec, vecLen,
113                      originalIncident, originalTarget, modifiedOriginal,
114                      targetNucleus, currentParticle, targetParticle,
115                      incidentHasChanged, targetHasChanged, quasiElastic );
116   
117    SetUpChange( vec, vecLen,
118                 currentParticle, targetParticle,
119                 incidentHasChanged );
120   
121  delete originalTarget;
122  return &theParticleChange;
123}
124
125 
126void G4RPGAntiSigmaPlusInelastic::Cascade(
127   G4FastVector<G4ReactionProduct,GHADLISTSIZE> &vec,
128   G4int& vecLen,
129   const G4HadProjectile *originalIncident,
130   G4ReactionProduct &currentParticle,
131   G4ReactionProduct &targetParticle,
132   G4bool &incidentHasChanged,
133   G4bool &targetHasChanged,
134   G4bool &quasiElastic )
135{
136  // Derived from H. Fesefeldt's original FORTRAN code CASASP
137  // AntiSigmaPlus 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->GetTotalEnergy()/MeV;
147  const G4double pOriginal = originalIncident->GetTotalMomentum()/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 numMulA = 400;
157    const G4int numSec = 60;
158    static G4double protmul[numMul], protnorm[numSec]; // proton constants
159    static G4double neutmul[numMul], neutnorm[numSec]; // neutron constants
160    static G4double protmulA[numMulA], protnormA[numSec]; // proton constants
161    static G4double neutmulA[numMulA], neutnormA[numSec]; // neutron constants
162    // np = number of pi+, nm = number of pi-, nz = number of pi0
163    G4int counter, nt=0, np=0, nm=0, nz=0;
164    G4double test;
165    const G4double c = 1.25;   
166    const G4double b[] = { 0.7, 0.7 };
167    if( first )       // compute normalization constants, this will only be Done once
168    {
169      first = false;
170      G4int i;
171      for( i=0; i<numMul; ++i )protmul[i] = 0.0;
172      for( i=0; i<numSec; ++i )protnorm[i] = 0.0;
173      counter = -1;
174      for( np=0; np<(numSec/3); ++np )
175      {
176        for( nm=std::max(0,np-1); nm<=(np+1); ++nm )
177        {
178          for( nz=0; nz<numSec/3; ++nz )
179          {
180            if( ++counter < numMul )
181            {
182              nt = np+nm+nz;
183              if( nt>0 && nt<=numSec )
184              {
185                protmul[counter] = Pmltpc(np,nm,nz,nt,b[0],c);
186                protnorm[nt-1] += protmul[counter];
187              }
188            }
189          }
190        }
191      }
192      for( i=0; i<numMul; ++i )neutmul[i] = 0.0;
193      for( i=0; i<numSec; ++i )neutnorm[i] = 0.0;
194      counter = -1;
195      for( np=0; np<numSec/3; ++np )
196      {
197        for( nm=np; nm<=(np+2); ++nm )
198        {
199          for( nz=0; nz<numSec/3; ++nz )
200          {
201            if( ++counter < numMul )
202            {
203              nt = np+nm+nz;
204              if( nt>0 && nt<=numSec )
205              {
206                neutmul[counter] = Pmltpc(np,nm,nz,nt,b[1],c);
207                neutnorm[nt-1] += neutmul[counter];
208              }
209            }
210          }
211        }
212      }
213      for( i=0; i<numSec; ++i )
214      {
215        if( protnorm[i] > 0.0 )protnorm[i] = 1.0/protnorm[i];
216        if( neutnorm[i] > 0.0 )neutnorm[i] = 1.0/neutnorm[i];
217      }
218      //
219      // do the same for annihilation channels
220      //
221      for( i=0; i<numMulA; ++i )protmulA[i] = 0.0;
222      for( i=0; i<numSec; ++i )protnormA[i] = 0.0;
223      counter = -1;
224      for( np=1; np<(numSec/3); ++np )
225      {
226        nm = np;
227        for( nz=0; nz<numSec/3; ++nz )
228        {
229          if( ++counter < numMulA )
230          {
231            nt = np+nm+nz;
232            if( nt>1 && nt<=numSec )
233            {
234              protmulA[counter] = Pmltpc(np,nm,nz,nt,b[0],c);
235              protnormA[nt-1] += protmulA[counter];
236            }
237          }
238        }
239      }
240      for( i=0; i<numMulA; ++i )neutmulA[i] = 0.0;
241      for( i=0; i<numSec; ++i )neutnormA[i] = 0.0;
242      counter = -1;
243       for( np=0; np<numSec/3; ++np )
244       {
245        nm = np+1;
246        for( nz=0; nz<numSec/3; ++nz )
247        {
248          if( ++counter < numMulA )
249          {
250            nt = np+nm+nz;
251            if( nt>1 && nt<=numSec )
252            {
253              neutmulA[counter] = Pmltpc(np,nm,nz,nt,b[1],c);
254              neutnormA[nt-1] += neutmulA[counter];
255            }
256          }
257        }
258      }
259      for( i=0; i<numSec; ++i )
260      {
261        if( protnormA[i] > 0.0 )protnormA[i] = 1.0/protnormA[i];
262        if( neutnormA[i] > 0.0 )neutnormA[i] = 1.0/neutnormA[i];
263      }
264    }   // end of initialization
265   
266    const G4double expxu = 82.;           // upper bound for arg. of exp
267    const G4double expxl = -expxu;        // lower bound for arg. of exp
268    G4ParticleDefinition *aNeutron = G4Neutron::Neutron();
269    G4ParticleDefinition *aProton = G4Proton::Proton();
270    G4ParticleDefinition *aPiPlus = G4PionPlus::PionPlus();
271    G4ParticleDefinition *anAntiLambda = G4AntiLambda::AntiLambda();
272    G4ParticleDefinition *aKaonMinus = G4KaonMinus::KaonMinus();
273    G4ParticleDefinition *aKaonPlus = G4KaonPlus::KaonPlus();
274    G4ParticleDefinition *aKaonZL = G4KaonZeroLong::KaonZeroLong();
275    G4ParticleDefinition *anAntiSigmaZero = G4AntiSigmaZero::AntiSigmaZero();
276    const G4double anhl[] = {1.00,1.00,1.00,1.00,1.00,1.00,1.00,1.00,0.97,0.88,
277                             0.85,0.81,0.75,0.64,0.64,0.55,0.55,0.45,0.47,0.40,
278                             0.39,0.36,0.33,0.10,0.01};
279    G4int iplab = G4int( pOriginal/GeV*10.0 );
280    if( iplab >  9 )iplab = G4int( (pOriginal/GeV- 1.0)*5.0  ) + 10;
281    if( iplab > 14 )iplab = G4int(  pOriginal/GeV- 2.0       ) + 15;
282    if( iplab > 22 )iplab = G4int( (pOriginal/GeV-10.0)/10.0 ) + 23;
283    if( iplab > 24 )iplab = 24;
284    if( G4UniformRand() > anhl[iplab] )
285    {
286      if( availableEnergy <= aPiPlus->GetPDGMass()/MeV )
287      {
288        quasiElastic = true;
289        return;
290      }
291      G4double n, anpn;
292      GetNormalizationConstant( availableEnergy, n, anpn );
293      G4double ran = G4UniformRand();
294      G4double dum, excs = 0.0;
295      if( targetParticle.GetDefinition() == aProton )
296      {
297        counter = -1;
298        for( np=0; np<numSec/3 && ran>=excs; ++np )
299        {
300          for( nm=std::max(0,np-1); nm<=(np+1) && ran>=excs; ++nm )
301          {
302            for( nz=0; nz<numSec/3 && ran>=excs; ++nz )
303            {
304              if( ++counter < numMul )
305              {
306                nt = np+nm+nz;
307                if( (nt>0) && (nt<=numSec) )
308                {
309                  test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
310                  dum = (pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n);
311                  if( std::fabs(dum) < 1.0 )
312                  {
313                    if( test >= 1.0e-10 )excs += dum*test;
314                  }
315                  else
316                    excs += dum*test;
317                }
318              }
319            }
320          }
321        }
322        if( ran >= excs )  // 3 previous loops continued to the end
323        {
324          quasiElastic = true;
325          return;
326        }
327        np--; nm--; nz--;
328        G4int ncht = std::min( 3, std::max( 1, np-nm+2 ) );
329        switch( ncht )
330        {
331         case 1:
332           if( G4UniformRand() < 0.5 )
333             currentParticle.SetDefinitionAndUpdateE( anAntiLambda );
334           else
335             currentParticle.SetDefinitionAndUpdateE( anAntiSigmaZero );
336           incidentHasChanged = true;
337           break;
338         case 2:
339           if( G4UniformRand() >= 0.5 )
340           {
341             if( G4UniformRand() < 0.5 )
342               currentParticle.SetDefinitionAndUpdateE( anAntiLambda );
343             else
344               currentParticle.SetDefinitionAndUpdateE( anAntiSigmaZero );
345             incidentHasChanged = true;
346           }             
347           targetParticle.SetDefinitionAndUpdateE( aNeutron );
348           targetHasChanged = true;
349           break;
350         case 3:
351           targetParticle.SetDefinitionAndUpdateE( aNeutron );
352           targetHasChanged = true;
353           break;
354        }
355      }
356      else  // target must be a neutron
357      {
358        counter = -1;
359        for( np=0; np<numSec/3 && ran>=excs; ++np )
360        {
361          for( nm=np; nm<=(np+2) && ran>=excs; ++nm )
362          {
363            for( nz=0; nz<numSec/3 && ran>=excs; ++nz )
364            {
365              if( ++counter < numMul )
366              {
367                nt = np+nm+nz;
368                if( (nt>0) && (nt<=numSec) )
369                {
370                  test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
371                  dum = (pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n);
372                  if( std::fabs(dum) < 1.0 )
373                  {
374                    if( test >= 1.0e-10 )excs += dum*test;
375                  }
376                  else
377                    excs += dum*test;
378                }
379              }
380            }
381          }
382        }
383        if( ran >= excs )  // 3 previous loops continued to the end
384        {
385          quasiElastic = true;
386          return;
387        }
388        np--; nm--; nz--;
389        G4int ncht = std::min( 3, std::max( 1, np-nm+3 ) );
390        switch( ncht )
391        {
392         case 1:
393           if( G4UniformRand() < 0.5 )
394             currentParticle.SetDefinitionAndUpdateE( anAntiLambda );
395           else
396             currentParticle.SetDefinitionAndUpdateE( anAntiSigmaZero );
397           incidentHasChanged = true;
398           targetParticle.SetDefinitionAndUpdateE( aProton );
399           targetHasChanged = true;
400           break;
401         case 2:
402           if( G4UniformRand() < 0.5 )
403           {
404             if( G4UniformRand() < 0.5 )
405             {
406               currentParticle.SetDefinitionAndUpdateE( anAntiLambda );
407               incidentHasChanged = true;
408             }
409             else
410             {
411               targetParticle.SetDefinitionAndUpdateE( aProton );
412               targetHasChanged = true;
413             }
414           }
415           else
416           {
417             if( G4UniformRand() < 0.5 )
418             {
419               currentParticle.SetDefinitionAndUpdateE( anAntiSigmaZero );
420               incidentHasChanged = true;
421             }
422             else
423             {
424               targetParticle.SetDefinitionAndUpdateE( aProton );
425               targetHasChanged = true;
426             }
427           }
428           break;
429         case 3:
430           break;
431        }
432      }
433    }
434    else  // random number <= anhl[iplab]
435    {
436      if( centerofmassEnergy <= aPiPlus->GetPDGMass()/MeV+aKaonPlus->GetPDGMass()/MeV )
437      {
438        quasiElastic = true;
439        return;
440      }
441      G4double n, anpn;
442      GetNormalizationConstant( -centerofmassEnergy, n, anpn );
443      G4double ran = G4UniformRand();
444      G4double dum, excs = 0.0;
445      if( targetParticle.GetDefinition() == aProton )
446      {
447        counter = -1;
448        for( np=1; np<numSec/3 && ran>=excs; ++np )
449        {
450          nm = np;
451          for( nz=0; nz<numSec/3 && ran>=excs; ++nz )
452          {
453            if( ++counter < numMulA )
454            {
455              nt = np+nm+nz;
456              if( nt>1 && nt<=numSec )
457              {
458                test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
459                dum = (pi/anpn)*nt*protmulA[counter]*protnormA[nt-1]/(2.0*n*n);
460                if( std::fabs(dum) < 1.0 )
461                {
462                  if( test >= 1.0e-10 )excs += dum*test;
463                }
464                else
465                  excs += dum*test;
466              }
467            }
468          }
469        }
470        if( ran >= excs )  // 3 previous loops continued to the end
471        {
472          quasiElastic = true;
473          return;
474        }
475        np--; nz--;
476      }
477      else  // target must be a neutron
478      {
479        counter = -1;
480        for( np=0; np<numSec/3 && ran>=excs; ++np )
481        {
482          nm = np+1;
483          for( nz=0; nz<numSec/3 && ran>=excs; ++nz )
484          {
485            if( ++counter < numMulA )
486            {
487              nt = np+nm+nz;
488              if( nt>1 && nt<=numSec )
489              {
490                test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
491                dum = (pi/anpn)*nt*neutmulA[counter]*neutnormA[nt-1]/(2.0*n*n);
492                if( std::fabs(dum) < 1.0 )
493                {
494                  if( test >= 1.0e-10 )excs += dum*test;
495                }
496                else
497                  excs += dum*test;
498              }
499            }
500          }
501        }
502        if( ran >= excs )  // 3 previous loops continued to the end
503        {
504          quasiElastic = true;
505          return;
506        }
507        np--; nz--;
508      }
509      if( nz > 0 )
510      {
511        if( nm > 0 )
512        {
513          if( G4UniformRand() < 0.5 )
514          {
515            vec.Initialize( 1 );
516            G4ReactionProduct *p= new G4ReactionProduct;
517            p->SetDefinition( aKaonMinus );
518            (G4UniformRand() < 0.5) ? p->SetSide( -1 ) : p->SetSide( 1 );
519            vec.SetElement( vecLen++, p );
520            --nm;
521          }
522          else
523          {
524            vec.Initialize( 1 );
525            G4ReactionProduct *p= new G4ReactionProduct ;
526            p->SetDefinition( aKaonZL );
527            (G4UniformRand() < 0.5) ? p->SetSide( -1 ) : p->SetSide( 1 );
528            vec.SetElement( vecLen++, p );
529            --nz;
530          }
531        }
532        else   // nm == 0
533        {
534          vec.Initialize( 1 );
535          G4ReactionProduct *p = new G4ReactionProduct;
536          p->SetDefinition( aKaonZL );
537          (G4UniformRand() < 0.5) ? p->SetSide( -1 ) : p->SetSide( 1 );
538          vec.SetElement( vecLen++, p );
539          --nz;
540        }
541      }
542      else    //  nz == 0
543      {
544        if( nm > 0 )
545        {
546          vec.Initialize( 1 );
547          G4ReactionProduct *p = new G4ReactionProduct;
548          p->SetDefinition( aKaonMinus );
549          (G4UniformRand() < 0.5) ? p->SetSide( -1 ) : p->SetSide( 1 );
550          vec.SetElement( vecLen++, p );
551          --nm;
552        }
553      }
554      currentParticle.SetMass( 0.0 );
555      targetParticle.SetMass( 0.0 );
556    }
557
558  SetUpPions( np, nm, nz, vec, vecLen );
559  return;
560}
561
562 /* end of file */
563 
Note: See TracBrowser for help on using the repository browser.