CMS 3D CMS Logo

PPSPixelDigiAnalyzer.cc
Go to the documentation of this file.
1 #ifndef SimPPS_RPIXDigiAnalyzer_h
2 #define SimPPS_RPIXDigiAnalyzer_h
3 
14 
24 
25 #include <iostream>
26 #include <string>
27 
28 #include "TH2D.h"
29 
30 #define SELECTED_PIXEL_ROW 89
31 #define SELECTED_PIXEL_COLUMN 23
32 #define SELECTED_UNITID 2014314496
33 #define TG184 0.332655724
34 
35 #define USE_MIDDLE_OF_PIXEL_2
36 #define CENTERX 1.05
37 #define CENTERY -8.475
38 
39 using namespace edm;
40 using namespace std;
41 
42 class PSimHit;
43 
44 namespace edm {
45  class ParameterSet;
46  class Event;
47  class EventSetup;
48 } // namespace edm
49 
50 class PPSPixelDigiAnalyzer : public edm::one::EDAnalyzer<edm::one::SharedResources> {
51 public:
53  ~PPSPixelDigiAnalyzer() override;
54  void endJob() override;
55  void beginJob() override;
56  void analyze(const edm::Event &event, const edm::EventSetup &eventSetup) override;
57 
58 private:
59  TH2D *hAllHits;
64  //TFile *file;
66 
70 
73  unsigned int cumulative_cluster_size_[3];
74 };
75 
77  : hAllHits(nullptr),
78  hOneHitperEvent(nullptr),
79  hOneHitperEvent2(nullptr),
80  hOneHitperEventCenter(nullptr),
81  hOneHitperEvent2Center(nullptr),
82  theRPixDetTopology_() {
83  label_ = pset.getUntrackedParameter<string>("label");
84  verbosity_ = pset.getParameter<int>("Verbosity");
86 #ifdef USE_MIDDLE_OF_PIXEL
87  hOneHitperEvent = file->make<TH2D>("OneHitperEvent", "One Hit per Event", 30, -8.511, -8.361, 20, 1, 1.1);
88  hOneHitperEvent2 = file->make<TH2D>("OneHitperEvent2", "One Hit per Event 2", 30, -8.511, -8.361, 20, 1, 1.1);
89 #else
90  hOneHitperEvent = file->make<TH2D>("OneHitperEvent", "One Hit per Event", 30, -8.55, -8.4, 20, 1, 1.1);
91  hOneHitperEvent2 = file->make<TH2D>("OneHitperEvent2", "One Hit per Event 2", 30, -8.55, -8.4, 20, 1, 1.1);
93  file->make<TH2D>("OneHitperEventCenter", "One Hit per Event Center", 30, -0.075, 0.075, 20, -0.05, 0.05);
95  file->make<TH2D>("OneHitperEvent2Center", "Cluster Size 2", 30, -0.075, 0.075, 20, -0.05, 0.05);
96 #endif
97  file->cd();
98  hAllHits = file->make<TH2D>("AllHits", "All Hits", 10, 1, 1.1, 10, -8.55, -8.4);
99 
100  psim_token = consumes<PSimHitContainer>(edm::InputTag("g4SimHits", "CTPPSPixelHits"));
101  pixel_token = consumes<edm::DetSetVector<CTPPSPixelDigi>>(edm::InputTag(label_, "")); //label=RPixDetDigitizer???
102 }
103 
105 
108  for (int a = 0; a < 3; a++)
110 }
112  edm::LogInfo("PPSPixelDigiAnalyzer") << "found_corresponding_digi_count_: " << found_corresponding_digi_count_;
113  edm::LogInfo("PPSPixelDigiAnalyzer") << "Cumulative cluster size (1,2,>2) = " << cumulative_cluster_size_[0] << ", "
115 }
116 
117 void PPSPixelDigiAnalyzer::analyze(const Event &event, const EventSetup &eventSetup) {
118  if (verbosity_ > 0)
119  edm::LogInfo("PPSPixelDigiAnalyzer") << "--- Run: " << event.id().run() << " Event: " << event.id().event();
120 
121  edm::LogInfo("PPSPixelDigiAnalyzer")
122  << " I do love Pixels ";
124  event.getByToken(psim_token, simHits);
125 
127  event.getByToken(pixel_token, CTPPSPixelDigis);
128 
129  if (verbosity_ > 0)
130  edm::LogInfo("PPSPixelDigiAnalyzer") << "\n=================== RPDA Starting SimHit access"
131  << " ===================";
132 
133  if (verbosity_ > 1)
134  edm::LogInfo("PPSPixelDigiAnalyzer") << simHits->size();
135 
136  double selected_pixel_lower_x;
137  double selected_pixel_lower_y;
138  double selected_pixel_upper_x;
139  double selected_pixel_upper_y;
140  double myX = 0;
141  double myY = 0;
142 
145  selected_pixel_lower_x,
146  selected_pixel_upper_x,
147  selected_pixel_lower_y,
148  selected_pixel_upper_y);
149 
150  double hit_inside_selected_pixel[2];
151  bool found_hit_inside_selected_pixel = false;
152 
153  for (vector<PSimHit>::const_iterator hit = simHits->begin(); hit != simHits->end(); hit++) {
154  LocalPoint entryP = hit->entryPoint();
155  LocalPoint exitP = hit->exitPoint();
156  LocalPoint midP((entryP.x() + exitP.x()) / 2., (entryP.y() + exitP.y()) / 2.);
157 
158 #ifdef USE_MIDDLE_OF_PIXEL
159  if (entryP.x() > selected_pixel_lower_x && entryP.x() < selected_pixel_upper_x &&
160  entryP.y() > (selected_pixel_lower_y + 0.115 * TG184) && entryP.y() < (selected_pixel_upper_y + 0.115 * TG184)
161 #else
163  if (midP.x() > selected_pixel_lower_x && midP.x() < selected_pixel_upper_x && midP.y() > selected_pixel_lower_y &&
164  midP.y() < selected_pixel_upper_y
165 #else
166  if (entryP.x() > selected_pixel_lower_x && entryP.x() < selected_pixel_upper_x &&
167  entryP.y() > selected_pixel_lower_y && entryP.y() < selected_pixel_upper_y
168 #endif
169 #endif
170  && hit->detUnitId() == SELECTED_UNITID) {
171  hit_inside_selected_pixel[0] = entryP.x();
172  hit_inside_selected_pixel[1] = entryP.y();
173  found_hit_inside_selected_pixel = true;
174 #ifdef USE_MIDDLE_OF_PIXEL_2
175  hAllHits->Fill(midP.x(), midP.y());
176  myX = midP.x();
177  myY = midP.y();
178 #else
179  hAllHits->Fill(entryP.x(), entryP.y());
180  myX = entryP.x();
181  myY = entryP.y();
182 #endif
183  if (verbosity_ > 2)
184  edm::LogInfo("PPSPixelDigiAnalyzer") << hit_inside_selected_pixel[0] << " " << hit_inside_selected_pixel[1];
185  }
186 
187  //--------------
188 
189  if (verbosity_ > 1)
190  if (hit->timeOfFlight() > 0) {
191  edm::LogInfo("PPSPixelDigiAnalyzer")
192  << "DetId: " << hit->detUnitId() << "PID: " << hit->particleType() << " TOF: " << hit->timeOfFlight()
193  << " Proc Type: " << hit->processType() << " p: " << hit->pabs() << " x = " << entryP.x()
194  << " y = " << entryP.y() << " z = " << entryP.z();
195  }
196  }
197 
198  if (verbosity_ > 0)
199  edm::LogInfo("PPSPixelDigiAnalyzer") << "\n=================== RPDA Starting Digi access"
200  << " ===================";
201  int numberOfDetUnits = 0;
202 
203  // Iterate on detector units
204  edm::DetSetVector<CTPPSPixelDigi>::const_iterator DSViter = CTPPSPixelDigis->begin();
205 
206  for (; DSViter != CTPPSPixelDigis->end(); DSViter++) {
207  ++numberOfDetUnits;
208 
209  DetId detIdObject(DSViter->detId());
210  if (verbosity_ > 1)
211  edm::LogInfo("PPSPixelDigiAnalyzer") << "DetId: " << DSViter->detId();
212 
213  bool found_corresponding_digi = false;
214  unsigned int corresponding_digi_cluster_size = 0;
215 
216  // looping over digis in a unit id
219 
220  if (verbosity_ > 2) {
221  edm::LogInfo("PPSPixelDigiAnalyzer") << "FF " << DSViter->detId();
222  for (edm::DetSet<CTPPSPixelDigi>::const_iterator di = begin; di != end; di++) {
223  edm::LogInfo("PPSPixelDigiAnalyzer") << " Digi row " << di->row() << ", col " << di->column();
224 
225  // reconvert the digi to local coordinates
226  double lx;
227  double ly;
228  double ux;
229  double uy;
230  unsigned int rr = di->row();
231  unsigned int cc = di->column();
232  theRPixDetTopology_.pixelRange(rr, cc, lx, ux, ly, uy);
233 
234  edm::LogInfo("PPSPixelDigiAnalyzer")
235  << " pixel boundaries x low up, y low up " << lx << " " << ux << " " << ly << " " << uy;
236  }
237  }
238  if (DSViter->detId() == SELECTED_UNITID && found_hit_inside_selected_pixel) {
239  for (edm::DetSet<CTPPSPixelDigi>::const_iterator di = begin; di != end; di++) {
240  if (verbosity_ > 1)
241  edm::LogInfo("PPSPixelDigiAnalyzer") << " Digi row " << di->row() << ", col " << di->column();
242 
243  if (di->row() == SELECTED_PIXEL_ROW && di->column() == SELECTED_PIXEL_COLUMN) {
245  found_corresponding_digi = true;
246  corresponding_digi_cluster_size = 1;
247  }
248  }
249  //if coresponding digi found, re-loop to look for adjacent pixels
250  if (found_corresponding_digi) {
251  for (edm::DetSet<CTPPSPixelDigi>::const_iterator di = begin; di != end; di++) {
252  if (verbosity_ > 1)
253  edm::LogInfo("PPSPixelDigiAnalyzer") << " Digi row " << di->row() << ", col " << di->column();
254 
255  if ((di->row() == SELECTED_PIXEL_ROW + 1 && di->column() == SELECTED_PIXEL_COLUMN) ||
256  (di->row() == SELECTED_PIXEL_ROW - 1 && di->column() == SELECTED_PIXEL_COLUMN) ||
257  (di->row() == SELECTED_PIXEL_ROW && di->column() == SELECTED_PIXEL_COLUMN + 1) ||
258  (di->row() == SELECTED_PIXEL_ROW && di->column() == SELECTED_PIXEL_COLUMN - 1)) {
259  corresponding_digi_cluster_size++;
260  edm::LogInfo("PPSPixelDigiAnalyzer") << " Digi row " << di->row() << ", col " << di->column();
261  }
262  }
263  }
264  }
265  if (corresponding_digi_cluster_size > 0) {
266  edm::LogInfo("PPSPixelDigiAnalyzer")
267  << "corresponding_digi_cluster_size in the event: " << corresponding_digi_cluster_size;
268  hOneHitperEvent->Fill(myY, myX);
269  hOneHitperEventCenter->Fill(myY - CENTERY, myX - CENTERX);
270  if (corresponding_digi_cluster_size < 3) {
271  cumulative_cluster_size_[corresponding_digi_cluster_size - 1]++;
272  if (corresponding_digi_cluster_size > 1) {
273  hOneHitperEvent2->Fill(myY, myX);
274  hOneHitperEvent2Center->Fill(myY - CENTERY, myX - CENTERX);
275  }
276  } else {
278  hOneHitperEvent2->Fill(myY, myX);
279  hOneHitperEvent2Center->Fill(myY - CENTERY, myX - CENTERX);
280  }
281  }
282  }
283 
284  if (verbosity_ > 1)
285  edm::LogInfo("PPSPixelDigiAnalyzer") << "numberOfDetUnits in the event: " << numberOfDetUnits;
286 }
287 
290 
292 
293 #endif
T getParameter(std::string const &) const
unsigned int cumulative_cluster_size_[3]
T getUntrackedParameter(std::string const &, T const &) const
CTPPSPixelSimTopology theRPixDetTopology_
#define nullptr
T y() const
Definition: PV3DBase.h:60
example_stream void analyze(const edm::Event &, const edm::EventSetup &) override
unsigned int found_corresponding_digi_count_
void beginJob()
Definition: Breakpoints.cc:14
#define SELECTED_PIXEL_ROW
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
#define TG184
#define CENTERX
T z() const
Definition: PV3DBase.h:61
edm::EDGetTokenT< edm::DetSetVector< CTPPSPixelDigi > > pixel_token
#define end
Definition: vmac.h:39
void pixelRange(unsigned int arow, unsigned int acol, double &lower_x, double &higher_x, double &lower_y, double &higher_y) const
#define CENTERY
Definition: DetId.h:17
#define SELECTED_PIXEL_COLUMN
edm::EDGetTokenT< edm::PSimHitContainer > psim_token
void analyze(const edm::Event &event, const edm::EventSetup &eventSetup) override
#define SELECTED_UNITID
#define USE_MIDDLE_OF_PIXEL_2
#define begin
Definition: vmac.h:32
HLT enums.
double a
Definition: hdecay.h:119
PPSPixelDigiAnalyzer(const edm::ParameterSet &pset)
collection_type::const_iterator const_iterator
Definition: DetSet.h:32
collection_type::const_iterator const_iterator
Definition: DetSetVector.h:102
T x() const
Definition: PV3DBase.h:59
Definition: event.py:1