Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4eLowEnergyLoss.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 //
27 // $Id$
28 // GEANT4 tag $Name: geant4-09-01-ref-00 $
29 //
30 // -------------------------------------------------------------------
31 
32 // Class description:
33 // Low Energy electromagnetic process, electron energy loss
34 // Further documentation available from http://www.ge.infn.it/geant4/lowE
35 
36 // -------------------------------------------------------------------
37 //
38 // This class is the implementation of the unified Energy Loss process.
39 // It calculates the continuous energy loss for e+/e-.
40 // The following processes give contributions to the continuous
41 // energy loss (by default) :
42 // --- ionisation (= cont.ion.loss + delta ray production)
43 // --- bremsstrahlung (= cont.loss due to soft brems+discrete bremsstrahlung)
44 // more can be added ..........
45 // This class creates static dE/dx and range tables for e+ and e-,
46 // which tables can be used by other processes , too.
47 // G4eLowEnergyLoss is the base class for the processes giving contribution
48 // to the (continuous) energy loss of e+/e- .
49 //
50 // History: first implementation, based on object model of
51 // 2nd December 1995, G.Cosmo
52 // ---------- G4eLowEnergyLoss physics process -----------
53 // by Laszlo Urban, 20 March 1997
54 //
55 // 27.05.98 OldGetRange removed + other corrs , L.Urban
56 // 10.09.98 cleanup
57 // 16.10.98 public method SetStepFunction() + messenger class
58 // 20.01.99 new data members , L.Urban
59 // 10.02.00 modifications, new e.m. structure , L.Urban
60 // 18.10.01 Revision to improve code quality and consistency with design
61 // 23.11.01 V.Ivanchenko Move static member-functions from header to source
62 // 28.03.02 V.Ivanchenko add fluorescence flag
63 // 21.01.03 V.Ivanchenko cut per region
64 // ------------------------------------------------------------
65 
66 #ifndef G4RDeLowEnergyLoss_h
67 #define G4RDeLowEnergyLoss_h 1
68 
69 #include "G4ios.hh"
70 #include "globals.hh"
71 #include "Randomize.hh"
72 #include "G4RDVeLowEnergyLoss.hh"
73 #include "G4Material.hh"
74 #include "G4Element.hh"
76 #include "globals.hh"
77 #include "G4Track.hh"
78 #include "G4Step.hh"
79 #include "G4Electron.hh"
80 #include "G4Positron.hh"
81 #include "G4PhysicsLogVector.hh"
82 #include "G4PhysicsLinearVector.hh"
83 #include "G4EnergyLossTables.hh"
84 
85 class G4EnergyLossMessenger;
86 
88 
89 {
90  public:
91 
92  G4eLowEnergyLoss(const G4String& );
93 
95 
97  // true for e+/e- , false otherwise
98 
99  void BuildDEDXTable(const G4ParticleDefinition& aParticleType);
100  // It builds dE/dx and range tables for aParticleType and
101  // for every material contained in the materialtable.
102 
104  G4double previousStepSize,
105  G4double currentMinimumStep,
106  G4double& currentSafety);
107  // Computes the steplimit due to the energy loss process.
108 
110  const G4Step& Step) ;
111  // Performs the computation of the (continuous) energy loss
112  // after the step (with fluctuation).
113 
114  virtual G4double GetMeanFreePath(const G4Track& track,
115  G4double previousStepSize,
117  // Virtual function to be overridden in the derived classes
118  // ( ionisation and bremsstrahlung) .
119 
120  virtual G4VParticleChange* PostStepDoIt(const G4Track& track,
121  const G4Step& step) = 0;
122  // Virtual function to be overridden in the derived classes
123  // ( ionisation and bremsstrahlung) .
124 
125  static void SetNbOfProcesses(G4int nb);
126  // Sets number of processes giving contribution to the energy loss
127 
128  static void PlusNbOfProcesses();
129  // Increases number of processes giving contribution to the energy loss
130 
131  static void MinusNbOfProcesses();
132  // Decreases number of processes giving contribution to the energy loss
133 
134  static G4int GetNbOfProcesses();
135  // Gets number of processes giving contribution to the energy loss
136  // ( default value = 2)
137 
138  static void SetLowerBoundEloss(G4double val);
139  static void SetUpperBoundEloss(G4double val);
140  static void SetNbinEloss(G4int nb);
141 
142  static G4double GetLowerBoundEloss();
143  static G4double GetUpperBoundEloss();
144  static G4int GetNbinEloss();
145 
146  void ActivateFluorescence(G4bool val);
147  // Set fluorescence flag on/off
148 
149  G4bool Fluorescence() const;
150  // Get flurescence flag
151 
152  protected:
153 
154  virtual std::vector<G4DynamicParticle*>* DeexciteAtom(const G4MaterialCutsCouple* ,
155  G4double, G4double) // incidentEnergy, eLoss
156  { return 0; };
157 
159 
160  G4double MinKineticEnergy ; // particle with kinetic energy
161  // smaller than MinKineticEnergy
162  // is stopped in AlongStepDoIt
163 
165 
166  //basic DEDX and Range tables
171 
172  //inverse tables of the range tables
175 
176  // lab and proper time tables
181 
182  //processes inherited from G4eLowEnergyLoss
183  //register themselves in the static array Recorder
184  //for electrons/positrons separately
185  //nb of contributing processes = NbOfProcesses
191 
192 
193  private:
194 
195  G4double GetConstraints(const G4DynamicParticle* aParticle,
196  const G4MaterialCutsCouple* couple);
197 
198  // hide assignment operator
200  G4eLowEnergyLoss & operator=(const G4eLowEnergyLoss &right);
201 
202 
203  G4PhysicsTable* theDEDXTable;
204 
205  G4int CounterOfProcess;
206  G4PhysicsTable** RecorderOfProcess;
207 
208  G4double fdEdx; // computed in GetConstraints
209  G4double fRangeNow; // computed in GetConstraints
210 
211  G4double linLossLimit ; // used in AlongStepDoIt
212 
213 
214  //New ParticleChange
215  G4ParticleChangeForLoss fParticleChange ;
216 
217  //
218  // static part of the class
219  //
220 
221  static G4int NbinEloss; // number of bins in table,
222  // calculated in BuildPhysicTable
223  static G4double LowerBoundEloss;
224  static G4double UpperBoundEloss;
225  static G4double RTable,LOGRTable; // LOGRTable=std::log(UpperBoundEloss-
226  // LowerBoundEloss)/NbinEloss
227  // RTable = std::exp(LOGRTable)
228 
229  //for interpolation within the tables
230  static G4PhysicsTable* theeRangeCoeffATable;
231  static G4PhysicsTable* theeRangeCoeffBTable;
232  static G4PhysicsTable* theeRangeCoeffCTable;
233  static G4PhysicsTable* thepRangeCoeffATable;
234  static G4PhysicsTable* thepRangeCoeffBTable;
235  static G4PhysicsTable* thepRangeCoeffCTable;
236 
237  static G4EnergyLossMessenger* eLossMessenger;
238 
239  G4bool theFluo; // Fluorescence flag
240 
241 };
242 
243 #include "G4eLowEnergyLoss.icc"
244 
245 #endif
246 
G4double condition(const G4ErrorSymMatrix &m)
static void SetLowerBoundEloss(G4double val)
static G4PhysicsTable * theInverseRangePositronTable
static G4PhysicsTable * theProperTimePositronTable
static G4int CounterOfElectronProcess
G4bool Fluorescence() const
static G4int NbOfProcesses
void ActivateFluorescence(G4bool val)
virtual G4double GetMeanFreePath(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)=0
static G4int GetNbOfProcesses()
static G4PhysicsTable * theRangeElectronTable
static G4double GetLowerBoundEloss()
static void SetUpperBoundEloss(G4double val)
G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &Step)
int G4int
Definition: G4Types.hh:78
G4bool IsApplicable(const G4ParticleDefinition &)
static G4PhysicsTable * theLabTimePositronTable
virtual G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &step)=0
bool G4bool
Definition: G4Types.hh:79
static G4PhysicsTable * theDEDXElectronTable
void BuildDEDXTable(const G4ParticleDefinition &aParticleType)
G4eLowEnergyLoss(const G4String &)
Definition: G4Step.hh:76
static G4PhysicsTable * theRangePositronTable
static void MinusNbOfProcesses()
static void SetNbinEloss(G4int nb)
static G4PhysicsTable ** RecorderOfPositronProcess
static G4PhysicsTable * theInverseRangeElectronTable
static G4double GetUpperBoundEloss()
static G4PhysicsTable ** RecorderOfElectronProcess
static G4PhysicsTable * theLabTimeElectronTable
G4double GetContinuousStepLimit(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety)
static G4int GetNbinEloss()
static G4PhysicsTable * theProperTimeElectronTable
Definition: Step.hh:41
static void PlusNbOfProcesses()
G4PhysicsTable * theLossTable
static G4int CounterOfPositronProcess
static void SetNbOfProcesses(G4int nb)
double G4double
Definition: G4Types.hh:76
static G4PhysicsTable * theDEDXPositronTable
G4ForceCondition
virtual std::vector< G4DynamicParticle * > * DeexciteAtom(const G4MaterialCutsCouple *, G4double, G4double)