G4LESigmaMinusInelastic Class Reference

#include <G4LESigmaMinusInelastic.hh>

Inheritance diagram for G4LESigmaMinusInelastic:

G4InelasticInteraction G4HadronicInteraction

Public Member Functions

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

Detailed Description

Definition at line 45 of file G4LESigmaMinusInelastic.hh.


Constructor & Destructor Documentation

G4LESigmaMinusInelastic::G4LESigmaMinusInelastic (  )  [inline]

Definition at line 49 of file G4LESigmaMinusInelastic.hh.

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

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

G4LESigmaMinusInelastic::~G4LESigmaMinusInelastic (  )  [inline]

Definition at line 57 of file G4LESigmaMinusInelastic.hh.

00057 {}


Member Function Documentation

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

Implements G4HadronicInteraction.

Definition at line 51 of file G4LESigmaMinusInelastic.cc.

References G4InelasticInteraction::CalculateMomenta(), G4Nucleus::Cinema(), G4InelasticInteraction::DoIsotopeCounting(), G4Nucleus::EvaporationEffects(), G4cout, G4endl, G4HadProjectile::Get4Momentum(), G4HadProjectile::GetDefinition(), G4DynamicParticle::GetDefinition(), G4HadProjectile::GetKineticEnergy(), G4HadProjectile::GetMaterial(), G4ReactionProduct::GetMomentum(), G4Material::GetName(), G4ParticleDefinition::GetParticleName(), G4ParticleDefinition::GetPDGMass(), 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.

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

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

Reimplemented from G4HadronicInteraction.

Definition at line 37 of file G4LESigmaMinusInelastic.cc.

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


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