[819] | 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 | // |
---|
[1340] | 27 | // $Id: G4TransitionRadiation.cc,v 1.10 2010/10/14 18:38:21 vnivanch Exp $ |
---|
| 28 | // GEANT4 tag $Name: xrays-V09-03-05 $ |
---|
[819] | 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 | |
---|
| 45 | #include "G4TransitionRadiation.hh" |
---|
| 46 | #include "G4Material.hh" |
---|
[1340] | 47 | #include "G4EmProcessSubType.hh" |
---|
[819] | 48 | |
---|
| 49 | // Local constants |
---|
| 50 | |
---|
| 51 | const G4int G4TransitionRadiation::fSympsonNumber = 100 ; |
---|
| 52 | const G4int G4TransitionRadiation::fGammaNumber = 15 ; |
---|
| 53 | const G4int G4TransitionRadiation::fPointNumber = 100 ; |
---|
| 54 | |
---|
| 55 | |
---|
| 56 | /////////////////////////////////////////////////////////////////////// |
---|
| 57 | // |
---|
| 58 | // Constructor for selected couple of materials |
---|
| 59 | // |
---|
| 60 | |
---|
| 61 | G4TransitionRadiation:: |
---|
| 62 | G4TransitionRadiation( const G4String& processName, G4ProcessType type ) |
---|
| 63 | : G4VDiscreteProcess(processName, type) |
---|
| 64 | { |
---|
[1340] | 65 | SetProcessSubType(fTransitionRadiation); |
---|
[819] | 66 | // fMatIndex1 = pMat1->GetIndex() ; |
---|
| 67 | // fMatIndex2 = pMat2->GetIndex() ; |
---|
| 68 | } |
---|
| 69 | |
---|
| 70 | ////////////////////////////////////////////////////////////////////// |
---|
| 71 | // |
---|
| 72 | // Destructor |
---|
| 73 | // |
---|
| 74 | |
---|
| 75 | G4TransitionRadiation::~G4TransitionRadiation() |
---|
[1340] | 76 | {} |
---|
[819] | 77 | |
---|
| 78 | |
---|
| 79 | /////////////////////////////////////////////////////////////////// |
---|
| 80 | // |
---|
| 81 | // Sympson integral of TR spectral-angle density over energy between |
---|
[1337] | 82 | // the limits energy 1 and energy2 at fixed varAngle = 1 - std::cos(Theta) |
---|
[819] | 83 | |
---|
| 84 | G4double |
---|
| 85 | G4TransitionRadiation::IntegralOverEnergy( G4double energy1, |
---|
| 86 | G4double energy2, |
---|
| 87 | G4double varAngle ) const |
---|
| 88 | { |
---|
| 89 | G4int i ; |
---|
| 90 | G4double h , sumEven = 0.0 , sumOdd = 0.0 ; |
---|
| 91 | h = 0.5*(energy2 - energy1)/fSympsonNumber ; |
---|
| 92 | for(i=1;i<fSympsonNumber;i++) |
---|
| 93 | { |
---|
| 94 | sumEven += SpectralAngleTRdensity(energy1 + 2*i*h,varAngle) ; |
---|
| 95 | sumOdd += SpectralAngleTRdensity(energy1 + (2*i - 1)*h,varAngle) ; |
---|
| 96 | } |
---|
| 97 | sumOdd += SpectralAngleTRdensity(energy1 + (2*fSympsonNumber - 1)*h,varAngle) ; |
---|
| 98 | return h*( SpectralAngleTRdensity(energy1,varAngle) |
---|
| 99 | + SpectralAngleTRdensity(energy2,varAngle) |
---|
| 100 | + 4.0*sumOdd + 2.0*sumEven )/3.0 ; |
---|
| 101 | } |
---|
| 102 | |
---|
| 103 | |
---|
| 104 | |
---|
| 105 | /////////////////////////////////////////////////////////////////// |
---|
| 106 | // |
---|
| 107 | // Sympson integral of TR spectral-angle density over energy between |
---|
| 108 | // the limits varAngle1 and varAngle2 at fixed energy |
---|
| 109 | |
---|
| 110 | G4double |
---|
| 111 | G4TransitionRadiation::IntegralOverAngle( G4double energy, |
---|
| 112 | G4double varAngle1, |
---|
| 113 | G4double varAngle2 ) const |
---|
| 114 | { |
---|
| 115 | G4int i ; |
---|
| 116 | G4double h , sumEven = 0.0 , sumOdd = 0.0 ; |
---|
| 117 | h = 0.5*(varAngle2 - varAngle1)/fSympsonNumber ; |
---|
| 118 | for(i=1;i<fSympsonNumber;i++) |
---|
| 119 | { |
---|
| 120 | sumEven += SpectralAngleTRdensity(energy,varAngle1 + 2*i*h) ; |
---|
| 121 | sumOdd += SpectralAngleTRdensity(energy,varAngle1 + (2*i - 1)*h) ; |
---|
| 122 | } |
---|
| 123 | sumOdd += SpectralAngleTRdensity(energy,varAngle1 + (2*fSympsonNumber - 1)*h) ; |
---|
| 124 | |
---|
| 125 | return h*( SpectralAngleTRdensity(energy,varAngle1) |
---|
| 126 | + SpectralAngleTRdensity(energy,varAngle2) |
---|
| 127 | + 4.0*sumOdd + 2.0*sumEven )/3.0 ; |
---|
| 128 | } |
---|
| 129 | |
---|
| 130 | /////////////////////////////////////////////////////////////////// |
---|
| 131 | // |
---|
| 132 | // The number of transition radiation photons generated in the |
---|
| 133 | // angle interval between varAngle1 and varAngle2 |
---|
| 134 | // |
---|
| 135 | |
---|
| 136 | G4double G4TransitionRadiation:: |
---|
| 137 | AngleIntegralDistribution( G4double varAngle1, |
---|
| 138 | G4double varAngle2 ) const |
---|
| 139 | { |
---|
| 140 | G4int i ; |
---|
| 141 | G4double h , sumEven = 0.0 , sumOdd = 0.0 ; |
---|
| 142 | h = 0.5*(varAngle2 - varAngle1)/fSympsonNumber ; |
---|
| 143 | for(i=1;i<fSympsonNumber;i++) |
---|
| 144 | { |
---|
| 145 | sumEven += IntegralOverEnergy(fMinEnergy, |
---|
| 146 | fMinEnergy +0.3*(fMaxEnergy-fMinEnergy), |
---|
| 147 | varAngle1 + 2*i*h) |
---|
| 148 | + IntegralOverEnergy(fMinEnergy + 0.3*(fMaxEnergy - fMinEnergy), |
---|
| 149 | fMaxEnergy, |
---|
| 150 | varAngle1 + 2*i*h); |
---|
| 151 | sumOdd += IntegralOverEnergy(fMinEnergy, |
---|
| 152 | fMinEnergy + 0.3*(fMaxEnergy - fMinEnergy), |
---|
| 153 | varAngle1 + (2*i - 1)*h) |
---|
| 154 | + IntegralOverEnergy(fMinEnergy + 0.3*(fMaxEnergy - fMinEnergy), |
---|
| 155 | fMaxEnergy, |
---|
| 156 | varAngle1 + (2*i - 1)*h) ; |
---|
| 157 | } |
---|
| 158 | sumOdd += IntegralOverEnergy(fMinEnergy, |
---|
| 159 | fMinEnergy + 0.3*(fMaxEnergy - fMinEnergy), |
---|
| 160 | varAngle1 + (2*fSympsonNumber - 1)*h) |
---|
| 161 | + IntegralOverEnergy(fMinEnergy + 0.3*(fMaxEnergy - fMinEnergy), |
---|
| 162 | fMaxEnergy, |
---|
| 163 | varAngle1 + (2*fSympsonNumber - 1)*h) ; |
---|
| 164 | |
---|
| 165 | return h*(IntegralOverEnergy(fMinEnergy, |
---|
| 166 | fMinEnergy + 0.3*(fMaxEnergy - fMinEnergy), |
---|
| 167 | varAngle1) |
---|
| 168 | + IntegralOverEnergy(fMinEnergy + 0.3*(fMaxEnergy - fMinEnergy), |
---|
| 169 | fMaxEnergy, |
---|
| 170 | varAngle1) |
---|
| 171 | + IntegralOverEnergy(fMinEnergy, |
---|
| 172 | fMinEnergy + 0.3*(fMaxEnergy - fMinEnergy), |
---|
| 173 | varAngle2) |
---|
| 174 | + IntegralOverEnergy(fMinEnergy + 0.3*(fMaxEnergy - fMinEnergy), |
---|
| 175 | fMaxEnergy, |
---|
| 176 | varAngle2) |
---|
| 177 | + 4.0*sumOdd + 2.0*sumEven )/3.0 ; |
---|
| 178 | } |
---|
| 179 | |
---|
| 180 | /////////////////////////////////////////////////////////////////// |
---|
| 181 | // |
---|
| 182 | // The number of transition radiation photons, generated in the |
---|
| 183 | // energy interval between energy1 and energy2 |
---|
| 184 | // |
---|
| 185 | |
---|
| 186 | G4double G4TransitionRadiation:: |
---|
| 187 | EnergyIntegralDistribution( G4double energy1, |
---|
| 188 | G4double energy2 ) const |
---|
| 189 | { |
---|
| 190 | G4int i ; |
---|
| 191 | G4double h , sumEven = 0.0 , sumOdd = 0.0 ; |
---|
| 192 | h = 0.5*(energy2 - energy1)/fSympsonNumber ; |
---|
| 193 | for(i=1;i<fSympsonNumber;i++) |
---|
| 194 | { |
---|
| 195 | sumEven += IntegralOverAngle(energy1 + 2*i*h,0.0,0.01*fMaxTheta ) |
---|
| 196 | + IntegralOverAngle(energy1 + 2*i*h,0.01*fMaxTheta,fMaxTheta); |
---|
| 197 | sumOdd += IntegralOverAngle(energy1 + (2*i - 1)*h,0.0,0.01*fMaxTheta) |
---|
| 198 | + IntegralOverAngle(energy1 + (2*i - 1)*h,0.01*fMaxTheta,fMaxTheta) ; |
---|
| 199 | } |
---|
| 200 | sumOdd += IntegralOverAngle(energy1 + (2*fSympsonNumber - 1)*h, |
---|
| 201 | 0.0,0.01*fMaxTheta) |
---|
| 202 | + IntegralOverAngle(energy1 + (2*fSympsonNumber - 1)*h, |
---|
| 203 | 0.01*fMaxTheta,fMaxTheta) ; |
---|
| 204 | |
---|
| 205 | return h*(IntegralOverAngle(energy1,0.0,0.01*fMaxTheta) |
---|
| 206 | + IntegralOverAngle(energy1,0.01*fMaxTheta,fMaxTheta) |
---|
| 207 | + IntegralOverAngle(energy2,0.0,0.01*fMaxTheta) |
---|
| 208 | + IntegralOverAngle(energy2,0.01*fMaxTheta,fMaxTheta) |
---|
| 209 | + 4.0*sumOdd + 2.0*sumEven )/3.0 ; |
---|
| 210 | } |
---|
| 211 | |
---|
| 212 | |
---|
| 213 | |
---|
| 214 | |
---|
| 215 | // end of G4TransitionRadiation implementation file -------------------------- |
---|