Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4RDVeLowEnergyLoss.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 // 3.4.2000 Veronique Lefebure:
31 // Move utils/include/G4VEnergyLoss.hh to
32 // lowenergy/include/G4RDVeLowEnergyLoss.hh
33 //
34 // ------------------------------------------------------------
35 // GEANT 4 class header file
36 //
37 //
38 // Class Description
39 //
40 // General service class for the energy loss classes
41 //
42 // It contains code needed to compute the range tables,
43 // time tables, the inverse range tables and some auxiliary
44 // tables.
45 // The energy loss fluctuation code is here,too.
46 //
47 // All the EnergyLoss classes are inherited from G4RDVeLowEnergyLoss
48 // class.
49 //
50 // -----------------------------------------------------------
51 // created on 28 January 2000 by L. Urban
52 // -----------------------------------------------------------
53 //
54 // Modifications:
55 // 20/09/00 V.Ivanchenko update fluctuations
56 // 23/11/01 V.Ivanchenko Move static member-functions from header to source
57 // 22/01/03 V.Ivanchenko Cut per range
58 //
59 // Class description:
60 // Abstract class for Low Energy Electromagnetic electron energy loss
61 // Further documentation available from http://www.ge.infn.it/geant4/lowE
62 
63 // -----------------------------------------------------------
64 
65 #ifndef G4RDVeLowEnergyLoss_h
66 #define G4RDVeLowEnergyLoss_h 1
67 
68 #include "globals.hh"
69 #include "G4ios.hh"
70 #include "Randomize.hh"
71 #include "G4Poisson.hh"
72 #include "G4Electron.hh"
74 #include "G4PhysicsLogVector.hh"
75 #include "G4PhysicsLinearVector.hh"
76 #include "G4MaterialCutsCouple.hh"
77 
79 {
80  public:
81 
83  G4ProcessType aType = fNotDefined );
85 
86  virtual ~G4RDVeLowEnergyLoss();
87 
88  virtual G4double GetContinuousStepLimit(const G4Track& track,
89  G4double previousStepSize,
90  G4double currentMinimumStep,
91  G4double& currentSafety) = 0 ;
92 
93  virtual G4VParticleChange* AlongStepDoIt(const G4Track& track,
94  const G4Step& Step) = 0 ;
95 
96  virtual G4double GetMeanFreePath(const G4Track& track,
97  G4double previousStepSize,
99 
100  virtual G4VParticleChange* PostStepDoIt(const G4Track& track,
101  const G4Step& Step) = 0;
102 
103 
104 
105  protected:// with description
106 
107  // code for the energy loss fluctuation
108 
110  const G4MaterialCutsCouple* couple,
111  G4double MeanLoss,
112  G4double step);
113 
114 
115  private:
116 
117  // hide default constructor and assignment operator as private
119  G4RDVeLowEnergyLoss & operator=(const G4RDVeLowEnergyLoss &right);
120 
121  protected:
122 
123  // data members to speed up the fluctuation calculation
128 
130 
131  // static part of the class
132 
133  public: // With description
134 
135  static void SetRndmStep (G4bool value);
136  // use / do not use randomisation in energy loss steplimit
137  // ( default = no randomisation)
138 
139  static void SetEnlossFluc (G4bool value);
140  // compute energy loss with/without fluctuation
141  // ( default : with fluctuation)
142 
143  static void SetStepFunction (G4double c1, G4double c2);
144  // sets values for data members used to compute the step limit:
145  // dRoverRange : max. relative range change in one step,
146  // finalRange : if range <= finalRange --> last step for the particle.
147 
148 
149  protected: // With description
150 
151  // Build range table starting from the DEDXtable
152  static G4PhysicsTable*
153  BuildRangeTable(G4PhysicsTable* theDEDXTable,
154  G4PhysicsTable* theRangeTable,
155  G4double Tmin,G4double Tmax,G4int nbin);
156 
157  // Build time tables starting from the DEDXtable
158  static G4PhysicsTable*
159  BuildLabTimeTable(G4PhysicsTable* theDEDXTable,
160  G4PhysicsTable* theLabTimeTable,
161  G4double Tmin,G4double Tmax,G4int nbin);
162 
163  static G4PhysicsTable*
164  BuildProperTimeTable(G4PhysicsTable* theDEDXTable,
165  G4PhysicsTable* ProperTimeTable,
166  G4double Tmin,G4double Tmax,G4int nbin);
167 
168  // Build tables of coefficients needed for inverting the range table
169  static G4PhysicsTable*
170  BuildRangeCoeffATable(G4PhysicsTable* theRangeTable,
171  G4PhysicsTable* theCoeffATable,
172  G4double Tmin,G4double Tmax,G4int nbin);
173  static G4PhysicsTable*
174  BuildRangeCoeffBTable(G4PhysicsTable* theRangeTable,
175  G4PhysicsTable* theCoeffBTable,
176  G4double Tmin,G4double Tmax,G4int nbin);
177  static G4PhysicsTable*
178  BuildRangeCoeffCTable(G4PhysicsTable* theRangeTable,
179  G4PhysicsTable* theCoeffCTable,
180  G4double Tmin,G4double Tmax,G4int nbin);
181 
182  // Invert range table
183  static G4PhysicsTable*
184  BuildInverseRangeTable(G4PhysicsTable* theRangeTable,
185  G4PhysicsTable* theRangeCoeffATable,
186  G4PhysicsTable* theRangeCoeffBTable,
187  G4PhysicsTable* theRangeCoeffCTable,
188  G4PhysicsTable* theInverseRangeTable,
189  G4double Tmin,G4double Tmax,G4int nbin);
190 
191  private:
192 
193  static void BuildRangeVector(G4PhysicsTable* theDEDXTable,
194  G4double Tmin,G4double Tmax,G4int nbin,
195  G4int materialIndex,G4PhysicsLogVector* rangeVector);
196 
197  static void BuildRangeVectorNew(const G4PhysicsTable*,G4int,
199 
200  static G4double RangeIntLin(G4PhysicsVector* physicsVector
201  ,G4int nbin);
202 
203  static G4double RangeIntLog(G4PhysicsVector* physicsVector
204  ,G4int nbin);
205 
206  static void BuildLabTimeVector(G4PhysicsTable* theDEDXTable,
207  G4double Tmin,G4double Tmax,G4int nbin,
208  G4int materialIndex,G4PhysicsLogVector* rangeVector);
209 
210  static void BuildProperTimeVector(G4PhysicsTable* theDEDXTable,
211  G4double Tmin,G4double Tmax,G4int nbin,
212  G4int materialIndex,G4PhysicsLogVector* rangeVector);
213 
214  static G4double LabTimeIntLog(G4PhysicsVector* physicsVector
215  ,G4int nbin);
216 
217  static G4double ProperTimeIntLog(G4PhysicsVector* physicsVector,
218  G4int nbin);
219 
220  static void InvertRangeVector(G4PhysicsTable* theRangeTable,
221  G4PhysicsTable* theRangeCoeffATable,
222  G4PhysicsTable* theRangeCoeffBTable,
223  G4PhysicsTable* theRangeCoeffCTable,
224  G4double Tmin,G4double Tmax,G4int nbin,
225  G4int materialIndex,G4PhysicsLogVector* rangeVector);
226 
227 
228  // data members
229  protected:
230 
231  // variables for the integration routines
233 
234 
235  static G4double dRoverRange; // dRoverRange is the maximum allowed
236  // deltarange/range in one Step
237  static G4double finalRange; // final step before stopping
238  static G4double c1lim,c2lim,c3lim ; // coeffs for computing steplimit
239 
240  static G4bool rndmStepFlag; // control the randomization of the step
241  static G4bool EnlossFlucFlag; // control the energy loss fluctuation
242 
243 
244 };
245 
246 #endif
247 
248 
249 
G4double condition(const G4ErrorSymMatrix &m)
static c2_factory< G4double > c2
static G4double finalRange
virtual G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &Step)=0
static G4double dRoverRange
static G4double ParticleMass
int G4int
Definition: G4Types.hh:78
static G4PhysicsTable * BuildRangeTable(G4PhysicsTable *theDEDXTable, G4PhysicsTable *theRangeTable, G4double Tmin, G4double Tmax, G4int nbin)
static G4PhysicsTable * BuildProperTimeTable(G4PhysicsTable *theDEDXTable, G4PhysicsTable *ProperTimeTable, G4double Tmin, G4double Tmax, G4int nbin)
virtual G4double GetMeanFreePath(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)=0
static G4PhysicsTable * BuildLabTimeTable(G4PhysicsTable *theDEDXTable, G4PhysicsTable *theLabTimeTable, G4double Tmin, G4double Tmax, G4int nbin)
const XML_Char int const XML_Char * value
Definition: expat.h:331
bool G4bool
Definition: G4Types.hh:79
static G4PhysicsTable * BuildRangeCoeffBTable(G4PhysicsTable *theRangeTable, G4PhysicsTable *theCoeffBTable, G4double Tmin, G4double Tmax, G4int nbin)
Definition: G4Step.hh:76
static void SetRndmStep(G4bool value)
static void SetStepFunction(G4double c1, G4double c2)
static G4PhysicsTable * BuildInverseRangeTable(G4PhysicsTable *theRangeTable, G4PhysicsTable *theRangeCoeffATable, G4PhysicsTable *theRangeCoeffBTable, G4PhysicsTable *theRangeCoeffCTable, G4PhysicsTable *theInverseRangeTable, G4double Tmin, G4double Tmax, G4int nbin)
Definition: Step.hh:41
double G4double
Definition: G4Types.hh:76
virtual G4double GetContinuousStepLimit(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety)=0
virtual G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &Step)=0
G4ForceCondition
G4double GetLossWithFluct(const G4DynamicParticle *aParticle, const G4MaterialCutsCouple *couple, G4double MeanLoss, G4double step)
static G4PhysicsTable * BuildRangeCoeffATable(G4PhysicsTable *theRangeTable, G4PhysicsTable *theCoeffATable, G4double Tmin, G4double Tmax, G4int nbin)
const G4Material * lastMaterial
tuple c1
Definition: plottest35.py:14
static void SetEnlossFluc(G4bool value)
static G4PhysicsTable * BuildRangeCoeffCTable(G4PhysicsTable *theRangeTable, G4PhysicsTable *theCoeffCTable, G4double Tmin, G4double Tmax, G4int nbin)
G4ProcessType