Geant4  10.02.p01
G4BraggIonModel.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: G4BraggIonModel.hh 93362 2015-10-19 13:45:19Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class header file
31 //
32 //
33 // File name: G4BraggIonModel
34 //
35 // Author: Vladimir Ivanchenko
36 //
37 // Creation date: 13.10.2004
38 //
39 // Modifications:
40 // 09-11-04 Migration to new interface of Store/Retrieve tables (V.Ivantchenko)
41 // 11-05-05 Major optimisation of internal interfaces (V.Ivantchenko)
42 // 15-02-06 ComputeCrossSectionPerElectron, ComputeCrossSectionPerAtom (mma)
43 // 25-04-06 Add stopping data from ASTAR (V.Ivanchenko)
44 // 12-08-08 Added methods GetParticleCharge, GetChargeSquareRatio,
45 // CorrectionsAlongStep needed for ions(V.Ivanchenko)
46 
47 //
48 // Class Description:
49 //
50 // Implementation of energy loss and delta-electron production
51 // by heavy slow charged particles using ICRU'49 and NIST evaluated data
52 // for He4 ions
53 
54 // -------------------------------------------------------------------
55 //
56 
57 #ifndef G4BraggIonModel_h
58 #define G4BraggIonModel_h 1
59 
60 #include <CLHEP/Units/PhysicalConstants.h>
61 
62 #include "G4VEmModel.hh"
63 #include "G4ASTARStopping.hh"
64 
66 class G4EmCorrections;
67 
69 {
70 
71 public:
72 
73  G4BraggIonModel(const G4ParticleDefinition* p = nullptr,
74  const G4String& nam = "BraggIon");
75 
76  virtual ~G4BraggIonModel();
77 
78  virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
79 
81  const G4MaterialCutsCouple* couple);
82 
84  const G4ParticleDefinition*,
85  G4double kineticEnergy,
86  G4double cutEnergy,
87  G4double maxEnergy);
88 
90  const G4ParticleDefinition*,
91  G4double kineticEnergy,
92  G4double Z, G4double A,
93  G4double cutEnergy,
94  G4double maxEnergy);
95 
97  const G4ParticleDefinition*,
98  G4double kineticEnergy,
99  G4double cutEnergy,
100  G4double maxEnergy);
101 
103  const G4ParticleDefinition*,
104  G4double kineticEnergy,
105  G4double cutEnergy);
106 
107  virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
108  const G4MaterialCutsCouple*,
109  const G4DynamicParticle*,
110  G4double tmin,
111  G4double maxEnergy);
112 
113  // Compute ion charge
115  const G4Material*,
116  G4double kineticEnergy);
117 
119  const G4Material* mat,
120  G4double kineticEnergy);
121 
122  // add correction to energy loss and ompute non-ionizing energy loss
123  virtual void CorrectionsAlongStep(const G4MaterialCutsCouple*,
124  const G4DynamicParticle*,
125  G4double& eloss,
126  G4double& niel,
127  G4double length);
128 
129 protected:
130 
132  G4double kinEnergy);
133 
134 private:
135 
136  void SetParticle(const G4ParticleDefinition* p);
137 
138  G4double HeEffChargeSquare(G4double z, G4double kinEnergyInMeV) const;
139 
140  // hide assignment operator
143 
144  G4bool HasMaterial(const G4Material* material);
145 
146  G4double StoppingPower(const G4Material* material,
147  G4double kineticEnergy);
148 
150  G4double kineticEnergy) const;
151 
152  G4double DEDX(const G4Material* material, G4double kineticEnergy);
153 
155 
159 
161 
163 
175 
176  G4int iMolecula; // index in the molecula's table
177  G4int iASTAR; // index in ASTAR
180 };
181 
182 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
183 
185 {
186  particle = p;
187  mass = particle->GetPDGMass();
188  spin = particle->GetPDGSpin();
190  chargeSquare = q*q;
191  massRate = mass/CLHEP::proton_mass_c2;
192  ratio = CLHEP::electron_mass_c2/mass;
193 }
194 
195 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
196 
197 #endif
static G4ASTARStopping * fASTAR
G4double lowestKinEnergy
G4double StoppingPower(const G4Material *material, G4double kineticEnergy)
G4bool HasMaterial(const G4Material *material)
G4double z
Definition: TRTMaterials.hh:39
virtual void CorrectionsAlongStep(const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double &eloss, G4double &niel, G4double length)
int G4int
Definition: G4Types.hh:78
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
G4double theZieglerFactor
const G4ParticleDefinition * particle
virtual G4double MinEnergyCut(const G4ParticleDefinition *, const G4MaterialCutsCouple *couple)
G4ParticleChangeForLoss * fParticleChange
virtual G4double GetParticleCharge(const G4ParticleDefinition *p, const G4Material *mat, G4double kineticEnergy)
double A(double temperature)
G4double HeEffChargeSquare(G4double z, G4double kinEnergyInMeV) const
bool G4bool
Definition: G4Types.hh:79
G4double DEDX(const G4Material *material, G4double kineticEnergy)
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
G4double ElectronicStoppingPower(G4double z, G4double kineticEnergy) const
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
G4BraggIonModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="BraggIon")
G4double GetPDGMass() const
G4BraggIonModel & operator=(const G4BraggIonModel &right)
void SetParticle(const G4ParticleDefinition *p)
G4ParticleDefinition * theElectron
G4double GetPDGSpin() const
virtual G4double GetChargeSquareRatio(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy)
double G4double
Definition: G4Types.hh:76
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
static const double eplus
Definition: G4SIunits.hh:196
virtual ~G4BraggIonModel()
G4double GetPDGCharge() const
const G4Material * currentMaterial
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy)
G4EmCorrections * corr
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy)