source: trunk/source/processes/electromagnetic/lowenergy/include/G4PenelopeCompton.hh @ 924

Last change on this file since 924 was 819, checked in by garnier, 16 years ago

import all except CVS

File size: 4.5 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: G4PenelopeCompton.hh,v 1.10 2006/06/29 19:36:21 gunter Exp $
27// GEANT4 tag $Name: geant4-09-01-patch-02 $
28//
29// Author: Luciano Pandola
30//
31// History:
32// -----------
33// 05 Dec 2002   L. Pandola   1st implementation
34// 12 Feb 2003   MG Pia       const argument in SelectRandomAtomForCompton
35// 14 Feb 2003   MG Pia       Modified some variables to lowercase initial
36// 26 Mar 2003   L.Pandola    Added fluorescence
37// 18 Mar 2004   L.Pandola    Use of std::map (code review)
38//
39// -------------------------------------------------------------------
40
41// Class description:
42// Low Energy Electromagnetic Physics, Compton Scattering
43// Penelope Model
44// -------------------------------------------------------------------
45
46#ifndef G4PENELOPECOMPTON_HH
47#define G4PENELOPECOMPTON_HH 1
48
49#include "globals.hh"
50#include "G4VDiscreteProcess.hh"
51
52class G4Track;
53class G4Step;
54class G4ParticleDefinition;
55class G4VParticleChange;
56class G4VEMDataSet;
57class G4VRangeTest;
58class G4Material;
59
60class G4PenelopeCompton : public G4VDiscreteProcess {
61
62public:
63 
64  G4PenelopeCompton(const G4String& processName ="PenCompton");
65 
66  ~G4PenelopeCompton();
67
68  G4bool IsApplicable(const G4ParticleDefinition& definition);
69 
70  void BuildPhysicsTable(const G4ParticleDefinition& photon);
71 
72  G4VParticleChange* PostStepDoIt(const G4Track& aTrack, const G4Step& aStep);
73 
74  // For testing purpose only
75  G4double DumpMeanFreePath(const G4Track& aTrack, 
76                            G4double previousStepSize, 
77                            G4ForceCondition* condition) 
78  { return GetMeanFreePath(aTrack, previousStepSize, condition); }
79
80protected:
81
82  G4double GetMeanFreePath(const G4Track& aTrack, 
83                           G4double previousStepSize, 
84                           G4ForceCondition* condition);
85
86private: 
87
88  // Hide copy constructor and assignment operator as private
89  G4PenelopeCompton& operator=(const G4PenelopeCompton& right);
90  G4PenelopeCompton(const G4PenelopeCompton& );
91
92  void ReadData();
93
94  G4double CrossSection(G4double energy,G4int Z);
95  G4double DifferentialCrossSection (G4double cdt);
96
97  G4double lowEnergyLimit;  // low energy limit  applied to the process
98  G4double highEnergyLimit; // high energy limit applied to the process
99
100  G4VEMDataSet* meanFreePathTable;
101  G4VRangeTest* rangeTest;
102
103  const G4double intrinsicLowEnergyLimit; // intrinsic validity range
104  const G4double intrinsicHighEnergyLimit;
105
106  G4double energyForIntegration; //for numerical integration of
107  G4int ZForIntegration;// analytical cross section
108 
109  //Parameters of atomic shells
110  std::map<G4int,G4DataVector*> *ionizationEnergy;
111  std::map<G4int,G4DataVector*> *hartreeFunction;
112  std::map<G4int,G4DataVector*> *occupationNumber;
113
114 
115  G4int SelectRandomAtomForCompton(const G4Material* material,G4double e) const;
116
117  const G4int nBins; //for building cross section table
118
119  std::vector<G4VEMDataSet*> *matCrossSections; //for random choice of atom
120  G4double cutForLowEnergySecondaryPhotons;
121};
122
123#endif
124
Note: See TracBrowser for help on using the repository browser.