CMS 3D CMS Logo

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

#include <L1TStage2uGT.h>

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

Public Member Functions

 L1TStage2uGT (const edm::ParameterSet &ps)
 
 ~L1TStage2uGT () 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 &, const edm::EventSetup &) override
 
void bookHistograms (DQMStore::IBooker &ibooker, edm::Run const &, edm::EventSetup const &) override
 
void dqmBeginRun (const edm::Run &, const edm::EventSetup &) override
 

Private Attributes

MonitorElementalgoBits_after_prescale_
 
MonitorElementalgoBits_after_prescale_bx_global_
 
MonitorElementalgoBits_after_prescale_bx_inEvt_
 
MonitorElementalgoBits_after_prescale_corr_
 
MonitorElementalgoBits_after_prescale_lumi_
 
MonitorElementalgoBits_before_bxmask_
 
MonitorElementalgoBits_before_bxmask_bx_global_
 
MonitorElementalgoBits_before_bxmask_bx_inEvt_
 
MonitorElementalgoBits_before_bxmask_corr_
 
MonitorElementalgoBits_before_bxmask_lumi_
 
MonitorElementalgoBits_before_prescale_
 
MonitorElementalgoBits_before_prescale_bx_global_
 
MonitorElementalgoBits_before_prescale_bx_inEvt_
 
MonitorElementalgoBits_before_prescale_corr_
 
MonitorElementalgoBits_before_prescale_lumi_
 
std::shared_ptr< l1t::L1TGlobalUtilgtUtil_
 
edm::EDGetTokenT< GlobalAlgBlkBxCollectionl1tStage2uGtSource_
 
std::string monitorDir_
 
int numAlgs_
 
MonitorElementprescaleFactorSet_
 
bool verbose_
 

Detailed Description

Description: DQM for L1 Micro Global Trigger.

Author
Mateusz Zarucki 2016
J. Berryhill, I. Mikulec
Vasile Mihai Ghete - HEPHY Vienna

Definition at line 44 of file L1TStage2uGT.h.

Constructor & Destructor Documentation

L1TStage2uGT::L1TStage2uGT ( const edm::ParameterSet ps)

Definition at line 15 of file L1TStage2uGT.cc.

15  :
16  l1tStage2uGtSource_(consumes<GlobalAlgBlkBxCollection>(params.getParameter<edm::InputTag>("l1tStage2uGtSource"))),
17  monitorDir_(params.getUntrackedParameter<std::string> ("monitorDir", "")),
18  verbose_(params.getUntrackedParameter<bool>("verbose", false)),
19  gtUtil_(new l1t::L1TGlobalUtil(params, consumesCollector(), *this, params.getParameter<edm::InputTag>("l1tStage2uGtSource"), params.getParameter<edm::InputTag>("l1tStage2uGtSource"))),
20  numAlgs_(0)
21 {
22 }
std::shared_ptr< l1t::L1TGlobalUtil > gtUtil_
Definition: L1TStage2uGT.h:65
edm::EDGetTokenT< GlobalAlgBlkBxCollection > l1tStage2uGtSource_
Definition: L1TStage2uGT.h:58
std::string monitorDir_
Definition: L1TStage2uGT.h:60
L1TStage2uGT::~L1TStage2uGT ( )
override

Definition at line 25 of file L1TStage2uGT.cc.

25 {}

Member Function Documentation

void L1TStage2uGT::analyze ( const edm::Event evt,
const edm::EventSetup evtSetup 
)
overrideprotected

Definition at line 114 of file L1TStage2uGT.cc.

References algoBits_after_prescale_, algoBits_after_prescale_bx_global_, algoBits_after_prescale_bx_inEvt_, algoBits_after_prescale_corr_, algoBits_after_prescale_lumi_, algoBits_before_bxmask_, algoBits_before_bxmask_bx_global_, algoBits_before_bxmask_bx_inEvt_, algoBits_before_bxmask_corr_, algoBits_before_bxmask_lumi_, algoBits_before_prescale_, algoBits_before_prescale_bx_global_, algoBits_before_prescale_bx_inEvt_, algoBits_before_prescale_corr_, algoBits_before_prescale_lumi_, BXVector< T >::begin(), edm::EventBase::bunchCrossing(), BXVector< T >::end(), MonitorElement::Fill(), edm::Event::getByToken(), BXVector< T >::getFirstBX(), BXVector< T >::getLastBX(), edm::HandleBase::isValid(), l1tStage2uGtSource_, edm::EventBase::luminosityBlock(), numAlgs_, prescaleFactorSet_, and verbose_.

