CMS 3D CMS Logo

BtlLocalRecoValidation.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: Validation/MtdValidation
4 // Class: BtlLocalRecoValidation
5 //
14 #include <string>
15 
20 
23 
30 
34 
39 
42 
44 
45 struct MTDHit {
46  float energy;
47  float time;
48  float x_local;
49  float y_local;
50  float z_local;
51 };
52 
54 public:
56  ~BtlLocalRecoValidation() override;
57 
58  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
59 
60 private:
61  void bookHistograms(DQMStore::IBooker&, edm::Run const&, edm::EventSetup const&) override;
62 
63  void analyze(const edm::Event&, const edm::EventSetup&) override;
64 
65  bool isSameCluster(const FTLCluster&, const FTLCluster&);
66 
67  // ------------ member data ------------
68 
70  const double hitMinEnergy_;
71  const bool LocalPosDebug_;
73  const double hitMinAmplitude_;
74 
80 
83 
84  // --- histograms declaration
85 
87 
91 
93 
94  //local position monitoring
99 
105 
115 
118 
122 
125 
135 
154 
156 
157  // --- UncalibratedRecHits histograms
158 
159  static constexpr int nBinsQ_ = 20;
160  static constexpr float binWidthQ_ = 30.;
161  static constexpr int nBinsQEta_ = 3;
162  static constexpr float binsQEta_[nBinsQEta_ + 1] = {0., 0.65, 1.15, 1.55};
163 
166 
167  static constexpr int nBinsEta_ = 31;
168  static constexpr float binWidthEta_ = 0.05;
169  static constexpr int nBinsEtaQ_ = 7;
170  static constexpr float binsEtaQ_[nBinsEtaQ_ + 1] = {0., 30., 60., 90., 120., 150., 360., 600.};
171 
174 };
175 
177  return clu1.id() == clu2.id() && clu1.size() == clu2.size() && clu1.x() == clu2.x() && clu1.y() == clu2.y() &&
178  clu1.time() == clu2.time();
179 }
180 
181 // ------------ constructor and destructor --------------
183  : folder_(iConfig.getParameter<std::string>("folder")),
184  hitMinEnergy_(iConfig.getParameter<double>("HitMinimumEnergy")),
185  LocalPosDebug_(iConfig.getParameter<bool>("LocalPositionDebug")),
186  uncalibRecHitsPlots_(iConfig.getParameter<bool>("UncalibRecHitsPlots")),
187  hitMinAmplitude_(iConfig.getParameter<double>("HitMinimumAmplitude")) {
188  btlRecHitsToken_ = consumes<FTLRecHitCollection>(iConfig.getParameter<edm::InputTag>("recHitsTag"));
191  consumes<FTLUncalibratedRecHitCollection>(iConfig.getParameter<edm::InputTag>("uncalibRecHitsTag"));
192  btlSimHitsToken_ = consumes<CrossingFrame<PSimHit> >(iConfig.getParameter<edm::InputTag>("simHitsTag"));
193  btlRecCluToken_ = consumes<FTLClusterCollection>(iConfig.getParameter<edm::InputTag>("recCluTag"));
194  mtdTrackingHitToken_ = consumes<MTDTrackingDetSetVector>(iConfig.getParameter<edm::InputTag>("trkHitTag"));
195 
196  mtdgeoToken_ = esConsumes<MTDGeometry, MTDDigiGeometryRecord>();
197  mtdtopoToken_ = esConsumes<MTDTopology, MTDTopologyRcd>();
198 }
199 
201 
202 // ------------ method called for each event ------------
204  using namespace edm;
205  using namespace std;
206  using namespace geant_units::operators;
207 
208  auto geometryHandle = iSetup.getTransientHandle(mtdgeoToken_);
209  const MTDGeometry* geom = geometryHandle.product();
210 
211  auto topologyHandle = iSetup.getTransientHandle(mtdtopoToken_);
212  const MTDTopology* topology = topologyHandle.product();
213 
214  auto btlRecHitsHandle = makeValid(iEvent.getHandle(btlRecHitsToken_));
215  auto btlSimHitsHandle = makeValid(iEvent.getHandle(btlSimHitsToken_));
216  auto btlRecCluHandle = makeValid(iEvent.getHandle(btlRecCluToken_));
217  auto mtdTrkHitHandle = makeValid(iEvent.getHandle(mtdTrackingHitToken_));
218  MixCollection<PSimHit> btlSimHits(btlSimHitsHandle.product());
219 
220  // --- Loop over the BTL SIM hits
221  std::unordered_map<uint32_t, MTDHit> m_btlSimHits;
222  for (auto const& simHit : btlSimHits) {
223  // --- Use only hits compatible with the in-time bunch-crossing
224  if (simHit.tof() < 0 || simHit.tof() > 25.)
225  continue;
226 
227  DetId id = simHit.detUnitId();
228 
229  auto simHitIt = m_btlSimHits.emplace(id.rawId(), MTDHit()).first;
230 
231  // --- Accumulate the energy (in MeV) of SIM hits in the same detector cell
232  (simHitIt->second).energy += convertUnitsTo(0.001_MeV, simHit.energyLoss());
233 
234  // --- Get the time of the first SIM hit in the cell
235  if ((simHitIt->second).time == 0 || simHit.tof() < (simHitIt->second).time) {
236  (simHitIt->second).time = simHit.tof();
237 
238  auto hit_pos = simHit.entryPoint();
239  (simHitIt->second).x_local = hit_pos.x();
240  (simHitIt->second).y_local = hit_pos.y();
241  (simHitIt->second).z_local = hit_pos.z();
242  }
243 
244  } // simHit loop
245 
246  // --- Loop over the BTL RECO hits
247  unsigned int n_reco_btl = 0;
248  for (const auto& recHit : *btlRecHitsHandle) {
249  BTLDetId detId = recHit.id();
251  const MTDGeomDet* thedet = geom->idToDet(geoId);
252  if (thedet == nullptr)
253  throw cms::Exception("BtlLocalRecoValidation") << "GeographicalID: " << std::hex << geoId.rawId() << " ("
254  << detId.rawId() << ") is invalid!" << std::dec << std::endl;
255  const ProxyMTDTopology& topoproxy = static_cast<const ProxyMTDTopology&>(thedet->topology());
256  const RectangularMTDTopology& topo = static_cast<const RectangularMTDTopology&>(topoproxy.specificTopology());
257 
258  Local3DPoint local_point(0., 0., 0.);
259  local_point = topo.pixelToModuleLocalPoint(local_point, detId.row(topo.nrows()), detId.column(topo.nrows()));
260  const auto& global_point = thedet->toGlobal(local_point);
261 
262  meHitEnergy_->Fill(recHit.energy());
263  meHitTime_->Fill(recHit.time());
264  meHitTimeError_->Fill(recHit.timeError());
265  meHitLongPos_->Fill(recHit.position());
266  meHitLongPosErr_->Fill(recHit.positionError());
267 
268  meOccupancy_->Fill(global_point.z(), global_point.phi());
269 
270  if (LocalPosDebug_) {
271  meLocalOccupancy_->Fill(local_point.x() + recHit.position(), local_point.y());
272  meHitXlocal_->Fill(local_point.x());
273  meHitYlocal_->Fill(local_point.y());
274  meHitZlocal_->Fill(local_point.z());
275  }
276  meHitX_->Fill(global_point.x());
277  meHitY_->Fill(global_point.y());
278  meHitZ_->Fill(global_point.z());
279  meHitPhi_->Fill(global_point.phi());
280  meHitEta_->Fill(global_point.eta());
281 
282  meHitTvsE_->Fill(recHit.energy(), recHit.time());
283  meHitEvsPhi_->Fill(global_point.phi(), recHit.energy());
284  meHitEvsEta_->Fill(global_point.eta(), recHit.energy());
285  meHitEvsZ_->Fill(global_point.z(), recHit.energy());
286  meHitTvsPhi_->Fill(global_point.phi(), recHit.time());
287  meHitTvsEta_->Fill(global_point.eta(), recHit.time());
288  meHitTvsZ_->Fill(global_point.z(), recHit.time());
289 
290  // Resolution histograms
291  if (m_btlSimHits.count(detId.rawId()) == 1 && m_btlSimHits[detId.rawId()].energy > hitMinEnergy_) {
292  float longpos_res = recHit.position() - convertMmToCm(m_btlSimHits[detId.rawId()].x_local);
293  float time_res = recHit.time() - m_btlSimHits[detId.rawId()].time;
294  float energy_res = recHit.energy() - m_btlSimHits[detId.rawId()].energy;
295 
296  Local3DPoint local_point_sim(convertMmToCm(m_btlSimHits[detId.rawId()].x_local),
297  convertMmToCm(m_btlSimHits[detId.rawId()].y_local),
298  convertMmToCm(m_btlSimHits[detId.rawId()].z_local));
299  local_point_sim =
300  topo.pixelToModuleLocalPoint(local_point_sim, detId.row(topo.nrows()), detId.column(topo.nrows()));
301  const auto& global_point_sim = thedet->toGlobal(local_point_sim);
302 
303  meTimeRes_->Fill(time_res);
304  meEnergyRes_->Fill(energy_res);
305 
306  meLongPosPull_->Fill(longpos_res / recHit.positionError());
307  meLongPosPullvsEta_->Fill(std::abs(global_point_sim.eta()), longpos_res / recHit.positionError());
308  meLongPosPullvsE_->Fill(m_btlSimHits[detId.rawId()].energy, longpos_res / recHit.positionError());
309 
310  meTPullvsEta_->Fill(std::abs(global_point_sim.eta()), time_res / recHit.timeError());
311  meTPullvsE_->Fill(m_btlSimHits[detId.rawId()].energy, time_res / recHit.timeError());
312  }
313 
314  n_reco_btl++;
315 
316  } // recHit loop
317 
318  if (n_reco_btl > 0)
319  meNhits_->Fill(log10(n_reco_btl));
320 
321  // --- Loop over the BTL RECO clusters ---
322  for (const auto& DetSetClu : *btlRecCluHandle) {
323  for (const auto& cluster : DetSetClu) {
324  if (cluster.energy() < hitMinEnergy_)
325  continue;
326  BTLDetId cluId = cluster.id();
327  DetId detIdObject(cluId);
328  const auto& genericDet = geom->idToDetUnit(detIdObject);
329  if (genericDet == nullptr) {
330  throw cms::Exception("BtlLocalRecoValidation")
331  << "GeographicalID: " << std::hex << cluId << " is invalid!" << std::dec << std::endl;
332  }
333 
334  const ProxyMTDTopology& topoproxy = static_cast<const ProxyMTDTopology&>(genericDet->topology());
335  const RectangularMTDTopology& topo = static_cast<const RectangularMTDTopology&>(topoproxy.specificTopology());
336 
337  // --- Cluster position in the module reference frame
338  Local3DPoint local_point(topo.localX(cluster.x()), topo.localY(cluster.y()), 0.);
339  const auto& global_point = genericDet->toGlobal(local_point);
340 
341  meCluEnergy_->Fill(cluster.energy());
342  meCluTime_->Fill(cluster.time());
343  meCluTimeError_->Fill(cluster.timeError());
344  meCluPhi_->Fill(global_point.phi());
345  meCluEta_->Fill(global_point.eta());
346  meCluZvsPhi_->Fill(global_point.z(), global_point.phi());
347  meCluHits_->Fill(cluster.size());
348 
349  // --- Get the SIM hits associated to the cluster and calculate
350  // the cluster SIM energy, time and position
351 
352  double cluEneSIM = 0.;
353  double cluTimeSIM = 0.;
354  double cluLocXSIM = 0.;
355  double cluLocYSIM = 0.;
356  double cluLocZSIM = 0.;
357 
358  for (int ihit = 0; ihit < cluster.size(); ++ihit) {
359  int hit_row = cluster.minHitRow() + cluster.hitOffset()[ihit * 2];
360  int hit_col = cluster.minHitCol() + cluster.hitOffset()[ihit * 2 + 1];
361 
362  // Match the RECO hit to the corresponding SIM hit
363  for (const auto& recHit : *btlRecHitsHandle) {
364  BTLDetId hitId(recHit.id().rawId());
365 
366  if (m_btlSimHits.count(hitId.rawId()) == 0)
367  continue;
368 
369  // Check the hit position
370  if (hitId.mtdSide() != cluId.mtdSide() || hitId.mtdRR() != cluId.mtdRR() || recHit.row() != hit_row ||
371  recHit.column() != hit_col)
372  continue;
373 
374  // Check the hit energy and time
375  if (recHit.energy() != cluster.hitENERGY()[ihit] || recHit.time() != cluster.hitTIME()[ihit])
376  continue;
377 
378  // SIM hit's position in the module reference frame
379  Local3DPoint local_point_sim(convertMmToCm(m_btlSimHits[recHit.id().rawId()].x_local),
380  convertMmToCm(m_btlSimHits[recHit.id().rawId()].y_local),
381  convertMmToCm(m_btlSimHits[recHit.id().rawId()].z_local));
382  local_point_sim =
383  topo.pixelToModuleLocalPoint(local_point_sim, hitId.row(topo.nrows()), hitId.column(topo.nrows()));
384 
385  // Calculate the SIM cluster's position in the module reference frame
386  cluLocXSIM += local_point_sim.x() * m_btlSimHits[recHit.id().rawId()].energy;
387  cluLocYSIM += local_point_sim.y() * m_btlSimHits[recHit.id().rawId()].energy;
388  cluLocZSIM += local_point_sim.z() * m_btlSimHits[recHit.id().rawId()].energy;
389 
390  // Calculate the SIM cluster energy and time
391  cluEneSIM += m_btlSimHits[recHit.id().rawId()].energy;
392  cluTimeSIM += m_btlSimHits[recHit.id().rawId()].time * m_btlSimHits[recHit.id().rawId()].energy;
393 
394  } // recHit loop
395 
396  } // ihit loop
397 
398  // Find the MTDTrackingRecHit corresponding to the cluster
399  MTDTrackingRecHit* comp(nullptr);
400  bool matchClu = false;
401  const auto& trkHits = (*mtdTrkHitHandle)[detIdObject];
402  for (const auto& trkHit : trkHits) {
403  if (isSameCluster(trkHit.mtdCluster(), cluster)) {
404  comp = trkHit.clone();
405  matchClu = true;
406  break;
407  }
408  }
409  if (!matchClu) {
410  edm::LogWarning("BtlLocalRecoValidation")
411  << "No valid TrackingRecHit corresponding to cluster, detId = " << detIdObject.rawId();
412  }
413 
414  // --- Fill the cluster resolution histograms
415  if (cluTimeSIM > 0. && cluEneSIM > 0.) {
416  cluTimeSIM /= cluEneSIM;
417 
418  Local3DPoint cluLocalPosSIM(cluLocXSIM / cluEneSIM, cluLocYSIM / cluEneSIM, cluLocZSIM / cluEneSIM);
419  const auto& cluGlobalPosSIM = genericDet->toGlobal(cluLocalPosSIM);
420 
421  float time_res = cluster.time() - cluTimeSIM;
422  float energy_res = cluster.energy() - cluEneSIM;
423  meCluTimeRes_->Fill(time_res);
424  meCluEnergyRes_->Fill(energy_res);
425 
426  float rho_res = global_point.perp() - cluGlobalPosSIM.perp();
427  float phi_res = global_point.phi() - cluGlobalPosSIM.phi();
428 
429  meCluRhoRes_->Fill(rho_res);
430  meCluPhiRes_->Fill(phi_res);
431 
432  float xlocal_res = local_point.x() - cluLocalPosSIM.x();
433  float ylocal_res = local_point.y() - cluLocalPosSIM.y();
434  float z_res = global_point.z() - cluGlobalPosSIM.z();
435 
436  meCluZRes_->Fill(z_res);
437 
438  if (LocalPosDebug_) {
439  if (matchClu && comp != nullptr) {
440  meCluLocalXRes_->Fill(xlocal_res);
441  meCluLocalYRes_->Fill(ylocal_res);
442  meCluLocalXPull_->Fill(xlocal_res / std::sqrt(comp->localPositionError().xx()));
443  meCluLocalYPull_->Fill(ylocal_res / std::sqrt(comp->localPositionError().yy()));
444  meCluZPull_->Fill(z_res / std::sqrt(comp->globalPositionError().czz()));
445  meCluXLocalErr_->Fill(std::sqrt(comp->localPositionError().xx()));
446  meCluYLocalErr_->Fill(std::sqrt(comp->localPositionError().yy()));
447  }
448 
449  meCluYXLocal_->Fill(local_point.x(), local_point.y());
450  meCluYXLocalSim_->Fill(cluLocalPosSIM.x(), cluLocalPosSIM.y());
451  }
452 
453  meCluEnergyvsEta_->Fill(std::abs(cluGlobalPosSIM.eta()), cluster.energy());
454  meCluHitsvsEta_->Fill(std::abs(cluGlobalPosSIM.eta()), cluster.size());
455 
456  meCluTResvsEta_->Fill(std::abs(cluGlobalPosSIM.eta()), time_res);
457  meCluTResvsE_->Fill(cluEneSIM, time_res);
458 
459  meCluTPullvsEta_->Fill(std::abs(cluGlobalPosSIM.eta()), time_res / cluster.timeError());
460  meCluTPullvsE_->Fill(cluEneSIM, time_res / cluster.timeError());
461 
462  } // if ( cluTimeSIM > 0. && cluEneSIM > 0. )
463  else {
464  meUnmatchedCluEnergy_->Fill(std::log10(cluster.energy()));
465  }
466 
467  } // cluster loop
468 
469  } // DetSetClu loop
470 
471  // --- Loop over the BTL Uncalibrated RECO hits
472  if (uncalibRecHitsPlots_) {
473  auto btlUncalibRecHitsHandle = makeValid(iEvent.getHandle(btlUncalibRecHitsToken_));
474 
475  for (const auto& uRecHit : *btlUncalibRecHitsHandle) {
476  BTLDetId detId = uRecHit.id();
477 
478  // --- Skip UncalibratedRecHits not matched to SimHits
479  if (m_btlSimHits.count(detId.rawId()) != 1)
480  continue;
481 
483  const MTDGeomDet* thedet = geom->idToDet(geoId);
484  if (thedet == nullptr)
485  throw cms::Exception("BtlLocalRecoValidation") << "GeographicalID: " << std::hex << geoId.rawId() << " ("
486  << detId.rawId() << ") is invalid!" << std::dec << std::endl;
487  const ProxyMTDTopology& topoproxy = static_cast<const ProxyMTDTopology&>(thedet->topology());
488  const RectangularMTDTopology& topo = static_cast<const RectangularMTDTopology&>(topoproxy.specificTopology());
489 
490  Local3DPoint local_point(0., 0., 0.);
491  local_point = topo.pixelToModuleLocalPoint(local_point, detId.row(topo.nrows()), detId.column(topo.nrows()));
492  const auto& global_point = thedet->toGlobal(local_point);
493 
494  // --- Combine the information from the left and right BTL cell sides
495 
496  float nHits = 0.;
497  float hit_amplitude = 0.;
498  float hit_time = 0.;
499 
500  // left side:
501  if (uRecHit.amplitude().first > 0.) {
502  hit_amplitude += uRecHit.amplitude().first;
503  hit_time += uRecHit.time().first;
504  nHits += 1.;
505  }
506  // right side:
507  if (uRecHit.amplitude().second > 0.) {
508  hit_amplitude += uRecHit.amplitude().second;
509  hit_time += uRecHit.time().second;
510  nHits += 1.;
511  }
512 
513  hit_amplitude /= nHits;
514  hit_time /= nHits;
515 
516  // --- Fill the histograms
517 
518  if (hit_amplitude < hitMinAmplitude_)
519  continue;
520 
521  float time_res = hit_time - m_btlSimHits[detId.rawId()].time;
522 
523  // amplitude histograms
524 
525  int qBin = (int)(hit_amplitude / binWidthQ_);
526  if (qBin > nBinsQ_ - 1)
527  qBin = nBinsQ_ - 1;
528 
529  meTimeResQ_[qBin]->Fill(time_res);
530 
531  int etaBin = 0;
532  for (int ibin = 1; ibin < nBinsQEta_; ++ibin)
533  if (fabs(global_point.eta()) >= binsQEta_[ibin] && fabs(global_point.eta()) < binsQEta_[ibin + 1])
534  etaBin = ibin;
535 
536  meTimeResQvsEta_[qBin][etaBin]->Fill(time_res);
537 
538  // eta histograms
539 
540  etaBin = (int)(fabs(global_point.eta()) / binWidthEta_);
541  if (etaBin > nBinsEta_ - 1)
542  etaBin = nBinsEta_ - 1;
543 
544  meTimeResEta_[etaBin]->Fill(time_res);
545 
546  qBin = 0;
547  for (int ibin = 1; ibin < nBinsEtaQ_; ++ibin)
548  if (hit_amplitude >= binsEtaQ_[ibin] && hit_amplitude < binsEtaQ_[ibin + 1])
549  qBin = ibin;
550 
551  meTimeResEtavsQ_[etaBin][qBin]->Fill(time_res);
552 
553  } // uRecHit loop
554  }
555 }
556 
557 // ------------ method for histogram booking ------------
559  edm::Run const& run,
560  edm::EventSetup const& iSetup) {
561  ibook.setCurrentFolder(folder_);
562 
563  // --- histograms booking
564 
565  meNhits_ = ibook.book1D("BtlNhits", "Number of BTL RECO hits;log_{10}(N_{RECO})", 100, 0., 5.25);
566 
567  meHitEnergy_ = ibook.book1D("BtlHitEnergy", "BTL RECO hits energy;E_{RECO} [MeV]", 100, 0., 20.);
568  meHitTime_ = ibook.book1D("BtlHitTime", "BTL RECO hits ToA;ToA_{RECO} [ns]", 100, 0., 25.);
569  meHitTimeError_ = ibook.book1D("BtlHitTimeError", "BTL RECO hits ToA error;#sigma^{ToA}_{RECO} [ns]", 50, 0., 0.1);
570  meOccupancy_ = ibook.book2D(
571  "BtlOccupancy", "BTL RECO hits occupancy;Z_{RECO} [cm]; #phi_{RECO} [rad]", 65, -260., 260., 126, -3.2, 3.2);
572  if (LocalPosDebug_) {
573  meLocalOccupancy_ = ibook.book2D(
574  "BtlLocalOccupancy", "BTL RECO hits local occupancy;X_{RECO} [cm]; Y_{RECO} [cm]", 100, 10., 10., 60, -3., 3.);
575  meHitXlocal_ = ibook.book1D("BtlHitXlocal", "BTL RECO local X;X_{RECO}^{LOC} [cm]", 100, -10., 10.);
576  meHitYlocal_ = ibook.book1D("BtlHitYlocal", "BTL RECO local Y;Y_{RECO}^{LOC} [cm]", 60, -3, 3);
577  meHitZlocal_ = ibook.book1D("BtlHitZlocal", "BTL RECO local z;z_{RECO}^{LOC} [cm]", 10, -1, 1);
578  }
579  meHitX_ = ibook.book1D("BtlHitX", "BTL RECO hits X;X_{RECO} [cm]", 60, -120., 120.);
580  meHitY_ = ibook.book1D("BtlHitY", "BTL RECO hits Y;Y_{RECO} [cm]", 60, -120., 120.);
581  meHitZ_ = ibook.book1D("BtlHitZ", "BTL RECO hits Z;Z_{RECO} [cm]", 100, -260., 260.);
582  meHitPhi_ = ibook.book1D("BtlHitPhi", "BTL RECO hits #phi;#phi_{RECO} [rad]", 126, -3.2, 3.2);
583  meHitEta_ = ibook.book1D("BtlHitEta", "BTL RECO hits #eta;#eta_{RECO}", 100, -1.55, 1.55);
584  meHitTvsE_ =
585  ibook.bookProfile("BtlHitTvsE", "BTL RECO ToA vs energy;E_{RECO} [MeV];ToA_{RECO} [ns]", 50, 0., 20., 0., 100.);
586  meHitEvsPhi_ = ibook.bookProfile(
587  "BtlHitEvsPhi", "BTL RECO energy vs #phi;#phi_{RECO} [rad];E_{RECO} [MeV]", 50, -3.2, 3.2, 0., 100.);
588  meHitEvsEta_ = ibook.bookProfile(
589  "BtlHitEvsEta", "BTL RECO energy vs #eta;#eta_{RECO};E_{RECO} [MeV]", 50, -1.55, 1.55, 0., 100.);
590  meHitEvsZ_ =
591  ibook.bookProfile("BtlHitEvsZ", "BTL RECO energy vs Z;Z_{RECO} [cm];E_{RECO} [MeV]", 50, -260., 260., 0., 100.);
592  meHitTvsPhi_ = ibook.bookProfile(
593  "BtlHitTvsPhi", "BTL RECO ToA vs #phi;#phi_{RECO} [rad];ToA_{RECO} [ns]", 50, -3.2, 3.2, 0., 100.);
594  meHitTvsEta_ =
595  ibook.bookProfile("BtlHitTvsEta", "BTL RECO ToA vs #eta;#eta_{RECO};ToA_{RECO} [ns]", 50, -1.6, 1.6, 0., 100.);
596  meHitTvsZ_ =
597  ibook.bookProfile("BtlHitTvsZ", "BTL RECO ToA vs Z;Z_{RECO} [cm];ToA_{RECO} [ns]", 50, -260., 260., 0., 100.);
598  meHitLongPos_ = ibook.book1D("BtlLongPos", "BTL RECO hits longitudinal position;long. pos._{RECO}", 100, -10, 10);
600  ibook.book1D("BtlLongPosErr", "BTL RECO hits longitudinal position error; long. pos. error_{RECO}", 100, -1, 1);
601  meTimeRes_ = ibook.book1D("BtlTimeRes", "BTL time resolution;T_{RECO}-T_{SIM}", 100, -0.5, 0.5);
602  meEnergyRes_ = ibook.book1D("BtlEnergyRes", "BTL energy resolution;E_{RECO}-E_{SIM}", 100, -0.5, 0.5);
603  meLongPosPull_ = ibook.book1D("BtlLongPosPull",
604  "BTL longitudinal position pull;X^{loc}_{RECO}-X^{loc}_{SIM}/#sigma_{xloc_{RECO}}",
605  100,
606  -5.,
607  5.);
609  "BtlLongposPullvsE",
610  "BTL longitudinal position pull vs E;E_{SIM} [MeV];X^{loc}_{RECO}-X^{loc}_{SIM}/#sigma_{xloc_{RECO}}",
611  20,
612  0.,
613  20.,
614  -5.,
615  5.,
616  "S");
618  "BtlLongposPullvsEta",
619  "BTL longitudinal position pull vs #eta;|#eta_{RECO}|;X^{loc}_{RECO}-X^{loc}_{SIM}/#sigma_{xloc_{RECO}}",
620  32,
621  0,
622  1.55,
623  -5.,
624  5.,
625  "S");
626  meTPullvsE_ = ibook.bookProfile(
627  "BtlTPullvsE", "BTL time pull vs E;E_{SIM} [MeV];(T_{RECO}-T_{SIM})/#sigma_{T_{RECO}}", 20, 0., 20., -5., 5., "S");
628  meTPullvsEta_ = ibook.bookProfile("BtlTPullvsEta",
629  "BTL time pull vs #eta;|#eta_{RECO}|;(T_{RECO}-T_{SIM})/#sigma_{T_{RECO}}",
630  30,
631  0,
632  1.55,
633  -5.,
634  5.,
635  "S");
636  meCluTime_ = ibook.book1D("BtlCluTime", "BTL cluster time ToA;ToA [ns]", 250, 0, 25);
637  meCluTimeError_ = ibook.book1D("BtlCluTimeError", "BTL cluster time error;#sigma_{t} [ns]", 100, 0, 0.1);
638  meCluEnergy_ = ibook.book1D("BtlCluEnergy", "BTL cluster energy;E_{RECO} [MeV]", 100, 0, 20);
639  meCluPhi_ = ibook.book1D("BtlCluPhi", "BTL cluster #phi;#phi_{RECO} [rad]", 144, -3.2, 3.2);
640  meCluEta_ = ibook.book1D("BtlCluEta", "BTL cluster #eta;#eta_{RECO}", 100, -1.55, 1.55);
641  meCluHits_ = ibook.book1D("BtlCluHitNumber", "BTL hits per cluster; Cluster size", 10, 0, 10);
642  meCluZvsPhi_ = ibook.book2D(
643  "BtlOccupancy", "BTL cluster Z vs #phi;Z_{RECO} [cm]; #phi_{RECO} [rad]", 144, -260., 260., 50, -3.2, 3.2);
645  "BtlCluEnergyVsEta", "BTL cluster energy vs #eta; |#eta_{RECO}|; E_{RECO} [cm]", 30, 0., 1.55, 0., 20., "S");
647  "BtlCluHitsVsEta", "BTL hits per cluster vs #eta; |#eta_{RECO}|;Cluster size", 30, 0., 1.55, 0., 10., "S");
648 
649  meCluTimeRes_ = ibook.book1D("BtlCluTimeRes", "BTL cluster time resolution;T_{RECO}-T_{SIM} [ns]", 100, -0.5, 0.5);
651  ibook.book1D("BtlCluEnergyRes", "BTL cluster energy resolution;E_{RECO}-E_{SIM} [MeV]", 100, -0.5, 0.5);
652  meCluTResvsE_ = ibook.bookProfile("BtlCluTResvsE",
653  "BTL cluster time resolution vs E;E_{SIM} [MeV];(T_{RECO}-T_{SIM}) [ns]",
654  20,
655  0.,
656  20.,
657  -0.5,
658  0.5,
659  "S");
660  meCluTResvsEta_ = ibook.bookProfile("BtlCluTResvsEta",
661  "BTL cluster time resolution vs #eta;|#eta_{RECO}|;(T_{RECO}-T_{SIM}) [ns]",
662  30,
663  0,
664  1.55,
665  -0.5,
666  0.5,
667  "S");
668  meCluTPullvsE_ = ibook.bookProfile("BtlCluTPullvsE",
669  "BTL cluster time pull vs E;E_{SIM} [MeV];(T_{RECO}-T_{SIM})/#sigma_{T_{RECO}}",
670  20,
671  0.,
672  20.,
673  -5.,
674  5.,
675  "S");
677  ibook.bookProfile("BtlCluTPullvsEta",
678  "BTL cluster time pull vs #eta;|#eta_{RECO}|;(T_{RECO}-T_{SIM})/#sigma_{T_{RECO}}",
679  30,
680  0,
681  1.55,
682  -5.,
683  5.,
684  "S");
685  meCluRhoRes_ =
686  ibook.book1D("BtlCluRhoRes", "BTL cluster #rho resolution;#rho_{RECO}-#rho_{SIM} [cm]", 100, -0.5, 0.5);
687  meCluPhiRes_ =
688  ibook.book1D("BtlCluPhiRes", "BTL cluster #phi resolution;#phi_{RECO}-#phi_{SIM} [rad]", 100, -0.03, 0.03);
689  meCluZRes_ = ibook.book1D("BtlCluZRes", "BTL cluster Z resolution;Z_{RECO}-Z_{SIM} [cm]", 100, -0.2, 0.2);
690  if (LocalPosDebug_) {
692  ibook.book1D("BtlCluLocalXRes", "BTL cluster local X resolution;X_{RECO}-X_{SIM} [cm]", 100, -3.1, 3.1);
694  ibook.book1D("BtlCluLocalYRes", "BTL cluster local Y resolution;Y_{RECO}-Y_{SIM} [cm]", 100, -3.1, 3.1);
696  ibook.book1D("BtlCluLocalXPull", "BTL cluster local X pull;X_{RECO}-X_{SIM}/sigmaX_[RECO]", 100, -5., 5.);
698  ibook.book1D("BtlCluLocalYPull", "BTL cluster local Y pull;Y_{RECO}-Y_{SIM}/sigmaY_[RECO]", 100, -5., 5.);
699  meCluZPull_ = ibook.book1D("BtlCluZPull", "BTL cluster Z pull;Z_{RECO}-Z_{SIM}/sigmaZ_[RECO]", 100, -5., 5.);
700  meCluYXLocal_ = ibook.book2D("BtlCluYXLocal",
701  "BTL cluster local Y vs X;X^{local}_{RECO} [cm];Y^{local}_{RECO} [cm]",
702  200,
703  -9.5,
704  9.5,
705  200,
706  -2.8,
707  2.8);
708  meCluYXLocalSim_ = ibook.book2D("BtlCluYXLocalSim",
709  "BTL cluster local Y vs X;X^{local}_{SIM} [cm];Y^{local}_{SIM} [cm]",
710  200,
711  -9.5,
712  9.5,
713  200,
714  -2.8,
715  2.8);
716  meCluXLocalErr_ = ibook.book1D("BtlCluXLocalErr", "BTL cluster X local error;sigmaX_{RECO,loc} [cm]", 30, 0., 3.);
717  meCluYLocalErr_ = ibook.book1D("BtlCluYLocalErr", "BTL cluster Y local error;sigmaY_{RECO,loc} [cm]", 30, 0., 0.9);
718  }
720  ibook.book1D("BtlUnmatchedCluEnergy", "BTL unmatched cluster log10(energy);log10(E_{RECO} [MeV])", 5, -3, 2);
721 
722  // --- UncalibratedRecHits histograms
723 
724  if (uncalibRecHitsPlots_) {
725  for (unsigned int ihistoQ = 0; ihistoQ < nBinsQ_; ++ihistoQ) {
726  std::string hname = Form("TimeResQ_%d", ihistoQ);
727  std::string htitle = Form("BTL time resolution (Q bin = %d);T_{RECO} - T_{SIM} [ns]", ihistoQ);
728  meTimeResQ_[ihistoQ] = ibook.book1D(hname, htitle, 200, -0.3, 0.7);
729 
730  for (unsigned int ihistoEta = 0; ihistoEta < nBinsQEta_; ++ihistoEta) {
731  hname = Form("TimeResQvsEta_%d_%d", ihistoQ, ihistoEta);
732  htitle = Form("BTL time resolution (Q bin = %d, |#eta| bin = %d);T_{RECO} - T_{SIM} [ns]", ihistoQ, ihistoEta);
733  meTimeResQvsEta_[ihistoQ][ihistoEta] = ibook.book1D(hname, htitle, 200, -0.3, 0.7);
734 
735  } // ihistoEta loop
736 
737  } // ihistoQ loop
738 
739  for (unsigned int ihistoEta = 0; ihistoEta < nBinsEta_; ++ihistoEta) {
740  std::string hname = Form("TimeResEta_%d", ihistoEta);
741  std::string htitle = Form("BTL time resolution (|#eta| bin = %d);T_{RECO} - T_{SIM} [ns]", ihistoEta);
742  meTimeResEta_[ihistoEta] = ibook.book1D(hname, htitle, 200, -0.3, 0.7);
743 
744  for (unsigned int ihistoQ = 0; ihistoQ < nBinsEtaQ_; ++ihistoQ) {
745  hname = Form("TimeResEtavsQ_%d_%d", ihistoEta, ihistoQ);
746  htitle = Form("BTL time resolution (|#eta| bin = %d, Q bin = %d);T_{RECO} - T_{SIM} [ns]", ihistoEta, ihistoQ);
747  meTimeResEtavsQ_[ihistoEta][ihistoQ] = ibook.book1D(hname, htitle, 200, -0.3, 0.7);
748 
749  } // ihistoQ loop
750 
751  } // ihistoEta loop
752  }
753 }
754 
755 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
758 
759  desc.add<std::string>("folder", "MTD/BTL/LocalReco");
760  desc.add<edm::InputTag>("recHitsTag", edm::InputTag("mtdRecHits", "FTLBarrel"));
761  desc.add<edm::InputTag>("uncalibRecHitsTag", edm::InputTag("mtdUncalibratedRecHits", "FTLBarrel"));
762  desc.add<edm::InputTag>("simHitsTag", edm::InputTag("mix", "g4SimHitsFastTimerHitsBarrel"));
763  desc.add<edm::InputTag>("recCluTag", edm::InputTag("mtdClusters", "FTLBarrel"));
764  desc.add<edm::InputTag>("trkHitTag", edm::InputTag("mtdTrackingRecHits"));
765  desc.add<double>("HitMinimumEnergy", 1.); // [MeV]
766  desc.add<bool>("LocalPositionDebug", false);
767  desc.add<bool>("UncalibRecHitsPlots", false);
768  desc.add<double>("HitMinimumAmplitude", 30.); // [pC]
769 
770  descriptions.add("btlLocalRecoValid", desc);
771 }
772 
int getMTDTopologyMode() const
Definition: MTDTopology.h:27
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
edm::ESGetToken< MTDGeometry, MTDDigiGeometryRecord > mtdgeoToken_
std::string folder_
float x() const
Definition: FTLCluster.h:120
int size() const
Definition: FTLCluster.h:141
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:32
T z() const
Definition: PV3DBase.h:61
virtual const Topology & topology() const
Definition: GeomDet.cc:67
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
virtual const PixelTopology & specificTopology() const
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
int mtdSide() const
Definition: MTDDetId.h:59
int row(unsigned nrows=16) const
Definition: BTLDetId.h:105
static constexpr float binWidthQ_
MonitorElement * meUnmatchedCluEnergy_
static constexpr int nBinsEtaQ_
LocalPoint pixelToModuleLocalPoint(const LocalPoint &plp, int row, int col) const
edm::EDGetTokenT< MTDTrackingDetSetVector > mtdTrackingHitToken_
MonitorElement * meLongPosPullvsEta_
constexpr NumType convertUnitsTo(double desiredUnits, NumType val)
Definition: GeantUnits.h:73
bool isSameCluster(const FTLCluster &, const FTLCluster &)
void Fill(long long x)
int column(unsigned nrows=16) const
Definition: BTLDetId.h:110
T x() const
Definition: PV3DBase.h:59
T y() const
Definition: PV3DBase.h:60
float y() const
Definition: FTLCluster.h:125
int iEvent
Definition: GenABIO.cc:224
float localX(const float mpX) const override
MonitorElement * bookProfile(TString const &name, TString const &title, int nchX, double lowX, double highX, int, double lowY, double highY, char const *option="s", FUNC onbooking=NOOP())
Definition: DQMStore.h:399
static constexpr float binWidthEta_
T sqrt(T t)
Definition: SSEVec.h:19
edm::EDGetTokenT< FTLClusterCollection > btlRecCluToken_
static constexpr float binsQEta_[nBinsQEta_+1]
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int mtdRR() const
Definition: MTDDetId.h:64
static constexpr int nBinsQEta_
int nrows() const override
edm::ESGetToken< MTDTopology, MTDTopologyRcd > mtdtopoToken_
MonitorElement * meLocalOccupancy_
MonitorElement * meTimeResEta_[nBinsEta_]
float localY(const float mpY) const override
BtlLocalRecoValidation(const edm::ParameterSet &)
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:49
edm::EDGetTokenT< FTLUncalibratedRecHitCollection > btlUncalibRecHitsToken_
Definition: DetId.h:17
const DetId & id() const
Definition: FTLCluster.h:178
constexpr NumType convertMmToCm(NumType millimeters)
Definition: angle_units.h:44
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
edm::EDGetTokenT< FTLRecHitCollection > btlRecHitsToken_
A 2D TrackerRecHit with time and time error information.
MonitorElement * book2D(TString const &name, TString const &title, int nchX, double lowX, double highX, int nchY, double lowY, double highY, FUNC onbooking=NOOP())
Definition: DQMStore.h:212
void add(std::string const &label, ParameterSetDescription const &psetDescription)
caConstants::TupleMultiplicity const CAHitNtupletGeneratorKernelsGPU::HitToTuple const cms::cuda::AtomicPairCounter GPUCACell const *__restrict__ uint32_t const *__restrict__ gpuPixelDoublets::CellNeighborsVector const gpuPixelDoublets::CellTracksVector const GPUCACell::OuterHitOfCell const int32_t nHits
void analyze(const edm::Event &, const edm::EventSetup &) override
static constexpr int nBinsQ_
MonitorElement * meTimeResQvsEta_[nBinsQ_][nBinsQEta_]
static constexpr float binsEtaQ_[nBinsEtaQ_+1]
HLT enums.
MonitorElement * meTimeResEtavsQ_[nBinsEta_][nBinsEtaQ_]
Detector identifier class for the Barrel Timing Layer. The crystal count must start from 0...
Definition: BTLDetId.h:18
ESTransientHandle< T > getTransientHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:162
static constexpr int nBinsEta_
float time() const
Definition: FTLCluster.h:130
BTLDetId::CrysLayout crysLayoutFromTopoMode(const int &topoMode)
Log< level::Warning, false > LogWarning
auto makeValid(const U &iOtherHandleType) noexcept(false)
Definition: ValidHandle.h:52
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
edm::EDGetTokenT< CrossingFrame< PSimHit > > btlSimHitsToken_
MonitorElement * meTimeResQ_[nBinsQ_]
int etaBin(const l1t::HGCalMulticluster *cl)
BTLDetId geographicalId(CrysLayout lay) const
Definition: BTLDetId.cc:171
Definition: Run.h:45