CMS 3D CMS Logo

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

#include <HLTInfo.h>

Public Member Functions

void analyze (const edm::Handle< edm::TriggerResults > &hltresults, const edm::Handle< GlobalAlgBlkBxCollection > &l1results, edm::EventSetup const &eventSetup, edm::Event const &iEvent, TTree *tree)
 
void beginRun (const edm::Run &, const edm::EventSetup &)
 
template<typename T >
 HLTInfo (edm::ParameterSet const &pset, edm::ConsumesCollector &&iC, T &module)
 
template<typename T >
 HLTInfo (edm::ParameterSet const &pset, edm::ConsumesCollector &iC, T &module)
 
void setup (const edm::ParameterSet &pSet, TTree *tree)
 

Private Member Functions

 HLTInfo ()
 

Private Attributes

bool _Debug
 
bool _OR_BXes
 
TString * algoBitToName
 
unsigned long long cache_id_
 
std::vector< std::string > dummyBranches_
 
int HltEvtCnt
 
float * hltpeta
 
float * hltppt
 
std::unique_ptr< HLTPrescaleProviderhltPrescaleProvider_
 
int L1EvtCnt
 
int * l1flag
 
int * l1flag5Bx
 
int * l1Prescl
 
int * l1techflag
 
int * l1techPrescl
 
edm::ESGetToken< L1TUtmTriggerMenu, L1TUtmTriggerMenuRcdl1tUtmTriggerMenuToken_
 
int nhltpart
 
std::string processName_
 
TString * techBitToName
 
int * trigflag
 
double * trigPrescl
 
int UnpackBxInEvent
 

Detailed Description

$Date: November 2006 $Revision:

Author
P. Bargassa - Rice U.

$Date: April 2016 $Revision:

Author
G. Karapostoli - ULB

Definition at line 53 of file HLTInfo.h.

Constructor & Destructor Documentation

◆ HLTInfo() [1/3]

template<typename T >
HLTInfo::HLTInfo ( edm::ParameterSet const &  pset,
edm::ConsumesCollector &&  iC,
T module 
)

Definition at line 112 of file HLTInfo.h.

◆ HLTInfo() [2/3]

template<typename T >
HLTInfo::HLTInfo ( edm::ParameterSet const &  pset,
edm::ConsumesCollector iC,
T module 
)

Definition at line 115 of file HLTInfo.h.

References edm::ConsumesCollector::esConsumes(), hltPrescaleProvider_, l1tUtmTriggerMenuToken_, callgraph::module, and muonDTDigis_cfi::pset.

115  : HLTInfo() {
117  hltPrescaleProvider_ = std::make_unique<HLTPrescaleProvider>(pset, iC, module);
118 }
HLTInfo()
Definition: HLTInfo.cc:22
edm::ESGetToken< L1TUtmTriggerMenu, L1TUtmTriggerMenuRcd > l1tUtmTriggerMenuToken_
Definition: HLTInfo.h:76
std::unique_ptr< HLTPrescaleProvider > hltPrescaleProvider_
Definition: HLTInfo.h:92

◆ HLTInfo() [3/3]

HLTInfo::HLTInfo ( )
private

Definition at line 22 of file HLTInfo.cc.

References _Debug, _OR_BXes, and UnpackBxInEvent.

22  {
23  //set parameter defaults
24  _Debug = false;
25  _OR_BXes = false;
26  UnpackBxInEvent = 1;
27 }
bool _Debug
Definition: HLTInfo.h:108
int UnpackBxInEvent
Definition: HLTInfo.h:96
bool _OR_BXes
Definition: HLTInfo.h:95

Member Function Documentation

◆ analyze()

void HLTInfo::analyze ( const edm::Handle< edm::TriggerResults > &  hltresults,
const edm::Handle< GlobalAlgBlkBxCollection > &  l1results,
edm::EventSetup const &  eventSetup,
edm::Event const &  iEvent,
TTree *  tree 
)

Analyze the Data

Definition at line 84 of file HLTInfo.cc.

References _Debug, accept(), algoBitToName, L1GtUtils::AlgorithmTrigger, DMR_cfg::cerr, gather_cfg::cout, dummyBranches_, options_cfi::eventSetup, HltEvtCnt, hltPrescaleProvider_, HLTBitAnalyser_cfi::hltresults, iEvent, L1EvtCnt, l1flag, l1Prescl, HLTBitAnalyser_cfi::l1results, l1tUtmTriggerMenuToken_, prototype_2023_v1_0_0::menu, L1GtUtils::prescaleFactor(), L1GtUtils::prescaleFactorSetIndex(), mps_fire::result, AlCaHLTBitMon_QueryRunRegistry::string, trigflag, L1TEGammaOffline_cfi::triggerNames, cscTnPEfficiencyTask_cfi::trigName, and trigPrescl.

Referenced by HLTBitAnalyzer::analyze().

88  {
89  // std::cout << " Beginning HLTInfo " << std::endl;
90 
92  if (hltresults.isValid()) {
93  int ntrigs = hltresults->size();
94  if (ntrigs == 0) {
95  std::cout << "%HLTInfo -- No trigger name given in TriggerResults of the input " << std::endl;
96  }
97 
98  edm::TriggerNames const& triggerNames = iEvent.triggerNames(*hltresults);
99 
100  // 1st event : Book as many branches as trigger paths provided in the input...
101  if (HltEvtCnt == 0) {
102  for (int itrig = 0; itrig != ntrigs; ++itrig) {
103  TString trigName = triggerNames.triggerName(itrig);
104  HltTree->Branch(trigName, trigflag + itrig, trigName + "/I");
105  HltTree->Branch(trigName + "_Prescl", trigPrescl + itrig, trigName + "_Prescl/I");
106  }
107 
108  int itdum = ntrigs;
109  for (auto& dummyBranche : dummyBranches_) {
110  TString trigName(dummyBranche.data());
111  bool addThisBranch = true;
112  for (int itrig = 0; itrig != ntrigs; ++itrig) {
113  TString realTrigName = triggerNames.triggerName(itrig);
114  if (trigName == realTrigName)
115  addThisBranch = false;
116  }
117  if (addThisBranch) {
118  HltTree->Branch(trigName, trigflag + itdum, trigName + "/I");
119  HltTree->Branch(trigName + "_Prescl", trigPrescl + itdum, trigName + "_Prescl/I");
120  trigflag[itdum] = 0;
121  trigPrescl[itdum] = 0;
122  ++itdum;
123  }
124  }
125 
126  HltEvtCnt++;
127  }
128  // ...Fill the corresponding accepts in branch-variables
129 
130  //std::cout << "Number of prescale sets: " << hltConfig_.prescaleSize() << std::endl;
131  //std::cout << "Number of HLT paths: " << hltConfig_.size() << std::endl;
132  //int presclSet = hltConfig_.prescaleSet(iEvent, eventSetup);
133  //std::cout<<"\tPrescale set number: "<< presclSet <<std::endl;
134 
135  for (int itrig = 0; itrig != ntrigs; ++itrig) {
136  const std::string& trigName = triggerNames.triggerName(itrig);
137  bool accept = hltresults->accept(itrig);
138 
139  trigPrescl[itrig] = hltPrescaleProvider_->prescaleValue<double>(iEvent, eventSetup, trigName);
140 
141  if (accept) {
142  trigflag[itrig] = 1;
143  } else {
144  trigflag[itrig] = 0;
145  }
146 
147  if (_Debug) {
148  if (_Debug)
149  std::cout << "%HLTInfo -- Number of HLT Triggers: " << ntrigs << std::endl;
150  std::cout << "%HLTInfo -- HLTTrigger(" << itrig << "): " << trigName << " = " << accept << std::endl;
151  }
152  }
153  } else {
154  if (_Debug)
155  std::cout << "%HLTInfo -- No Trigger Result" << std::endl;
156  }
157 
158  //==============L1 information=======================================
159 
160  // L1 Triggers from Menu
161  L1GtUtils const& l1GtUtils = hltPrescaleProvider_->l1GtUtils();
162 
163  // m_l1GtUtils.retrieveL1EventSetup(eventSetup);
164  //m_l1GtUtils.getL1GtRunCache(iEvent,eventSetup,useL1EventSetup,useL1GtTriggerMenuLite);
165  /*
166  unsigned long long id = eventSetup.get<L1TUtmTriggerMenuRcd>().cacheIdentifier();
167 
168  if (id != cache_id_) {
169  cache_id_ = id;
170  */
171  auto const& menu = eventSetup.getHandle(l1tUtmTriggerMenuToken_);
172 
173  //std::map<std::string, L1TUtmAlgorithm> const & algorithmMap_ = &(menu->getAlgorithmMap());
174  /*
175  // get the bit/name association
176  for (auto const & keyval: menu->getAlgorithmMap()) {
177  std::string const & name = keyval.second.getName();
178  unsigned int index = keyval.second.getIndex();
179  std::cerr << "bit: " << index << "\tname: " << name << std::endl;
180  }
181  */
182  //} // end get menu
183 
184  int iErrorCode = -1;
186  const int pfSetIndexAlgorithmTrigger = l1GtUtils.prescaleFactorSetIndex(iEvent, trigCategory, iErrorCode);
187  if (iErrorCode == 0) {
188  if (_Debug)
189  std::cout << "%Prescale set index: " << pfSetIndexAlgorithmTrigger << std::endl;
190  } else {
191  std::cout << "%Could not extract Prescale set index from event record. Error code: " << iErrorCode << std::endl;
192  }
193 
194  // 1st event : Book as many branches as trigger paths provided in the input...
195  if (l1results.isValid()) {
196  int ntrigs = l1results->size();
197  if (ntrigs == 0) {
198  std::cout << "%L1Results -- No trigger name given in TriggerResults of the input " << std::endl;
199  }
200  /*
201  edm::TriggerNames const& triggerNames = iEvent.triggerNames(&results);
202  // 1st event : Book as many branches as trigger paths provided in the input...
203  */
204  if (L1EvtCnt == 0) {
205  // get the bit/name association
206  for (auto const& keyval : menu->getAlgorithmMap()) {
207  std::string const& trigName = keyval.second.getName();
208  unsigned int index = keyval.second.getIndex();
209  if (_Debug)
210  std::cerr << "bit: " << index << "\tname: " << trigName << std::endl;
211 
212  int itrig = index;
213  algoBitToName[itrig] = TString(trigName);
214 
215  TString l1trigName = static_cast<const char*>(algoBitToName[itrig]);
216  std::string l1triggername = static_cast<const char*>(algoBitToName[itrig]);
217 
218  HltTree->Branch(l1trigName, l1flag + itrig, l1trigName + "/I");
219  HltTree->Branch(l1trigName + "_Prescl", l1Prescl + itrig, l1trigName + "_Prescl/I");
220 
221  } // end algo Map
222 
223  L1EvtCnt++;
224  } // end l1evtCnt=0
225 
226  GlobalAlgBlk const& result = l1results->at(0, 0);
227 
228  // get the individual decisions from the GlobalAlgBlk
229  for (unsigned int itrig = 0; itrig < result.maxPhysicsTriggers; ++itrig) {
230  // std::cerr << "bit: " << itrig << "\tresult: " << results.getAlgoDecisionFinal(itrig) << std::endl;
231 
232  bool myflag = result.getAlgoDecisionFinal(itrig);
233  if (myflag) {
234  l1flag[itrig] = 1;
235  } else {
236  l1flag[itrig] = 0;
237  }
238 
239  std::string l1triggername = static_cast<const char*>(algoBitToName[itrig]);
240  l1Prescl[itrig] = l1GtUtils.prescaleFactor(iEvent, l1triggername, iErrorCode);
241 
242  if (_Debug)
243  std::cout << "L1 TD: " << itrig << " " << algoBitToName[itrig] << " " << l1flag[itrig] << " " << l1Prescl[itrig]
244  << std::endl;
245  }
246 
247  // L1EvtCnt++;
248  if (_Debug)
249  std::cout << "%L1Info -- Done with routine" << std::endl;
250 
251  } // l1results.isValid
252  else {
253  if (_Debug)
254  std::cout << "%L1Results -- No Trigger Result" << std::endl;
255  }
256 }
const int prescaleFactorSetIndex(const edm::Event &iEvent, const TriggerCategory &trigCategory, int &errorCode) const
Definition: L1GtUtils.cc:1269
bool _Debug
Definition: HLTInfo.h:108
int L1EvtCnt
Definition: HLTInfo.h:80
double * trigPrescl
Definition: HLTInfo.h:83
std::vector< std::string > dummyBranches_
Definition: HLTInfo.h:88
edm::ESGetToken< L1TUtmTriggerMenu, L1TUtmTriggerMenuRcd > l1tUtmTriggerMenuToken_
Definition: HLTInfo.h:76
const int prescaleFactor(const edm::Event &iEvent, const std::string &nameAlgoTechTrig, int &errorCode) const
return prescale factor for a given algorithm or technical trigger
Definition: L1GtUtils.cc:1078
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:31
int iEvent
Definition: GenABIO.cc:224
TString * algoBitToName
Definition: HLTInfo.h:86
int HltEvtCnt
Definition: HLTInfo.h:80
int * l1flag
Definition: HLTInfo.h:82
int * l1Prescl
Definition: HLTInfo.h:84
int * trigflag
Definition: HLTInfo.h:82
std::unique_ptr< HLTPrescaleProvider > hltPrescaleProvider_
Definition: HLTInfo.h:92

