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
virtual void clearHits()=0
~SensitiveDetector() override
G4bool ProcessHits(G4Step *step, G4TouchableHistory *tHistory) override=0
SensitiveDetector(const std::string &iname, const SensitiveDetectorCatalog &, bool calo)
virtual uint32_t setDetUnitId(const G4Step *step)=0
void AssignSD(const std::string &vname)
Local3DPoint FinalStepPosition(const G4Step *step, coordinates) const
std::vector< std::string > m_namesOfSD
Local3DPoint ConvertToLocal3DPoint(const G4ThreeVector &point) const
Local3DPoint InitialStepPosition(const G4Step *step, coordinates) const
bool isCaloSD() const
TrackInformation * cmsTrackInformation(const G4Track *aTrack)
Local3DPoint LocalPreStepPosition(const G4Step *step) const
void EndOfEvent(G4HCofThisEvent *eventHC) override
void setNames(const std::vector< std::string > &)
const std::vector< std::string > & getNames() const
Local3DPoint LocalPostStepPosition(const G4Step *step) const
step
Definition: StallMonitor.cc:98
void Initialize(G4HCofThisEvent *eventHC) override
void NaNTrap(const G4Step *step) const
Definition: Common.h:9
*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
Point3DBase< float, LocalTag > Local3DPoint
Definition: LocalPoint.h:9