BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
Loading...
Searching...
No Matches
BDSSDEnergyDepositionGlobal.hh
1/*
2Beam Delivery Simulation (BDSIM) Copyright (C) Royal Holloway,
3University of London 2001 - 2023.
4
5This file is part of BDSIM.
6
7BDSIM is free software: you can redistribute it and/or modify
8it under the terms of the GNU General Public License as published
9by the Free Software Foundation version 3 of the License.
10
11BDSIM is distributed in the hope that it will be useful, but
12WITHOUT ANY WARRANTY; without even the implied warranty of
13MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14GNU General Public License for more details.
15
16You should have received a copy of the GNU General Public License
17along with BDSIM. If not, see <http://www.gnu.org/licenses/>.
18*/
19#ifndef BDSSDENERGYDEPOSITIONGLOBAL_H
20#define BDSSDENERGYDEPOSITIONGLOBAL_H
21
22#include "BDSHitEnergyDepositionGlobal.hh"
23#include "BDSSensitiveDetector.hh"
24
25class G4HCofThisEvent;
26class G4Step;
27class G4TouchableHistory;
28class G4Track;
29
43{
44public:
45 explicit BDSSDEnergyDepositionGlobal(const G4String& name,
46 G4bool killedParticleMassAddedToElossIn = false);
48
49 virtual void Initialize(G4HCofThisEvent* HCE);
50
54 virtual G4bool ProcessHits(G4Step* aStep,
55 G4TouchableHistory* th);
56
61 virtual G4bool ProcessHitsTrack(const G4Track* track,
62 G4TouchableHistory* th);
63
65 virtual G4VHit* last() const;
66
67private:
72
73 G4bool killedParticleMassAddedToEloss;
74 G4String colName;
76 G4int HCIDe;
77
79 G4double energy;
82 G4double stepLength;
83 G4double weight;
84 G4double X,Y,Z; // Global coordinates.
85 G4double globalTime; // Time since start of event.
86 G4int pdgID;
87 G4int trackID;
88 G4int parentID;
91};
92
93#endif
94
Generates BDSHitsEnergyDepositionGlobal from step information.
G4double globalTime
Per hit variable.
G4double stepLength
Per hit variable.
G4double postStepKineticEnergy
Per hit variable.
virtual G4VHit * last() const
Provide access to last hit.
virtual G4bool ProcessHits(G4Step *aStep, G4TouchableHistory *th)
virtual G4bool ProcessHitsTrack(const G4Track *track, G4TouchableHistory *th)
G4double preStepKineticEnergy
Per hit variable.
BDSSDEnergyDepositionGlobal & operator=(const BDSSDEnergyDepositionGlobal &)
assignment and copy constructor not implemented nor used
Virtual class to define interface for ordered multi-sensitive detector.