Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DNADingfelderChargeIncreaseModel.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#ifndef G4DNADingfelderChargeIncreaseModel_h
29#define G4DNADingfelderChargeIncreaseModel_h 1
30
31#include "G4VEmModel.hh"
34#include "G4Electron.hh"
35#include "G4Proton.hh"
37#include "G4NistManager.hh"
38
40{
41public:
42
44 const G4String& nam = "DNADingfelderChargeIncreaseModel");
46
49
50 void Initialise(const G4ParticleDefinition*, const G4DataVector&) override;
52 const G4ParticleDefinition* p,
53 G4double ekin,
54 G4double emin,
55 G4double emax) override;
56 void SampleSecondaries(std::vector<G4DynamicParticle*>*,
58 const G4DynamicParticle*,
59 G4double tmin,
60 G4double maxEnergy) override;
61 inline void SelectStationary(G4bool input);
62protected:
64private:
65 // Water density table
66 const std::vector<G4double>* fpMolWaterDensity = nullptr;
67
68 std::map<G4String,G4double,std::less<G4String> > lowEnergyLimit;
69 std::map<G4String,G4double,std::less<G4String> > highEnergyLimit;
70
71 G4bool isInitialised = false, statCode = false;
72 // Verbosity scale:
73 // 0 = nothing
74 // 1 = warning for energy non-conservation
75 // 2 = details of energy budget
76 // 3 = calculation of cross sections, file openings, sampling of atoms
77 // 4 = entering in methods
78 G4int verboseLevel = 0;
79
80 // Partial cross section
81
82 G4double PartialCrossSection(const G4double& energy, const G4int& level, const G4ParticleDefinition* particle);
83
84 G4double Sum(const G4double& energy, const G4ParticleDefinition* particle);
85
86 G4int RandomSelect(const G4double& energy, const G4ParticleDefinition* particle);
87
88 G4int numberOfPartialCrossSections[2] = {0}; // 2 is the particle type index
89 G4double f0[2][2] = {{0, 0},{0, 0}};
90 G4double a0[2][2] = {{0, 0},{0, 0}};
91 G4double a1[2][2] = {{0, 0},{0, 0}};
92 G4double b0[2][2] = {{0, 0},{0, 0}};
93 G4double b1[2][2] = {{0, 0},{0, 0}};
94 G4double c0[2][2] = {{0, 0},{0, 0}};
95 G4double d0[2][2] = {{0, 0},{0, 0}};
96 G4double x0[2][2] = {{0, 0},{0, 0}};
97 G4double x1[2][2] = {{0, 0},{0, 0}};
98
99 // Final state
100 G4int NumberOfFinalStates(G4ParticleDefinition* particleDefinition, G4int finalStateIndex);
101 G4ParticleDefinition* OutgoingParticleDefinition(G4ParticleDefinition* particleDefinition, G4int finalStateIndex);
102 G4double IncomingParticleBindingEnergyConstant(G4ParticleDefinition* particleDefinition, G4int finalStateIndex);
103 // Reusable particle definitions
104 G4ParticleDefinition* hydrogenDef = nullptr;
105 G4ParticleDefinition* alphaPlusPlusDef = nullptr;
106 G4ParticleDefinition* alphaPlusDef = nullptr;
107 G4ParticleDefinition* heliumDef = nullptr;
108
109};
110//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
111
113{
114 statCode = input;
115}
116
117//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
118
119#endif
120
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax) override
G4DNADingfelderChargeIncreaseModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="DNADingfelderChargeIncreaseModel")
G4DNADingfelderChargeIncreaseModel & operator=(const G4DNADingfelderChargeIncreaseModel &right)=delete
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
G4DNADingfelderChargeIncreaseModel(const G4DNADingfelderChargeIncreaseModel &)=delete
~G4DNADingfelderChargeIncreaseModel() override=default