G4LEAntiSigmaMinusInelastic Class Reference

#include <G4LEAntiSigmaMinusInelastic.hh>

Inheritance diagram for G4LEAntiSigmaMinusInelastic:

G4InelasticInteraction G4HadronicInteraction

Public Member Functions

 G4LEAntiSigmaMinusInelastic ()
 ~G4LEAntiSigmaMinusInelastic ()
G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
virtual void ModelDescription (std::ostream &outFile) const

Detailed Description

Definition at line 45 of file G4LEAntiSigmaMinusInelastic.hh.


Constructor & Destructor Documentation

G4LEAntiSigmaMinusInelastic::G4LEAntiSigmaMinusInelastic (  )  [inline]

Definition at line 49 of file G4LEAntiSigmaMinusInelastic.hh.

References G4cout, G4endl, G4HadronicInteraction::SetMaxEnergy(), and G4HadronicInteraction::SetMinEnergy().

00049                                   : G4InelasticInteraction("G4LEAntiSigmaMinusInelastic")
00050     {
00051       SetMinEnergy(0.0);
00052       SetMaxEnergy(25.*CLHEP::GeV);
00053       G4cout << "WARNING: model G4LEAntiSigmaMinusInelastic is being deprecated and will\n"
00054              << "disappear in Geant4 version 10.0"  << G4endl;
00055     }

G4LEAntiSigmaMinusInelastic::~G4LEAntiSigmaMinusInelastic (  )  [inline]

Definition at line 57 of file G4LEAntiSigmaMinusInelastic.hh.

00057 { }


Member Function Documentation

G4HadFinalState * G4LEAntiSigmaMinusInelastic::ApplyYourself ( const G4HadProjectile aTrack,
G4Nucleus targetNucleus 
) [virtual]

Implements G4HadronicInteraction.

Definition at line 53 of file G4LEAntiSigmaMinusInelastic.cc.

References G4InelasticInteraction::CalculateMomenta(), G4Nucleus::Cinema(), G4InelasticInteraction::DoIsotopeCounting(), G4Nucleus::EvaporationEffects(), G4cout, G4endl, G4UniformRand, G4HadProjectile::Get4Momentum(), G4HadProjectile::GetDefinition(), G4DynamicParticle::GetDefinition(), G4ReactionProduct::GetKineticEnergy(), G4HadProjectile::GetKineticEnergy(), G4HadProjectile::GetMaterial(), G4ReactionProduct::GetMomentum(), G4Material::GetName(), G4ParticleDefinition::GetParticleName(), G4ParticleDefinition::GetPDGMass(), G4ReactionProduct::GetTotalMomentum(), G4FastVector< Type, N >::Initialize(), isAlive, G4InelasticInteraction::isotopeProduction, G4InuclParticleNames::pp, G4Nucleus::ReturnTargetParticle(), G4HadFinalState::SetEnergyChange(), G4ReactionProduct::SetKineticEnergy(), G4ReactionProduct::SetMomentum(), G4HadFinalState::SetMomentumChange(), G4ReactionProduct::SetSide(), G4HadFinalState::SetStatusChange(), G4InelasticInteraction::SetUpChange(), G4HadronicInteraction::theParticleChange, and G4HadronicInteraction::verboseLevel.

00055 { 
00056   const G4HadProjectile *originalIncident = &aTrack;
00057   if (originalIncident->GetKineticEnergy()<= 0.1*MeV) {
00058     theParticleChange.SetStatusChange(isAlive);
00059     theParticleChange.SetEnergyChange(aTrack.GetKineticEnergy());
00060     theParticleChange.SetMomentumChange(aTrack.Get4Momentum().vect().unit()); 
00061     return &theParticleChange;      
00062   }
00063 
00064   // create the target particle
00065   G4DynamicParticle* originalTarget = targetNucleus.ReturnTargetParticle();
00066     
00067   if (verboseLevel > 1) {
00068     const G4Material *targetMaterial = aTrack.GetMaterial();
00069     G4cout << "G4LEAntiSigmaMinusInelastic::ApplyYourself called" << G4endl;
00070     G4cout << "kinetic energy = " << originalIncident->GetKineticEnergy()/MeV << "MeV, ";
00071     G4cout << "target material = " << targetMaterial->GetName() << ", ";
00072     G4cout << "target particle = " << originalTarget->GetDefinition()->GetParticleName()
00073            << G4endl;
00074   }
00075 
00076   // Fermi motion and evaporation
00077   // As of Geant3, the Fermi energy calculation had not been Done
00078   G4double ek = originalIncident->GetKineticEnergy()/MeV;
00079   G4double amas = originalIncident->GetDefinition()->GetPDGMass()/MeV;
00080   G4ReactionProduct modifiedOriginal;
00081   modifiedOriginal = *originalIncident;
00082     
00083   G4double tkin = targetNucleus.Cinema( ek );
00084   ek += tkin;
00085   modifiedOriginal.SetKineticEnergy( ek*MeV );
00086   G4double et = ek + amas;
00087   G4double p = std::sqrt( std::abs((et-amas)*(et+amas)) );
00088   G4double pp = modifiedOriginal.GetMomentum().mag()/MeV;
00089   if (pp > 0.0) {
00090     G4ThreeVector momentum = modifiedOriginal.GetMomentum();
00091     modifiedOriginal.SetMomentum( momentum * (p/pp) );
00092   }
00093 
00094   // calculate black track energies
00095   tkin = targetNucleus.EvaporationEffects( ek );
00096   ek -= tkin;
00097   modifiedOriginal.SetKineticEnergy( ek*MeV );
00098   et = ek + amas;
00099   p = std::sqrt( std::abs((et-amas)*(et+amas)) );
00100   pp = modifiedOriginal.GetMomentum().mag()/MeV;
00101   if (pp > 0.0) {
00102     G4ThreeVector momentum = modifiedOriginal.GetMomentum();
00103     modifiedOriginal.SetMomentum( momentum * (p/pp) );
00104   }
00105   G4ReactionProduct currentParticle = modifiedOriginal;
00106   G4ReactionProduct targetParticle;
00107   targetParticle = *originalTarget;
00108   currentParticle.SetSide( 1 ); // incident always goes in forward hemisphere
00109   targetParticle.SetSide( -1 );  // target always goes in backward hemisphere
00110   G4bool incidentHasChanged = false;
00111   G4bool targetHasChanged = false;
00112   G4bool quasiElastic = false;
00113   G4FastVector<G4ReactionProduct,GHADLISTSIZE> vec;  // vec will contain the secondary particles
00114   G4int vecLen = 0;
00115   vec.Initialize(0);
00116     
00117   const G4double cutOff = 0.1;
00118   const G4double anni = std::min( 1.3*currentParticle.GetTotalMomentum()/GeV, 0.4 );
00119   if ((currentParticle.GetKineticEnergy()/MeV > cutOff) ||
00120       (G4UniformRand() > anni) )
00121     Cascade(vec, vecLen, originalIncident, currentParticle, targetParticle,
00122             incidentHasChanged, targetHasChanged, quasiElastic);
00123     
00124   CalculateMomenta(vec, vecLen, originalIncident, originalTarget,
00125                    modifiedOriginal, targetNucleus, currentParticle,
00126                    targetParticle, incidentHasChanged, targetHasChanged,
00127                    quasiElastic);
00128     
00129   SetUpChange(vec, vecLen, currentParticle, targetParticle, incidentHasChanged);
00130 
00131   if (isotopeProduction) DoIsotopeCounting(originalIncident, targetNucleus);
00132 
00133   delete originalTarget;
00134   return &theParticleChange;
00135 }

void G4LEAntiSigmaMinusInelastic::ModelDescription ( std::ostream &  outFile  )  const [virtual]

Reimplemented from G4HadronicInteraction.

Definition at line 38 of file G4LEAntiSigmaMinusInelastic.cc.

00039 {
00040   outFile << "G4LEAntiSigmaMinusInelastic is one of the Low Energy\n"
00041           << "Parameterized (LEP) models used to implement inelastic\n"
00042           << "antiSigma- scattering from nuclei.  It is a re-engineered\n"
00043           << "version of the GHEISHA code of H. Fesefeldt.  It divides the\n"
00044           << "initial collision products into backward- and forward-going\n"
00045           << "clusters which are then decayed into final state hadrons.  The\n"
00046           << "model does not conserve energy on an event-by-event basis.  It\n"
00047           << "may be applied to antiSigma- with initial energies between 0\n"
00048           << "and 25 GeV.\n";
00049 }


The documentation for this class was generated from the following files:
Generated on Mon May 27 17:52:22 2013 for Geant4 by  doxygen 1.4.7