Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4LEProtonInelastic Class Reference

#include <G4LEProtonInelastic.hh>

+ Inheritance diagram for G4LEProtonInelastic:

Public Member Functions

 G4LEProtonInelastic (const G4String &name="G4LEProtonInelastic")
 
 ~G4LEProtonInelastic ()
 
G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
 
virtual void ModelDescription (std::ostream &outFile) const
 
- Public Member Functions inherited from G4InelasticInteraction
 G4InelasticInteraction (const G4String &name="LEInelastic")
 
virtual ~G4InelasticInteraction ()
 
void RegisterIsotopeProductionModel (G4VIsotopeProduction *aModel)
 
void TurnOnIsotopeProduction ()
 
virtual const std::pair< G4double, G4doubleGetFatalEnergyCheckLevels () const
 
- Public Member Functions inherited from G4HadronicInteraction
 G4HadronicInteraction (const G4String &modelName="HadronicModel")
 
virtual ~G4HadronicInteraction ()
 
virtual G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)=0
 
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
 
G4int GetVerboseLevel () const
 
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, G4doubleGetFatalEnergyCheckLevels () const
 
virtual std::pair< G4double, G4doubleGetEnergyMomentumCheckLevels () const
 
void SetEnergyMomentumCheckLevels (G4double relativeLevel, G4double absoluteLevel)
 
virtual void ModelDescription (std::ostream &outFile) const
 

Additional Inherited Members

- Static Public Member Functions inherited from G4InelasticInteraction
static G4IsoParticleChangeGetIsotopeProductionInfo ()
 
- Protected Member Functions inherited from G4InelasticInteraction
G4double Pmltpc (G4int np, G4int nm, G4int nz, G4int n, G4double b, G4double c)
 
G4bool MarkLeadingStrangeParticle (const G4ReactionProduct &currentParticle, const G4ReactionProduct &targetParticle, G4ReactionProduct &leadParticle)
 
void SetUpPions (const G4int np, const G4int nm, const G4int nz, G4FastVector< G4ReactionProduct, GHADLISTSIZE > &vec, G4int &vecLen)
 
void Rotate (G4FastVector< G4ReactionProduct, GHADLISTSIZE > &vec, G4int &vecLen)
 
void GetNormalizationConstant (const G4double availableEnergy, G4double &n, G4double &anpn)
 
void CalculateMomenta (G4FastVector< G4ReactionProduct, GHADLISTSIZE > &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, GHADLISTSIZE > &vec, G4int &vecLen, G4ReactionProduct &currentParticle, G4ReactionProduct &targetParticle, G4bool &incidentHasChanged)
 
void DoIsotopeCounting (const G4HadProjectile *theProjectile, const G4Nucleus &aNucleus)
 
G4IsoResultExtractResidualNucleus (const G4Nucleus &aNucleus)
 
- Protected Member Functions inherited from G4HadronicInteraction
void SetModelName (const G4String &nam)
 
G4bool IsBlocked () const
 
void Block ()
 
- Protected Attributes inherited from G4InelasticInteraction
G4bool isotopeProduction
 
G4ReactionDynamics theReactionDynamics
 
- Protected Attributes inherited from G4HadronicInteraction
G4HadFinalState theParticleChange
 
G4int verboseLevel
 
G4double theMinEnergy
 
G4double theMaxEnergy
 
G4bool isBlocked
 

Detailed Description

Definition at line 46 of file G4LEProtonInelastic.hh.

Constructor & Destructor Documentation

◆ G4LEProtonInelastic()

G4LEProtonInelastic::G4LEProtonInelastic ( const G4String name = "G4LEProtonInelastic")

Definition at line 36 of file G4LEProtonInelastic.cc.

38{
39 SetMinEnergy(0.0);
40 SetMaxEnergy(55.*GeV);
41 G4cout << "WARNING: model G4LEProtonInelastic is being deprecated and will\n"
42 << "disappear in Geant4 version 10.0" << G4endl;
43}
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
void SetMinEnergy(G4double anEnergy)
void SetMaxEnergy(const G4double anEnergy)

