CMS 3D CMS Logo

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

#include <L1TStage2uGTTiming.h>

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

Public Member Functions

 L1TStage2uGTTiming (const edm::ParameterSet &ps)
 
 ~L1TStage2uGTTiming () 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

int algoBitFirstBxInTrain_
 
int algoBitIsoBx_
 
int algoBitLastBxInTrain_
 
const std::string algoNameFirstBxInTrain_
 
const std::string algoNameIsoBx_
 
const std::string algoNameLastBxInTrain_
 
MonitorElementden_first_collision_in_train_
 
MonitorElementden_first_collision_in_train_minus1_
 
MonitorElementden_first_collision_in_train_minus2_
 
MonitorElementden_isolated_collision_
 
MonitorElementden_last_collision_in_train_
 
MonitorElementden_prescaled_algo_first_collision_in_train_
 
MonitorElementden_prescaled_algo_isolated_collision_
 
MonitorElementden_prescaled_algo_last_collision_in_train_
 
MonitorElementden_unprescaled_algo_first_collision_in_train_
 
MonitorElementden_unprescaled_algo_isolated_collision_
 
MonitorElementden_unprescaled_algo_last_collision_in_train_
 
MonitorElementfirst_collision_in_train_
 
MonitorElementfirst_collision_in_train_minus1_
 
MonitorElementfirst_collision_in_train_minus2_
 
std::shared_ptr< l1t::L1TGlobalUtilgtUtil_
 
MonitorElementisolated_collision_
 
edm::EDGetTokenT< GlobalAlgBlkBxCollectionl1tStage2uGtSource_
 
MonitorElementlast_collision_in_train_
 
std::string monitorDir_
 
int numAlgs_
 
MonitorElementprescaled_algo_first_collision_in_train_
 
MonitorElementprescaled_algo_isolated_collision_
 
MonitorElementprescaled_algo_last_collision_in_train_
 
std::vector< std::pair< std::string, int > > prescaledAlgoBitName_
 
std::vector< std::string > prescaledAlgoShortList_
 
MonitorElementunprescaled_algo_first_collision_in_train_
 
MonitorElementunprescaled_algo_isolated_collision_
 
MonitorElementunprescaled_algo_last_collision_in_train_
 
std::vector< std::pair< std::string, int > > unprescaledAlgoBitName_
 
std::vector< std::string > unprescaledAlgoShortList_
 
unsigned int useAlgoDecision_
 
bool verbose_
 

Detailed Description

Description: DQM for L1 Micro Global Trigger timing.

Definition at line 37 of file L1TStage2uGTTiming.h.

Constructor & Destructor Documentation

L1TStage2uGTTiming::L1TStage2uGTTiming ( const edm::ParameterSet ps)

Definition at line 11 of file L1TStage2uGTTiming.cc.

References edm::ParameterSet::getUntrackedParameter(), AlCaHLTBitMon_QueryRunRegistry::string, and useAlgoDecision_.

11  :
12  l1tStage2uGtSource_(consumes<GlobalAlgBlkBxCollection>(params.getParameter<edm::InputTag>("l1tStage2uGtSource"))),
13  monitorDir_(params.getUntrackedParameter<std::string> ("monitorDir", "")),
14  verbose_(params.getUntrackedParameter<bool>("verbose", false)),
15  gtUtil_(new l1t::L1TGlobalUtil(params, consumesCollector(), *this, params.getParameter<edm::InputTag>("l1tStage2uGtSource"), params.getParameter<edm::InputTag>("l1tStage2uGtSource"))),
16  numAlgs_(0),
19  algoBitIsoBx_(-1),
20  algoNameFirstBxInTrain_(params.getUntrackedParameter<std::string>("firstBXInTrainAlgo", "")),
21  algoNameLastBxInTrain_(params.getUntrackedParameter<std::string>("lastBXInTrainAlgo", "")),
22  algoNameIsoBx_(params.getUntrackedParameter<std::string>("isoBXAlgo", "")),
23  unprescaledAlgoShortList_(params.getUntrackedParameter<std::vector<std::string>> ("unprescaledAlgoShortList")),
24  prescaledAlgoShortList_(params.getUntrackedParameter<std::vector<std::string>> ("prescaledAlgoShortList"))
25 {
26  if (params.getUntrackedParameter<std::string>("useAlgoDecision").find("final") == 0) {
27  useAlgoDecision_ = 2;
28  } else if (params.getUntrackedParameter<std::string>("useAlgoDecision").find("intermediate") == 0) {
29  useAlgoDecision_ = 1;
30  } else {
31  useAlgoDecision_ = 0;
32  }
33 }
edm::EDGetTokenT< GlobalAlgBlkBxCollection > l1tStage2uGtSource_
const std::string algoNameLastBxInTrain_
unsigned int useAlgoDecision_
std::shared_ptr< l1t::L1TGlobalUtil > gtUtil_
std::vector< std::string > unprescaledAlgoShortList_
std::vector< std::string > prescaledAlgoShortList_
const std::string algoNameIsoBx_
const std::string algoNameFirstBxInTrain_
L1TStage2uGTTiming::~L1TStage2uGTTiming ( )
override

Definition at line 36 of file L1TStage2uGTTiming.cc.

36 {}

Member Function Documentation

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

Definition at line 167 of file L1TStage2uGTTiming.cc.

