Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4LEPTSExcitationModel.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//
27
28
29//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
31 : G4VLEPTSModel( modelName )
32{
34} // constructor
35
36
37//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
39}
40
41
42//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
44 const G4DataVector&)
45{
46 Init();
47 BuildPhysicsTable( *aParticle );
48
49 fParticleChangeForGamma = GetParticleChangeForGamma();
50
51 LowestExcitationEnergy = 0;
52 LowestNeutralDisociationEnergy = 0;
53
54
55}
56
57
58std::map<G4int,std::vector<G4double> > G4LEPTSExcitationModel::ReadIXS(G4String fileTXS, const G4Material* aMaterial)
59{
60 std::map<G4int,std::vector<G4double> > integralXS = G4VLEPTSModel::ReadIXS( fileTXS, aMaterial);
61
62 if( integralXS.size() == 0 ) return integralXS;
63
64 for (G4int jj=theNXSdat[aMaterial]; jj>=0; jj--) {
65 if( integralXS[XSDissociation][jj] > 0.001) LowestExcitationEnergy = integralXS[XSTotal][jj-1];
66 if( integralXS[XSVibration][jj] > 0.001) LowestNeutralDisociationEnergy = integralXS[XSTotal][jj-1]*CLHEP::eV;
67 // if( txs[5][j] > 0.001) LowestExcitationEnergy = txs[0][j-1];
68 // if( txs[6][j] > 0.001) LowestNeutralDisociationEnergy = txs[0][j-1]*CLHEP::eV;
69 }
70
71 if( verboseLevel >= 1) G4cout << " LowestExcitationEnergy: " << LowestExcitationEnergy << G4endl
72 << "LowestNeutralDisociationEnergy: " << LowestNeutralDisociationEnergy/CLHEP::eV
73 << G4endl;
74
75 return integralXS;
76}
77
78//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
80 const G4ParticleDefinition* aParticle,
81 G4double kineticEnergy,
84{
85 return 1./GetMeanFreePath( mate, aParticle, kineticEnergy );
86
87}
88
89
90//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
91void G4LEPTSExcitationModel::SampleSecondaries(std::vector<G4DynamicParticle*>*,
92 const G4MaterialCutsCouple* mateCuts,
93 const G4DynamicParticle* aDynamicParticle,
96{
97 G4double P0KinEn = aDynamicParticle->GetKineticEnergy();
98
99 G4double Edep=0;
100 G4double Energylost=0;
101 G4ThreeVector P0Dir = aDynamicParticle->GetMomentumDirection();
102
103 G4double eMin = 0.0;
104 const G4Material* aMaterial = mateCuts->GetMaterial();
105 G4double eMax = std::min(theIonisPot[aMaterial], P0KinEn);
106 Energylost = SampleEnergyLoss(aMaterial, eMin, eMax);
107
108 Edep = Energylost;
109
110 G4ThreeVector P1Dir = SampleNewDirection(aMaterial, P0Dir, P0KinEn/CLHEP::eV, Energylost/CLHEP::eV);
111 G4double P1KinEn = P0KinEn - Edep;
112
113 fParticleChangeForGamma->ProposeMomentumDirection( P1Dir);
114 fParticleChangeForGamma->SetProposedKineticEnergy( P1KinEn);
115 fParticleChangeForGamma->ProposeLocalEnergyDeposit( Edep);
116
117}
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
@ XSExcitation
@ XSDissociation
@ XSTotal
@ XSVibration
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
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 std::map< G4int, std::vector< G4double > > ReadIXS(G4String, const G4Material *aMaterial)
G4LEPTSExcitationModel(const G4String &modelName="G4LEPTSExcitationModel")
const G4Material * GetMaterial() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(const G4ThreeVector &Pfinal)
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:124
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 > theIonisPot
void BuildPhysicsTable(const G4ParticleDefinition &aParticleType)
virtual std::map< G4int, std::vector< G4double > > ReadIXS(G4String fileName, const G4Material *aMaterial)
std::map< const G4Material *, G4int > theNXSdat
void ProposeLocalEnergyDeposit(G4double anEnergyPart)