G4LELambdaInelastic Class Reference

#include <G4LELambdaInelastic.hh>

Inheritance diagram for G4LELambdaInelastic:

G4InelasticInteraction G4HadronicInteraction

Public Member Functions

 G4LELambdaInelastic (const G4String &name="G4LELambdaInelastic")
 ~G4LELambdaInelastic ()
G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
virtual void ModelDescription (std::ostream &outFile) const

Detailed Description

Definition at line 44 of file G4LELambdaInelastic.hh.


Constructor & Destructor Documentation

G4LELambdaInelastic::G4LELambdaInelastic ( const G4String name = "G4LELambdaInelastic"  ) 

Definition at line 39 of file G4LELambdaInelastic.cc.

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

00040  :G4InelasticInteraction(name)
00041 {
00042   SetMinEnergy(0.0);
00043   SetMaxEnergy(25.*GeV);
00044   G4cout << "WARNING: model G4LELambdaInelastic is being deprecated and will\n"
00045          << "disappear in Geant4 version 10.0"  << G4endl;
00046 }

G4LELambdaInelastic::~G4LELambdaInelastic (  )  [inline]

Definition at line 50 of file G4LELambdaInelastic.hh.

00050 { }


Member Function Documentation

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

Implements G4HadronicInteraction.

Definition at line 64 of file G4LELambdaInelastic.cc.

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

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

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

Reimplemented from G4HadronicInteraction.

Definition at line 49 of file G4LELambdaInelastic.cc.

00050 {
00051   outFile << "G4LELambdaInelastic is one of the Low Energy Parameterized\n"
00052           << "(LEP) models used to implement inelastic Lambda scattering\n"
00053           << "from nuclei.  It is a re-engineered version of the GHEISHA\n"
00054           << "code of H. Fesefeldt.  It divides the initial collision\n"
00055           << "products into backward- and forward-going clusters which are\n"
00056           << "then decayed into final state hadrons.  The model does not\n"
00057           << "conserve energy on an event-by-event basis.  It may be\n"
00058           << "applied to lambdas with initial energies between 0 and 25\n"
00059           << "GeV.\n";
00060 }


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