Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4BraggModel.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: G4BraggModel.hh 96934 2016-05-18 09:10:41Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class header file
31 //
32 //
33 // File name: G4BraggModel
34 //
35 // Author: Vladimir Ivanchenko
36 //
37 // Creation date: 03.01.2002
38 //
39 // Modifications:
40 // 23-12-02 V.Ivanchenko change interface in order to moveto cut per region
41 // 24-01-03 Make models region aware (V.Ivanchenko)
42 // 13-02-03 Add name (V.Ivanchenko)
43 // 12-11-03 Fix for GenericIons (V.Ivanchenko)
44 // 11-04-05 Major optimisation of internal interfaces (V.Ivantchenko)
45 // 15-02-06 ComputeCrossSectionPerElectron, ComputeCrossSectionPerAtom (mma)
46 // 25-04-06 Added stopping data from PSTAR (V.Ivanchenko)
47 // 12-08-08 Added methods GetParticleCharge, GetChargeSquareRatio,
48 // CorrectionsAlongStep needed for ions(V.Ivanchenko)
49 
50 //
51 // Class Description:
52 //
53 // Implementation of energy loss and delta-electron production
54 // by heavy slow charged particles using ICRU'49 and NIST evaluated data
55 // for protons
56 
57 // -------------------------------------------------------------------
58 //
59 
60 #ifndef G4BraggModel_h
61 #define G4BraggModel_h 1
62 
64 
65 #include "G4VEmModel.hh"
66 #include "G4PSTARStopping.hh"
67 
69 class G4EmCorrections;
70 
71 class G4BraggModel : public G4VEmModel
72 {
73 
74 public:
75 
76  explicit G4BraggModel(const G4ParticleDefinition* p = nullptr,
77  const G4String& nam = "Bragg");
78 
79  virtual ~G4BraggModel();
80 
81  virtual void Initialise(const G4ParticleDefinition*,
82  const G4DataVector&) override;
83 
85  const G4ParticleDefinition*,
86  G4double kineticEnergy,
87  G4double cutEnergy,
88  G4double maxEnergy);
89 
91  const G4ParticleDefinition*,
92  G4double kineticEnergy,
94  G4double cutEnergy,
95  G4double maxEnergy) override;
96 
98  const G4ParticleDefinition*,
99  G4double kineticEnergy,
100  G4double cutEnergy,
101  G4double maxEnergy) override;
102 
104  const G4ParticleDefinition*,
105  G4double kineticEnergy,
106  G4double cutEnergy) override;
107 
108  virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
109  const G4MaterialCutsCouple*,
110  const G4DynamicParticle*,
111  G4double tmin,
112  G4double maxEnergy) override;
113 
114  // Compute ion charge
116  const G4Material*,
117  G4double kineticEnergy) override;
118 
120  const G4Material* mat,
121  G4double kineticEnergy) override;
122 
123 protected:
124 
126  G4double kinEnergy) final;
127 
128  inline G4double GetChargeSquareRatio() const;
129 
130  inline void SetChargeSquareRatio(G4double val);
131 
132 private:
133 
134  inline void SetParticle(const G4ParticleDefinition* p);
135 
136  G4bool HasMaterial(const G4Material* material);
137 
138  G4double StoppingPower(const G4Material* material,
139  G4double kineticEnergy);
140 
141  G4double ElectronicStoppingPower(G4double z,
142  G4double kineticEnergy) const;
143 
144  G4double DEDX(const G4Material* material, G4double kineticEnergy);
145 
146  G4bool MolecIsInZiegler1988(const G4Material* material);
147 
148  G4double ChemicalFactor(G4double kineticEnergy, G4double eloss125) const;
149 
150  // hide assignment operator
151  G4BraggModel & operator=(const G4BraggModel &right) = delete;
152  G4BraggModel(const G4BraggModel&) = delete;
153 
154  G4EmCorrections* corr;
155 
156  const G4ParticleDefinition* particle;
157  G4ParticleDefinition* theElectron;
158  G4ParticleChangeForLoss* fParticleChange;
159  static G4PSTARStopping* fPSTAR;
160 
161  const G4Material* currentMaterial;
162 
163  G4double mass;
164  G4double spin;
165  G4double chargeSquare;
166  G4double massRate;
167  G4double ratio;
168  G4double lowestKinEnergy;
169  G4double protonMassAMU;
170  G4double theZieglerFactor;
171  G4double expStopPower125; // Experimental Stopping power at 125keV
172 
173  G4int iMolecula; // index in the molecula's table
174  G4int iPSTAR; // index in PSTAR
175  G4bool isIon;
176 };
177 
178 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
179 
180 inline void G4BraggModel::SetParticle(const G4ParticleDefinition* p)
181 {
182  particle = p;
183  mass = particle->GetPDGMass();
184  spin = particle->GetPDGSpin();
185  G4double q = particle->GetPDGCharge()/CLHEP::eplus;
186  chargeSquare = q*q;
187  massRate = mass/CLHEP::proton_mass_c2;
188  ratio = CLHEP::electron_mass_c2/mass;
189 }
190 
192 {
193  return chargeSquare;
194 }
195 
197 {
198  chargeSquare = val;
199 }
200 
201 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
202 
203 #endif
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
static constexpr double proton_mass_c2
const char * p
Definition: xmltok.h:285
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy) final
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy) override
G4double GetChargeSquareRatio() const
int G4int
Definition: G4Types.hh:78
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy) override
static constexpr double electron_mass_c2
double A(double temperature)
bool G4bool
Definition: G4Types.hh:79
static constexpr double eplus
virtual ~G4BraggModel()
G4double GetPDGMass() const
virtual G4double GetParticleCharge(const G4ParticleDefinition *p, const G4Material *mat, G4double kineticEnergy) override
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
G4double GetPDGSpin() const
G4BraggModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="Bragg")
Definition: G4BraggModel.cc:84
void SetChargeSquareRatio(G4double val)
double G4double
Definition: G4Types.hh:76
G4double GetPDGCharge() const