G4LEKaonZeroInelastic Class Reference

#include <G4LEKaonZeroInelastic.hh>

Inheritance diagram for G4LEKaonZeroInelastic:

G4InelasticInteraction G4HadronicInteraction

Public Member Functions

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

Detailed Description

Definition at line 45 of file G4LEKaonZeroInelastic.hh.


Constructor & Destructor Documentation

G4LEKaonZeroInelastic::G4LEKaonZeroInelastic (  )  [inline]

Definition at line 49 of file G4LEKaonZeroInelastic.hh.

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

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

G4LEKaonZeroInelastic::~G4LEKaonZeroInelastic (  )  [inline]

Definition at line 57 of file G4LEKaonZeroInelastic.hh.

00058     { }


Member Function Documentation

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

Implements G4HadronicInteraction.

Definition at line 53 of file G4LEKaonZeroInelastic.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.

Referenced by G4LEKaonZeroSInelastic::ApplyYourself(), and G4LEKaonZeroLInelastic::ApplyYourself().

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

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

Reimplemented from G4HadronicInteraction.

Definition at line 39 of file G4LEKaonZeroInelastic.cc.

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