114  {
115  if (verbose_) {
116  edm::LogInfo("L1TStage2uGT") << "L1TStage2uGT DQM: Analyzing.." << std::endl;
117  }
118 
119  // Get standard event parameters
120  int lumi = evt.luminosityBlock();
121  int bx = evt.bunchCrossing();
122 
123  // Open uGT readout record
125  evt.getByToken(l1tStage2uGtSource_, uGtAlgs);
126 
127  if (!uGtAlgs.isValid()) {
128  edm::LogInfo("L1TStage2uGT") << "Cannot find uGT readout record.";
129  return;
130  }
131 
132  //algoBits_->Fill(-1.); // fill underflow to normalize // FIXME: needed?
133  for (int ibx = uGtAlgs->getFirstBX(); ibx <= uGtAlgs->getLastBX(); ++ibx) {
134  for (auto itr = uGtAlgs->begin(ibx); itr != uGtAlgs->end(ibx); ++itr) { // FIXME: redundant loop?
135 
136  // Fills prescale factor set histogram
137  prescaleFactorSet_->Fill(lumi, itr->getPreScColumn());
138 
139  // Fills algorithm bits histograms
140  for(int algoBit = 0; algoBit < numAlgs_; ++algoBit) {
141 
142  // Algorithm bits before AlgoBX mask
143  if(itr->getAlgoDecisionInitial(algoBit)) {
144  algoBits_before_bxmask_->Fill(algoBit);
145  algoBits_before_bxmask_lumi_->Fill(lumi, algoBit);
146  algoBits_before_bxmask_bx_inEvt_->Fill(ibx, algoBit); // FIXME: or itr->getbxInEventNr()/getbxNr()?
147  algoBits_before_bxmask_bx_global_->Fill(bx + ibx, algoBit);
148 
149  for(int algoBit2 = 0; algoBit2 < numAlgs_; ++algoBit2) {
150  if(itr->getAlgoDecisionInitial(algoBit2)) {
151  algoBits_before_bxmask_corr_->Fill(algoBit, algoBit2);
152  }
153  }
154  }
155 
156  // Algorithm bits before prescale
157  if(itr->getAlgoDecisionInterm(algoBit)) {
159  algoBits_before_prescale_lumi_->Fill(lumi, algoBit);
161  algoBits_before_prescale_bx_global_->Fill(bx + ibx, algoBit);
162 
163  for(int algoBit2 = 0; algoBit2 < numAlgs_; ++algoBit2) {
164  if(itr->getAlgoDecisionInterm(algoBit2)) {
165  algoBits_before_prescale_corr_->Fill(algoBit, algoBit2);
166  }
167  }
168  }
169 
170  // Algorithm bits after prescale
171  if(itr->getAlgoDecisionFinal(algoBit)) {
172  algoBits_after_prescale_->Fill(algoBit);
173  algoBits_after_prescale_lumi_->Fill(lumi, algoBit);
175  algoBits_after_prescale_bx_global_->Fill(bx + ibx, algoBit);
176 
177  for(int algoBit2 = 0; algoBit2 < numAlgs_; ++algoBit2) {
178  if(itr->getAlgoDecisionFinal(algoBit2)) {
179  algoBits_after_prescale_corr_->Fill(algoBit, algoBit2);
180  }
181  }
182  }
183  }
184  }
185  }
186 }
const_iterator end(int bx) const
MonitorElement * algoBits_after_prescale_bx_global_
Definition: L1TStage2uGT.h:81
MonitorElement * algoBits_before_prescale_corr_
Definition: L1TStage2uGT.h:75
MonitorElement * algoBits_before_prescale_bx_global_
Definition: L1TStage2uGT.h:80
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
int bunchCrossing() const
Definition: EventBase.h:64
edm::LuminosityBlockNumber_t luminosityBlock() const
Definition: EventBase.h:61
MonitorElement * algoBits_before_bxmask_bx_global_
Definition: L1TStage2uGT.h:79
MonitorElement * algoBits_before_prescale_lumi_
Definition: L1TStage2uGT.h:90
void Fill(long long x)
MonitorElement * algoBits_before_bxmask_corr_
Definition: L1TStage2uGT.h:74
MonitorElement * algoBits_after_prescale_corr_
Definition: L1TStage2uGT.h:76
MonitorElement * algoBits_after_prescale_bx_inEvt_
Definition: L1TStage2uGT.h:86
MonitorElement * algoBits_before_bxmask_lumi_
Definition: L1TStage2uGT.h:89
MonitorElement * algoBits_before_prescale_
Definition: L1TStage2uGT.h:70
bool isValid() const
Definition: HandleBase.h:74
MonitorElement * prescaleFactorSet_
Definition: L1TStage2uGT.h:94
int getFirstBX() const
edm::EDGetTokenT< GlobalAlgBlkBxCollection > l1tStage2uGtSource_
Definition: L1TStage2uGT.h:58
MonitorElement * algoBits_before_prescale_bx_inEvt_
Definition: L1TStage2uGT.h:85
int getLastBX() const
MonitorElement * algoBits_before_bxmask_
Definition: L1TStage2uGT.h:69
MonitorElement * algoBits_after_prescale_
Definition: L1TStage2uGT.h:71
const_iterator begin(int bx) const
MonitorElement * algoBits_before_bxmask_bx_inEvt_
Definition: L1TStage2uGT.h:84
MonitorElement * algoBits_after_prescale_lumi_
Definition: L1TStage2uGT.h:91
void L1TStage2uGT::bookHistograms ( DQMStore::IBooker ibooker,
edm::Run const &  ,
edm::EventSetup const &  evtSetup 
)
overrideprotected

