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

Last change on this file since 1201 was 1196, checked in by garnier, 16 years ago

update CVS release candidate geant4.9.3.01

File size: 4.7 KB
RevLine 
[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//
[961]26// $Id: G4PenelopeCompton.hh,v 1.11 2008/03/26 15:30:19 pandola Exp $
[1196]27// GEANT4 tag $Name: geant4-09-03-cand-01 $
[819]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)
[961]38// 26 Mar 2008 L.Pandola Add boolean flag to control atomic de-excitation
[819]39//
40// -------------------------------------------------------------------
41
42// Class description:
43// Low Energy Electromagnetic Physics, Compton Scattering
44// Penelope Model
45// -------------------------------------------------------------------
46
47#ifndef G4PENELOPECOMPTON_HH
48#define G4PENELOPECOMPTON_HH 1
49
50#include "globals.hh"
51#include "G4VDiscreteProcess.hh"
52
53class G4Track;
54class G4Step;
55class G4ParticleDefinition;
56class G4VParticleChange;
57class G4VEMDataSet;
58class G4VRangeTest;
59class G4Material;
60
61class G4PenelopeCompton : public G4VDiscreteProcess {
62
63public:
64
65 G4PenelopeCompton(const G4String& processName ="PenCompton");
66
67 ~G4PenelopeCompton();
68
69 G4bool IsApplicable(const G4ParticleDefinition& definition);
70
71 void BuildPhysicsTable(const G4ParticleDefinition& photon);
72
73 G4VParticleChange* PostStepDoIt(const G4Track& aTrack, const G4Step& aStep);
74
75 // For testing purpose only
76 G4double DumpMeanFreePath(const G4Track& aTrack,
77 G4double previousStepSize,
78 G4ForceCondition* condition)
79 { return GetMeanFreePath(aTrack, previousStepSize, condition); }
80
[961]81 void SetUseAtomicDeexcitation(G4bool value){fUseAtomicDeexcitation = value;};
82
83 G4bool GetUseAtomicDeexcitation(){return fUseAtomicDeexcitation;};
84
[819]85protected:
86
87 G4double GetMeanFreePath(const G4Track& aTrack,
88 G4double previousStepSize,
89 G4ForceCondition* condition);
90
91private:
92
93 // Hide copy constructor and assignment operator as private
94 G4PenelopeCompton& operator=(const G4PenelopeCompton& right);
95 G4PenelopeCompton(const G4PenelopeCompton& );
96
97 void ReadData();
98
99 G4double CrossSection(G4double energy,G4int Z);
100 G4double DifferentialCrossSection (G4double cdt);
101
102 G4double lowEnergyLimit; // low energy limit applied to the process
103 G4double highEnergyLimit; // high energy limit applied to the process
104
105 G4VEMDataSet* meanFreePathTable;
106 G4VRangeTest* rangeTest;
107
108 const G4double intrinsicLowEnergyLimit; // intrinsic validity range
109 const G4double intrinsicHighEnergyLimit;
110
111 G4double energyForIntegration; //for numerical integration of
112 G4int ZForIntegration;// analytical cross section
113
114 //Parameters of atomic shells
115 std::map<G4int,G4DataVector*> *ionizationEnergy;
116 std::map<G4int,G4DataVector*> *hartreeFunction;
117 std::map<G4int,G4DataVector*> *occupationNumber;
118
119
120 G4int SelectRandomAtomForCompton(const G4Material* material,G4double e) const;
121
122 const G4int nBins; //for building cross section table
123
124 std::vector<G4VEMDataSet*> *matCrossSections; //for random choice of atom
[961]125
[819]126 G4double cutForLowEnergySecondaryPhotons;
[961]127 G4bool fUseAtomicDeexcitation;
[819]128};
129
130#endif
131
Note: See TracBrowser for help on using the repository browser.