CMS 3D CMS Logo

SensitiveDetector.h
Go to the documentation of this file.
1 #ifndef SimG4Core_SensitiveDetector_H
2 #define SimG4Core_SensitiveDetector_H
3 
7 
8 #include "G4VSensitiveDetector.hh"
9 
10 #include <string>
11 #include <cstdint>
12 
13 class G4Track;
14 class G4Step;
15 class G4HCofThisEvent;
16 class G4TouchableHistory;
17 class G4VPhysicalVolume;
18 
19 class SensitiveDetector : public G4VSensitiveDetector {
20 public:
21  explicit SensitiveDetector(const std::string& iname, const SensitiveDetectorCatalog&, bool calo);
22 
23  ~SensitiveDetector() override;
24 
25  void Initialize(G4HCofThisEvent* eventHC) override;
26  G4bool ProcessHits(G4Step* step, G4TouchableHistory* tHistory) override = 0;
27  void EndOfEvent(G4HCofThisEvent* eventHC) override;
28 
29  virtual uint32_t setDetUnitId(const G4Step* step) = 0;
30  virtual void clearHits() = 0;
31 
32  inline const std::vector<std::string>& getNames() const { return m_namesOfSD; }
33 
34  inline bool isCaloSD() const { return m_isCalo; }
35 
36 protected:
37  // generic geometry methods, all coordinates in mm
39  Local3DPoint InitialStepPosition(const G4Step* step, coordinates) const;
40  Local3DPoint FinalStepPosition(const G4Step* step, coordinates) const;
41 
42  // local coordinates to the volume, in which the step is done
43  Local3DPoint LocalPreStepPosition(const G4Step* step) const;
44  Local3DPoint LocalPostStepPosition(const G4Step* step) const;
45 
46  inline Local3DPoint ConvertToLocal3DPoint(const G4ThreeVector& point) const {
47  return Local3DPoint(point.x(), point.y(), point.z());
48  }
49 
50  void setNames(const std::vector<std::string>&);
51  void NaNTrap(const G4Step* step) const;
52 
53  TrackInformation* cmsTrackInformation(const G4Track* aTrack);
54 
55 private:
56  void AssignSD(const std::string& vname);
57 
58  std::vector<std::string> m_namesOfSD;
59  bool m_isCalo;
60 };
61 
62 #endif
SensitiveDetector::SensitiveDetector
SensitiveDetector(const std::string &iname, const SensitiveDetectorCatalog &, bool calo)
Definition: SensitiveDetector.cc:18
SensitiveDetector::clearHits
virtual void clearHits()=0
SensitiveDetector::~SensitiveDetector
~SensitiveDetector() override
Definition: SensitiveDetector.cc:39
step
step
Definition: StallMonitor.cc:94
SensitiveDetector::NaNTrap
void NaNTrap(const G4Step *step) const
Definition: SensitiveDetector.cc:106
SensitiveDetector::ProcessHits
G4bool ProcessHits(G4Step *step, G4TouchableHistory *tHistory) override=0
SensitiveDetector::AssignSD
void AssignSD(const std::string &vname)
Definition: SensitiveDetector.cc:45
SensitiveDetector::coordinates
coordinates
Definition: SensitiveDetector.h:38
SensitiveDetector
Definition: SensitiveDetector.h:19
SensitiveDetector::m_namesOfSD
std::vector< std::string > m_namesOfSD
Definition: SensitiveDetector.h:58
SensitiveDetector::setDetUnitId
virtual uint32_t setDetUnitId(const G4Step *step)=0
SensitiveDetectorCatalog
Definition: SensitiveDetectorCatalog.h:10
Point3DBase< float, LocalTag >
SensitiveDetector::FinalStepPosition
Local3DPoint FinalStepPosition(const G4Step *step, coordinates) const
Definition: SensitiveDetector.cc:66
SensitiveDetector::LocalPreStepPosition
Local3DPoint LocalPreStepPosition(const G4Step *step) const
Definition: SensitiveDetector.cc:78
SensitiveDetector::InitialStepPosition
Local3DPoint InitialStepPosition(const G4Step *step, coordinates) const
Definition: SensitiveDetector.cc:54
SensitiveDetectorCatalog.h
SensitiveDetector::ConvertToLocal3DPoint
Local3DPoint ConvertToLocal3DPoint(const G4ThreeVector &point) const
Definition: SensitiveDetector.h:46
TrackInformation
Definition: TrackInformation.h:12
SensitiveDetector::LocalCoordinates
Definition: SensitiveDetector.h:38
SensitiveDetector::isCaloSD
bool isCaloSD() const
Definition: SensitiveDetector.h:34
TrackInformation.h
SensitiveDetector::m_isCalo
bool m_isCalo
Definition: SensitiveDetector.h:59
AlCaHLTBitMon_QueryRunRegistry.string
string string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
SensitiveDetector::getNames
const std::vector< std::string > & getNames() const
Definition: SensitiveDetector.h:32
SensitiveDetector::setNames
void setNames(const std::vector< std::string > &)
Definition: SensitiveDetector.cc:101
SensitiveDetector::cmsTrackInformation
TrackInformation * cmsTrackInformation(const G4Track *aTrack)
Definition: SensitiveDetector.cc:92
SensitiveDetector::Initialize
void Initialize(G4HCofThisEvent *eventHC) override
Definition: SensitiveDetector.cc:41
LocalPoint.h
SensitiveDetector::LocalPostStepPosition
Local3DPoint LocalPostStepPosition(const G4Step *step) const
Definition: SensitiveDetector.cc:85
SensitiveDetector::WorldCoordinates
Definition: SensitiveDetector.h:38
SensitiveDetector::EndOfEvent
void EndOfEvent(G4HCofThisEvent *eventHC) override
Definition: SensitiveDetector.cc:43
calo
Definition: Common.h:9
point
*vegas h *****************************************************used in the default bin number in original ***version of VEGAS is ***a higher bin number might help to derive a more precise ***grade subtle point
Definition: invegas.h:5
Local3DPoint
Point3DBase< float, LocalTag > Local3DPoint
Definition: LocalPoint.h:9