[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 | // |
---|
| 27 | // $Id: G4PAIxSection.hh,v 1.12 2006/06/29 19:50:44 gunter Exp $ |
---|
| 28 | // GEANT4 tag $Name: $ |
---|
| 29 | // |
---|
| 30 | // |
---|
| 31 | // G4PAIxSection.hh -- header file |
---|
| 32 | // |
---|
| 33 | // GEANT 4 class header file --- Copyright CERN 1995 |
---|
| 34 | // CERB Geneva Switzerland |
---|
| 35 | // |
---|
| 36 | // for information related to this code, please, contact |
---|
| 37 | // CERN, CN Division, ASD Group |
---|
| 38 | // |
---|
| 39 | // Preparation of ionizing collision cross section according to Photo Absorption |
---|
| 40 | // Ionization (PAI) model for simulation of ionization energy losses in very thin |
---|
| 41 | // absorbers. Author: Vladimir.Grichine@cern.ch |
---|
| 42 | // |
---|
| 43 | // History: |
---|
| 44 | // |
---|
| 45 | // 19.10.03, V. Grichine: Integral dEdx was added for G4PAIModel class |
---|
| 46 | // |
---|
| 47 | // 13.05.03, V. Grichine: Numerical instability was fixed in SumOverInterval/Border |
---|
| 48 | // functions |
---|
| 49 | // 10.02.02, V. Grichine: New functions and arrays/gets for Cerenkov and |
---|
| 50 | // plasmon collisions dN/dx |
---|
| 51 | // 27.10.99, V. Grichine: Bug fixed in constructors, 3rd constructor and |
---|
| 52 | // GetStepEnergyLoss(step) were added, fDelta = 0.005 |
---|
| 53 | // 30.11.97, V. Grichine: 2nd version |
---|
| 54 | // 11.06.97, V. Grichine: 1st version |
---|
| 55 | |
---|
| 56 | #ifndef G4PAIXSECTION_HH |
---|
| 57 | #define G4PAIXSECTION_HH |
---|
| 58 | |
---|
| 59 | #include "G4ios.hh" |
---|
| 60 | #include "globals.hh" |
---|
| 61 | #include "Randomize.hh" |
---|
| 62 | |
---|
| 63 | #include"G4SandiaTable.hh" |
---|
| 64 | |
---|
| 65 | class G4MaterialCutsCouple; |
---|
| 66 | class G4Sandiatable; |
---|
| 67 | |
---|
| 68 | |
---|
| 69 | class G4PAIxSection |
---|
| 70 | { |
---|
| 71 | public: |
---|
| 72 | // Constructors |
---|
| 73 | G4PAIxSection( G4MaterialCutsCouple* matCC); |
---|
| 74 | |
---|
| 75 | G4PAIxSection( G4int materialIndex, |
---|
| 76 | G4double maxEnergyTransfer ) ; |
---|
| 77 | |
---|
| 78 | G4PAIxSection( G4int materialIndex, // for proton loss table |
---|
| 79 | G4double maxEnergyTransfer, |
---|
| 80 | G4double betaGammaSq , |
---|
| 81 | G4double** photoAbsCof, G4int intNumber ) ; |
---|
| 82 | |
---|
| 83 | G4PAIxSection( G4int materialIndex, // test constructor |
---|
| 84 | G4double maxEnergyTransfer, |
---|
| 85 | G4double betaGammaSq ) ; |
---|
| 86 | |
---|
| 87 | // G4PAIxSection(const G4PAIxSection& right) ; |
---|
| 88 | |
---|
| 89 | // Destructor |
---|
| 90 | |
---|
| 91 | ~G4PAIxSection() ; |
---|
| 92 | |
---|
| 93 | // Operators |
---|
| 94 | // G4PAIxSection& operator=(const G4PAIxSection& right) ; |
---|
| 95 | // G4int operator==(const G4PAIxSection& right)const ; |
---|
| 96 | // G4int operator!=(const G4PAIxSection& right)const ; |
---|
| 97 | |
---|
| 98 | // Methods |
---|
| 99 | |
---|
| 100 | // General control functions |
---|
| 101 | |
---|
| 102 | void InitPAI() ; |
---|
| 103 | |
---|
| 104 | void NormShift( G4double betaGammaSq ) ; |
---|
| 105 | |
---|
| 106 | void SplainPAI( G4double betaGammaSq ) ; |
---|
| 107 | |
---|
| 108 | // Physical methods |
---|
| 109 | |
---|
| 110 | |
---|
| 111 | G4double RutherfordIntegral( G4int intervalNumber, |
---|
| 112 | G4double limitLow, |
---|
| 113 | G4double limitHigh ) ; |
---|
| 114 | |
---|
| 115 | G4double ImPartDielectricConst( G4int intervalNumber, |
---|
| 116 | G4double energy ) ; |
---|
| 117 | |
---|
| 118 | G4double RePartDielectricConst(G4double energy) ; |
---|
| 119 | |
---|
| 120 | G4double DifPAIxSection( G4int intervalNumber, |
---|
| 121 | G4double betaGammaSq ) ; |
---|
| 122 | |
---|
| 123 | G4double PAIdNdxCerenkov( G4int intervalNumber, |
---|
| 124 | G4double betaGammaSq ) ; |
---|
| 125 | |
---|
| 126 | G4double PAIdNdxPlasmon( G4int intervalNumber, |
---|
| 127 | G4double betaGammaSq ) ; |
---|
| 128 | |
---|
| 129 | void IntegralPAIxSection() ; |
---|
| 130 | void IntegralCerenkov() ; |
---|
| 131 | void IntegralPlasmon() ; |
---|
| 132 | |
---|
| 133 | G4double SumOverInterval(G4int intervalNumber) ; |
---|
| 134 | G4double SumOverIntervaldEdx(G4int intervalNumber) ; |
---|
| 135 | G4double SumOverInterCerenkov(G4int intervalNumber) ; |
---|
| 136 | G4double SumOverInterPlasmon(G4int intervalNumber) ; |
---|
| 137 | |
---|
| 138 | G4double SumOverBorder( G4int intervalNumber, |
---|
| 139 | G4double energy ) ; |
---|
| 140 | G4double SumOverBorderdEdx( G4int intervalNumber, |
---|
| 141 | G4double energy ) ; |
---|
| 142 | G4double SumOverBordCerenkov( G4int intervalNumber, |
---|
| 143 | G4double energy ) ; |
---|
| 144 | G4double SumOverBordPlasmon( G4int intervalNumber, |
---|
| 145 | G4double energy ) ; |
---|
| 146 | |
---|
| 147 | G4double GetStepEnergyLoss( G4double step ) ; |
---|
| 148 | G4double GetStepCerenkovLoss( G4double step ) ; |
---|
| 149 | G4double GetStepPlasmonLoss( G4double step ) ; |
---|
| 150 | |
---|
| 151 | // Inline access functions |
---|
| 152 | |
---|
| 153 | G4int GetNumberOfGammas() const { return fNumberOfGammas ; } |
---|
| 154 | |
---|
| 155 | G4int GetSplineSize() const { return fSplineNumber ; } |
---|
| 156 | |
---|
| 157 | G4int GetIntervalNumber() const { return fIntervalNumber ; } |
---|
| 158 | |
---|
| 159 | G4double GetEnergyInterval(G4int i){ return fEnergyInterval[i] ; } |
---|
| 160 | |
---|
| 161 | G4double GetDifPAIxSection(G4int i){ return fDifPAIxSection[i] ; } |
---|
| 162 | G4double GetPAIdNdxCrenkov(G4int i){ return fdNdxCerenkov[i] ; } |
---|
| 163 | G4double GetPAIdNdxPlasmon(G4int i){ return fdNdxPlasmon[i] ; } |
---|
| 164 | |
---|
| 165 | G4double GetMeanEnergyLoss() const {return fIntegralPAIxSection[0] ; } |
---|
| 166 | G4double GetMeanCerenkovLoss() const {return fIntegralCerenkov[0] ; } |
---|
| 167 | G4double GetMeanPlasmonLoss() const {return fIntegralPlasmon[0] ; } |
---|
| 168 | |
---|
| 169 | G4double GetNormalizationCof() const { return fNormalizationCof ; } |
---|
| 170 | |
---|
| 171 | inline G4double GetPAItable(G4int i,G4int j) const ; |
---|
| 172 | |
---|
| 173 | inline G4double GetLorentzFactor(G4int i) const ; |
---|
| 174 | |
---|
| 175 | inline G4double GetSplineEnergy(G4int i) const ; |
---|
| 176 | |
---|
| 177 | inline G4double GetIntegralPAIxSection(G4int i) const ; |
---|
| 178 | inline G4double GetIntegralPAIdEdx(G4int i) const ; |
---|
| 179 | inline G4double GetIntegralCerenkov(G4int i) const ; |
---|
| 180 | inline G4double GetIntegralPlasmon(G4int i) const ; |
---|
| 181 | |
---|
| 182 | protected : |
---|
| 183 | |
---|
| 184 | private : |
---|
| 185 | |
---|
| 186 | // Local class constants |
---|
| 187 | |
---|
| 188 | static const G4double fDelta ; // energy shift from interval border = 0.001 |
---|
| 189 | static const G4double fError ; // error in lin-log approximation = 0.005 |
---|
| 190 | |
---|
| 191 | static G4int fNumberOfGammas ; // = 111 ; |
---|
| 192 | static const G4double fLorentzFactor[112] ; // static gamma array |
---|
| 193 | |
---|
| 194 | static |
---|
| 195 | const G4int fRefGammaNumber ; // The number of gamma for creation of spline (15) |
---|
| 196 | |
---|
| 197 | G4int fIntervalNumber ; // The number of energy intervals |
---|
| 198 | G4double fNormalizationCof ; // Normalization cof for PhotoAbsorptionXsection |
---|
| 199 | |
---|
| 200 | // G4double fBetaGammaSq ; // (beta*gamma)^2 |
---|
| 201 | |
---|
| 202 | G4double fDensity ; // Current density |
---|
| 203 | G4double fElectronDensity ; // Current electron (number) density |
---|
| 204 | G4int fSplineNumber ; // Current size of spline |
---|
| 205 | |
---|
| 206 | // Arrays of Sandia coefficients |
---|
| 207 | |
---|
| 208 | G4OrderedTable* fMatSandiaMatrix; |
---|
| 209 | G4SandiaTable* fSandia; |
---|
| 210 | |
---|
| 211 | G4double* fEnergyInterval ; |
---|
| 212 | G4double* fA1 ; |
---|
| 213 | G4double* fA2 ; |
---|
| 214 | G4double* fA3 ; |
---|
| 215 | G4double* fA4 ; |
---|
| 216 | |
---|
| 217 | static |
---|
| 218 | const G4int fMaxSplineSize ; // Max size of output splain arrays = 500 |
---|
| 219 | |
---|
| 220 | /* ****************** |
---|
| 221 | G4double* fSplineEnergy ; // energy points of splain |
---|
| 222 | G4double* fRePartDielectricConst ; // Real part of dielectric const |
---|
| 223 | G4double* fImPartDielectricConst ; // Imaginary part of dielectric const |
---|
| 224 | G4double* fIntegralTerm ; // Integral term in PAI cross section |
---|
| 225 | G4double* fDifPAIxSection ; // Differential PAI cross section |
---|
| 226 | G4double* fIntegralPAIxSection ; // Integral PAI cross section ? |
---|
| 227 | */ /////////////// |
---|
| 228 | |
---|
| 229 | |
---|
| 230 | G4double fSplineEnergy[500] ; // energy points of splain |
---|
| 231 | G4double fRePartDielectricConst[500] ; // Real part of dielectric const |
---|
| 232 | G4double fImPartDielectricConst[500] ; // Imaginary part of dielectric const |
---|
| 233 | G4double fIntegralTerm[500] ; // Integral term in PAI cross section |
---|
| 234 | G4double fDifPAIxSection[500] ; // Differential PAI cross section |
---|
| 235 | G4double fdNdxCerenkov[500] ; // dNdx of Cerenkov collisions |
---|
| 236 | G4double fdNdxPlasmon[500] ; // dNdx of Plasmon collisions |
---|
| 237 | |
---|
| 238 | G4double fIntegralPAIxSection[500] ; // Integral PAI cross section ? |
---|
| 239 | G4double fIntegralPAIdEdx[500] ; // Integral PAI dEdx ? |
---|
| 240 | G4double fIntegralCerenkov[500] ; // Integral Cerenkov N>omega ? |
---|
| 241 | G4double fIntegralPlasmon[500] ; // Integral Plasmon N>omega ? |
---|
| 242 | |
---|
| 243 | G4double fPAItable[500][112] ; // Output array |
---|
| 244 | |
---|
| 245 | } ; |
---|
| 246 | |
---|
| 247 | //////////////// Inline methods ////////////////////////////////// |
---|
| 248 | // |
---|
| 249 | |
---|
| 250 | |
---|
| 251 | inline G4double G4PAIxSection::GetPAItable(G4int i, G4int j) const |
---|
| 252 | { |
---|
| 253 | return fPAItable[i][j] ; |
---|
| 254 | } |
---|
| 255 | |
---|
| 256 | inline G4double G4PAIxSection::GetLorentzFactor(G4int j) const |
---|
| 257 | { |
---|
| 258 | return fLorentzFactor[j] ; |
---|
| 259 | } |
---|
| 260 | |
---|
| 261 | inline G4double G4PAIxSection::GetSplineEnergy(G4int i) const |
---|
| 262 | { |
---|
| 263 | if(i < 1 || i > fSplineNumber) |
---|
| 264 | { |
---|
| 265 | G4Exception("Invalid argument in G4PAIxSection::GetSplineEnergy"); |
---|
| 266 | } |
---|
| 267 | return fSplineEnergy[i] ; |
---|
| 268 | } |
---|
| 269 | |
---|
| 270 | inline G4double G4PAIxSection::GetIntegralPAIxSection(G4int i) const |
---|
| 271 | { |
---|
| 272 | if(i < 1 || i > fSplineNumber) |
---|
| 273 | { |
---|
| 274 | G4Exception("Invalid argument in G4PAIxSection::GetIntegralPAIxSection"); |
---|
| 275 | } |
---|
| 276 | return fIntegralPAIxSection[i] ; |
---|
| 277 | } |
---|
| 278 | |
---|
| 279 | inline G4double G4PAIxSection::GetIntegralPAIdEdx(G4int i) const |
---|
| 280 | { |
---|
| 281 | if(i < 1 || i > fSplineNumber) |
---|
| 282 | { |
---|
| 283 | G4Exception("Invalid argument in G4PAIxSection::GetIntegralPAIxSection"); |
---|
| 284 | } |
---|
| 285 | return fIntegralPAIdEdx[i] ; |
---|
| 286 | } |
---|
| 287 | |
---|
| 288 | inline G4double G4PAIxSection::GetIntegralCerenkov(G4int i) const |
---|
| 289 | { |
---|
| 290 | if(i < 1 || i > fSplineNumber) |
---|
| 291 | { |
---|
| 292 | G4Exception("Invalid argument in G4PAIxSection::GetIntegralCerenkov"); |
---|
| 293 | } |
---|
| 294 | return fIntegralCerenkov[i] ; |
---|
| 295 | } |
---|
| 296 | |
---|
| 297 | inline G4double G4PAIxSection::GetIntegralPlasmon(G4int i) const |
---|
| 298 | { |
---|
| 299 | if(i < 1 || i > fSplineNumber) |
---|
| 300 | { |
---|
| 301 | G4Exception("Invalid argument in G4PAIxSection::GetIntegralPlasmon"); |
---|
| 302 | } |
---|
| 303 | return fIntegralPlasmon[i] ; |
---|
| 304 | } |
---|
| 305 | |
---|
| 306 | #endif |
---|
| 307 | |
---|
| 308 | // ----------------- end of G4PAIxSection header file ------------------- |
---|