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

#include <G4LEPTSIonisationModel.hh>

+ Inheritance diagram for G4LEPTSIonisationModel:

Public Member Functions

 G4LEPTSIonisationModel (const G4String &modelName="G4LEPTSIonisationModel")
 
 ~G4LEPTSIonisationModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double tmax=DBL_MAX)
 
virtual G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
- Public Member Functions inherited from G4VLEPTSModel
 G4VLEPTSModel (const G4String &processName)
 
 ~G4VLEPTSModel ()
 
void BuildPhysicsTable (const G4ParticleDefinition &aParticleType)
 
G4double GetMeanFreePath (const G4Material *mate, const G4ParticleDefinition *aParticle, G4double kineticEnergy)
 
G4ThreeVector SampleNewDirection (const G4Material *aMaterial, G4ThreeVector Dir, G4double e, G4double el)
 
G4double SampleAngle (const G4Material *aMaterial, G4double e, G4double el)
 
G4ThreeVector SampleNewDirection (G4ThreeVector Dir, G4double ang)
 
G4VLEPTSModeloperator= (const G4VLEPTSModel &right)
 
 G4VLEPTSModel (const G4VLEPTSModel &)
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)=0
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double tmax=DBL_MAX)=0
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel)
 
virtual void InitialiseForMaterial (const G4ParticleDefinition *, const G4Material *)
 
