Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
G4RPGKZeroInelastic Class Reference

#include <G4RPGKZeroInelastic.hh>

Inheritance diagram for G4RPGKZeroInelastic:
G4RPGInelastic G4HadronicInteraction

Public Member Functions

 G4RPGKZeroInelastic ()
 
 ~G4RPGKZeroInelastic ()
 
G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
 
- Public Member Functions inherited from G4RPGInelastic
 G4RPGInelastic (const G4String &modelName="RPGInelastic")
 
virtual ~G4RPGInelastic ()
 
- Public Member Functions inherited from G4HadronicInteraction
 G4HadronicInteraction (const G4String &modelName="HadronicModel")
 
virtual ~G4HadronicInteraction ()
 
virtual G4double SampleInvariantT (const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A)
 
virtual G4bool IsApplicable (const G4HadProjectile &, G4Nucleus &)
 
G4double GetMinEnergy () const
 
G4double GetMinEnergy (const G4Material *aMaterial, const G4Element *anElement) const
 
void SetMinEnergy (G4double anEnergy)
 
void SetMinEnergy (G4double anEnergy, const G4Element *anElement)
 
void SetMinEnergy (G4double anEnergy, const G4Material *aMaterial)
 
G4double GetMaxEnergy () const
 
G4double GetMaxEnergy (const G4Material *aMaterial, const G4Element *anElement) const
 
void SetMaxEnergy (const G4double anEnergy)
 
void SetMaxEnergy (G4double anEnergy, const G4Element *anElement)
 
void SetMaxEnergy (G4double anEnergy, const G4Material *aMaterial)
 
const G4HadronicInteractionGetMyPointer () const
 
virtual G4int GetVerboseLevel () const
 
virtual void SetVerboseLevel (G4int value)
 
const G4StringGetModelName () const
 
void DeActivateFor (const G4Material *aMaterial)
 
void ActivateFor (const G4Material *aMaterial)
 
void DeActivateFor (const G4Element *anElement)
 
void ActivateFor (const G4Element *anElement)
 
G4bool IsBlocked (const G4Material *aMaterial) const
 
G4bool IsBlocked (const G4Element *anElement) const
 
void SetRecoilEnergyThreshold (G4double val)
 
G4double GetRecoilEnergyThreshold () const
 
G4bool operator== (const G4HadronicInteraction &right) const
 
G4bool operator!= (const G4HadronicInteraction &right) const
 
virtual const std::pair
< G4double, G4double
GetFatalEnergyCheckLevels () const
 
virtual std::pair< G4double,
G4double
GetEnergyMomentumCheckLevels () const
 
void SetEnergyMomentumCheckLevels (G4double relativeLevel, G4double absoluteLevel)
 
virtual void ModelDescription (std::ostream &outFile) const
 

Additional Inherited Members

- Protected Types inherited from G4RPGInelastic
enum  {
  pi0, pip, pim, kp,
  km, k0, k0b, pro,
  neu, lam, sp, s0,
  sm, xi0, xim, om,
  ap, an
}
 
- Protected Member Functions inherited from G4RPGInelastic
G4double Pmltpc (G4int np, G4int nm, G4int nz, G4int n, G4double b, G4double c)
 
G4int Factorial (G4int n)
 
G4bool MarkLeadingStrangeParticle (const G4ReactionProduct &currentParticle, const G4ReactionProduct &targetParticle, G4ReactionProduct &leadParticle)
 
void SetUpPions (const G4int np, const G4int nm, const G4int nz, G4FastVector< G4ReactionProduct, 256 > &vec, G4int &vecLen)
 
void GetNormalizationConstant (const G4double availableEnergy, G4double &n, G4double &anpn)
 
void CalculateMomenta (G4FastVector< G4ReactionProduct, 256 > &vec, G4int &vecLen, const G4HadProjectile *originalIncident, const G4DynamicParticle *originalTarget, G4ReactionProduct &modifiedOriginal, G4Nucleus &targetNucleus, G4ReactionProduct &currentParticle, G4ReactionProduct &targetParticle, G4bool &incidentHasChanged, G4bool &targetHasChanged, G4bool quasiElastic)
 
void SetUpChange (G4FastVector< G4ReactionProduct, 256 > &vec, G4int &vecLen, G4ReactionProduct &currentParticle, G4ReactionProduct &targetParticle, G4bool &incidentHasChanged)
 
std::pair< G4int, G4doubleinterpolateEnergy (G4double ke) const
 
