CMS 3D CMS Logo

MuonSensitiveDetector.h
Go to the documentation of this file.
1 #ifndef SimG4CMS_Muon_MuonSensitiveDetector_h
2 #define SimG4CMS_Muon_MuonSensitiveDetector_h
3 
22 
23 #include <string>
24 
25 class MuonSlaveSD;
27 class MuonFrameRotation;
28 class UpdatablePSimHit;
29 class MuonSubDetector;
30 class MuonG4Numbering;
31 class SimHitPrinter;
32 class G4Step;
34 class SimTrackManager;
35 
36 class MuonSensitiveDetector : public SensitiveTkDetector, public Observer<const BeginOfEvent*> {
37 public:
38  explicit MuonSensitiveDetector(const std::string&,
39  const DDCompactView&,
41  edm::ParameterSet const&,
42  const SimTrackManager*);
43  ~MuonSensitiveDetector() override;
44  G4bool ProcessHits(G4Step*, G4TouchableHistory*) override;
45  uint32_t setDetUnitId(const G4Step*) override;
46  void EndOfEvent(G4HCofThisEvent*) override;
47 
48  void fillHits(edm::PSimHitContainer&, const std::string&) override;
49  void clearHits() override;
50 
51  const MuonSlaveSD* GetSlaveMuon() const { return slaveMuon; }
52 
53 protected:
54  void update(const BeginOfEvent*) override;
55 
56 private:
57  inline Local3DPoint cmsUnits(const Local3DPoint& v) { return Local3DPoint(v.x() * 0.1, v.y() * 0.1, v.z() * 0.1); }
58 
64 
65  bool newHit(const G4Step*);
66  void createHit(const G4Step*);
67  void updateHit(const G4Step*);
68  void saveHit();
69 
75  Local3DPoint InitialStepPositionVsParent(const G4Step* currentStep, G4int levelsUp);
76  Local3DPoint FinalStepPositionVsParent(const G4Step* currentStep, G4int levelsUp);
77 
78  const G4VPhysicalVolume* thePV;
80  uint32_t theDetUnitId;
81  uint32_t newDetUnitId;
83 
84  bool printHits;
86 
87  //--- SimTracks cuts
90 
92 
94 };
95 
96 #endif // MuonSensitiveDetector_h
const G4VPhysicalVolume * thePV
Local3DPoint FinalStepPositionVsParent(const G4Step *currentStep, G4int levelsUp)
G4bool ProcessHits(G4Step *, G4TouchableHistory *) override
void updateHit(const G4Step *)
T y() const
Definition: PV3DBase.h:63
Compact representation of the geometrical detector hierarchy.
Definition: DDCompactView.h:80
MuonSubDetector * detector
void fillHits(edm::PSimHitContainer &, const std::string &) override
Local3DPoint cmsUnits(const Local3DPoint &v)
void createHit(const G4Step *)
G4ProcessTypeEnumerator * theG4ProcessTypeEnumerator
T z() const
Definition: PV3DBase.h:64
bool newHit(const G4Step *)
const SimTrackManager * theManager
MuonFrameRotation * theRotation
MuonSensitiveDetector(const std::string &, const DDCompactView &, const SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
const MuonSlaveSD * GetSlaveMuon() const
uint32_t setDetUnitId(const G4Step *) override
Point3DBase< float, LocalTag > Local3DPoint
Definition: LocalPoint.h:9
void EndOfEvent(G4HCofThisEvent *) override
MuonSimHitNumberingScheme * numbering
void update(const BeginOfEvent *) override
This routine will be called when the appropriate signal arrives.
std::vector< PSimHit > PSimHitContainer
T x() const
Definition: PV3DBase.h:62
Local3DPoint InitialStepPositionVsParent(const G4Step *currentStep, G4int levelsUp)
UpdatablePSimHit * theHit
MuonG4Numbering * g4numbering