Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4TransitionRadiation.cc
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 // $Id: G4TransitionRadiation.cc 68037 2013-03-13 14:15:08Z gcosmo $
27 //
28 // G4TransitionRadiation class -- implementation file
29 
30 // GEANT 4 class implementation file --- Copyright CERN 1995
31 // CERN Geneva Switzerland
32 
33 // For information related to this code, please, contact
34 // CERN, CN Division, ASD Group
35 // History:
36 // 1st version 11.09.97 V. Grichine (Vladimir.Grichine@cern.ch )
37 // 2nd version 16.12.97 V. Grichine
38 // 3rd version 28.07.05, P.Gumplinger add G4ProcessType to constructor
39 
40 
41 #include <cmath>
42 
43 #include "G4TransitionRadiation.hh"
44 #include "G4Material.hh"
45 #include "G4EmProcessSubType.hh"
46 
47 // Local constants
48 
52 
53 
55 //
56 // Constructor for selected couple of materials
57 //
58 
60 G4TransitionRadiation( const G4String& processName, G4ProcessType type )
61  : G4VDiscreteProcess(processName, type)
62 {
64  fMatIndex1 = fMatIndex2 = 0;
65 
67 }
68 
70 //
71 // Destructor
72 //
73 
75 {}
76 
77 G4bool
79 {
80  return ( aParticleType.GetPDGCharge() != 0.0 );
81 }
82 
84  G4double,
86 {
87  *condition = Forced;
88  return DBL_MAX; // so TR doesn't limit mean free path
89 }
90 
92  const G4Step&)
93 {
95  return &aParticleChange;
96 }
97 
99 //
100 // Sympson integral of TR spectral-angle density over energy between
101 // the limits energy 1 and energy2 at fixed varAngle = 1 - std::cos(Theta)
102 
103 G4double
105  G4double energy2,
106  G4double varAngle ) const
107 {
108  G4int i ;
109  G4double h , sumEven = 0.0 , sumOdd = 0.0 ;
110  h = 0.5*(energy2 - energy1)/fSympsonNumber ;
111  for(i=1;i<fSympsonNumber;i++)
112  {
113  sumEven += SpectralAngleTRdensity(energy1 + 2*i*h,varAngle) ;
114  sumOdd += SpectralAngleTRdensity(energy1 + (2*i - 1)*h,varAngle) ;
115  }
116  sumOdd += SpectralAngleTRdensity(energy1 + (2*fSympsonNumber - 1)*h,varAngle) ;
117  return h*( SpectralAngleTRdensity(energy1,varAngle)
118  + SpectralAngleTRdensity(energy2,varAngle)
119  + 4.0*sumOdd + 2.0*sumEven )/3.0 ;
120 }
121 
122 
123 
125 //
126 // Sympson integral of TR spectral-angle density over energy between
127 // the limits varAngle1 and varAngle2 at fixed energy
128 
129 G4double
131  G4double varAngle1,
132  G4double varAngle2 ) const
133 {
134  G4int i ;
135  G4double h , sumEven = 0.0 , sumOdd = 0.0 ;
136  h = 0.5*(varAngle2 - varAngle1)/fSympsonNumber ;
137  for(i=1;i<fSympsonNumber;i++)
138  {
139  sumEven += SpectralAngleTRdensity(energy,varAngle1 + 2*i*h) ;
140  sumOdd += SpectralAngleTRdensity(energy,varAngle1 + (2*i - 1)*h) ;
141  }
142  sumOdd += SpectralAngleTRdensity(energy,varAngle1 + (2*fSympsonNumber - 1)*h) ;
143 
144  return h*( SpectralAngleTRdensity(energy,varAngle1)
145  + SpectralAngleTRdensity(energy,varAngle2)
146  + 4.0*sumOdd + 2.0*sumEven )/3.0 ;
147 }
148 
150 //
151 // The number of transition radiation photons generated in the
152 // angle interval between varAngle1 and varAngle2
153 //
154 
157  G4double varAngle2 ) const
158 {
159  G4int i ;
160  G4double h , sumEven = 0.0 , sumOdd = 0.0 ;
161  h = 0.5*(varAngle2 - varAngle1)/fSympsonNumber ;
162  for(i=1;i<fSympsonNumber;i++)
163  {
164  sumEven += IntegralOverEnergy(fMinEnergy,
166  varAngle1 + 2*i*h)
168  fMaxEnergy,
169  varAngle1 + 2*i*h);
170  sumOdd += IntegralOverEnergy(fMinEnergy,
171  fMinEnergy + 0.3*(fMaxEnergy - fMinEnergy),
172  varAngle1 + (2*i - 1)*h)
174  fMaxEnergy,
175  varAngle1 + (2*i - 1)*h) ;
176  }
177  sumOdd += IntegralOverEnergy(fMinEnergy,
178  fMinEnergy + 0.3*(fMaxEnergy - fMinEnergy),
179  varAngle1 + (2*fSympsonNumber - 1)*h)
181  fMaxEnergy,
182  varAngle1 + (2*fSympsonNumber - 1)*h) ;
183 
184  return h*(IntegralOverEnergy(fMinEnergy,
185  fMinEnergy + 0.3*(fMaxEnergy - fMinEnergy),
186  varAngle1)
188  fMaxEnergy,
189  varAngle1)
191  fMinEnergy + 0.3*(fMaxEnergy - fMinEnergy),
192  varAngle2)
194  fMaxEnergy,
195  varAngle2)
196  + 4.0*sumOdd + 2.0*sumEven )/3.0 ;
197 }
198 
200 //
201 // The number of transition radiation photons, generated in the
202 // energy interval between energy1 and energy2
203 //
204 
207  G4double energy2 ) const
208 {
209  G4int i ;
210  G4double h , sumEven = 0.0 , sumOdd = 0.0 ;
211  h = 0.5*(energy2 - energy1)/fSympsonNumber ;
212  for(i=1;i<fSympsonNumber;i++)
213  {
214  sumEven += IntegralOverAngle(energy1 + 2*i*h,0.0,0.01*fMaxTheta )
215  + IntegralOverAngle(energy1 + 2*i*h,0.01*fMaxTheta,fMaxTheta);
216  sumOdd += IntegralOverAngle(energy1 + (2*i - 1)*h,0.0,0.01*fMaxTheta)
217  + IntegralOverAngle(energy1 + (2*i - 1)*h,0.01*fMaxTheta,fMaxTheta) ;
218  }
219  sumOdd += IntegralOverAngle(energy1 + (2*fSympsonNumber - 1)*h,
220  0.0,0.01*fMaxTheta)
221  + IntegralOverAngle(energy1 + (2*fSympsonNumber - 1)*h,
222  0.01*fMaxTheta,fMaxTheta) ;
223 
224  return h*(IntegralOverAngle(energy1,0.0,0.01*fMaxTheta)
225  + IntegralOverAngle(energy1,0.01*fMaxTheta,fMaxTheta)
226  + IntegralOverAngle(energy2,0.0,0.01*fMaxTheta)
227  + IntegralOverAngle(energy2,0.01*fMaxTheta,fMaxTheta)
228  + 4.0*sumOdd + 2.0*sumEven )/3.0 ;
229 }
230 
231 
232 
233 
234 // end of G4TransitionRadiation implementation file --------------------------
G4double condition(const G4ErrorSymMatrix &m)
G4double AngleIntegralDistribution(G4double varAngle1, G4double varAngle2) const
static const G4int fSympsonNumber
G4TransitionRadiation(const G4String &processName="TR", G4ProcessType type=fElectromagnetic)
void ClearNumberOfInteractionLengthLeft()
Definition: G4VProcess.hh:447
int G4int
Definition: G4Types.hh:78
virtual G4double SpectralAngleTRdensity(G4double energy, G4double varAngle) const =0
bool G4bool
Definition: G4Types.hh:79
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
Definition: G4Step.hh:76
virtual G4double GetMeanFreePath(const G4Track &, G4double, G4ForceCondition *condition) override
G4bool IsApplicable(const G4ParticleDefinition &aParticleType) override
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &) override
G4double energy(const ThreeVector &p, const G4double m)
G4double IntegralOverAngle(G4double energy, G4double varAngle1, G4double varAngle2) const
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:289
G4double EnergyIntegralDistribution(G4double energy1, G4double energy2) const
double G4double
Definition: G4Types.hh:76
G4ForceCondition
G4double IntegralOverEnergy(G4double energy1, G4double energy2, G4double varAngle) const
G4double GetPDGCharge() const
#define DBL_MAX
Definition: templates.hh:83
static const G4int fGammaNumber
static const G4int fPointNumber
G4ProcessType