G4LEKaonMinusInelastic Class Reference

#include <G4LEKaonMinusInelastic.hh>

Inheritance diagram for G4LEKaonMinusInelastic:

G4InelasticInteraction G4HadronicInteraction

Public Member Functions

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

Detailed Description

Definition at line 43 of file G4LEKaonMinusInelastic.hh.


Constructor & Destructor Documentation

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

Definition at line 41 of file G4LEKaonMinusInelastic.cc.

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

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

G4LEKaonMinusInelastic::~G4LEKaonMinusInelastic (  )  [inline]

Definition at line 49 of file G4LEKaonMinusInelastic.hh.

00049 {}


Member Function Documentation

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

Implements G4HadronicInteraction.

Definition at line 66 of file G4LEKaonMinusInelastic.cc.

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

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

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

Reimplemented from G4HadronicInteraction.

Definition at line 51 of file G4LEKaonMinusInelastic.cc.

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


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