Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PAIPhotonModel.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$
27//
28// -------------------------------------------------------------------
29//
30// GEANT4 Class header file
31//
32//
33// File name: G4PAIPhotonModel
34//
35// Author: V. Grichine based on Vladimir Ivanchenko code
36//
37// Creation date: 05.10.2003
38//
39// Modifications:
40// 11.04.05 Major optimisation of internal interfaces (V.Ivantchenko)
41// 26-09-07 Fixed tmax computation (V.Ivantchenko)
42//
43//
44// Class Description:
45//
46// Implementation of PAI model of energy loss and
47// delta-electron production by heavy charged particles
48
49// -------------------------------------------------------------------
50//
51
52#ifndef G4PAIPhotonModel_h
53#define G4PAIPhotonModel_h 1
54
55#include <vector>
56#include "G4VEmModel.hh"
57#include "globals.hh"
59
61class G4PhysicsTable;
62class G4Region;
65
67{
68
69public:
70
71 G4PAIPhotonModel(const G4ParticleDefinition* p = 0, const G4String& nam = "PAIPhoton");
72
73 virtual ~G4PAIPhotonModel();
74
75 virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
76
77 virtual void InitialiseMe(const G4ParticleDefinition*);
78
81 G4double kineticEnergy,
82 G4double cutEnergy);
83
86 G4double kineticEnergy,
87 G4double cutEnergy,
88 G4double maxEnergy);
89
90 virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
92 const G4DynamicParticle*,
93 G4double tmin,
94 G4double maxEnergy);
95
97 const G4DynamicParticle*,
98 G4double&,
99 G4double&,
100 G4double&);
101
102 virtual G4double Dispersion( const G4Material*,
103 const G4DynamicParticle*,
104 G4double&,
105 G4double&);
106
107 void DefineForRegion(const G4Region* r) ;
110 void BuildLambdaVector(const G4MaterialCutsCouple* matCutsCouple);
111
112 G4double GetdNdxCut( G4int iPlace, G4double transferCut);
113 G4double GetdNdxPhotonCut( G4int iPlace, G4double transferCut);
114 G4double GetdNdxPlasmonCut( G4int iPlace, G4double transferCut);
115
116 G4double GetdEdxCut( G4int iPlace, G4double transferCut);
117
119 G4int iPlace, G4double scaledTkin );
121 G4int iPlace, G4double scaledTkin,G4double step, G4double cof );
123 G4double position, G4int iTransfer );
124
125protected:
126
128 G4double kinEnergy);
129
130private:
131
132 void SetParticle(const G4ParticleDefinition* p);
133
134 // hide assignment operator
135 G4PAIPhotonModel & operator=(const G4PAIPhotonModel &right);
137
138 // The vector over proton kinetic energies: the range of gammas
139
140 G4double fLowestKineticEnergy;
141 G4double fHighestKineticEnergy;
142 G4int fTotBin;
143 G4int fMeanNumber;
144 G4int fVerbose;
145 G4PhysicsLogVector* fProtonEnergyVector ;
146
147 // vectors
148
149 G4PhysicsTable* fPAItransferTable;
150 std::vector<G4PhysicsTable*> fPAIxscBank;
151
152 G4PhysicsTable* fPAIphotonTable;
153 std::vector<G4PhysicsTable*> fPAIphotonBank;
154
155 G4PhysicsTable* fPAIplasmonTable;
156 std::vector<G4PhysicsTable*> fPAIplasmonBank;
157
158 G4PhysicsTable* fPAIdEdxTable;
159 std::vector<G4PhysicsTable*> fPAIdEdxBank;
160
161 std::vector<const G4MaterialCutsCouple*> fMaterialCutsCoupleVector;
162 std::vector<const G4Region*> fPAIRegionVector;
163
164 size_t fMatIndex ;
165 G4double** fSandiaPhotoAbsCof ;
166 G4int fSandiaIntervalNumber ;
167
168 G4PhysicsLogVector* fdEdxVector ;
169 std::vector<G4PhysicsLogVector*> fdEdxTable ;
170
171 G4PhysicsLogVector* fLambdaVector ;
172 std::vector<G4PhysicsLogVector*> fLambdaTable ;
173
174 G4PhysicsLogVector* fdNdxCutVector ;
175 std::vector<G4PhysicsLogVector*> fdNdxCutTable ;
176
177 G4PhysicsLogVector* fdNdxCutPhotonVector ;
178 std::vector<G4PhysicsLogVector*> fdNdxCutPhotonTable ;
179
180 G4PhysicsLogVector* fdNdxCutPlasmonVector ;
181 std::vector<G4PhysicsLogVector*> fdNdxCutPlasmonTable ;
182
183
184 const G4ParticleDefinition* fParticle;
185 const G4ParticleDefinition* fElectron;
186 const G4ParticleDefinition* fPositron;
187 G4ParticleChangeForLoss* fParticleChange;
188
189 G4double fMass;
190 G4double fSpin;
191 G4double fChargeSquare;
192 G4double fRatio;
193 G4double fHighKinEnergy;
194 G4double fLowKinEnergy;
195 G4double fTwoln10;
196 G4double fBg2lim;
197 G4double fTaulim;
198 G4double fQc;
199
200 G4bool isInitialised;
201};
202
203#endif
204
205
206
207
208
209
210
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
void BuildLambdaVector(const G4MaterialCutsCouple *matCutsCouple)
G4double GetdNdxCut(G4int iPlace, G4double transferCut)
void DefineForRegion(const G4Region *r)
G4double GetdNdxPlasmonCut(G4int iPlace, G4double transferCut)
virtual G4double SampleFluctuations(const G4Material *, const G4DynamicParticle *, G4double &, G4double &, G4double &)
G4double GetdNdxPhotonCut(G4int iPlace, G4double transferCut)
G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
virtual void InitialiseMe(const G4ParticleDefinition *)
virtual ~G4PAIPhotonModel()
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
G4double GetPostStepTransfer(G4PhysicsTable *, G4PhysicsLogVector *, G4int iPlace, G4double scaledTkin)
G4double GetAlongStepTransfer(G4PhysicsTable *, G4PhysicsLogVector *, G4int iPlace, G4double scaledTkin, G4double step, G4double cof)
G4double GetEnergyTransfer(G4PhysicsTable *, G4int iPlace, G4double position, G4int iTransfer)
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
virtual G4double Dispersion(const G4Material *, const G4DynamicParticle *, G4double &, G4double &)
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy)
G4double GetdEdxCut(G4int iPlace, G4double transferCut)