G4RPGProtonInelastic Class Reference

#include <G4RPGProtonInelastic.hh>

Inheritance diagram for G4RPGProtonInelastic:

G4RPGNucleonInelastic G4RPGInelastic G4HadronicInteraction

Public Member Functions

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

Detailed Description

Definition at line 43 of file G4RPGProtonInelastic.hh.


Constructor & Destructor Documentation

G4RPGProtonInelastic::G4RPGProtonInelastic (  )  [inline]

Definition at line 47 of file G4RPGProtonInelastic.hh.

00047                           : G4RPGNucleonInelastic("RPGProtonInelastic")
00048    {}

G4RPGProtonInelastic::~G4RPGProtonInelastic (  )  [inline]

Definition at line 50 of file G4RPGProtonInelastic.hh.

00051    {}


Member Function Documentation

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

Implements G4HadronicInteraction.

Definition at line 34 of file G4RPGProtonInelastic.cc.

References G4RPGInelastic::CalculateMomenta(), G4Nucleus::Cinema(), G4HadFinalState::Clear(), G4Nucleus::EvaporationEffects(), G4UniformRand, G4HadProjectile::Get4Momentum(), G4HadProjectile::GetDefinition(), G4ReactionProduct::GetKineticEnergy(), G4HadProjectile::GetKineticEnergy(), G4ReactionProduct::GetMomentum(), G4ParticleDefinition::GetPDGMass(), G4FastVector< Type, N >::Initialize(), isAlive, G4InuclParticleNames::pp, G4Nucleus::ReturnTargetParticle(), G4HadFinalState::SetEnergyChange(), G4ReactionProduct::SetKineticEnergy(), G4ReactionProduct::SetMomentum(), G4HadFinalState::SetMomentumChange(), G4ReactionProduct::SetSide(), G4HadFinalState::SetStatusChange(), G4RPGInelastic::SetUpChange(), and G4HadronicInteraction::theParticleChange.

00036 {
00037   theParticleChange.Clear();
00038   const G4HadProjectile *originalIncident = &aTrack;
00039   if (originalIncident->GetKineticEnergy()<= 0.1) 
00040   {
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 (originalIncident->GetKineticEnergy()/GeV < 0.01+2.*G4UniformRand()/9. )
00053   {
00054     SlowProton( originalIncident, targetNucleus );
00055     delete originalTarget;
00056     return &theParticleChange;
00057   }
00058 
00059   // Fermi motion and evaporation
00060   // As of Geant3, the Fermi energy calculation had not been Done
00061 
00062   G4double ek = originalIncident->GetKineticEnergy();
00063   G4double amas = originalIncident->GetDefinition()->GetPDGMass();
00064   G4ReactionProduct modifiedOriginal;
00065   modifiedOriginal = *originalIncident;
00066     
00067   G4double tkin = targetNucleus.Cinema( ek );
00068   ek += tkin;
00069   modifiedOriginal.SetKineticEnergy(ek);
00070   G4double et = ek + amas;
00071   G4double p = std::sqrt( std::abs((et-amas)*(et+amas)) );
00072   G4double pp = modifiedOriginal.GetMomentum().mag();
00073   if (pp > 0.0) {
00074     G4ThreeVector momentum = modifiedOriginal.GetMomentum();
00075     modifiedOriginal.SetMomentum( momentum * (p/pp) );
00076   }
00077   //
00078   // calculate black track energies
00079   //
00080   tkin = targetNucleus.EvaporationEffects(ek);
00081   ek -= tkin;
00082   modifiedOriginal.SetKineticEnergy(ek);
00083   et = ek + amas;
00084   p = std::sqrt( std::abs((et-amas)*(et+amas)) );
00085   pp = modifiedOriginal.GetMomentum().mag();
00086   if (pp > 0.0) {
00087     G4ThreeVector momentum = modifiedOriginal.GetMomentum();
00088     modifiedOriginal.SetMomentum( momentum * (p/pp) );
00089   }
00090   const G4double cutOff = 0.1;
00091   if (modifiedOriginal.GetKineticEnergy() < cutOff) {
00092     SlowProton( originalIncident, targetNucleus );
00093     delete originalTarget;
00094     return &theParticleChange;
00095   }
00096 
00097   G4ReactionProduct currentParticle = modifiedOriginal;
00098   G4ReactionProduct targetParticle;
00099   targetParticle = *originalTarget;
00100   currentParticle.SetSide( 1 ); // incident always goes in forward hemisphere
00101   targetParticle.SetSide( -1 );  // target always goes in backward hemisphere
00102   G4bool incidentHasChanged = false;
00103   G4bool targetHasChanged = false;
00104   G4bool quasiElastic = false;
00105   G4FastVector<G4ReactionProduct,256> vec;  // vec will contain the sec. particles
00106   G4int vecLen = 0;
00107   vec.Initialize( 0 );
00108 
00109   InitialCollision(vec, vecLen, currentParticle, targetParticle,
00110                    incidentHasChanged, targetHasChanged);
00111 
00112   CalculateMomenta(vec, vecLen,
00113                    originalIncident, originalTarget, modifiedOriginal,
00114                    targetNucleus, currentParticle, targetParticle,
00115                    incidentHasChanged, targetHasChanged, quasiElastic);
00116 
00117   SetUpChange( vec, vecLen,
00118                currentParticle, targetParticle,
00119                incidentHasChanged );
00120 
00121   delete originalTarget;
00122   return &theParticleChange;
00123 }


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