source: trunk/source/processes/electromagnetic/xrays/src/G4TransitionRadiation.cc @ 961

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

update processes

File size: 8.4 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: G4TransitionRadiation.cc,v 1.7 2006/06/29 19:56:19 gunter Exp $
28// GEANT4 tag $Name: geant4-09-02-ref-02 $
29//
30// G4TransitionRadiation class -- implementation file
31
32// GEANT 4 class implementation file --- Copyright CERN 1995
33// CERN Geneva Switzerland
34
35// For information related to this code, please, contact
36// CERN, CN Division, ASD Group
37// History:
38// 1st version 11.09.97 V. Grichine (Vladimir.Grichine@cern.ch )
39// 2nd version 16.12.97 V. Grichine
40// 3rd version 28.07.05, P.Gumplinger add G4ProcessType to constructor
41
42
43#include <cmath>
44// #include "G4ios.hh"
45// #include <fstream.h>
46// #include <stdlib.h>
47
48#include "G4TransitionRadiation.hh"
49#include "G4Material.hh"
50
51// Init gamma array
52
53
54// Local constants
55
56const G4int   G4TransitionRadiation::fSympsonNumber = 100 ;
57const G4int   G4TransitionRadiation::fGammaNumber = 15 ;
58const G4int   G4TransitionRadiation::fPointNumber = 100 ;
59
60using namespace std;
61
62///////////////////////////////////////////////////////////////////////
63//
64// Constructor for selected couple of materials
65//
66
67G4TransitionRadiation::
68G4TransitionRadiation( const G4String& processName, G4ProcessType type )
69  : G4VDiscreteProcess(processName, type)
70{
71  //  fMatIndex1 = pMat1->GetIndex() ;
72  //  fMatIndex2 = pMat2->GetIndex() ;
73}
74
75//////////////////////////////////////////////////////////////////////
76//
77// Destructor
78//
79
80G4TransitionRadiation::~G4TransitionRadiation()
81{
82        ;
83}
84
85
86///////////////////////////////////////////////////////////////////
87//
88// Sympson integral of TR spectral-angle density over energy between
89// the limits energy 1 and energy2 at fixed varAngle = 1 - cos(Theta)
90
91G4double
92G4TransitionRadiation::IntegralOverEnergy( G4double energy1,
93                                           G4double energy2,
94                                           G4double varAngle     )  const
95{
96  G4int i ;
97  G4double h , sumEven = 0.0 , sumOdd = 0.0 ;
98  h = 0.5*(energy2 - energy1)/fSympsonNumber ;
99  for(i=1;i<fSympsonNumber;i++)
100  {
101    sumEven += SpectralAngleTRdensity(energy1 + 2*i*h,varAngle)  ;
102    sumOdd  += SpectralAngleTRdensity(energy1 + (2*i - 1)*h,varAngle) ;
103  }
104  sumOdd += SpectralAngleTRdensity(energy1 + (2*fSympsonNumber - 1)*h,varAngle) ;
105  return h*(  SpectralAngleTRdensity(energy1,varAngle)
106            + SpectralAngleTRdensity(energy2,varAngle)
107            + 4.0*sumOdd + 2.0*sumEven    )/3.0 ;
108}
109
110
111
112///////////////////////////////////////////////////////////////////
113//
114// Sympson integral of TR spectral-angle density over energy between
115// the limits varAngle1 and varAngle2 at fixed energy
116
117G4double
118G4TransitionRadiation::IntegralOverAngle( G4double energy,
119                                          G4double varAngle1,
120                                          G4double varAngle2     ) const
121{
122  G4int i ;
123  G4double h , sumEven = 0.0 , sumOdd = 0.0 ;
124  h = 0.5*(varAngle2 - varAngle1)/fSympsonNumber ;
125  for(i=1;i<fSympsonNumber;i++)
126  {
127    sumEven += SpectralAngleTRdensity(energy,varAngle1 + 2*i*h)  ;
128    sumOdd  += SpectralAngleTRdensity(energy,varAngle1 + (2*i - 1)*h) ;
129  }
130  sumOdd += SpectralAngleTRdensity(energy,varAngle1 + (2*fSympsonNumber - 1)*h) ;
131
132  return h*(  SpectralAngleTRdensity(energy,varAngle1)
133            + SpectralAngleTRdensity(energy,varAngle2)
134            + 4.0*sumOdd + 2.0*sumEven    )/3.0 ;
135}
136
137///////////////////////////////////////////////////////////////////
138//
139// The number of transition radiation photons generated in the
140// angle interval between varAngle1 and varAngle2
141//
142
143G4double G4TransitionRadiation::
144AngleIntegralDistribution( G4double varAngle1,
145                           G4double varAngle2     )   const
146{
147  G4int i ;
148  G4double h , sumEven = 0.0 , sumOdd = 0.0 ;
149  h = 0.5*(varAngle2 - varAngle1)/fSympsonNumber ;
150  for(i=1;i<fSympsonNumber;i++)
151  {
152   sumEven += IntegralOverEnergy(fMinEnergy,
153                                 fMinEnergy +0.3*(fMaxEnergy-fMinEnergy),
154                                 varAngle1 + 2*i*h)
155            + IntegralOverEnergy(fMinEnergy + 0.3*(fMaxEnergy - fMinEnergy),
156                                 fMaxEnergy,
157                                 varAngle1 + 2*i*h);
158   sumOdd  += IntegralOverEnergy(fMinEnergy,
159                                 fMinEnergy + 0.3*(fMaxEnergy - fMinEnergy),
160                                 varAngle1 + (2*i - 1)*h)
161            + IntegralOverEnergy(fMinEnergy + 0.3*(fMaxEnergy - fMinEnergy),
162                                 fMaxEnergy,
163                                 varAngle1 + (2*i - 1)*h) ;
164  }
165  sumOdd += IntegralOverEnergy(fMinEnergy,
166                               fMinEnergy + 0.3*(fMaxEnergy - fMinEnergy),
167                               varAngle1 + (2*fSympsonNumber - 1)*h)
168          + IntegralOverEnergy(fMinEnergy + 0.3*(fMaxEnergy - fMinEnergy),
169                               fMaxEnergy,
170                               varAngle1 + (2*fSympsonNumber - 1)*h) ;
171
172  return h*(IntegralOverEnergy(fMinEnergy,
173                               fMinEnergy + 0.3*(fMaxEnergy - fMinEnergy),
174                               varAngle1)
175          + IntegralOverEnergy(fMinEnergy + 0.3*(fMaxEnergy - fMinEnergy),
176                               fMaxEnergy,
177                               varAngle1)
178          + IntegralOverEnergy(fMinEnergy,
179                               fMinEnergy + 0.3*(fMaxEnergy - fMinEnergy),
180                               varAngle2)
181          + IntegralOverEnergy(fMinEnergy + 0.3*(fMaxEnergy - fMinEnergy),
182                               fMaxEnergy,
183                               varAngle2)
184            + 4.0*sumOdd + 2.0*sumEven    )/3.0 ;
185}
186
187///////////////////////////////////////////////////////////////////
188//
189// The number of transition radiation photons, generated in the
190// energy interval between energy1 and energy2
191//
192
193G4double G4TransitionRadiation::
194EnergyIntegralDistribution( G4double energy1,
195                            G4double energy2     )  const
196{
197  G4int i ;
198  G4double h , sumEven = 0.0 , sumOdd = 0.0 ;
199  h = 0.5*(energy2 - energy1)/fSympsonNumber ;
200  for(i=1;i<fSympsonNumber;i++)
201  {
202   sumEven += IntegralOverAngle(energy1 + 2*i*h,0.0,0.01*fMaxTheta )
203            + IntegralOverAngle(energy1 + 2*i*h,0.01*fMaxTheta,fMaxTheta);
204   sumOdd  += IntegralOverAngle(energy1 + (2*i - 1)*h,0.0,0.01*fMaxTheta)
205            + IntegralOverAngle(energy1 + (2*i - 1)*h,0.01*fMaxTheta,fMaxTheta) ;
206  }
207  sumOdd += IntegralOverAngle(energy1 + (2*fSympsonNumber - 1)*h,
208                              0.0,0.01*fMaxTheta)
209          + IntegralOverAngle(energy1 + (2*fSympsonNumber - 1)*h,
210                              0.01*fMaxTheta,fMaxTheta) ;
211
212  return h*(IntegralOverAngle(energy1,0.0,0.01*fMaxTheta)
213          + IntegralOverAngle(energy1,0.01*fMaxTheta,fMaxTheta)
214          + IntegralOverAngle(energy2,0.0,0.01*fMaxTheta)
215          + IntegralOverAngle(energy2,0.01*fMaxTheta,fMaxTheta)
216            + 4.0*sumOdd + 2.0*sumEven    )/3.0 ;
217}
218
219
220
221
222// end of G4TransitionRadiation implementation file --------------------------
Note: See TracBrowser for help on using the repository browser.