Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4VLowEnergyDiscretePhotonProcess.hh
Go to the documentation of this file.
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$
27 //
28 // --------------------------------------------------------------
29 //
30 // File name: G4VLowEnergyDiscretePhotonProcess.hh
31 //
32 // Author: Capra Riccardo
33 //
34 // Creation date: May 2005
35 //
36 // History:
37 // -----------
38 // 02 May 2005 R. Capra 1st implementation
39 // 20 May 2005 MGP Modified protected member functions to return a const pointer
40 //
41 //----------------------------------------------------------------
42 
43 
44 
45 #ifndef G4VLowEnergyDiscretePhotonProcess_hh
46 #define G4VLowEnergyDiscretePhotonProcess_hh
47 
48 // Base class
50 
51 // Forward declaration
52 class G4String;
54 class G4DynamicParticle;
55 class G4Track;
57 class G4VEMDataSet;
59 
60 // G4VLowEnergyDiscretePhotonProcess
61 // A common class for Rayleigh and Compton processes
63 {
64 public:
65  // Class constructor
66  //
67  // Creates crossSectionHandler and scatterFunctionData
68  // scatterFunctionData are loaded from the scatterFile file, and are interpolated with scatterInterpolation algorithm
69  // processName The name of the process
70  // aCrossSectionFileName The name of the cross-section data file
71  // aScatterFileName The name of the scatter function data
72  // aScatterInterpolation The interpolation algorithm
73  // aLowEnergyLimit The lower energy limit of validity against data
74  // aHighEnergyLimit The higher energy limit of validity against data
75  G4VLowEnergyDiscretePhotonProcess(const G4String& processName,
76  const G4String& aCrossSectionFileName,
77  const G4String& aScatterFileName,
78  G4VDataSetAlgorithm* aScatterInterpolation,
79  G4double aLowEnergyLimit,
80  G4double aHighEnergyLimit);
81 
82  // Class destructor
83  //
84  // Deletes crossSectionHandler, meanFreePathTable and scatterFunctionData
86 
87 
88 
89  // Checks if the process is applicable to the particle
90  //
91  // For processes inheriting from this class the only applicable particle is the photon
92  // particleDefinition Is the particle to be checked for the applicability of the process
93  // true only if the particle is a photon
94  virtual G4bool IsApplicable(const G4ParticleDefinition& particleDefinition);
95 
96  // Updates the crossSectionHandler and builds the meanFreePathTable from it
97  //
98  // crossSectionHandler data is loaded from crossSectionFile file
99  // photon particle is always a photon
100  virtual void BuildPhysicsTable(const G4ParticleDefinition& photon);
101 
102 
103 
104  // Returns the low energy limit of the process
105  // The low energy limit of the process
106  inline G4double GetLowEnergyLimit(void) const {return lowEnergyLimit;}
107 
108  // Returns the high energy limit of the process
109  // The high energy limit of the process
110  inline G4double GetHighEnergyLimit(void) const {return highEnergyLimit;}
111 
112 
113 
114 protected:
115  // Evaluates the process mean free path
116  //
117  // Mean free path evaluation is based on meanFreePathTable generated by the crossSectionHandler with data taken from the crossSectionFile.
118  // The method uses lowEnergyLimit and highEnergyLimit
119  // aTrack the particle momentum for which the mean free path must be evaluates (a photon)
120  // previousStepSize the size of the prevous step (not used by this implementation)
121  // condition the consition to be updated (not used by this implementation)
122  // The mean free path for the process
123  virtual G4double GetMeanFreePath(const G4Track &aTrack,
124  G4double previousStepSize,
126 
127  // Returns the cross-section handler
128  // The cross-section handler
129  inline const G4VCrossSectionHandler* GetCrossSectionHandler(void) const {return crossSectionHandler;}
130 
131  // Returns the mean free path table
132  // The mean free path table
133  inline const G4VEMDataSet* GetMeanFreePathTable(void) const {return meanFreePathTable;}
134 
135  // Returns the scatter function data
136  // The scatter function data
137  inline const G4VEMDataSet* GetScatterFunctionData(void) const {return scatterFunctionData;}
138 
139 
140 
141  // Verifies if the polarization vector is orthonormal to the photon direction
142  //
143  // This method is used by polarized processes. It returns always a vector orthonormal to the photon direction.
144  // When G4DynamicParticle::GetPolarization() is well defined the returned vector matches this vector.
145  // When G4DynamicParticle::GetPolarization() is sensibly different from a null-vector, the returned vector is
146  // orthonormal to the photon direction and on the plane made with G4DynamicParticle::GetPolarization()
147  // When G4DynamicParticle::GetPolarization() is almost a null-vector, the returned vector is a random
148  // vector orthonormal to the photon direction
149  // photon The incoming photon
150  // s Always a vector ortogonal to the photon direction.
152 
153 
154 
155 private:
156 
157  // Hides copy constructor
159 
160  // Hides assignment operator
162 
163  // Low energy limit
164  G4double lowEnergyLimit;
165 
166  // High energy limit
167  G4double highEnergyLimit;
168 
169  // Cross-section file name
170  G4String crossSectionFileName;
171 
172  // Cross-section handler
173  G4VCrossSectionHandler* crossSectionHandler;
174 
175  // The mean free path table based on cross-section data
176  G4VEMDataSet* meanFreePathTable;
177 
178  // Scatter function table
179  G4VEMDataSet* scatterFunctionData;
180 };
181 
182 #endif /* G4VLowEnergyDiscretePhotonProcess_hh */