G4RPGPiPlusInelastic Class Reference

#include <G4RPGPiPlusInelastic.hh>

Inheritance diagram for G4RPGPiPlusInelastic:

G4RPGPionInelastic G4RPGInelastic G4HadronicInteraction

Public Member Functions

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

Detailed Description

Definition at line 41 of file G4RPGPiPlusInelastic.hh.


Constructor & Destructor Documentation

G4RPGPiPlusInelastic::G4RPGPiPlusInelastic (  )  [inline]

Definition at line 45 of file G4RPGPiPlusInelastic.hh.

00045                           : G4RPGPionInelastic("RPGPiPlusInelastic")
00046    {}

G4RPGPiPlusInelastic::~G4RPGPiPlusInelastic (  )  [inline]

Definition at line 48 of file G4RPGPiPlusInelastic.hh.

00049    {}


Member Function Documentation

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

Implements G4HadronicInteraction.

Definition at line 34 of file G4RPGPiPlusInelastic.cc.

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

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


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