Loading [MathJax]/extensions/tex2jax.js
 |
ATLAS Offline Software
|
Go to the documentation of this file.
5 #ifndef CaloG4_SimulationEnergies_H
6 #define CaloG4_SimulationEnergies_H
8 #include "G4TrackStatus.hh"
9 #include "G4ParticleDefinition.hh"
10 #include "G4DynamicParticle.hh"
11 #include "G4ThreeVector.hh"
20 class G4ParticleDefinition;
56 void Energies(
const G4Step* , std::vector<G4double> & )
const;
82 std::map<eEnergyCategory, G4double>
energy;
110 G4double totalEnergy,
111 G4double kineticEnergy)
const;
115 G4double totalEnergy,
116 G4double kineticEnergy)
const;
123 std::unique_ptr<G4Step>
CreateFakeStep(G4Track* a_track, G4double a_energy)
const;
129 #endif // CaloG4_SimulationEnergies_H
void Energies(const G4Step *, std::vector< G4double > &) const
The simple method to call from a calibration calculator: Examine the G4Step and return the energies r...
@ kNumberOfEnergyCategories
ClassifyResult_t Classify(const G4Step *, const G4bool processEscaped=true) const
This method performs the actual function of this class: It classifies the components of the energy de...
G4double measurableEnergyV2(const G4ParticleDefinition *particleDef, G4int PDGEncoding, G4double totalEnergy, G4double kineticEnergy) const
std::unique_ptr< G4Step > CreateFakeStep(G4Track *a_track, G4double a_energy) const
G4bool ParticleIsNeutrino(G4ParticleDefinition *particle) const
G4double measurableEnergy(const G4ParticleDefinition *particleDef, G4int PDGEncoding, G4double totalEnergy, G4double kineticEnergy) const
Some private methods for internal calculations:
This defines the results returned by the energy classification; these detailed results are mostly use...
std::map< eEnergyCategory, G4double > energy
virtual ~SimulationEnergies()
G4bool ProcessEscapedEnergy(G4Step *fakeStep) const
eEnergyCategory
Accessing detailed information: