CMS 3D CMS Logo

List of all members | Public Member Functions | Protected Member Functions | Private Attributes
L1TPUM Class Reference

#include <L1TPUM.h>

Inheritance diagram for L1TPUM:
one::DQMEDAnalyzer< T > one::dqmimplementation::DQMBaseClass< T... >

Public Member Functions

 L1TPUM (const edm::ParameterSet &ps)
 
 ~L1TPUM () override
 
- Public Member Functions inherited from one::DQMEDAnalyzer< T >
 DQMEDAnalyzer ()=default
 
 DQMEDAnalyzer (DQMEDAnalyzer< T... > const &)=delete
 
 DQMEDAnalyzer (DQMEDAnalyzer< T... > &&)=delete
 
 ~DQMEDAnalyzer () override=default
 

Protected Member Functions

void analyze (const edm::Event &e, const edm::EventSetup &c) override
 
void bookHistograms (DQMStore::IBooker &ibooker, const edm::Run &, const edm::EventSetup &) override
 
void dqmBeginRun (const edm::Run &, const edm::EventSetup &) override
 

Private Attributes

std::string histFolder_
 
MonitorElementnonZeroRegionsBx0_
 
MonitorElementnonZeroRegionsBxM2_
 
MonitorElementnonZeroRegionsBxP2_
 
MonitorElementregionBxEtSum_
 
MonitorElementregionBxPopulation_
 
MonitorElementregionsAvgEtBx0_
 
MonitorElementregionsAvgEtBxM2_
 
MonitorElementregionsAvgEtBxP2_
 
MonitorElementregionsAvgNonZeroEtBx0_
 
MonitorElementregionsAvgNonZeroEtBxM2_
 
MonitorElementregionsAvgNonZeroEtBxP2_
 
edm::EDGetTokenT< L1CaloRegionCollectionregionSource_
 
std::vector< MonitorElement * > regionsPUMEtaBx0_
 
std::vector< MonitorElement * > regionsPUMEtaBxM2_
 
std::vector< MonitorElement * > regionsPUMEtaBxP2_
 
MonitorElementregionsTotalEtBx0_
 
MonitorElementregionsTotalEtBxM2_
 
MonitorElementregionsTotalEtBxP2_
 

Detailed Description

Definition at line 25 of file L1TPUM.h.

Constructor & Destructor Documentation

L1TPUM::L1TPUM ( const edm::ParameterSet ps)

Definition at line 30 of file L1TPUM.cc.

30  :
31  regionSource_(consumes<L1CaloRegionCollection>(ps.getParameter<edm::InputTag>("regionSource"))),
32  histFolder_(ps.getParameter<std::string>("histFolder"))
33 {
34 }
T getParameter(std::string const &) const
edm::EDGetTokenT< L1CaloRegionCollection > regionSource_
Definition: L1TPUM.h:36
std::string histFolder_
Definition: L1TPUM.h:37
L1TPUM::~L1TPUM ( )
override

Definition at line 36 of file L1TPUM.cc.

37 {
38 }

Member Function Documentation

void L1TPUM::analyze ( const edm::Event e,
const edm::EventSetup c 
)
overrideprotected

Definition at line 44 of file L1TPUM.cc.

References conversionPostprocessing_cfi::etaBin, MonitorElement::Fill(), HcalObjRepresent::Fill(), nonZeroRegionsBx0_, nonZeroRegionsBxM2_, nonZeroRegionsBxP2_, regionBxEtSum_, regionBxPopulation_, regionsAvgEtBx0_, regionsAvgEtBxM2_, regionsAvgEtBxP2_, regionsAvgNonZeroEtBx0_, regionsAvgNonZeroEtBxM2_, regionsAvgNonZeroEtBxP2_, regionSource_, regionsPUMEtaBx0_, regionsPUMEtaBxM2_, regionsPUMEtaBxP2_, regionsTotalEtBx0_, regionsTotalEtBxM2_, and regionsTotalEtBxP2_.

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 }
MonitorElement * regionsTotalEtBxP2_
Definition: L1TPUM.h:39
MonitorElement * regionsAvgNonZeroEtBxM2_
Definition: L1TPUM.h:49
MonitorElement * regionsAvgEtBxM2_
Definition: L1TPUM.h:45
edm::EDGetTokenT< L1CaloRegionCollection > regionSource_
Definition: L1TPUM.h:36
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)
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
MonitorElement * regionsAvgEtBx0_
Definition: L1TPUM.h:44
MonitorElement * regionsTotalEtBx0_
Definition: L1TPUM.h:40
MonitorElement * regionsAvgNonZeroEtBx0_
Definition: L1TPUM.h:48
MonitorElement * regionsTotalEtBxM2_
Definition: L1TPUM.h:41
MonitorElement * regionBxEtSum_
Definition: L1TPUM.h:56
void L1TPUM::bookHistograms ( DQMStore::IBooker ibooker,
const edm::Run run,
const edm::EventSetup es 
)
overrideprotected