References patPFMETCorrections_cff::algo, algoBitFirstBxInTrain_, algoBitIsoBx_, algoBitLastBxInTrain_, BXVector< T >::begin(), den_first_collision_in_train_, den_first_collision_in_train_minus1_, den_first_collision_in_train_minus2_, den_isolated_collision_, den_last_collision_in_train_, den_prescaled_algo_first_collision_in_train_, den_prescaled_algo_isolated_collision_, den_prescaled_algo_last_collision_in_train_, den_unprescaled_algo_first_collision_in_train_, den_unprescaled_algo_isolated_collision_, den_unprescaled_algo_last_collision_in_train_, BXVector< T >::end(), MonitorElement::Fill(), first_collision_in_train_, first_collision_in_train_minus1_, first_collision_in_train_minus2_, edm::Event::getByToken(), BXVector< T >::getFirstBX(), BXVector< T >::getLastBX(), isolated_collision_, edm::HandleBase::isValid(), l1tStage2uGtSource_, last_collision_in_train_, SiStripPI::max, min(), prescaled_algo_first_collision_in_train_, prescaled_algo_isolated_collision_, prescaled_algo_last_collision_in_train_, prescaledAlgoBitName_, unprescaled_algo_first_collision_in_train_, unprescaled_algo_isolated_collision_, unprescaled_algo_last_collision_in_train_, unprescaledAlgoBitName_, useAlgoDecision_, and verbose_.

