CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
PPSPixelSD.h
Go to the documentation of this file.
1 #ifndef _PPSPixelSD_h
2 #define _PPSPixelSD_h
3 // -*- C++ -*-
4 //
5 // Package: PPS
6 // Class : PPSPixelSD
7 //
16 //
17 // Original Author:
18 // Created: Tue May 16 10:14:34 CEST 2006
19 //
20 
21 // system include files
22 
23 // user include files
24 
29 
34 
37 
38 #include <string>
39 
40 class TrackingSlaveSD;
41 class SimTrackManager;
42 class G4Step;
43 class G4StepPoint;
44 class G4Track;
45 
47  public Observer<const BeginOfEvent*>,
48  public Observer<const EndOfEvent*> {
49 public:
51  ~PPSPixelSD() override;
52 
53  bool ProcessHits(G4Step*, G4TouchableHistory*) override;
54  uint32_t setDetUnitId(const G4Step*) override;
55 
56  void Initialize(G4HCofThisEvent* HCE) override;
57  void EndOfEvent(G4HCofThisEvent* eventHC) override;
58  void clear() override;
59  void DrawAll() override;
60  void PrintAll() override;
61 
62  void fillHits(edm::PSimHitContainer&, const std::string&) override;
63 
64 private:
65  void update(const BeginOfEvent*) override;
66  void update(const ::EndOfEvent*) override;
67  void clearHits() override;
68 
69 private:
70  static constexpr unsigned int maxPixelHits_ = 15000;
71  G4ThreeVector setToLocal(const G4ThreeVector& globalPoint);
72  void stepInfo(const G4Step* aStep);
73  bool hitExists();
74  void createNewHit();
75  void updateHit();
76  void storeHit(PPSPixelG4Hit*);
77  void resetForNewPrimary();
78  void summarize();
79 
80 private:
81  std::unique_ptr<TrackingSlaveSD> slave_;
82  std::unique_ptr<PPSVDetectorOrganization> numberingScheme_;
83 
84  // Data relative to primary particle (the one which triggers a shower)
85  // These data are common to all Hits of a given shower.
86  // One shower is made of several hits which differ by the
87  // unit ID (cristal/fiber/scintillator) and the Time slice ID.
88 
89  G4ThreeVector entrancePoint_;
91  G4int primID_; //@@ ID of the primary particle.
92 
94  G4int hcID_;
97 
98  int tsID_;
100  G4Track* theTrack_;
101  G4VPhysicalVolume* currentPV_;
104  double tSlice_;
105 
106  G4StepPoint* preStepPoint_;
107  G4StepPoint* postStepPoint_;
108  float edeposit_;
109  G4ThreeVector hitPoint_;
110 
111  G4ThreeVector position_;
112  G4ThreeVector theEntryPoint_;
113  G4ThreeVector theExitPoint_;
114  float Pabs_;
115  float Tof_;
116  float Eloss_;
118 
120  float PhiAtEntry_;
121 
123  float Vx_, Vy_, Vz_;
124 
125  int eventno_;
126 };
127 
128 #endif
void DrawAll() override
Definition: PPSPixelSD.cc:131
PPSPixelSD(const std::string &, const SensitiveDetectorCatalog &, edm::ParameterSet const &, SimTrackManager const *)
Definition: PPSPixelSD.cc:37
uint32_t previousUnitID_
Definition: PPSPixelSD.h:102
G4ThreeVector setToLocal(const G4ThreeVector &globalPoint)
Definition: PPSPixelSD.cc:154
G4ThreeVector theEntryPoint_
Definition: PPSPixelSD.h:112
uint32_t setDetUnitId(const G4Step *) override
Definition: PPSPixelSD.cc:83
void stepInfo(const G4Step *aStep)
Definition: PPSPixelSD.cc:161
float Pabs_
Definition: PPSPixelSD.h:114
G4int hcID_
Definition: PPSPixelSD.h:94
void clearHits() override
Definition: PPSPixelSD.cc:152
float Tof_
Definition: PPSPixelSD.h:115
void PrintAll() override
Definition: PPSPixelSD.cc:133
G4VPhysicalVolume * currentPV_
Definition: PPSPixelSD.h:101
float Vy_
Definition: PPSPixelSD.h:123
const SimTrackManager * theManager_
Definition: PPSPixelSD.h:96
G4StepPoint * preStepPoint_
Definition: PPSPixelSD.h:106
short ParticleType_
Definition: PPSPixelSD.h:117
void createNewHit()
Definition: PPSPixelSD.cc:246
void storeHit(PPSPixelG4Hit *)
Definition: PPSPixelSD.cc:307
int ParentId_
Definition: PPSPixelSD.h:122
bool ProcessHits(G4Step *, G4TouchableHistory *) override
Definition: PPSPixelSD.cc:73
float Vx_
Definition: PPSPixelSD.h:123
void update(const BeginOfEvent *) override
This routine will be called when the appropriate signal arrives.
Definition: PPSPixelSD.cc:144
int eventno_
Definition: PPSPixelSD.h:125
int primaryID_
Definition: PPSPixelSD.h:103
PPSPixelG4Hit * currentHit_
Definition: PPSPixelSD.h:99
G4int primID_
Definition: PPSPixelSD.h:91
void resetForNewPrimary()
Definition: PPSPixelSD.cc:318
G4Track * theTrack_
Definition: PPSPixelSD.h:100
bool hitExists()
Definition: PPSPixelSD.cc:211
void Initialize(G4HCofThisEvent *HCE) override
Definition: PPSPixelSD.cc:87
void clear() override
Definition: PPSPixelSD.cc:129
float incidentEnergy_
Definition: PPSPixelSD.h:90
std::unique_ptr< TrackingSlaveSD > slave_
Definition: PPSPixelSD.h:81
std::string name_
Definition: PPSPixelSD.h:93
PPSPixelG4HitCollection * theHC_
Definition: PPSPixelSD.h:95
double tSlice_
Definition: PPSPixelSD.h:104
float Eloss_
Definition: PPSPixelSD.h:116
float PhiAtEntry_
Definition: PPSPixelSD.h:120
std::unique_ptr< PPSVDetectorOrganization > numberingScheme_
Definition: PPSPixelSD.h:82
G4ThreeVector theExitPoint_
Definition: PPSPixelSD.h:113
void fillHits(edm::PSimHitContainer &, const std::string &) override
Definition: PPSPixelSD.cc:138
~PPSPixelSD() override
Definition: PPSPixelSD.cc:71
G4ThreeVector position_
Definition: PPSPixelSD.h:111
int tSliceID_
Definition: PPSPixelSD.h:103
G4ThreeVector hitPoint_
Definition: PPSPixelSD.h:109
static constexpr unsigned int maxPixelHits_
Definition: PPSPixelSD.h:70
G4StepPoint * postStepPoint_
Definition: PPSPixelSD.h:107
std::vector< PSimHit > PSimHitContainer
uint32_t unitID_
Definition: PPSPixelSD.h:102
void EndOfEvent(G4HCofThisEvent *eventHC) override
Definition: PPSPixelSD.cc:101
float edeposit_
Definition: PPSPixelSD.h:108
void updateHit()
Definition: PPSPixelSD.cc:293
G4ThreeVector entrancePoint_
Definition: PPSPixelSD.h:89
float ThetaAtEntry_
Definition: PPSPixelSD.h:119
void summarize()
Definition: PPSPixelSD.cc:324
float Vz_
Definition: PPSPixelSD.h:123