Geant4 11.2.2
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// base class for 'fast' parametrisation model describing X-ray transition
29// created in some G4Envelope. Angular distribuiton is very rough !!! (see DoIt
30// method
31//
32// History:
33// 06.10.05 V. Grichine first step to discrete process
34// 15.01.02 V. Grichine first version
35// 28.07.05, P.Gumplinger add G4ProcessType to constructor
36// 28.09.07, V.Ivanchenko general cleanup without change of algorithms
37// 19.09.21, V. Grichine, set/get functions for angle anf energy ranges and number of bins
38
39#ifndef G4VXTRenergyLoss_h
40#define G4VXTRenergyLoss_h 1
41
42#include "globals.hh"
43#include "G4Gamma.hh"
44#include "G4LogicalVolume.hh"
45#include "G4Material.hh"
46#include "G4ParticleChange.hh"
47#include "G4PhysicsTable.hh"
48#include "G4Step.hh"
49#include "G4Track.hh"
50#include "G4VDiscreteProcess.hh"
51
52class G4SandiaTable;
57
59{
60 public:
61 explicit G4VXTRenergyLoss(G4LogicalVolume* anEnvelope, G4Material*,
63 const G4String& processName = "XTRenergyLoss",
65 virtual ~G4VXTRenergyLoss();
66
67 virtual void ProcessDescription(std::ostream&) const override;
68 virtual void DumpInfo() const override { ProcessDescription(G4cout); };
69
72
73 // Virtual methods to be implemented in inherited particular TR radiators
74 virtual G4double GetStackFactor(G4double energy, G4double gamma,
75 G4double varAngle);
76
77 virtual G4bool IsApplicable(const G4ParticleDefinition&) override;
78
79 virtual G4VParticleChange* PostStepDoIt(const G4Track& aTrack,
80 const G4Step& aStep) override;
81
82 virtual G4double GetMeanFreePath(const G4Track& aTrack,
83 G4double previousStepSize,
84 G4ForceCondition* condition) override;
85
86 virtual void BuildPhysicsTable(const G4ParticleDefinition&) override;
87 void BuildEnergyTable();
89
90 void BuildTable(){};
91 void BuildAngleTable();
93
95 G4double varAngle);
96
98
99 virtual G4double SpectralXTRdEdx(G4double energy);
100
102
104
106 G4double varAngle) const;
107
108 // for photon energy distribution tables
111
112 // for photon angle distribution tables
115
116 void GetNumberOfPhotons();
117
118 // Auxiliary functions for plate/gas material parameters
123 void GetPlateZmuProduct();
125
130 void GetGasZmuProduct();
132
136
137 G4double GetXTRrandomEnergy(G4double scaledTkin, G4int iTkin);
138 G4double GetXTRenergy(G4int iPlace, G4double position, G4int iTransfer);
139
140 G4double GetRandomAngle(G4double energyXTR, G4int iTkin);
142
143 // set/get methods for class fields
144
145 void SetGamma(G4double gamma) { fGamma = gamma; };
146 G4double GetGamma() { return fGamma; };
147 void SetEnergy(G4double energy) { fEnergy = energy; };
149 void SetVarAngle(G4double varAngle) { fVarAngle = varAngle; };
151 void SetCompton(G4bool pC) { fCompton = pC; };
153
155 void SetKrange( G4int kk ){ fKrange = kk;};
156
157
162
167
168 void SetMinEnergyTR(G4double minetr){ fMinEnergyTR = minetr;};
170 void SetMaxEnergyTR(G4double maxetr){ fMaxEnergyTR = maxetr;};
172
173 void SetTheMinAngle(G4double minang){ fTheMinAngle = minang;};
175 void SetTheMaxAngle(G4double maxang){ fTheMaxAngle = maxang;};
177
178 void SetMinThetaTR(G4double minatr){ fMinThetaTR = minatr;};
180 void SetMaxThetaTR(G4double maxatr){ fMaxThetaTR = maxatr;};
182
183 // modes of XTR angle distribution
184
185 void SetFastAngle(G4bool fatr){ fFastAngle = fatr;};
189
190
191
192
194 G4int GetTotBin() { return fTotBin; };
196
197 protected:
198 // min TR energy
200 // max TR energy
202 G4double fTheMinAngle; // min theta of TR quanta
203 G4double fTheMaxAngle; // 1.e-4; // max theta of TR quanta
204
205 // static const members
206
207 // min Tkin of proton in tables
208 static constexpr G4double fMinProtonTkin = 100. * CLHEP::GeV;
209 // max Tkin of proton in tables
210 static constexpr G4double fMaxProtonTkin = 100. * CLHEP::TeV;
211 // physical constants for plasma energy
212 static constexpr G4double fPlasmaCof =
213 4. * CLHEP::pi * CLHEP::fine_structure_const * CLHEP::hbarc * CLHEP::hbarc *
214 CLHEP::hbarc / CLHEP::electron_mass_c2;
215 static constexpr G4double fCofTR = CLHEP::fine_structure_const / CLHEP::pi;
216
217 G4int fTotBin; // number of bins in log-gamma scale
218 G4int fBinTR; // number of bins in TR energy-angle vectors
220 G4ParticleDefinition* fPtrGamma; // pointer to TR photon
221
222 G4double* fGammaCutInKineticEnergy; // TR photon cut in energy array
231
233 std::vector<G4PhysicsTable*> fAngleBank;
234
235 G4double fGammaTkinCut; // Tkin cut of TR photon in current mat.
236 G4double fMinEnergyTR; // min TR energy in material
237 G4double fMaxEnergyTR; // max TR energy in material
238 G4double fMinThetaTR, fMaxThetaTR; // min-max theta of TR quanta
244 G4double fGamma; // current Lorentz factor
245 G4double fEnergy; // energy and
246 G4double fVarAngle; // angle squared!
249 G4double fSigma2; // plasma energy Sq of matter1/2
250
254
258
259 G4int secID = -1; // creator modelID
260};
261
262#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
G4GLOB_DLL std::ostream G4cout
G4double GetPlateLinearPhotoAbs(G4double)
G4double AngleSpectralXTRdEdx(G4double energy)
G4double GetTheMinEnergyTR()
void SetEnergy(G4double energy)
G4double XTRNAngleDensity(G4double varAngle)
G4PhysicsTable * fEnergyDistrTable
void SetAlphaPlate(G4double ap)
void SetFastAngle(G4bool fatr)
G4VXTRenergyLoss & operator=(const G4VXTRenergyLoss &right)=delete
void SetTheMinAngle(G4double minang)
G4double GetAngleXTR(G4int iTR, G4double position, G4int iAngle)
virtual G4bool IsApplicable(const G4ParticleDefinition &) override
void SetKrange(G4int kk)
G4PhysicsTable * fAngleForEnergyTable
static constexpr G4double fMaxProtonTkin
G4PhysicsLogVector * fProtonEnergyVector
G4PhysicsLogVector * fXTREnergyVector
void SetAlphaGas(G4double ag)
void SetTheMaxEnergyTR(G4double maxetr)
G4double GetGasFormationZone(G4double, G4double, G4double)
G4SandiaTable * fPlatePhotoAbsCof
G4double OneBoundaryXTRNdensity(G4double energy, G4double gamma, G4double varAngle) const
G4PhysicsFreeVector * GetAngleVector(G4double energy, G4int n)
void SetMinEnergyTR(G4double minetr)
G4VXTRenergyLoss(G4VXTRenergyLoss &)=delete
void SetMaxEnergyTR(G4double maxetr)
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)
void SetMinThetaTR(G4double minatr)
G4VXTRenergyLoss(G4LogicalVolume *anEnvelope, G4Material *, G4Material *, G4double, G4double, G4int, const G4String &processName="XTRenergyLoss", G4ProcessType type=fElectromagnetic)
G4double XTRNAngleSpectralDensity(G4double energy)
void SetGamma(G4double gamma)
G4double SpectralAngleXTRdEdx(G4double varAngle)
G4double GetXTRrandomEnergy(G4double scaledTkin, G4int iTkin)
G4SandiaTable * fGasPhotoAbsCof
void SetCompton(G4bool pC)
static constexpr G4double fCofTR
G4LogicalVolume * fEnvelope
std::vector< G4PhysicsTable * > fAngleBank
G4PhysicsLogVector * GetProtonVector()
void SetTheMaxAngle(G4double maxang)
void SetTheMinEnergyTR(G4double minetr)
static constexpr G4double fPlasmaCof
virtual G4double SpectralXTRdEdx(G4double energy)
void SetVarAngle(G4double varAngle)
G4double GetComptonPerAtom(G4double, G4double)
virtual void DumpInfo() const override
void SetMaxThetaTR(G4double maxatr)
void SetAngleRadDistr(G4bool fatr)
static constexpr G4double fMinProtonTkin
G4double GetGasCompton(G4double)
G4double XTRNSpectralAngleDensity(G4double varAngle)
virtual void ProcessDescription(std::ostream &) const override
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 GetTheMaxEnergyTR()
G4double AngleXTRdEdx(G4double varAngle)