Definition at line 35 of file L1TStage2uGT.cc.

References algoBits_after_prescale_, algoBits_after_prescale_bx_global_, algoBits_after_prescale_bx_inEvt_, algoBits_after_prescale_corr_, algoBits_after_prescale_lumi_, algoBits_before_bxmask_, algoBits_before_bxmask_bx_global_, algoBits_before_bxmask_bx_inEvt_, algoBits_before_bxmask_corr_, algoBits_before_bxmask_lumi_, algoBits_before_prescale_, algoBits_before_prescale_bx_global_, algoBits_before_prescale_bx_inEvt_, algoBits_before_prescale_corr_, algoBits_before_prescale_lumi_, DQMStore::IBooker::book1D(), DQMStore::IBooker::book2D(), monitorDir_, numAlgs_, prescaleFactorSet_, MonitorElement::setAxisTitle(), and DQMStore::IBooker::setCurrentFolder().

35  {
36 
37  // Book histograms
38  const int numLS = 2000;
39  const double numLS_d = static_cast<double>(numLS);
40  const double numAlgs_d = static_cast<double>(numAlgs_);
41  const int numBx = 3564;
42  const double numBx_d = static_cast<double>(numBx);
43 
45 
46  // Algorithm bits
47  algoBits_before_bxmask_ = ibooker.book1D("algoBits_before_bxmask", "uGT: Algorithm Trigger Bits (before AlgoBX mask)", numAlgs_, -0.5, numAlgs_d-0.5);
48  algoBits_before_bxmask_->setAxisTitle("Algorithm Trigger Bits (before AlgoBX mask)", 1);
49 
50  algoBits_before_prescale_ = ibooker.book1D("algoBits_before_prescale", "uGT: Algorithm Trigger Bits (before prescale)", numAlgs_, -0.5, numAlgs_d-0.5);
51  algoBits_before_prescale_->setAxisTitle("Algorithm Trigger Bits (before prescale)", 1);
52 
53  algoBits_after_prescale_ = ibooker.book1D("algoBits_after_prescale", "uGT: Algorithm Trigger Bits (after prescale)", numAlgs_, -0.5, numAlgs_d-0.5);
54  algoBits_after_prescale_->setAxisTitle("Algorithm Trigger Bits (after prescale)", 1);
55 
56  // Algorithm bits correlation
57  algoBits_before_bxmask_corr_ = ibooker.book2D("algoBits_before_bxmask_corr","uGT: Algorithm Trigger Bit Correlation (before AlgoBX mask)", numAlgs_, -0.5, numAlgs_d-0.5, numAlgs_, -0.5, numAlgs_d-0.5);
58  algoBits_before_bxmask_corr_->setAxisTitle("Algorithm Trigger Bits (before AlgoBX mask)", 1);
59  algoBits_before_bxmask_corr_->setAxisTitle("Algorithm Trigger Bits (before AlgoBX mask)", 2);
60 
61  algoBits_before_prescale_corr_ = ibooker.book2D("algoBits_before_prescale_corr","uGT: Algorithm Trigger Bit Correlation (before prescale)", numAlgs_, -0.5, numAlgs_d-0.5, numAlgs_, -0.5, numAlgs_d-0.5);
62  algoBits_before_prescale_corr_->setAxisTitle("Algorithm Trigger Bits (before prescale)", 1);
63  algoBits_before_prescale_corr_->setAxisTitle("Algorithm Trigger Bits (before prescale)", 2);
64 
65  algoBits_after_prescale_corr_ = ibooker.book2D("algoBits_after_prescale_corr","uGT: Algorithm Trigger Bit Correlation (after prescale)", numAlgs_, -0.5, numAlgs_d-0.5, numAlgs_, -0.5, numAlgs_d-0.5);
66  algoBits_after_prescale_corr_->setAxisTitle("Algorithm Trigger Bits (after prescale)", 1);
67  algoBits_after_prescale_corr_->setAxisTitle("Algorithm Trigger Bits (after prescale)", 2);
68 
69  // Algorithm bits vs global BX number
70  algoBits_before_bxmask_bx_global_ = ibooker.book2D("algoBits_before_bxmask_bx_global", "uGT: Algorithm Trigger Bits (before AlgoBX mask) vs. Global BX Number", numBx, 0.5, numBx_d + 0.5, numAlgs_, -0.5, numAlgs_d-0.5);
71  algoBits_before_bxmask_bx_global_->setAxisTitle("Global Bunch Crossing Number", 1);
72  algoBits_before_bxmask_bx_global_->setAxisTitle("Algorithm Trigger Bits (before AlgoBX mask)", 2);
73 
74  algoBits_before_prescale_bx_global_ = ibooker.book2D("algoBits_before_prescale_bx_global", "uGT: Algorithm Trigger Bits (before prescale) vs. Global BX Number", numBx, 0.5, numBx_d + 0.5, numAlgs_, -0.5, numAlgs_d-0.5);
75  algoBits_before_prescale_bx_global_->setAxisTitle("Global Bunch Crossing Number", 1);
76  algoBits_before_prescale_bx_global_->setAxisTitle("Algorithm Trigger Bits (before prescale)", 2);
77 
78  algoBits_after_prescale_bx_global_ = ibooker.book2D("algoBits_after_prescale_bx_global", "uGT: Algorithm Trigger Bits (after prescale) vs. Global BX Number", numBx, 0.5, numBx_d + 0.5, numAlgs_, -0.5, numAlgs_d-0.5);
79  algoBits_after_prescale_bx_global_->setAxisTitle("Global Bunch Crossing Number", 1);
80  algoBits_after_prescale_bx_global_->setAxisTitle("Algorithm Trigger Bits (after prescale)", 2);
81 
82  // Algorithm bits vs BX number in event
83  algoBits_before_bxmask_bx_inEvt_ = ibooker.book2D("algoBits_before_bxmask_bx_inEvt", "uGT: Algorithm Trigger Bits (before AlgoBX mask) vs. BX Number in Event", 5, -2.5, 2.5, numAlgs_, -0.5, numAlgs_d-0.5);
84  algoBits_before_bxmask_bx_inEvt_->setAxisTitle("Bunch Crossing Number in Event", 1);
85  algoBits_before_bxmask_bx_inEvt_->setAxisTitle("Algorithm Trigger Bits (before AlgoBX mask)", 2);
86 
87  algoBits_before_prescale_bx_inEvt_ = ibooker.book2D("algoBits_before_prescale_bx_inEvt", "uGT: Algorithm Trigger Bits (before prescale) vs. BX Number in Event", 5, -2.5, 2.5, numAlgs_, -0.5, numAlgs_d-0.5);
88  algoBits_before_prescale_bx_inEvt_->setAxisTitle("Bunch Crossing Number in Event", 1);
89  algoBits_before_prescale_bx_inEvt_->setAxisTitle("Algorithm Trigger Bits (before prescale)", 2);
90 
91  algoBits_after_prescale_bx_inEvt_ = ibooker.book2D("algoBits_after_prescale_bx_inEvt", "uGT: Algorithm Trigger Bits (after prescale) vs. BX Number in Event", 5, -2.5, 2.5, numAlgs_, -0.5, numAlgs_d-0.5);
92  algoBits_after_prescale_bx_inEvt_->setAxisTitle("Bunch Crossing Number in Event", 1);
93  algoBits_after_prescale_bx_inEvt_->setAxisTitle("Algorithm Trigger Bits (after prescale)", 2);
94 
95  // Algorithm bits vs LS
96  algoBits_before_bxmask_lumi_ = ibooker.book2D("algoBits_before_bxmask_lumi","uGT: Algorithm Trigger Bits (before AlgoBX mask) vs. LS", numLS, 0., numLS_d, numAlgs_, -0.5, numAlgs_d-0.5);
97  algoBits_before_bxmask_lumi_->setAxisTitle("Luminosity Segment", 1);
98  algoBits_before_bxmask_lumi_->setAxisTitle("Algorithm Trigger Bits (before AlgoBX mask)", 2);
99 
100  algoBits_before_prescale_lumi_ = ibooker.book2D("algoBits_before_prescale_lumi","uGT: Algorithm Trigger Bits (before prescale) vs. LS", numLS, 0., numLS_d, numAlgs_, -0.5, numAlgs_d-0.5);
101  algoBits_before_prescale_lumi_->setAxisTitle("Luminosity Segment", 1);
102  algoBits_before_prescale_lumi_->setAxisTitle("Algorithm Trigger Bits (before prescale)", 2);
103 
104  algoBits_after_prescale_lumi_ = ibooker.book2D("algoBits_after_prescale_lumi","uGT: Algorithm Trigger Bits (after prescale) vs. LS", numLS, 0., numLS_d, numAlgs_, -0.5, numAlgs_d-0.5);
105  algoBits_after_prescale_lumi_->setAxisTitle("Luminosity Segment", 1);
106  algoBits_after_prescale_lumi_->setAxisTitle("Algorithm Trigger Bits (after prescale)", 2);
107 
108  // Prescale factor index
109  prescaleFactorSet_ = ibooker.book2D("prescaleFactorSet", "uGT: Index of Prescale Factor Set vs. LS", numLS, 0., numLS_d, 25, 0., 25.);
110  prescaleFactorSet_->setAxisTitle("Luminosity Segment", 1);
111  prescaleFactorSet_->setAxisTitle("Prescale Factor Set Index", 2);
112 }
MonitorElement * algoBits_after_prescale_bx_global_
Definition: L1TStage2uGT.h:81
MonitorElement * algoBits_before_prescale_corr_
Definition: L1TStage2uGT.h:75
MonitorElement * algoBits_before_prescale_bx_global_
Definition: L1TStage2uGT.h:80
MonitorElement * algoBits_before_bxmask_bx_global_
Definition: L1TStage2uGT.h:79
MonitorElement * algoBits_before_prescale_lumi_
Definition: L1TStage2uGT.h:90
MonitorElement * algoBits_before_bxmask_corr_
Definition: L1TStage2uGT.h:74
MonitorElement * algoBits_after_prescale_bx_inEvt_
Definition: L1TStage2uGT.h:86
MonitorElement * algoBits_after_prescale_corr_
Definition: L1TStage2uGT.h:76
MonitorElement * algoBits_before_bxmask_lumi_
Definition: L1TStage2uGT.h:89
void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:268
MonitorElement * algoBits_before_prescale_
Definition: L1TStage2uGT.h:70
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:106
MonitorElement * prescaleFactorSet_
Definition: L1TStage2uGT.h:94
MonitorElement * book2D(Args &&...args)
Definition: DQMStore.h:109
std::string monitorDir_
Definition: L1TStage2uGT.h:60
MonitorElement * algoBits_before_prescale_bx_inEvt_
Definition: L1TStage2uGT.h:85
MonitorElement * algoBits_before_bxmask_
Definition: L1TStage2uGT.h:69
MonitorElement * algoBits_after_prescale_
Definition: L1TStage2uGT.h:71
void setAxisTitle(const std::string &title, int axis=1)
set x-, y- or z-axis title (axis=1, 2, 3 respectively)
MonitorElement * algoBits_before_bxmask_bx_inEvt_
Definition: L1TStage2uGT.h:84
MonitorElement * algoBits_after_prescale_lumi_
Definition: L1TStage2uGT.h:91
void L1TStage2uGT::dqmBeginRun ( const edm::Run ,
const edm::EventSetup  
)
overrideprotected