◆ beginRun()

void HLTInfo::beginRun ( const edm::Run run,
const edm::EventSetup c 
)

Definition at line 29 of file HLTInfo.cc.

References HltBtagPostValidation_cff::c, gather_cfg::cout, hltPrescaleProvider_, processName_, and writedatasetfile::run.

Referenced by HLTBitAnalyzer::beginRun().

29  {
30  bool changed(true);
31  if (hltPrescaleProvider_->init(run, c, processName_, changed)) {
32  // if init returns TRUE, initialisation has succeeded!
33  if (changed) {
34  // The HLT config has actually changed wrt the previous Run, hence rebook your
35  // histograms or do anything else dependent on the revised HLT config
36  std::cout << "Initalizing HLTConfigProvider" << std::endl;
37  }
38  } else {
39  // if init returns FALSE, initialisation has NOT succeeded, which indicates a problem
40  // with the file and/or code and needs to be investigated!
41  std::cout << " HLT config extraction failure with process name " << processName_ << std::endl;
42  // In this case, all access methods will return empty values!
43  }
44 }
std::string processName_
Definition: HLTInfo.h:93
std::unique_ptr< HLTPrescaleProvider > hltPrescaleProvider_
Definition: HLTInfo.h:92

◆ setup()

void HLTInfo::setup ( const edm::ParameterSet pSet,
TTree *  tree 
)

Definition at line 47 of file HLTInfo.cc.

References _Debug, algoBitToName, dummyBranches_, edm::ParameterSet::getParameter(), edm::ParameterSet::getParameterNames(), edm::ParameterSet::getUntrackedParameter(), HltEvtCnt, hltpeta, hltppt, L1EvtCnt, l1flag, l1flag5Bx, l1Prescl, l1techflag, l1techPrescl, processName_, AlCaHLTBitMon_QueryRunRegistry::string, techBitToName, trigflag, and trigPrescl.

Referenced by HLTBitAnalyzer::HLTBitAnalyzer().

47  {
48  processName_ = pSet.getParameter<std::string>("HLTProcessName");
49 
50  edm::ParameterSet myHltParams = pSet.getParameter<edm::ParameterSet>("RunParameters");
51  std::vector<std::string> parameterNames = myHltParams.getParameterNames();
52 
53  for (auto& parameterName : parameterNames) {
54  if (parameterName == "Debug")
55  _Debug = myHltParams.getParameter<bool>(parameterName);
56  }
57 
58  dummyBranches_ = pSet.getUntrackedParameter<std::vector<std::string> >("dummyBranches", std::vector<std::string>(0));
59 
60  HltEvtCnt = 0;
61  const int kMaxTrigFlag = 10000;
62  trigflag = new int[kMaxTrigFlag];
63  trigPrescl = new double[kMaxTrigFlag];
64 
65  L1EvtCnt = 0;
66  const int kMaxL1Flag = 10000;
67  l1flag = new int[kMaxL1Flag];
68  l1flag5Bx = new int[kMaxTrigFlag];
69  l1Prescl = new int[kMaxL1Flag];
70 
71  l1techflag = new int[kMaxL1Flag];
72  // l1techflag5Bx = new int[kMaxTrigFlag];
73  l1techPrescl = new int[kMaxTrigFlag];
74 
75  const int kMaxHLTPart = 10000;
76  hltppt = new float[kMaxHLTPart];
77  hltpeta = new float[kMaxHLTPart];
78 
79  algoBitToName = new TString[512];
80  techBitToName = new TString[512];
81 }
bool _Debug
Definition: HLTInfo.h:108
int L1EvtCnt
Definition: HLTInfo.h:80
double * trigPrescl
Definition: HLTInfo.h:83
std::vector< std::string > dummyBranches_
Definition: HLTInfo.h:88
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
int * l1techPrescl
Definition: HLTInfo.h:84
float * hltpeta
Definition: HLTInfo.h:79
T getUntrackedParameter(std::string const &, T const &) const
std::string processName_
Definition: HLTInfo.h:93
TString * algoBitToName
Definition: HLTInfo.h:86
int HltEvtCnt
Definition: HLTInfo.h:80
int * l1flag
Definition: HLTInfo.h:82
TString * techBitToName
Definition: HLTInfo.h:87
int * l1Prescl
Definition: HLTInfo.h:84
int * l1techflag
Definition: HLTInfo.h:82
int * trigflag
Definition: HLTInfo.h:82
int * l1flag5Bx
Definition: HLTInfo.h:82
float * hltppt
Definition: HLTInfo.h:79
std::vector< std::string > getParameterNames() const

