G4RPGKPlusInelastic Class Reference

#include <G4RPGKPlusInelastic.hh>

Inheritance diagram for G4RPGKPlusInelastic:

G4RPGInelastic G4HadronicInteraction

Public Member Functions

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

Detailed Description

Definition at line 43 of file G4RPGKPlusInelastic.hh.


Constructor & Destructor Documentation

G4RPGKPlusInelastic::G4RPGKPlusInelastic (  )  [inline]

Definition at line 47 of file G4RPGKPlusInelastic.hh.

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

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

G4RPGKPlusInelastic::~G4RPGKPlusInelastic (  )  [inline]

Definition at line 53 of file G4RPGKPlusInelastic.hh.

00054     { }


Member Function Documentation

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

Implements G4HadronicInteraction.

Definition at line 36 of file G4RPGKPlusInelastic.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.

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


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