Definition at line 27 of file L1TStage2uGT.cc.

References gtUtil_, and numAlgs_.

27  {
28  // Get the trigger menu information
29  gtUtil_->retrieveL1Setup(evtSetup);
30  // Find the number of algos defined
31  numAlgs_ = static_cast<int>(gtUtil_->decisionsInitial().size());
32 }
std::shared_ptr< l1t::L1TGlobalUtil > gtUtil_
Definition: L1TStage2uGT.h:65

Member Data Documentation

MonitorElement* L1TStage2uGT::algoBits_after_prescale_
private

Definition at line 71 of file L1TStage2uGT.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* L1TStage2uGT::algoBits_after_prescale_bx_global_
private

Definition at line 81 of file L1TStage2uGT.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* L1TStage2uGT::algoBits_after_prescale_bx_inEvt_
private

Definition at line 86 of file L1TStage2uGT.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* L1TStage2uGT::algoBits_after_prescale_corr_
private

Definition at line 76 of file L1TStage2uGT.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* L1TStage2uGT::algoBits_after_prescale_lumi_
private

Definition at line 91 of file L1TStage2uGT.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* L1TStage2uGT::algoBits_before_bxmask_
private

Definition at line 69 of file L1TStage2uGT.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* L1TStage2uGT::algoBits_before_bxmask_bx_global_
private

