G4LEKaonPlusInelastic Class Reference

#include <G4LEKaonPlusInelastic.hh>

Inheritance diagram for G4LEKaonPlusInelastic:

G4InelasticInteraction G4HadronicInteraction

Public Member Functions

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

Detailed Description

Definition at line 44 of file G4LEKaonPlusInelastic.hh.


Constructor & Destructor Documentation

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

Definition at line 39 of file G4LEKaonPlusInelastic.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 G4LEKaonPlusInelastic is being deprecated and will\n"
00045          << "disappear in Geant4 version 10.0"  << G4endl;
00046 }

G4LEKaonPlusInelastic::~G4LEKaonPlusInelastic (  )  [inline]

Definition at line 50 of file G4LEKaonPlusInelastic.hh.

00050 {}


Member Function Documentation

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

Implements G4HadronicInteraction.

Definition at line 64 of file G4LEKaonPlusInelastic.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.

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

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

Reimplemented from G4HadronicInteraction.

Definition at line 49 of file G4LEKaonPlusInelastic.cc.

00050 {
00051   outFile << "G4LEKaonPlusInelastic is one of the Low Energy Parameterized\n"
00052           << "(LEP) models used to implement inelastic K+ 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 kaons 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