CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
SiStripGainsPCLWorker.cc
Go to the documentation of this file.
4 
5 #include <iostream>
6 #include <sstream>
7 
8 //********************************************************************************//
10  MinTrackMomentum = iConfig.getUntrackedParameter<double>("minTrackMomentum", 3.0);
11  MaxTrackMomentum = iConfig.getUntrackedParameter<double>("maxTrackMomentum", 99999.0);
12  MinTrackEta = iConfig.getUntrackedParameter<double>("minTrackEta", -5.0);
13  MaxTrackEta = iConfig.getUntrackedParameter<double>("maxTrackEta", 5.0);
14  MaxNrStrips = iConfig.getUntrackedParameter<unsigned>("maxNrStrips", 2);
15  MinTrackHits = iConfig.getUntrackedParameter<unsigned>("MinTrackHits", 8);
16  MaxTrackChiOverNdf = iConfig.getUntrackedParameter<double>("MaxTrackChiOverNdf", 3);
17  MaxTrackingIteration = iConfig.getUntrackedParameter<int>("MaxTrackingIteration", 7);
18  AllowSaturation = iConfig.getUntrackedParameter<bool>("AllowSaturation", false);
19  FirstSetOfConstants = iConfig.getUntrackedParameter<bool>("FirstSetOfConstants", true);
20  Validation = iConfig.getUntrackedParameter<bool>("Validation", false);
21  OldGainRemoving = iConfig.getUntrackedParameter<bool>("OldGainRemoving", false);
22  useCalibration = iConfig.getUntrackedParameter<bool>("UseCalibration", false);
23  doChargeMonitorPerPlane = iConfig.getUntrackedParameter<bool>("doChargeMonitorPerPlane", false);
24  m_DQMdir = iConfig.getUntrackedParameter<std::string>("DQMdir", "AlCaReco/SiStripGains");
25  m_calibrationMode = iConfig.getUntrackedParameter<std::string>("calibrationMode", "StdBunch");
26  VChargeHisto = iConfig.getUntrackedParameter<std::vector<std::string>>("ChargeHisto");
27 
28  // fill in the mapping between the histogram indices and the (id,side,plane) tuple
29  std::vector<std::pair<std::string, std::string>> hnames =
30  APVGain::monHnames(VChargeHisto, doChargeMonitorPerPlane, "");
31  for (unsigned int i = 0; i < hnames.size(); i++) {
32  int id = APVGain::subdetectorId((hnames[i]).first);
33  int side = APVGain::subdetectorSide((hnames[i]).first);
34  int plane = APVGain::subdetectorPlane((hnames[i]).first);
35  int thick = APVGain::thickness((hnames[i]).first);
36  std::string s = hnames[i].first;
37 
38  auto loc = APVloc(thick, id, side, plane, s);
39  theTopologyMap.insert(std::make_pair(i, loc));
40  }
41 
42  //Set the monitoring element tag and store
43  dqm_tag_.reserve(7);
44  dqm_tag_.clear();
45  dqm_tag_.push_back("StdBunch"); // statistic collection from Standard Collision Bunch @ 3.8 T
46  dqm_tag_.push_back("StdBunch0T"); // statistic collection from Standard Collision Bunch @ 0 T
47  dqm_tag_.push_back("AagBunch"); // statistic collection from First Collision After Abort Gap @ 3.8 T
48  dqm_tag_.push_back("AagBunch0T"); // statistic collection from First Collision After Abort Gap @ 0 T
49  dqm_tag_.push_back("IsoMuon"); // statistic collection from Isolated Muon @ 3.8 T
50  dqm_tag_.push_back("IsoMuon0T"); // statistic collection from Isolated Muon @ 0 T
51  dqm_tag_.push_back("Harvest"); // statistic collection: Harvest
52 
53  m_tracks_token = consumes<edm::View<reco::Track>>(iConfig.getParameter<edm::InputTag>("tracks"));
54  m_association_token = consumes<TrajTrackAssociationCollection>(iConfig.getParameter<edm::InputTag>("tracks"));
55 
57  tkGeomTokenBR_ = esConsumes<edm::Transition::BeginRun>();
58  tkGeomToken_ = esConsumes<>();
59  gainToken_ = esConsumes<edm::Transition::BeginRun>();
60  qualityToken_ = esConsumes<edm::Transition::BeginRun>();
61 }
62 
63 //********************************************************************************//
65  edm::EventSetup const& iSetup,
67  using namespace edm;
68  static constexpr float defaultGainTick = 690. / 640.;
69 
70  // fills the APV collections at each begin run
71  const TrackerGeometry* bareTkGeomPtr = &iSetup.getData(tkGeomTokenBR_);
72  checkBookAPVColls(bareTkGeomPtr, histograms);
73 
74  const auto gainHandle = iSetup.getHandle(gainToken_);
75  if (!gainHandle.isValid()) {
76  edm::LogError("SiStripGainPCLWorker") << "gainHandle is not valid\n";
77  exit(0);
78  }
79 
80  const auto& siStripQuality = iSetup.getData(qualityToken_);
81 
82  for (unsigned int a = 0; a < histograms.APVsCollOrdered.size(); a++) {
83  std::shared_ptr<stAPVGain> APV = histograms.APVsCollOrdered[a];
84 
85  if (APV->SubDet == PixelSubdetector::PixelBarrel || APV->SubDet == PixelSubdetector::PixelEndcap)
86  continue;
87 
88  APV->isMasked = siStripQuality.IsApvBad(APV->DetId, APV->APVId);
89 
90  if (gainHandle->getNumberOfTags() != 2) {
91  edm::LogError("SiStripGainPCLWorker") << "NUMBER OF GAIN TAG IS EXPECTED TO BE 2\n";
92  fflush(stdout);
93  exit(0);
94  };
95  float newPreviousGain = gainHandle->getApvGain(APV->APVId, gainHandle->getRange(APV->DetId, 1), 1);
96  if (APV->PreviousGain != 1 and newPreviousGain != APV->PreviousGain)
97  edm::LogWarning("SiStripGainPCLWorker") << "WARNING: ParticleGain in the global tag changed\n";
98  APV->PreviousGain = newPreviousGain;
99 
100  float newPreviousGainTick =
101  APV->isMasked ? defaultGainTick : gainHandle->getApvGain(APV->APVId, gainHandle->getRange(APV->DetId, 0), 0);
102  if (APV->PreviousGainTick != 1 and newPreviousGainTick != APV->PreviousGainTick) {
103  edm::LogWarning("SiStripGainPCLWorker")
104  << "WARNING: TickMarkGain in the global tag changed\n"
105  << std::endl
106  << " APV->SubDet: " << APV->SubDet << " APV->APVId:" << APV->APVId << std::endl
107  << " APV->PreviousGainTick: " << APV->PreviousGainTick << " newPreviousGainTick: " << newPreviousGainTick
108  << std::endl;
109  }
110  APV->PreviousGainTick = newPreviousGainTick;
111  }
112 }
113 
114 namespace {
115  struct HitCluster {
116  uint32_t det;
117  const SiStripCluster* strip;
118  const SiPixelCluster* pixel;
119  HitCluster(uint32_t detId, const SiStripCluster* strip, const SiPixelCluster* pixel)
120  : det(detId), strip(strip), pixel(pixel) {}
121  };
122  std::vector<HitCluster> getClusters(const TrackingRecHit* hit) {
123  const auto simple1d = dynamic_cast<const SiStripRecHit1D*>(hit);
124  const auto simple = dynamic_cast<const SiStripRecHit2D*>(hit);
125  const auto matched = dynamic_cast<const SiStripMatchedRecHit2D*>(hit);
126  const auto pixel = dynamic_cast<const SiPixelRecHit*>(hit);
127  std::vector<HitCluster> clusters;
128  if (matched) {
129  clusters.emplace_back(matched->monoId(), &matched->monoCluster(), nullptr);
130  clusters.emplace_back(matched->stereoId(), &matched->stereoCluster(), nullptr);
131  } else if (simple) {
132  clusters.emplace_back(simple->geographicalId().rawId(), simple->cluster().get(), nullptr);
133  } else if (simple1d) {
134  clusters.emplace_back(simple1d->geographicalId().rawId(), simple1d->cluster().get(), nullptr);
135  } else if (pixel) {
136  clusters.emplace_back(pixel->geographicalId().rawId(), nullptr, pixel->cluster().get());
137  }
138  return clusters;
139  }
140 
141  bool isFarFromBorder(const TrajectoryStateOnSurface& trajState, uint32_t detId, const TrackerGeometry* tGeom) {
142  const auto gdu = tGeom->idToDetUnit(detId);
143  if ((!dynamic_cast<const StripGeomDetUnit*>(gdu)) && (!dynamic_cast<const PixelGeomDetUnit*>(gdu))) {
144  edm::LogWarning("SiStripGainCalibTableProducer")
145  << "DetId " << detId << " does not seem to belong to the tracker";
146  return false;
147  }
148  const auto plane = gdu->surface();
149  const auto trapBounds = dynamic_cast<const TrapezoidalPlaneBounds*>(&plane.bounds());
150  const auto rectBounds = dynamic_cast<const RectangularPlaneBounds*>(&plane.bounds());
151 
152  static constexpr double distFromBorder = 1.0;
153  double halfLength = 0.;
154  if (trapBounds) {
155  halfLength = trapBounds->parameters()[3];
156  } else if (rectBounds) {
157  halfLength = .5 * gdu->surface().bounds().length();
158  } else {
159  return false;
160  }
161 
162  const auto pos = trajState.localPosition();
163  const auto posError = trajState.localError().positionError();
164  if (std::abs(pos.y()) + posError.yy() >= (halfLength - distFromBorder))
165  return false;
166 
167  return true;
168  }
169 } // namespace
170 
171 //********************************************************************************//
172 // ------------ method called for each event ------------
174  edm::EventSetup const& iSetup,
175  APVGain::APVGainHistograms const& histograms) const {
176  using namespace edm;
177 
178  unsigned int eventnumber = iEvent.id().event();
179  unsigned int runnumber = iEvent.id().run();
180 
181  edm::LogInfo("SiStripGainsPCLWorker") << "Processing run " << runnumber << " and event " << eventnumber << std::endl;
182 
183  const TrackerTopology* topo = &iSetup.getData(tTopoToken_);
184  const TrackerGeometry* tGeom = &iSetup.getData(tkGeomToken_);
185 
186  // Event data handles
189  edm::Handle<TrajTrackAssociationCollection> trajTrackAssociations;
190  iEvent.getByToken(m_association_token, trajTrackAssociations);
191 
192  for (const auto& elem : theTopologyMap) {
193  LogDebug("SiStripGainsPCLWorker") << elem.first << " - " << elem.second.m_string << " "
194  << elem.second.m_subdetectorId << " " << elem.second.m_subdetectorSide << " "
195  << elem.second.m_subdetectorPlane << std::endl;
196  }
197 
198  LogDebug("SiStripGainsPCLWorker") << "for mode" << m_calibrationMode << std::endl;
199 
200  int elepos = statCollectionFromMode(m_calibrationMode.c_str());
201 
202  std::size_t nStoredClusters{0};
203  for (const auto& assoc : *trajTrackAssociations) {
204  const auto traj = assoc.key.get();
205  const auto track = assoc.val.get();
206 
207  if ((track->eta() < MinTrackEta) || (track->eta() > MaxTrackEta) || (track->p() < MinTrackMomentum) ||
208  (track->p() > MaxTrackMomentum) || (track->numberOfValidHits() < MinTrackHits) ||
209  ((track->chi2() / track->ndof()) > MaxTrackChiOverNdf) || (track->algo() > MaxTrackingIteration))
210  continue;
211 
212  int iCluster{-1};
213  for (const auto& meas : traj->measurements()) {
214  const auto& trajState = meas.updatedState();
215  if (!trajState.isValid())
216  continue;
217 
218  // there can be 2 (stereo module), 1 (no stereo module), or 0 (no pixel or strip hit) clusters
219  auto clusters = getClusters(meas.recHit()->hit());
220  for (const auto hitCluster : clusters) {
221  ++iCluster;
222  bool saturation = false;
223  bool overlapping = false;
224  unsigned int charge = 0;
225  int firstStrip = 0;
226  unsigned int nStrips = 0;
227  if (hitCluster.strip) {
228  const auto& ampls = hitCluster.strip->amplitudes();
229  firstStrip = hitCluster.strip->firstStrip();
230  nStrips = ampls.size();
231  charge = hitCluster.strip->charge();
232  saturation = std::any_of(ampls.begin(), ampls.end(), [](uint8_t amp) { return amp >= 254; });
233 
234  overlapping = (((firstStrip % 128) == 0) || ((firstStrip / 128) != ((firstStrip + int(nStrips)) / 128)));
235  } else if (hitCluster.pixel) {
236  const auto& ampls = hitCluster.pixel->pixelADC();
237  const int firstRow = hitCluster.pixel->minPixelRow();
238  const int firstCol = hitCluster.pixel->minPixelCol();
239  firstStrip = ((firstRow / 80) << 3 | (firstCol / 52)) * 128; //Hack to save the APVId
240  nStrips = 0;
241  for (const auto amp : ampls) {
242  charge += amp;
243  if (amp >= 254)
244  saturation = true;
245  }
246  }
247  // works for both strip and pixel thanks to firstStrip encoding for pixel above, as in the calibTree
248  std::shared_ptr<stAPVGain> APV = histograms.APVsColl.at((hitCluster.det << 4) | (firstStrip / 128));
249 
250  const auto farFromEdge = (hitCluster.strip ? isFarFromBorder(trajState, hitCluster.det, tGeom) : true);
251  if ((APV->SubDet > 2) &&
252  ((!farFromEdge) || overlapping || (saturation && !AllowSaturation) || (nStrips > MaxNrStrips)))
253  continue;
254 
255  int clusterCharge = 0;
256  if (APV->SubDet > 2) { // strip
258  saturation = false;
259  for (const auto origCharge : hitCluster.strip->amplitudes()) {
260  int stripCharge;
261  if (useCalibration) {
262  if (FirstSetOfConstants) {
263  stripCharge = int(origCharge / APV->CalibGain);
264  } else {
265  stripCharge = int(origCharge * (APV->PreviousGain / APV->CalibGain));
266  }
267  } else {
268  if (FirstSetOfConstants) {
269  stripCharge = origCharge;
270  } else {
271  stripCharge = int(origCharge * APV->PreviousGain);
272  }
273  }
274  if (stripCharge > 1024) {
275  stripCharge = 255;
276  saturation = true;
277  } else if (stripCharge > 254) {
278  stripCharge = 254;
279  saturation = true;
280  }
281  clusterCharge += stripCharge;
282  }
283  if (saturation && !AllowSaturation)
284  continue;
285  } else {
286  clusterCharge = charge;
287  }
288  } else { // pixel
289  clusterCharge = charge / 265.0; //expected scale factor between pixel and strip charge
290  }
291 
292  const auto trackDir = trajState.localDirection();
293  const auto path = (10. * APV->Thickness) / std::abs(trackDir.z() / trackDir.mag());
294  double ClusterChargeOverPath = ((double)clusterCharge) / path;
295  if (APV->SubDet > 2) {
296  if (Validation) {
297  ClusterChargeOverPath /= APV->PreviousGain;
298  }
299  if (OldGainRemoving) {
300  ClusterChargeOverPath *= APV->PreviousGain;
301  }
302  } else {
303  // keep processing of pixel cluster charge until here
304  continue;
305  }
306  ++nStoredClusters;
307 
308  // real histogram for calibration
309  histograms.Charge_Vs_Index[elepos]->Fill(APV->Index, ClusterChargeOverPath);
310  LogDebug("SiStripGainsPCLWorker")
311  << " for mode " << m_calibrationMode << "\n"
312  << " i " << iCluster << " useCalibration " << useCalibration << " FirstSetOfConstants "
313  << FirstSetOfConstants << " APV->PreviousGain " << APV->PreviousGain << " APV->CalibGain " << APV->CalibGain
314  << " APV->DetId " << APV->DetId << " APV->Index " << APV->Index << " Charge " << clusterCharge << " Path "
315  << path << " ClusterChargeOverPath " << ClusterChargeOverPath << std::endl;
316 
317  // Fill monitoring histograms
318  int mCharge1 = 0;
319  for (const auto sCharge : hitCluster.strip->amplitudes()) {
320  if (sCharge > 254) {
321  mCharge1 += 254;
322  } else {
323  mCharge1 += sCharge;
324  }
325  }
326  // Revome gains for monitoring
327  int mCharge2 = mCharge1 * APV->PreviousGain; // remove G2
328  int mCharge3 = mCharge1 * APV->PreviousGainTick; // remove G1
329  int mCharge4 = mCharge1 * APV->PreviousGain * APV->PreviousGainTick; // remove G1 and G2
330 
331  LogDebug("SiStripGainsPCLWorker") << " full charge " << mCharge1 << " remove G2 " << mCharge2 << " remove G1 "
332  << mCharge3 << " remove G1*G2 " << mCharge4 << std::endl;
333 
334  auto indices = APVGain::FetchIndices(theTopologyMap, hitCluster.det, topo);
335 
336  for (auto m : indices)
337  histograms.Charge_1[elepos][m]->Fill(((double)mCharge1) / path);
338  for (auto m : indices)
339  histograms.Charge_2[elepos][m]->Fill(((double)mCharge2) / path);
340  for (auto m : indices)
341  histograms.Charge_3[elepos][m]->Fill(((double)mCharge3) / path);
342  for (auto m : indices)
343  histograms.Charge_4[elepos][m]->Fill(((double)mCharge4) / path);
344 
345  if (APV->SubDet == StripSubdetector::TIB) {
346  histograms.Charge_Vs_PathlengthTIB[elepos]->Fill(path, clusterCharge); // TIB
347 
348  } else if (APV->SubDet == StripSubdetector::TOB) {
349  histograms.Charge_Vs_PathlengthTOB[elepos]->Fill(path, clusterCharge); // TOB
350 
351  } else if (APV->SubDet == StripSubdetector::TID) {
352  if (APV->Eta < 0) {
353  histograms.Charge_Vs_PathlengthTIDM[elepos]->Fill(path, clusterCharge);
354  } // TID minus
355  else if (APV->Eta > 0) {
356  histograms.Charge_Vs_PathlengthTIDP[elepos]->Fill(path, clusterCharge);
357  } // TID plus
358 
359  } else if (APV->SubDet == StripSubdetector::TEC) {
360  if (APV->Eta < 0) {
361  if (APV->Thickness < 0.04) {
362  histograms.Charge_Vs_PathlengthTECM1[elepos]->Fill(path, clusterCharge);
363  } // TEC minus, type 1
364  else if (APV->Thickness > 0.04) {
365  histograms.Charge_Vs_PathlengthTECM2[elepos]->Fill(path, clusterCharge);
366  } // TEC minus, type 2
367  } else if (APV->Eta > 0) {
368  if (APV->Thickness < 0.04) {
369  histograms.Charge_Vs_PathlengthTECP1[elepos]->Fill(path, clusterCharge);
370  } // TEC plus, type 1
371  else if (APV->Thickness > 0.04) {
372  histograms.Charge_Vs_PathlengthTECP2[elepos]->Fill(path, clusterCharge);
373  } // TEC plus, type 2
374  }
375  }
376  }
377  }
378  }
379 
380  histograms.EventStats->Fill(0., 0., 1);
381  histograms.EventStats->Fill(1., 0., tracks->size());
382  histograms.EventStats->Fill(2., 0., nStoredClusters);
383 
384  //LogDebug("SiStripGainsPCLWorker")<<" for mode"<< m_calibrationMode
385  // <<" entries in histogram:"<< histograms.Charge_Vs_Index[elepos].getEntries()
386  // <<std::endl;
387 }
388 
389 //********************************************************************************//
390 // ------------ method called once each job just before starting event loop ------------
393  if (bareTkGeomPtr) { // pointer not yet set: called the first time => fill the APVColls
394  auto const& Det = bareTkGeomPtr->dets();
395 
396  edm::LogInfo("SiStripGainsPCLWorker") << " Resetting APV struct" << std::endl;
397 
398  unsigned int Index = 0;
399 
400  for (unsigned int i = 0; i < Det.size(); i++) {
401  DetId Detid = Det[i]->geographicalId();
402  int SubDet = Detid.subdetId();
403 
404  if (SubDet == StripSubdetector::TIB || SubDet == StripSubdetector::TID || SubDet == StripSubdetector::TOB ||
405  SubDet == StripSubdetector::TEC) {
406  auto DetUnit = dynamic_cast<const StripGeomDetUnit*>(Det[i]);
407  if (!DetUnit)
408  continue;
409 
410  const StripTopology& Topo = DetUnit->specificTopology();
411  unsigned int NAPV = Topo.nstrips() / 128;
412 
413  for (unsigned int j = 0; j < NAPV; j++) {
414  auto APV = std::make_shared<stAPVGain>();
415  APV->Index = Index;
416  APV->Bin = -1;
417  APV->DetId = Detid.rawId();
418  APV->APVId = j;
419  APV->SubDet = SubDet;
420  APV->FitMPV = -1;
421  APV->FitMPVErr = -1;
422  APV->FitWidth = -1;
423  APV->FitWidthErr = -1;
424  APV->FitChi2 = -1;
425  APV->FitNorm = -1;
426  APV->Gain = -1;
427  APV->PreviousGain = 1;
428  APV->PreviousGainTick = 1;
429  APV->x = DetUnit->position().basicVector().x();
430  APV->y = DetUnit->position().basicVector().y();
431  APV->z = DetUnit->position().basicVector().z();
432  APV->Eta = DetUnit->position().basicVector().eta();
433  APV->Phi = DetUnit->position().basicVector().phi();
434  APV->R = DetUnit->position().basicVector().transverse();
435  APV->Thickness = DetUnit->surface().bounds().thickness();
436  APV->NEntries = 0;
437  APV->isMasked = false;
438 
439  histograms.APVsCollOrdered.push_back(APV);
440  histograms.APVsColl[(APV->DetId << 4) | APV->APVId] = APV;
441  Index++;
442  histograms.NStripAPVs++;
443  } // loop on APVs
444  } // if is Strips
445  } // loop on dets
446 
447  for (unsigned int i = 0; i < Det.size();
448  i++) { //Make two loop such that the Pixel information is added at the end --> make transition simpler
449  DetId Detid = Det[i]->geographicalId();
450  int SubDet = Detid.subdetId();
452  auto DetUnit = dynamic_cast<const PixelGeomDetUnit*>(Det[i]);
453  if (!DetUnit)
454  continue;
455 
456  const PixelTopology& Topo = DetUnit->specificTopology();
457  unsigned int NROCRow = Topo.nrows() / (80.);
458  unsigned int NROCCol = Topo.ncolumns() / (52.);
459 
460  for (unsigned int j = 0; j < NROCRow; j++) {
461  for (unsigned int i = 0; i < NROCCol; i++) {
462  auto APV = std::make_shared<stAPVGain>();
463  APV->Index = Index;
464  APV->Bin = -1;
465  APV->DetId = Detid.rawId();
466  APV->APVId = (j << 3 | i);
467  APV->SubDet = SubDet;
468  APV->FitMPV = -1;
469  APV->FitMPVErr = -1;
470  APV->FitWidth = -1;
471  APV->FitWidthErr = -1;
472  APV->FitChi2 = -1;
473  APV->Gain = -1;
474  APV->PreviousGain = 1;
475  APV->PreviousGainTick = 1;
476  APV->x = DetUnit->position().basicVector().x();
477  APV->y = DetUnit->position().basicVector().y();
478  APV->z = DetUnit->position().basicVector().z();
479  APV->Eta = DetUnit->position().basicVector().eta();
480  APV->Phi = DetUnit->position().basicVector().phi();
481  APV->R = DetUnit->position().basicVector().transverse();
482  APV->Thickness = DetUnit->surface().bounds().thickness();
483  APV->isMasked = false; //SiPixelQuality_->IsModuleBad(Detid.rawId());
484  APV->NEntries = 0;
485 
486  histograms.APVsCollOrdered.push_back(APV);
487  histograms.APVsColl[(APV->DetId << 4) | APV->APVId] = APV;
488  Index++;
489  histograms.NPixelDets++;
490 
491  } // loop on ROC cols
492  } // loop on ROC rows
493  } // if Pixel
494  } // loop on Dets
495  } //if (!bareTkGeomPtr_) ...
496 }
497 
498 //********************************************************************************//
501  desc.setUnknown();
502  descriptions.addDefault(desc);
503 }
504 
505 //********************************************************************************//
507  edm::Run const& run,
508  edm::EventSetup const& setup,
510  ibooker.cd();
511  std::string dqm_dir = m_DQMdir;
512  const char* tag = dqm_tag_[statCollectionFromMode(m_calibrationMode.c_str())].c_str();
513 
514  edm::LogInfo("SiStripGainsPCLWorker") << "Setting " << dqm_dir << " in DQM and booking histograms for tag " << tag
515  << std::endl;
516 
517  ibooker.setCurrentFolder(dqm_dir);
518 
519  // this MonitorElement is created to log the number of events / tracks and clusters used
520  // by the calibration algorithm
521 
522  histograms.EventStats = ibooker.book2S("EventStats", "Statistics", 3, -0.5, 2.5, 1, 0, 1);
523  histograms.EventStats->setBinLabel(1, "events count", 1);
524  histograms.EventStats->setBinLabel(2, "tracks count", 1);
525  histograms.EventStats->setBinLabel(3, "clusters count", 1);
526 
527  std::string stag(tag);
528  if (!stag.empty() && stag[0] != '_')
529  stag.insert(0, 1, '_');
530 
531  std::string cvi = std::string("Charge_Vs_Index") + stag;
532  std::string cvpTIB = std::string("Charge_Vs_PathlengthTIB") + stag;
533  std::string cvpTOB = std::string("Charge_Vs_PathlengthTOB") + stag;
534  std::string cvpTIDP = std::string("Charge_Vs_PathlengthTIDP") + stag;
535  std::string cvpTIDM = std::string("Charge_Vs_PathlengthTIDM") + stag;
536  std::string cvpTECP1 = std::string("Charge_Vs_PathlengthTECP1") + stag;
537  std::string cvpTECP2 = std::string("Charge_Vs_PathlengthTECP2") + stag;
538  std::string cvpTECM1 = std::string("Charge_Vs_PathlengthTECM1") + stag;
539  std::string cvpTECM2 = std::string("Charge_Vs_PathlengthTECM2") + stag;
540 
541  int elepos = statCollectionFromMode(tag);
542 
543  histograms.Charge_Vs_Index.reserve(dqm_tag_.size());
544  histograms.Charge_Vs_PathlengthTIB.reserve(dqm_tag_.size());
545  histograms.Charge_Vs_PathlengthTOB.reserve(dqm_tag_.size());
546  histograms.Charge_Vs_PathlengthTIDP.reserve(dqm_tag_.size());
547  histograms.Charge_Vs_PathlengthTIDM.reserve(dqm_tag_.size());
548  histograms.Charge_Vs_PathlengthTECP1.reserve(dqm_tag_.size());
549  histograms.Charge_Vs_PathlengthTECP2.reserve(dqm_tag_.size());
550  histograms.Charge_Vs_PathlengthTECM1.reserve(dqm_tag_.size());
551  histograms.Charge_Vs_PathlengthTECM2.reserve(dqm_tag_.size());
552 
553  // The cluster charge is stored by exploiting a non uniform binning in order
554  // reduce the histogram memory size. The bin width is relaxed with a falling
555  // exponential function and the bin boundaries are stored in the binYarray.
556  // The binXarray is used to provide as many bins as the APVs.
557  //
558  // More details about this implementations are here:
559  // https://indico.cern.ch/event/649344/contributions/2672267/attachments/1498323/2332518/OptimizeChHisto.pdf
560 
561  std::vector<float> binXarray;
562  binXarray.reserve(histograms.NStripAPVs + 1);
563  for (unsigned int a = 0; a <= histograms.NStripAPVs; a++) {
564  binXarray.push_back((float)a);
565  }
566 
567  std::array<float, 688> binYarray;
568  double p0 = 5.445;
569  double p1 = 0.002113;
570  double p2 = 69.01576;
571  double y = 0.;
572  for (int b = 0; b < 687; b++) {
573  binYarray[b] = y;
574  if (y <= 902.)
575  y = y + 2.;
576  else
577  y = (p0 - log(exp(p0 - p1 * y) - p2 * p1)) / p1;
578  }
579  binYarray[687] = 4000.;
580 
581  histograms.Charge_1[elepos].clear();
582  histograms.Charge_2[elepos].clear();
583  histograms.Charge_3[elepos].clear();
584  histograms.Charge_4[elepos].clear();
585 
586  auto it = histograms.Charge_Vs_Index.begin();
587  histograms.Charge_Vs_Index.insert(
588  it + elepos,
589  ibooker.book2S(cvi.c_str(), cvi.c_str(), histograms.NStripAPVs, &binXarray[0], 687, binYarray.data()));
590 
591  it = histograms.Charge_Vs_PathlengthTIB.begin();
592  histograms.Charge_Vs_PathlengthTIB.insert(it + elepos,
593  ibooker.book2S(cvpTIB.c_str(), cvpTIB.c_str(), 20, 0.3, 1.3, 250, 0, 2000));
594 
595  it = histograms.Charge_Vs_PathlengthTOB.begin();
596  histograms.Charge_Vs_PathlengthTOB.insert(it + elepos,
597  ibooker.book2S(cvpTOB.c_str(), cvpTOB.c_str(), 20, 0.3, 1.3, 250, 0, 2000));
598 
599  it = histograms.Charge_Vs_PathlengthTIDP.begin();
600  histograms.Charge_Vs_PathlengthTIDP.insert(
601  it + elepos, ibooker.book2S(cvpTIDP.c_str(), cvpTIDP.c_str(), 20, 0.3, 1.3, 250, 0, 2000));
602 
603  it = histograms.Charge_Vs_PathlengthTIDM.begin();
604  histograms.Charge_Vs_PathlengthTIDM.insert(
605  it + elepos, ibooker.book2S(cvpTIDM.c_str(), cvpTIDM.c_str(), 20, 0.3, 1.3, 250, 0, 2000));
606 
607  it = histograms.Charge_Vs_PathlengthTECP1.begin();
608  histograms.Charge_Vs_PathlengthTECP1.insert(
609  it + elepos, ibooker.book2S(cvpTECP1.c_str(), cvpTECP1.c_str(), 20, 0.3, 1.3, 250, 0, 2000));
610 
611  it = histograms.Charge_Vs_PathlengthTECP2.begin();
612  histograms.Charge_Vs_PathlengthTECP2.insert(
613  it + elepos, ibooker.book2S(cvpTECP2.c_str(), cvpTECP2.c_str(), 20, 0.3, 1.3, 250, 0, 2000));
614 
615  it = histograms.Charge_Vs_PathlengthTECM1.begin();
616  histograms.Charge_Vs_PathlengthTECM1.insert(
617  it + elepos, ibooker.book2S(cvpTECM1.c_str(), cvpTECM1.c_str(), 20, 0.3, 1.3, 250, 0, 2000));
618 
619  it = histograms.Charge_Vs_PathlengthTECM2.begin();
620  histograms.Charge_Vs_PathlengthTECM2.insert(
621  it + elepos, ibooker.book2S(cvpTECM2.c_str(), cvpTECM2.c_str(), 20, 0.3, 1.3, 250, 0, 2000));
622 
623  std::vector<std::pair<std::string, std::string>> hnames =
625  for (unsigned int i = 0; i < hnames.size(); i++) {
626  std::string htag = (hnames[i]).first + stag;
627  histograms.Charge_1[elepos].push_back(ibooker.book1DD(htag.c_str(), (hnames[i]).second.c_str(), 100, 0., 1000.));
628  }
629 
631  for (unsigned int i = 0; i < hnames.size(); i++) {
632  std::string htag = (hnames[i]).first + stag;
633  histograms.Charge_2[elepos].push_back(ibooker.book1DD(htag.c_str(), (hnames[i]).second.c_str(), 100, 0., 1000.));
634  }
635 
637  for (unsigned int i = 0; i < hnames.size(); i++) {
638  std::string htag = (hnames[i]).first + stag;
639  histograms.Charge_3[elepos].push_back(ibooker.book1DD(htag.c_str(), (hnames[i]).second.c_str(), 100, 0., 1000.));
640  }
641 
643  for (unsigned int i = 0; i < hnames.size(); i++) {
644  std::string htag = (hnames[i]).first + stag;
645  histograms.Charge_4[elepos].push_back(ibooker.book1DD(htag.c_str(), (hnames[i]).second.c_str(), 100, 0., 1000.));
646  }
647 }
RunNumber_t run() const
Definition: EventID.h:38
EventNumber_t event() const
Definition: EventID.h:40
static constexpr auto TEC
virtual int nstrips() const =0
T getUntrackedParameter(std::string const &, T const &) const
std::array< std::vector< dqm::reco::MonitorElement * >, 7 > Charge_4
std::vector< unsigned int > FetchIndices(std::map< unsigned int, APVloc >, uint32_t, const TrackerTopology *topo=nullptr)
static std::vector< std::string > checklist log
const TrackerGeomDet * idToDetUnit(DetId) const override
Return the pointer to the GeomDetUnit corresponding to a given DetId.
MonitorElement * book2S(TString const &name, TString const &title, int nchX, double lowX, double highX, int nchY, double lowY, double highY, FUNC onbooking=NOOP())
Definition: DQMStore.h:219
edm::ESGetToken< SiStripQuality, SiStripQualityRcd > qualityToken_
tuple SubDet
Definition: ntupleEnum.py:14
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:32
virtual int ncolumns() const =0
const TString p2
Definition: fwPaths.cc:13
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
std::vector< dqm::reco::MonitorElement * > Charge_Vs_Index
LocalVector localDirection() const
void dqmBeginRun(edm::Run const &, edm::EventSetup const &, APVGain::APVGainHistograms &) const override
SiStripGainsPCLWorker(const edm::ParameterSet &)
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
std::vector< dqm::reco::MonitorElement * > Charge_Vs_PathlengthTOB
virtual int nrows() const =0
int subdetectorPlane(uint32_t, const TrackerTopology *)
std::vector< dqm::reco::MonitorElement * > Charge_Vs_PathlengthTIDP
auto const & tracks
cannot be loose
const DetContainer & dets() const override
Returm a vector of all GeomDet (including all GeomDetUnits)
Exp< T >::type exp(const T &t)
Definition: Exp.h:22
Log< level::Error, false > LogError
edm::EDGetTokenT< edm::View< reco::Track > > m_tracks_token
float length() const override
Lenght along local Y.
tuple nStrips
1.2 is to make the matching window safely the two nearest strips 0.35 is the size of an ME0 chamber i...
LocalError positionError() const
edm::EDGetTokenT< TrajTrackAssociationCollection > m_association_token
std::vector< dqm::reco::MonitorElement * > Charge_Vs_PathlengthTIDM
std::vector< std::string > dqm_tag_
void Fill(long long x)
bool getData(T &iHolder) const
Definition: EventSetup.h:122
U second(std::pair< T, U > const &p)
void dqmAnalyze(edm::Event const &, edm::EventSetup const &, APVGain::APVGainHistograms const &) const override
MonitorElement * book1DD(TString const &name, TString const &title, int nchX, double lowX, double highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:155
std::vector< dqm::reco::MonitorElement * > Charge_Vs_PathlengthTECP2
int iEvent
Definition: GenABIO.cc:224
void addDefault(ParameterSetDescription const &psetDescription)
std::atomic< unsigned int > NPixelDets
edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > tkGeomTokenBR_
edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > tkGeomToken_
std::array< std::vector< dqm::reco::MonitorElement * >, 7 > Charge_3
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:48
std::atomic< unsigned int > NStripAPVs
const TString p1
Definition: fwPaths.cc:12
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::map< unsigned int, APVloc > theTopologyMap
const LocalTrajectoryError & localError() const
static constexpr auto TOB
virtual void setBinLabel(int bin, const std::string &label, int axis=1)
set bin label for x, y or z axis (axis=1, 2, 3 respectively)
Log< level::Info, false > LogInfo
edm::ESGetToken< SiStripGain, SiStripGainRcd > gainToken_
Definition: DetId.h:17
static constexpr auto TIB
void checkBookAPVColls(const TrackerGeometry *bareTkGeomPtr, APVGain::APVGainHistograms &histograms) const
std::vector< std::shared_ptr< stAPVGain > > APVsCollOrdered
int subdetectorSide(uint32_t, const TrackerTopology *)
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
double b
Definition: hdecay.h:118
std::vector< std::string > VChargeHisto
int thickness(uint32_t)
std::unordered_map< unsigned int, std::shared_ptr< stAPVGain > > APVsColl
edm::EventID id() const
Definition: EventBase.h:59
Pixel cluster – collection of neighboring pixels above threshold.
dqm::reco::MonitorElement * EventStats
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &, APVGain::APVGainHistograms &) const override
double a
Definition: hdecay.h:119
std::array< std::vector< dqm::reco::MonitorElement * >, 7 > Charge_1
std::vector< dqm::reco::MonitorElement * > Charge_Vs_PathlengthTECM2
int statCollectionFromMode(const char *tag) const
edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > tTopoToken_
std::array< std::vector< dqm::reco::MonitorElement * >, 7 > Charge_2
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:151
Log< level::Warning, false > LogWarning
std::vector< dqm::reco::MonitorElement * > Charge_Vs_PathlengthTECM1
list indices
Definition: dqmdumpme.py:50
std::vector< std::pair< std::string, std::string > > monHnames(std::vector< std::string >, bool, const char *tag)
static constexpr auto TID
std::vector< dqm::reco::MonitorElement * > Charge_Vs_PathlengthTIB
int subdetectorId(uint32_t)
std::vector< dqm::reco::MonitorElement * > Charge_Vs_PathlengthTECP1
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: Run.h:45
Our base class.
Definition: SiPixelRecHit.h:23
#define LogDebug(id)