Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VXTRenergyLoss.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///////////////////////////////////////////////////////////////////////////
30//
31// base class for 'fast' parametrisation model describing X-ray transition
32// created in some G4Envelope. Anglur distribuiton is very rough !!! (see DoIt
33// method
34//
35// History:
36// 06.10.05 V. Grichine first step to discrete process
37// 15.01.02 V. Grichine first version
38// 28.07.05, P.Gumplinger add G4ProcessType to constructor
39// 28.09.07, V.Ivanchenko general cleanup without change of algorithms
40//
41
42#ifndef G4VXTRenergyLoss_h
43#define G4VXTRenergyLoss_h 1
44
45#include <complex>
46#include "globals.hh"
47#include "Randomize.hh"
48
49#include "G4LogicalVolume.hh"
50
51#include "G4PhysicsTable.hh"
52#include "G4PhysicsLogVector.hh"
53#include "G4Gamma.hh"
54#include "G4ThreeVector.hh"
55#include "G4ParticleMomentum.hh"
56#include "G4Step.hh"
57#include "G4Track.hh"
59#include "G4VDiscreteProcess.hh"
60#include "G4DynamicParticle.hh"
61#include "G4Material.hh"
62#include "G4PhysicsTable.hh"
65#include "G4Integrator.hh"
66#include "G4ParticleChange.hh"
67
68class G4SandiaTable;
72
73class G4VXTRenergyLoss : public G4VDiscreteProcess // G4VContinuousProcess
74{
75public:
76
77 explicit G4VXTRenergyLoss (G4LogicalVolume *anEnvelope,G4Material*,
79 const G4String & processName = "XTRenergyLoss",
81 virtual ~G4VXTRenergyLoss ();
82
83 // These virtual has to be implemented in inherited particular TR radiators
84
85 virtual G4double GetStackFactor( G4double energy, G4double gamma,
86 G4double varAngle );
87
88 virtual G4bool IsApplicable(const G4ParticleDefinition&) override;
89
90 virtual G4VParticleChange* PostStepDoIt(const G4Track& aTrack,
91 const G4Step& aStep) override;
92
93 virtual G4double GetMeanFreePath(const G4Track& aTrack,
94 G4double previousStepSize,
95 G4ForceCondition* condition) override;
96
97 virtual void BuildPhysicsTable(const G4ParticleDefinition&) override;
98 void BuildEnergyTable() ;
100
101 void BuildTable(){} ;
102 void BuildAngleTable() ;
103 void BuildGlobalAngleTable() ;
104
106 G4double gamma,
107 G4double varAngle ) ;
108
110
111 virtual G4double SpectralXTRdEdx(G4double energy) ;
112
114
115 G4double AngleXTRdEdx(G4double varAngle) ;
116
117
118 /////////////////////////////////////////////////////////////
119
121 G4double gamma,
122 G4double varAngle ) const ;
123
124
125 // for photon energy distribution tables
126
129
130 // for photon angle distribution tables
131
134
135 void GetNumberOfPhotons() ;
136
137 // Auxiliary functions for plate/gas material parameters
138
143 void GetPlateZmuProduct() ;
145
150 void GetGasZmuProduct();
152
156
157 G4double GetXTRrandomEnergy( G4double scaledTkin, G4int iTkin );
158 G4double GetXTRenergy( G4int iPlace, G4double position, G4int iTransfer );
159
160 G4double GetRandomAngle( G4double energyXTR, G4int iTkin );
162
166
167 void SetGamma(G4double gamma) {fGamma = gamma;};
168 void SetEnergy(G4double energy) {fEnergy = energy;};
169 void SetVarAngle(G4double varAngle){fVarAngle = varAngle;};
170 void SetAngleRadDistr(G4bool pAngleRadDistr){fAngleRadDistr=pAngleRadDistr;};
171 void SetCompton(G4bool pC){fCompton=pC;};
172
176
177protected:
178
179 G4ParticleDefinition* fPtrGamma ; // pointer to TR photon
180
181 G4double* fGammaCutInKineticEnergy ; // TR photon cut in energy array
182
183 G4double fGammaTkinCut ; // Tkin cut of TR photon in current mat.
187
190
191 G4double fTheMinEnergyTR; // min TR energy
192 G4double fTheMaxEnergyTR; // max TR energy
193 G4double fMinEnergyTR; // min TR energy in material
194 G4double fMaxEnergyTR; // max TR energy in material
195 G4double fTheMaxAngle; // max theta of TR quanta
196 G4double fTheMinAngle; // max theta of TR quanta
197 G4double fMaxThetaTR; // max theta of TR quanta
198 G4int fBinTR; // number of bins in TR vectors
199
200 G4double fMinProtonTkin; // min Tkin of proton in tables
201 G4double fMaxProtonTkin; // max Tkin of proton in tables
202 G4int fTotBin; // number of bins in log scale
203 G4double fGamma; // current Lorentz factor
204 G4double fEnergy; // energy and
205 G4double fVarAngle; // angle squared
207
208 G4double fPlasmaCof ; // physical consts for plasma energy
210
215 G4double fSigma2; // plasma energy Sq of matter1/2
216
220
226
228
230
232
234 std::vector<G4PhysicsTable*> fAngleBank;
235
236private:
237
238 // copy constructor and hide assignment operator
240 G4VXTRenergyLoss & operator=(const G4VXTRenergyLoss &right) = delete;
241
242};
243
244#endif
G4double condition(const G4ErrorSymMatrix &m)
G4ForceCondition
G4ProcessType
@ fElectromagnetic
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
std::complex< G4double > G4complex
Definition: G4Types.hh:88
int G4int
Definition: G4Types.hh:85
Definition: G4Step.hh:62
G4double GetPlateLinearPhotoAbs(G4double)
G4double AngleSpectralXTRdEdx(G4double energy)
void SetEnergy(G4double energy)
G4double XTRNAngleDensity(G4double varAngle)
G4PhysicsTable * fEnergyDistrTable
G4double GetAngleXTR(G4int iTR, G4double position, G4int iAngle)
virtual G4bool IsApplicable(const G4ParticleDefinition &) override
G4PhysicsTable * fAngleForEnergyTable
G4PhysicsLogVector * fProtonEnergyVector
G4PhysicsLogVector * fXTREnergyVector
G4double GetGasFormationZone(G4double, G4double, G4double)
G4SandiaTable * fPlatePhotoAbsCof
G4double OneBoundaryXTRNdensity(G4double energy, G4double gamma, G4double varAngle) const
G4PhysicsFreeVector * GetAngleVector(G4double energy, G4int n)
G4PhysicsTable * fAngleDistrTable
G4double GetRandomAngle(G4double energyXTR, G4int iTkin)
G4complex GetPlateComplexFZ(G4double, G4double, G4double)
virtual void BuildPhysicsTable(const G4ParticleDefinition &) override
G4complex OneInterfaceXTRdEdx(G4double energy, G4double gamma, G4double varAngle)
virtual G4double GetMeanFreePath(const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition) override
G4double GetXTRenergy(G4int iPlace, G4double position, G4int iTransfer)
G4double XTRNAngleSpectralDensity(G4double energy)
void SetGamma(G4double gamma)
G4double SpectralAngleXTRdEdx(G4double varAngle)
G4double GetXTRrandomEnergy(G4double scaledTkin, G4int iTkin)
void SetAngleRadDistr(G4bool pAngleRadDistr)
G4SandiaTable * fGasPhotoAbsCof
void SetCompton(G4bool pC)
G4LogicalVolume * fEnvelope
std::vector< G4PhysicsTable * > fAngleBank
G4PhysicsLogVector * GetProtonVector()
virtual ~G4VXTRenergyLoss()
virtual G4double SpectralXTRdEdx(G4double energy)
void SetVarAngle(G4double varAngle)
G4double GetComptonPerAtom(G4double, G4double)
G4double GetGasCompton(G4double)
G4double XTRNSpectralAngleDensity(G4double varAngle)
G4double GetPlateFormationZone(G4double, G4double, G4double)
G4double GetGasLinearPhotoAbs(G4double)
G4double XTRNSpectralDensity(G4double energy)
virtual G4double GetStackFactor(G4double energy, G4double gamma, G4double varAngle)
G4ParticleDefinition * fPtrGamma
virtual G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep) override
G4ParticleChange fParticleChange
G4double * fGammaCutInKineticEnergy
G4complex GetGasComplexFZ(G4double, G4double, G4double)
G4double GetPlateCompton(G4double)
G4double AngleXTRdEdx(G4double varAngle)