Geant4  10.01
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 73607 2013-09-02 10:04:03Z 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 
57 #include <CLHEP/Units/PhysicalConstants.h>
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  G4PAIModel(const G4ParticleDefinition* p = 0, const G4String& nam = "PAI");
75 
76  virtual ~G4PAIModel();
77 
78  virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
79 
80  virtual void InitialiseLocal(const G4ParticleDefinition*,
81  G4VEmModel* masterModel);
82 
84  const G4ParticleDefinition*,
85  G4double kineticEnergy,
86  G4double cutEnergy);
87 
89  const G4ParticleDefinition*,
90  G4double kineticEnergy,
91  G4double cutEnergy,
92  G4double maxEnergy);
93 
94  virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
95  const G4MaterialCutsCouple*,
96  const G4DynamicParticle*,
97  G4double tmin,
98  G4double maxEnergy);
99 
101  const G4DynamicParticle*,
102  G4double,
103  G4double,
104  G4double);
105 
106  virtual G4double Dispersion( const G4Material*,
107  const G4DynamicParticle*,
108  G4double,
109  G4double);
110 
111  void DefineForRegion(const G4Region* r);
112 
114 
115  inline G4double ComputeMaxEnergy(G4double scaledEnergy);
116 
117  inline void SetVerboseLevel(G4int verbose);
118 
119 protected:
120 
122  G4double kinEnergy);
123 
124 private:
125 
127 
128  inline void SetParticle(const G4ParticleDefinition* p);
129 
130  // hide assignment operator
132  G4PAIModel(const G4PAIModel&);
133 
135 
137 
138  std::vector<const G4MaterialCutsCouple*> fMaterialCutsCoupleVector;
139  std::vector<const G4Region*> fPAIRegionVector;
140 
145 
149 
151 };
152 
154 {
155  return fModelData;
156 }
157 
159 {
160  return MaxSecondaryEnergy(fParticle, scaledEnergy/fRatio);
161 }
162 
163 inline void G4PAIModel::SetVerboseLevel(G4int verbose)
164 {
165  fVerbose=verbose;
166 }
167 
169 {
170  G4int idx = -1;
171  size_t jMatMax = fMaterialCutsCoupleVector.size();
172  for(size_t jMat = 0;jMat < jMatMax; ++jMat) {
173  if(couple == fMaterialCutsCoupleVector[jMat]) {
174  idx = jMat;
175  break;
176  }
177  }
178  return idx;
179 }
180 
182 {
183  if(fParticle != p) {
184  fParticle = p;
186  fRatio = CLHEP::proton_mass_c2/fMass;
188  fChargeSquare = q*q;
189  }
190 }
191 
192 #endif
193 
194 
195 
196 
197 
198 
199 
G4PAIModelData * GetPAIModelData()
Definition: G4PAIModel.hh:153
std::vector< const G4MaterialCutsCouple * > fMaterialCutsCoupleVector
Definition: G4PAIModel.hh:138
const G4ParticleDefinition * fPositron
Definition: G4PAIModel.hh:143
G4int FindCoupleIndex(const G4MaterialCutsCouple *)
Definition: G4PAIModel.hh:168
G4double fMass
Definition: G4PAIModel.hh:146
std::vector< const G4Region * > fPAIRegionVector
Definition: G4PAIModel.hh:139
G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy)
Definition: G4PAIModel.cc:371
G4double fRatio
Definition: G4PAIModel.hh:147
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4PAIModel.cc:111
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel)
Definition: G4PAIModel.cc:192
G4int fVerbose
Definition: G4PAIModel.hh:134
G4bool isInitialised
Definition: G4PAIModel.hh:150
int G4int
Definition: G4Types.hh:78
virtual G4double Dispersion(const G4Material *, const G4DynamicParticle *, G4double, G4double)
Definition: G4PAIModel.cc:340
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
Definition: G4PAIModel.cc:225
G4PAIModel & operator=(const G4PAIModel &right)
void SetParticle(const G4ParticleDefinition *p)
Definition: G4PAIModel.hh:181
const G4ParticleDefinition * fParticle
Definition: G4PAIModel.hh:141
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
Definition: G4PAIModel.cc:249
G4double fChargeSquare
Definition: G4PAIModel.hh:148
G4PAIModelData * fModelData
Definition: G4PAIModel.hh:136
bool G4bool
Definition: G4Types.hh:79
const G4ParticleDefinition * fElectron
Definition: G4PAIModel.hh:142
G4ParticleChangeForLoss * fParticleChange
Definition: G4PAIModel.hh:144
virtual ~G4PAIModel()
Definition: G4PAIModel.cc:104
G4double GetPDGMass() const
G4PAIModel(const G4ParticleDefinition *p=0, const G4String &nam="PAI")
Definition: G4PAIModel.cc:82
virtual G4double SampleFluctuations(const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double, G4double, G4double)
Definition: G4PAIModel.cc:305
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy)
Definition: G4PAIModel.cc:201
double G4double
Definition: G4Types.hh:76
void DefineForRegion(const G4Region *r)
Definition: G4PAIModel.cc:388
static const double eplus
Definition: G4SIunits.hh:178
G4double GetPDGCharge() const
G4double ComputeMaxEnergy(G4double scaledEnergy)
Definition: G4PAIModel.hh:158
void SetVerboseLevel(G4int verbose)
Definition: G4PAIModel.hh:163