source: trunk/source/processes/electromagnetic/standard/include/G4BetheBlochModel.hh @ 1007

Last change on this file since 1007 was 1007, checked in by garnier, 15 years ago

update to geant4.9.2

File size: 5.9 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// $Id: G4BetheBlochModel.hh,v 1.16 2008/10/22 16:00:57 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-02 $
28//
29// -------------------------------------------------------------------
30//
31// GEANT4 Class header file
32//
33//
34// File name:     G4BetheBlochModel
35//
36// Author:        Vladimir Ivanchenko on base of Laszlo Urban code
37//
38// Creation date: 03.01.2002
39//
40// Modifications:
41//
42// 23-12-02 V.Ivanchenko change interface in order to moveto cut per region
43// 24-01-03 Make models region aware (V.Ivanchenko)
44// 13-02-03 Add name (V.Ivanchenko)
45// 12-11-03 Fix for GenericIons (V.Ivanchenko)
46// 24-03-05 Add G4EmCorrections (V.Ivanchenko)
47// 11-04-05 Major optimisation of internal interfaces (V.Ivanchenko)
48// 11-04-04 Move MaxSecondaryEnergy to models (V.Ivanchenko)
49// 11-02-06 ComputeCrossSectionPerElectron, ComputeCrossSectionPerAtom (mma)
50// 12-08-08 Added methods GetParticleCharge, GetChargeSquareRatio,
51//          CorrectionsAlongStep needed for ions(V.Ivanchenko)
52
53//
54// Class Description:
55//
56// Implementation of Bethe-Bloch model of energy loss and
57// delta-electron production by heavy charged particles
58
59// -------------------------------------------------------------------
60//
61
62#ifndef G4BetheBlochModel_h
63#define G4BetheBlochModel_h 1
64
65#include "G4VEmModel.hh"
66
67class G4EmCorrections;
68class G4ParticleChangeForLoss;
69class G4NistManager;
70
71
72class G4BetheBlochModel : public G4VEmModel
73{
74
75public:
76
77  G4BetheBlochModel(const G4ParticleDefinition* p = 0,
78                    const G4String& nam = "BetheBloch");
79
80  virtual ~G4BetheBlochModel();
81
82  virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
83
84  G4double MinEnergyCut(const G4ParticleDefinition*,
85                        const G4MaterialCutsCouple*);
86                       
87  virtual G4double ComputeCrossSectionPerElectron(
88                                 const G4ParticleDefinition*,
89                                 G4double kineticEnergy,
90                                 G4double cutEnergy,
91                                 G4double maxEnergy);
92                                 
93  virtual G4double ComputeCrossSectionPerAtom(
94                                 const G4ParticleDefinition*,
95                                 G4double kineticEnergy,
96                                 G4double Z, G4double A,
97                                 G4double cutEnergy,
98                                 G4double maxEnergy);
99                                                                 
100  virtual G4double CrossSectionPerVolume(const G4Material*,
101                                 const G4ParticleDefinition*,
102                                 G4double kineticEnergy,
103                                 G4double cutEnergy,
104                                 G4double maxEnergy);
105                                 
106  virtual G4double ComputeDEDXPerVolume(const G4Material*,
107                                        const G4ParticleDefinition*,
108                                        G4double kineticEnergy,
109                                        G4double cutEnergy);
110
111  virtual G4double GetChargeSquareRatio(const G4ParticleDefinition* p,
112                                        const G4Material* mat,
113                                        G4double kineticEnergy);
114
115  virtual G4double GetParticleCharge(const G4ParticleDefinition* p,
116                                     const G4Material* mat,
117                                     G4double kineticEnergy);
118
119  virtual void CorrectionsAlongStep(const G4MaterialCutsCouple*,
120                                    const G4DynamicParticle*,
121                                    G4double& eloss,
122                                    G4double& niel,
123                                    G4double length);
124
125  virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
126                                 const G4MaterialCutsCouple*,
127                                 const G4DynamicParticle*,
128                                 G4double tmin,
129                                 G4double maxEnergy);
130
131protected:
132
133  G4double MaxSecondaryEnergy(const G4ParticleDefinition*,
134                              G4double kinEnergy);
135
136private:
137
138  void SetParticle(const G4ParticleDefinition* p);
139
140  // hide assignment operator
141  G4BetheBlochModel & operator=(const  G4BetheBlochModel &right);
142  G4BetheBlochModel(const  G4BetheBlochModel&);
143
144  const G4ParticleDefinition* particle;
145  const G4Material*           currentMaterial;
146  G4ParticleDefinition*       theElectron;
147  G4EmCorrections*            corr;
148  G4ParticleChangeForLoss*    fParticleChange;
149  G4NistManager*              nist;
150
151  G4double mass;
152  G4double tlimit;
153  G4double spin;
154  G4double magMoment2;
155  G4double chargeSquare;
156  G4double ratio;
157  G4double formfact;
158  G4double twoln10;
159  G4double bg2lim;
160  G4double taulim;
161  G4double corrFactor;
162  G4bool   isIon;
163  G4bool   isInitialised;
164};
165
166//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
167
168inline G4double G4BetheBlochModel::MaxSecondaryEnergy(
169          const G4ParticleDefinition* pd,
170                G4double kinEnergy) 
171{
172  if(isIon) SetParticle(pd);
173  G4double tau  = kinEnergy/mass;
174  G4double tmax = 2.0*electron_mass_c2*tau*(tau + 2.) /
175                  (1. + 2.0*(tau + 1.)*ratio + ratio*ratio);
176  return std::min(tmax,tlimit);
177}
178
179//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
180
181#endif
Note: See TracBrowser for help on using the repository browser.