167  {
168  if (verbose_) {
169  edm::LogInfo("L1TStage2uGTTiming") << "L1TStage2uGTTiming DQM: Analyzing.." << std::endl;
170  }
171 
172  // Open uGT readout record
174  evt.getByToken(l1tStage2uGtSource_, uGtAlgs);
175 
176  if (!uGtAlgs.isValid()) {
177  edm::LogInfo("L1TStage2uGTTiming") << "Cannot find uGT readout record.";
178  return;
179  }
180 
181  // Find out in which BX the first collision in train, isolated bunch, and last collision in train have fired.
182  // In case of pre firing it will be in BX 1 or BX 2 and this will determine the BX shift that
183  // will be applied to the timing histogram later.
184  int bxShiftFirst = -999;
185  int bxShiftLast = -999;
186  int bxShiftIso = -999;
187  for (int bx = uGtAlgs->getFirstBX(); bx <= uGtAlgs->getLastBX(); ++bx) {
188  for (GlobalAlgBlkBxCollection::const_iterator itr = uGtAlgs->begin(bx); itr != uGtAlgs->end(bx); ++itr) {
189  // first bunch in train
190  if (algoBitFirstBxInTrain_ != -1) {
191  bool bit = false;
192  switch (useAlgoDecision_) {
193  case 0:
194  bit = itr->getAlgoDecisionInitial(algoBitFirstBxInTrain_);
195  break;
196  case 1:
197  bit = itr->getAlgoDecisionInterm(algoBitFirstBxInTrain_);
198  break;
199  case 2:
200  bit = itr->getAlgoDecisionFinal(algoBitFirstBxInTrain_);
201  break;
202  }
203  if (bit) {
204  bxShiftFirst = bx;
205  }
206  }
207  // last bunch in train
208  if(algoBitLastBxInTrain_ != -1) {
209  bool bit = false;
210  switch (useAlgoDecision_) {
211  case 0:
212  bit = itr->getAlgoDecisionInitial(algoBitLastBxInTrain_);
213  break;
214  case 1:
215  bit = itr->getAlgoDecisionInterm(algoBitLastBxInTrain_);
216  break;
217  case 2:
218  bit = itr->getAlgoDecisionFinal(algoBitLastBxInTrain_);
219  break;
220  }
221  if (bit) {
222  bxShiftLast = bx;
223  }
224  }
225  // isolated bunch
226  if (algoBitIsoBx_ != -1) {
227  bool bit = false;
228  switch (useAlgoDecision_) {
229  case 0:
230  bit = itr->getAlgoDecisionInitial(algoBitIsoBx_);
231  break;
232  case 1:
233  bit = itr->getAlgoDecisionInterm(algoBitIsoBx_);
234  break;
235  case 2:
236  bit = itr->getAlgoDecisionFinal(algoBitIsoBx_);
237  break;
238  }
239  if (bit) {
240  bxShiftIso = bx;
241  }
242  }
243  }
244  }
245 
246  // fill the first bunch in train maps
247  if (bxShiftFirst > -999) {
248  auto minBx = std::max(uGtAlgs->getFirstBX(), uGtAlgs->getFirstBX() + bxShiftFirst);
249  auto maxBx = std::min(uGtAlgs->getLastBX(), uGtAlgs->getLastBX() + bxShiftFirst);
250 
251  for (GlobalAlgBlkBxCollection::const_iterator itr = uGtAlgs->begin(bxShiftFirst); itr != uGtAlgs->end(bxShiftFirst); ++itr) {
252  for (int ibx = minBx; ibx <= maxBx; ++ibx) {
253  for (auto itr2 = uGtAlgs->begin(ibx); itr2 != uGtAlgs->end(ibx); ++itr2) {
254  auto algoBits = itr2->getAlgoDecisionInitial();
255  for (size_t algo = 0; algo < algoBits.size(); ++algo) {
256  if (algoBits.at(algo)) {
257  first_collision_in_train_->Fill(ibx - bxShiftFirst, algo);
258  for (int ibx2 = minBx; ibx2 <= maxBx; ++ibx2) {
259  den_first_collision_in_train_->Fill(ibx2 - bxShiftFirst, algo);
260  }
261  }
262  }
263  for (unsigned int algo = 0; algo < prescaledAlgoBitName_.size(); ++algo) {
264  if (itr2->getAlgoDecisionInitial(prescaledAlgoBitName_.at(algo).second)) {
265  prescaled_algo_first_collision_in_train_->Fill(ibx - bxShiftFirst, algo);
266  for (int ibx2 = minBx; ibx2 <= maxBx; ++ibx2) {
268  }
269  }
270  }
271  for (unsigned int algo = 0; algo < unprescaledAlgoBitName_.size(); ++algo) {
272  if (itr2->getAlgoDecisionInitial(unprescaledAlgoBitName_.at(algo).second)) {
274  for (int ibx2 = minBx; ibx2 <= maxBx; ++ibx2) {
276  }
277  }
278  }
279  }
280  }
281  }
282  }
283 
284  // fill the last bunch in train maps
285  if (bxShiftLast > -999) {
286  auto minBx = std::max(uGtAlgs->getFirstBX(), uGtAlgs->getFirstBX() + bxShiftLast);
287  auto maxBx = std::min(uGtAlgs->getLastBX(), uGtAlgs->getLastBX() + bxShiftLast);
288 
289  for (GlobalAlgBlkBxCollection::const_iterator itr = uGtAlgs->begin(bxShiftLast); itr != uGtAlgs->end(bxShiftLast); ++itr) {
290  for (int ibx = minBx; ibx <= maxBx; ++ibx) {
291  for (auto itr2 = uGtAlgs->begin(ibx); itr2 != uGtAlgs->end(ibx); ++itr2) {
292  auto algoBits = itr2->getAlgoDecisionInitial();
293  for (size_t algo = 0; algo < algoBits.size(); ++algo) {
294  if (algoBits.at(algo)) {
295  last_collision_in_train_->Fill(ibx - bxShiftLast, algo);
296  for (int ibx2 = minBx; ibx2 <= maxBx; ++ibx2) {
297  den_last_collision_in_train_->Fill(ibx2 - bxShiftLast, algo);
298  }
299  }
300  }
301  for (unsigned int algo = 0; algo < prescaledAlgoBitName_.size(); ++algo) {
302  if (itr2->getAlgoDecisionInitial(prescaledAlgoBitName_.at(algo).second)) {
304  for (int ibx2 = minBx; ibx2 <= maxBx; ++ibx2) {
306  }
307  }
308  }
309  for (unsigned int algo = 0; algo < unprescaledAlgoBitName_.size(); ++algo) {
310  if (itr2->getAlgoDecisionInitial(unprescaledAlgoBitName_.at(algo).second)) {
312  for (int ibx2 = minBx; ibx2 <= maxBx; ++ibx2) {
314  }
315  }
316  }
317  }
318  }
319  }
320  }
321 
322  // fill the isolated bunch maps
323  if (bxShiftIso > -999) {
324  auto minBx = std::max(uGtAlgs->getFirstBX(), uGtAlgs->getFirstBX() + bxShiftIso);
325  auto maxBx = std::min(uGtAlgs->getLastBX(), uGtAlgs->getLastBX() + bxShiftIso);
326 
327  for (GlobalAlgBlkBxCollection::const_iterator itr = uGtAlgs->begin(bxShiftIso); itr != uGtAlgs->end(bxShiftIso); ++itr) {
328  for (int ibx = minBx; ibx <= maxBx; ++ibx) {
329  for (auto itr2 = uGtAlgs->begin(ibx); itr2 != uGtAlgs->end(ibx); ++itr2) {
330  auto algoBits = itr2->getAlgoDecisionInitial();
331  for (size_t algo = 0; algo < algoBits.size(); ++algo) {
332  if (algoBits.at(algo)) {
333  isolated_collision_->Fill(ibx - bxShiftIso, algo);
334  for (int ibx2 = minBx; ibx2 <= maxBx; ++ibx2) {
335  den_isolated_collision_->Fill(ibx2 - bxShiftIso, algo);
336  }
337  }
338  }
339  for (unsigned int algo = 0; algo < prescaledAlgoBitName_.size(); ++algo) {
340  if (itr2->getAlgoDecisionInitial(prescaledAlgoBitName_.at(algo).second)) {
341  prescaled_algo_isolated_collision_->Fill(ibx - bxShiftIso, algo);
342  for (int ibx2 = minBx; ibx2 <= maxBx; ++ibx2) {
343  den_prescaled_algo_isolated_collision_->Fill(ibx2 - bxShiftIso, algo);
344  }
345  }
346  }
347  for (unsigned int algo = 0; algo < unprescaledAlgoBitName_.size(); ++algo) {
348  if (itr2->getAlgoDecisionInitial(unprescaledAlgoBitName_.at(algo).second)) {
349  unprescaled_algo_isolated_collision_->Fill(ibx - bxShiftIso, algo);
350  for (int ibx2 = minBx; ibx2 <= maxBx; ++ibx2) {
352  }
353  }
354  }
355  }
356  }
357  }
358  }
359 
360  // If algoBitFirstBxInTrain_ fired in L1A BX 2 something else must have prefired in the actual BX -2 before the first bunch crossing in the train
361  if (uGtAlgs->getLastBX() >= 2) {
362  for(auto itr = uGtAlgs->begin(2); itr != uGtAlgs->end(2); ++itr) {
363  if(algoBitFirstBxInTrain_ != -1 && itr->getAlgoDecisionInitial(algoBitFirstBxInTrain_)) {
364  for(int ibx = uGtAlgs->getFirstBX(); ibx <= uGtAlgs->getLastBX(); ++ibx) {
365  for(auto itr2 = uGtAlgs->begin(ibx); itr2 != uGtAlgs->end(ibx); ++itr2) {
366  auto algoBits = itr2->getAlgoDecisionInitial();
367  for(size_t algo = 0; algo < algoBits.size(); ++algo) {
368  if(algoBits.at(algo)) {
370  for(int ibx2 = uGtAlgs->getFirstBX(); ibx2 <= uGtAlgs->getLastBX(); ++ibx2) {
372  }
373  }
374  }
375  }
376  }
377  }
378  }
379  }
380 
381  // If algoBitFirstBxInTrain_ fired in L1A BX 1 something else must have prefired in the actual BX -1 before the first bunch crossing in the train
382  if (uGtAlgs->getLastBX() >= 1) {
383  for(auto itr = uGtAlgs->begin(1); itr != uGtAlgs->end(1); ++itr) {
384  if(algoBitFirstBxInTrain_ != -1 && itr->getAlgoDecisionInitial(algoBitFirstBxInTrain_)) {
385  for(int ibx = uGtAlgs->getFirstBX(); ibx <= uGtAlgs->getLastBX(); ++ibx) {
386  for(auto itr2 = uGtAlgs->begin(ibx); itr2 != uGtAlgs->end(ibx); ++itr2) {
387  auto algoBits = itr2->getAlgoDecisionInitial();
388  for(size_t algo = 0; algo < algoBits.size(); ++algo) {
389  if(algoBits.at(algo)) {
391  for(int ibx2 = uGtAlgs->getFirstBX(); ibx2 <= uGtAlgs->getLastBX(); ++ibx2) {
393  }
394  }
395  }
396  }
397  }
398  }
399  }
400  }
401 }
MonitorElement * prescaled_algo_first_collision_in_train_
const_iterator end(int bx) const
MonitorElement * den_unprescaled_algo_first_collision_in_train_
std::vector< std::pair< std::string, int > > prescaledAlgoBitName_
edm::EDGetTokenT< GlobalAlgBlkBxCollection > l1tStage2uGtSource_
MonitorElement * den_prescaled_algo_first_collision_in_train_
MonitorElement * den_prescaled_algo_last_collision_in_train_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
MonitorElement * unprescaled_algo_isolated_collision_
MonitorElement * isolated_collision_
MonitorElement * den_first_collision_in_train_
MonitorElement * den_unprescaled_algo_last_collision_in_train_
std::vector< std::pair< std::string, int > > unprescaledAlgoBitName_
MonitorElement * den_first_collision_in_train_minus2_
MonitorElement * last_collision_in_train_
MonitorElement * first_collision_in_train_minus2_
void Fill(long long x)
unsigned int useAlgoDecision_
MonitorElement * unprescaled_algo_last_collision_in_train_
MonitorElement * first_collision_in_train_minus1_
MonitorElement * prescaled_algo_isolated_collision_
MonitorElement * den_last_collision_in_train_
MonitorElement * den_first_collision_in_train_minus1_
T min(T a, T b)
Definition: MathUtil.h:58
bool isValid() const
Definition: HandleBase.h:74
int getFirstBX() const
MonitorElement * den_prescaled_algo_isolated_collision_
MonitorElement * den_isolated_collision_
MonitorElement * den_unprescaled_algo_isolated_collision_
MonitorElement * prescaled_algo_last_collision_in_train_
int getLastBX() const
MonitorElement * unprescaled_algo_first_collision_in_train_
const_iterator begin(int bx) const
MonitorElement * first_collision_in_train_
std::vector< T >::const_iterator const_iterator
Definition: BXVector.h:20
void L1TStage2uGTTiming::bookHistograms ( DQMStore::IBooker ibooker,
edm::Run const &  ,
edm::EventSetup const &  evtSetup 
)
overrideprotected

Definition at line 79 of file L1TStage2uGTTiming.cc.

References patPFMETCorrections_cff::algo, DQMStore::IBooker::book2D(), den_first_collision_in_train_, den_first_collision_in_train_minus1_, den_first_collision_in_train_minus2_, den_isolated_collision_, den_last_collision_in_train_, den_prescaled_algo_first_collision_in_train_, den_prescaled_algo_isolated_collision_, den_prescaled_algo_last_collision_in_train_, den_unprescaled_algo_first_collision_in_train_, den_unprescaled_algo_isolated_collision_, den_unprescaled_algo_last_collision_in_train_, first_collision_in_train_, first_collision_in_train_minus1_, first_collision_in_train_minus2_, isolated_collision_, last_collision_in_train_, monitorDir_, numAlgs_, prescaled_algo_first_collision_in_train_, prescaled_algo_isolated_collision_, prescaled_algo_last_collision_in_train_, prescaledAlgoBitName_, MonitorElement::setBinLabel(), DQMStore::IBooker::setCurrentFolder(), unprescaled_algo_first_collision_in_train_, unprescaled_algo_isolated_collision_, unprescaled_algo_last_collision_in_train_, and unprescaledAlgoBitName_.

79  {
80  // Book histograms
81  const auto numAlgs_d = static_cast<double>(numAlgs_);
82  const auto preAlgs_d = static_cast<double>(prescaledAlgoBitName_.size());
83  const auto unpreAlgs_d = static_cast<double>(unprescaledAlgoBitName_.size());
84 
86 
87  first_collision_in_train_minus2_ = ibooker.book2D("first_bunch_in_train_minus2", "uGT: Algorithm Trigger Bits (first bunch in train minus 2) vs. BX Number In Event;Bunch Crossing Number In Event;Algorithm Trigger Bits", 5, -2.5, 2.5, numAlgs_, -0.5, numAlgs_d-0.5);
88  den_first_collision_in_train_minus2_ = ibooker.book2D("den_first_bunch_in_train_minus2", "uGT: Algorithm Trigger Bits (all entries for each trigget bit first bunch in train minus 2) vs. BX Number In Event;Bunch Crossing Number In Event;Algorithm Trigger Bits", 5, -2.5, 2.5, numAlgs_, -0.5, numAlgs_d-0.5);
89 
90  first_collision_in_train_minus1_ = ibooker.book2D("first_bunch_in_train_minus1", "uGT: Algorithm Trigger Bits (first bunch in train minus 1) vs. BX Number In Event;Bunch Crossing Number In Event;Algorithm Trigger Bits", 5, -2.5, 2.5, numAlgs_, -0.5, numAlgs_d-0.5);
91  den_first_collision_in_train_minus1_ = ibooker.book2D("den_first_bunch_in_train_minus1", "uGT: Algorithm Trigger Bits (all entries for each trigget bit first bunch in train minus 1) vs. BX Number In Event;Bunch Crossing Number In Event;Algorithm Trigger Bits", 5, -2.5, 2.5, numAlgs_, -0.5, numAlgs_d-0.5);
92 
93  first_collision_in_train_ = ibooker.book2D("first_bunch_in_train", "uGT: Algorithm Trigger Bits (first bunch in train) vs. BX Number In Event;Bunch Crossing Number In Event;Algorithm Trigger Bits (first bunch in train)", 5, -2.5, 2.5, numAlgs_, -0.5, numAlgs_d-0.5);
94  den_first_collision_in_train_ = ibooker.book2D("den_first_bunch_in_train", "uGT: Algorithm Trigger Bits (all entries for each trigget bit first bunch in train) vs. BX Number In Event;Bunch Crossing Number In Event;Algorithm Trigger Bits (first bunch in train)", 5, -2.5, 2.5, numAlgs_, -0.5, numAlgs_d-0.5);
95 
96  last_collision_in_train_ = ibooker.book2D("last_bunch_in_train", "uGT: Algorithm Trigger Bits (last bunch in train) vs. BX Number In Event;Bunch Crossing Number In Event;Algorithm Trigger Bits (last bunch in train)", 5, -2.5, 2.5, numAlgs_, -0.5, numAlgs_d-0.5);
97  den_last_collision_in_train_ = ibooker.book2D("den_last_bunch_in_train", "uGT: Algorithm Trigger Bits (all entries for each trigget bit last bunch in train) vs. BX Number In Event;Bunch Crossing Number In Event;Algorithm Trigger Bits (last bunch in train)", 5, -2.5, 2.5, numAlgs_, -0.5, numAlgs_d-0.5);
98 
99  isolated_collision_ = ibooker.book2D("isolated_bunch", "uGT: Algorithm Trigger Bits vs. BX Number In Event;Bunch Crossing Number In Event;Algorithm Trigger Bits", 5, -2.5, 2.5, numAlgs_, -0.5, numAlgs_d-0.5);
100  den_isolated_collision_ = ibooker.book2D("den_isolated_bunch_in_train", "uGT: Algorithm Trigger Bits (all entries for each trigget bit isolated bunch in train) vs. BX Number In Event;Bunch Crossing Number In Event;Algorithm Trigger Bits (isolated bunch in train)", 5, -2.5, 2.5, numAlgs_, -0.5, numAlgs_d-0.5);
101 
102  // Prescaled and Unprescaled Algo Trigger Bits
103  // First bunch in train
104  prescaled_algo_first_collision_in_train_ = ibooker.book2D("prescaled_algo_first_collision_in_train_", "uGT: Prescaled Algorithm Trigger Bits vs. BX Number In Event;Bunch Crossing Number In Event;Algorithm Trigger Names + Bits", 5, -2.5, 2.5, prescaledAlgoBitName_.size(), -0.5, preAlgs_d-0.5);
105  for(unsigned int algo=0; algo<prescaledAlgoBitName_.size(); ++algo) {
107  }
108 
109  den_prescaled_algo_first_collision_in_train_ = ibooker.book2D("den_prescaled_algo_first_collision_in_train_", "uGT: Prescaled Algorithm Trigger Bits Deno vs. BX Number In Event;Bunch Crossing Number In Event;Algorithm Trigger Names + Bits", 5, -2.5, 2.5, prescaledAlgoBitName_.size(), -0.5, preAlgs_d-0.5);
110  for(unsigned int algo=0; algo<prescaledAlgoBitName_.size(); ++algo) {
112  }
113 
114  unprescaled_algo_first_collision_in_train_ = ibooker.book2D("unprescaled_algo_first_collision_in_train_", "uGT: Unprescaled Algorithm Trigger Bits vs. BX Number In Event;Bunch Crossing Number In Event;Algorithm Trigger Names + Bits", 5, -2.5, 2.5, unprescaledAlgoBitName_.size(), -0.5, unpreAlgs_d-0.5);
115  for(unsigned int algo=0; algo<unprescaledAlgoBitName_.size(); ++algo) {
117  }
118 
119  den_unprescaled_algo_first_collision_in_train_ = ibooker.book2D("den_unprescaled_algo_first_collision_in_train_", "uGT: Unprescaled Algorithm Trigger Bits Deno vs. BX Number In Event;Bunch Crossing Number In Event;Algorithm Trigger Names + Bits", 5, -2.5, 2.5, unprescaledAlgoBitName_.size(), -0.5, unpreAlgs_d-0.5);
120  for(unsigned int algo=0; algo<unprescaledAlgoBitName_.size(); ++algo) {
122  }
123 
124  // Isolated bunch
125  prescaled_algo_isolated_collision_ = ibooker.book2D("prescaled_algo_isolated_collision_", "uGT: Prescaled Algorithm Trigger Bits vs. BX Number In Event;Bunch Crossing Number In Event;Algorithm Trigger Names + Bits", 5, -2.5, 2.5, prescaledAlgoBitName_.size(), -0.5, preAlgs_d-0.5);
126  for(unsigned int algo=0; algo<prescaledAlgoBitName_.size(); ++algo) {
127  prescaled_algo_isolated_collision_->setBinLabel(algo+1, prescaledAlgoBitName_.at(algo).first+" ("+std::to_string(prescaledAlgoBitName_.at(algo).second)+")", 2);
128  }
129 
130  den_prescaled_algo_isolated_collision_ = ibooker.book2D("den_prescaled_algo_isolated_collision_", "uGT: Prescaled Algorithm Trigger Bits Deno vs. BX Number In Event;Bunch Crossing Number In Event;Algorithm Trigger Names + Bits", 5, -2.5, 2.5, prescaledAlgoBitName_.size(), -0.5, preAlgs_d-0.5);
131  for(unsigned int algo=0; algo<prescaledAlgoBitName_.size(); ++algo) {
133  }
134 
135  unprescaled_algo_isolated_collision_ = ibooker.book2D("unprescaled_algo_isolated_collision_", "uGT: Unprescaled Algorithm Trigger Bits vs. BX Number In Event;Bunch Crossing Number In Event;Algorithm Trigger Names + Bits", 5, -2.5, 2.5, unprescaledAlgoBitName_.size(), -0.5, unpreAlgs_d-0.5);
136  for(unsigned int algo=0; algo<unprescaledAlgoBitName_.size(); ++algo) {
138  }
139 
140  den_unprescaled_algo_isolated_collision_ = ibooker.book2D("den_unprescaled_algo_isolated_collision_", "uGT: Unprescaled Algorithm Trigger Bits Deno vs. BX Number In Event;Bunch Crossing Number In Event;Algorithm Trigger Names + Bits", 5, -2.5, 2.5, unprescaledAlgoBitName_.size(), -0.5, unpreAlgs_d-0.5);
141  for(unsigned int algo=0; algo<unprescaledAlgoBitName_.size(); ++algo) {
143  }
144 
145  // Last bunch in train
146  prescaled_algo_last_collision_in_train_ = ibooker.book2D("prescaled_algo_last_collision_in_train_", "uGT: Prescaled Algorithm Trigger Bits vs. BX Number In Event;Bunch Crossing Number In Event;Algorithm Trigger Names + Bits", 5, -2.5, 2.5, prescaledAlgoBitName_.size(), -0.5, preAlgs_d-0.5);
147  for(unsigned int algo=0; algo<prescaledAlgoBitName_.size(); ++algo) {
149  }
150 
151  den_prescaled_algo_last_collision_in_train_ = ibooker.book2D("den_prescaled_algo_last_collision_in_train_", "uGT: Prescaled Algorithm Trigger Bits Deno vs. BX Number In Event;Bunch Crossing Number In Event;Algorithm Trigger Names + Bits", 5, -2.5, 2.5, prescaledAlgoBitName_.size(), -0.5, preAlgs_d-0.5);
152  for(unsigned int algo=0; algo<prescaledAlgoBitName_.size(); ++algo) {
154  }
155 
156  unprescaled_algo_last_collision_in_train_ = ibooker.book2D("unprescaled_algo_last_collision_in_train_", "uGT: Unprescaled Algorithm Trigger Bits vs. BX Number In Event;Bunch Crossing Number In Event;Algorithm Trigger Names + Bits", 5, -2.5, 2.5, unprescaledAlgoBitName_.size(), -0.5, unpreAlgs_d-0.5);
157  for(unsigned int algo=0; algo<unprescaledAlgoBitName_.size(); ++algo) {
159  }
160 
161  den_unprescaled_algo_last_collision_in_train_ = ibooker.book2D("den_unprescaled_algo_last_collision_in_train_", "uGT: Unprescaled Algorithm Trigger Bits Deno vs. BX Number In Event;Bunch Crossing Number In Event;Algorithm Trigger Names + Bits", 5, -2.5, 2.5, unprescaledAlgoBitName_.size(), -0.5, unpreAlgs_d-0.5);
162  for(unsigned int algo=0; algo<unprescaledAlgoBitName_.size(); ++algo) {
164  }
165 }
MonitorElement * prescaled_algo_first_collision_in_train_
MonitorElement * den_unprescaled_algo_first_collision_in_train_
std::vector< std::pair< std::string, int > > prescaledAlgoBitName_
MonitorElement * den_prescaled_algo_first_collision_in_train_
MonitorElement * den_prescaled_algo_last_collision_in_train_
MonitorElement * unprescaled_algo_isolated_collision_
MonitorElement * isolated_collision_
MonitorElement * den_first_collision_in_train_
MonitorElement * den_unprescaled_algo_last_collision_in_train_
std::vector< std::pair< std::string, int > > unprescaledAlgoBitName_
void setBinLabel(int bin, const std::string &label, int axis=1)
set bin label for x, y or z axis (axis=1, 2, 3 respectively)
MonitorElement * den_first_collision_in_train_minus2_
MonitorElement * last_collision_in_train_
MonitorElement * first_collision_in_train_minus2_
void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:268
MonitorElement * unprescaled_algo_last_collision_in_train_
MonitorElement * prescaled_algo_isolated_collision_
MonitorElement * first_collision_in_train_minus1_
MonitorElement * den_last_collision_in_train_
MonitorElement * den_first_collision_in_train_minus1_
MonitorElement * book2D(Args &&...args)
Definition: DQMStore.h:109
MonitorElement * den_prescaled_algo_isolated_collision_
MonitorElement * den_isolated_collision_
MonitorElement * den_unprescaled_algo_isolated_collision_
MonitorElement * prescaled_algo_last_collision_in_train_
MonitorElement * unprescaled_algo_first_collision_in_train_
MonitorElement * first_collision_in_train_
void L1TStage2uGTTiming::dqmBeginRun ( const edm::Run ,
const edm::EventSetup  
)
overrideprotected

