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

#include <G4LEPionPlusInelastic.hh>

+ Inheritance diagram for G4LEPionPlusInelastic:

Public Member Functions

 G4LEPionPlusInelastic (const G4String &name="G4LEPionPlusInelastic")
 
 ~G4LEPionPlusInelastic ()
 
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 43 of file G4LEPionPlusInelastic.hh.

Constructor & Destructor Documentation

◆ G4LEPionPlusInelastic()

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

Definition at line 41 of file G4LEPionPlusInelastic.cc.

43{
44 SetMinEnergy(0.0);
45 SetMaxEnergy(55.*GeV);
46 G4cout << "WARNING: model G4LEPionPlusInelastic is being deprecated and will\n"
47 << "disappear in Geant4 version 10.0" << G4endl;
48}
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
void SetMinEnergy(G4double anEnergy)
void SetMaxEnergy(const G4double anEnergy)

◆ ~G4LEPionPlusInelastic()

G4LEPionPlusInelastic::~G4LEPionPlusInelastic ( )
inline

Definition at line 49 of file G4LEPionPlusInelastic.hh.

49{ }

Member Function Documentation

◆ ApplyYourself()

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

Implements G4HadronicInteraction.

Definition at line 66 of file G4LEPionPlusInelastic.cc.

68{
69 const G4HadProjectile *originalIncident = &aTrack;
70 if (originalIncident->GetKineticEnergy()<= 0.1*MeV) {
74 return &theParticleChange;
75 }
76
77 // create the target particle
78
79 G4DynamicParticle *originalTarget = targetNucleus.ReturnTargetParticle();
80 G4ReactionProduct targetParticle( originalTarget->GetDefinition() );
81
82 if (verboseLevel > 1) {
83 const G4Material* targetMaterial = aTrack.GetMaterial();
84 G4cout << "G4LEPionPlusInelastic::ApplyYourself called" << G4endl;
85 G4cout << "kinetic energy = " << originalIncident->GetKineticEnergy() << "MeV, ";
86 G4cout << "target material = " << targetMaterial->GetName() << ", ";
87 G4cout << "target particle = " << originalTarget->GetDefinition()->GetParticleName()
88 << G4endl;
89 }
90 G4ReactionProduct currentParticle(
91 const_cast<G4ParticleDefinition *>(originalIncident->GetDefinition() ) );
92 currentParticle.SetMomentum( originalIncident->Get4Momentum().vect() );
93 currentParticle.SetKineticEnergy( originalIncident->GetKineticEnergy() );
94
95 // Fermi motion and evaporation
96 // As of Geant3, the Fermi energy calculation had not been Done
97 G4double ek = originalIncident->GetKineticEnergy();
98 G4double amas = originalIncident->GetDefinition()->GetPDGMass();
99
100 G4double tkin = targetNucleus.Cinema(ek);
101 ek += tkin;
102 currentParticle.SetKineticEnergy(ek);
103 G4double et = ek + amas;
104 G4double p = std::sqrt( std::abs((et-amas)*(et+amas)) );
105 G4double pp = currentParticle.GetMomentum().mag();
106 if (pp > 0.0) {
107 G4ThreeVector momentum = currentParticle.GetMomentum();
108 currentParticle.SetMomentum( momentum * (p/pp) );
109 }
110
111 // calculate black track energies
112 tkin = targetNucleus.EvaporationEffects(ek);
113 ek -= tkin;
114 currentParticle.SetKineticEnergy(ek);
115 et = ek + amas;
116 p = std::sqrt( std::abs((et-amas)*(et+amas)) );
117 pp = currentParticle.GetMomentum().mag();
118 if (pp > 0.0) {
119 G4ThreeVector momentum = currentParticle.GetMomentum();
120 currentParticle.SetMomentum( momentum * (p/pp) );
121 }
122
123 G4ReactionProduct modifiedOriginal = currentParticle;
124
125 currentParticle.SetSide(1); // incident always goes in forward hemisphere
126 targetParticle.SetSide(-1); // target always goes in backward hemisphere
127 G4bool incidentHasChanged = false;
128 G4bool targetHasChanged = false;
129 G4bool quasiElastic = false;
130 G4FastVector<G4ReactionProduct,GHADLISTSIZE> vec; // vec will contain the secondary particles
131 G4int vecLen = 0;
132 vec.Initialize(0);
133
134 const G4double cutOff = 0.1*MeV;
135 if (currentParticle.GetKineticEnergy() > cutOff)
136 Cascade(vec, vecLen, originalIncident, currentParticle, targetParticle,
137 incidentHasChanged, targetHasChanged, quasiElastic);
138
139 CalculateMomenta(vec, vecLen,
140 originalIncident, originalTarget, modifiedOriginal,
141 targetNucleus, currentParticle, targetParticle,
142 incidentHasChanged, targetHasChanged, quasiElastic);
143
144 SetUpChange(vec, vecLen, currentParticle, targetParticle, incidentHasChanged);
145
146 if (isotopeProduction) DoIsotopeCounting(originalIncident, targetNucleus);
147
148 delete originalTarget;
149 return &theParticleChange;
150}
@ isAlive
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
Hep3Vector unit() 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 SetSide(const G4int sid)

◆ ModelDescription()

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

Reimplemented from G4HadronicInteraction.

Definition at line 51 of file G4LEPionPlusInelastic.cc.

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

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