Geant4  10.03.p01
 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 // $Id: G4PAIModel.hh 96934 2016-05-18 09:10:41Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class header file
31 //
32 //
33 // File name: G4PAIModel
34 //
35 // Author: V. Grichine based on Vladimir Ivanchenko code
36 //
37 // Creation date: 05.10.2003
38 //
39 // Modifications:
40 // 08-04-05 Major optimisation of internal interfaces (V.Ivantchenko)
41 // 26-09-07 Fixed tmax computation (V.Ivantchenko)
42 // 19.08.13 V.Ivanchenko extract data handling to G4PAIModelData class
43 // added sharing of internal data between threads (MT migration)
44 //
45 //
46 // Class Description:
47 //
48 // Implementation of PAI model of energy loss and
49 // delta-electron production by heavy charged particles
50 
51 // -------------------------------------------------------------------
52 //
53 
54 #ifndef G4PAIModel_h
55 #define G4PAIModel_h 1
56 
58 
59 #include "G4VEmModel.hh"
60 #include "G4VEmFluctuationModel.hh"
61 #include "globals.hh"
62 #include <vector>
63 
64 class G4Region;
67 class G4PAIModelData;
68 
70 {
71 
72 public:
73 
74  explicit G4PAIModel(const G4ParticleDefinition* p = nullptr,
75  const G4String& nam = "PAI");
76 
77  virtual ~G4PAIModel();
78 
79  virtual void Initialise(const G4ParticleDefinition*,
80  const G4DataVector&) final;
81 
82  virtual void InitialiseLocal(const G4ParticleDefinition*,
83  G4VEmModel* masterModel) final;
84 
86  const G4ParticleDefinition*,
87  G4double kineticEnergy,
88  G4double cutEnergy) final;
89 
91  const G4ParticleDefinition*,
92  G4double kineticEnergy,
93  G4double cutEnergy,
94  G4double maxEnergy) final;
95 
96  virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
97  const G4MaterialCutsCouple*,
98  const G4DynamicParticle*,
99  G4double tmin,
100  G4double maxEnergy) final;
101 
103  const G4DynamicParticle*,
104  G4double,
105  G4double,
106  G4double) final;
107 
108  virtual G4double Dispersion( const G4Material*,
109  const G4DynamicParticle*,
110  G4double,
111  G4double) final;
112 
113  virtual void DefineForRegion(const G4Region* r) final;
114 
116 
117  inline const std::vector<const G4MaterialCutsCouple*>& GetVectorOfCouples();
118 
119  inline G4double ComputeMaxEnergy(G4double scaledEnergy);
120 
121  inline void SetVerboseLevel(G4int verbose);
122 
123 protected:
124 
126  G4double kinEnergy) final;
127 
128 private:
129 
130  inline G4int FindCoupleIndex(const G4MaterialCutsCouple*);
131 
132  inline void SetParticle(const G4ParticleDefinition* p);
133 
134  // hide assignment operator
135  G4PAIModel & operator=(const G4PAIModel &right) = delete;
136  G4PAIModel(const G4PAIModel&) = delete;
137 
138  G4int fVerbose;
139 
140  G4PAIModelData* fModelData;
141 
142  std::vector<const G4MaterialCutsCouple*> fMaterialCutsCoupleVector;
143  std::vector<const G4Region*> fPAIRegionVector;
144 
145  const G4ParticleDefinition* fParticle;
146  const G4ParticleDefinition* fElectron;
147  const G4ParticleDefinition* fPositron;
148  G4ParticleChangeForLoss* fParticleChange;
149 
150  G4double fMass;
151  G4double fRatio;
152  G4double fChargeSquare;
153 };
154 
156 {
157  return fModelData;
158 }
159 
160 inline const std::vector<const G4MaterialCutsCouple*>&
162 {
163  return fMaterialCutsCoupleVector;
164 }
165 
167 {
168  return MaxSecondaryEnergy(fParticle, scaledEnergy/fRatio);
169 }
170 
171 inline void G4PAIModel::SetVerboseLevel(G4int verbose)
172 {
173  fVerbose=verbose;
174 }
175 
176 inline G4int G4PAIModel::FindCoupleIndex(const G4MaterialCutsCouple* couple)
177 {
178  G4int idx = -1;
179  size_t jMatMax = fMaterialCutsCoupleVector.size();
180  for(size_t jMat = 0;jMat < jMatMax; ++jMat) {
181  if(couple == fMaterialCutsCoupleVector[jMat]) {
182  idx = jMat;
183  break;
184  }
185  }
186  return idx;
187 }
188 
189 inline 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
201 
202 
203 
204 
205 
206 
207 
const std::vector< const G4MaterialCutsCouple * > & GetVectorOfCouples()
Definition: G4PAIModel.hh:161
G4PAIModelData * GetPAIModelData()
Definition: G4PAIModel.hh:155
G4PAIModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="PAI")
Definition: G4PAIModel.cc:77
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) final
Definition: G4PAIModel.cc:257
static constexpr double proton_mass_c2
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) final
Definition: G4PAIModel.cc:104
const char * p
Definition: xmltok.h:285
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy) final
Definition: G4PAIModel.cc:227
virtual G4double Dispersion(const G4Material *, const G4DynamicParticle *, G4double, G4double) final
Definition: G4PAIModel.cc:350
int G4int
Definition: G4Types.hh:78
virtual G4double SampleFluctuations(const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double, G4double, G4double) final
Definition: G4PAIModel.cc:315
virtual void DefineForRegion(const G4Region *r) final
Definition: G4PAIModel.cc:386
static constexpr double eplus
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel) final
Definition: G4PAIModel.cc:192
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) final
Definition: G4PAIModel.cc:204
virtual ~G4PAIModel()
Definition: G4PAIModel.cc:97
G4double GetPDGMass() const
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy) final
Definition: G4PAIModel.cc:369
double G4double
Definition: G4Types.hh:76
G4double GetPDGCharge() const
G4double ComputeMaxEnergy(G4double scaledEnergy)
Definition: G4PAIModel.hh:166
void SetVerboseLevel(G4int verbose)
Definition: G4PAIModel.hh:171