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
 
int nhltpart
 
std::string processName_
 
TString * techBitToName
 
int * trigflag
 
int * trigPrescl
 
int UnpackBxInEvent
 

Detailed Description

$Date: November 2006 $Revision:

Author
P. Bargassa - Rice U. $Date: April 2016 $Revision:
G. Karapostoli - ULB

Definition at line 52 of file HLTInfo.h.

Constructor & Destructor Documentation

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

Definition at line 108 of file HLTInfo.h.

108 : HLTInfo(pset, iC, module) {}
HLTInfo()
Definition: HLTInfo.cc:22
Definition: vlib.h:198
template<typename T >
HLTInfo::HLTInfo ( edm::ParameterSet const &  pset,
edm::ConsumesCollector iC,
T module 
)

Definition at line 111 of file HLTInfo.h.

References hltPrescaleProvider_.

111  : HLTInfo() {
113 }
HLTInfo()
Definition: HLTInfo.cc:22
Definition: vlib.h:198
std::unique_ptr< HLTPrescaleProvider > hltPrescaleProvider_
Definition: HLTInfo.h:88
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:104
int UnpackBxInEvent
Definition: HLTInfo.h:92
bool _OR_BXes
Definition: HLTInfo.h:91

Member Function Documentation

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(), edm::HLTGlobalStatus::accept(), algoBitToName, L1GtUtils::AlgorithmTrigger, BXVector< T >::at(), beam_dqm_sourceclient-live_cfg::cerr, gather_cfg::cout, dummyBranches_, edm::EventSetup::get(), GlobalAlgBlk::getAlgoDecisionFinal(), L1TUtmTriggerMenu::getAlgorithmMap(), HltEvtCnt, hltPrescaleProvider_, edm::HandleBase::isValid(), L1EvtCnt, l1flag, l1Prescl, GlobalAlgBlk::maxPhysicsTriggers, optionsL1T::menu, L1GtUtils::prescaleFactor(), L1GtUtils::prescaleFactorSetIndex(), mps_fire::result, edm::HLTGlobalStatus::size(), BXVector< T >::size(), AlCaHLTBitMon_QueryRunRegistry::string, trigflag, edm::TriggerNames::triggerName(), L1TEGammaOffline_cfi::triggerNames, edm::Event::triggerNames(), EgHLTOffTrigSelection_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] = hltConfig_.prescaleValue(iEvent, eventSetup, trigName);
140  trigPrescl[itrig] = hltPrescaleProvider_->prescaleValue(iEvent, eventSetup, trigName);
141 
142  if (accept) {
143  trigflag[itrig] = 1;
144  } else {
145  trigflag[itrig] = 0;
146  }
147 
148  if (_Debug) {
149  if (_Debug)
150  std::cout << "%HLTInfo -- Number of HLT Triggers: " << ntrigs << std::endl;
151  std::cout << "%HLTInfo -- HLTTrigger(" << itrig << "): " << trigName << " = " << accept << std::endl;
152  }
153  }
154  } else {
155  if (_Debug)
156  std::cout << "%HLTInfo -- No Trigger Result" << std::endl;
157  }
158 
159  //==============L1 information=======================================
160 
161  // L1 Triggers from Menu
162  L1GtUtils const& l1GtUtils = hltPrescaleProvider_->l1GtUtils();
163 
164  // m_l1GtUtils.retrieveL1EventSetup(eventSetup);
165  //m_l1GtUtils.getL1GtRunCache(iEvent,eventSetup,useL1EventSetup,useL1GtTriggerMenuLite);
166  /*
167  unsigned long long id = eventSetup.get<L1TUtmTriggerMenuRcd>().cacheIdentifier();
168 
169  if (id != cache_id_) {
170  cache_id_ = id;
171  */
173  eventSetup.get<L1TUtmTriggerMenuRcd>().get(menu);
174  //std::map<std::string, L1TUtmAlgorithm> const & algorithmMap_ = &(menu->getAlgorithmMap());
175  /*
176  // get the bit/name association
177  for (auto const & keyval: menu->getAlgorithmMap()) {
178  std::string const & name = keyval.second.getName();
179  unsigned int index = keyval.second.getIndex();
180  std::cerr << "bit: " << index << "\tname: " << name << std::endl;
181  }
182  */
183  //} // end get menu
184 
185  int iErrorCode = -1;
187  const int pfSetIndexAlgorithmTrigger = l1GtUtils.prescaleFactorSetIndex(iEvent, trigCategory, iErrorCode);
188  if (iErrorCode == 0) {
189  if (_Debug)
190  std::cout << "%Prescale set index: " << pfSetIndexAlgorithmTrigger << std::endl;
191  } else {
192  std::cout << "%Could not extract Prescale set index from event record. Error code: " << iErrorCode << std::endl;
193  }
194 
195  // 1st event : Book as many branches as trigger paths provided in the input...
196  if (l1results.isValid()) {
197  int ntrigs = l1results->size();
198  if (ntrigs == 0) {
199  std::cout << "%L1Results -- No trigger name given in TriggerResults of the input " << std::endl;
200  }
201  /*
202  edm::TriggerNames const& triggerNames = iEvent.triggerNames(&results);
203  // 1st event : Book as many branches as trigger paths provided in the input...
204  */
205  if (L1EvtCnt == 0) {
206  // get the bit/name association
207  for (auto const& keyval : menu->getAlgorithmMap()) {
208  std::string const& trigName = keyval.second.getName();
209  unsigned int index = keyval.second.getIndex();
210  if (_Debug)
211  std::cerr << "bit: " << index << "\tname: " << trigName << std::endl;
212 
213  int itrig = index;
214  algoBitToName[itrig] = TString(trigName);
215 
216  TString l1trigName = static_cast<const char*>(algoBitToName[itrig]);
217  std::string l1triggername = static_cast<const char*>(algoBitToName[itrig]);
218 
219  HltTree->Branch(l1trigName, l1flag + itrig, l1trigName + "/I");
220  HltTree->Branch(l1trigName + "_Prescl", l1Prescl + itrig, l1trigName + "_Prescl/I");
221 
222  } // end algo Map
223 
224  L1EvtCnt++;
225  } // end l1evtCnt=0
226 
227  GlobalAlgBlk const& result = l1results->at(0, 0);
228 
229  // get the individual decisions from the GlobalAlgBlk
230  for (unsigned int itrig = 0; itrig < result.maxPhysicsTriggers; ++itrig) {
231  // std::cerr << "bit: " << itrig << "\tresult: " << results.getAlgoDecisionFinal(itrig) << std::endl;
232 
233  bool myflag = result.getAlgoDecisionFinal(itrig);
234  if (myflag) {
235  l1flag[itrig] = 1;
236  } else {
237  l1flag[itrig] = 0;
238  }
239 
240  std::string l1triggername = static_cast<const char*>(algoBitToName[itrig]);
241  l1Prescl[itrig] = l1GtUtils.prescaleFactor(iEvent, l1triggername, iErrorCode);
242 
243  if (_Debug)
244  std::cout << "L1 TD: " << itrig << " " << algoBitToName[itrig] << " " << l1flag[itrig] << " " << l1Prescl[itrig]
245  << std::endl;
246  }
247 
248  // L1EvtCnt++;
249  if (_Debug)
250  std::cout << "%L1Info -- Done with routine" << std::endl;
251 
252  } // l1results.isValid
253  else {
254  if (_Debug)
255  std::cout << "%L1Results -- No Trigger Result" << std::endl;
256  }
257 }
bool _Debug
Definition: HLTInfo.h:104
int L1EvtCnt
Definition: HLTInfo.h:77
std::vector< std::string > dummyBranches_
Definition: HLTInfo.h:84
unsigned size(int bx) const
std::vector< bool > const & getAlgoDecisionFinal() const
Definition: GlobalAlgBlk.h:84
int * l1Prescl
Definition: HLTInfo.h:80
const int prescaleFactorSetIndex(const edm::Event &iEvent, const TriggerCategory &trigCategory, int &errorCode) const
Definition: L1GtUtils.cc:1227
TriggerCategory
Definition: L1GtUtils.h:96
bool accept() const
Has at least one path accepted the event?
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:30
int iEvent
Definition: GenABIO.cc:224
TString * algoBitToName
Definition: HLTInfo.h:82
unsigned int size() const
Get number of paths stored.
int HltEvtCnt
Definition: HLTInfo.h:77
bool isValid() const
Definition: HandleBase.h:70
int * l1flag
Definition: HLTInfo.h:79
int * trigPrescl
Definition: HLTInfo.h:80
static unsigned int maxPhysicsTriggers
Definition: GlobalAlgBlk.h:52
std::string const & triggerName(unsigned int index) const
Definition: TriggerNames.cc:22
const std::map< std::string, L1TUtmAlgorithm > & getAlgorithmMap() const
int * trigflag
Definition: HLTInfo.h:79
std::unique_ptr< HLTPrescaleProvider > hltPrescaleProvider_
Definition: HLTInfo.h:88
const T & at(int bx, unsigned i) const
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:1036
void HLTInfo::beginRun ( const edm::Run run,
const edm::EventSetup c 
)

Definition at line 29 of file HLTInfo.cc.

References gather_cfg::cout, hltPrescaleProvider_, and processName_.

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:89
std::unique_ptr< HLTPrescaleProvider > hltPrescaleProvider_
Definition: HLTInfo.h:88
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 int[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:104
int L1EvtCnt
Definition: HLTInfo.h:77
T getParameter(std::string const &) const
std::vector< std::string > dummyBranches_
Definition: HLTInfo.h:84
T getUntrackedParameter(std::string const &, T const &) const
int * l1Prescl
Definition: HLTInfo.h:80
int * l1techPrescl
Definition: HLTInfo.h:80
float * hltpeta
Definition: HLTInfo.h:76
std::string processName_
Definition: HLTInfo.h:89
TString * algoBitToName
Definition: HLTInfo.h:82
int HltEvtCnt
Definition: HLTInfo.h:77
std::vector< std::string > getParameterNames() const
int * l1flag
Definition: HLTInfo.h:79
TString * techBitToName
Definition: HLTInfo.h:83
int * trigPrescl
Definition: HLTInfo.h:80
int * l1techflag
Definition: HLTInfo.h:79
int * trigflag
Definition: HLTInfo.h:79
int * l1flag5Bx
Definition: HLTInfo.h:79
float * hltppt
Definition: HLTInfo.h:76

Member Data Documentation

bool HLTInfo::_Debug
private

Definition at line 104 of file HLTInfo.h.

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

bool HLTInfo::_OR_BXes
private

Definition at line 91 of file HLTInfo.h.

Referenced by HLTInfo().

TString* HLTInfo::algoBitToName
private

Definition at line 82 of file HLTInfo.h.

Referenced by analyze(), and setup().

unsigned long long HLTInfo::cache_id_
private

Definition at line 97 of file HLTInfo.h.

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

Definition at line 84 of file HLTInfo.h.

Referenced by analyze(), and setup().

int HLTInfo::HltEvtCnt
private

Definition at line 77 of file HLTInfo.h.

Referenced by analyze(), and setup().

float * HLTInfo::hltpeta
private

Definition at line 76 of file HLTInfo.h.

Referenced by setup().

float* HLTInfo::hltppt
private

Definition at line 76 of file HLTInfo.h.

Referenced by setup().

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

Definition at line 88 of file HLTInfo.h.

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

int HLTInfo::L1EvtCnt
private

Definition at line 77 of file HLTInfo.h.

Referenced by analyze(), and setup().

int * HLTInfo::l1flag
private

Definition at line 79 of file HLTInfo.h.

Referenced by analyze(), and setup().

int * HLTInfo::l1flag5Bx
private

Definition at line 79 of file HLTInfo.h.

Referenced by setup().

int * HLTInfo::l1Prescl
private

Definition at line 80 of file HLTInfo.h.

Referenced by analyze(), and setup().

int * HLTInfo::l1techflag
private

Definition at line 79 of file HLTInfo.h.

Referenced by setup().

int * HLTInfo::l1techPrescl
private

Definition at line 80 of file HLTInfo.h.

Referenced by setup().

int HLTInfo::nhltpart
private

Definition at line 77 of file HLTInfo.h.

std::string HLTInfo::processName_
private

Definition at line 89 of file HLTInfo.h.

Referenced by beginRun(), and setup().

TString* HLTInfo::techBitToName
private

Definition at line 83 of file HLTInfo.h.

Referenced by setup().

int* HLTInfo::trigflag
private

Definition at line 79 of file HLTInfo.h.

Referenced by analyze(), and setup().

int* HLTInfo::trigPrescl
private

Definition at line 80 of file HLTInfo.h.

Referenced by analyze(), and setup().

int HLTInfo::UnpackBxInEvent
private

Definition at line 92 of file HLTInfo.h.

Referenced by HLTInfo().