◆ ~G4LEProtonInelastic()

G4LEProtonInelastic::~G4LEProtonInelastic ( )
inline

Definition at line 52 of file G4LEProtonInelastic.hh.

52{}

Member Function Documentation

◆ ApplyYourself()

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

Implements G4HadronicInteraction.

Definition at line 61 of file G4LEProtonInelastic.cc.

63{
65 const G4HadProjectile *originalIncident = &aTrack;
66 if (originalIncident->GetKineticEnergy()<= 0.1*MeV) {
70 return &theParticleChange;
71 }
72
73 // Create the target particle
74
75 G4DynamicParticle *originalTarget = targetNucleus.ReturnTargetParticle();
76 if (verboseLevel > 1) {
77 const G4Material *targetMaterial = aTrack.GetMaterial();
78 G4cout << "G4LEProtonInelastic::ApplyYourself called" << G4endl;
79 G4cout << "kinetic energy = " << originalIncident->GetKineticEnergy()/MeV << "MeV, ";
80 G4cout << "target material = " << targetMaterial->GetName() << ", ";
81 G4cout << "target particle = " << originalTarget->GetDefinition()->GetParticleName()
82 << G4endl;
83 }
84
85 if (originalIncident->GetKineticEnergy()/GeV < 0.01+2.*G4UniformRand()/9.) {
86 SlowProton(originalIncident, targetNucleus);
87 if (isotopeProduction) DoIsotopeCounting(originalIncident, targetNucleus);
88 delete originalTarget;
89 return &theParticleChange;
90 }
91
92 // Fermi motion and evaporation
93 // As of Geant3, the Fermi energy calculation had not been Done
94
95 G4double ek = originalIncident->GetKineticEnergy()/MeV;
96 G4double amas = originalIncident->GetDefinition()->GetPDGMass()/MeV;
97 G4ReactionProduct modifiedOriginal;
98 modifiedOriginal = *originalIncident;
99
100 G4double tkin = targetNucleus.Cinema( ek );
101 ek += tkin;
102 modifiedOriginal.SetKineticEnergy( ek*MeV );
103 G4double et = ek + amas;
104 G4double p = std::sqrt( std::abs((et-amas)*(et+amas)) );
105 G4double pp = modifiedOriginal.GetMomentum().mag()/MeV;
106 if (pp > 0.0) {
107 G4ThreeVector momentum = modifiedOriginal.GetMomentum();
108 modifiedOriginal.SetMomentum( momentum * (p/pp) );
109 }
110
111 // calculate black track energies
112
113 tkin = targetNucleus.EvaporationEffects( ek );
114 ek -= tkin;
115 modifiedOriginal.SetKineticEnergy( ek*MeV );
116 et = ek + amas;
117 p = std::sqrt( std::abs((et-amas)*(et+amas)) );
118 pp = modifiedOriginal.GetMomentum().mag()/MeV;
119 if (pp > 0.0) {
120 G4ThreeVector momentum = modifiedOriginal.GetMomentum();
121 modifiedOriginal.SetMomentum( momentum * (p/pp) );
122 }
123 const G4double cutOff = 0.1;
124 if (modifiedOriginal.GetKineticEnergy()/MeV <= cutOff) {
125 SlowProton( originalIncident, targetNucleus);
126 if (isotopeProduction) DoIsotopeCounting(originalIncident, targetNucleus);
127 delete originalTarget;
128 return &theParticleChange;
129 }
130 G4ReactionProduct currentParticle = modifiedOriginal;
131 G4ReactionProduct targetParticle;
132 targetParticle = *originalTarget;
133 currentParticle.SetSide( 1 ); // incident always goes in forward hemisphere
134 targetParticle.SetSide( -1 ); // target always goes in backward hemisphere
135 G4bool incidentHasChanged = false;
136 G4bool targetHasChanged = false;
137 G4bool quasiElastic = false;
138 G4FastVector<G4ReactionProduct,GHADLISTSIZE> vec; // vec will contain the sec. particles
139 G4int vecLen = 0;
140 vec.Initialize( 0 );
141
142 Cascade(vec, vecLen,
143 originalIncident, currentParticle, targetParticle,
144 incidentHasChanged, targetHasChanged, quasiElastic);
145
146 CalculateMomenta(vec, vecLen,
147 originalIncident, originalTarget, modifiedOriginal,
148 targetNucleus, currentParticle, targetParticle,
149 incidentHasChanged, targetHasChanged, quasiElastic);
150
151 SetUpChange(vec, vecLen,
152 currentParticle, targetParticle,
153 incidentHasChanged);
154
155 if (isotopeProduction) DoIsotopeCounting(originalIncident, targetNucleus);
156 delete originalTarget;
157 return &theParticleChange;
158}
@ isAlive
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
#define G4UniformRand()
Definition: Randomize.hh:53
Hep3Vector unit() const
double mag() const
Hep3Vector vect() const
G4ParticleDefinition * GetDefinition() const
void Initialize(G4int items)
Definition: G4FastVector.hh:63
void SetStatusChange(G4HadFinalStateStatus aS)
void SetEnergyChange(G4double anEnergy)
void SetMomentumChange(const G4ThreeVector &aV)
const G4Material * GetMaterial() const
const G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
const G4LorentzVector & Get4Momentum() const
void CalculateMomenta(G4FastVector< G4ReactionProduct, GHADLISTSIZE > &vec, G4int &vecLen, const G4HadProjectile *originalIncident, const G4DynamicParticle *originalTarget, G4ReactionProduct &modifiedOriginal, G4Nucleus &targetNucleus, G4ReactionProduct &currentParticle, G4ReactionProduct &targetParticle, G4bool &incidentHasChanged, G4bool &targetHasChanged, G4bool quasiElastic)
void DoIsotopeCounting(const G4HadProjectile *theProjectile, const G4Nucleus &aNucleus)
void SetUpChange(G4FastVector< G4ReactionProduct, GHADLISTSIZE > &vec, G4int &vecLen, G4ReactionProduct &currentParticle, G4ReactionProduct &targetParticle, G4bool &incidentHasChanged)
const G4String & GetName() const
Definition: G4Material.hh:177
G4double EvaporationEffects(G4double kineticEnergy)
Definition: G4Nucleus.cc:264
G4double Cinema(G4double kineticEnergy)
Definition: G4Nucleus.cc:368
G4DynamicParticle * ReturnTargetParticle() const
Definition: G4Nucleus.cc:227
const G4String & GetParticleName() const
void SetMomentum(const G4double x, const G4double y, const G4double z)
G4double GetKineticEnergy() const
G4ThreeVector GetMomentum() const
void SetSide(const G4int sid)
void SetKineticEnergy(const G4double en)

◆ ModelDescription()

void G4LEProtonInelastic::ModelDescription ( std::ostream &  outFile) const
virtual

Reimplemented from G4HadronicInteraction.

Definition at line 46 of file G4LEProtonInelastic.cc.

47{
48 outFile << "G4LEProtonInelastic is one of the Low Energy Parameterized\n"
49 << "(LEP) models used to implement inelastic proton scattering\n"
50 << "from nuclei. It is a re-engineered version of the GHEISHA\n"
51 << "code of H. Fesefeldt. It divides the initial collision\n"
52 << "products into backward- and forward-going clusters which are\n"
53 << "then decayed into final state hadrons. The model does not\n"
54 << "conserve energy on an event-by-event basis. It may be\n"
55 << "applied to protons with initial energies between 0 and 25\n"
56 << "GeV.\n";
57}

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