CMS 3D CMS Logo

L1TPUM.cc
Go to the documentation of this file.
1 /*
2  * \file L1TPUM.cc
3  *
4  * N. Smith <nick.smith@cern.ch>
5  */
6 
8 
11 
14 
15 // TODO: move to configuration?
16 namespace {
17  const unsigned int R10BINS = 1024;
18  const float R10MIN = -0.5;
19  const float R10MAX = 1023.5;
20 
21  const unsigned int PUMETABINS = 22;
22  const unsigned int PUMNORMALIZE = 22;
23 
24  const unsigned int PUMBINS = 18;
25  const float PUMMIN = -0.5;
26  const float PUMMAX = 17.5;
27 } // namespace
28 
30  : regionSource_(consumes<L1CaloRegionCollection>(ps.getParameter<edm::InputTag>("regionSource"))),
31  histFolder_(ps.getParameter<std::string>("histFolder")) {}
32 
34 
36 
38  edm::Handle<L1CaloRegionCollection> regionCollection;
39  event.getByToken(regionSource_, regionCollection);
40 
41  int nonzeroRegionsBxP2{0};
42  int nonzeroRegionsBx0{0};
43  int nonzeroRegionsBxM2{0};
44 
45  float regionsTotalEtBxP2{0.};
46  float regionsTotalEtBx0{0.};
47  float regionsTotalEtBxM2{0.};
48 
49  for (const auto& region : *regionCollection) {
50  if (region.et() > 0) {
51  if (region.bx() == 0) {
52  nonzeroRegionsBx0++;
53  regionsTotalEtBx0 += region.et();
54  } else if (region.bx() == 2) {
55  nonzeroRegionsBxP2++;
56  regionsTotalEtBxP2 += region.et();
57  } else if (region.bx() == -2) {
58  nonzeroRegionsBxM2++;
59  regionsTotalEtBxM2 += region.et();
60  }
61  }
62  }
63 
64  regionsTotalEtBxP2_->Fill(regionsTotalEtBxP2);
65  regionsTotalEtBx0_->Fill(regionsTotalEtBx0);
66  regionsTotalEtBxM2_->Fill(regionsTotalEtBxM2);
67 
68  regionsAvgEtBxP2_->Fill(regionsTotalEtBxP2 / 396.);
69  regionsAvgEtBx0_->Fill(regionsTotalEtBx0 / 396.);
70  regionsAvgEtBxM2_->Fill(regionsTotalEtBxM2 / 396.);
71 
72  regionsAvgNonZeroEtBxP2_->Fill(regionsTotalEtBxP2 / nonzeroRegionsBxP2);
73  regionsAvgNonZeroEtBx0_->Fill(regionsTotalEtBx0 / nonzeroRegionsBx0);
74  regionsAvgNonZeroEtBxM2_->Fill(regionsTotalEtBxM2 / nonzeroRegionsBxM2);
75 
76  nonZeroRegionsBxP2_->Fill(nonzeroRegionsBxP2);
77  nonZeroRegionsBx0_->Fill(nonzeroRegionsBx0);
78  nonZeroRegionsBxM2_->Fill(nonzeroRegionsBxM2);
79 
80  for (const auto& region : *regionCollection) {
81  size_t etaBin = region.gctEta();
82  regionBxPopulation_->Fill(etaBin * 18 + region.gctPhi(), region.bx());
83  regionBxEtSum_->Fill(etaBin * 18 + region.gctPhi(), region.bx(), region.et());
84  if (region.bx() == 0)
85  regionsPUMEtaBx0_[etaBin]->Fill(nonzeroRegionsBx0 / PUMNORMALIZE, region.et());
86  else if (region.bx() == 2)
87  regionsPUMEtaBxP2_[etaBin]->Fill(nonzeroRegionsBxP2 / PUMNORMALIZE, region.et());
88  else if (region.bx() == -2)
89  regionsPUMEtaBxM2_[etaBin]->Fill(nonzeroRegionsBxM2 / PUMNORMALIZE, region.et());
90  }
91 }
92 
94  ibooker.setCurrentFolder(histFolder_ + "/BX0");
95  regionsPUMEtaBx0_.resize(PUMETABINS);
96  for (size_t ieta = 0; ieta < PUMETABINS; ++ieta) {
97  regionsTotalEtBx0_ = ibooker.book1D("regionsTotalEt", "Total ET distribution;Sum Rank;Counts", 200, 0, 2000);
98  regionsAvgEtBx0_ = ibooker.book1D("regionsAvgEt", "Average Rank;Average Rank;Counts", R10BINS, R10MIN, R10MAX);
100  ibooker.book1D("regionsAvgNonZeroEt", "Average Rank >0;Average Rank Regions>0;Counts", R10BINS, R10MIN, R10MAX);
101  nonZeroRegionsBx0_ = ibooker.book1D("nonZeroRegions", "Nonzero regions;Number Regions >0;Counts", 397, -0.5, 396.5);
102  regionsPUMEtaBx0_[ieta] = ibooker.book2D("regionsPUMEta" + std::to_string(ieta),
103  "PUM Bin rank distribution;PU bin;Rank",
104  PUMBINS,
105  PUMMIN,
106  PUMMAX,
107  R10BINS,
108  R10MIN,
109  R10MAX);
110  }
111 
112  ibooker.setCurrentFolder(histFolder_ + "/BXP2");
113  regionsPUMEtaBxP2_.resize(PUMETABINS);
114  for (size_t ieta = 0; ieta < PUMETABINS; ++ieta) {
115  regionsTotalEtBxP2_ = ibooker.book1D("regionsTotalEt", "Total ET distribution;Sum Rank;Counts", 200, 0, 2000);
116  regionsAvgEtBxP2_ = ibooker.book1D("regionsAvgEt", "Average Rank;Average Rank;Counts", R10BINS, R10MIN, R10MAX);
118  ibooker.book1D("regionsAvgNonZeroEt", "Average Rank >0;Average Rank Regions>0;Counts", R10BINS, R10MIN, R10MAX);
120  ibooker.book1D("nonZeroRegions", "Nonzero regions;Number Regions >0;Counts", 397, -0.5, 396.5);
121  regionsPUMEtaBxP2_[ieta] = ibooker.book2D("regionsPUMEta" + std::to_string(ieta),
122  "PUM Bin rank distribution;PU bin;Rank",
123  PUMBINS,
124  PUMMIN,
125  PUMMAX,
126  R10BINS,
127  R10MIN,
128  R10MAX);
129  }
130 
131  ibooker.setCurrentFolder(histFolder_ + "/BXM2");
132  regionsPUMEtaBxM2_.resize(PUMETABINS);
133  for (size_t ieta = 0; ieta < PUMETABINS; ++ieta) {
134  regionsTotalEtBxM2_ = ibooker.book1D("regionsTotalEt", "Total ET distribution;Sum Rank;Counts", 200, 0, 2000);
135  regionsAvgEtBxM2_ = ibooker.book1D("regionsAvgEt", "Average Rank;Average Rank;Counts", R10BINS, R10MIN, R10MAX);
137  ibooker.book1D("regionsAvgNonZeroEt", "Average Rank >0;Average Rank Regions>0;Counts", R10BINS, R10MIN, R10MAX);
139  ibooker.book1D("nonZeroRegions", "Nonzero regions;Number Regions >0;Counts", 397, -0.5, 396.5);
140  regionsPUMEtaBxM2_[ieta] = ibooker.book2D("regionsPUMEta" + std::to_string(ieta),
141  "PUM Bin rank distribution;PU bin;Rank",
142  PUMBINS,
143  PUMMIN,
144  PUMMAX,
145  R10BINS,
146  R10MIN,
147  R10MAX);
148  }
149 
150  ibooker.setCurrentFolder(histFolder_ + "/RegionBxInfo");
152  ibooker.book2D("regionBxPopulation",
153  "Event counts per region per bunch crossing;Region index (18*eta+phi);BX index;Counts",
154  396,
155  -0.5,
156  395.5,
157  5,
158  -2.5,
159  2.5);
160  regionBxEtSum_ = ibooker.book2D("regionBxEtSum",
161  "Et per region per bunch crossing;Region index (18*eta+phi);BX index;Counts*et",
162  396,
163  -0.5,
164  395.5,
165  5,
166  -2.5,
167  2.5);
168 }
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX)
Definition: DQMStore.cc:239
MonitorElement * regionsTotalEtBxP2_
Definition: L1TPUM.h:39
MonitorElement * regionsAvgNonZeroEtBxM2_
Definition: L1TPUM.h:49
MonitorElement * regionsAvgEtBxM2_
Definition: L1TPUM.h:45
L1TPUM(const edm::ParameterSet &ps)
Definition: L1TPUM.cc:29
edm::EDGetTokenT< L1CaloRegionCollection > regionSource_
Definition: L1TPUM.h:36
void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:418
MonitorElement * regionsAvgEtBxP2_
Definition: L1TPUM.h:43
std::vector< MonitorElement * > regionsPUMEtaBxP2_
Definition: L1TPUM.h:58
MonitorElement * nonZeroRegionsBxM2_
Definition: L1TPUM.h:53
MonitorElement * regionBxPopulation_
Definition: L1TPUM.h:55
void Fill(long long x)
MonitorElement * nonZeroRegionsBx0_
Definition: L1TPUM.h:52
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
void analyze(const edm::Event &e, const edm::EventSetup &c) override
Definition: L1TPUM.cc:37
MonitorElement * nonZeroRegionsBxP2_
Definition: L1TPUM.h:51
MonitorElement * regionsAvgNonZeroEtBxP2_
Definition: L1TPUM.h:47
std::vector< MonitorElement * > regionsPUMEtaBxM2_
Definition: L1TPUM.h:60
std::vector< MonitorElement * > regionsPUMEtaBx0_
Definition: L1TPUM.h:59
std::string histFolder_
Definition: L1TPUM.h:37
void bookHistograms(DQMStore::IBooker &ibooker, const edm::Run &, const edm::EventSetup &) override
Definition: L1TPUM.cc:93
HLT enums.
MonitorElement * book2D(TString const &name, TString const &title, int nchX, double lowX, double highX, int nchY, double lowY, double highY)
Definition: DQMStore.cc:266
MonitorElement * regionsAvgEtBx0_
Definition: L1TPUM.h:44
MonitorElement * regionsTotalEtBx0_
Definition: L1TPUM.h:40
MonitorElement * regionsAvgNonZeroEtBx0_
Definition: L1TPUM.h:48
std::vector< L1CaloRegion > L1CaloRegionCollection
~L1TPUM() override
Definition: L1TPUM.cc:33
MonitorElement * regionsTotalEtBxM2_
Definition: L1TPUM.h:41
MonitorElement * regionBxEtSum_
Definition: L1TPUM.h:56
void dqmBeginRun(const edm::Run &, const edm::EventSetup &) override
Definition: L1TPUM.cc:35
Definition: event.py:1
Definition: Run.h:45