Definition at line 38 of file L1TStage2uGTTiming.cc.

References algoBitFirstBxInTrain_, algoBitIsoBx_, algoBitLastBxInTrain_, algoNameFirstBxInTrain_, algoNameIsoBx_, algoNameLastBxInTrain_, gtUtil_, mps_fire::i, numAlgs_, prescaledAlgoBitName_, prescaledAlgoShortList_, unprescaledAlgoBitName_, and unprescaledAlgoShortList_.

38  {
39  // Get the trigger menu information
40  gtUtil_->retrieveL1Setup(evtSetup);
41 
42  // Find the number of algos defined
43  numAlgs_ = static_cast<int>(gtUtil_->decisionsInitial().size());
44 
45  // Get the algo bits needed for the timing histograms
46  if (!gtUtil_->getAlgBitFromName(algoNameFirstBxInTrain_, algoBitFirstBxInTrain_)) {
47  edm::LogWarning("L1TStage2uGTTiming") << "Algo \"" << algoNameFirstBxInTrain_ << "\" not found in the trigger menu " << gtUtil_->gtTriggerMenuName() << ". Could not retrieve algo bit number.";
48  }
49 
50  if (!gtUtil_->getAlgBitFromName(algoNameLastBxInTrain_, algoBitLastBxInTrain_)) {
51  edm::LogWarning("L1TStage2uGTTiming") << "Algo \"" << algoNameLastBxInTrain_ << "\" not found in the trigger menu " << gtUtil_->gtTriggerMenuName() << ". Could not retrieve algo bit number.";
52  }
53 
54  if (!gtUtil_->getAlgBitFromName(algoNameIsoBx_, algoBitIsoBx_)) {
55  edm::LogWarning("L1TStage2uGTTiming") << "Algo \"" << algoNameIsoBx_ << "\" not found in the trigger menu " << gtUtil_->gtTriggerMenuName() << ". Could not retrieve algo bit number.";
56  }
57 
58  int algoBitUnpre_=-1;
59  for(unsigned int i=0;i<unprescaledAlgoShortList_.size();i++){
60  if (gtUtil_->getAlgBitFromName(unprescaledAlgoShortList_.at(i), algoBitUnpre_)) {
61  unprescaledAlgoBitName_.emplace_back(unprescaledAlgoShortList_.at(i), algoBitUnpre_);
62  }
63  else {
64  edm::LogWarning("L1TStage2uGTTiming") << "Algo \"" << unprescaledAlgoShortList_.at(i) << "\" not found in the trigger menu " << gtUtil_->gtTriggerMenuName() << ". Could not retrieve algo bit number.";
65  }
66  }
67 
68  int algoBitPre_=-1;
69  for(unsigned int i=0;i<prescaledAlgoShortList_.size();i++){
70  if ((gtUtil_->getAlgBitFromName(prescaledAlgoShortList_.at(i), algoBitPre_))) {
71  prescaledAlgoBitName_.emplace_back(prescaledAlgoShortList_.at(i), algoBitPre_);
72  }
73  else {
74  edm::LogWarning("L1TStage2uGTTiming") << "Algo \"" << prescaledAlgoShortList_.at(i) << "\" not found in the trigger menu " << gtUtil_->gtTriggerMenuName() << ". Could not retrieve algo bit number.";
75  }
76  }
77 }
std::vector< std::pair< std::string, int > > prescaledAlgoBitName_
std::vector< std::pair< std::string, int > > unprescaledAlgoBitName_
const std::string algoNameLastBxInTrain_
std::shared_ptr< l1t::L1TGlobalUtil > gtUtil_
std::vector< std::string > unprescaledAlgoShortList_
std::vector< std::string > prescaledAlgoShortList_
const std::string algoNameIsoBx_
const std::string algoNameFirstBxInTrain_

