Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4LEPTSIonisationModel.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
28
29//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
31 : G4VLEPTSModel( modelName )
32{
34 fParticleChangeForGamma = 0;
36
37} // constructor
38
39//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
41{
42}
43
44
45//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
47 const G4DataVector&)
48{
49 Init();
50 BuildPhysicsTable( *aParticle );
51 fParticleChangeForGamma = GetParticleChangeForGamma();
52
53}
54
55
56//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
58 const G4ParticleDefinition* aParticle,
59 G4double kineticEnergy,
62{
63 return 1./GetMeanFreePath( mate, aParticle, kineticEnergy );
64
65}
66
67
68void G4LEPTSIonisationModel::SampleSecondaries(std::vector<G4DynamicParticle*>* fvect,
69 const G4MaterialCutsCouple* mateCuts,
70 const G4DynamicParticle* aDynamicParticle,
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
@ XSIonisation
#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
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
G4LEPTSIonisationModel(const G4String &modelName="G4LEPTSIonisationModel")
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double tmax=DBL_MAX)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
const G4Material * GetMaterial() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(const G4ThreeVector &Pfinal)
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:124
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:802
G4ThreeVector SampleNewDirection(const G4Material *aMaterial, G4ThreeVector Dir, G4double e, G4double el)
G4double SampleEnergyLoss(const G4Material *aMaterial, G4double eMin, G4double eMax)
G4double GetMeanFreePath(const G4Material *mate, const G4ParticleDefinition *aParticle, G4double kineticEnergy)
std::map< const G4Material *, G4double > theIonisPotInt
std::map< const G4Material *, G4double > theIonisPot
void BuildPhysicsTable(const G4ParticleDefinition &aParticleType)
G4double theLowestEnergyLimit
void ProposeLocalEnergyDeposit(G4double anEnergyPart)