DESY Hbb Analysis Framework
TriggerAccepts.cc
Go to the documentation of this file.
1 
8 //
9 // Original Author: Roberval Walsh Bastos Rangel
10 // Created: Mon, 20 Oct 2014 14:24:08 GMT
11 //
12 //
13 
14 // system include files
15 #include <iostream>
16 #include <boost/algorithm/string.hpp>
17 
18 //
19 // user include files
20 #include "FWCore/Framework/interface/Frameworkfwd.h"
21 #include "FWCore/Framework/interface/EDProducer.h"
22 #include "FWCore/Framework/interface/Event.h"
23 #include "FWCore/ParameterSet/interface/ParameterSet.h"
24 
26 
27 //
28 // class declaration
29 //
30 
31 using namespace analysis;
32 using namespace analysis::ntuple;
33 
34 //
35 // constructors and destructor
36 //
38 {
39  // default constructor
40 }
41 
42 TriggerAccepts::TriggerAccepts(const edm::InputTag& tag, TTree* tree, const std::vector<std::string>& paths, const std::vector<std::string>& seeds, const std::shared_ptr<HLTPrescaleProvider> hltPrescale)
43 {
44  hlt_prescale_ = hltPrescale;
45  input_collection_ = tag;
46  tree_ = tree;
47  paths_.clear();
48  seeds_.clear();
49  paths_ = paths;
50  seeds_ = seeds;
51 
52  // remove duplicates of paths
53  sort( paths_.begin(), paths_.end() );
54  paths_.erase( unique( paths_.begin(), paths_.end() ), paths_.end() );
55  // remove duplicates of seeds
56  sort( seeds_.begin(), seeds_.end() );
57  seeds_.erase( unique( seeds_.begin(), seeds_.end() ), seeds_.end() );
58 
59  first_ = true;
60  psinfo_ = true;
61 }
62 
64 {
65  // do anything here that needs to be done at desctruction time
66  // (e.g. close files, deallocate resources etc.)
67 }
68 
69 
70 //
71 // member functions
72 //
73 
74 // ------------ method called for each event ------------
75 void TriggerAccepts::Fill(const edm::Event& event, const edm::EventSetup & setup)
76 {
77  using namespace edm;
78 
79  // reset trigger accepts and prescales to default -1
80  for (size_t i = 0; i < paths_.size() ; ++i )
81  {
82  accept_[i] = false;
83  pshlt_[i] = -1;
84  }
85  std::map<std::string, bool> l1done; // L1 prescale only once per event
86  for (size_t i = 0; i < seeds_.size() ; ++i )
87  {
88  psl1_[i] = -1;
89  l1accept_[i] = false;
90  l1done[seeds_[i]] = false;
91  }
92 
93  Handle<TriggerResults> handler;
94  event.getByLabel(input_collection_, handler);
95  const TriggerResults & triggers = *(handler.product());
96 
97  // l1 accept
98  for ( size_t j = 0 ; j < hlt_config_.size() ; ++j )
99  {
100  for (size_t i = 0; i < paths_.size() ; ++i )
101  {
102  if ( hlt_config_.triggerName(j).find(paths_[i]) == 0 )
103  {
104  // trigger accepted?
105  accept_[i] = triggers.accept(j);
106 
107  // get prescale info if requested
108  if ( psinfo_ )
109  {
110  const std::pair<std::vector<std::pair<std::string,int> >,int> ps = hlt_prescale_->prescaleValuesInDetail(event,setup,hlt_config_.triggerName(j));
111  // HLT prescale
112  pshlt_[i] = ps.second;
113  // Get L1 prescale of all seeds of the path
114  for ( size_t k = 0; k < ps.first.size(); ++k ) // loop over seeds of the path
115  {
116  for ( size_t l = 0; l < seeds_.size(); ++l ) // loop over seeds passed by python config
117  {
118  if ( ! l1done[seeds_[l]] && ps.first[k].first == seeds_[l] ) // if prescale of L1 seed not read and seed is in path
119  {
120  psl1_[l] = ps.first[k].second;
121  l1done[seeds_[l]] = true;
122  hlt_prescale_->l1tGlobalUtil().getFinalDecisionByName (seeds_[l], l1accept_[l]);
123  break;
124  }
125  }
126  }
127  }
128  else
129  {
130  std::vector<std::string> l1seeds = hlt_config_.hltL1TSeeds(hlt_config_.triggerName(j));
131  for ( size_t l = 0; l < seeds_.size(); ++l ) // loop over seeds passed by python config
132  {
133  for ( auto & l1 : l1seeds )
134  {
135  if ( l1.find(seeds_[l]) == 0 && ! l1done[seeds_[l]] )
136  {
137  hlt_prescale_->l1tGlobalUtil().getFinalDecisionByName (seeds_[l], l1accept_[l]);
138  l1done[seeds_[l]] = true;
139  break;
140  }
141  }
142  }
143  }
144  }
145  }
146  }
147 
148  tree_ -> Fill();
149 
150 }
151 
152 // ------------ method called once each job just before starting event loop ------------
154 {
155  // two loops for separation of accepts and prescales(?)
156  for (size_t i = 0; i < paths_.size() ; ++i )
157  {
158  tree_->Branch(paths_[i].c_str(), &accept_[i], (paths_[i]+"/O").c_str());
159  }
160  for (size_t i = 0; i < seeds_.size() ; ++i )
161  {
162  tree_->Branch(seeds_[i].c_str(), &l1accept_[i], (seeds_[i]+"/O").c_str());
163  }
164  for (size_t i = 0; i < paths_.size() ; ++i )
165  {
166  tree_->Branch(("ps_"+paths_[i]).c_str(), &pshlt_[i], ("ps_"+paths_[i]+"/I").c_str());
167  }
168  for (size_t i = 0; i < seeds_.size() ; ++i )
169  {
170  tree_->Branch(("ps_"+seeds_[i]).c_str(), &psl1_[i], ("ps_"+seeds_[i]+"/I").c_str());
171  }
172  // std::cout << "TriggerAccepts Branches ok" << std::endl;
173 }
174 
175 void TriggerAccepts::Run(edm::Run const & run, edm::EventSetup const& setup)
176 {
177  bool changed;
178  hlt_prescale_->init(run, setup, input_collection_.process(), changed);
179  hlt_config_ = hlt_prescale_->hltConfigProvider();
180 
181 }
182 void TriggerAccepts::ReadPrescaleInfo(const bool & ok)
183 {
184  psinfo_ = ok;
185 }
187 {
188  return psinfo_;
189 }
190 
192 {
193  Branches();
194 }
void Run(edm::Run const &, edm::EventSetup const &)
std::vector< std::string > paths_
std::shared_ptr< HLTPrescaleProvider > hlt_prescale_
void Fill(const edm::Event &event, const edm::EventSetup &setup)
std::vector< std::string > seeds_