Member Data Documentation

int L1TStage2uGTTiming::algoBitFirstBxInTrain_
private

Definition at line 63 of file L1TStage2uGTTiming.h.

Referenced by analyze(), and dqmBeginRun().

int L1TStage2uGTTiming::algoBitIsoBx_
private

Definition at line 65 of file L1TStage2uGTTiming.h.

Referenced by analyze(), and dqmBeginRun().

int L1TStage2uGTTiming::algoBitLastBxInTrain_
private

Definition at line 64 of file L1TStage2uGTTiming.h.

Referenced by analyze(), and dqmBeginRun().

const std::string L1TStage2uGTTiming::algoNameFirstBxInTrain_
private

Definition at line 66 of file L1TStage2uGTTiming.h.

Referenced by dqmBeginRun().

const std::string L1TStage2uGTTiming::algoNameIsoBx_
private

Definition at line 68 of file L1TStage2uGTTiming.h.

Referenced by dqmBeginRun().

const std::string L1TStage2uGTTiming::algoNameLastBxInTrain_
private

Definition at line 67 of file L1TStage2uGTTiming.h.

Referenced by dqmBeginRun().

MonitorElement* L1TStage2uGTTiming::den_first_collision_in_train_
private

Definition at line 87 of file L1TStage2uGTTiming.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* L1TStage2uGTTiming::den_first_collision_in_train_minus1_
private

