Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/electromagnetic/lowenergy/src/G4IonChuFluctuationModel.cc

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  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 // -------------------------------------------------------------------
 28 //
 29 // GEANT4 Class file
 30 //
 31 //
 32 // File name:     G4IonChuFluctuationModel
 33 //
 34 // Author:        V.Ivanchenko (Vladimir.Ivanchenko@cern.ch)
 35 // 
 36 // Creation date: 18 August 2000
 37 //
 38 // Modifications: 
 39 // 18/08/2000  V.Ivanchenko First implementation
 40 // 04/09/2000  V.Ivanchenko Rename fluctuations            
 41 // 03/10/2000  V.Ivanchenko CodeWizard clean up
 42 // 10/05/2001  V.Ivanchenko Clean up againist Linux compilation with -Wall
 43 //
 44 // -------------------------------------------------------------------
 45 //
 46 // Class Description: 
 47 //
 48 // The aproximation of additional ion energy loss fluctuations 
 49 // W.K.Chu, In: Ion Beam Handbook for Material Analysis.
 50 // eds. J.W. Mayer and E. Rimini (Academic Press, New York, 1977).
 51 // Q.Yang et al., NIM B61(1991)149-155.
 52 //
 53 // Class Description: End 
 54 //
 55 // -------------------------------------------------------------------
 56 //
 57 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
 58 
 59 #include "G4IonChuFluctuationModel.hh"
 60 
 61 #include "globals.hh"
 62 #include "G4PhysicalConstants.hh"
 63 #include "G4SystemOfUnits.hh"
 64 #include "G4DynamicParticle.hh"
 65 #include "G4ParticleDefinition.hh"
 66 #include "G4Material.hh"
 67 
 68 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
 69 
 70 G4IonChuFluctuationModel::G4IonChuFluctuationModel(const G4String& name)
 71   : G4VLowEnergyModel(name) 
 72 {;}
 73 
 74 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
 75 
 76 G4IonChuFluctuationModel::~G4IonChuFluctuationModel() 
 77 {;}
 78 
 79 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
 80 
 81 G4double G4IonChuFluctuationModel::TheValue(const G4DynamicParticle* particle,
 82                                       const G4Material* material) 
 83 {
 84   G4double energy = particle->GetKineticEnergy() ;
 85   G4double particleMass = particle->GetMass() ;
 86 
 87   G4double q = ChuFluctuationModel(material,energy,particleMass) ;
 88 
 89   return q ;
 90 }
 91 
 92 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
 93 
 94 G4double G4IonChuFluctuationModel::TheValue(
 95                                    const G4ParticleDefinition* aParticle,
 96                              const G4Material* material,
 97                                          G4double kineticEnergy) 
 98 {
 99   G4double particleMass = aParticle->GetPDGMass() ;
100 
101   G4double q = ChuFluctuationModel(material,kineticEnergy,particleMass);
102 
103   return q ;
104 }
105 
106 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
107 
108 G4double G4IonChuFluctuationModel::HighEnergyLimit(const G4ParticleDefinition*,
109                const G4Material* ) const
110 {
111   return 1.0*TeV ;
112 }
113 
114 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
115 
116 G4double G4IonChuFluctuationModel::LowEnergyLimit(const G4ParticleDefinition*,
117               const G4Material*) const
118 {
119   return 0.0 ;
120 }
121 
122 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
123 
124 G4double G4IonChuFluctuationModel::HighEnergyLimit(const G4ParticleDefinition*) const
125 {
126   return 1.0*TeV ;
127 }
128 
129 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
130 
131 G4double G4IonChuFluctuationModel::LowEnergyLimit(const G4ParticleDefinition*) const
132 {
133   return 0.0 ;
134 }
135 
136 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
137  
138 G4bool G4IonChuFluctuationModel::IsInCharge(const G4DynamicParticle*,
139               const G4Material*) const
140 {
141   return true ;
142 }
143 
144 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
145  
146 G4bool G4IonChuFluctuationModel::IsInCharge(const G4ParticleDefinition*,
147               const G4Material*) const
148 {
149   return true ;
150 }
151 
152 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
153 
154 G4double G4IonChuFluctuationModel::ChuFluctuationModel(const G4Material* material,
155                    G4double kineticEnergy,
156                    G4double particleMass) const
157 {
158   // The aproximation of energy loss fluctuations 
159   // Q.Yang et al., NIM B61(1991)149-155.
160 
161   // Reduced energy in MeV/AMU
162   G4double energy = kineticEnergy * amu_c2/(particleMass*MeV) ;
163 
164   G4double zeff = (material->GetElectronDensity())/
165                   (material->GetTotNbOfAtomsPerVolume()) ;
166 
167   static const G4double a[96][4] = {
168  {-0.3291, -0.8312,  0.2460, -1.0220},
169  {-0.5615, -0.5898,  0.5205, -0.7258},
170  {-0.5280, -0.4981,  0.5519, -0.5865},
171  {-0.5125, -0.4625,  0.5660, -0.5190},
172  {-0.5127, -0.8595,  0.5626, -0.8721},
173  {-0.5174, -1.1930,  0.5565, -1.1980},
174  {-0.5179, -1.1850,  0.5560, -1.2070},
175  {-0.5209, -0.9355,  0.5590, -1.0250},
176  {-0.5255, -0.7766,  0.5720, -0.9412},
177 
178  {-0.5776, -0.6665,  0.6598, -0.8484},
179  {-0.6013, -0.6045,  0.7321, -0.7671},
180  {-0.5781, -0.5518,  0.7605, -0.6919},
181  {-0.5587, -0.4981,  0.7835, -0.6195},
182  {-0.5466, -0.4656,  0.7978, -0.5771},
183  {-0.5406, -0.4690,  0.8031, -0.5718},
184  {-0.5391, -0.5061,  0.8024, -0.5974},
185  {-0.5380, -0.6483,  0.7962, -0.6970},
186  {-0.5355, -0.7722,  0.7962, -0.7839},
187  {-0.5329, -0.7720,  0.7988, -0.7846},
188 
189  {-0.5335, -0.7671,  0.7984, -0.7933},
190  {-0.5324, -0.7612,  0.7998, -0.8031},
191  {-0.5305, -0.7300,  0.8031, -0.7990},
192  {-0.5307, -0.7178,  0.8049, -0.8216},
193  {-0.5248, -0.6621,  0.8165, -0.7919},
194  {-0.5180, -0.6502,  0.8266, -0.7986},
195  {-0.5084, -0.6408,  0.8396, -0.8048},
196  {-0.4967, -0.6331,  0.8549, -0.8093},
197  {-0.4861, -0.6508,  0.8712, -0.8432},
198  {-0.4700, -0.6186,  0.8961, -0.8132},
199 
200  {-0.4545, -0.5720,  0.9227, -0.7710},
201  {-0.4404, -0.5226,  0.9481, -0.7254},
202  {-0.4288, -0.4778,  0.9701, -0.6850},
203  {-0.4199, -0.4425,  0.9874, -0.6539},
204  {-0.4131, -0.4188,  0.9998, -0.6332},
205  {-0.4089, -0.4057,  1.0070, -0.6218},
206  {-0.4039, -0.3913,  1.0150, -0.6107},
207  {-0.3987, -0.3698,  1.0240, -0.5938},
208  {-0.3977, -0.3608,  1.0260, -0.5852},
209  {-0.3972, -0.3600,  1.0260, -0.5842},
210 
211  {-0.3985, -0.3803,  1.0200, -0.6013},
212  {-0.3985, -0.3979,  1.0150, -0.6168},
213  {-0.3968, -0.3990,  1.0160, -0.6195},
214  {-0.3971, -0.4432,  1.0050, -0.6591},
215  {-0.3944, -0.4665,  1.0010, -0.6825},
216  {-0.3924, -0.5109,  0.9921, -0.7235},
217  {-0.3882, -0.5158,  0.9947, -0.7343},
218  {-0.3838, -0.5125,  0.9999, -0.7370},
219  {-0.3786, -0.4976,  1.0090, -0.7310},
220  {-0.3741, -0.4738,  1.0200, -0.7155},
221 
222  {-0.3969, -0.4496,  1.0320, -0.6982},
223  {-0.3663, -0.4297,  1.0430, -0.6828},
224  {-0.3630, -0.4120,  1.0530, -0.6689},
225  {-0.3597, -0.3964,  1.0620, -0.6564},
226  {-0.3555, -0.3809,  1.0720, -0.6454},
227  {-0.3525, -0.3607,  1.0820, -0.6289},
228  {-0.3505, -0.3465,  1.0900, -0.6171},
229  {-0.3397, -0.3570,  1.1020, -0.6384},
230  {-0.3314, -0.3552,  1.1130, -0.6441},
231  {-0.3235, -0.3531,  1.1230, -0.6498},
232 
233  {-0.3150, -0.3483,  1.1360, -0.6539},
234  {-0.3060, -0.3441,  1.1490, -0.6593},
235  {-0.2968, -0.3396,  1.1630, -0.6649},
236  {-0.2935, -0.3225,  1.1760, -0.6527},
237  {-0.2797, -0.3262,  1.1940, -0.6722},
238  {-0.2704, -0.3202,  1.2100, -0.6770},
239  {-0.2815, -0.3227,  1.2480, -0.6775},
240  {-0.2880, -0.3245,  1.2810, -0.6801},
241  {-0.3034, -0.3263,  1.3270, -0.6778},
242  {-0.2936, -0.3215,  1.3430, -0.6835},
243 
244  {-0.3282, -0.3200,  1.3980, -0.6650},
245  {-0.3260, -0.3070,  1.4090, -0.6552},
246  {-0.3511, -0.3074,  1.4470, -0.6442},
247  {-0.3501, -0.3064,  1.4500, -0.6442},
248  {-0.3490, -0.3027,  1.4550, -0.6418},
249  {-0.3487, -0.3048,  1.4570, -0.6447},
250  {-0.3478, -0.3074,  1.4600, -0.6483},
251  {-0.3501, -0.3283,  1.4540, -0.6669},
252  {-0.3494, -0.3373,  1.4550, -0.6765},
253  {-0.3485, -0.3373,  1.4570, -0.6774},
254 
255  {-0.3462, -0.3300,  1.4630, -0.6728},
256  {-0.3462, -0.3225,  1.4690, -0.6662},
257  {-0.3453, -0.3094,  1.4790, -0.6553},
258  {-0.3844, -0.3134,  1.5240, -0.6412},
259  {-0.3848, -0.3018,  1.5310, -0.6303},
260  {-0.3862, -0.2955,  1.5360, -0.6237},
261  {-0.4262, -0.2991,  1.5860, -0.6115},
262  {-0.4278, -0.2910,  1.5900, -0.6029},
263  {-0.4303, -0.2817,  1.5940, -0.5927},
264  {-0.4315, -0.2719,  1.6010, -0.5829},
265 
266  {-0.4359, -0.2914,  1.6050, -0.6010},
267  {-0.4365, -0.2982,  1.6080, -0.6080},
268  {-0.4253, -0.3037,  1.6120, -0.6150},
269  {-0.4335, -0.3245,  1.6160, -0.6377},
270  {-0.4307, -0.3292,  1.6210, -0.6447},
271  {-0.4284, -0.3204,  1.6290, -0.6380},
272  {-0.4227, -0.3217,  1.6360, -0.6438}
273   } ;
274 
275   G4int iz = (G4int)zeff - 2 ;
276   if( 0 > iz ) iz = 0 ;
277   if(95 < iz ) iz = 95 ;
278 
279   G4double q = 1.0 / (1.0 + a[iz][0]*std::pow(energy,a[iz][1])+
280                           + a[iz][2]*std::pow(energy,a[iz][3])) ; 
281 
282   return q ;    
283 }
284 
285 
286