Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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 //
27 // $Id: G4ForwardXrayTR.hh 97385 2016-06-02 09:59:53Z gcosmo $
28 //
29 // G4ForwardXrayTR
30 //
31 // Class for description
32 //
33 // Class for forward X-ray transition radiation generated
34 // by relativistic charged particle crossed interface between material 1
35 // and material 2 (1 -> 2)
36 
37 // History:
38 // 22.09.97, V. Grichine (Vladimir.Grichine@cern.ch)
39 // 26.01.00, V.Grichine, new constructor and protected DM for fast sim. models
40 // 10.03.03, V.Ivanchenko migrade to "cut per region"
41 // 03.06.03, V.Ivanchenko fix compilation warnings
42 
43 #ifndef G4FORWARDXRAYTR_H
44 #define G4FORWARDXRAYTR_H
45 
46 
47 #include "globals.hh"
48 #include "templates.hh"
49 #include "geomdefs.hh"
50 #include "Randomize.hh"
51 #include "G4Step.hh"
52 #include "G4VDiscreteProcess.hh"
53 #include "G4DynamicParticle.hh"
54 #include "G4Material.hh"
56 #include "G4LogicalSkinSurface.hh"
57 #include "G4OpticalSurface.hh"
58 #include "G4OpticalPhoton.hh"
60 
61 #include "G4TransitionRadiation.hh"
62 #include "G4PhysicsTable.hh"
63 #include "G4Gamma.hh"
64 #include "G4PhysicsLogVector.hh"
65 
67 {
68 public:
69 
70  // Constructors
71 
72  explicit G4ForwardXrayTR( const G4String& matName1, // G4Material* pMat1,
73  const G4String& matName2, // G4Material* pMat2,
74  const G4String& processName="XrayTR" );
75 
76  explicit G4ForwardXrayTR( const G4String& processName="XrayTR" );
77 
78  // Destructor // virtual
79 
80  virtual ~G4ForwardXrayTR();
81 
83 
84  void BuildXrayTRtables();
85 
87  G4ForceCondition* condition) override;
88 
89  G4VParticleChange* PostStepDoIt( const G4Track& aTrack,
90  const G4Step& aStep ) override;
91 
92  G4double GetEnergyTR(G4int iMat, G4int jMat, G4int iTkin) const;
93 
94  G4double GetThetaTR(G4int iMat, G4int jMat, G4int iTkin) const;
95 
96 
98 //
99 
101  G4double varAngle ) const override;
102 
104  G4double varAngle ) const;
105 
107  G4double energy2,
108  G4double varAngle ) const;
109 
110 G4double AngleSum( G4double varAngle1,
111  G4double varAngle2 ) const;
112 
114 
116  G4double x ) const;
117 
119  G4double varAngle1,
120  G4double varAngle2 ) const;
121 
122 G4double EnergySum( G4double energy1,
123  G4double energy2 ) const;
124 
125 
127 
130 
131  static G4int GetSympsonNumber() { return fSympsonNumber; };
132  static G4int GetBinTR() { return fBinTR; };
133 
136  static G4int GetTotBin() { return fTotBin; };
137 
138 
139 protected: // for access from X-ray TR fast simulation models
140 
141  // private : /////////////// Data members ///////////////////////////
142 
143 G4ParticleDefinition* fPtrGamma; // pointer to TR photon
144 
145 const std::vector<G4double>* fGammaCutInKineticEnergy;
146  // TR photon cut in energy array
147 G4double fGammaTkinCut; // Tkin cut of TR photon in current mat.
148 
151 
153 
154 static G4int fSympsonNumber; // Accuracy of Sympson integration
155 
156 static G4double fTheMinEnergyTR; // static min TR energy
157 static G4double fTheMaxEnergyTR; // static max TR energy
158  G4double fMinEnergyTR; // min TR energy in material
159  G4double fMaxEnergyTR; // max TR energy in material
160 static G4double fTheMaxAngle; // max theta of TR quanta
161 static G4double fTheMinAngle; // max theta of TR quanta
162  G4double fMaxThetaTR; // max theta of TR quanta
163 static G4int fBinTR; // number of bins in TR vectors
164 
165 static G4double fMinProtonTkin; // min Tkin of proton in tables
166 static G4double fMaxProtonTkin; // max Tkin of proton in tables
167 static G4int fTotBin; // number of bins in log scale
168  G4double fGamma; // current Lorentz factor
169 
170 static G4double fPlasmaCof; // physical consts for plasma energy
172 
173 G4double fSigma1; // plasma energy Sq of matter1
174 G4double fSigma2; // plasma energy Sq of matter2
175 
176 private:
177  // Operators
178 
179  G4ForwardXrayTR(const G4ForwardXrayTR& right) = delete;
180 
181  G4ForwardXrayTR& operator=(const G4ForwardXrayTR& right) = delete;
182 
183  // G4int operator==(const G4ForwardXrayTR& right)const;
184  // G4int operator!=(const G4ForwardXrayTR& right)const;
185 
186 }; // end of G4ForwardXrayTR class ---------------------------
187 
188 #endif // G4FORWARDXRAYTR_H
189 
190 
191 
static G4double fMaxProtonTkin
G4double AngleDensity(G4double energy, G4double varAngle) const
G4PhysicsTable * GetAngleDistrTable()
G4double condition(const G4ErrorSymMatrix &m)
G4double AngleInterval(G4double energy, G4double varAngle1, G4double varAngle2) const
static G4int GetBinTR()
static G4int fBinTR
static G4int GetSympsonNumber()
static G4double fTheMaxEnergyTR
G4double GetEnergyTR(G4int iMat, G4int jMat, G4int iTkin) const
static G4double fMinProtonTkin
G4ForwardXrayTR(const G4String &matName1, const G4String &matName2, const G4String &processName="XrayTR")
G4double SpectralDensity(G4double energy, G4double x) const
G4double AngleSum(G4double varAngle1, G4double varAngle2) const
G4double GetThetaTR(G4int iMat, G4int jMat, G4int iTkin) const
tuple x
Definition: test.py:50
int G4int
Definition: G4Types.hh:78
G4PhysicsTable * fAngleDistrTable
static G4double fPlasmaCof
static G4double fTheMinAngle
static G4int fSympsonNumber
virtual ~G4ForwardXrayTR()
static G4int fTotBin
Definition: G4Step.hh:76
G4PhysicsTable * fEnergyDistrTable
G4double SpectralAngleTRdensity(G4double energy, G4double varAngle) const override
G4PhysicsTable * GetEnergyDistrTable()
G4ParticleDefinition * fPtrGamma
G4double energy(const ThreeVector &p, const G4double m)
G4double EnergySum(G4double energy1, G4double energy2) const
static G4int GetTotBin()
static G4double fCofTR
static G4double fTheMinEnergyTR
G4double GetMeanFreePath(const G4Track &, G4double, G4ForceCondition *condition) override
static G4double fTheMaxAngle
double G4double
Definition: G4Types.hh:76
G4ForceCondition
G4double EnergyInterval(G4double energy1, G4double energy2, G4double varAngle) const
G4PhysicsLogVector * fProtonEnergyVector
const std::vector< G4double > * fGammaCutInKineticEnergy
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep) override
static G4double GetMaxProtonTkin()
static G4double GetMinProtonTkin()