Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VRangeToEnergyConverter.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// G4VRangeToEnergyConverter
27//
28// Class description:
29//
30// Base class for Range to Energy Converters.
31// Cut in energy corresponding to given cut value in range
32// is calculated for a material by using Convert() method.
33
34// Author: H.Kurashige, 05 October 2002 - First implementation
35// --------------------------------------------------------------------
36#ifndef G4VRangeToEnergyConverter_hh
37#define G4VRangeToEnergyConverter_hh 1
38
39#include <vector>
40
41#include "globals.hh"
43#include "G4Material.hh"
44
46{
47public:
48
50
52
53 // operators are not used
56 (const G4VRangeToEnergyConverter &r) = delete;
59
60 // Calculate energy cut from given range cut for the material
61 virtual G4double Convert(const G4double rangeCut, const G4Material* material);
62
63 // Set energy range for all particle type
64 // if highedge > 10 GeV, highedge value is not changed
65 static void SetEnergyRange(const G4double lowedge, const G4double highedge);
66
67 // Get energy range for all particle type
70
71 // Get/set max cut energy for all particle type
72 // No check on the value
74 static void SetMaxEnergyCut(const G4double value);
75
76 // Return pointer to the particle type which this converter takes care of
77 inline const G4ParticleDefinition* GetParticleType() const;
78
79 inline void SetVerboseLevel(G4int value);
80 inline G4int GetVerboseLevel() const;
81 // control flag for output message
82 // 0: Silent
83 // 1: Warning message
84 // 2: More
85
86protected:
87
88 virtual G4double ComputeValue(const G4int Z, const G4double kinEnergy) = 0;
89
90private:
91
92 static void FillEnergyVector(const G4double emin, const G4double emax);
93
94 G4double ConvertForGamma(const G4double rangeCut, const G4Material* material);
95
96 G4double ConvertForElectron(const G4double rangeCut,
97 const G4Material* material);
98
99 inline G4double LiniearInterpolation(const G4double e1, const G4double e2,
100 const G4double r1, const G4double r2,
101 const G4double r);
102
103protected:
104
107
108private:
109
110 static G4double sEmin;
111 static G4double sEmax;
112 static std::vector<G4double>* sEnergy;
113 static G4int sNbinPerDecade;
114 static G4int sNbin;
115
116 G4int verboseLevel = 1;
117 G4bool isFirstInstance = false;
118};
119
120// ------------------
121// Inline methods
122// ------------------
123
124inline
126{
127 verboseLevel = value;
128}
129
130inline
132{
133 return verboseLevel;
134}
135
136inline
141
142inline G4double G4VRangeToEnergyConverter::LiniearInterpolation(
143 const G4double e1, const G4double e2,
144 const G4double r1, const G4double r2, const G4double r)
145{
146 return (r1 == r2) ? e1 : e1 + (e2 - e1)*(r - r1)/(r2 - r1);
147}
148
149#endif
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
static void SetMaxEnergyCut(const G4double value)
virtual G4double ComputeValue(const G4int Z, const G4double kinEnergy)=0
const G4ParticleDefinition * GetParticleType() const
virtual G4double Convert(const G4double rangeCut, const G4Material *material)
G4bool operator==(const G4VRangeToEnergyConverter &r) const =delete
static void SetEnergyRange(const G4double lowedge, const G4double highedge)
G4VRangeToEnergyConverter(const G4VRangeToEnergyConverter &r)=delete
const G4ParticleDefinition * theParticle
G4bool operator!=(const G4VRangeToEnergyConverter &r) const =delete