Definition at line 86 of file L1TStage2uGTTiming.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* L1TStage2uGTTiming::den_first_collision_in_train_minus2_
private

Definition at line 85 of file L1TStage2uGTTiming.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* L1TStage2uGTTiming::den_isolated_collision_
private

Definition at line 89 of file L1TStage2uGTTiming.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* L1TStage2uGTTiming::den_last_collision_in_train_
private

Definition at line 88 of file L1TStage2uGTTiming.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* L1TStage2uGTTiming::den_prescaled_algo_first_collision_in_train_
private

Definition at line 99 of file L1TStage2uGTTiming.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* L1TStage2uGTTiming::den_prescaled_algo_isolated_collision_
private

Definition at line 101 of file L1TStage2uGTTiming.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* L1TStage2uGTTiming::den_prescaled_algo_last_collision_in_train_
private

Definition at line 103 of file L1TStage2uGTTiming.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* L1TStage2uGTTiming::den_unprescaled_algo_first_collision_in_train_
private

Definition at line 100 of file L1TStage2uGTTiming.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* L1TStage2uGTTiming::den_unprescaled_algo_isolated_collision_
private

Definition at line 102 of file L1TStage2uGTTiming.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* L1TStage2uGTTiming::den_unprescaled_algo_last_collision_in_train_
private