G4int sampleFlat (std::vector< G4double > sigma) const
 
void CheckQnums (G4FastVector< G4ReactionProduct, 256 > &vec, G4int &vecLen, G4ReactionProduct &currentParticle, G4ReactionProduct &targetParticle, G4double Q, G4double B, G4double S)
 
- Protected Member Functions inherited from G4HadronicInteraction
void SetModelName (const G4String &nam)
 
G4bool IsBlocked () const
 
void Block ()
 
- Protected Attributes inherited from G4RPGInelastic
G4RPGFragmentation fragmentation
 
G4RPGTwoCluster twoCluster
 
G4RPGPionSuppression pionSuppression
 
G4RPGStrangeProduction strangeProduction
 
G4RPGTwoBody twoBody
 
G4ParticleDefinitionparticleDef [18]
 
- Protected Attributes inherited from G4HadronicInteraction
G4HadFinalState theParticleChange
 
G4int verboseLevel
 
G4double theMinEnergy
 
G4double theMaxEnergy
 
G4bool isBlocked
 

Detailed Description

Definition at line 43 of file G4RPGKZeroInelastic.hh.

Constructor & Destructor Documentation

G4RPGKZeroInelastic::G4RPGKZeroInelastic ( )
inline

Definition at line 47 of file G4RPGKZeroInelastic.hh.

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

47  : G4RPGInelastic("G4RPGKZeroInelastic")
48  {
49  SetMinEnergy( 0.0 );
50  SetMaxEnergy( 25.*CLHEP::GeV );
51  }
void SetMinEnergy(G4double anEnergy)
G4RPGInelastic(const G4String &modelName="RPGInelastic")
void SetMaxEnergy(const G4double anEnergy)
G4RPGKZeroInelastic::~G4RPGKZeroInelastic ( )
inline

Definition at line 53 of file G4RPGKZeroInelastic.hh.

54  { }

Member Function Documentation

G4HadFinalState * G4RPGKZeroInelastic::ApplyYourself ( const G4HadProjectile aTrack,
G4Nucleus targetNucleus 
)
virtual

Implements G4HadronicInteraction.

Definition at line 35 of file G4RPGKZeroInelastic.cc.

References G4RPGInelastic::CalculateMomenta(), G4Nucleus::Cinema(), G4Nucleus::EvaporationEffects(), G4cout, G4endl, G4HadProjectile::GetDefinition(), G4DynamicParticle::GetDefinition(), G4HadProjectile::GetKineticEnergy(), G4ReactionProduct::GetKineticEnergy(), G4HadProjectile::GetMaterial(), G4ReactionProduct::GetMomentum(), G4Material::GetName(), G4ParticleDefinition::GetParticleName(), G4ParticleDefinition::GetPDGMass(), G4FastVector< Type, N >::Initialize(), CLHEP::Hep3Vector::mag(), python.hepunit::MeV, G4InuclParticleNames::pp, G4Nucleus::ReturnTargetParticle(), G4ReactionProduct::SetKineticEnergy(), G4ReactionProduct::SetMomentum(), G4ReactionProduct::SetSide(), G4RPGInelastic::SetUpChange(), G4HadronicInteraction::theParticleChange, and G4HadronicInteraction::verboseLevel.

Referenced by G4RPGKShortInelastic::ApplyYourself(), and G4RPGKLongInelastic::ApplyYourself().

