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  edm::Handle<L1CaloRegionCollection> regionCollection;
37  event.getByToken(regionSource_, regionCollection);
38 
39  int nonzeroRegionsBxP2{0};
40  int nonzeroRegionsBx0{0};
41  int nonzeroRegionsBxM2{0};
42 
43  float regionsTotalEtBxP2{0.};
44  float regionsTotalEtBx0{0.};
45  float regionsTotalEtBxM2{0.};
46 
47  for (const auto& region : *regionCollection) {
48  if (region.et() > 0) {
49  if (region.bx() == 0) {
50  nonzeroRegionsBx0++;
51  regionsTotalEtBx0 += region.et();
52  } else if (region.bx() == 2) {
53  nonzeroRegionsBxP2++;
54  regionsTotalEtBxP2 += region.et();
55  } else if (region.bx() == -2) {
56  nonzeroRegionsBxM2++;
57  regionsTotalEtBxM2 += region.et();
58  }
59  }
60  }
61 
62  regionsTotalEtBxP2_->Fill(regionsTotalEtBxP2);
63  regionsTotalEtBx0_->Fill(regionsTotalEtBx0);
64  regionsTotalEtBxM2_->Fill(regionsTotalEtBxM2);
65 
66  regionsAvgEtBxP2_->Fill(regionsTotalEtBxP2 / 396.);
67  regionsAvgEtBx0_->Fill(regionsTotalEtBx0 / 396.);
68  regionsAvgEtBxM2_->Fill(regionsTotalEtBxM2 / 396.);
69 
70  regionsAvgNonZeroEtBxP2_->Fill(regionsTotalEtBxP2 / nonzeroRegionsBxP2);
71  regionsAvgNonZeroEtBx0_->Fill(regionsTotalEtBx0 / nonzeroRegionsBx0);
72  regionsAvgNonZeroEtBxM2_->Fill(regionsTotalEtBxM2 / nonzeroRegionsBxM2);
73 
74  nonZeroRegionsBxP2_->Fill(nonzeroRegionsBxP2);
75  nonZeroRegionsBx0_->Fill(nonzeroRegionsBx0);
76  nonZeroRegionsBxM2_->Fill(nonzeroRegionsBxM2);
77 
78  for (const auto& region : *regionCollection) {
79  size_t etaBin = region.gctEta();
80  regionBxPopulation_->Fill(etaBin * 18 + region.gctPhi(), region.bx());
81  regionBxEtSum_->Fill(etaBin * 18 + region.gctPhi(), region.bx(), region.et());
82  if (region.bx() == 0)
83  regionsPUMEtaBx0_[etaBin]->Fill(nonzeroRegionsBx0 / PUMNORMALIZE, region.et());
84  else if (region.bx() == 2)
85  regionsPUMEtaBxP2_[etaBin]->Fill(nonzeroRegionsBxP2 / PUMNORMALIZE, region.et());
86  else if (region.bx() == -2)
87  regionsPUMEtaBxM2_[etaBin]->Fill(nonzeroRegionsBxM2 / PUMNORMALIZE, region.et());
88  }
89 }
90 
92  ibooker.setCurrentFolder(histFolder_ + "/BX0");
93  regionsPUMEtaBx0_.resize(PUMETABINS);
94  for (size_t ieta = 0; ieta < PUMETABINS; ++ieta) {
95  regionsTotalEtBx0_ = ibooker.book1D("regionsTotalEt", "Total ET distribution;Sum Rank;Counts", 200, 0, 2000);
96  regionsAvgEtBx0_ = ibooker.book1D("regionsAvgEt", "Average Rank;Average Rank;Counts", R10BINS, R10MIN, R10MAX);
98  ibooker.book1D("regionsAvgNonZeroEt", "Average Rank >0;Average Rank Regions>0;Counts", R10BINS, R10MIN, R10MAX);
99  nonZeroRegionsBx0_ = ibooker.book1D("nonZeroRegions", "Nonzero regions;Number Regions >0;Counts", 397, -0.5, 396.5);
100  regionsPUMEtaBx0_[ieta] = ibooker.book2D("regionsPUMEta" + std::to_string(ieta),
101  "PUM Bin rank distribution;PU bin;Rank",
102  PUMBINS,
103  PUMMIN,
104  PUMMAX,
105  R10BINS,
106  R10MIN,
107  R10MAX);
108  }
109 
110  ibooker.setCurrentFolder(histFolder_ + "/BXP2");
111  regionsPUMEtaBxP2_.resize(PUMETABINS);
112  for (size_t ieta = 0; ieta < PUMETABINS; ++ieta) {
113  regionsTotalEtBxP2_ = ibooker.book1D("regionsTotalEt", "Total ET distribution;Sum Rank;Counts", 200, 0, 2000);
114  regionsAvgEtBxP2_ = ibooker.book1D("regionsAvgEt", "Average Rank;Average Rank;Counts", R10BINS, R10MIN, R10MAX);
116  ibooker.book1D("regionsAvgNonZeroEt", "Average Rank >0;Average Rank Regions>0;Counts", R10BINS, R10MIN, R10MAX);
118  ibooker.book1D("nonZeroRegions", "Nonzero regions;Number Regions >0;Counts", 397, -0.5, 396.5);
119  regionsPUMEtaBxP2_[ieta] = ibooker.book2D("regionsPUMEta" + std::to_string(ieta),
120  "PUM Bin rank distribution;PU bin;Rank",
121  PUMBINS,
122  PUMMIN,
123  PUMMAX,
124  R10BINS,
125  R10MIN,
126  R10MAX);
127  }
128 
129  ibooker.setCurrentFolder(histFolder_ + "/BXM2");
130  regionsPUMEtaBxM2_.resize(PUMETABINS);
131  for (size_t ieta = 0; ieta < PUMETABINS; ++ieta) {
132  regionsTotalEtBxM2_ = ibooker.book1D("regionsTotalEt", "Total ET distribution;Sum Rank;Counts", 200, 0, 2000);
133  regionsAvgEtBxM2_ = ibooker.book1D("regionsAvgEt", "Average Rank;Average Rank;Counts", R10BINS, R10MIN, R10MAX);
135  ibooker.book1D("regionsAvgNonZeroEt", "Average Rank >0;Average Rank Regions>0;Counts", R10BINS, R10MIN, R10MAX);
137  ibooker.book1D("nonZeroRegions", "Nonzero regions;Number Regions >0;Counts", 397, -0.5, 396.5);
138  regionsPUMEtaBxM2_[ieta] = ibooker.book2D("regionsPUMEta" + std::to_string(ieta),
139  "PUM Bin rank distribution;PU bin;Rank",
140  PUMBINS,
141  PUMMIN,
142  PUMMAX,
143  R10BINS,
144  R10MIN,
145  R10MAX);
146  }
147 
148  ibooker.setCurrentFolder(histFolder_ + "/RegionBxInfo");
150  ibooker.book2D("regionBxPopulation",
151  "Event counts per region per bunch crossing;Region index (18*eta+phi);BX index;Counts",
152  396,
153  -0.5,
154  395.5,
155  5,
156  -2.5,
157  2.5);
158  regionBxEtSum_ = ibooker.book2D("regionBxEtSum",
159  "Et per region per bunch crossing;Region index (18*eta+phi);BX index;Counts*et",
160  396,
161  -0.5,
162  395.5,
163  5,
164  -2.5,
165  2.5);
166 }
L1TPUM::regionsAvgEtBxP2_
MonitorElement * regionsAvgEtBxP2_
Definition: L1TPUM.h:42
L1TPUM::regionsAvgEtBx0_
MonitorElement * regionsAvgEtBx0_
Definition: L1TPUM.h:43
L1TPUM::~L1TPUM
~L1TPUM() override
Definition: L1TPUM.cc:33
RunSummaryRcd.h
L1TPUM::regionBxEtSum_
MonitorElement * regionBxEtSum_
Definition: L1TPUM.h:55
ESHandle.h
etaBin
int etaBin(const l1t::HGCalMulticluster *cl)
Definition: L1EGammaEEProducer.cc:19
edm::Run
Definition: Run.h:45
edm
HLT enums.
Definition: AlignableModifier.h:19
dqm::implementation::NavigatorBase::setCurrentFolder
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:32
L1TPUM::regionsAvgEtBxM2_
MonitorElement * regionsAvgEtBxM2_
Definition: L1TPUM.h:44
L1TPUM.h
L1TPUM::analyze
void analyze(const edm::Event &e, const edm::EventSetup &c) override
Definition: L1TPUM.cc:35
edm::Handle
Definition: AssociativeIterator.h:50
L1TPUM::nonZeroRegionsBxP2_
MonitorElement * nonZeroRegionsBxP2_
Definition: L1TPUM.h:50
L1TPUM::nonZeroRegionsBxM2_
MonitorElement * nonZeroRegionsBxM2_
Definition: L1TPUM.h:52
dqm::impl::MonitorElement::Fill
void Fill(long long x)
Definition: MonitorElement.h:290
L1TPUM::regionsAvgNonZeroEtBxP2_
MonitorElement * regionsAvgNonZeroEtBxP2_
Definition: L1TPUM.h:46
L1TPUM::regionsPUMEtaBxP2_
std::vector< MonitorElement * > regionsPUMEtaBxP2_
Definition: L1TPUM.h:57
L1TPUM::regionBxPopulation_
MonitorElement * regionBxPopulation_
Definition: L1TPUM.h:54
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
L1TPUM::histFolder_
std::string histFolder_
Definition: L1TPUM.h:36
LEDCalibrationChannels.ieta
ieta
Definition: LEDCalibrationChannels.py:63
HLT_2018_cff.InputTag
InputTag
Definition: HLT_2018_cff.py:79016
edm::ParameterSet
Definition: ParameterSet.h:36
L1TPUM::regionsPUMEtaBx0_
std::vector< MonitorElement * > regionsPUMEtaBx0_
Definition: L1TPUM.h:58
L1TPUM::nonZeroRegionsBx0_
MonitorElement * nonZeroRegionsBx0_
Definition: L1TPUM.h:51
L1TPUM::regionsTotalEtBx0_
MonitorElement * regionsTotalEtBx0_
Definition: L1TPUM.h:39
edm::EventSetup
Definition: EventSetup.h:57
L1TPUM::regionsAvgNonZeroEtBx0_
MonitorElement * regionsAvgNonZeroEtBx0_
Definition: L1TPUM.h:47
HcalObjRepresent::Fill
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
Definition: HcalObjRepresent.h:1053
L1TPUM::regionsTotalEtBxM2_
MonitorElement * regionsTotalEtBxM2_
Definition: L1TPUM.h:40
RunInfo.h
std
Definition: JetResolutionObject.h:76
writedatasetfile.run
run
Definition: writedatasetfile.py:27
L1CaloRegionCollection
std::vector< L1CaloRegion > L1CaloRegionCollection
Definition: L1CaloCollections.h:11
dqm::implementation::IBooker::book2D
MonitorElement * book2D(TString const &name, TString const &title, int nchX, double lowX, double highX, int nchY, double lowY, double highY, FUNC onbooking=NOOP())
Definition: DQMStore.h:177
EventAuxiliary.h
L1TPUM::regionsPUMEtaBxM2_
std::vector< MonitorElement * > regionsPUMEtaBxM2_
Definition: L1TPUM.h:59
HLT_2018_cff.region
region
Definition: HLT_2018_cff.py:81479
L1TPUM::regionsAvgNonZeroEtBxM2_
MonitorElement * regionsAvgNonZeroEtBxM2_
Definition: L1TPUM.h:48
L1TPUM::L1TPUM
L1TPUM(const edm::ParameterSet &ps)
Definition: L1TPUM.cc:29
dqm::implementation::IBooker
Definition: DQMStore.h:43
L1TPUM::regionsTotalEtBxP2_
MonitorElement * regionsTotalEtBxP2_
Definition: L1TPUM.h:38
L1TPUM::bookHistograms
void bookHistograms(DQMStore::IBooker &ibooker, const edm::Run &, const edm::EventSetup &) override
Definition: L1TPUM.cc:91
event
Definition: event.py:1
edm::Event
Definition: Event.h:73
L1TPUM::regionSource_
edm::EDGetTokenT< L1CaloRegionCollection > regionSource_
Definition: L1TPUM.h:35
dqm::implementation::IBooker::book1D
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98