#include <G4LEAntiProtonInelastic.hh>
Inheritance diagram for G4LEAntiProtonInelastic:
Public Member Functions | |
G4LEAntiProtonInelastic (const G4String &name="G4LEAntiProtonInelastic") | |
~G4LEAntiProtonInelastic () | |
G4HadFinalState * | ApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus) |
virtual void | ModelDescription (std::ostream &outFile) const |
Definition at line 45 of file G4LEAntiProtonInelastic.hh.
G4LEAntiProtonInelastic::G4LEAntiProtonInelastic | ( | const G4String & | name = "G4LEAntiProtonInelastic" |
) |
Definition at line 42 of file G4LEAntiProtonInelastic.cc.
References G4cout, G4endl, G4HadronicInteraction::SetMaxEnergy(), and G4HadronicInteraction::SetMinEnergy().
00043 :G4InelasticInteraction(name) 00044 { 00045 SetMinEnergy(0.0); 00046 SetMaxEnergy(25.*GeV); 00047 G4cout << "WARNING: model G4LEAntiProtonInelastic is being deprecated and will\n" 00048 << "disappear in Geant4 version 10.0" << G4endl; 00049 }
G4LEAntiProtonInelastic::~G4LEAntiProtonInelastic | ( | ) | [inline] |
G4HadFinalState * G4LEAntiProtonInelastic::ApplyYourself | ( | const G4HadProjectile & | aTrack, | |
G4Nucleus & | targetNucleus | |||
) | [virtual] |
Implements G4HadronicInteraction.
Definition at line 67 of file G4LEAntiProtonInelastic.cc.
References G4InelasticInteraction::CalculateMomenta(), G4Nucleus::Cinema(), 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(), G4HadProjectile::GetTotalMomentum(), 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.
00069 { 00070 const G4HadProjectile* originalIncident = &aTrack; 00071 00072 if (originalIncident->GetKineticEnergy() <= 0.1*MeV) { 00073 theParticleChange.SetStatusChange(isAlive); 00074 theParticleChange.SetEnergyChange(aTrack.GetKineticEnergy()); 00075 theParticleChange.SetMomentumChange(aTrack.Get4Momentum().vect().unit()); 00076 return &theParticleChange; 00077 } 00078 00079 // create the target particle 00080 00081 G4DynamicParticle *originalTarget = targetNucleus.ReturnTargetParticle(); 00082 00083 if (verboseLevel > 1) { 00084 const G4Material *targetMaterial = aTrack.GetMaterial(); 00085 G4cout << "G4LEAntiProtonInelastic::ApplyYourself called" << G4endl; 00086 G4cout << "kinetic energy = " << originalIncident->GetKineticEnergy()/MeV << "MeV, "; 00087 G4cout << "target material = " << targetMaterial->GetName() << ", "; 00088 G4cout << "target particle = " << originalTarget->GetDefinition()->GetParticleName() 00089 << G4endl; 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 tkin = targetNucleus.EvaporationEffects( ek ); 00113 ek -= tkin; 00114 modifiedOriginal.SetKineticEnergy( ek*MeV ); 00115 et = ek + amas; 00116 p = std::sqrt( std::abs((et-amas)*(et+amas)) ); 00117 pp = modifiedOriginal.GetMomentum().mag()/MeV; 00118 if (pp > 0.0) { 00119 G4ThreeVector momentum = modifiedOriginal.GetMomentum(); 00120 modifiedOriginal.SetMomentum( momentum * (p/pp) ); 00121 } 00122 G4ReactionProduct currentParticle = modifiedOriginal; 00123 G4ReactionProduct targetParticle; 00124 targetParticle = *originalTarget; 00125 currentParticle.SetSide( 1 ); // incident always goes in forward hemisphere 00126 targetParticle.SetSide( -1 ); // target always goes in backward hemisphere 00127 G4bool incidentHasChanged = false; 00128 G4bool targetHasChanged = false; 00129 G4bool quasiElastic = false; 00130 G4FastVector<G4ReactionProduct,GHADLISTSIZE> vec; // vec will contain the secondary particles 00131 G4int vecLen = 0; 00132 vec.Initialize( 0 ); 00133 00134 const G4double cutOff = 0.1; 00135 const G4double anni = std::min( 1.3*originalIncident->GetTotalMomentum()/GeV, 0.4 ); 00136 00137 if ((currentParticle.GetKineticEnergy()/MeV > cutOff) || 00138 (G4UniformRand() > anni) ) 00139 Cascade(vec, vecLen, originalIncident, currentParticle, targetParticle, 00140 incidentHasChanged, targetHasChanged, quasiElastic); 00141 else 00142 quasiElastic = true; 00143 00144 CalculateMomenta(vec, vecLen, originalIncident, originalTarget, 00145 modifiedOriginal, targetNucleus, currentParticle, 00146 targetParticle, incidentHasChanged, targetHasChanged, 00147 quasiElastic); 00148 00149 SetUpChange(vec, vecLen, currentParticle, targetParticle, incidentHasChanged); 00150 00151 if (isotopeProduction) DoIsotopeCounting(originalIncident, targetNucleus); 00152 00153 delete originalTarget; 00154 return &theParticleChange; 00155 }
void G4LEAntiProtonInelastic::ModelDescription | ( | std::ostream & | outFile | ) | const [virtual] |
Reimplemented from G4HadronicInteraction.
Definition at line 52 of file G4LEAntiProtonInelastic.cc.
00053 { 00054 outFile << "G4LEAntiProtonInelastic is one of the Low Energy Parameterized\n" 00055 << "(LEP) models used to implement inelastic anti-proton scattering\n" 00056 << "from nuclei. It is a re-engineered version of the GHEISHA\n" 00057 << "code of H. Fesefeldt. It divides the initial collision\n" 00058 << "products into backward- and forward-going clusters which are\n" 00059 << "then decayed into final state hadrons. The model does not\n" 00060 << "conserve energy on an event-by-event basis. It may be\n" 00061 << "applied to anti-protons with initial energies between 0 and 25\n" 00062 << "GeV.\n"; 00063 }