G4RPGSigmaPlusInelastic Class Reference

#include <G4RPGSigmaPlusInelastic.hh>

Inheritance diagram for G4RPGSigmaPlusInelastic:

G4RPGInelastic G4HadronicInteraction

Public Member Functions

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

Detailed Description

Definition at line 42 of file G4RPGSigmaPlusInelastic.hh.


Constructor & Destructor Documentation

G4RPGSigmaPlusInelastic::G4RPGSigmaPlusInelastic (  )  [inline]

Definition at line 46 of file G4RPGSigmaPlusInelastic.hh.

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

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

G4RPGSigmaPlusInelastic::~G4RPGSigmaPlusInelastic (  )  [inline]

Definition at line 52 of file G4RPGSigmaPlusInelastic.hh.

00052 { }


Member Function Documentation

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

Implements G4HadronicInteraction.

Definition at line 35 of file G4RPGSigmaPlusInelastic.cc.

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