CMS 3D CMS Logo

TkAccumulatingSensitiveDetector.h
Go to the documentation of this file.
1 #ifndef SimG4CMS_TkAccumulatingSensitiveDetector_H
2 #define SimG4CMS_TkAccumulatingSensitiveDetector_H
3 
14 
15 #include "G4Step.hh"
16 #include "G4Track.hh"
17 
18 #include <string>
19 
20 class TrackInformation;
21 class SimTrackManager;
22 class TrackingSlaveSD;
23 class FrameRotation;
24 class UpdatablePSimHit;
27 
29 public SensitiveTkDetector,
30 public Observer<const BeginOfEvent*>,
31 public Observer<const BeginOfTrack*>,
32 public Observer<const BeginOfJob*>
33 {
34 public:
37  edm::ParameterSet const &,
38  const SimTrackManager*);
40  bool ProcessHits(G4Step *,G4TouchableHistory *) override;
41  uint32_t setDetUnitId(const G4Step*) override;
42  void EndOfEvent(G4HCofThisEvent*) override;
43 
44  void fillHits(edm::PSimHitContainer&, const std::string&) override;
45  void clearHits() override;
46 
47 private:
48  void createHit(const G4Step *);
49  void sendHit();
50  void updateHit(const G4Step *);
51  bool newHit (const G4Step *);
52  bool closeHit (const G4Step *);
53 
54 protected:
55  void update(const BeginOfEvent *) override;
56  void update(const BeginOfTrack *) override;
57  void update(const BeginOfJob *) override;
58 
59 private:
60 
61  // data members initialised before run
63  std::unique_ptr<TrackingSlaveSD> slaveLowTof;
64  std::unique_ptr<TrackingSlaveSD> slaveHighTof;
65  std::unique_ptr<FrameRotation> theRotation;
66  std::unique_ptr<const G4ProcessTypeEnumerator> theG4ProcTypeEnumerator;
69  bool printHits;
71  double rTracker2; // tracker volume R^2
72  double rTracker; // tracker volume R
73  double zTracker; // tracker volume Z
74  float theTofLimit;
75  float energyCut;
77 
78  // run time cache
80  uint32_t lastId;
81  int lastTrack;
82 
83  // cache stuff for debugging and printout
86  const G4VPhysicalVolume* oldVolume;
87  float px,py,pz;
88  int eventno;
90 };
91 
92 #endif
void update(const BeginOfEvent *) override
This routine will be called when the appropriate signal arrives.
std::unique_ptr< TrackingSlaveSD > slaveHighTof
std::unique_ptr< FrameRotation > theRotation
std::unique_ptr< const G4ProcessTypeEnumerator > theG4ProcTypeEnumerator
TkAccumulatingSensitiveDetector(const std::string &, const DDCompactView &, const SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
bool ProcessHits(G4Step *, G4TouchableHistory *) override
uint32_t setDetUnitId(const G4Step *) override
Compact representation of the geometrical detector hierarchy.
Definition: DDCompactView.h:80
std::unique_ptr< TrackingSlaveSD > slaveLowTof
void fillHits(edm::PSimHitContainer &, const std::string &) override
std::vector< PSimHit > PSimHitContainer
void EndOfEvent(G4HCofThisEvent *) override
TrackerG4SimHitNumberingScheme * theNumberingScheme