test
CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PixelLumiDQM.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 
3 // Package: PixelLumiDQM
4 // Class: PixelLumiDQM
5 
6 // Author: Amita Raval
7 // Based on Jeroen Hegeman's code for Pixel Cluster Count Luminosity
8 
9 #include "PixelLumiDQM.h"
10 
37 
38 #include <fstream>
39 #include <map>
40 #include <string>
41 #include <sys/time.h>
42 #include <time.h>
43 #include <vector>
44 
45 const unsigned int PixelLumiDQM::lastBunchCrossing;
46 
47 // Constructors and destructor.
49  fPixelClusterLabel(consumes<edmNew::DetSetVector<SiPixelCluster> >(iConfig.getUntrackedParameter<edm::InputTag>("pixelClusterLabel",
50  edm::InputTag("siPixelClusters")))),
51  fIncludePixelClusterInfo(iConfig.getUntrackedParameter<bool>("includePixelClusterInfo", true)),
52  fIncludePixelQualCheckHistos(iConfig.getUntrackedParameter<bool>("includePixelQualCheckHistos", true)),
53  fResetIntervalInLumiSections(iConfig.getUntrackedParameter<int>("resetEveryNLumiSections", 1)),
54  fDeadModules(iConfig.getUntrackedParameter<std::vector<uint32_t> >("deadModules", std::vector<uint32_t>())),
55  fMinPixelsPerCluster(iConfig.getUntrackedParameter<int>("minNumPixelsPerCluster", 0)),
56  fMinClusterCharge(iConfig.getUntrackedParameter<double>("minChargePerCluster", 0)),
57  bunchTriggerMask(lastBunchCrossing+1,false),
58  filledAndUnmaskedBunches(0),
59  useInnerBarrelLayer(iConfig.getUntrackedParameter<bool>("useInnerBarrelLayer", false)),
60  fLogFileName_(iConfig.getUntrackedParameter<std::string>("logFileName","/tmp/pixel_lumi.txt"))
61 {
62  edm::LogInfo("Configuration")
63  << "PixelLumiDQM looking for pixel clusters in '"
64  << iConfig.getUntrackedParameter<edm::InputTag>("pixelClusterLabel",
65  edm::InputTag("siPixelClusters")) << "'";
66  edm::LogInfo("Configuration")
67  << "PixelLumiDQM storing pixel cluster info? "
69  edm::LogInfo("Configuration")
70  << "PixelLumiDQM storing pixel cluster quality check histograms? "
72 
73  if (!fDeadModules.size()) {
74  edm::LogInfo("Configuration")
75  << "No pixel modules specified to be ignored";
76  } else {
77  edm::LogInfo("Configuration")
78  << fDeadModules.size() << " pixel modules specified to be ignored:";
79  for (std::vector<uint32_t>::const_iterator it = fDeadModules.begin();
80  it != fDeadModules.end(); ++it) {
81  edm::LogInfo("Configuration")
82  << " " << *it;
83  }
84  }
85  edm::LogInfo("Configuration")
86  << "Ignoring pixel clusters with less than "
87  << fMinPixelsPerCluster << " pixels";
88  edm::LogInfo("Configuration")
89  << "Ignoring pixel clusters with charge less than "
91 }
92 
93 
95 {
96 }
97 
98 
99 void
102  desc.setUnknown();
103  descriptions.addDefault(desc);
104 }
105 
106 
107 void
109  const edm::EventSetup& iSetup)
110 {
111  // Collect all bookkeeping information.
112  fRunNo = iEvent.id().run();
113  fEvtNo = iEvent.id().event();
115  fBXNo = iEvent.bunchCrossing();
116  fTimestamp = iEvent.time().unixTime();
117  fHistBunchCrossings->Fill(float(fBXNo));
119  // This serves as event counter to compute luminosity from cluster counts.
120  std::map<int, PixelClusterCount>::iterator it = fNumPixelClusters.find(fBXNo);
122 
124  // Find tracker geometry.
126  iSetup.get<TrackerDigiGeometryRecord>().get(trackerGeo);
127 
128  // Find pixel clusters.
130  iEvent.getByToken(fPixelClusterLabel, pixelClusters);
131 
132  // Loop over entire tracker geometry.
133  for (TrackerGeometry::DetContainer::const_iterator
134  i = trackerGeo->dets().begin();
135  i != trackerGeo->dets().end(); ++i) {
136  // See if this is a pixel unit(?).
137 
138  if ( GeomDetEnumerators::isTrackerPixel((*i)->subDetector())) {
139  DetId detId = (*i)->geographicalId();
140  // Find all clusters on this detector module.
142  pixelClusters->find(detId);
143  if (iSearch != pixelClusters->end()) {
144 
145  // Count the number of clusters with at least a minimum
146  // number of pixels per cluster and at least a minimum charge.
147  size_t numClusters = 0;
149  itClus = iSearch->begin();
150  itClus != iSearch->end(); ++itClus) {
151  if ((itClus->size() >= fMinPixelsPerCluster) &&
152  (itClus->charge() >= fMinClusterCharge)) {
153  ++numClusters;
154  }
155  }
156  // DEBUG DEBUG DEBUG
157  assert(numClusters <= iSearch->size());
158  // DEBUG DEBUG DEBUG end
159 
160  // Add up the cluster count based on the position of this detector element.
161  if (detId.subdetId() == PixelSubdetector::PixelBarrel) {
162  PixelBarrelName detName = PixelBarrelName(detId);
163  int layer = detName.layerName();
164  fNumPixelClusters[fBXNo].numB.at(layer - 1) += numClusters;
165  fNumPixelClusters[fBXNo].dnumB.at(layer - 1) += sqrt(numClusters);
166  } else {
167  // DEBUG DEBUG DEBUG
169  // DEBUG DEBUG DEBUG end
170 
171  PixelEndcapName detName = PixelEndcapName(detId);
172  PixelEndcapName::HalfCylinder halfCylinder = detName.halfCylinder();
173  int disk = detName.diskName();
174  switch (halfCylinder) {
175  case PixelEndcapName::mO:
176  case PixelEndcapName::mI:
177  fNumPixelClusters[fBXNo].numFM.at(disk - 1) += numClusters;
178  fNumPixelClusters[fBXNo].dnumFM.at(disk - 1) += sqrt(numClusters);
179  break;
180  case PixelEndcapName::pO:
181  case PixelEndcapName::pI:
182  fNumPixelClusters[fBXNo].numFP.at(disk - 1) += numClusters;
183  fNumPixelClusters[fBXNo].dnumFP.at(disk - 1) += sqrt(numClusters);
184  break;
185  default:
186  assert(false);
187  break;
188  }
189  }
190  }
191  }
192  }
193  }
194  // ----------
195 
196  // Fill some pixel cluster quality check histograms if requested.
198 
199  // Find tracker geometry.
201  iSetup.get<TrackerDigiGeometryRecord>().get(trackerGeo);
202 
203  // Find pixel clusters.
205  iEvent.getByToken(fPixelClusterLabel, pixelClusters);
206 
207  bool filterDeadModules = (fDeadModules.size() > 0);
208  std::vector<uint32_t>::const_iterator deadModulesBegin = fDeadModules.begin();
209  std::vector<uint32_t>::const_iterator deadModulesEnd = fDeadModules.end();
210 
211  // Loop over entire tracker geometry.
212  for (TrackerGeometry::DetContainer::const_iterator
213  i = trackerGeo->dets().begin();
214  i != trackerGeo->dets().end(); ++i) {
215 
216  // See if this is a pixel module.
217  if ( GeomDetEnumerators::isTrackerPixel((*i)->subDetector())) {
218  DetId detId = (*i)->geographicalId();
219 
220  // Skip this module if it's on the list of modules to be ignored.
221  if (filterDeadModules &&
222  find(deadModulesBegin, deadModulesEnd, detId()) != deadModulesEnd) {
223  continue;
224  }
225 
226  // Find all clusters in this module.
228  pixelClusters->find(detId);
229 
230  // Loop over all clusters in this module.
231  if (iSearch != pixelClusters->end()) {
233  clus != iSearch->end(); ++clus) {
234 
235  if ((clus->size() >= fMinPixelsPerCluster) &&
236  (clus->charge() >= fMinClusterCharge)) {
237 
238  PixelGeomDetUnit const* theGeomDet =
239  dynamic_cast<PixelGeomDetUnit const*>(trackerGeo->idToDet(detId));
240  PixelTopology const* topol = &(theGeomDet->specificTopology());
241  double x = clus->x();
242  double y = clus->y();
243  LocalPoint clustLP = topol->localPosition(MeasurementPoint(x, y));
244  GlobalPoint clustGP = theGeomDet->surface().toGlobal(clustLP);
245  double charge = clus->charge() / 1.e3;
246  int size = clus->size();
247 
248  if (detId.subdetId() == PixelSubdetector::PixelBarrel) {
249  PixelBarrelName detName = PixelBarrelName(detId);
250  int layer = detName.layerName();
251  switch (layer) {
252  case 1:
253  fHistContainerThisRun["clusPosBarrel1"]->Fill(clustGP.z(), clustGP.phi());
254  fHistContainerThisRun["clusChargeBarrel1"]->Fill(iEvent.bunchCrossing(), charge);
255  fHistContainerThisRun["clusSizeBarrel1"]->Fill(iEvent.bunchCrossing(), size);
256  break;
257  case 2:
258  fHistContainerThisRun["clusPosBarrel2"]->Fill(clustGP.z(), clustGP.phi());
259  fHistContainerThisRun["clusChargeBarrel2"]->Fill(iEvent.bunchCrossing(), charge);
260  fHistContainerThisRun["clusSizeBarrel2"]->Fill(iEvent.bunchCrossing(), size);
261  break;
262  case 3:
263  fHistContainerThisRun["clusPosBarrel3"]->Fill(clustGP.z(), clustGP.phi());
264  fHistContainerThisRun["clusChargeBarrel3"]->Fill(iEvent.bunchCrossing(), charge);
265  fHistContainerThisRun["clusSizeBarrel3"]->Fill(iEvent.bunchCrossing(), size);
266  break;
267  default:
268  assert(false);
269  break;
270  }
271  } else {
272 
273  // DEBUG DEBUG DEBUG
275  // DEBUG DEBUG DEBUG end
276 
277  PixelEndcapName detName = PixelEndcapName(detId);
278  PixelEndcapName::HalfCylinder halfCylinder = detName.halfCylinder();
279  int disk = detName.diskName();
280  switch (halfCylinder) {
281  case PixelEndcapName::mO:
282  case PixelEndcapName::mI:
283  switch (disk) {
284  case 1:
285  fHistContainerThisRun["clusPosEndCapM1"]->Fill(clustGP.x(), clustGP.y());
286  fHistContainerThisRun["clusChargeEndCapM1"]->Fill(iEvent.bunchCrossing(), charge);
287  fHistContainerThisRun["clusSizeEndCapM1"]->Fill(iEvent.bunchCrossing(), size);
288  break;
289  case 2:
290  fHistContainerThisRun["clusPosEndCapM2"]->Fill(clustGP.x(), clustGP.y());
291  fHistContainerThisRun["clusChargeEndCapM2"]->Fill(iEvent.bunchCrossing(), charge);
292  fHistContainerThisRun["clusSizeEndCapM2"]->Fill(iEvent.bunchCrossing(), size);
293  break;
294  default:
295  assert(false);
296  break;
297  }
298  break;
299  case PixelEndcapName::pO:
300  case PixelEndcapName::pI:
301  switch (disk) {
302  case 1:
303  fHistContainerThisRun["clusPosEndCapP1"]->Fill(clustGP.x(), clustGP.y());
304  fHistContainerThisRun["clusChargeEndCapP1"]->Fill(iEvent.bunchCrossing(), charge);
305  fHistContainerThisRun["clusSizeEndCapP1"]->Fill(iEvent.bunchCrossing(), size);
306  break;
307  case 2:
308  fHistContainerThisRun["clusPosEndCapP2"]->Fill(clustGP.x(), clustGP.y());
309  fHistContainerThisRun["clusChargeEndCapP2"]->Fill(iEvent.bunchCrossing(), charge);
310  fHistContainerThisRun["clusSizeEndCapP2"]->Fill(iEvent.bunchCrossing(), size);
311  break;
312  default:
313  assert(false);
314  break;
315  }
316  break;
317  default:
318  assert(false);
319  break;
320  }
321  }
322  }
323  }
324  }
325  }
326  }
327  }
328 }
329 
330 void
332  edm::Run const & run,
333  edm::EventSetup const & /* iSetup */)
334 {
335  edm::LogInfo("Status") << "Starting processing of run #" << run.id().run();
336 
337  // Top folder containing high-level information about pixel and HF lumi.
338  std::string folder = "PixelLumi/";
339  ibooker.setCurrentFolder(folder);
340 
341  fHistTotalRecordedLumiByLS = ibooker.book1D("totalPixelLumiByLS","Pixel Lumi in nb vs LS",8000,0.5,8000.5);
342  fHistRecordedByBxCumulative = ibooker.book1D("PXLumiByBXsum","Pixel Lumi in nb by BX Cumulative vs LS",lastBunchCrossing,
343  0.5,float(lastBunchCrossing)+0.5);
344 
345  std::string subfolder = folder + "lastLS/";
346  ibooker.setCurrentFolder(subfolder);
347  fHistRecordedByBxLastLumi = ibooker.book1D("PXByBXLastLumi","Pixel By BX Last Lumi",lastBunchCrossing+1,
348  -0.5,float(lastBunchCrossing)+0.5);
349 
350  subfolder = folder+"ClusterCountingDetails/";
351  ibooker.setCurrentFolder(subfolder);
352 
353  fHistnBClusVsLS[0] = ibooker.book1D("nBClusVsLS_0","Fraction of Clusters vs LS Barrel layer 0",8000,0.5,8000.5);
354  fHistnBClusVsLS[1] = ibooker.book1D("nBClusVsLS_1","Fraction of Clusters vs LS Barrel layer 1",8000,0.5,8000.5);
355  fHistnBClusVsLS[2] = ibooker.book1D("nBClusVsLS_2","Fraction of Clusters vs LS Barrel layer 2",8000,0.5,8000.5);
356  fHistnFPClusVsLS[0] = ibooker.book1D("nFPClusVsLS_0","Fraction of Clusters vs LS Barrel layer 0",8000,0.5,8000.5);
357  fHistnFPClusVsLS[1] = ibooker.book1D("nFPClusVsLS_1","Fraction of Clusters vs LS Barrel layer 1",8000,0.5,8000.5);
358  fHistnFMClusVsLS[0] = ibooker.book1D("nFMClusVsLS_0","Fraction of Clusters vs LS Barrel layer 0",8000,0.5,8000.5);
359  fHistnFMClusVsLS[1] = ibooker.book1D("nFMClusVsLS_1","Fraction of Clusters vs LS Barrel layer 1",8000,0.5,8000.5);
360  fHistBunchCrossings = ibooker.book1D("BunchCrossings","Cumulative Bunch Crossings",lastBunchCrossing,
361  0.5,float(lastBunchCrossing)+0.5);
362  fHistBunchCrossingsLastLumi = ibooker.book1D("BunchCrossingsLL","Bunch Crossings Last Lumi",lastBunchCrossing,
363  0.5,float(lastBunchCrossing)+0.5);
364  fHistClusterCountByBxLastLumi = ibooker.book1D("ClusterCountByBxLL","Cluster Count by BX Last Lumi",lastBunchCrossing,
365  0.5,float(lastBunchCrossing)+0.5);
366  fHistClusterCountByBxCumulative = ibooker.book1D("ClusterCountByBxSum","Cluster Count by BX Cumulative",lastBunchCrossing,
367  0.5,float(lastBunchCrossing)+0.5);
368  fHistClusByLS = ibooker.book1D("totalClusByLS","Number of Clusters all dets vs LS",8000,0.5,8000.5);
369 
370  // Add some pixel cluster quality check histograms (in a subfolder).
371  subfolder = folder+"qualityChecks/";
372  ibooker.setCurrentFolder(subfolder);
373 
375  // Create histograms for this run if not already present in our list.
376  edm::LogInfo("Status") << "Creating histograms for run #" << run.id().run();
377 
378  // Pixel cluster positions in the barrel - (z, phi).
379  for (size_t i = 1; i <= kNumLayers; ++i) {
380  std::stringstream key;
381  key << "clusPosBarrel" << i;
382  std::stringstream name;
383  name << key.str() << "_" << run.run();
384  std::stringstream title;
385  title << "Pixel cluster position - barrel layer " << i;
386  fHistContainerThisRun[key.str()] = ibooker.book2D(name.str().c_str(),
387  title.str().c_str(),
388  100, -30., 30.,
389  64, -Geom::pi(), Geom::pi());
390  }
391 
392  // Pixel cluster positions in the endcaps (x, y).
393  std::vector<std::string> sides;
394  sides.push_back("M");
395  sides.push_back("P");
396  for (std::vector<std::string>::const_iterator side = sides.begin();
397  side != sides.end(); ++side) {
398  for (size_t i = 1; i <= kNumDisks; ++i) {
399  std::stringstream key;
400  key << "clusPosEndCap" << *side << i;
401  std::stringstream name;
402  name << key.str() << "_" << run.run();
403  std::stringstream title;
404  title << "Pixel cluster position - endcap disk " << i;
405  fHistContainerThisRun[key.str()] = ibooker.book2D(name.str().c_str(),
406  title.str().c_str(),
407  100, -20., 20.,
408  100, -20., 20.);
409  }
410  }
411 
412  // Pixel cluster charge in the barrel, per bx.
413  for (size_t i = 1; i <= kNumLayers; ++i) {
414  std::stringstream key;
415  key << "clusChargeBarrel" << i;
416  std::stringstream name;
417  name << key.str() << "_" << run.run();
418  std::stringstream title;
419  title << "Pixel cluster charge - barrel layer " << i;
420  fHistContainerThisRun[key.str()] = ibooker.book2D(name.str().c_str(),
421  title.str().c_str(),
422  3564, .5, 3564.5,
423  100, 0., 100.);
424  }
425 
426  // Pixel cluster charge in the endcaps, per bx.
427  for (std::vector<std::string>::const_iterator side = sides.begin();
428  side != sides.end(); ++side) {
429  for (size_t i = 1; i <= kNumDisks; ++i) {
430  std::stringstream key;
431  key << "clusChargeEndCap" << *side << i;
432  std::stringstream name;
433  name << key.str() << "_" << run.run();
434  std::stringstream title;
435  title << "Pixel cluster charge - endcap disk " << i;
436  fHistContainerThisRun[key.str()] = ibooker.book2D(name.str().c_str(),
437  title.str().c_str(),
438  3564, .5, 3564.5,
439  100, 0., 100.);
440  }
441  }
442 
443  // Pixel cluster size in the barrel, per bx.
444  for (size_t i = 1; i <= kNumLayers; ++i) {
445  std::stringstream key;
446  key << "clusSizeBarrel" << i;
447  std::stringstream name;
448  name << key.str() << "_" << run.run();
449  std::stringstream title;
450  title << "Pixel cluster size - barrel layer " << i;
451  fHistContainerThisRun[key.str()] = ibooker.book2D(name.str().c_str(),
452  title.str().c_str(),
453  3564, .5, 3564.5,
454  100, 0., 100.);
455  }
456 
457  // Pixel cluster size in the endcaps, per bx.
458  for (std::vector<std::string>::const_iterator side = sides.begin();
459  side != sides.end(); ++side) {
460  for (size_t i = 1; i <= kNumDisks; ++i) {
461  std::stringstream key;
462  key << "clusSizeEndCap" << *side << i;
463  std::stringstream name;
464  name << key.str() << "_" << run.run();
465  std::stringstream title;
466  title << "Pixel cluster size - endcap disk " << i;
467  fHistContainerThisRun[key.str()] = ibooker.book2D(name.str().c_str(),
468  title.str().c_str(),
469  3564, .5, 3564.5,
470  100, 0., 100.);
471  }
472  }
473  }
474 }
475 
476 // ------------ Method called when ending the processing of a run. ------------
477 void
479 {
480 }
481 
482 // ------------ Method called when ending the processing of a run. ------------
483 void
485 {
486 }
487 
488 
489 // ------------ Method called when starting to process a luminosity block. ------------
490 void
492  edm::EventSetup const&)
493 {
494  // Only reset and fill every fResetIntervalInLumiSections (default is 1 LS)
495  // Return unless the PREVIOUS LS was at the right modulo value
496  // (e.g. is resetinterval = 5 the rest will only be executed at LS=6
497  // NB: reset is done here so the histograms by LS are sent before resetting.
498  // NB: not being used for now since default is 1 LS. There is a bug here.
499 
500  unsigned int ls = lumiBlock.luminosityBlockAuxiliary().luminosityBlock();
501 
502  if((ls-1)%fResetIntervalInLumiSections==0){
506  }
507 }
508 
509 
510 // ------------ Method called when ending the processing of a luminosity block. ------------
511 void
513  edm::EventSetup const&es)
514 {
515 
516  unsigned int ls = lumiBlock.luminosityBlockAuxiliary().luminosityBlock();
517 
518  // Only fill every fResetIntervalInLumiSections (default is 1 LS)
519  if(ls%fResetIntervalInLumiSections!=0) return;
520 
521  printf("Lumi Block = %d\n",ls);
522 
523  if((ls-1)%fResetIntervalInLumiSections==0){
524  }
525 
526  unsigned int nBClus[3] = {0,0,0};
527  unsigned int nFPClus[2] = {0, 0};
528  unsigned int nFMClus[2] = {0, 0};
529 
530  double total_recorded = 0.;
531  double total_recorded_unc_square = 0.;
532 
533  // Obtain bunch-by-bunch cluster counts and compute totals for lumi calculation.
534  double totalcounts = 0.0;
535  double etotalcounts = 0.0;
536  double totalevents = 0.0;
537  double lumi_factor_per_bx = 0.0;
540  else
542 
543  for(std::map<int, PixelClusterCount>::iterator it = fNumPixelClusters.begin();
544  it != fNumPixelClusters.end(); it++) {
545 
546  // Sum all clusters for this BX.
547  unsigned int total = (*it).second.numB.at(1)+
548  (*it).second.numB.at(2)+(*it).second.numFP.at(0)+(*it).second.numFP.at(1)+
549  (*it).second.numFM.at(0)+(*it).second.numFM.at(1);
550  if(useInnerBarrelLayer) total+=(*it).second.numB.at(0);
551  totalcounts += total;
552  double etotal = (*it).second.dnumB.at(1)+
553  (*it).second.dnumB.at(2)+(*it).second.dnumFP.at(0)+(*it).second.dnumFP.at(1)+
554  (*it).second.dnumFM.at(0)+(*it).second.dnumFM.at(1);
555  if(useInnerBarrelLayer) etotal = (*it).second.dnumB.at(0);
556  etotalcounts += etotal;
557  etotal = sqrt(etotal);
558 
559  fHistClusterCountByBxLastLumi->setBinContent((*it).first,total);
560  fHistClusterCountByBxLastLumi->setBinError((*it).first,etotal);
562 
563  unsigned int events_per_bx = fHistBunchCrossingsLastLumi->getBinContent((*it).first);
564  totalevents += events_per_bx;
565  double average_cluster_count = events_per_bx !=0 ? double(total)/events_per_bx : 0.;
566  double average_cluster_count_unc = events_per_bx!=0 ? etotal/events_per_bx : 0.;
567  double pixel_bx_lumi_per_ls = lumi_factor_per_bx * average_cluster_count / CM2_TO_NANOBARN ;
568  double pixel_bx_lumi_per_ls_unc = 0.0;
570  pixel_bx_lumi_per_ls_unc = sqrt(lumi_factor_per_bx*lumi_factor_per_bx *
571  (average_cluster_count_unc*average_cluster_count_unc +
572  (average_cluster_count* XSEC_PIXEL_CLUSTER_UNC /
574  (average_cluster_count* XSEC_PIXEL_CLUSTER /
576  ) / CM2_TO_NANOBARN ;
577  else
578  pixel_bx_lumi_per_ls_unc = sqrt(lumi_factor_per_bx*lumi_factor_per_bx *
579  (average_cluster_count_unc*average_cluster_count_unc +
580  (average_cluster_count* rXSEC_PIXEL_CLUSTER_UNC /
582  (average_cluster_count* rXSEC_PIXEL_CLUSTER /
584  ) / CM2_TO_NANOBARN ;
585 
586  fHistRecordedByBxLastLumi->setBinContent((*it).first,pixel_bx_lumi_per_ls);
587  fHistRecordedByBxLastLumi->setBinError((*it).first,pixel_bx_lumi_per_ls_unc);
588 
591  pixel_bx_lumi_per_ls);
592 
593  /*
594  if(fHistRecordedByBxLastLumi->getBinContent((*it).first)!=0.)
595  fHistRecordedByBxLastLumi->getBinContent((*it).first));
596  if(fHistRecordedByBxCumulative->getBinContent((*it).first)!=0.)
597  fHistRecordedByBxCumulative->getBinContent((*it).first));
598  */
599 
600  nBClus[0] +=(*it).second.numB.at(0);
601  nBClus[1] +=(*it).second.numB.at(1);
602  nBClus[2] +=(*it).second.numB.at(2);
603  nFPClus[0] +=(*it).second.numFP.at(0);
604  nFPClus[1] +=(*it).second.numFP.at(1);
605  nFMClus[0] +=(*it).second.numFM.at(0);
606  nFMClus[1] +=(*it).second.numFM.at(1);
607 
608  // Reset counters
609  (*it).second.Reset();
610 
611  // std::cout << "bx="<< (*it).first << " clusters=" << (*it).second.numB.at(0)) << std::endl;
612  }
613 
615  for(unsigned int i = 0; i<= lastBunchCrossing; i++){
616  if(bunchTriggerMask[i]){
617  double err = fHistRecordedByBxLastLumi->getBinError(i);
618  total_recorded += fHistRecordedByBxLastLumi->getBinContent(i);
619  total_recorded_unc_square += err*err;
620  }
621  }
622 
623  // Replace the total obtained by summing over BXs with the average per BX from the total cluster count and rescale
624 
625  if(totalevents > 10){
626  total_recorded = lumi_factor_per_bx * totalcounts / totalevents / CM2_TO_NANOBARN ;
627  }
628  else total_recorded = 0.0;
629 
630  std::cout << " Total recorded " << total_recorded << std::endl;
631  fHistTotalRecordedLumiByLS->setBinContent(ls,total_recorded);
633  sqrt(total_recorded_unc_square));
634  }
635  // fill cluster counts by detector regions for sanity checks
636  unsigned int all_detectors_counts = 0;
637  for(unsigned int i = 0; i < 3; i++){
638  all_detectors_counts+=nBClus[i];
639  }
640  for(unsigned int i = 0; i < 2; i++){
641  all_detectors_counts+=nFPClus[i];
642  }
643  for(unsigned int i = 0; i < 2; i++){
644  all_detectors_counts+=nFMClus[i];
645  }
646 
647  fHistClusByLS->setBinContent(ls, all_detectors_counts);
648 
649  for(unsigned int i = 0; i < 3; i++){
651  float(nBClus[i])/float(all_detectors_counts));
652  }
653  for(unsigned int i = 0; i < 2; i++){
655  float(nFPClus[i])/float(all_detectors_counts));
656  }
657  for(unsigned int i = 0; i < 2; i++){
659  float(nFMClus[i])/float(all_detectors_counts));
660  }
661 
663 
664  timeval tv;
665  gettimeofday(&tv,0);
666  tm *ts = gmtime(&tv.tv_sec);
667  char datestring[256];
668  strftime(datestring, sizeof(datestring),"%Y.%m.%d %T GMT %s",ts);
669  logFile_ << "RunNumber "<< fRunNo << std::endl;
670  logFile_ << "EndTimeOfFit " << datestring << std::endl;
671  logFile_ << "LumiRange "<< ls << "-" << ls << std::endl;
672  logFile_ << "Fill "<< -99 << std::endl;
673  logFile_ << "ActiveBunchCrossings "<< filledAndUnmaskedBunches << std::endl;
674  logFile_ << "PixelLumi "<< fHistTotalRecordedLumiByLS->getBinContent(ls) * 0.98 << std::endl;
675  logFile_ << "HFLumi "<< -99 << std::endl;
676  logFile_ << "Ratio " << -99 << std::endl;
677  logFile_.close();
678 }
679 
680 unsigned int PixelLumiDQM::calculateBunchMask(MonitorElement *e, std::vector<bool> &mask){
681  unsigned int nbins = e->getNbinsX();
682  std::vector<float> ar(nbins,0.);
683  for(unsigned int i = 1; i<= nbins; i++){
684  ar[i] = e->getBinContent(i);
685  }
686  return calculateBunchMask(ar,nbins,mask);
687 }
688 unsigned int PixelLumiDQM::calculateBunchMask(std::vector<float> &e, unsigned int nbins, std::vector<bool> &mask){
689  // Take the cumulative cluster count histogram and find max and average of non-empty bins.
690  unsigned int active_count = 0;
691  double maxc = 0.0;
692  double ave = 0.0; // Average of non-empty bins
693  unsigned int non_empty_bins = 0;
694 
695  for(unsigned int i = 1; i<= nbins; i++){
696  double bin = e[i];
697  if(bin !=0.0){
698  if(maxc<bin) maxc = bin;
699  ave += bin;
700  non_empty_bins++;
701  }
702  }
703 
704  ave /= non_empty_bins;
705  std::cout << "Bunch mask finder - non empty bins " << non_empty_bins
706  << " average of non empty bins " << ave
707  << " max content of one bin " << maxc << std::endl;
708  double mean = 0.;
709  double sigma = 0.;
710  if(non_empty_bins < 50){
711  mean = maxc;
712  sigma = (maxc)/20;
713  }
714  else{
715  TH1F dist("dist","dist",500,0.,maxc+(maxc/500.)*20.);
716  for(unsigned int i = 1; i<= nbins; i++){
717  double bin = e[i];
718  dist.Fill(bin);
719  }
720  dist.Fit("gaus","","",fmax(0.,ave-(maxc-ave)/5.),maxc);
721  TF1 *fit = dist.GetFunction("gaus");
722  mean = fit->GetParameter("Mean");
723  sigma = fit->GetParameter("Sigma");
724  }
725  std::cout << "Bunch mask will use mean" << mean << " sigma " << sigma << std::endl;
726  // Active BX defined as those which have nclus within fixed standard deviations of peak.
727  for(unsigned int i = 1; i<= nbins; i++){
728  double bin = e[i];
729  if(bin>0. && abs(bin-mean)<5.*(sigma)){ mask[i]=true; active_count++;}
730  }
731  std::cout << "Bunch mask finds " << active_count << " active bunch crossings " << std::endl;
732  // std::cout << "this is the full bx mask " ;
733  // for(unsigned int i = 1; i<= nbins; i++)
734  // std::cout << ((mask[i]) ? 1:0);
735  // std::cout << std::endl;
736  return active_count;
737 }
738 // Define this as a CMSSW plug-in.
RunNumber_t run() const
Definition: EventID.h:39
bool fIncludePixelQualCheckHistos
Definition: PixelLumiDQM.h:134
GlobalPoint toGlobal(const Point2DBase< Scalar, LocalTag > lp) const
Definition: Surface.h:106
EventNumber_t event() const
Definition: EventID.h:41
T getUntrackedParameter(std::string const &, T const &) const
static double SECONDS_PER_LS
Definition: PixelLumiDQM.h:43
int i
Definition: DBlmapReader.cc:9
unsigned int filledAndUnmaskedBunches
Definition: PixelLumiDQM.h:167
boost::transform_iterator< IterHelp, const_IdIter > const_iterator
void setBinContent(int binx, double content)
set content of bin (1-D)
const_iterator end(bool update=false) const
RunID const & id() const
Definition: RunBase.h:39
RunNumber_t run() const
Definition: RunBase.h:40
virtual LocalPoint localPosition(const MeasurementPoint &) const =0
static const unsigned int lastBunchCrossing
Definition: PixelLumiDQM.h:53
RunNumber_t run() const
Definition: RunID.h:39
virtual void analyze(const edm::Event &, const edm::EventSetup &)
std::vector< uint32_t > fDeadModules
Definition: PixelLumiDQM.h:141
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:462
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
assert(m_qm.get())
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
def ls
Definition: eostools.py:348
T y() const
Definition: PV3DBase.h:63
int bunchCrossing() const
Definition: EventBase.h:65
int fResetIntervalInLumiSections
Definition: PixelLumiDQM.h:135
data_type const * const_iterator
Definition: DetSetNew.h:30
virtual void endRun(edm::Run const &, edm::EventSetup const &)
MonitorElement * fHistRecordedByBxCumulative
Definition: PixelLumiDQM.h:164
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:42
bool useInnerBarrelLayer
Definition: PixelLumiDQM.h:168
void Fill(long long x)
LuminosityBlockNumber_t luminosityBlock() const
MonitorElement * fHistClusByLS
Definition: PixelLumiDQM.h:161
Measurement2DPoint MeasurementPoint
Measurement points are two-dimensional by default.
int iEvent
Definition: GenABIO.cc:230
void addDefault(ParameterSetDescription const &psetDescription)
MonitorElement * fHistnFMClusVsLS[2]
Definition: PixelLumiDQM.h:156
edm::EDGetTokenT< edmNew::DetSetVector< SiPixelCluster > > fPixelClusterLabel
Definition: PixelLumiDQM.h:115
unsigned int calculateBunchMask(MonitorElement *, std::vector< bool > &)
static double rXSEC_PIXEL_CLUSTER
Definition: PixelLumiDQM.h:50
T sqrt(T t)
Definition: SSEVec.h:18
T z() const
Definition: PV3DBase.h:64
virtual void endLuminosityBlock(edm::LuminosityBlock const &, edm::EventSetup const &)
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:115
MonitorElement * fHistClusterCountByBxLastLumi
Definition: PixelLumiDQM.h:159
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
MonitorElement * fHistTotalRecordedLumiByLS
Definition: PixelLumiDQM.h:162
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
unsigned int unixTime() const
Time in seconds since January 1, 1970.
Definition: Timestamp.h:46
virtual void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
static double FREQ_ORBIT
Definition: PixelLumiDQM.h:42
static double XSEC_PIXEL_CLUSTER_UNC
Definition: PixelLumiDQM.h:47
LuminosityBlock const & getLuminosityBlock() const
Definition: Event.h:85
void setBinError(int binx, double error)
set uncertainty on content of bin (1-D)
MonitorElement * fHistClusterCountByBxCumulative
Definition: PixelLumiDQM.h:160
MonitorElement * fHistnBClusVsLS[3]
Definition: PixelLumiDQM.h:154
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:37
LuminosityBlockNumber_t luminosityBlock() const
string key
FastSim: produces sample of signal events, overlayed with premixed minbias events.
double fMinClusterCharge
Definition: PixelLumiDQM.h:149
std::vector< bool > bunchTriggerMask
Definition: PixelLumiDQM.h:166
Definition: DetId.h:18
static double CM2_TO_NANOBARN
Definition: PixelLumiDQM.h:52
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:276
LuminosityBlockAuxiliary const & luminosityBlockAuxiliary() const
bool fIncludePixelClusterInfo
Definition: PixelLumiDQM.h:133
static double XSEC_PIXEL_CLUSTER
Definition: PixelLumiDQM.h:46
MonitorElement * book2D(Args &&...args)
Definition: DQMStore.h:133
MonitorElement * fHistBunchCrossings
Definition: PixelLumiDQM.h:157
std::map< std::string, MonitorElement * > fHistContainerThisRun
Definition: PixelLumiDQM.h:137
double getBinError(int binx) const
get uncertainty on content of bin (1-D) - See TH1::GetBinError for details
MonitorElement * fHistBunchCrossingsLastLumi
Definition: PixelLumiDQM.h:158
const T & get() const
Definition: EventSetup.h:56
virtual const PixelTopology & specificTopology() const
Returns a reference to the pixel proxy topology.
int layerName() const
layer id
UInt_t fTimestamp
Definition: PixelLumiDQM.h:121
string const
Definition: compareJSON.py:14
bool isTrackerPixel(const GeomDetEnumerators::SubDetector m)
std::string fLogFileName_
Definition: PixelLumiDQM.h:170
PixelLumiDQM(const edm::ParameterSet &)
Definition: PixelLumiDQM.cc:48
static double rXSEC_PIXEL_CLUSTER_UNC
Definition: PixelLumiDQM.h:51
virtual void dqmBeginRun(edm::Run const &, edm::EventSetup const &)
double getBinContent(int binx) const
get content of bin (1-D)
edm::EventID id() const
Definition: EventBase.h:59
Pixel cluster – collection of neighboring pixels above threshold.
std::map< int, PixelClusterCount > fNumPixelClusters
Definition: PixelLumiDQM.h:130
MonitorElement * fHistRecordedByBxLastLumi
Definition: PixelLumiDQM.h:163
int fMinPixelsPerCluster
Definition: PixelLumiDQM.h:146
int getNbinsX(void) const
get # of bins in X-axis
virtual void beginLuminosityBlock(edm::LuminosityBlock const &, edm::EventSetup const &)
tuple cout
Definition: gather_cfg.py:145
volatile std::atomic< bool > shutdown_flag false
constexpr double pi()
Definition: Pi.h:31
T x() const
Definition: PV3DBase.h:62
static size_t const kNumDisks
Definition: PixelLumiDQM.h:77
void Reset(void)
reset ME (ie. contents, errors, etc)
MonitorElement * fHistnFPClusVsLS[2]
Definition: PixelLumiDQM.h:155
tuple size
Write out results.
edm::Timestamp time() const
Definition: EventBase.h:60
const_iterator begin(bool update=false) const
tuple folder
Histograms Source for live online DQM in P5
Definition: Run.h:43
std::ofstream logFile_
Definition: PixelLumiDQM.h:172
static size_t const kNumLayers
Definition: PixelLumiDQM.h:76