CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
SUSY_HLT_alphaT.cc
Go to the documentation of this file.
7 //#include "HLTriggerOffline/SUSYBSM/interface/AlphaT.h"
8 
9 typedef ROOT::Math::LorentzVector<ROOT::Math::PtEtaPhiM4D<double>> LorentzV;
10 
12  edm::LogInfo("SUSY_HLT_alphaT") << "Constructor SUSY_HLT_alphaT::SUSY_HLT_alphaT " << std::endl;
13  // Get parameters from configuration file
14  theTrigSummary_ = consumes<trigger::TriggerEvent>(ps.getParameter<edm::InputTag>("trigSummary"));
15  thePfJetCollection_ = consumes<reco::PFJetCollection>(ps.getParameter<edm::InputTag>("pfJetCollection"));
16  // theCaloJetCollection_ =
17  // consumes<reco::CaloJetCollection>(ps.getParameter<edm::InputTag>("caloJetCollection"));
18  triggerResults_ = consumes<edm::TriggerResults>(ps.getParameter<edm::InputTag>("TriggerResults"));
19  HLTProcess_ = ps.getParameter<std::string>("HLTProcess");
20  triggerPath_ = ps.getParameter<std::string>("TriggerPath");
21  triggerPathAuxiliaryForHadronic_ = ps.getParameter<std::string>("TriggerPathAuxiliaryForHadronic");
22  triggerFilter_ = ps.getParameter<edm::InputTag>("TriggerFilter");
23  triggerPreFilter_ = ps.getParameter<edm::InputTag>("TriggerPreFilter");
24  ptThrJet_ = ps.getUntrackedParameter<double>("PtThrJet");
25  etaThrJet_ = ps.getUntrackedParameter<double>("EtaThrJet");
26  // caloAlphaTThrTurnon_ =
27  // ps.getUntrackedParameter<double>("caloAlphaTThrTurnon"); caloHtThrTurnon_ =
28  // ps.getUntrackedParameter<double>("caloHtThrTurnon");
29  pfAlphaTThrTurnon_ = ps.getUntrackedParameter<double>("pfAlphaTThrTurnon");
30  pfHtThrTurnon_ = ps.getUntrackedParameter<double>("pfHtThrTurnon");
31 }
32 
34  edm::LogInfo("SUSY_HLT_alphaT") << "Destructor SUSY_HLT_alphaT::~SUSY_HLT_alphaT " << std::endl;
35 }
36 
38  bool changed;
39 
40  if (!fHltConfig.init(run, e, HLTProcess_, changed)) {
41  edm::LogError("SUSY_HLT_alphaT") << "Initialization of HLTConfigProvider failed!!";
42  return;
43  }
44 
45  bool pathFound = false;
46  const std::vector<std::string> allTrigNames = fHltConfig.triggerNames();
47  for (size_t j = 0; j < allTrigNames.size(); ++j) {
48  if (allTrigNames[j].find(triggerPath_) != std::string::npos) {
49  pathFound = true;
50  }
51  }
52 
53  if (!pathFound) {
54  LogDebug("SUSY_HLT_alphaT") << "Path not found"
55  << "\n";
56  return;
57  }
58  // std::vector<std::string> filtertags = fHltConfig.moduleLabels( triggerPath_
59  // ); triggerFilter_ =
60  // edm::InputTag(filtertags[filtertags.size()-1],"",fHltConfig.processName());
61  // triggerFilter_ = edm::InputTag("hltPFMET120Mu5L3PreFiltered", "",
62  // fHltConfig.processName());
63 
64  edm::LogInfo("SUSY_HLT_alphaT") << "SUSY_HLT_alphaT::beginRun" << std::endl;
65 }
66 
68  edm::LogInfo("SUSY_HLT_alphaT") << "SUSY_HLT_alphaT::bookHistograms" << std::endl;
69  // book at beginRun
70  bookHistos(ibooker_);
71 }
72 
74  edm::LogInfo("SUSY_HLT_alphaT") << "SUSY_HLT_alphaT::analyze" << std::endl;
75 
76  //-------------------------------
77  //--- Trigger
78  //-------------------------------
80  e.getByToken(triggerResults_, hltresults);
81  if (!hltresults.isValid()) {
82  edm::LogWarning("SUSY_HLT_alphaT") << "invalid collection: TriggerResults"
83  << "\n";
84  return;
85  }
87  e.getByToken(theTrigSummary_, triggerSummary);
88  if (!triggerSummary.isValid()) {
89  edm::LogWarning("SUSY_HLT_alphaT") << "invalid collection: TriggerSummary"
90  << "\n";
91  return;
92  }
93 
94  //-------------------------------
95  //--- Jets
96  //-------------------------------
98  e.getByToken(thePfJetCollection_, pfJetCollection);
99  if (!pfJetCollection.isValid()) {
100  edm::LogWarning("SUSY_HLT_alphaT") << "invalid collection: PFJets"
101  << "\n";
102  return;
103  }
104  // edm::Handle<reco::CaloJetCollection> caloJetCollection;
105  // e.getByToken (theCaloJetCollection_,caloJetCollection);
106  // if ( !caloJetCollection.isValid() ){
107  // edm::LogWarning ("SUSY_HLT_alphaT") << "invalid collection: CaloJets"
108  // << "\n"; return;
109  // }
110 
111  // get online objects
112  // For now just get the jets and recalculate ht and alphaT
113  size_t filterIndex = triggerSummary->filterIndex(triggerFilter_);
114  // size_t preFilterIndex = triggerSummary->filterIndex( triggerPreFilter_ );
115  trigger::TriggerObjectCollection triggerObjects = triggerSummary->getObjects();
116 
117  // Get the PF objects from the filter
118  double hltPfHt = 0.;
119  std::vector<LorentzV> hltPfJets;
120  if (!(filterIndex >= triggerSummary->sizeFilters())) {
121  const trigger::Keys &keys = triggerSummary->filterKeys(filterIndex);
122 
123  for (size_t j = 0; j < keys.size(); ++j) {
124  trigger::TriggerObject foundObject = triggerObjects[keys[j]];
125 
126  // if(foundObject.id() == 85){ //It's a jet
127  if (foundObject.pt() > ptThrJet_ && fabs(foundObject.eta()) < etaThrJet_) {
128  hltPfHt += foundObject.pt();
129  LorentzV JetLVec(foundObject.pt(), foundObject.eta(), foundObject.phi(), foundObject.mass());
130  hltPfJets.push_back(JetLVec);
131  }
132  // }
133  }
134  }
135 
136  // //Get the Calo objects from the prefilter
137  // double hltCaloHt=0.;
138  // std::vector<LorentzV> hltCaloJets;
139  // if( !(preFilterIndex >= triggerSummary->sizeFilters()) ){
140  // const trigger::Keys& keys = triggerSummary->filterKeys( preFilterIndex
141  // );
142 
143  // for( size_t j = 0; j < keys.size(); ++j ){
144  // trigger::TriggerObject foundObject = triggerObjects[keys[j]];
145 
146  // // if(foundObject.id() == 85){ //It's a jet
147  // if(foundObject.pt()>ptThrJet_ && fabs(foundObject.eta()) <
148  // etaThrJet_){
149  // hltCaloHt += foundObject.pt();
150  // LorentzV
151  // JetLVec(foundObject.pt(),foundObject.eta(),foundObject.phi(),foundObject.mass());
152  // hltCaloJets.push_back(JetLVec);
153  // }
154  // // }
155  // }
156  // }
157 
158  // Fill the alphaT and HT histograms
159  if (!hltPfJets.empty()) {
160  double hltPfAlphaT = AlphaT(hltPfJets, true).value();
161  h_triggerPfAlphaT->Fill(hltPfAlphaT);
162  h_triggerPfHt->Fill(hltPfHt);
163  h_triggerPfAlphaT_triggerPfHt->Fill(hltPfHt, hltPfAlphaT);
164  }
165 
166  // if(hltCaloJets.size()>0){
167  // double hltCaloAlphaT = AlphaT(hltCaloJets,true).value();
168  // h_triggerCaloAlphaT->Fill(hltCaloAlphaT);
169  // h_triggerCaloHt->Fill(hltCaloHt);
170  // h_triggerCaloAlphaT_triggerCaloHt->Fill(hltCaloHt, hltCaloAlphaT);
171  // }
172 
173  bool hasFired = false;
174  bool hasFiredAuxiliaryForHadronicLeg = false;
175  const edm::TriggerNames &trigNames = e.triggerNames(*hltresults);
176  unsigned int numTriggers = trigNames.size();
177  for (unsigned int hltIndex = 0; hltIndex < numTriggers; ++hltIndex) {
178  if (trigNames.triggerName(hltIndex).find(triggerPath_) != std::string::npos && hltresults->wasrun(hltIndex) &&
179  hltresults->accept(hltIndex))
180  hasFired = true;
181  if (trigNames.triggerName(hltIndex).find(triggerPathAuxiliaryForHadronic_) != std::string::npos &&
182  hltresults->wasrun(hltIndex) && hltresults->accept(hltIndex))
183  hasFiredAuxiliaryForHadronicLeg = true;
184  }
185 
186  if (hasFiredAuxiliaryForHadronicLeg) {
187  float pfHT = 0.0;
188  std::vector<LorentzV> pfJets;
189  for (reco::PFJetCollection::const_iterator i_pfjet = pfJetCollection->begin(); i_pfjet != pfJetCollection->end();
190  ++i_pfjet) {
191  if (i_pfjet->pt() < ptThrJet_)
192  continue;
193  if (fabs(i_pfjet->eta()) > etaThrJet_)
194  continue;
195  pfHT += i_pfjet->pt();
196  LorentzV JetLVec(i_pfjet->pt(), i_pfjet->eta(), i_pfjet->phi(), i_pfjet->mass());
197  pfJets.push_back(JetLVec);
198  }
199 
200  // //Make the gen Calo HT and AlphaT
201  // float caloHT = 0.0;
202  // std::vector<LorentzV> caloJets;
203  // for (reco::CaloJetCollection::const_iterator i_calojet =
204  // caloJetCollection->begin(); i_calojet != caloJetCollection->end();
205  // ++i_calojet){
206  // if (i_calojet->pt() < ptThrJet_) continue;
207  // if (fabs(i_calojet->eta()) > etaThrJet_) continue;
208  // caloHT += i_calojet->pt();
209  // LorentzV
210  // JetLVec(i_calojet->pt(),i_calojet->eta(),i_calojet->phi(),i_calojet->mass());
211  // caloJets.push_back(JetLVec);
212  // }
213 
214  // AlphaT aT = AlphaT(jets);
215  // double caloAlphaT = AlphaT(caloJets).value();
216  double pfAlphaT = AlphaT(pfJets).value();
217 
218  // Fill the turnons
219  if (hasFired) {
220  // if(caloHT>caloHtThrTurnon_) h_caloAlphaTTurnOn_num-> Fill(caloAlphaT);
221  // if(caloAlphaT>caloAlphaTThrTurnon_) h_caloHtTurnOn_num-> Fill(caloHT);
222 
223  if (pfHT > pfHtThrTurnon_)
224  h_pfAlphaTTurnOn_num->Fill(pfAlphaT);
225  if (pfAlphaT > pfAlphaTThrTurnon_)
226  h_pfHtTurnOn_num->Fill(pfHT);
227  }
228  // if(caloHT>caloHtThrTurnon_) h_caloAlphaTTurnOn_den-> Fill(caloAlphaT);
229  // if(caloAlphaT>caloAlphaTThrTurnon_) h_caloHtTurnOn_den-> Fill(caloHT);
230 
231  if (pfHT > pfHtThrTurnon_)
232  h_pfAlphaTTurnOn_den->Fill(pfAlphaT);
233  if (pfAlphaT > pfAlphaTThrTurnon_)
234  h_pfHtTurnOn_den->Fill(pfHT);
235  }
236 }
237 
239  ibooker_.cd();
240  ibooker_.setCurrentFolder("HLT/SUSYBSM/" + triggerPath_);
241 
242  // offline quantities
243 
244  // online quantities
245  // h_triggerCaloHt = ibooker_.book1D("triggerCaloHt", "Trigger Calo Ht; HT
246  // (GeV)", 60, 0.0, 1500.0); h_triggerCaloAlphaT =
247  // ibooker_.book1D("triggerCaloAlphaT", "Trigger Calo AlphaT; AlphaT", 80,
248  // 0., 1.0); h_triggerCaloAlphaT_triggerCaloHt =
249  // ibooker_.book2D("triggerCaloAlphaT_triggerCaloHt","Trigger Calo HT vs
250  // Trigger Calo AlphaT; HT (GeV); AlphaT", 60,0.0,1500.,80,0.,1.0);
251  h_triggerPfHt = ibooker_.book1D("triggerPfHt", "Trigger PF Ht; HT (GeV)", 60, 0.0, 1500.0);
252  h_triggerPfAlphaT = ibooker_.book1D("triggerPfAlphaT", "Trigger PF AlphaT; AlphaT", 80, 0., 1.0);
253  h_triggerPfAlphaT_triggerPfHt = ibooker_.book2D("triggerPfAlphaT_triggerPfHt",
254  "Trigger PF HT vs Trigger PF AlphaT; HT (GeV); AlphaT",
255  60,
256  0.0,
257  1500.,
258  80,
259  0.,
260  1.0);
261 
262  // num and den hists to be divided in harvesting step to make turn on curves
263  // h_caloAlphaTTurnOn_num = ibooker_.book1D("caloAlphaTTurnOn_num", "Calo
264  // AlphaT Turn On Numerator; AlphaT", 40, 0.0, 1.0 ); h_caloAlphaTTurnOn_den =
265  // ibooker_.book1D("caloAlphaTTurnOn_den", "Calo AlphaT Turn OnDenominator;
266  // AlphaT", 40, 0.0, 1.0 ); h_caloHtTurnOn_num =
267  // ibooker_.book1D("caloHtTurnOn_num", "Calo HT Turn On Numerator; HT (GeV)",
268  // 30, 0.0, 1500.0 ); h_caloHtTurnOn_den = ibooker_.book1D("caloHtTurnOn_den",
269  // "Calo HT Turn On Denominator; HT (GeV)", 30, 0.0, 1500.0 );
270 
271  h_pfAlphaTTurnOn_num = ibooker_.book1D("pfAlphaTTurnOn_num", "PF AlphaT Turn On Numerator; AlphaT", 40, 0.0, 1.0);
272  h_pfAlphaTTurnOn_den = ibooker_.book1D("pfAlphaTTurnOn_den", "PF AlphaT Turn OnDenominator; AlphaT", 40, 0.0, 1.0);
273  h_pfHtTurnOn_num = ibooker_.book1D("pfHtTurnOn_num", "PF HT Turn On Numerator; HT (GeV)", 30, 0.0, 1500.0);
274  h_pfHtTurnOn_den = ibooker_.book1D("pfHtTurnOn_den", "PF HT Turn On Denominator; HT (GeV)", 30, 0.0, 1500.0);
275 
276  ibooker_.cd();
277 }
278 
279 // define this as a plug-in
std::size_t size() const
Definition: TriggerNames.cc:59
void dqmBeginRun(edm::Run const &, edm::EventSetup const &) override
edm::InputTag triggerFilter_
T getUntrackedParameter(std::string const &, T const &) const
Definition: AlphaT.h:8
MonitorElement * h_triggerPfHt
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:32
edm::EDGetTokenT< trigger::TriggerEvent > theTrigSummary_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
float phi() const
Definition: TriggerObject.h:54
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
const std::vector< std::string > & triggerNames() const
names of trigger paths
MonitorElement * h_triggerPfAlphaT_triggerPfHt
float eta() const
Definition: TriggerObject.h:53
Log< level::Error, false > LogError
~SUSY_HLT_alphaT() override
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
edm::TriggerNames const & triggerNames(edm::TriggerResults const &triggerResults) const override
Definition: Event.cc:275
void Fill(long long x)
Single trigger physics object (e.g., an isolated muon)
Definition: TriggerObject.h:21
tuple pfJetCollection
void bookHistos(DQMStore::IBooker &)
MonitorElement * h_pfAlphaTTurnOn_num
edm::EDGetTokenT< edm::TriggerResults > triggerResults_
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
bool isValid() const
Definition: HandleBase.h:70
ROOT::Math::LorentzVector< ROOT::Math::PtEtaPhiM4D< double > > LorentzV
HLTConfigProvider fHltConfig
std::vector< TriggerObject > TriggerObjectCollection
collection of trigger physics objects (e.g., all isolated muons)
Definition: TriggerObject.h:75
Log< level::Info, false > LogInfo
static const char *const trigNames[]
Definition: EcalDumpRaw.cc:57
void analyze(edm::Event const &e, edm::EventSetup const &eSetup) override
std::string triggerPathAuxiliaryForHadronic_
std::string const & triggerName(unsigned int index) const
Definition: TriggerNames.cc:50
std::vector< size_type > Keys
MonitorElement * h_pfHtTurnOn_den
bool init(const edm::Run &iRun, const edm::EventSetup &iSetup, const std::string &processName, bool &changed)
d&#39;tor
MonitorElement * book2D(TString const &name, TString const &title, int nchX, double lowX, double highX, int nchY, double lowY, double highY, FUNC onbooking=NOOP())
Definition: DQMStore.h:212
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
MonitorElement * h_triggerPfAlphaT
MonitorElement * h_pfAlphaTTurnOn_den
std::string HLTProcess_
std::string triggerPath_
MonitorElement * h_pfHtTurnOn_num
edm::EDGetTokenT< reco::PFJetCollection > thePfJetCollection_
Log< level::Warning, false > LogWarning
float mass() const
Definition: TriggerObject.h:55
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
edm::InputTag triggerPreFilter_
tuple pfJets
Definition: pfJets_cff.py:8
SUSY_HLT_alphaT(const edm::ParameterSet &ps)
double pfAlphaTThrTurnon_
double value(void) const
Definition: AlphaT.h:40
Definition: Run.h:45
#define LogDebug(id)