Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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 "G4BraggModel.hh"
60
61class G4ASTARStopping;
62
64{
65
66public:
67
68 explicit G4BraggIonModel(const G4ParticleDefinition* p = nullptr,
69 const G4String& nam = "BraggIon");
70
71 ~G4BraggIonModel() override;
72
74 const G4DataVector&) override;
75
78 G4double kineticEnergy,
80 G4double cutEnergy,
81 G4double maxEnergy) override;
82
85 G4double kineticEnergy,
86 G4double cutEnergy,
87 G4double maxEnergy) override;
88
91 G4double kineticEnergy,
92 G4double cutEnergy) override;
93
94 // Compute ion charge not applied to alpha
96 const G4Material*,
97 G4double kineticEnergy) override;
98
99 // add correction to energy loss and ompute non-ionizing energy loss
101 const G4DynamicParticle*,
102 const G4double& length,
103 G4double& eloss) override;
104
105 // hide assignment operator
106 G4BraggIonModel & operator=(const G4BraggIonModel &right) = delete;
108
109private:
110
111 G4double HeEffChargeSquare(const G4double z,
112 const G4double kinEnergyInMeV) const;
113
114 G4int HasMaterialForHe(const G4Material* material) const;
115
116 G4double HeStoppingPower(const G4double kinEnergy) const;
117
118 G4double HeElectronicStoppingPower(const G4int z, const G4double kinEnergy) const;
119
120 G4double HeDEDX(const G4Material* material, const G4double kinEnergy);
121
122 static G4ASTARStopping* fASTAR;
123
124 G4double heChargeSquare = 4.0;
125 G4double HeMass;
126 G4double massFactor;
127
128 G4int iASTAR = -1; // index in ASTAR
129 G4bool isAlpha = false;
130 G4bool isFirstAlpha = false;
131};
132
133//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
134
135#endif
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
const G4double A[17]
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
G4BraggIonModel(const G4BraggIonModel &)=delete
G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
G4double GetChargeSquareRatio(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy) override
~G4BraggIonModel() override
G4BraggIonModel & operator=(const G4BraggIonModel &right)=delete
G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy) override
G4BraggIonModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="BraggIon")
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