Definition at line 103 of file L1TPUM.cc.

References DQMStore::IBooker::book1D(), DQMStore::IBooker::book2D(), histFolder_, nonZeroRegionsBx0_, nonZeroRegionsBxM2_, nonZeroRegionsBxP2_, regionBxEtSum_, regionBxPopulation_, regionsAvgEtBx0_, regionsAvgEtBxM2_, regionsAvgEtBxP2_, regionsAvgNonZeroEtBx0_, regionsAvgNonZeroEtBxM2_, regionsAvgNonZeroEtBxP2_, regionsPUMEtaBx0_, regionsPUMEtaBxM2_, regionsPUMEtaBxP2_, regionsTotalEtBx0_, regionsTotalEtBxM2_, regionsTotalEtBxP2_, and DQMStore::IBooker::setCurrentFolder().

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 }
MonitorElement * regionsTotalEtBxP2_
Definition: L1TPUM.h:39
MonitorElement * regionsAvgNonZeroEtBxM2_
Definition: L1TPUM.h:49
MonitorElement * regionsAvgEtBxM2_
Definition: L1TPUM.h:45
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
MonitorElement * nonZeroRegionsBx0_
Definition: L1TPUM.h:52
void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:268
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:106
MonitorElement * nonZeroRegionsBxP2_
Definition: L1TPUM.h:51
MonitorElement * regionsAvgNonZeroEtBxP2_
Definition: L1TPUM.h:47
std::vector< MonitorElement * > regionsPUMEtaBxM2_
Definition: L1TPUM.h:60
MonitorElement * book2D(Args &&...args)
Definition: DQMStore.h:109
std::vector< MonitorElement * > regionsPUMEtaBx0_
Definition: L1TPUM.h:59
std::string histFolder_
Definition: L1TPUM.h:37
MonitorElement * regionsAvgEtBx0_
Definition: L1TPUM.h:44
MonitorElement * regionsTotalEtBx0_
Definition: L1TPUM.h:40
MonitorElement * regionsAvgNonZeroEtBx0_
Definition: L1TPUM.h:48
MonitorElement * regionsTotalEtBxM2_
Definition: L1TPUM.h:41
MonitorElement * regionBxEtSum_
Definition: L1TPUM.h:56
void L1TPUM::dqmBeginRun ( const edm::Run ,
const edm::EventSetup  
)
overrideprotected

Definition at line 40 of file L1TPUM.cc.

41 {
42 }

Member Data Documentation

std::string L1TPUM::histFolder_
private

Definition at line 37 of file L1TPUM.h.

Referenced by bookHistograms().

MonitorElement* L1TPUM::nonZeroRegionsBx0_
private

Definition at line 52 of file L1TPUM.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* L1TPUM::nonZeroRegionsBxM2_
private

Definition at line 53 of file L1TPUM.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* L1TPUM::nonZeroRegionsBxP2_
private

Definition at line 51 of file L1TPUM.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* L1TPUM::regionBxEtSum_
private

Definition at line 56 of file L1TPUM.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* L1TPUM::regionBxPopulation_
private

Definition at line 55 of file L1TPUM.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* L1TPUM::regionsAvgEtBx0_
private

Definition at line 44 of file L1TPUM.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* L1TPUM::regionsAvgEtBxM2_
private

Definition at line 45 of file L1TPUM.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* L1TPUM::regionsAvgEtBxP2_
private

Definition at line 43 of file L1TPUM.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* L1TPUM::regionsAvgNonZeroEtBx0_
private

Definition at line 48 of file L1TPUM.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* L1TPUM::regionsAvgNonZeroEtBxM2_
private

Definition at line 49 of file L1TPUM.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* L1TPUM::regionsAvgNonZeroEtBxP2_
private

Definition at line 47 of file L1TPUM.h.

Referenced by analyze(), and bookHistograms().

edm::EDGetTokenT<L1CaloRegionCollection> L1TPUM::regionSource_
private

Definition at line 36 of file L1TPUM.h.

Referenced by analyze().

std::vector<MonitorElement*> L1TPUM::regionsPUMEtaBx0_
private

Definition at line 59 of file L1TPUM.h.

Referenced by analyze(), and bookHistograms().

std::vector<MonitorElement*> L1TPUM::regionsPUMEtaBxM2_
private

Definition at line 60 of file L1TPUM.h.

Referenced by analyze(), and bookHistograms().

std::vector<MonitorElement*> L1TPUM::regionsPUMEtaBxP2_
private

Definition at line 58 of file L1TPUM.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* L1TPUM::regionsTotalEtBx0_
private

Definition at line 40 of file L1TPUM.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* L1TPUM::regionsTotalEtBxM2_
private

Definition at line 41 of file L1TPUM.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* L1TPUM::regionsTotalEtBxP2_
private

Definition at line 39 of file L1TPUM.h.

Referenced by analyze(), and bookHistograms().