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