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
G4PAIModel.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: G4PAIModel
33//
34// Author: V. Grichine based on Vladimir Ivanchenko code
35//
36// Creation date: 05.10.2003
37//
38// Modifications:
39// 08-04-05 Major optimisation of internal interfaces (V.Ivantchenko)
40// 26-09-07 Fixed tmax computation (V.Ivantchenko)
41// 19.08.13 V.Ivanchenko extract data handling to G4PAIModelData class
42// added sharing of internal data between threads (MT migration)
43//
44//
45// Class Description:
46//
47// Implementation of PAI model of energy loss and
48// delta-electron production by heavy charged particles
49
50// -------------------------------------------------------------------
51//
52
53#ifndef G4PAIModel_h
54#define G4PAIModel_h 1
55
57
58#include "G4VEmModel.hh"
60#include "globals.hh"
61#include <vector>
62
63class G4Region;
66class G4PAIModelData;
67
68class G4PAIModel final : public G4VEmModel, public G4VEmFluctuationModel
69{
70
71public:
72
73 explicit G4PAIModel(const G4ParticleDefinition* p = nullptr,
74 const G4String& nam = "PAI");
75
76 ~G4PAIModel() final;
77
78 void Initialise(const G4ParticleDefinition*, const G4DataVector&) final;
79
81 G4VEmModel* masterModel) final;
82
84 const G4MaterialCutsCouple* couple) final;
85
88 G4double kineticEnergy,
89 G4double cutEnergy) final;
90
93 G4double kineticEnergy,
94 G4double cutEnergy,
95 G4double maxEnergy) final;
96
97 void SampleSecondaries(std::vector<G4DynamicParticle*>*,
99 const G4DynamicParticle*,
100 G4double tmin,
101 G4double maxEnergy) final;
102
104 const G4DynamicParticle*,
105 const G4double, const G4double,
106 const G4double, const G4double) final;
107
109 const G4double, const G4double,
110 const G4double) final;
111
112 void DefineForRegion(const G4Region* r) final;
113
115
116 inline const std::vector<const G4MaterialCutsCouple*>& GetVectorOfCouples();
117
118 inline G4double ComputeMaxEnergy(G4double scaledEnergy);
119
120 inline void SetVerboseLevel(G4int verbose);
121
122protected:
123
125 G4double kinEnergy) final;
126
127 // hide assignment operator
128 G4PAIModel & operator=(const G4PAIModel &right) = delete;
129 G4PAIModel(const G4PAIModel&) = delete;
130
131private:
132
133 inline G4int FindCoupleIndex(const G4MaterialCutsCouple*);
134
135 inline void SetParticle(const G4ParticleDefinition* p);
136
137 G4int fVerbose;
138
139 G4PAIModelData* fModelData;
140
141 std::vector<const G4MaterialCutsCouple*> fMaterialCutsCoupleVector;
142 std::vector<const G4Region*> fPAIRegionVector;
143
144 const G4ParticleDefinition* fParticle;
145 const G4ParticleDefinition* fElectron;
146 const G4ParticleDefinition* fPositron;
147 G4ParticleChangeForLoss* fParticleChange;
148
149 G4double fMass;
150 G4double fRatio;
151 G4double fChargeSquare;
152 G4double fLowestTcut;
153};
154
156{
157 return fModelData;
158}
159
160inline const std::vector<const G4MaterialCutsCouple*>&
162{
163 return fMaterialCutsCoupleVector;
164}
165
167{
168 return MaxSecondaryEnergy(fParticle, scaledEnergy/fRatio);
169}
170
172{
173 fVerbose=verbose;
174}
175
176inline G4int G4PAIModel::FindCoupleIndex(const G4MaterialCutsCouple* couple)
177{
178 G4int idx = -1;
179 G4int jMatMax = (G4int)fMaterialCutsCoupleVector.size();
180 for(G4int jMat = 0;jMat < jMatMax; ++jMat) {
181 if(couple == fMaterialCutsCoupleVector[jMat]) {
182 idx = jMat;
183 break;
184 }
185 }
186 return idx;
187}
188
189inline void G4PAIModel::SetParticle(const G4ParticleDefinition* p)
190{
191 if(fParticle != p) {
192 fParticle = p;
193 fMass = fParticle->GetPDGMass();
194 fRatio = CLHEP::proton_mass_c2/fMass;
195 G4double q = fParticle->GetPDGCharge()/CLHEP::eplus;
196 fChargeSquare = q*q;
197 }
198}
199
200#endif
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
G4PAIModel & operator=(const G4PAIModel &right)=delete
void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel) final
Definition: G4PAIModel.cc:192
G4PAIModelData * GetPAIModelData()
Definition: G4PAIModel.hh:155
const std::vector< const G4MaterialCutsCouple * > & GetVectorOfCouples()
Definition: G4PAIModel.hh:161
G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) final
Definition: G4PAIModel.cc:212
void DefineForRegion(const G4Region *r) final
Definition: G4PAIModel.cc:399
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) final
Definition: G4PAIModel.cc:265
G4double ComputeMaxEnergy(G4double scaledEnergy)
Definition: G4PAIModel.hh:166
G4double SampleFluctuations(const G4MaterialCutsCouple *, const G4DynamicParticle *, const G4double, const G4double, const G4double, const G4double) final
Definition: G4PAIModel.cc:325
void Initialise(const G4ParticleDefinition *, const G4DataVector &) final
Definition: G4PAIModel.cc:104
G4double Dispersion(const G4Material *, const G4DynamicParticle *, const G4double, const G4double, const G4double) final
Definition: G4PAIModel.cc:362
G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy) final
Definition: G4PAIModel.cc:382
G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy) final
Definition: G4PAIModel.cc:235
G4double MinEnergyCut(const G4ParticleDefinition *, const G4MaterialCutsCouple *couple) final
Definition: G4PAIModel.cc:204
~G4PAIModel() final
Definition: G4PAIModel.cc:97
void SetVerboseLevel(G4int verbose)
Definition: G4PAIModel.hh:171
G4PAIModel(const G4PAIModel &)=delete
G4double GetPDGCharge() const