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 
7 
9 
12 
15 
16 // TODO: move to configuration?
17 namespace {
18  const unsigned int R10BINS = 1024;
19  const float R10MIN = -0.5;
20  const float R10MAX = 1023.5;
21 
22  const unsigned int PUMETABINS = 22;
23  const unsigned int PUMNORMALIZE = 22;
24 
25  const unsigned int PUMBINS = 18;
26  const float PUMMIN = -0.5;
27  const float PUMMAX = 17.5;
28 }
29 
31  regionSource_(consumes<L1CaloRegionCollection>(ps.getParameter<edm::InputTag>("regionSource"))),
32  histFolder_(ps.getParameter<std::string>("histFolder"))
33 {
34 }
35 
37 {
38 }
39 
41 {
42 }
43 
45 {
46  edm::Handle<L1CaloRegionCollection> regionCollection;
47  event.getByToken(regionSource_, regionCollection);
48 
49  int nonzeroRegionsBxP2{0};
50  int nonzeroRegionsBx0{0};
51  int nonzeroRegionsBxM2{0};
52 
53  float regionsTotalEtBxP2{0.};
54  float regionsTotalEtBx0{0.};
55  float regionsTotalEtBxM2{0.};
56 
57  for (const auto& region : *regionCollection) {
58  if ( region.et() > 0 ) {
59  if ( region.bx() == 0 ) {
60  nonzeroRegionsBx0++;
61  regionsTotalEtBx0 += region.et();
62  }
63  else if ( region.bx() == 2 ) {
64  nonzeroRegionsBxP2++;
65  regionsTotalEtBxP2 += region.et();
66  }
67  else if ( region.bx() == -2 ) {
68  nonzeroRegionsBxM2++;
69  regionsTotalEtBxM2 += region.et();
70  }
71  }
72  }
73 
74  regionsTotalEtBxP2_->Fill(regionsTotalEtBxP2);
75  regionsTotalEtBx0_->Fill(regionsTotalEtBx0);
76  regionsTotalEtBxM2_->Fill(regionsTotalEtBxM2);
77 
78  regionsAvgEtBxP2_->Fill(regionsTotalEtBxP2/396.);
79  regionsAvgEtBx0_->Fill(regionsTotalEtBx0/396.);
80  regionsAvgEtBxM2_->Fill(regionsTotalEtBxM2/396.);
81 
82  regionsAvgNonZeroEtBxP2_->Fill(regionsTotalEtBxP2/nonzeroRegionsBxP2);
83  regionsAvgNonZeroEtBx0_->Fill(regionsTotalEtBx0/nonzeroRegionsBx0);
84  regionsAvgNonZeroEtBxM2_->Fill(regionsTotalEtBxM2/nonzeroRegionsBxM2);
85 
86  nonZeroRegionsBxP2_->Fill(nonzeroRegionsBxP2);
87  nonZeroRegionsBx0_->Fill(nonzeroRegionsBx0);
88  nonZeroRegionsBxM2_->Fill(nonzeroRegionsBxM2);
89 
90  for (const auto& region : *regionCollection) {
91  size_t etaBin = region.gctEta();
92  regionBxPopulation_->Fill(etaBin*18+region.gctPhi(), region.bx());
93  regionBxEtSum_->Fill(etaBin*18+region.gctPhi(), region.bx(), region.et());
94  if ( region.bx() == 0 )
95  regionsPUMEtaBx0_[etaBin]->Fill(nonzeroRegionsBx0/PUMNORMALIZE, region.et());
96  else if ( region.bx() == 2 )
97  regionsPUMEtaBxP2_[etaBin]->Fill(nonzeroRegionsBxP2/PUMNORMALIZE, region.et());
98  else if ( region.bx() == -2 )
99  regionsPUMEtaBxM2_[etaBin]->Fill(nonzeroRegionsBxM2/PUMNORMALIZE, region.et());
100  }
101 }
102 
104 {
105  ibooker.setCurrentFolder(histFolder_+"/BX0");
106  regionsPUMEtaBx0_.resize(PUMETABINS);
107  for (size_t ieta=0; ieta<PUMETABINS; ++ieta) {
108  regionsTotalEtBx0_ = ibooker.book1D("regionsTotalEt", "Total ET distribution;Sum Rank;Counts", 200, 0, 2000);
109  regionsAvgEtBx0_ = ibooker.book1D("regionsAvgEt", "Average Rank;Average Rank;Counts", R10BINS, R10MIN, R10MAX);
110  regionsAvgNonZeroEtBx0_ = ibooker.book1D("regionsAvgNonZeroEt", "Average Rank >0;Average Rank Regions>0;Counts", R10BINS, R10MIN, R10MAX);
111  nonZeroRegionsBx0_ = ibooker.book1D("nonZeroRegions", "Nonzero regions;Number Regions >0;Counts", 397, -0.5, 396.5);
112  regionsPUMEtaBx0_[ieta] = ibooker.book2D("regionsPUMEta"+std::to_string(ieta), "PUM Bin rank distribution;PU bin;Rank", PUMBINS, PUMMIN, PUMMAX, R10BINS, R10MIN, R10MAX);
113  }
114 
115  ibooker.setCurrentFolder(histFolder_+"/BXP2");
116  regionsPUMEtaBxP2_.resize(PUMETABINS);
117  for (size_t ieta=0; ieta<PUMETABINS; ++ieta) {
118  regionsTotalEtBxP2_ = ibooker.book1D("regionsTotalEt", "Total ET distribution;Sum Rank;Counts", 200, 0, 2000);
119  regionsAvgEtBxP2_ = ibooker.book1D("regionsAvgEt", "Average Rank;Average Rank;Counts", R10BINS, R10MIN, R10MAX);
120  regionsAvgNonZeroEtBxP2_ = ibooker.book1D("regionsAvgNonZeroEt", "Average Rank >0;Average Rank Regions>0;Counts", R10BINS, R10MIN, R10MAX);
121  nonZeroRegionsBxP2_ = ibooker.book1D("nonZeroRegions", "Nonzero regions;Number Regions >0;Counts", 397, -0.5, 396.5);
122  regionsPUMEtaBxP2_[ieta] = ibooker.book2D("regionsPUMEta"+std::to_string(ieta), "PUM Bin rank distribution;PU bin;Rank", PUMBINS, PUMMIN, PUMMAX, R10BINS, R10MIN, R10MAX);
123  }
124 
125  ibooker.setCurrentFolder(histFolder_+"/BXM2");
126  regionsPUMEtaBxM2_.resize(PUMETABINS);
127  for (size_t ieta=0; ieta<PUMETABINS; ++ieta) {
128  regionsTotalEtBxM2_ = ibooker.book1D("regionsTotalEt", "Total ET distribution;Sum Rank;Counts", 200, 0, 2000);
129  regionsAvgEtBxM2_ = ibooker.book1D("regionsAvgEt", "Average Rank;Average Rank;Counts", R10BINS, R10MIN, R10MAX);
130  regionsAvgNonZeroEtBxM2_ = ibooker.book1D("regionsAvgNonZeroEt", "Average Rank >0;Average Rank Regions>0;Counts", R10BINS, R10MIN, R10MAX);
131  nonZeroRegionsBxM2_ = ibooker.book1D("nonZeroRegions", "Nonzero regions;Number Regions >0;Counts", 397, -0.5, 396.5);
132  regionsPUMEtaBxM2_[ieta] = ibooker.book2D("regionsPUMEta"+std::to_string(ieta), "PUM Bin rank distribution;PU bin;Rank", PUMBINS, PUMMIN, PUMMAX, R10BINS, R10MIN, R10MAX);
133  }
134 
135  ibooker.setCurrentFolder(histFolder_+"/RegionBxInfo");
136  regionBxPopulation_ = ibooker.book2D("regionBxPopulation", "Event counts per region per bunch crossing;Region index (18*eta+phi);BX index;Counts", 396, -0.5, 395.5, 5, -2.5, 2.5);
137  regionBxEtSum_ = ibooker.book2D("regionBxEtSum", "Et per region per bunch crossing;Region index (18*eta+phi);BX index;Counts*et", 396, -0.5, 395.5, 5, -2.5, 2.5);
138 }
139 
141 {
142 }
143 
MonitorElement * regionsTotalEtBxP2_
Definition: L1TPUM.h:40
MonitorElement * regionsAvgNonZeroEtBxM2_
Definition: L1TPUM.h:50
MonitorElement * regionsAvgEtBxM2_
Definition: L1TPUM.h:46
L1TPUM(const edm::ParameterSet &ps)
Definition: L1TPUM.cc:30
const unsigned int R10BINS
Definition: L1TRCT.cc:27
edm::EDGetTokenT< L1CaloRegionCollection > regionSource_
Definition: L1TPUM.h:37
const float R10MIN
Definition: L1TRCT.cc:28
MonitorElement * regionsAvgEtBxP2_
Definition: L1TPUM.h:44
std::vector< MonitorElement * > regionsPUMEtaBxP2_
Definition: L1TPUM.h:59
MonitorElement * nonZeroRegionsBxM2_
Definition: L1TPUM.h:54
MonitorElement * regionBxPopulation_
Definition: L1TPUM.h:56
void Fill(long long x)
MonitorElement * nonZeroRegionsBx0_
Definition: L1TPUM.h:53
void beginLuminosityBlock(const edm::LuminosityBlock &, const edm::EventSetup &) override
Definition: L1TPUM.cc:140
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
void analyze(const edm::Event &e, const edm::EventSetup &c) override
Definition: L1TPUM.cc:44
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:118
MonitorElement * nonZeroRegionsBxP2_
Definition: L1TPUM.h:52
MonitorElement * regionsAvgNonZeroEtBxP2_
Definition: L1TPUM.h:48
def ls(path, rec=False)
Definition: eostools.py:348
std::vector< MonitorElement * > regionsPUMEtaBxM2_
Definition: L1TPUM.h:61
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:279
MonitorElement * book2D(Args &&...args)
Definition: DQMStore.h:136
std::vector< MonitorElement * > regionsPUMEtaBx0_
Definition: L1TPUM.h:60
std::string histFolder_
Definition: L1TPUM.h:38
void bookHistograms(DQMStore::IBooker &ibooker, const edm::Run &, const edm::EventSetup &) override
Definition: L1TPUM.cc:103
HLT enums.
MonitorElement * regionsAvgEtBx0_
Definition: L1TPUM.h:45
MonitorElement * regionsTotalEtBx0_
Definition: L1TPUM.h:41
MonitorElement * regionsAvgNonZeroEtBx0_
Definition: L1TPUM.h:49
std::vector< L1CaloRegion > L1CaloRegionCollection
const float R10MAX
Definition: L1TRCT.cc:29
~L1TPUM() override
Definition: L1TPUM.cc:36
MonitorElement * regionsTotalEtBxM2_
Definition: L1TPUM.h:42
MonitorElement * regionBxEtSum_
Definition: L1TPUM.h:57
void dqmBeginRun(const edm::Run &, const edm::EventSetup &) override
Definition: L1TPUM.cc:40
Definition: event.py:1
Definition: Run.h:43