Definition at line 104 of file L1TStage2uGTTiming.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* L1TStage2uGTTiming::first_collision_in_train_
private

Definition at line 81 of file L1TStage2uGTTiming.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* L1TStage2uGTTiming::first_collision_in_train_minus1_
private

Definition at line 80 of file L1TStage2uGTTiming.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* L1TStage2uGTTiming::first_collision_in_train_minus2_
private

Definition at line 79 of file L1TStage2uGTTiming.h.

Referenced by analyze(), and bookHistograms().

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

Definition at line 58 of file L1TStage2uGTTiming.h.

Referenced by dqmBeginRun().

MonitorElement* L1TStage2uGTTiming::isolated_collision_
private

Definition at line 83 of file L1TStage2uGTTiming.h.

Referenced by analyze(), and bookHistograms().

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

Definition at line 51 of file L1TStage2uGTTiming.h.

Referenced by analyze().

MonitorElement* L1TStage2uGTTiming::last_collision_in_train_
private

Definition at line 82 of file L1TStage2uGTTiming.h.

Referenced by analyze(), and bookHistograms().

std::string L1TStage2uGTTiming::monitorDir_
private

Definition at line 53 of file L1TStage2uGTTiming.h.

Referenced by bookHistograms().

int L1TStage2uGTTiming::numAlgs_
private

Definition at line 60 of file L1TStage2uGTTiming.h.

Referenced by bookHistograms(), and dqmBeginRun().

MonitorElement* L1TStage2uGTTiming::prescaled_algo_first_collision_in_train_
private

Definition at line 92 of file L1TStage2uGTTiming.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* L1TStage2uGTTiming::prescaled_algo_isolated_collision_
private

Definition at line 94 of file L1TStage2uGTTiming.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* L1TStage2uGTTiming::prescaled_algo_last_collision_in_train_
private

Definition at line 96 of file L1TStage2uGTTiming.h.

Referenced by analyze(), and bookHistograms().

std::vector<std::pair<std::string,int> > L1TStage2uGTTiming::prescaledAlgoBitName_
private

Definition at line 76 of file L1TStage2uGTTiming.h.

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

std::vector<std::string> L1TStage2uGTTiming::prescaledAlgoShortList_
private

Definition at line 73 of file L1TStage2uGTTiming.h.

Referenced by dqmBeginRun().

MonitorElement* L1TStage2uGTTiming::unprescaled_algo_first_collision_in_train_
private

Definition at line 93 of file L1TStage2uGTTiming.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* L1TStage2uGTTiming::unprescaled_algo_isolated_collision_
private

Definition at line 95 of file L1TStage2uGTTiming.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* L1TStage2uGTTiming::unprescaled_algo_last_collision_in_train_
private

Definition at line 97 of file L1TStage2uGTTiming.h.

Referenced by analyze(), and bookHistograms().

std::vector<std::pair<std::string,int> > L1TStage2uGTTiming::unprescaledAlgoBitName_
private

Definition at line 75 of file L1TStage2uGTTiming.h.

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

std::vector<std::string> L1TStage2uGTTiming::unprescaledAlgoShortList_
private

Definition at line 72 of file L1TStage2uGTTiming.h.

Referenced by dqmBeginRun().

unsigned int L1TStage2uGTTiming::useAlgoDecision_
private

Definition at line 70 of file L1TStage2uGTTiming.h.

Referenced by analyze(), and L1TStage2uGTTiming().

bool L1TStage2uGTTiming::verbose_
private

Definition at line 55 of file L1TStage2uGTTiming.h.

Referenced by analyze().