37 {
38  const G4HadProjectile *originalIncident = &aTrack;
39 
40  // create the target particle
41 
42  G4DynamicParticle *originalTarget = targetNucleus.ReturnTargetParticle();
43 
44  if( verboseLevel > 1 )
45  {
46  const G4Material *targetMaterial = aTrack.GetMaterial();
47  G4cout << "G4RPGKZeroInelastic::ApplyYourself called" << G4endl;
48  G4cout << "kinetic energy = " << originalIncident->GetKineticEnergy()/MeV << "MeV, ";
49  G4cout << "target material = " << targetMaterial->GetName() << ", ";
50  G4cout << "target particle = " << originalTarget->GetDefinition()->GetParticleName()
51  << G4endl;
52  }
53 
54  // Fermi motion and evaporation
55  // As of Geant3, the Fermi energy calculation had not been Done
56 
57  G4double ek = originalIncident->GetKineticEnergy()/MeV;
58  G4double amas = originalIncident->GetDefinition()->GetPDGMass()/MeV;
59  G4ReactionProduct modifiedOriginal;
60  modifiedOriginal = *originalIncident;
61 
62  G4double tkin = targetNucleus.Cinema( ek );
63  ek += tkin;
64  modifiedOriginal.SetKineticEnergy( ek*MeV );
65  G4double et = ek + amas;
66  G4double p = std::sqrt( std::abs((et-amas)*(et+amas)) );
67  G4double pp = modifiedOriginal.GetMomentum().mag()/MeV;
68  if( pp > 0.0 )
69  {
70  G4ThreeVector momentum = modifiedOriginal.GetMomentum();
71  modifiedOriginal.SetMomentum( momentum * (p/pp) );
72  }
73  //
74  // calculate black track energies
75  //
76  tkin = targetNucleus.EvaporationEffects( ek );
77  ek -= tkin;
78  modifiedOriginal.SetKineticEnergy( ek*MeV );
79  et = ek + amas;
80  p = std::sqrt( std::abs((et-amas)*(et+amas)) );
81  pp = modifiedOriginal.GetMomentum().mag()/MeV;
82  if( pp > 0.0 )
83  {
84  G4ThreeVector momentum = modifiedOriginal.GetMomentum();
85  modifiedOriginal.SetMomentum( momentum * (p/pp) );
86  }
87  G4ReactionProduct currentParticle = modifiedOriginal;
88  G4ReactionProduct targetParticle;
89  targetParticle = *originalTarget;
90  currentParticle.SetSide( 1 ); // incident always goes in forward hemisphere
91  targetParticle.SetSide( -1 ); // target always goes in backward hemisphere
92  G4bool incidentHasChanged = false;
93  G4bool targetHasChanged = false;
94  G4bool quasiElastic = false;
95  G4FastVector<G4ReactionProduct,GHADLISTSIZE> vec; // vec will contain the secondary particles
96  G4int vecLen = 0;
97  vec.Initialize( 0 );
98 
99  const G4double cutOff = 0.1;
100  if( currentParticle.GetKineticEnergy()/MeV > cutOff )
101  Cascade( vec, vecLen,
102  originalIncident, currentParticle, targetParticle,
103  incidentHasChanged, targetHasChanged, quasiElastic );
104 
105  CalculateMomenta( vec, vecLen,
106  originalIncident, originalTarget, modifiedOriginal,
107  targetNucleus, currentParticle, targetParticle,
108  incidentHasChanged, targetHasChanged, quasiElastic );
109 
110  SetUpChange( vec, vecLen,
111  currentParticle, targetParticle,
112  incidentHasChanged );
113 
114  delete originalTarget;
115  return &theParticleChange;
116 }
G4double EvaporationEffects(G4double kineticEnergy)
Definition: G4Nucleus.cc:264
void SetUpChange(G4FastVector< G4ReactionProduct, 256 > &vec, G4int &vecLen, G4ReactionProduct &currentParticle, G4ReactionProduct &targetParticle, G4bool &incidentHasChanged)
void SetKineticEnergy(const G4double en)
void SetMomentum(const G4double x, const G4double y, const G4double z)
const char * p
Definition: xmltok.h:285
const G4String & GetName() const
Definition: G4Material.hh:176
void SetSide(const G4int sid)
void CalculateMomenta(G4FastVector< G4ReactionProduct, 256 > &vec, G4int &vecLen, const G4HadProjectile *originalIncident, const G4DynamicParticle *originalTarget, G4ReactionProduct &modifiedOriginal, G4Nucleus &targetNucleus, G4ReactionProduct &currentParticle, G4ReactionProduct &targetParticle, G4bool &incidentHasChanged, G4bool &targetHasChanged, G4bool quasiElastic)
G4ParticleDefinition * GetDefinition() const
void Initialize(G4int items)
Definition: G4FastVector.hh:63
int G4int
Definition: G4Types.hh:78
G4DynamicParticle * ReturnTargetParticle() const
Definition: G4Nucleus.cc:227
const G4String & GetParticleName() const
G4GLOB_DLL std::ostream G4cout
const G4ParticleDefinition * GetDefinition() const
bool G4bool
Definition: G4Types.hh:79
G4double GetKineticEnergy() const
G4double GetKineticEnergy() const
G4double GetPDGMass() const
G4double Cinema(G4double kineticEnergy)
Definition: G4Nucleus.cc:368
G4ThreeVector GetMomentum() const
#define G4endl
Definition: G4ios.hh:61
const G4Material * GetMaterial() const
double G4double
Definition: G4Types.hh:76
double mag() const

The documentation for this class was generated from the following files: