Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ForwardXrayTR.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// G4ForwardXrayTR
27//
28// Class for description
29//
30// Class for forward X-ray transition radiation generated
31// by relativistic charged particle crossed interface between material 1
32// and material 2 (1 -> 2)
33
34// History:
35// 22.09.97, V. Grichine ([email protected])
36// 26.01.00, V.Grichine, new constructor and protected DM for fast sim. models
37// 10.03.03, V.Ivanchenko migrade to "cut per region"
38// 03.06.03, V.Ivanchenko fix compilation warnings
39
40#ifndef G4FORWARDXRAYTR_H
41#define G4FORWARDXRAYTR_H
42
43#include "globals.hh"
44#include "G4Track.hh"
46#include "G4VParticleChange.hh"
47
49class G4PhysicsTable;
51
53{
54 public:
55 explicit G4ForwardXrayTR(const G4String& matName1, const G4String& matName2,
56 const G4String& processName = "XrayTR");
57
58 explicit G4ForwardXrayTR(const G4String& processName = "XrayTR");
59
61
62 G4ForwardXrayTR(const G4ForwardXrayTR& right) = delete;
63 G4ForwardXrayTR& operator=(const G4ForwardXrayTR& right) = delete;
64
65 /////////////////////// Methods /////////////////////////////////
66
67 void ProcessDescription(std::ostream&) const override;
68 void DumpInfo() const override { ProcessDescription(G4cout); };
69
70 void BuildXrayTRtables();
71
73 G4ForceCondition* condition) override;
74
76 const G4Step& aStep) override;
77
78 G4double GetEnergyTR(G4int iMat, G4int jMat, G4int iTkin) const;
79
80 G4double GetThetaTR(G4int iMat, G4int jMat, G4int iTkin) const;
81
82 ///////////////////// Angle distribution /////////////////////////////
83
85 G4double varAngle) const override;
86
87 G4double AngleDensity(G4double energy, G4double varAngle) const;
88
90 G4double varAngle) const;
91
92 G4double AngleSum(G4double varAngle1, G4double varAngle2) const;
93
94 ///////////////////////// Energy distribution ///////////////////////////////
95
97
98 G4double AngleInterval(G4double energy, G4double varAngle1,
99 G4double varAngle2) const;
100
101 G4double EnergySum(G4double energy1, G4double energy2) const;
102
103 /////////////////////////// Access functions ////////////////////////////
104
107
108 static G4int GetSympsonNumber();
109 static G4int GetBinTR();
110
111 static G4double GetMinProtonTkin();
112 static G4double GetMaxProtonTkin();
113 static G4int GetTotBin();
114
115 protected: // for access from X-ray TR fast simulation models
116 static constexpr G4double fTheMinEnergyTR =
117 1. * CLHEP::keV; // static min TR energy
118 static constexpr G4double fTheMaxEnergyTR =
119 100. * CLHEP::keV; // static max TR energy
120 static constexpr G4double fTheMaxAngle = 1.0e-3; // max theta of TR quanta
121 static constexpr G4double fTheMinAngle = 5.0e-6; // min theta of TR quanta
122 static constexpr G4double fMinProtonTkin =
123 100. * CLHEP::GeV; // min Tkin of proton in tables
124 static constexpr G4double fMaxProtonTkin =
125 100. * CLHEP::TeV; // max Tkin of proton in tables
126 static constexpr G4double fPlasmaCof =
127 4.0 * CLHEP::pi * CLHEP::fine_structure_const * CLHEP::hbarc *
128 CLHEP::hbarc * CLHEP::hbarc /
129 CLHEP::electron_mass_c2; // physical consts for plasma energy
130 static constexpr G4double fCofTR = CLHEP::fine_structure_const / CLHEP::pi;
131
132 static constexpr G4int fSympsonNumber =
133 100; // Accuracy of Sympson integration
134 static constexpr G4int fBinTR = 50; // number of bins in TR vectors
135 static constexpr G4int fTotBin = 50; // number of bins in log scale
136
137 const std::vector<G4double>* fGammaCutInKineticEnergy;
138 // TR photon cut in energy array
139
140 G4ParticleDefinition* fPtrGamma; // pointer to TR photon
141
144
146
147 G4double fMinEnergyTR; // min TR energy in material
148 G4double fMaxEnergyTR; // max TR energy in material
149 G4double fMaxThetaTR; // max theta of TR quanta
150 G4double fGamma; // current Lorentz factor
151 G4double fGammaTkinCut; // Tkin cut of TR photon in current mat.
152 G4double fSigma1; // plasma energy Sq of matter1
153 G4double fSigma2; // plasma energy Sq of matter2
154
155 G4int secID = -1; // creator modelID
156};
157
158#endif // G4FORWARDXRAYTR_H
G4double condition(const G4ErrorSymMatrix &m)
G4ForceCondition
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
G4GLOB_DLL std::ostream G4cout
static constexpr G4double fPlasmaCof
G4double SpectralDensity(G4double energy, G4double x) const
static constexpr G4int fTotBin
G4ForwardXrayTR(const G4String &matName1, const G4String &matName2, const G4String &processName="XrayTR")
static G4double GetMinProtonTkin()
G4double EnergyInterval(G4double energy1, G4double energy2, G4double varAngle) const
static constexpr G4double fTheMaxEnergyTR
G4double EnergySum(G4double energy1, G4double energy2) const
static G4int GetSympsonNumber()
static constexpr G4int fBinTR
G4ForwardXrayTR & operator=(const G4ForwardXrayTR &right)=delete
static constexpr G4double fCofTR
static G4double GetMaxProtonTkin()
static constexpr G4double fTheMinAngle
static G4int GetTotBin()
static G4int GetBinTR()
const std::vector< G4double > * fGammaCutInKineticEnergy
G4PhysicsTable * fEnergyDistrTable
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep) override
G4PhysicsTable * fAngleDistrTable
G4double AngleDensity(G4double energy, G4double varAngle) const
static constexpr G4double fTheMinEnergyTR
static constexpr G4int fSympsonNumber
G4double GetEnergyTR(G4int iMat, G4int jMat, G4int iTkin) const
static constexpr G4double fMaxProtonTkin
static constexpr G4double fTheMaxAngle
G4PhysicsTable * GetAngleDistrTable()
static constexpr G4double fMinProtonTkin
G4double AngleInterval(G4double energy, G4double varAngle1, G4double varAngle2) const
G4double GetThetaTR(G4int iMat, G4int jMat, G4int iTkin) const
G4PhysicsLogVector * fProtonEnergyVector
G4ParticleDefinition * fPtrGamma
G4ForwardXrayTR(const G4ForwardXrayTR &right)=delete
G4double SpectralAngleTRdensity(G4double energy, G4double varAngle) const override
void ProcessDescription(std::ostream &) const override
G4double AngleSum(G4double varAngle1, G4double varAngle2) const
G4PhysicsTable * GetEnergyDistrTable()
void DumpInfo() const override
G4double GetMeanFreePath(const G4Track &, G4double, G4ForceCondition *condition) override