Geant4_10
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 74790 2013-10-22 07:31:37Z 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 
77  const G4String& nam = "Bragg");
78 
79  virtual ~G4BraggModel();
80 
81  virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
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 protected:
123 
125  G4double kinEnergy);
126 
127  inline G4double GetChargeSquareRatio() const;
128 
129  inline void SetChargeSquareRatio(G4double val);
130 
131 private:
132 
133  inline void SetParticle(const G4ParticleDefinition* p);
134 
135  G4bool HasMaterial(const G4Material* material);
136 
137  G4double StoppingPower(const G4Material* material,
138  G4double kineticEnergy);
139 
140  G4double ElectronicStoppingPower(G4double z,
141  G4double kineticEnergy) const;
142 
143  G4double DEDX(const G4Material* material, G4double kineticEnergy);
144 
145  G4bool MolecIsInZiegler1988(const G4Material* material);
146 
147  G4double ChemicalFactor(G4double kineticEnergy, G4double eloss125) const;
148 
149  void SetExpStopPower125(G4double value) {expStopPower125 = value;};
150 
151  // hide assignment operator
152  G4BraggModel & operator=(const G4BraggModel &right);
153  G4BraggModel(const G4BraggModel&);
154 
155  G4EmCorrections* corr;
156 
157  const G4ParticleDefinition* particle;
158  G4ParticleDefinition* theElectron;
159  G4ParticleChangeForLoss* fParticleChange;
160  G4PSTARStopping pstar;
161 
162  const G4Material* currentMaterial;
163 
164  G4double mass;
165  G4double spin;
166  G4double chargeSquare;
167  G4double massRate;
168  G4double ratio;
169  G4double lowestKinEnergy;
170  G4double protonMassAMU;
171  G4double theZieglerFactor;
172  G4double expStopPower125; // Experimental Stopping power at 125keV
173 
174  G4int iMolecula; // index in the molecula's table
175  G4int iPSTAR; // index in PSTAR
176  G4bool isIon;
177  G4bool isInitialised;
178 };
179 
180 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
181 
182 inline void G4BraggModel::SetParticle(const G4ParticleDefinition* p)
183 {
184  particle = p;
185  mass = particle->GetPDGMass();
186  spin = particle->GetPDGSpin();
187  G4double q = particle->GetPDGCharge()/CLHEP::eplus;
188  chargeSquare = q*q;
189  massRate = mass/CLHEP::proton_mass_c2;
190  ratio = CLHEP::electron_mass_c2/mass;
191 }
192 
194 {
195  return chargeSquare;
196 }
197 
199 {
200  chargeSquare = val;
201 }
202 
203 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
204 
205 #endif
G4BraggModel(const G4ParticleDefinition *p=0, const G4String &nam="Bragg")
Definition: G4BraggModel.cc:84
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy)
const char * p
Definition: xmltok.h:285
G4double GetChargeSquareRatio() const
int G4int
Definition: G4Types.hh:78
Float_t mat
Definition: plot.C:40
string material
Definition: eplot.py:19
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy)
Float_t Z
Definition: plot.C:39
bool G4bool
Definition: G4Types.hh:79
virtual G4double GetParticleCharge(const G4ParticleDefinition *p, const G4Material *mat, G4double kineticEnergy)
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy)
virtual ~G4BraggModel()
G4double GetPDGMass() const
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
tuple z
Definition: test.py:28
G4double GetPDGSpin() const
const XML_Char int const XML_Char * value
Definition: expat.h:331
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
void SetChargeSquareRatio(G4double val)
double G4double
Definition: G4Types.hh:76
G4double GetPDGCharge() const
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)