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