source: trunk/source/processes/electromagnetic/lowenergy/include/G4eLowEnergyLoss.icc@ 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.6 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//
27// $Id: G4eLowEnergyLoss.icc,v 1.8 2006/06/29 19:37:44 gunter Exp $
28// GEANT4 tag $Name: geant4-09-03-cand-01 $
29//
30// ---------------------------------------------------------------
31// GEANT 4 class inlined methods file
32//
33// History: based on object model of
34// 2nd December 1995, G.Cosmo
35// ------------ G4eLowEnergyLoss physics process ------------
36// by Laszlo Urban, 20 March 1997
37// ***************************************************************
38// It is the first implementation of the NEW UNIFIED ENERGY LOSS PROCESS.
39// It calculates the energy loss of e+/e-.
40// -------------------------------------------------------------
41//
42// 08-09-98 cleanup
43// 28-03-02 V.Ivanchenko add fluorescence flag
44// 21-01-03 V.Ivanchenko cut per region
45// 18-04-03 V.Ivanchenko finalRange redefinition
46//
47// ---------------------------------------------------------------
48
49inline G4double G4eLowEnergyLoss::GetConstraints(
50 const G4DynamicParticle* aParticle,
51 const G4MaterialCutsCouple* couple)
52{
53 G4double StepLimit;
54 // returns the Step limit
55 // dRoverRange is the max. allowed relative range loss in one Step
56 // it calculates dEdx and the range as well....
57
58 const G4ParticleDefinition* ParticleType=aParticle->GetDefinition();
59
60 Charge = aParticle->GetDefinition()->GetPDGCharge();
61 if(Charge != lastCharge) lastCharge = Charge ;
62
63 G4double KineticEnergy = aParticle->GetKineticEnergy();
64
65 fdEdx = G4EnergyLossTables::GetDEDX(ParticleType,KineticEnergy,couple);
66 fRangeNow =
67 G4EnergyLossTables::GetRange(ParticleType,KineticEnergy,couple);
68
69 G4double r = std::min(finalRange, couple->GetProductionCuts()
70 ->GetProductionCut(idxG4ElectronCut));
71
72 if (fRangeNow > r) {
73 StepLimit = dRoverRange*fRangeNow + r*(1.0 - dRoverRange)*(2.0 - r/fRangeNow);
74 //randomise this value
75 if (rndmStepFlag) StepLimit = finalRange + (StepLimit-finalRange)*G4UniformRand();
76 if (StepLimit > fRangeNow) StepLimit = fRangeNow;
77 }
78 else StepLimit = fRangeNow;
79
80 return StepLimit;
81}
82
83//
84
85inline G4double G4eLowEnergyLoss::GetContinuousStepLimit(
86 const G4Track& track,
87 G4double,
88 G4double currentMinimumStep,
89 G4double&)
90{
91
92 G4double Step =
93 GetConstraints(track.GetDynamicParticle(),track.GetMaterialCutsCouple());
94
95 if ((Step>0.0)&&(Step<currentMinimumStep)) currentMinimumStep = Step;
96
97 return Step ;
98}
99
100//
101inline G4bool G4eLowEnergyLoss::IsApplicable(const G4ParticleDefinition&
102 particle)
103{
104 return( (&particle == G4Electron::Electron())
105 ||(&particle == G4Positron::Positron()) );
106}
107
108inline void G4eLowEnergyLoss::ActivateFluorescence(G4bool val)
109{
110 theFluo = val;
111}
112
113inline G4bool G4eLowEnergyLoss::Fluorescence() const
114{
115 return theFluo;
116}
117
118//
119
Note: See TracBrowser for help on using the repository browser.