virtual void InitialiseForElement (const G4ParticleDefinition *, G4int Z)
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
virtual G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double GetPartialCrossSection (const G4Material *, G4int level, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ComputeCrossSectionPerShell (const G4ParticleDefinition *, G4int Z, G4int shellIdx, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ChargeSquareRatio (const G4Track &)
 
virtual G4double GetChargeSquareRatio (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual G4double GetParticleCharge (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void StartTracking (G4Track *)
 
virtual void CorrectionsAlongStep (const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double &eloss, G4double &niel, G4double length)
 
virtual G4double Value (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *, G4double cut=0.0)
 
virtual G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
virtual void ModelDescription (std::ostream &outFile) const
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
std::vector< G4EmElementSelector * > * GetElementSelectors ()
 
void SetElementSelectors (std::vector< G4EmElementSelector * > *)
 
virtual G4double ComputeDEDX (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
G4double CrossSection (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeMeanFreePath (const G4ParticleDefinition *, G4double kineticEnergy, const G4Material *, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, const G4Element *, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectTargetAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double logKineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4int SelectRandomAtomNumber (const G4Material *)
 
G4int SelectIsotopeNumber (const G4Element *)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
G4VEmModelGetTripletModel ()
 
void SetTripletModel (G4VEmModel *)
 
void SetAngularDistribution (G4VEmAngularDistribution *)
 
G4double HighEnergyLimit () const
 
G4double LowEnergyLimit () const
 
G4double HighEnergyActivationLimit () const
 
G4double LowEnergyActivationLimit () const
 
G4double PolarAngleLimit () const
 
G4double SecondaryThreshold () const
 
G4bool LPMFlag () const
 
G4bool DeexcitationFlag () const
 
G4bool ForceBuildTableFlag () const
 
G4bool UseAngularGeneratorFlag () const
 
void SetAngularGeneratorFlag (G4bool)
 
void SetHighEnergyLimit (G4double)
 
void SetLowEnergyLimit (G4double)
 
void SetActivationHighEnergyLimit (G4double)
 
void SetActivationLowEnergyLimit (G4double)
 
G4bool IsActive (G4double kinEnergy) const
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetLPMFlag (G4bool val)
 
void SetDeexcitationFlag (G4bool val)
 
void SetForceBuildTable (G4bool val)
 
void SetFluctuationFlag (G4bool val)
 
void SetMasterThread (G4bool val)
 
G4bool IsMaster () const
 
void SetUseBaseMaterials (G4bool val)
 
G4bool UseBaseMaterials () const
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
const G4ElementGetCurrentElement () const
 
const G4IsotopeGetCurrentIsotope () const
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

Additional Inherited Members

- Protected Member Functions inherited from G4VLEPTSModel
void Init ()
 
G4bool ReadParam (G4String fileName, const G4Material *aMaterial)
 
virtual std::map< G4int, std::vector< G4double > > ReadIXS (G4String fileName, const G4Material *aMaterial)
 
G4double SampleEnergyLoss (const G4Material *aMaterial, G4double eMin, G4double eMax)
 
void BuildMeanFreePathTable (const G4Material *aMaterial, std::map< G4int, std::vector< G4double > > &integralXS)
 
- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 
- Protected Attributes inherited from G4VLEPTSModel
G4PhysicsTabletheMeanFreePathTable
 
G4double theLowestEnergyLimit
 
G4double theHighestEnergyLimit
 
G4int theNumbBinTable
 
std::map< const G4Material *, G4doubletheIonisPot
 
std::map< const G4Material *, G4doubletheIonisPotInt
 
std::map< const G4Material *, G4doubletheMolecularMass
 
std::map< const G4Material *, G4LEPTSDiffXS * > theDiffXS
 
std::map< const G4Material *, G4LEPTSDistribution * > theRMTDistr
 
std::map< const G4Material *, G4LEPTSElossDistr * > theElostDistr
 
std::map< const G4Material *, G4LEPTSDistribution * > theElostDistr2
 
std::map< const G4Material *, G4inttheNXSdat
 
std::map< const G4Material *, G4inttheNXSsub
 
G4bool isInitialised
 
XSType theXSType
 
G4int verboseLevel
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const G4MaterialpBaseMaterial
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 
G4bool lossFlucFlag
 
G4double inveplus
 
G4double pFactor
 

Detailed Description

Definition at line 31 of file G4LEPTSIonisationModel.hh.

Constructor & Destructor Documentation

◆ G4LEPTSIonisationModel()

G4LEPTSIonisationModel::G4LEPTSIonisationModel ( const G4String modelName = "G4LEPTSIonisationModel")

Definition at line 30 of file G4LEPTSIonisationModel.cc.

31 : G4VLEPTSModel( modelName )
32{
34 fParticleChangeForGamma = 0;
36
37} // constructor
@ XSIonisation
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:813

◆ ~G4LEPTSIonisationModel()

G4LEPTSIonisationModel::~G4LEPTSIonisationModel ( )

Definition at line 40 of file G4LEPTSIonisationModel.cc.

41{
42}

Member Function Documentation

◆ CrossSectionPerVolume()

G4double G4LEPTSIonisationModel::CrossSectionPerVolume ( const G4Material mate,
const G4ParticleDefinition aParticle,
G4double  kineticEnergy,
G4double  cutEnergy = 0.0,
G4double  maxEnergy = DBL_MAX 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 57 of file G4LEPTSIonisationModel.cc.

62{
63 return 1./GetMeanFreePath( mate, aParticle, kineticEnergy );
64
65}
G4double GetMeanFreePath(const G4Material *mate, const G4ParticleDefinition *aParticle, G4double kineticEnergy)

◆ Initialise()

void G4LEPTSIonisationModel::Initialise ( const G4ParticleDefinition aParticle,
const G4DataVector  
)
virtual

Implements G4VEmModel.

Definition at line 46 of file G4LEPTSIonisationModel.cc.

48{
49 Init();
50 BuildPhysicsTable( *aParticle );
51 fParticleChangeForGamma = GetParticleChangeForGamma();
52
53}
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:133
void BuildPhysicsTable(const G4ParticleDefinition &aParticleType)

◆ SampleSecondaries()

void G4LEPTSIonisationModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  fvect,
const G4MaterialCutsCouple mateCuts,
const G4DynamicParticle aDynamicParticle,
G4double  tmin = 0.0,
G4double  tmax = DBL_MAX 
)
virtual

Implements G4VEmModel.

Definition at line 68 of file G4LEPTSIonisationModel.cc.

73{
74 G4double P0KinEn = aDynamicParticle->GetKineticEnergy();
75
76 G4double Edep=0;
77 G4double Energylost=0;
78 G4ThreeVector P0Dir = aDynamicParticle->GetMomentumDirection();
79
80 const G4Material* aMaterial = mateCuts->GetMaterial();
81 if(P0KinEn < theIonisPot[aMaterial]) {
82 theIonisPot[aMaterial] = P0KinEn;
83 }
84 Energylost = SampleEnergyLoss(aMaterial, theIonisPot[aMaterial], P0KinEn);
85 G4ThreeVector P1Dir = SampleNewDirection(aMaterial, P0Dir, P0KinEn/CLHEP::eV, Energylost/CLHEP::eV);
86 G4double P1KinEn = std::max(0., P0KinEn - Energylost);
87 fParticleChangeForGamma->ProposeMomentumDirection( P1Dir);
88 fParticleChangeForGamma->SetProposedKineticEnergy( P1KinEn);
89#ifdef DEBUG_LEPTS
90 G4cout << " G4LEPTSIonisationModel::SampleSecondaries( SetProposedKineticEnergy " << P1KinEn << " " << P0KinEn << " - " << Energylost << G4endl;
91#endif
92
93 G4double P2KinEn;
94
95 if( Energylost < theIonisPotInt[aMaterial]) { // External Ionisation
96 //- SetModelName("Ionisation");
97 Edep = theIonisPot[aMaterial];
98 P2KinEn = std::max(0.001*CLHEP::eV, (Energylost - theIonisPot[aMaterial]) );
99 }
100 else { // Auger
101 //- SetModelName("IonisAuger");
102 Edep = 35*CLHEP::eV;
103 P2KinEn = std::max(0.0, (Energylost - theIonisPotInt[aMaterial]) );
104 G4double P3KinEn = std::max(0.0, theIonisPotInt[aMaterial] - Edep);
105
106 G4ThreeVector P3Dir;
107 P3Dir.setX( G4UniformRand() );
108 P3Dir.setY( G4UniformRand() );
109 P3Dir.setZ( G4UniformRand() );
110 P3Dir /= P3Dir.mag();
111
112 G4DynamicParticle* e3 = new G4DynamicParticle(G4Electron::Electron(), P3Dir, P3KinEn);
113 fvect->push_back(e3);
114 }
115
116 fParticleChangeForGamma->ProposeLocalEnergyDeposit(Edep);
117
118 if( P2KinEn > theLowestEnergyLimit) {
119 G4double cp0 = std::sqrt(P0KinEn*(P0KinEn + 2.*CLHEP::electron_mass_c2));
120 G4double cp1 = std::sqrt(P1KinEn*(P1KinEn + 2.*CLHEP::electron_mass_c2));
121 G4ThreeVector P2Momentum = cp0*P0Dir -cp1*P1Dir;
122 G4ThreeVector P2Dir = P2Momentum / P2Momentum.mag();
123 P2Dir.rotateUz(P0Dir);
124 G4DynamicParticle* e2 = new G4DynamicParticle(G4Electron::Electron(), P2Dir, P2KinEn);
125 fvect->push_back(e2);
126 }
127
128}
double G4double
Definition: G4Types.hh:83
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
void setY(double)
void setZ(double)
double mag() const
void setX(double)
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:33
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
static G4Electron * Electron()
Definition: G4Electron.cc:93
const G4Material * GetMaterial() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
G4ThreeVector SampleNewDirection(const G4Material *aMaterial, G4ThreeVector Dir, G4double e, G4double el)
G4double SampleEnergyLoss(const G4Material *aMaterial, G4double eMin, G4double eMax)
std::map< const G4Material *, G4double > theIonisPotInt
std::map< const G4Material *, G4double > theIonisPot
G4double theLowestEnergyLimit
void ProposeLocalEnergyDeposit(G4double anEnergyPart)

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