#include <G4LEOmegaMinusInelastic.hh>
Inheritance diagram for G4LEOmegaMinusInelastic:
Public Member Functions | |
G4LEOmegaMinusInelastic () | |
~G4LEOmegaMinusInelastic () | |
G4HadFinalState * | ApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus) |
virtual void | ModelDescription (std::ostream &outFile) const |
Definition at line 45 of file G4LEOmegaMinusInelastic.hh.
G4LEOmegaMinusInelastic::G4LEOmegaMinusInelastic | ( | ) | [inline] |
Definition at line 49 of file G4LEOmegaMinusInelastic.hh.
References G4cout, G4endl, G4HadronicInteraction::SetMaxEnergy(), and G4HadronicInteraction::SetMinEnergy().
00049 : G4InelasticInteraction("G4LEOmegaMinusInelastic") 00050 { 00051 SetMinEnergy(0.0); 00052 SetMaxEnergy(25.*CLHEP::GeV); 00053 G4cout << "WARNING: model G4LEOmegaMinusInelastic is being deprecated and will\n" 00054 << "disappear in Geant4 version 10.0" << G4endl; 00055 }
G4LEOmegaMinusInelastic::~G4LEOmegaMinusInelastic | ( | ) | [inline] |
G4HadFinalState * G4LEOmegaMinusInelastic::ApplyYourself | ( | const G4HadProjectile & | aTrack, | |
G4Nucleus & | targetNucleus | |||
) | [virtual] |
Implements G4HadronicInteraction.
Definition at line 51 of file G4LEOmegaMinusInelastic.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.
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 G4ReactionProduct targetParticle( originalTarget->GetDefinition() ); 00065 00066 if (verboseLevel > 1) { 00067 const G4Material *targetMaterial = aTrack.GetMaterial(); 00068 G4cout << "G4LEOmegaMinusInelastic::ApplyYourself called" << G4endl; 00069 G4cout << "kinetic energy = " << originalIncident->GetKineticEnergy() << "MeV, "; 00070 G4cout << "target material = " << targetMaterial->GetName() << ", "; 00071 G4cout << "target particle = " << originalTarget->GetDefinition()->GetParticleName() 00072 << G4endl; 00073 } 00074 G4ReactionProduct currentParticle( const_cast<G4ParticleDefinition *>(originalIncident->GetDefinition() )); 00075 currentParticle.SetMomentum( originalIncident->Get4Momentum().vect() ); 00076 currentParticle.SetKineticEnergy( originalIncident->GetKineticEnergy() ); 00077 00078 // Fermi motion and evaporation 00079 // As of Geant3, the Fermi energy calculation had not been Done 00080 G4double ek = originalIncident->GetKineticEnergy(); 00081 G4double amas = originalIncident->GetDefinition()->GetPDGMass(); 00082 00083 G4double tkin = targetNucleus.Cinema(ek); 00084 ek += tkin; 00085 currentParticle.SetKineticEnergy(ek); 00086 G4double et = ek + amas; 00087 G4double p = std::sqrt( std::abs((et-amas)*(et+amas)) ); 00088 G4double pp = currentParticle.GetMomentum().mag(); 00089 if (pp > 0.0) { 00090 G4ThreeVector momentum = currentParticle.GetMomentum(); 00091 currentParticle.SetMomentum( momentum * (p/pp) ); 00092 } 00093 00094 // calculate black track energies 00095 tkin = targetNucleus.EvaporationEffects( ek ); 00096 ek -= tkin; 00097 currentParticle.SetKineticEnergy( ek ); 00098 et = ek + amas; 00099 p = std::sqrt( std::abs((et-amas)*(et+amas)) ); 00100 pp = currentParticle.GetMomentum().mag(); 00101 if (pp > 0.0) { 00102 G4ThreeVector momentum = currentParticle.GetMomentum(); 00103 currentParticle.SetMomentum( momentum * (p/pp) ); 00104 } 00105 00106 G4ReactionProduct modifiedOriginal = currentParticle; 00107 00108 currentParticle.SetSide(1); // incident always goes in forward hemisphere 00109 targetParticle.SetSide(-1); // target always goes in backward hemisphere 00110 G4bool incidentHasChanged = false; 00111 G4bool targetHasChanged = false; 00112 G4bool quasiElastic = false; 00113 G4FastVector<G4ReactionProduct,GHADLISTSIZE> vec; // vec will contain the secondary particles 00114 G4int vecLen = 0; 00115 vec.Initialize( 0 ); 00116 00117 const G4double cutOff = 0.1*MeV; 00118 if (currentParticle.GetKineticEnergy() > cutOff) 00119 Cascade(vec, vecLen, originalIncident, currentParticle, targetParticle, 00120 incidentHasChanged, targetHasChanged, quasiElastic); 00121 00122 CalculateMomenta(vec, vecLen, originalIncident, originalTarget, 00123 modifiedOriginal, targetNucleus, currentParticle, 00124 targetParticle, incidentHasChanged, targetHasChanged, 00125 quasiElastic); 00126 00127 SetUpChange(vec, vecLen, currentParticle, targetParticle, incidentHasChanged); 00128 00129 if (isotopeProduction) DoIsotopeCounting(originalIncident, targetNucleus); 00130 00131 delete originalTarget; 00132 return &theParticleChange; 00133 }
void G4LEOmegaMinusInelastic::ModelDescription | ( | std::ostream & | outFile | ) | const [virtual] |
Reimplemented from G4HadronicInteraction.
Definition at line 37 of file G4LEOmegaMinusInelastic.cc.
00038 { 00039 outFile << "G4LEOmegaMinusInelastic is one of the Low Energy Parameterized\n" 00040 << "(LEP) models used to implement inelastic Omega- scattering from\n" 00041 << "nuclei. It is a re-engineered version of the GHEISHA code of\n" 00042 << "H. Fesefeldt. It divides the initial collision products into\n" 00043 << "backward- and forward-going clusters which are then decayed\n" 00044 << "into final state hadrons. The model does not conserve energy\n" 00045 << "on an event-by-event basis. It may be applied to Omega- with\n" 00046 << "initial energies between 0 and 25 GeV.\n"; 00047 }