Geant4-11
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//
27// -------------------------------------------------------------------
28//
29// GEANT4 Class header file
30//
31//
32// File name: G4BraggIonModel
33//
34// Author: Vladimir Ivanchenko
35//
36// Creation date: 13.10.2004
37//
38// Modifications:
39// 09-11-04 Migration to new interface of Store/Retrieve tables (V.Ivantchenko)
40// 11-05-05 Major optimisation of internal interfaces (V.Ivantchenko)
41// 15-02-06 ComputeCrossSectionPerElectron, ComputeCrossSectionPerAtom (mma)
42// 25-04-06 Add stopping data from ASTAR (V.Ivanchenko)
43// 12-08-08 Added methods GetParticleCharge, GetChargeSquareRatio,
44// CorrectionsAlongStep needed for ions(V.Ivanchenko)
45
46//
47// Class Description:
48//
49// Implementation of energy loss and delta-electron production
50// by heavy slow charged particles using ICRU'49, NIST, and ICRU90
51// evaluated data for alpha and protons
52
53// -------------------------------------------------------------------
54//
55
56#ifndef G4BraggIonModel_h
57#define G4BraggIonModel_h 1
58
59#include "G4VEmModel.hh"
60#include "G4ASTARStopping.hh"
61
63class G4EmCorrections;
65
67{
68
69public:
70
71 explicit G4BraggIonModel(const G4ParticleDefinition* p = nullptr,
72 const G4String& nam = "BraggIon");
73
74 ~G4BraggIonModel() override;
75
77 const G4DataVector&) override;
78
80 const G4MaterialCutsCouple* couple) override;
81
84 G4double kineticEnergy,
85 G4double cutEnergy,
86 G4double maxEnergy);
87
90 G4double kineticEnergy,
92 G4double cutEnergy,
93 G4double maxEnergy) override;
94
97 G4double kineticEnergy,
98 G4double cutEnergy,
99 G4double maxEnergy) override;
100
103 G4double kineticEnergy,
104 G4double cutEnergy) override;
105
106 void SampleSecondaries(std::vector<G4DynamicParticle*>*,
108 const G4DynamicParticle*,
109 G4double tmin,
110 G4double maxEnergy) override;
111
112 // Compute ion charge not applied to alpha
114 const G4Material*,
115 G4double kineticEnergy) override;
116
118 const G4Material* mat,
119 G4double kineticEnergy) override;
120
121 // add correction to energy loss and ompute non-ionizing energy loss
123 const G4DynamicParticle*,
124 const G4double& length,
125 G4double& eloss) override;
126
127 // hide assignment operator
128 G4BraggIonModel & operator=(const G4BraggIonModel &right) = delete;
130
131protected:
132
134 G4double kinEnergy) final;
135
136private:
137
138 void SetParticle(const G4ParticleDefinition* p);
139
141 const G4double kinEnergyInMeV) const;
142
143 G4int HasMaterial(const G4Material* material) const;
144
146 const G4double kinEnergy) const;
147
149 const G4double kinEnergy) const;
150
151 G4double DEDX(const G4Material* material, const G4double kinEnergy);
152
158
159 const G4Material* currentMaterial = nullptr;
160 const G4Material* baseMaterial = nullptr;
161
163
175
176 G4int iMolecula = -1; // index in the molecula's table
177 G4int iASTAR = -1; // index in ASTAR
179 G4bool isIon = false;
181};
182
183//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
184
185#endif
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4double A[17]
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
G4BraggIonModel(const G4BraggIonModel &)=delete
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
G4double MinEnergyCut(const G4ParticleDefinition *, const G4MaterialCutsCouple *couple) override
void SetParticle(const G4ParticleDefinition *p)
G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
static G4ASTARStopping * fASTAR
G4double GetChargeSquareRatio(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy) override
~G4BraggIonModel() override
G4double DEDX(const G4Material *material, const G4double kinEnergy)
G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy) final
const G4ParticleDefinition * theElectron
G4double effChargeSquare
G4BraggIonModel & operator=(const G4BraggIonModel &right)=delete
G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy) override
G4double StoppingPower(const G4Material *material, const G4double kinEnergy) const
G4double HeEffChargeSquare(const G4double z, const G4double kinEnergyInMeV) const
const G4Material * currentMaterial
G4double GetParticleCharge(const G4ParticleDefinition *p, const G4Material *mat, G4double kineticEnergy) override
G4BraggIonModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="BraggIon")
const G4ParticleDefinition * particle
G4int HasMaterial(const G4Material *material) const
G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy) override
void CorrectionsAlongStep(const G4MaterialCutsCouple *, const G4DynamicParticle *, const G4double &length, G4double &eloss) override
G4EmCorrections * corr
G4double theZieglerFactor
G4double lowestKinEnergy
G4double ElectronicStoppingPower(const G4double z, const G4double kinEnergy) const
G4ICRU90StoppingData * fICRU90
G4ParticleChangeForLoss * fParticleChange
const G4Material * baseMaterial
string material
Definition: eplot.py:19