source: trunk/source/processes/electromagnetic/xrays/include/G4VXTRenergyLoss.hh @ 1196

Last change on this file since 1196 was 1196, checked in by garnier, 15 years ago

update CVS release candidate geant4.9.3.01

File size: 8.6 KB
Line 
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,v 1.24 2007/09/29 17:49:34 vnivanch Exp $
28// GEANT4 tag $Name: geant4-09-03-cand-01 $
29//
30//
31///////////////////////////////////////////////////////////////////////////
32//
33// base class for 'fast' parametrisation model describing X-ray transition
34// created in some G4Envelope. Anglur distribuiton is very rough !!! (see DoIt
35// method
36//
37// History:
38// 06.10.05 V. Grichine first step to discrete process
39// 15.01.02 V. Grichine first version
40// 28.07.05, P.Gumplinger add G4ProcessType to constructor
41// 28.09.07, V.Ivanchenko general cleanup without change of algorithms
42//
43
44#ifndef G4VXTRenergyLoss_h
45#define G4VXTRenergyLoss_h 1
46
47#include <complex>
48#include "globals.hh"
49#include "Randomize.hh"
50
51#include "G4LogicalVolume.hh"
52
53#include "G4PhysicsTable.hh"
54#include "G4PhysicsLogVector.hh"
55#include "G4Gamma.hh"
56#include "G4ThreeVector.hh"
57#include "G4ParticleMomentum.hh"
58#include "G4Step.hh"
59#include "G4Track.hh"
60#include "G4VContinuousProcess.hh"
61#include "G4VDiscreteProcess.hh"
62#include "G4DynamicParticle.hh"
63#include "G4Material.hh"
64#include "G4PhysicsTable.hh"
65#include "G4MaterialPropertiesTable.hh"
66#include "G4PhysicsOrderedFreeVector.hh"
67#include "G4Integrator.hh"
68#include "G4ParticleChange.hh"
69
70class G4SandiaTable;
71class G4VParticleChange;
72class G4PhysicsFreeVector;
73
74class G4VXTRenergyLoss : public G4VDiscreteProcess  // G4VContinuousProcess
75{
76public:
77
78  G4VXTRenergyLoss (G4LogicalVolume *anEnvelope,G4Material*,G4Material*,
79                    G4double,G4double,G4int,
80                    const G4String & processName = "XTRenergyLoss",
81                    G4ProcessType type = fElectromagnetic);
82  virtual  ~G4VXTRenergyLoss ();
83
84  // These virtual has to be implemented in inherited particular TR radiators
85 
86  virtual  G4double GetStackFactor( G4double energy, G4double gamma,
87                                                     G4double varAngle );
88
89  G4bool IsApplicable(const G4ParticleDefinition&);
90
91  G4VParticleChange* PostStepDoIt(const G4Track& aTrack, 
92                                   const G4Step&  aStep);
93
94  G4double GetMeanFreePath(const G4Track& aTrack,
95                           G4double previousStepSize,
96                           G4ForceCondition* condition);
97
98  void BuildPhysicsTable(const G4ParticleDefinition&);
99  void BuildTable() ;
100  void BuildEnergyTable() ;
101  void BuildAngleTable() ;
102  void BuildGlobalAngleTable() ;
103
104  G4complex OneInterfaceXTRdEdx( G4double energy, 
105                                G4double gamma,
106                                G4double varAngle ) ;
107
108  G4double SpectralAngleXTRdEdx(G4double varAngle) ;
109
110  virtual  G4double SpectralXTRdEdx(G4double energy) ;
111
112  G4double AngleSpectralXTRdEdx(G4double energy) ;
113
114  G4double AngleXTRdEdx(G4double varAngle) ;
115
116
117  /////////////////////////////////////////////////////////////
118
119  G4double OneBoundaryXTRNdensity( G4double energy,
120                                   G4double gamma,
121                                   G4double varAngle ) const ;
122
123
124  // for photon energy distribution tables
125
126  G4double XTRNSpectralAngleDensity(G4double varAngle) ;
127  G4double XTRNSpectralDensity(G4double energy) ;
128 
129  // for photon angle distribution tables
130
131  G4double XTRNAngleSpectralDensity(G4double energy) ;
132  G4double XTRNAngleDensity(G4double varAngle) ;
133
134  void GetNumberOfPhotons() ; 
135
136  // Auxiliary functions for plate/gas material parameters
137
138  G4double  GetPlateFormationZone(G4double,G4double,G4double);
139  G4complex GetPlateComplexFZ(G4double,G4double,G4double);
140  void      ComputePlatePhotoAbsCof();
141  G4double  GetPlateLinearPhotoAbs(G4double);
142  void      GetPlateZmuProduct() ;
143  G4double  GetPlateZmuProduct(G4double,G4double,G4double);
144
145  G4double  GetGasFormationZone(G4double,G4double,G4double);
146  G4complex GetGasComplexFZ(G4double,G4double,G4double);
147  void      ComputeGasPhotoAbsCof();
148  G4double  GetGasLinearPhotoAbs(G4double);
149  void      GetGasZmuProduct();
150  G4double  GetGasZmuProduct(G4double,G4double,G4double);
151
152  G4double  GetPlateCompton(G4double);
153  G4double  GetGasCompton(G4double);
154  G4double  GetComptonPerAtom(G4double,G4double);
155
156  G4double  GetXTRrandomEnergy( G4double scaledTkin, G4int iTkin );
157  G4double  GetXTRenergy( G4int iPlace, G4double position, G4int iTransfer  );
158
159  G4double  GetRandomAngle( G4double energyXTR, G4int iTkin );
160  G4double  GetAngleXTR(G4int iTR,G4double position,G4int iAngle);
161
162  G4double  GetGamma()   {return fGamma;}; 
163  G4double  GetEnergy()  {return fEnergy;};               
164  G4double  GetVarAngle(){return fVarAngle;};
165               
166  void SetGamma(G4double gamma)      {fGamma    = gamma;}; 
167  void SetEnergy(G4double energy)    {fEnergy   = energy;};               
168  void SetVarAngle(G4double varAngle){fVarAngle = varAngle;};               
169  void SetAngleRadDistr(G4bool pAngleRadDistr){fAngleRadDistr=pAngleRadDistr;};               
170  void SetCompton(G4bool pC){fCompton=pC;};               
171
172  G4PhysicsLogVector* GetProtonVector(){ return fProtonEnergyVector;};
173  G4int GetTotBin(){return fTotBin;};           
174  G4PhysicsFreeVector* GetAngleVector(G4double energy, G4int n);
175
176protected:
177
178  G4ParticleDefinition* fPtrGamma ;    // pointer to TR photon
179
180  G4double* fGammaCutInKineticEnergy ; // TR photon cut in energy array
181
182  G4double         fGammaTkinCut ;     // Tkin cut of TR photon in current mat.
183  G4LogicalVolume* fEnvelope ;
184  G4PhysicsTable*  fAngleDistrTable ;
185  G4PhysicsTable*  fEnergyDistrTable ;
186
187  G4PhysicsLogVector* fProtonEnergyVector ;
188  G4PhysicsLogVector* fXTREnergyVector ;
189
190  G4double fTheMinEnergyTR;            //   min TR energy
191  G4double fTheMaxEnergyTR;            //   max TR energy
192  G4double fMinEnergyTR;               //  min TR energy in material
193  G4double fMaxEnergyTR;               //  max TR energy in material
194  G4double fTheMaxAngle;               //  max theta of TR quanta
195  G4double fTheMinAngle;               //  max theta of TR quanta
196  G4double fMaxThetaTR;                //  max theta of TR quanta
197  G4int    fBinTR;                     //  number of bins in TR vectors
198
199  G4double fMinProtonTkin;             // min Tkin of proton in tables
200  G4double fMaxProtonTkin;             // max Tkin of proton in tables
201  G4int    fTotBin;                    // number of bins in log scale
202  G4double fGamma;                     // current Lorentz factor
203  G4double fEnergy;                    // energy and
204  G4double fVarAngle;                  // angle squared
205  G4double fLambda;
206
207  G4double fPlasmaCof ;                // physical consts for plasma energy
208  G4double fCofTR ; 
209
210  G4bool   fExitFlux;
211  G4bool   fAngleRadDistr;
212  G4bool   fCompton;
213  G4double fSigma1; 
214  G4double fSigma2;                    // plasma energy Sq of matter1/2
215
216  G4int    fMatIndex1;
217  G4int    fMatIndex2;
218  G4int    fPlateNumber;
219
220  G4double fTotalDist;
221  G4double fPlateThick;
222  G4double fGasThick;     
223  G4double fAlphaPlate;
224  G4double fAlphaGas ;
225
226  G4SandiaTable* fPlatePhotoAbsCof;
227 
228  G4SandiaTable* fGasPhotoAbsCof;
229
230  G4ParticleChange fParticleChange;
231
232  G4PhysicsTable*                    fAngleForEnergyTable;
233  std::vector<G4PhysicsTable*>       fAngleBank;
234
235};
236
237#endif
Note: See TracBrowser for help on using the repository browser.