Definition at line 79 of file L1TStage2uGT.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* L1TStage2uGT::algoBits_before_bxmask_bx_inEvt_
private

Definition at line 84 of file L1TStage2uGT.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* L1TStage2uGT::algoBits_before_bxmask_corr_
private

Definition at line 74 of file L1TStage2uGT.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* L1TStage2uGT::algoBits_before_bxmask_lumi_
private

Definition at line 89 of file L1TStage2uGT.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* L1TStage2uGT::algoBits_before_prescale_
private

Definition at line 70 of file L1TStage2uGT.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* L1TStage2uGT::algoBits_before_prescale_bx_global_
private

Definition at line 80 of file L1TStage2uGT.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* L1TStage2uGT::algoBits_before_prescale_bx_inEvt_
private

Definition at line 85 of file L1TStage2uGT.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* L1TStage2uGT::algoBits_before_prescale_corr_
private

Definition at line 75 of file L1TStage2uGT.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* L1TStage2uGT::algoBits_before_prescale_lumi_
private

Definition at line 90 of file L1TStage2uGT.h.

Referenced by analyze(), and bookHistograms().

std::shared_ptr<l1t::L1TGlobalUtil> L1TStage2uGT::gtUtil_
private

Definition at line 65 of file L1TStage2uGT.h.

Referenced by dqmBeginRun().

edm::EDGetTokenT<GlobalAlgBlkBxCollection> L1TStage2uGT::l1tStage2uGtSource_
private

Definition at line 58 of file L1TStage2uGT.h.

Referenced by analyze().

std::string L1TStage2uGT::monitorDir_
private

Definition at line 60 of file L1TStage2uGT.h.

Referenced by bookHistograms().

int L1TStage2uGT::numAlgs_
private

Definition at line 66 of file L1TStage2uGT.h.

Referenced by analyze(), bookHistograms(), and dqmBeginRun().

MonitorElement* L1TStage2uGT::prescaleFactorSet_
private

Definition at line 94 of file L1TStage2uGT.h.

Referenced by analyze(), and bookHistograms().

bool L1TStage2uGT::verbose_
private

Definition at line 62 of file L1TStage2uGT.h.

Referenced by analyze().