source: trunk/source/processes/electromagnetic/standard/src/G4ComptonScattering52.cc @ 1228

Last change on this file since 1228 was 1228, checked in by garnier, 14 years ago

update geant4.9.3 tag

File size: 18.1 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// $Id: G4ComptonScattering52.cc,v 1.7 2008/10/15 17:53:44 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-03 $
28//
29//
30//------------ G4ComptonScattering52 physics process -----------------------------
31//                   by Michel Maire, April 1996
32//
33// 28-05-96, DoIt() small change in ElecDirection, by M.Maire
34// 10-06-96, simplification in ComputeMicroscopicCrossSection(), by M.Maire
35// 21-06-96, SetCuts implementation, M.Maire
36// 13-09-96, small changes in DoIt for better efficiency. Thanks to P.Urban
37// 06-01-97, crossection table + meanfreepath table, M.Maire
38// 05-03-97, new Physics scheme, M.Maire
39// 28-03-97, protection in BuildPhysicsTable, M.Maire
40// 07-04-98, remove 'tracking cut' of the scattered gamma, MMa
41// 04-06-98, in DoIt, secondary production condition:
42//                                     range>std::min(threshold,safety)
43// 13-08-98, new methods SetBining()  PrintInfo()
44// 15-12-98, cross section=0 below 10 keV
45// 28-05-01, V.Ivanchenko minor changes to provide ANSI -wall compilation
46// 13-07-01, DoIt: suppression of production cut for the electron (mma)
47// 03-08-01, new methods Store/Retrieve PhysicsTable (mma)
48// 06-08-01, BuildThePhysicsTable() called from constructor (mma)
49// 17-09-01, migration of Materials to pure STL (mma)
50// 20-09-01, DoIt: fminimalEnergy = 1*eV (mma)
51// 01-10-01, come back to BuildPhysicsTable(const G4ParticleDefinition&)
52// 17-04-02, LowestEnergyLimit = 1*keV     
53// 26-05-04, cross section parametrization improved for low energy :
54//           Egamma <~ 15 keV (Laszlo)
55// 08-11-04, Remove Store/Retrieve tables (V.Ivanchenko)
56// 04-05-05, Add 52 to class name (V.Ivanchenko)
57// -----------------------------------------------------------------------------
58
59#include "G4ComptonScattering52.hh"
60#include "G4UnitsTable.hh"
61#include "G4PhysicsTableHelper.hh"
62
63//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
64
65using namespace std;
66
67G4ComptonScattering52::G4ComptonScattering52(const G4String& processName,
68    G4ProcessType type):G4VDiscreteProcess (processName, type),
69    theCrossSectionTable(NULL),
70    theMeanFreePathTable(NULL),
71    LowestEnergyLimit (  1*keV),
72    HighestEnergyLimit(100*GeV),
73    NumbBinTable(80),
74    fminimalEnergy(1*eV)
75{
76  SetProcessSubType(13);
77  G4cout << "!!! G4ComptonScattering52 is the obsolete process class and will be removed soon !!!"
78         << G4endl;
79}
80
81//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
82 
83// destructor
84 
85G4ComptonScattering52::~G4ComptonScattering52()
86{
87   if (theCrossSectionTable) {
88      theCrossSectionTable->clearAndDestroy();
89      delete theCrossSectionTable;
90   }
91
92   if (theMeanFreePathTable) {
93      theMeanFreePathTable->clearAndDestroy();
94      delete theMeanFreePathTable;
95   }
96}
97
98//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
99
100G4bool G4ComptonScattering52::IsApplicable( const G4ParticleDefinition& particle)
101{
102   return ( &particle == G4Gamma::Gamma() ); 
103}
104 
105//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
106
107void G4ComptonScattering52::SetPhysicsTableBining(
108                                   G4double lowE, G4double highE, G4int nBins)
109{
110  LowestEnergyLimit = lowE; HighestEnergyLimit = highE; NumbBinTable = nBins;
111} 
112
113//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
114 
115void G4ComptonScattering52::BuildPhysicsTable(const G4ParticleDefinition&)
116// Build cross section and mean free path tables
117{
118   G4double LowEdgeEnergy, Value;
119   G4PhysicsLogVector* ptrVector;
120
121// Build cross section per atom tables for the Compton Scattering process
122
123   if (theCrossSectionTable) {
124       theCrossSectionTable->clearAndDestroy(); delete theCrossSectionTable;}
125
126   theCrossSectionTable = new G4PhysicsTable(G4Element::GetNumberOfElements());
127   const G4ElementTable* theElementTable = G4Element::GetElementTable();
128   G4double AtomicNumber;
129   size_t J;
130
131   for ( J=0 ; J < G4Element::GetNumberOfElements(); J++ )
132       {
133        //create physics vector then fill it ....
134        ptrVector = new G4PhysicsLogVector(LowestEnergyLimit,HighestEnergyLimit,
135                                           NumbBinTable );
136        AtomicNumber = (*theElementTable)[J]->GetZ();
137
138        for ( G4int i = 0 ; i < NumbBinTable ; i++ )
139           {
140             LowEdgeEnergy = ptrVector->GetLowEdgeEnergy(i);
141             Value = ComputeCrossSectionPerAtom(LowEdgeEnergy, AtomicNumber);
142             ptrVector->PutValue(i,Value);
143           }
144
145        theCrossSectionTable->insertAt( J , ptrVector ) ;
146
147      }
148
149// Build mean free path table for the Compton Scattering process
150
151   if (theMeanFreePathTable) {
152       theMeanFreePathTable->clearAndDestroy(); delete theMeanFreePathTable;}
153 
154   theMeanFreePathTable= new G4PhysicsTable(G4Material::GetNumberOfMaterials());
155   const G4MaterialTable* theMaterialTable = G4Material::GetMaterialTable();
156   G4Material* material;
157
158   for ( J=0 ; J < G4Material::GetNumberOfMaterials(); J++ ) 
159       { 
160        //create physics vector then fill it ....
161        ptrVector = new G4PhysicsLogVector(LowestEnergyLimit,HighestEnergyLimit,
162                                           NumbBinTable ) ;
163        material = (*theMaterialTable)[J];
164
165        for ( G4int i = 0 ; i < NumbBinTable ; i++ )     
166           {
167             LowEdgeEnergy = ptrVector->GetLowEdgeEnergy( i ) ;
168             Value = ComputeMeanFreePath( LowEdgeEnergy, material); 
169             ptrVector->PutValue( i , Value ) ;
170           }
171
172        theMeanFreePathTable->insertAt( J , ptrVector ) ;
173      }
174
175    PrintInfoDefinition(); 
176
177}
178
179//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
180
181G4double G4ComptonScattering52::ComputeCrossSectionPerAtom
182                              (G4double GammaEnergy, G4double Z)
183 
184// Calculates the cross section per atom in GEANT4 internal units.
185// A parametrized formula from L. Urban is used to estimate
186// the total cross section.
187// It gives a good description of the data from 10 keV to 100/Z GeV.
188// lower limit 1 keV now with a correction for low energy
189 
190{
191 G4double CrossSection = 0.0 ;
192 if ( Z < 1. )                     return CrossSection;
193 if ( GammaEnergy < 1.*keV       ) return CrossSection;
194 if ( GammaEnergy > (100.*GeV/Z) ) return CrossSection;
195
196 static const G4double a = 20.0 , b = 230.0 , c = 440.0;
197 
198 static const G4double
199 d1= 2.7965e-1*barn, d2=-1.8300e-1*barn, d3= 6.7527   *barn, d4=-1.9798e+1*barn,
200 e1= 1.9756e-5*barn, e2=-1.0205e-2*barn, e3=-7.3913e-2*barn, e4= 2.7079e-2*barn,
201 f1=-3.9178e-7*barn, f2= 6.8241e-5*barn, f3= 6.0480e-5*barn, f4= 3.0274e-4*barn;
202     
203 G4double p1Z = Z*(d1 + e1*Z + f1*Z*Z), p2Z = Z*(d2 + e2*Z + f2*Z*Z),
204          p3Z = Z*(d3 + e3*Z + f3*Z*Z), p4Z = Z*(d4 + e4*Z + f4*Z*Z);
205
206 G4double T0 = 15*keV; if (Z == 1.) T0 = 40*keV; 
207
208 G4double X = max(GammaEnergy, T0) / electron_mass_c2;
209 CrossSection = p1Z*log(1.+2*X)/X
210                + (p2Z + p3Z*X + p4Z*X*X)/(1. + a*X + b*X*X + c*X*X*X);
211               
212 //  modification for low energy. (special case for Hydrogen)
213 if (GammaEnergy < T0) {
214   G4double dT0 = 1.*keV;
215   X = (T0+dT0) / electron_mass_c2 ;
216   G4double sigma = p1Z*log(1.+2*X)/X
217                    + (p2Z + p3Z*X + p4Z*X*X)/(1. + a*X + b*X*X + c*X*X*X);
218   G4double c1 = -T0*(sigma-CrossSection)/(CrossSection*dT0);             
219   G4double c2 = 0.150; if (Z > 1.) c2 = 0.375-0.0556*log(Z);
220   G4double  y = log(GammaEnergy/T0);
221   CrossSection *= exp(-y*(c1+c2*y));         
222   }
223
224 return CrossSection;
225}
226
227//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
228
229G4double G4ComptonScattering52::ComputeMeanFreePath(G4double GammaEnergy,
230                                                         G4Material* aMaterial)
231
232// returns the gamma mean free path in GEANT4 internal units
233
234{
235  const G4ElementVector* theElementVector = aMaterial->GetElementVector() ;
236  const G4double* NbOfAtomsPerVolume = aMaterial->GetVecNbOfAtomsPerVolume();
237
238  G4double SIGMA = 0.;
239
240  for ( size_t elm=0 ; elm < aMaterial->GetNumberOfElements() ; elm++ )
241      {
242        SIGMA += NbOfAtomsPerVolume[elm] * 
243                 ComputeCrossSectionPerAtom(GammaEnergy,
244                                            (*theElementVector)[elm]->GetZ());
245      }       
246  return SIGMA > DBL_MIN ? 1./SIGMA : DBL_MAX;
247}
248
249 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
250
251G4double G4ComptonScattering52::GetCrossSectionPerAtom(
252                                 G4DynamicParticle* aDynamicGamma,
253                                 G4Element*         anElement)
254 
255// gives the microscopic total cross section in GEANT4 internal units
256
257{
258   G4double crossSection;
259   G4double GammaEnergy = aDynamicGamma->GetKineticEnergy();
260   G4bool isOutRange ;
261   if (GammaEnergy < LowestEnergyLimit || GammaEnergy > HighestEnergyLimit)
262      crossSection = 0.;
263   else
264      crossSection = (*theCrossSectionTable)(anElement->GetIndex())->
265                                       GetValue(GammaEnergy, isOutRange);
266
267   return crossSection;
268}
269
270//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
271
272
273G4double G4ComptonScattering52::GetMeanFreePath(const G4Track& aTrack,
274                                              G4double,
275                                              G4ForceCondition*)
276
277// returns the gamma mean free path in GEANT4 internal units
278
279{
280   const G4DynamicParticle* aDynamicGamma = aTrack.GetDynamicParticle();
281   G4double GammaEnergy = aDynamicGamma->GetKineticEnergy();
282   G4Material* aMaterial = aTrack.GetMaterial();
283
284   G4double MeanFreePath;
285   G4bool isOutRange;
286
287   if (GammaEnergy > HighestEnergyLimit || GammaEnergy < LowestEnergyLimit)
288     MeanFreePath = DBL_MAX;
289   else
290     MeanFreePath = (*theMeanFreePathTable)(aMaterial->GetIndex())->
291                                       GetValue(GammaEnergy, isOutRange);
292   return MeanFreePath;
293}
294
295//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
296
297G4VParticleChange* G4ComptonScattering52::PostStepDoIt(const G4Track& aTrack,
298                                                     const G4Step&  aStep)
299//
300// The scattered gamma energy is sampled according to Klein - Nishina formula.
301// The random number techniques of Butcher & Messel are used
302// (Nuc Phys 20(1960),15).
303// GEANT4 internal units
304//
305// Note : Effects due to binding of atomic electrons are negliged.
306 
307{
308   aParticleChange.Initialize(aTrack);
309
310   const G4DynamicParticle* aDynamicGamma = aTrack.GetDynamicParticle();
311   G4double GammaEnergy0 = aDynamicGamma->GetKineticEnergy();
312   G4double E0_m = GammaEnergy0 / electron_mass_c2 ;
313
314   G4ParticleMomentum GammaDirection0 = aDynamicGamma->GetMomentumDirection();
315
316   //
317   // sample the energy rate of the scattered gamma
318   //
319
320   G4double epsilon, epsilonsq, onecost, sint2, greject ;
321
322   G4double epsilon0 = 1./(1. + 2*E0_m) , epsilon0sq = epsilon0*epsilon0;
323   G4double alpha1   = - log(epsilon0)  , alpha2 = 0.5*(1.- epsilon0sq);
324
325   do {
326       if ( alpha1/(alpha1+alpha2) > G4UniformRand() )
327            { epsilon   = exp(-alpha1*G4UniformRand());   // epsilon0**r
328              epsilonsq = epsilon*epsilon; }
329       else {
330             epsilonsq = epsilon0sq + (1.- epsilon0sq)*G4UniformRand();
331             epsilon   = sqrt(epsilonsq);
332       };
333       onecost = (1.- epsilon)/(epsilon*E0_m);
334       sint2   = onecost*(2.-onecost);
335       greject = 1. - epsilon*sint2/(1.+ epsilonsq);
336   } while (greject < G4UniformRand());
337 
338   //
339   // scattered gamma angles. ( Z - axis along the parent gamma)
340   //
341
342   G4double cosTeta = 1. - onecost , sinTeta = sqrt (sint2);
343   G4double Phi     = twopi * G4UniformRand();
344   G4double dirx = sinTeta*cos(Phi), diry = sinTeta*sin(Phi), dirz = cosTeta;
345
346   //
347   // update G4VParticleChange for the scattered gamma
348   //
349   
350   G4ThreeVector GammaDirection1 ( dirx,diry,dirz );
351   GammaDirection1.rotateUz(GammaDirection0);
352   aParticleChange.ProposeMomentumDirection( GammaDirection1 );
353   G4double GammaEnergy1 = epsilon*GammaEnergy0;
354   G4double localEnergyDeposit = 0.;
355   
356   if (GammaEnergy1 > fminimalEnergy)
357     {
358       aParticleChange.ProposeEnergy( GammaEnergy1 );
359     }
360   else
361     {
362       localEnergyDeposit += GammaEnergy1;   
363       aParticleChange.ProposeEnergy(0.) ;
364       aParticleChange.ProposeTrackStatus(fStopAndKill);
365     }
366       
367   //
368   // kinematic of the scattered electron
369   //
370
371   G4double ElecKineEnergy = GammaEnergy0 - GammaEnergy1;
372
373    if (ElecKineEnergy > fminimalEnergy)       
374      {
375        G4double ElecMomentum = sqrt(ElecKineEnergy*
376                                    (ElecKineEnergy+2.*electron_mass_c2));
377        G4ThreeVector ElecDirection (
378        (GammaEnergy0*GammaDirection0 - GammaEnergy1*GammaDirection1)
379        *(1./ElecMomentum) );
380
381        // create G4DynamicParticle object for the electron.
382        G4DynamicParticle* aElectron= new G4DynamicParticle(
383                           G4Electron::Electron(),ElecDirection,ElecKineEnergy);
384
385        aParticleChange.SetNumberOfSecondaries(1);
386        aParticleChange.AddSecondary( aElectron );
387      }
388    else
389      {
390        aParticleChange.SetNumberOfSecondaries(0);
391        localEnergyDeposit += ElecKineEnergy;
392      }
393    aParticleChange.ProposeLocalEnergyDeposit (localEnergyDeposit);
394
395   //  Reset NbOfInteractionLengthLeft and return aParticleChange
396   return G4VDiscreteProcess::PostStepDoIt( aTrack, aStep);
397}
398
399//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
400
401G4bool G4ComptonScattering52::StorePhysicsTable(const G4ParticleDefinition* particle,
402                                              const G4String& directory,
403                                                    G4bool          ascii)
404{
405  G4String filename;
406
407  // store cross section table
408  filename = GetPhysicsTableFileName(particle,directory,"CrossSection",ascii);
409  if ( !theCrossSectionTable->StorePhysicsTable(filename, ascii) ){
410    G4cout << " FAIL theCrossSectionTable->StorePhysicsTable in " << filename
411           << G4endl;
412    return false;
413  }
414
415  // store mean free path table
416  filename = GetPhysicsTableFileName(particle,directory,"MeanFreePath",ascii);
417  if ( !theMeanFreePathTable->StorePhysicsTable(filename, ascii) ){
418    G4cout << " FAIL theMeanFreePathTable->StorePhysicsTable in " << filename
419           << G4endl;
420    return false;
421  }
422
423  G4cout << GetProcessName() << " for " << particle->GetParticleName()
424         << ": Success to store the PhysicsTables in "
425         << directory << G4endl;
426  return true;
427}
428
429//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
430/*
431G4bool G4ComptonScattering52::RetrievePhysicsTable(const G4ParticleDefinition* particle,
432                                                 const G4String& directory,
433                                                 G4bool          ascii)
434{
435  // delete theCrossSectionTable and theMeanFreePathTable
436  if (theCrossSectionTable != 0) {
437    theCrossSectionTable->clearAndDestroy();
438    delete theCrossSectionTable;
439  }
440  if (theMeanFreePathTable != 0) {
441    theMeanFreePathTable->clearAndDestroy();
442    delete theMeanFreePathTable;
443  }
444
445  G4String filename;
446
447  // retreive cross section table
448  filename = GetPhysicsTableFileName(particle,directory,"CrossSection",ascii);
449  theCrossSectionTable = new G4PhysicsTable(G4Element::GetNumberOfElements());
450  if ( !G4PhysicsTableHelper::RetrievePhysicsTable(filename, ascii) ){
451    G4cout << " FAIL theCrossSectionTable->RetrievePhysicsTable in " << filename
452           << G4endl;
453    return false;
454  }
455
456  // retreive mean free path table
457  filename = GetPhysicsTableFileName(particle,directory,"MeanFreePath",ascii);
458  theMeanFreePathTable = new G4PhysicsTable(G4Material::GetNumberOfMaterials());
459  if ( !G4PhysicsTableHelper::RetrievePhysicsTable(filename, ascii) ){
460    G4cout << " FAIL theMeanFreePathTable->RetrievePhysicsTable in " << filename
461           << G4endl;
462    return false;
463  }
464
465  G4cout << GetProcessName() << " for " << particle->GetParticleName()
466         << ": Success to retrieve the PhysicsTables from "
467         << directory << G4endl;
468  return true;
469}
470*/
471//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
472
473void G4ComptonScattering52::PrintInfoDefinition()
474{
475  G4String comments = "Total cross sections from a parametrisation. ";
476           comments += "Good description from 10 KeV to (100/Z) GeV. \n";
477           comments += "       Scattered gamma energy according Klein-Nishina.";
478                     
479  G4cout << G4endl << GetProcessName() << ":  " << comments
480         << "\n        PhysicsTables from "
481                   << G4BestUnit(LowestEnergyLimit,"Energy")
482         << " to " << G4BestUnit(HighestEnergyLimit,"Energy") 
483         << " in " << NumbBinTable << " bins. \n";
484  G4cout << "        WARNING: This process is obsolete and will be soon removed" 
485         << G4endl;
486}         
487
488//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
Note: See TracBrowser for help on using the repository browser.