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