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;
28 
30 public SensitiveTkDetector,
31 public Observer<const BeginOfEvent*>,
32 public Observer<const BeginOfTrack*>,
33 public Observer<const BeginOfJob*>
34 {
35 public:
38  edm::ParameterSet const &,
39  const SimTrackManager*);
41  bool ProcessHits(G4Step *,G4TouchableHistory *) override;
42  uint32_t setDetUnitId(const G4Step*) override;
43  void EndOfEvent(G4HCofThisEvent*) override;
44 
45  void fillHits(edm::PSimHitContainer&, const std::string&) override;
46  void clearHits() override;
47 
48 private:
49  void createHit(const G4Step *);
50  void sendHit();
51  void updateHit(const G4Step *);
52  bool newHit (const G4Step *);
53  bool closeHit (const G4Step *);
54 
55  void update(const BeginOfEvent *) override;
56  void update(const BeginOfTrack *) override;
57  void update(const BeginOfJob *) override;
58 
60  TrackInformation* getTrackInformation(const G4Track *);
61 
62  // data members initialised before run
64  std::unique_ptr<TrackingSlaveSD> slaveLowTof;
65  std::unique_ptr<TrackingSlaveSD> slaveHighTof;
66  std::unique_ptr<FrameRotation> theRotation;
67  std::unique_ptr<const G4ProcessTypeEnumerator> theG4ProcTypeEnumerator;
68  std::unique_ptr<const G4TrackToParticleID> theG4TrackToParticleID;
69  std::unique_ptr<TrackerG4SimHitNumberingScheme> theNumberingScheme;
71  bool printHits;
73  double rTracker2; // tracker volume R^2
74  double rTracker; // tracker volume R
75  double zTracker; // tracker volume Z
76  float theTofLimit;
77  float energyCut;
79 
80  // run time cache
82  uint32_t lastId;
83  int lastTrack;
84 
85  // cache stuff for debugging and printout
88  const G4VPhysicalVolume* oldVolume;
89  float px,py,pz;
90  int eventno;
92 };
93 
94 #endif
void update(const BeginOfEvent *) override
This routine will be called when the appropriate signal arrives.
std::unique_ptr< TrackingSlaveSD > slaveHighTof
TrackInformation * getTrackInformation(const G4Track *)
std::unique_ptr< FrameRotation > theRotation
std::unique_ptr< const G4ProcessTypeEnumerator > theG4ProcTypeEnumerator
TkAccumulatingSensitiveDetector(const std::string &, const DDCompactView &, const SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
std::unique_ptr< const G4TrackToParticleID > theG4TrackToParticleID
bool ProcessHits(G4Step *, G4TouchableHistory *) override
uint32_t setDetUnitId(const G4Step *) override
type of data representation of DDCompactView
Definition: DDCompactView.h:90
std::unique_ptr< TrackingSlaveSD > slaveLowTof
void fillHits(edm::PSimHitContainer &, const std::string &) override
Local3DPoint RotateAndScale(const Local3DPoint &)
std::unique_ptr< TrackerG4SimHitNumberingScheme > theNumberingScheme
std::vector< PSimHit > PSimHitContainer
void EndOfEvent(G4HCofThisEvent *) override