Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PAIPhotModel.hh
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//
26//
27// -------------------------------------------------------------------
28//
29// GEANT4 Class header file
30//
31//
32// File name: G4PAIPhotModel
33//
34// Author: V. Grichine based on G4PAIModel code for MT
35//
36// Creation date: 07.10.2013
37//
38//
39// Class Description:
40//
41// Implementation of PAI photon model of energy loss and
42// delta-electron or X-ray production by charged particles in MT framework
43//
44// -------------------------------------------------------------------
45//
46
47#ifndef G4PAIPhotModel_h
48#define G4PAIPhotModel_h 1
49
51
52#include "G4VEmModel.hh"
54#include "globals.hh"
55#include <vector>
56
57class G4Region;
60
61class G4PAIPhotData;
62
64{
65
66public:
67
68 explicit G4PAIPhotModel(const G4ParticleDefinition* p = nullptr,
69 const G4String& nam = "PAI");
70
71 ~G4PAIPhotModel() final;
72
73 void Initialise(const G4ParticleDefinition*, const G4DataVector&) final;
74
76 G4VEmModel* masterModel) final;
77
79 const G4MaterialCutsCouple* couple) final;
80
83 G4double kineticEnergy,
84 G4double cutEnergy) final;
85
88 G4double kineticEnergy,
89 G4double cutEnergy,
90 G4double maxEnergy) final;
91
92 void SampleSecondaries(std::vector<G4DynamicParticle*>*,
94 const G4DynamicParticle*,
95 G4double tmin,
96 G4double maxEnergy) final;
97
99 const G4DynamicParticle*,
100 G4double, G4double, G4double) final;
101
103 G4double, G4double) final;
104
105 void DefineForRegion(const G4Region* r) final;
106
108
109 inline const std::vector<const G4MaterialCutsCouple*>& GetVectorOfCouples();
110
111 inline G4double ComputeMaxEnergy(G4double scaledEnergy);
112
113 inline void SetVerboseLevel(G4int verbose);
114
115protected:
116
118 G4double kinEnergy) final;
119
120private:
121
122 inline G4int FindCoupleIndex(const G4MaterialCutsCouple*);
123
124 inline void SetParticle(const G4ParticleDefinition* p);
125
126 // hide assignment operator
127 G4PAIPhotModel & operator=(const G4PAIPhotModel &right) = delete;
128 G4PAIPhotModel(const G4PAIPhotModel&) = delete;
129
130 G4int fVerbose;
131
132 G4PAIPhotData* fModelData;
133
134 std::vector<const G4MaterialCutsCouple*> fMaterialCutsCoupleVector;
135 std::vector<const G4Region*> fPAIRegionVector;
136
137 const G4ParticleDefinition* fParticle;
138 const G4ParticleDefinition* fElectron;
139 const G4ParticleDefinition* fPositron;
140 G4ParticleChangeForLoss* fParticleChange;
141
142 G4double fMass;
143 G4double fRatio;
144 G4double fChargeSquare;
145 G4double fLowestTcut;
146};
147
149{
150 return fModelData;
151}
152
153inline const std::vector<const G4MaterialCutsCouple*>&
155{
156 return fMaterialCutsCoupleVector;
157}
158
160{
161 return MaxSecondaryEnergy(fParticle, scaledEnergy/fRatio);
162}
163
165{
166 fVerbose=verbose;
167}
168
169inline G4int G4PAIPhotModel::FindCoupleIndex(const G4MaterialCutsCouple* couple)
170{
171 G4int idx = -1;
172 size_t jMatMax = fMaterialCutsCoupleVector.size();
173 for(size_t jMat = 0;jMat < jMatMax; ++jMat) {
174 if(couple == fMaterialCutsCoupleVector[jMat]) {
175 idx = jMat;
176 break;
177 }
178 }
179 return idx;
180}
181
182inline void G4PAIPhotModel::SetParticle(const G4ParticleDefinition* p)
183{
184 if(fParticle != p) {
185 fParticle = p;
186 fMass = fParticle->GetPDGMass();
187 fRatio = CLHEP::proton_mass_c2/fMass;
188 G4double q = fParticle->GetPDGCharge()/CLHEP::eplus;
189 fChargeSquare = q*q;
190 }
191}
192
193#endif
194
195
196
197
198
199
200
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
~G4PAIPhotModel() final
G4double MinEnergyCut(const G4ParticleDefinition *, const G4MaterialCutsCouple *couple) final
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) final
G4double Dispersion(const G4Material *, const G4DynamicParticle *, G4double, G4double) final
G4PAIPhotData * GetPAIPhotData()
const std::vector< const G4MaterialCutsCouple * > & GetVectorOfCouples()
void DefineForRegion(const G4Region *r) final
void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel) final
G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy) final
void Initialise(const G4ParticleDefinition *, const G4DataVector &) final
G4double ComputeMaxEnergy(G4double scaledEnergy)
G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy) final
G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) final
G4double SampleFluctuations(const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double, G4double, G4double) final
void SetVerboseLevel(G4int verbose)
G4double GetPDGCharge() const