CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
SiPixelLorentzAnglePCLHarvester.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: CalibTracker/SiPixelLorentzAnglePCLHarvester
4 // Class: SiPixelLorentzAnglePCLHarvester
5 //
12 //
13 // Original Author: mmusich
14 // Created: Sat, 29 May 2021 14:46:19 GMT
15 //
16 //
17 
18 // system includes
19 #include <fmt/format.h>
20 #include <fmt/printf.h>
21 #include <fstream>
22 
23 // user includes
41 
42 /*
43  * Auxilliary struct to store fit results
44  */
45 namespace SiPixelLAHarvest {
46  struct fitResults {
47  public:
49  // set all parameters to default
50  p0 = p1 = p2 = p3 = p4 = p5 = 0.;
51  e0 = e1 = e2 = e3 = e4 = e5 = 0.;
52  chi2 = prob = redChi2 = tan_LA = error_LA = -9999.;
53  ndf = -999;
54  };
55 
56  double p0;
57  double e0;
58  double p1;
59  double e1;
60  double p2;
61  double e2;
62  double p3;
63  double e3;
64  double p4;
65  double e4;
66  double p5;
67  double e5;
68  double chi2;
69  int ndf;
70  double prob;
71  double redChi2;
72  double tan_LA;
73  double error_LA;
74  };
75 } // namespace SiPixelLAHarvest
76 
77 //------------------------------------------------------------------------------
79 public:
81  ~SiPixelLorentzAnglePCLHarvester() override = default;
82  void beginRun(const edm::Run&, const edm::EventSetup&) override;
83 
85 
86 private:
88  void findMean(MonitorElement* h_drift_depth_adc_slice_, int i, int i_ring);
89  SiPixelLAHarvest::fitResults fitAndStore(std::shared_ptr<SiPixelLorentzAngle> theLA, int i_idx, int i_lay, int i_mod);
90 
91  // es tokens
96 
97  std::vector<std::string> newmodulelist_;
99  const double fitChi2Cut_;
100  const int minHitsCut_;
102  std::unique_ptr<TF1> f1;
103  float width_;
104 
108 };
109 
110 //------------------------------------------------------------------------------
112  : geomEsToken_(esConsumes<edm::Transition::BeginRun>()),
113  topoEsToken_(esConsumes<edm::Transition::BeginRun>()),
114  siPixelLAEsToken_(esConsumes<edm::Transition::BeginRun>()),
115  magneticFieldToken_(esConsumes<edm::Transition::BeginRun>()),
116  newmodulelist_(iConfig.getParameter<std::vector<std::string>>("newmodulelist")),
117  dqmDir_(iConfig.getParameter<std::string>("dqmDir")),
118  fitChi2Cut_(iConfig.getParameter<double>("fitChi2Cut")),
119  minHitsCut_(iConfig.getParameter<int>("minHitsCut")),
120  recordName_(iConfig.getParameter<std::string>("record")) {
121  // first ensure DB output service is available
123  if (!poolDbService.isAvailable())
124  throw cms::Exception("SiPixelLorentzAnglePCLHarvester") << "PoolDBService required";
125 }
126 
127 //------------------------------------------------------------------------------
129  // geometry
130  const TrackerGeometry* geom = &iSetup.getData(geomEsToken_);
131  const TrackerTopology* tTopo = &iSetup.getData(topoEsToken_);
132 
135 
136  PixelTopologyMap map = PixelTopologyMap(geom, tTopo);
138  hists.nModules_.resize(hists.nlay);
139  for (int i = 0; i < hists.nlay; i++) {
140  hists.nModules_[i] = map.getPXBModules(i + 1);
141  }
142 
143  // list of modules already filled, then return (we already entered here)
144  if (!hists.BPixnewDetIds_.empty() || !hists.FPixnewDetIds_.empty())
145  return;
146 
147  if (!newmodulelist_.empty()) {
148  for (auto const& modulename : newmodulelist_) {
149  if (modulename.find("BPix_") != std::string::npos) {
150  PixelBarrelName bn(modulename, true);
151  const auto& detId = bn.getDetId(tTopo);
152  hists.BPixnewmodulename_.push_back(modulename);
153  hists.BPixnewDetIds_.push_back(detId.rawId());
154  hists.BPixnewModule_.push_back(bn.moduleName());
155  hists.BPixnewLayer_.push_back(bn.layerName());
156  } else if (modulename.find("FPix_") != std::string::npos) {
157  PixelEndcapName en(modulename, true);
158  const auto& detId = en.getDetId(tTopo);
159  hists.FPixnewmodulename_.push_back(modulename);
160  hists.FPixnewDetIds_.push_back(detId.rawId());
161  hists.FPixnewDisk_.push_back(en.diskName());
162  hists.FPixnewBlade_.push_back(en.bladeName());
163  }
164  }
165  }
166 
167  uint count = 0;
168  for (const auto& id : hists.BPixnewDetIds_) {
169  LogDebug("SiPixelLorentzAnglePCLHarvester") << id;
170  count++;
171  }
172  LogDebug("SiPixelLorentzAnglePCLHarvester") << "Stored a total of " << count << " new detIds.";
173 
174  // list of modules already filled, return (we already entered here)
175  if (!hists.detIdsList.empty())
176  return;
177 
178  std::vector<uint32_t> treatedIndices;
179 
180  for (auto det : geom->detsPXB()) {
181  const PixelGeomDetUnit* pixelDet = dynamic_cast<const PixelGeomDetUnit*>(det);
182  width_ = pixelDet->surface().bounds().thickness();
183  const auto& layer = tTopo->pxbLayer(pixelDet->geographicalId());
184  const auto& module = tTopo->pxbModule(pixelDet->geographicalId());
185  int i_index = module + (layer - 1) * hists.nModules_[layer - 1];
186 
187  uint32_t rawId = pixelDet->geographicalId().rawId();
188 
189  // if the detId is already accounted for in the special class, do not attach it
190  if (std::find(hists.BPixnewDetIds_.begin(), hists.BPixnewDetIds_.end(), rawId) != hists.BPixnewDetIds_.end())
191  continue;
192 
193  if (std::find(treatedIndices.begin(), treatedIndices.end(), i_index) != treatedIndices.end()) {
194  hists.detIdsList.at(i_index).push_back(rawId);
195  } else {
196  hists.detIdsList.insert(std::pair<uint32_t, std::vector<uint32_t>>(i_index, {rawId}));
197  treatedIndices.push_back(i_index);
198  }
199  }
200 
201  count = 0;
202  for (const auto& i : treatedIndices) {
203  for (const auto& id : hists.detIdsList.at(i)) {
204  LogDebug("SiPixelLorentzAnglePCLHarvester") << id;
205  count++;
206  };
207  }
208  LogDebug("SiPixelLorentzAnglePCLHarvester") << "Stored a total of " << count << " detIds.";
209 }
210 
211 //------------------------------------------------------------------------------
213  // go in the right directory
214  iGetter.cd();
215  iGetter.setCurrentFolder(dqmDir_);
216 
217  // fetch the 2D histograms
218  for (int i_layer = 1; i_layer <= hists.nlay; i_layer++) {
219  const auto& prefix_ = fmt::sprintf("%s/BPix/BPixLayer%i", dqmDir_, i_layer);
220  for (int i_module = 1; i_module <= hists.nModules_[i_layer - 1]; i_module++) {
221  int i_index = i_module + (i_layer - 1) * hists.nModules_[i_layer - 1];
222 
223  hists.h_drift_depth_[i_index] =
224  iGetter.get(fmt::format("{}/h_drift_depth_layer{}_module{}", prefix_, i_layer, i_module));
225 
226  if (hists.h_drift_depth_[i_index] == nullptr) {
227  edm::LogError("SiPixelLorentzAnglePCLHarvester::dqmEndJob")
228  << "Failed to retrieve electron drift over depth for layer " << i_layer << ", module " << i_module << ".";
229  continue;
230  }
231 
232  hists.h_drift_depth_adc_[i_index] =
233  iGetter.get(fmt::format("{}/h_drift_depth_adc_layer{}_module{}", prefix_, i_layer, i_module));
234 
235  hists.h_drift_depth_adc2_[i_index] =
236  iGetter.get(fmt::format("{}/h_drift_depth_adc2_layer{}_module{}", prefix_, i_layer, i_module));
237 
238  hists.h_drift_depth_noadc_[i_index] =
239  iGetter.get(fmt::format("{}/h_drift_depth_noadc_layer{}_module{}", prefix_, i_layer, i_module));
240 
241  hists.h_mean_[i_index] = iGetter.get(fmt::format("{}/h_mean_layer{}_module{}", dqmDir_, i_layer, i_module));
242 
243  hists.h_drift_depth_[i_index]->divide(
244  hists.h_drift_depth_adc_[i_index], hists.h_drift_depth_noadc_[i_index], 1., 1., "");
245  }
246  }
247 
248  // fetch the new modules 2D histograms
249  for (int i = 0; i < (int)hists.BPixnewDetIds_.size(); i++) {
250  int new_index = i + 1 + hists.nModules_[hists.nlay - 1] + (hists.nlay - 1) * hists.nModules_[hists.nlay - 1];
251 
252  hists.h_drift_depth_[new_index] = iGetter.get(
253  fmt::format("{}/h_BPixnew_drift_depth_{}", dqmDir_ + "/BPix/NewModules", hists.BPixnewmodulename_[i]));
254 
255  if (hists.h_drift_depth_[new_index] == nullptr) {
256  edm::LogError("SiPixelLorentzAnglePCLHarvester")
257  << "Failed to retrieve electron drift over depth for new module " << hists.BPixnewmodulename_[i] << ".";
258  continue;
259  }
260 
261  hists.h_drift_depth_adc_[new_index] = iGetter.get(
262  fmt::format("{}/h_BPixnew_drift_depth_adc_{}", dqmDir_ + "/BPix/NewModules", hists.BPixnewmodulename_[i]));
263 
264  hists.h_drift_depth_adc2_[new_index] = iGetter.get(
265  fmt::format("{}/h_BPixnew_drift_depth_adc2_{}", dqmDir_ + "/BPix/NewModules", hists.BPixnewmodulename_[i]));
266 
267  hists.h_drift_depth_noadc_[new_index] = iGetter.get(
268  fmt::format("{}/h_BPixnew_drift_depth_noadc_{}", dqmDir_ + "/BPix/NewModules", hists.BPixnewmodulename_[i]));
269 
270  hists.h_mean_[new_index] = iGetter.get(fmt::format("{}/h_BPixnew_mean_{}", dqmDir_, hists.BPixnewmodulename_[i]));
271 
272  hists.h_drift_depth_[new_index]->divide(
273  hists.h_drift_depth_adc_[new_index], hists.h_drift_depth_noadc_[new_index], 1., 1., "");
274  }
275 
276  hists.h_bySectOccupancy_ = iGetter.get(fmt::format("{}/h_bySectorOccupancy", dqmDir_ + "/SectorMonitoring"));
277  if (hists.h_bySectOccupancy_ == nullptr) {
278  edm::LogError("SiPixelLorentzAnglePCLHarvester") << "Failed to retrieve the hit on track occupancy.";
279  return;
280  }
281 
282  int hist_drift_;
283  int hist_depth_;
284  double min_drift_;
285  double max_drift_;
286 
287  if (hists.h_drift_depth_adc_[1] != nullptr) {
288  hist_drift_ = hists.h_drift_depth_adc_[1]->getNbinsX();
289  hist_depth_ = hists.h_drift_depth_adc_[1]->getNbinsY();
290  min_drift_ = hists.h_drift_depth_adc_[1]->getAxisMin(1);
291  max_drift_ = hists.h_drift_depth_adc_[1]->getAxisMax(1);
292  } else {
293  hist_drift_ = 100;
294  hist_depth_ = 50;
295  min_drift_ = -500.;
296  max_drift_ = 500.;
297  }
298 
299  iBooker.setCurrentFolder(fmt::format("{}Harvesting", dqmDir_));
300  MonitorElement* h_drift_depth_adc_slice_ =
301  iBooker.book1D("h_drift_depth_adc_slice", "slice of adc histogram", hist_drift_, min_drift_, max_drift_);
302 
303  // book histogram of differences
304  MonitorElement* h_diffLA = iBooker.book1D(
305  "h_diffLA", "difference in #mu_{H}; #Delta #mu_{H}/#mu_{H} (old-new)/old [%];n. modules", 100, -3., 3.);
306 
307  // retrieve the number of bins from the other monitoring histogram
308  const auto& maxSect = hists.h_bySectOccupancy_->getNbinsX();
309  const double lo = -0.5;
310  const double hi = maxSect + 0.5;
311 
312  // this will be booked in the Harvesting folder
313  iBooker.setCurrentFolder(fmt::format("{}Harvesting/SectorMonitoring", dqmDir_));
314  std::string repText = "%s tan#theta_{LA}/B by sector;pixel sector;%s tan(#theta_{LA})/B [1/T]";
316  iBooker.book1D("h_LAbySector_Measured", fmt::sprintf(repText, "measured", "measured"), maxSect, lo, hi);
318  iBooker.book1D("h_LAbySector_Accepted", fmt::sprintf(repText, "accepted", "accepted"), maxSect, lo, hi);
320  iBooker.book1D("h_LAbySector_Rejected", fmt::sprintf(repText, "rejected", "rejected"), maxSect, lo, hi);
321  hists.h_bySectLA_ = iBooker.book1D("h_LAbySector", fmt::sprintf(repText, "payload", "payload"), maxSect, lo, hi);
323  iBooker.book1D("h_deltaLAbySector", fmt::sprintf(repText, "#Delta", "#Delta"), maxSect, lo, hi);
325  iBooker.book1D("h_bySectorChi2", "Fit #chi^{2}/ndf by sector;pixel sector; fit #chi^{2}/ndf", maxSect, lo, hi);
326 
327  // copy the bin labels from the occupancy histogram
328  for (int bin = 1; bin < maxSect; bin++) {
329  const auto& binName = hists.h_bySectOccupancy_->getTH1()->GetXaxis()->GetBinLabel(bin);
331  hists.h_bySectSetLA_->setBinLabel(bin, binName);
333  hists.h_bySectLA_->setBinLabel(bin, binName);
335  hists.h_bySectChi2_->setBinLabel(bin, binName);
336  }
337 
338  // clang-format off
339  edm::LogPrint("LorentzAngle") << "module" << "\t" << "layer" << "\t"
340  << "offset" << "\t" << "e0" << "\t"
341  << "slope" << "\t" << "e1" << "\t"
342  << "rel.err" << "\t" << "pull" << "\t"
343  << "p2" << "\t" << "e2" << "\t"
344  << "p3" << "\t" << "e3" << "\t"
345  << "p4" << "\t" << "e4" << "\t"
346  << "p5" << "\t" << "e5" << "\t"
347  << "chi2" << "\t" << "prob" << "\t"
348  << "newDetId" << "\t" << "tan(LA)" << "\t"
349  << "Error(LA)" ;
350  // clang-format on
351 
352  // payload to be written out
353  std::shared_ptr<SiPixelLorentzAngle> LorentzAngle = std::make_shared<SiPixelLorentzAngle>();
354 
355  // fill the map of simulation values
356  double p1_simul_newmodule = 0.294044;
357  double p1_simul[hists.nlay + 1][hists.nModules_[hists.nlay - 1]];
358  for (int i_layer = 1; i_layer <= hists.nlay; i_layer++) {
359  for (int i_module = 1; i_module <= hists.nModules_[i_layer - 1]; i_module++) {
360  if (i_layer == 1)
361  p1_simul[i_layer - 1][i_module - 1] = 0.436848;
362  else if (i_layer == 2)
363  p1_simul[i_layer - 1][i_module - 1] = 0.25802;
364  else if (i_layer == 3 && i_module <= 4)
365  p1_simul[i_layer - 1][i_module - 1] = 0.29374;
366  else if (i_layer == 3 && i_module >= 5)
367  p1_simul[i_layer - 1][i_module - 1] = 0.31084;
368  else if (i_layer == 4 && i_module <= 4)
369  p1_simul[i_layer - 1][i_module - 1] = 0.29944;
370  else
371  p1_simul[i_layer - 1][i_module - 1] = 0.31426;
372  }
373  }
374  // fictitious n-th layer to store the values of new modules
375  for (int i_module = 1; i_module <= hists.nModules_[hists.nlay - 1]; i_module++) {
376  p1_simul[hists.nlay][i_module - 1] = p1_simul_newmodule;
377  }
378 
379  // loop over "new" BPix modules
380  for (int j = 0; j < (int)hists.BPixnewDetIds_.size(); j++) {
381  //uint32_t rawId = hists.BPixnewDetIds_[j];
382  int new_index = j + 1 + hists.nModules_[hists.nlay - 1] + (hists.nlay - 1) * hists.nModules_[hists.nlay - 1];
383  if (hists.h_drift_depth_adc_[new_index] == nullptr)
384  continue;
385  for (int i = 1; i <= hist_depth_; i++) {
386  findMean(h_drift_depth_adc_slice_, i, new_index);
387  }
388 
389  // fit the distributions and store the LA in the payload
390  const auto& res = fitAndStore(LorentzAngle, new_index, hists.BPixnewLayer_[j], hists.BPixnewModule_[j]);
391 
392  edm::LogPrint("SiPixelLorentzAngle") << std::setprecision(4) << hists.BPixnewModule_[j] << "\t"
393  << hists.BPixnewLayer_[j] << "\t" << res.p0 << "\t" << res.e0 << "\t" << res.p1
394  << std::setprecision(3) << "\t" << res.e1 << "\t" << res.e1 / res.p1 * 100.
395  << "\t" << (res.p1 - p1_simul[hists.nlay][0]) / res.e1 << "\t" << res.p2
396  << "\t" << res.e2 << "\t" << res.p3 << "\t" << res.e3 << "\t" << res.p4 << "\t"
397  << res.e4 << "\t" << res.p5 << "\t" << res.e5 << "\t" << res.chi2 << "\t"
398  << res.prob << "\t" << hists.BPixnewDetIds_[j] << "\t" << res.tan_LA << "\t"
399  << res.error_LA;
400  } // loop on BPix new modules
401 
402  //loop over modules and layers to fit the lorentz angle
403  for (int i_layer = 1; i_layer <= hists.nlay; i_layer++) {
404  for (int i_module = 1; i_module <= hists.nModules_[i_layer - 1]; i_module++) {
405  int i_index = i_module + (i_layer - 1) * hists.nModules_[i_layer - 1];
406  if (hists.h_drift_depth_adc_[i_index] == nullptr)
407  continue;
408  //loop over bins in depth (z-local-coordinate) (in order to fit slices)
409  for (int i = 1; i <= hist_depth_; i++) {
410  findMean(h_drift_depth_adc_slice_, i, i_index);
411  } // end loop over bins in depth
412 
413  // fit the distributions and store the LA in the payload
414  const auto& res = fitAndStore(LorentzAngle, i_index, i_layer, i_module);
415 
416  edm::LogPrint("SiPixelLorentzAngle")
417  << std::setprecision(4) << i_module << "\t" << i_layer << "\t" << res.p0 << "\t" << res.e0 << "\t" << res.p1
418  << std::setprecision(3) << "\t" << res.e1 << "\t" << res.e1 / res.p1 * 100. << "\t"
419  << (res.p1 - p1_simul[i_layer - 1][i_module - 1]) / res.e1 << "\t" << res.p2 << "\t" << res.e2 << "\t"
420  << res.p3 << "\t" << res.e3 << "\t" << res.p4 << "\t" << res.e4 << "\t" << res.p5 << "\t" << res.e5 << "\t"
421  << res.chi2 << "\t" << res.prob << "\t"
422  << "null"
423  << "\t" << res.tan_LA << "\t" << res.error_LA;
424  }
425  } // end loop over modules and layers
426 
427  // fill the rest of DetIds not filled above (for the moment FPix)
428  const auto& currentLAMap = currentLorentzAngle->getLorentzAngles();
429  const auto& newLAMap = LorentzAngle->getLorentzAngles();
430  std::vector<unsigned int> currentLADets;
431  std::vector<unsigned int> newLADets;
432 
433  std::transform(currentLAMap.begin(),
434  currentLAMap.end(),
435  std::back_inserter(currentLADets),
436  [](const std::map<unsigned int, float>::value_type& pair) { return pair.first; });
437 
438  std::transform(newLAMap.begin(),
439  newLAMap.end(),
440  std::back_inserter(newLADets),
441  [](const std::map<unsigned int, float>::value_type& pair) { return pair.first; });
442 
443  std::vector<unsigned int> notCommon;
444  std::set_symmetric_difference(
445  currentLADets.begin(), currentLADets.end(), newLADets.begin(), newLADets.end(), std::back_inserter(notCommon));
446 
447  for (const auto& id : notCommon) {
448  float fPixLorentzAnglePerTesla_ = currentLorentzAngle->getLorentzAngle(id);
449  if (!LorentzAngle->putLorentzAngle(id, fPixLorentzAnglePerTesla_)) {
450  edm::LogError("SiPixelLorentzAnglePCLHarvester")
451  << "[SiPixelLorentzAnglePCLHarvester::dqmEndRun] filling rest of payload: detid already exists";
452  }
453  }
454 
455  for (const auto& id : newLADets) {
456  float deltaMuHoverMuH = (currentLorentzAngle->getLorentzAngle(id) - LorentzAngle->getLorentzAngle(id)) /
458  h_diffLA->Fill(deltaMuHoverMuH);
459  }
460 
461  // fill the DB object record
463  if (mydbservice.isAvailable()) {
464  try {
465  mydbservice->writeOneIOV(*LorentzAngle, mydbservice->currentTime(), recordName_);
466  } catch (const cond::Exception& er) {
467  edm::LogError("SiPixelLorentzAngleDB") << er.what();
468  } catch (const std::exception& er) {
469  edm::LogError("SiPixelLorentzAngleDB") << "caught std::exception " << er.what();
470  }
471  } else {
472  edm::LogError("SiPixelLorentzAngleDB") << "Service is unavailable";
473  }
474 }
475 
476 //------------------------------------------------------------------------------
477 void SiPixelLorentzAnglePCLHarvester::findMean(MonitorElement* h_drift_depth_adc_slice_, int i, int i_ring) {
478  double nentries = 0;
479  h_drift_depth_adc_slice_->Reset();
480  int hist_drift_ = h_drift_depth_adc_slice_->getNbinsX();
481 
482  // determine sigma and sigma^2 of the adc counts and average adc counts
483  //loop over bins in drift width
484  for (int j = 1; j <= hist_drift_; j++) {
485  if (hists.h_drift_depth_noadc_[i_ring]->getBinContent(j, i) >= 1) {
486  double adc_error2 = (hists.h_drift_depth_adc2_[i_ring]->getBinContent(j, i) -
487  hists.h_drift_depth_adc_[i_ring]->getBinContent(j, i) *
488  hists.h_drift_depth_adc_[i_ring]->getBinContent(j, i) /
489  hists.h_drift_depth_noadc_[i_ring]->getBinContent(j, i)) /
490  hists.h_drift_depth_noadc_[i_ring]->getBinContent(j, i);
491 
492  hists.h_drift_depth_adc_[i_ring]->setBinError(j, i, sqrt(adc_error2));
493  double error2 = adc_error2 / (hists.h_drift_depth_noadc_[i_ring]->getBinContent(j, i) - 1.);
494  hists.h_drift_depth_[i_ring]->setBinError(j, i, sqrt(error2));
495  } else {
496  hists.h_drift_depth_[i_ring]->setBinError(j, i, 0);
497  hists.h_drift_depth_adc_[i_ring]->setBinError(j, i, 0);
498  }
499  h_drift_depth_adc_slice_->setBinContent(j, hists.h_drift_depth_adc_[i_ring]->getBinContent(j, i));
500  h_drift_depth_adc_slice_->setBinError(j, hists.h_drift_depth_adc_[i_ring]->getBinError(j, i));
501  nentries += hists.h_drift_depth_noadc_[i_ring]->getBinContent(j, i);
502  } // end loop over bins in drift width
503 
504  double mean = h_drift_depth_adc_slice_->getMean(1);
505  double error = 0;
506  if (nentries != 0) {
507  error = h_drift_depth_adc_slice_->getRMS(1) / std::sqrt(nentries);
508  }
509  hists.h_mean_[i_ring]->setBinContent(i, mean);
510  hists.h_mean_[i_ring]->setBinError(i, error);
511 }
512 
513 //------------------------------------------------------------------------------
515  std::shared_ptr<SiPixelLorentzAngle> theLAPayload, int i_index, int i_layer, int i_module) {
516  // output results
518 
519  // B-field value
520  // nominalValue returns the magnetic field value in kgauss (1T = 10 kgauss)
521  float theMagField = magField->nominalValue() / 10.;
522 
523  double half_width = width_ * 10000 / 2; // pixel half thickness in units of micro meter
524 
525  f1 = std::make_unique<TF1>("f1", "[0] + [1]*x + [2]*x*x + [3]*x*x*x + [4]*x*x*x*x + [5]*x*x*x*x*x", 5., 280.);
526  f1->SetParName(0, "offset");
527  f1->SetParName(1, "tan#theta_{LA}");
528  f1->SetParName(2, "quad term");
529  f1->SetParName(3, "cubic term");
530  f1->SetParName(4, "quartic term");
531  f1->SetParName(5, "quintic term");
532 
533  f1->SetParameter(0, 0);
534  f1->SetParError(0, 0);
535  f1->SetParameter(1, 0.4);
536  f1->SetParError(1, 0);
537  f1->SetParameter(2, 0.0);
538  f1->SetParError(2, 0);
539  f1->SetParameter(3, 0.0);
540  f1->SetParError(3, 0);
541  f1->SetParameter(4, 0.0);
542  f1->SetParError(4, 0);
543  f1->SetParameter(5, 0.0);
544  f1->SetParError(5, 0);
545  f1->SetChisquare(0);
546 
547  hists.h_mean_[i_index]->getTH1()->Fit(f1.get(), "ERQ");
548 
549  res.p0 = f1->GetParameter(0);
550  res.e0 = f1->GetParError(0);
551  res.p1 = f1->GetParameter(1);
552  res.e1 = f1->GetParError(1);
553  res.p2 = f1->GetParameter(2);
554  res.e2 = f1->GetParError(2);
555  res.p3 = f1->GetParameter(3);
556  res.e3 = f1->GetParError(3);
557  res.p4 = f1->GetParameter(4);
558  res.e4 = f1->GetParError(4);
559  res.p5 = f1->GetParameter(5);
560  res.e5 = f1->GetParError(5);
561  res.chi2 = f1->GetChisquare();
562  res.ndf = f1->GetNDF();
563  res.prob = f1->GetProb();
564  res.redChi2 = res.ndf > 0. ? res.chi2 / res.ndf : 0.;
565 
566  double f1_halfwidth = res.p0 + res.p1 * half_width + res.p2 * pow(half_width, 2) + res.p3 * pow(half_width, 3) +
567  res.p4 * pow(half_width, 4) + res.p5 * pow(half_width, 5);
568 
569  double f1_zerowidth = res.p0;
570 
571  // tan_LA = (f1(x = half_width) - f1(x = 0)) / (half_width - 0)
572  res.tan_LA = (f1_halfwidth - f1_zerowidth) / half_width;
573  double errsq_LA =
574  (pow(res.e1, 2) + pow((half_width * res.e2), 2) + pow((half_width * half_width * res.e3), 2) +
575  pow((half_width * half_width * half_width * res.e4), 2) +
576  pow((half_width * half_width * half_width * half_width * res.e5), 2)); // Propagation of uncertainty
577  res.error_LA = sqrt(errsq_LA);
578 
579  hists.h_bySectMeasLA_->setBinContent(i_index, (res.tan_LA / theMagField));
580  hists.h_bySectMeasLA_->setBinError(i_index, (res.error_LA / theMagField));
581  hists.h_bySectChi2_->setBinContent(i_index, res.redChi2);
582 
583  int nentries = hists.h_bySectOccupancy_->getBinContent(i_index); // number of on track hits in that sector
584 
585  bool isNew = (i_index > hists.nlay * hists.nModules_[hists.nlay - 1]);
586  int shiftIdx = i_index - hists.nlay * hists.nModules_[hists.nlay - 1] - 1;
587 
588  LogDebug("SiPixelLorentzAnglePCLHarvester")
589  << " isNew: " << isNew << " i_index: " << i_index << " shift index: " << shiftIdx;
590 
591  const auto& detIdsToFill =
592  isNew ? std::vector<unsigned int>({hists.BPixnewDetIds_[shiftIdx]}) : hists.detIdsList.at(i_index);
593 
594  LogDebug("SiPixelLorentzAnglePCLHarvester")
595  << "index: " << i_index << " i_module: " << i_module << " i_layer: " << i_layer;
596  for (const auto& id : detIdsToFill) {
597  LogDebug("SiPixelLorentzAnglePCLHarvester") << id << ",";
598  }
599 
600  float LorentzAnglePerTesla_;
601  float currentLA = currentLorentzAngle->getLorentzAngle(detIdsToFill.front());
602  // if the fit quality is OK
603  if ((res.redChi2 != 0.) && (res.redChi2 < fitChi2Cut_) && (nentries > minHitsCut_)) {
604  LorentzAnglePerTesla_ = res.tan_LA / theMagField;
605  // fill the LA actually written to payload
606  hists.h_bySectSetLA_->setBinContent(i_index, LorentzAnglePerTesla_);
607  hists.h_bySectRejectLA_->setBinContent(i_index, 0.);
608  hists.h_bySectLA_->setBinContent(i_index, LorentzAnglePerTesla_);
609 
610  const auto& deltaLA = (LorentzAnglePerTesla_ - currentLA);
611  hists.h_bySectDeltaLA_->setBinContent(i_index, deltaLA);
612 
613  for (const auto& id : detIdsToFill) {
614  if (!theLAPayload->putLorentzAngle(id, LorentzAnglePerTesla_)) {
615  edm::LogError("SiPixelLorentzAnglePCLHarvester") << "[SiPixelLorentzAnglePCLHarvester::fitAndStore]: detid ("
616  << i_layer << "," << i_module << ") already exists";
617  }
618  }
619  } else {
620  // just copy the values from the existing payload
621  hists.h_bySectSetLA_->setBinContent(i_index, 0.);
622  hists.h_bySectRejectLA_->setBinContent(i_index, (res.tan_LA / theMagField));
623  hists.h_bySectLA_->setBinContent(i_index, currentLA);
624 
625  hists.h_bySectDeltaLA_->setBinContent(i_index, 0.);
626 
627  for (const auto& id : detIdsToFill) {
628  LorentzAnglePerTesla_ = currentLorentzAngle->getLorentzAngle(id);
629  if (!theLAPayload->putLorentzAngle(id, LorentzAnglePerTesla_)) {
630  edm::LogError("SiPixelLorentzAnglePCLHarvester") << "[SiPixelLorentzAnglePCLHarvester::fitAndStore]: detid ("
631  << i_layer << "," << i_module << ") already exists";
632  }
633  }
634  }
635  // return the struct of fit details
636  return res;
637 }
638 
639 //------------------------------------------------------------------------------
642  desc.setComment("Harvester module of the SiPixel Lorentz Angle PCL monitoring workflow");
643  desc.add<std::vector<std::string>>("newmodulelist", {})->setComment("the list of DetIds for new sensors");
644  desc.add<std::string>("dqmDir", "AlCaReco/SiPixelLorentzAngle")->setComment("the directory of PCL Worker output");
645  desc.add<double>("fitChi2Cut", 20.)->setComment("cut on fit chi2/ndof to accept measurement");
646  desc.add<int>("minHitsCut", 10000)->setComment("cut on minimum number of on-track hits to accept measurement");
647  desc.add<std::string>("record", "SiPixelLorentzAngleRcd")->setComment("target DB record");
648  descriptions.addWithDefaultLabel(desc);
649 }
650 
void setComment(std::string const &value)
const std::map< unsigned int, float > & getLorentzAngles() const
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
Base exception class for the object to relational access.
Definition: Exception.h:11
uint16_t *__restrict__ id
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:32
int moduleName() const
module id (index in z)
int nominalValue() const
The nominal field value for this map in kGauss.
Definition: MagneticField.h:49
void beginRun(const edm::Run &, const edm::EventSetup &) override
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void findMean(MonitorElement *h_drift_depth_adc_slice_, int i, int i_ring)
SiPixelLorentzAnglePCLHarvester(const edm::ParameterSet &)
SiPixelLAHarvest::fitResults fitAndStore(std::shared_ptr< SiPixelLorentzAngle > theLA, int i_idx, int i_lay, int i_mod)
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
const Bounds & bounds() const
Definition: Surface.h:87
unsigned int pxbModule(const DetId &id) const
SiPixelLorentzAngleCalibrationHistograms hists
edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > topoEsToken_
Log< level::Error, false > LogError
static void fillDescriptions(edm::ConfigurationDescriptions &)
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:37
unsigned int numberOfLayers(int subdet) const
constexpr std::array< uint8_t, layerIndexSize > layer
void Fill(long long x)
virtual float thickness() const =0
bool getData(T &iHolder) const
Definition: EventSetup.h:128
~SiPixelLorentzAnglePCLHarvester() override=default
virtual void Reset()
Remove all data from the ME, keept the empty histogram with all its settings.
edm::ESGetToken< SiPixelLorentzAngle, SiPixelLorentzAngleRcd > siPixelLAEsToken_
int bladeName() const
blade id
void setComment(std::string const &value)
virtual int getNbinsX() const
get # of bins in X-axis
const unsigned getPXBModules(unsigned int lay) const
T sqrt(T t)
Definition: SSEVec.h:19
std::unordered_map< uint32_t, std::vector< uint32_t > > detIdsList
if(conf_.getParameter< bool >("UseStripCablingDB"))
virtual MonitorElement * get(std::string const &fullpath) const
Definition: DQMStore.cc:673
const DetContainer & detsPXB() const
bool isAvailable() const
Definition: Service.h:40
edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > geomEsToken_
char const * what() const noexceptoverride
Definition: Exception.cc:103
Hash writeOneIOV(const T &payload, Time_t time, const std::string &recordName)
Transition
Definition: Transition.h:12
DetId geographicalId() const
The label of this GeomDet.
Definition: GeomDet.h:64
ParameterDescriptionBase * add(U const &iLabel, T const &value)
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)
virtual double getBinContent(int binx) const
get content of bin (1-D)
Log< level::Warning, true > LogPrint
PXFDetId getDetId()
return DetId
unsigned int pxbLayer(const DetId &id) const
virtual double getMean(int axis=1) const
get mean value of histogram along x, y or z axis (axis=1, 2, 3 respectively)
virtual void setBinContent(int binx, double content)
set content of bin (1-D)
virtual double getRMS(int axis=1) const
get RMS of histogram along x, y or z axis (axis=1, 2, 3 respectively)
int layerName() const
layer id
float getLorentzAngle(const uint32_t &) const
int diskName() const
disk id
virtual void setBinError(int binx, double error)
set uncertainty on content of bin (1-D)
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::ESGetToken< MagneticField, IdealMagneticFieldRecord > magneticFieldToken_
PXBDetId getDetId()
return the DetId
ESGetTokenH3DDVariant esConsumes(std::string const &Reccord, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
Definition: Run.h:45
tuple module
Definition: callgraph.py:69
void dqmEndJob(DQMStore::IBooker &, DQMStore::IGetter &) override
#define LogDebug(id)
unsigned transform(const HcalDetId &id, unsigned transformCode)