Geant4  10.03.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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 // $Id: G4VXTRenergyLoss.hh 97385 2016-06-02 09:59:53Z gcosmo $
28 //
29 //
31 //
32 // base class for 'fast' parametrisation model describing X-ray transition
33 // created in some G4Envelope. Anglur distribuiton is very rough !!! (see DoIt
34 // method
35 //
36 // History:
37 // 06.10.05 V. Grichine first step to discrete process
38 // 15.01.02 V. Grichine first version
39 // 28.07.05, P.Gumplinger add G4ProcessType to constructor
40 // 28.09.07, V.Ivanchenko general cleanup without change of algorithms
41 //
42 
43 #ifndef G4VXTRenergyLoss_h
44 #define G4VXTRenergyLoss_h 1
45 
46 #include <complex>
47 #include "globals.hh"
48 #include "Randomize.hh"
49 
50 #include "G4LogicalVolume.hh"
51 
52 #include "G4PhysicsTable.hh"
53 #include "G4PhysicsLogVector.hh"
54 #include "G4Gamma.hh"
55 #include "G4ThreeVector.hh"
56 #include "G4ParticleMomentum.hh"
57 #include "G4Step.hh"
58 #include "G4Track.hh"
59 #include "G4VContinuousProcess.hh"
60 #include "G4VDiscreteProcess.hh"
61 #include "G4DynamicParticle.hh"
62 #include "G4Material.hh"
63 #include "G4PhysicsTable.hh"
66 #include "G4Integrator.hh"
67 #include "G4ParticleChange.hh"
68 
69 class G4SandiaTable;
70 class G4VParticleChange;
73 
74 class G4VXTRenergyLoss : public G4VDiscreteProcess // G4VContinuousProcess
75 {
76 public:
77 
78  explicit G4VXTRenergyLoss (G4LogicalVolume *anEnvelope,G4Material*,
80  const G4String & processName = "XTRenergyLoss",
82  virtual ~G4VXTRenergyLoss ();
83 
84  // These virtual has to be implemented in inherited particular TR radiators
85 
87  G4double varAngle );
88 
89  virtual G4bool IsApplicable(const G4ParticleDefinition&) override;
90 
91  virtual G4VParticleChange* PostStepDoIt(const G4Track& aTrack,
92  const G4Step& aStep) override;
93 
94  virtual G4double GetMeanFreePath(const G4Track& aTrack,
95  G4double previousStepSize,
96  G4ForceCondition* condition) override;
97 
98  virtual void BuildPhysicsTable(const G4ParticleDefinition&) override;
99  void BuildEnergyTable() ;
100  void BuildAngleForEnergyBank() ;
101 
102  void BuildTable(){} ;
103  void BuildAngleTable() ;
104  void BuildGlobalAngleTable() ;
105 
107  G4double gamma,
108  G4double varAngle ) ;
109 
111 
112  virtual G4double SpectralXTRdEdx(G4double energy) ;
113 
115 
116  G4double AngleXTRdEdx(G4double varAngle) ;
117 
118 
120 
122  G4double gamma,
123  G4double varAngle ) const ;
124 
125 
126  // for photon energy distribution tables
127 
130 
131  // for photon angle distribution tables
132 
135 
136  void GetNumberOfPhotons() ;
137 
138  // Auxiliary functions for plate/gas material parameters
139 
144  void GetPlateZmuProduct() ;
146 
149  void ComputeGasPhotoAbsCof();
151  void GetGasZmuProduct();
153 
157 
158  G4double GetXTRrandomEnergy( G4double scaledTkin, G4int iTkin );
159  G4double GetXTRenergy( G4int iPlace, G4double position, G4int iTransfer );
160 
161  G4double GetRandomAngle( G4double energyXTR, G4int iTkin );
163 
164  G4double GetGamma() {return fGamma;};
165  G4double GetEnergy() {return fEnergy;};
167 
168  void SetGamma(G4double gamma) {fGamma = gamma;};
169  void SetEnergy(G4double energy) {fEnergy = energy;};
170  void SetVarAngle(G4double varAngle){fVarAngle = varAngle;};
171  void SetAngleRadDistr(G4bool pAngleRadDistr){fAngleRadDistr=pAngleRadDistr;};
172  void SetCompton(G4bool pC){fCompton=pC;};
173 
175  G4int GetTotBin(){return fTotBin;};
177 
178 protected:
179 
180  G4ParticleDefinition* fPtrGamma ; // pointer to TR photon
181 
182  G4double* fGammaCutInKineticEnergy ; // TR photon cut in energy array
183 
184  G4double fGammaTkinCut ; // Tkin cut of TR photon in current mat.
188 
191 
192  G4double fTheMinEnergyTR; // min TR energy
193  G4double fTheMaxEnergyTR; // max TR energy
194  G4double fMinEnergyTR; // min TR energy in material
195  G4double fMaxEnergyTR; // max TR energy in material
196  G4double fTheMaxAngle; // max theta of TR quanta
197  G4double fTheMinAngle; // max theta of TR quanta
198  G4double fMaxThetaTR; // max theta of TR quanta
199  G4int fBinTR; // number of bins in TR vectors
200 
201  G4double fMinProtonTkin; // min Tkin of proton in tables
202  G4double fMaxProtonTkin; // max Tkin of proton in tables
203  G4int fTotBin; // number of bins in log scale
204  G4double fGamma; // current Lorentz factor
205  G4double fEnergy; // energy and
206  G4double fVarAngle; // angle squared
208 
209  G4double fPlasmaCof ; // physical consts for plasma energy
211 
216  G4double fSigma2; // plasma energy Sq of matter1/2
217 
221 
227 
229 
231 
233 
235  std::vector<G4PhysicsTable*> fAngleBank;
236 
237 private:
238 
239  // copy constructor and hide assignment operator
241  G4VXTRenergyLoss & operator=(const G4VXTRenergyLoss &right) = delete;
242 
243 };
244 
245 #endif
G4double condition(const G4ErrorSymMatrix &m)
G4double GetXTRrandomEnergy(G4double scaledTkin, G4int iTkin)
G4PhysicsTable * fEnergyDistrTable
G4PhysicsLogVector * GetProtonVector()
void SetEnergy(G4double energy)
G4double GetGasFormationZone(G4double, G4double, G4double)
G4double GetPlateLinearPhotoAbs(G4double)
void SetCompton(G4bool pC)
G4double XTRNSpectralDensity(G4double energy)
void SetGamma(G4double gamma)
G4LogicalVolume * fEnvelope
virtual G4double GetMeanFreePath(const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition) override
G4double GetPlateCompton(G4double)
G4double GetGasLinearPhotoAbs(G4double)
void SetAngleRadDistr(G4bool pAngleRadDistr)
int G4int
Definition: G4Types.hh:78
G4double GetPlateFormationZone(G4double, G4double, G4double)
G4double GetComptonPerAtom(G4double, G4double)
G4complex GetPlateComplexFZ(G4double, G4double, G4double)
std::complex< G4double > G4complex
Definition: G4Types.hh:81
virtual G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep) override
G4PhysicsLogVector * fProtonEnergyVector
G4double GetXTRenergy(G4int iPlace, G4double position, G4int iTransfer)
bool G4bool
Definition: G4Types.hh:79
G4double OneBoundaryXTRNdensity(G4double energy, G4double gamma, G4double varAngle) const
Definition: G4Step.hh:76
G4double XTRNSpectralAngleDensity(G4double varAngle)
const G4int n
G4PhysicsFreeVector * GetAngleVector(G4double energy, G4int n)
G4SandiaTable * fGasPhotoAbsCof
virtual void BuildPhysicsTable(const G4ParticleDefinition &) override
void SetVarAngle(G4double varAngle)
G4double GetAngleXTR(G4int iTR, G4double position, G4int iAngle)
G4double energy(const ThreeVector &p, const G4double m)
virtual ~G4VXTRenergyLoss()
G4double * fGammaCutInKineticEnergy
G4double SpectralAngleXTRdEdx(G4double varAngle)
virtual G4double SpectralXTRdEdx(G4double energy)
G4PhysicsLogVector * fXTREnergyVector
G4VXTRenergyLoss(G4LogicalVolume *anEnvelope, G4Material *, G4Material *, G4double, G4double, G4int, const G4String &processName="XTRenergyLoss", G4ProcessType type=fElectromagnetic)
std::vector< G4PhysicsTable * > fAngleBank
G4double GetRandomAngle(G4double energyXTR, G4int iTkin)
G4complex OneInterfaceXTRdEdx(G4double energy, G4double gamma, G4double varAngle)
G4PhysicsTable * fAngleForEnergyTable
G4double GetGasCompton(G4double)
G4PhysicsTable * fAngleDistrTable
G4SandiaTable * fPlatePhotoAbsCof
G4double XTRNAngleDensity(G4double varAngle)
G4complex GetGasComplexFZ(G4double, G4double, G4double)
double G4double
Definition: G4Types.hh:76
G4ParticleDefinition * fPtrGamma
G4ParticleChange fParticleChange
G4ForceCondition
G4double XTRNAngleSpectralDensity(G4double energy)
virtual G4bool IsApplicable(const G4ParticleDefinition &) override
G4double AngleSpectralXTRdEdx(G4double energy)
virtual G4double GetStackFactor(G4double energy, G4double gamma, G4double varAngle)
G4ProcessType
G4double AngleXTRdEdx(G4double varAngle)