CMS 3D CMS Logo

SensitiveDetector.h
Go to the documentation of this file.
1 #ifndef SimG4Core_SensitiveDetector_H
2 #define SimG4Core_SensitiveDetector_H
3 
5 
7 //#include "DataFormats/GeometryVector/interface/LocalVector.h"
9 
10 #include "G4VSensitiveDetector.hh"
11 
12 #include <string>
13 
14 class G4Step;
15 class G4HCofThisEvent;
16 class G4TouchableHistory;
17 class G4VPhysicalVolume;
18 class DDCompactView;
19 
20 class SensitiveDetector : public G4VSensitiveDetector
21 {
22 public:
23  explicit SensitiveDetector(const std::string & iname,
24  const DDCompactView & cpv,
25  const SensitiveDetectorCatalog & ,
26  edm::ParameterSet const & p);
27  ~SensitiveDetector() override;
28 
29  void Initialize(G4HCofThisEvent * eventHC) override;
30  G4bool ProcessHits(G4Step * step ,G4TouchableHistory * tHistory) override = 0;
31  void EndOfEvent(G4HCofThisEvent * eventHC) override;
32 
33  virtual uint32_t setDetUnitId(const G4Step * step) = 0;
34  virtual void clearHits() = 0;
35 
36  inline const std::vector<std::string>& getNames() const { return namesOfSD; }
37 
38 protected:
39 
41  Local3DPoint InitialStepPosition(const G4Step * step, coordinates) const;
42  Local3DPoint FinalStepPosition(const G4Step * step, coordinates) const;
43 
44  Local3DPoint LocalPreStepPosition(const G4Step * step) const;
45  Local3DPoint LocalPostStepPosition(const G4Step * step) const;
46 
47  inline Local3DPoint ConvertToLocal3DPoint(const G4ThreeVector& point) const
48  { return Local3DPoint(point.x(),point.y(),point.z()); }
49 
50  void setNames(const std::vector<std::string>&);
51  void NaNTrap(const G4Step* step ) const;
52 
53 private:
54 
55  void AssignSD(const std::string & vname);
56 
57  std::vector<std::string> namesOfSD;
58 };
59 
60 #endif
SensitiveDetector(const std::string &iname, const DDCompactView &cpv, const SensitiveDetectorCatalog &, edm::ParameterSet const &p)
Local3DPoint ConvertToLocal3DPoint(const G4ThreeVector &point) const
virtual void clearHits()=0
~SensitiveDetector() override
G4bool ProcessHits(G4Step *step, G4TouchableHistory *tHistory) override=0
virtual uint32_t setDetUnitId(const G4Step *step)=0
type of data representation of DDCompactView
Definition: DDCompactView.h:90
void AssignSD(const std::string &vname)
const std::vector< std::string > & getNames() const
Local3DPoint FinalStepPosition(const G4Step *step, coordinates) const
std::vector< std::string > namesOfSD
Point3DBase< float, LocalTag > Local3DPoint
Definition: LocalPoint.h:9
void EndOfEvent(G4HCofThisEvent *eventHC) override
Local3DPoint LocalPostStepPosition(const G4Step *step) const
void setNames(const std::vector< std::string > &)
Local3DPoint LocalPreStepPosition(const G4Step *step) const
step
void Initialize(G4HCofThisEvent *eventHC) override
*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
void NaNTrap(const G4Step *step) const
Local3DPoint InitialStepPosition(const G4Step *step, coordinates) const