G4RPGAntiProtonInelastic Class Reference

#include <G4RPGAntiProtonInelastic.hh>

Inheritance diagram for G4RPGAntiProtonInelastic:

G4RPGInelastic G4HadronicInteraction

Public Member Functions

 G4RPGAntiProtonInelastic ()
 ~G4RPGAntiProtonInelastic ()
G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)

Detailed Description

Definition at line 42 of file G4RPGAntiProtonInelastic.hh.


Constructor & Destructor Documentation

G4RPGAntiProtonInelastic::G4RPGAntiProtonInelastic (  )  [inline]

Definition at line 46 of file G4RPGAntiProtonInelastic.hh.

References G4HadronicInteraction::SetMaxEnergy(), and G4HadronicInteraction::SetMinEnergy().

00046                               :G4RPGInelastic("G4RPGAntiProtonInelastic")
00047     {
00048       SetMinEnergy( 0.0 );
00049       SetMaxEnergy( 25.*CLHEP::GeV );
00050     }

G4RPGAntiProtonInelastic::~G4RPGAntiProtonInelastic (  )  [inline]

Definition at line 52 of file G4RPGAntiProtonInelastic.hh.

00052 { }


Member Function Documentation

G4HadFinalState * G4RPGAntiProtonInelastic::ApplyYourself ( const G4HadProjectile aTrack,
G4Nucleus targetNucleus 
) [virtual]

Implements G4HadronicInteraction.

Definition at line 35 of file G4RPGAntiProtonInelastic.cc.

References G4RPGInelastic::CalculateMomenta(), G4Nucleus::Cinema(), 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, G4InuclParticleNames::pp, G4Nucleus::ReturnTargetParticle(), G4HadFinalState::SetEnergyChange(), G4ReactionProduct::SetKineticEnergy(), G4ReactionProduct::SetMomentum(), G4HadFinalState::SetMomentumChange(), G4ReactionProduct::SetSide(), G4HadFinalState::SetStatusChange(), G4RPGInelastic::SetUpChange(), G4HadronicInteraction::theParticleChange, and G4HadronicInteraction::verboseLevel.

00037 { 
00038   const G4HadProjectile* originalIncident = &aTrack;
00039 
00040   if (originalIncident->GetKineticEnergy() <= 0.1*MeV) {
00041     theParticleChange.SetStatusChange(isAlive);
00042     theParticleChange.SetEnergyChange(aTrack.GetKineticEnergy());
00043     theParticleChange.SetMomentumChange(aTrack.Get4Momentum().vect().unit());
00044     return &theParticleChange;
00045   }
00046 
00047     //
00048     // create the target particle
00049     //
00050     G4DynamicParticle *originalTarget = targetNucleus.ReturnTargetParticle();
00051     
00052     if( verboseLevel > 1 )
00053     {
00054       const G4Material *targetMaterial = aTrack.GetMaterial();
00055       G4cout << "G4RPGAntiProtonInelastic::ApplyYourself called" << G4endl;
00056       G4cout << "kinetic energy = " << originalIncident->GetKineticEnergy()/MeV << "MeV, ";
00057       G4cout << "target material = " << targetMaterial->GetName() << ", ";
00058       G4cout << "target particle = " << originalTarget->GetDefinition()->GetParticleName()
00059            << G4endl;
00060     }
00061     //
00062     // Fermi motion and evaporation
00063     // As of Geant3, the Fermi energy calculation had not been Done
00064     //
00065     G4double ek = originalIncident->GetKineticEnergy()/MeV;
00066     G4double amas = originalIncident->GetDefinition()->GetPDGMass()/MeV;
00067     G4ReactionProduct modifiedOriginal;
00068     modifiedOriginal = *originalIncident;
00069     
00070     G4double tkin = targetNucleus.Cinema( ek );
00071     ek += tkin;
00072     modifiedOriginal.SetKineticEnergy( ek*MeV );
00073     G4double et = ek + amas;
00074     G4double p = std::sqrt( std::abs((et-amas)*(et+amas)) );
00075     G4double pp = modifiedOriginal.GetMomentum().mag()/MeV;
00076     if( pp > 0.0 )
00077     {
00078       G4ThreeVector momentum = modifiedOriginal.GetMomentum();
00079       modifiedOriginal.SetMomentum( momentum * (p/pp) );
00080     }
00081     //
00082     // calculate black track energies
00083     //
00084     tkin = targetNucleus.EvaporationEffects( ek );
00085     ek -= tkin;
00086     modifiedOriginal.SetKineticEnergy( ek*MeV );
00087     et = ek + amas;
00088     p = std::sqrt( std::abs((et-amas)*(et+amas)) );
00089     pp = modifiedOriginal.GetMomentum().mag()/MeV;
00090     if( pp > 0.0 )
00091     {
00092       G4ThreeVector momentum = modifiedOriginal.GetMomentum();
00093       modifiedOriginal.SetMomentum( momentum * (p/pp) );
00094     }
00095     G4ReactionProduct currentParticle = modifiedOriginal;
00096     G4ReactionProduct targetParticle;
00097     targetParticle = *originalTarget;
00098     currentParticle.SetSide( 1 ); // incident always goes in forward hemisphere
00099     targetParticle.SetSide( -1 );  // target always goes in backward hemisphere
00100     G4bool incidentHasChanged = false;
00101     G4bool targetHasChanged = false;
00102     G4bool quasiElastic = false;
00103     G4FastVector<G4ReactionProduct,GHADLISTSIZE> vec;  // vec will contain the secondary particles
00104     G4int vecLen = 0;
00105     vec.Initialize( 0 );
00106     
00107     const G4double cutOff = 0.1;
00108     const G4double anni = std::min( 1.3*originalIncident->GetTotalMomentum()/GeV, 0.4 );
00109     
00110     if( (currentParticle.GetKineticEnergy()/MeV > cutOff) ||
00111         (G4UniformRand() > anni) )
00112       Cascade( vec, vecLen,
00113                originalIncident, currentParticle, targetParticle,
00114                incidentHasChanged, targetHasChanged, quasiElastic );
00115     else
00116       quasiElastic = true;
00117     
00118     CalculateMomenta( vec, vecLen,
00119                       originalIncident, originalTarget, modifiedOriginal,
00120                       targetNucleus, currentParticle, targetParticle,
00121                       incidentHasChanged, targetHasChanged, quasiElastic );
00122     
00123     SetUpChange( vec, vecLen,
00124                  currentParticle, targetParticle,
00125                  incidentHasChanged );
00126     
00127   delete originalTarget;
00128   return &theParticleChange;
00129 }


The documentation for this class was generated from the following files:
Generated on Mon May 27 17:53:18 2013 for Geant4 by  doxygen 1.4.7