Member Data Documentation

◆ _Debug

bool HLTInfo::_Debug
private

Definition at line 108 of file HLTInfo.h.

Referenced by analyze(), HLTInfo(), and setup().

◆ _OR_BXes

bool HLTInfo::_OR_BXes
private

Definition at line 95 of file HLTInfo.h.

Referenced by HLTInfo().

◆ algoBitToName

TString* HLTInfo::algoBitToName
private

Definition at line 86 of file HLTInfo.h.

Referenced by analyze(), and setup().

◆ cache_id_

unsigned long long HLTInfo::cache_id_
private

Definition at line 101 of file HLTInfo.h.

◆ dummyBranches_

std::vector<std::string> HLTInfo::dummyBranches_
private

Definition at line 88 of file HLTInfo.h.

Referenced by analyze(), and setup().

◆ HltEvtCnt

int HLTInfo::HltEvtCnt
private

Definition at line 80 of file HLTInfo.h.

Referenced by analyze(), and setup().

◆ hltpeta

float * HLTInfo::hltpeta
private

Definition at line 79 of file HLTInfo.h.

Referenced by setup().

◆ hltppt

float* HLTInfo::hltppt
private

Definition at line 79 of file HLTInfo.h.

Referenced by setup().

◆ hltPrescaleProvider_

std::unique_ptr<HLTPrescaleProvider> HLTInfo::hltPrescaleProvider_
private

Definition at line 92 of file HLTInfo.h.

Referenced by analyze(), beginRun(), and HLTInfo().

◆ L1EvtCnt

int HLTInfo::L1EvtCnt
private

Definition at line 80 of file HLTInfo.h.

Referenced by analyze(), and setup().

◆ l1flag

int * HLTInfo::l1flag
private

Definition at line 82 of file HLTInfo.h.

Referenced by analyze(), and setup().

◆ l1flag5Bx

int * HLTInfo::l1flag5Bx
private

Definition at line 82 of file HLTInfo.h.

Referenced by setup().

◆ l1Prescl

int* HLTInfo::l1Prescl
private

Definition at line 84 of file HLTInfo.h.

Referenced by analyze(), and setup().

◆ l1techflag

int * HLTInfo::l1techflag
private

Definition at line 82 of file HLTInfo.h.

Referenced by setup().

◆ l1techPrescl

int * HLTInfo::l1techPrescl
private

Definition at line 84 of file HLTInfo.h.

Referenced by setup().

◆ l1tUtmTriggerMenuToken_

edm::ESGetToken<L1TUtmTriggerMenu, L1TUtmTriggerMenuRcd> HLTInfo::l1tUtmTriggerMenuToken_
private

Definition at line 76 of file HLTInfo.h.

Referenced by analyze(), and HLTInfo().

◆ nhltpart

int HLTInfo::nhltpart
private

Definition at line 80 of file HLTInfo.h.

◆ processName_

std::string HLTInfo::processName_
private

Definition at line 93 of file HLTInfo.h.

Referenced by beginRun(), and setup().

◆ techBitToName

TString* HLTInfo::techBitToName
private

Definition at line 87 of file HLTInfo.h.

Referenced by setup().

◆ trigflag

int* HLTInfo::trigflag
private

Definition at line 82 of file HLTInfo.h.

Referenced by analyze(), and setup().

◆ trigPrescl

double* HLTInfo::trigPrescl
private

Definition at line 83 of file HLTInfo.h.

Referenced by analyze(), and setup().

◆ UnpackBxInEvent

int HLTInfo::UnpackBxInEvent
private

Definition at line 96 of file HLTInfo.h.

Referenced by HLTInfo().