ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/VHbbAnalysis/HbbAnalyzer/plugins/HLTInfoDumperGeneral.cc
Revision: 1.3
Committed: Mon Sep 19 14:46:12 2011 UTC (13 years, 7 months ago) by tboccali
Content type: text/plain
Branch: MAIN
CVS Tags: EDMV42_Step2_V6, EDMV42_Step2_V5a, EDMV42_Step2_V5, tauCandV42, hbbsubstructDev_11, hbbsubstructDev_10, hbbsubstructDev_9, hbbsubstructDev_8, hbbsubstructDev_7, hbbsubstructDev_6, hbbsubstructDev_5, hbbsubstructDev_4, hbbsubstructDev_3, hbbsubstructDev_2, hbbsubstructDev_1, hbbsubstructDev, V21TauCand_0, TauCandidates_0, EDMV42_Step2_V4a, EDMV42_Step2_V4, EDMV42_Step2_V3, EDMV42_Step2_V2, EDMV42_Step2_V1, EdmV42, EdmV41alpha1, EdmV40alpha1, EdmV40alpha, V21emuCand, EdmV33Jun12v2_consistent, Step2ForV33_v2, Step2ForV33_v1, EdmV33Jun12v2, EdmV33Jun12v1, EdmV33Jun12v0, Step2ForV32_v2, Step2ForV32_v0, Step2ForV31_v0, EdmV32May24v0, EdmV31May21v1, EdmV31May17v0, May14thStep2, EdmV30Apr10, EdmV21Apr10v2, EdmV22May9, EdmV21Apr06, EdmV21Apr10, EdmV21Apr04, EdmV21Apr03, EdmV21Apr2, EdmV21Mar30, EdmV20Mar12, EdmV11Oct2011_fixMET, EdmV11Oct2011, EdmV10Oct2011, EdmV9Sept2011, Sept19th2011_2, HEAD
Branch point for: V42TauCandidate, hbbsubstructDevPostHCP, V21TauCand, TauCandidatesV21, V21emuCandidate
Changes since 1.2: +1 -1 lines
Log Message:
add souvik s code

File Contents

# User Rev Content
1 tboccali 1.1 #ifndef HLT_FROM_SOUVIK_H
2     #define HLT_FROM_SOUVIK_H
3     #include "FWCore/Framework/interface/Frameworkfwd.h"
4     #include "FWCore/Framework/interface/Event.h"
5     #include "FWCore/Framework/interface/EventSetup.h"
6     #include "FWCore/Framework/interface/EDProducer.h"
7     #include "FWCore/Utilities/interface/InputTag.h"
8 tboccali 1.3 #include "FWCore/ParameterSet/interface/ParameterSet.h"
9 tboccali 1.1 #include "DataFormats/Common/interface/Handle.h"
10    
11     #include "FWCore/ServiceRegistry/interface/Service.h"
12     #include "CommonTools/UtilAlgos/interface/TFileService.h"
13    
14     #include "FWCore/Common/interface/TriggerNames.h"
15     #include "DataFormats/Common/interface/TriggerResults.h"
16     #include "DataFormats/HLTReco/interface/TriggerEvent.h"
17    
18     #include <algorithm>
19     #include <iostream>
20     #include <cmath>
21     #include <vector>
22     #include <string>
23    
24     using namespace edm;
25     using namespace std;
26     using namespace reco;
27    
28    
29     class HLTInfoDumperGeneral: public edm::EDProducer
30     {
31     public:
32     HLTInfoDumperGeneral(const edm::ParameterSet&);
33    
34     private:
35     void produce(edm::Event&, const edm::EventSetup&);
36     std::string hltPath_;
37     edm::InputTag trigTag_, trigEvt_;
38     int nFilters_;
39     typedef std::vector<edm::InputTag> VInputTag;
40     VInputTag filterNames_;
41     };
42    
43     HLTInfoDumperGeneral::HLTInfoDumperGeneral(const ParameterSet& cfg)
44     {
45     hltPath_=cfg.getUntrackedParameter<std::string>("HLTPath");
46     trigTag_=cfg.getUntrackedParameter<edm::InputTag>("TriggerResults");
47     trigEvt_=cfg.getUntrackedParameter<edm::InputTag>("TriggerEvent");
48     filterNames_=cfg.getUntrackedParameter<VInputTag>("FilterNames");
49    
50     produces<unsigned int>("Passed").setBranchAlias("Passed");
51     nFilters_=filterNames_.size();
52     for (int i=0; i<nFilters_; ++i)
53     {
54     std::string filterName=filterNames_.at(i).label();
55     std::string pT_string=filterName+"PT";
56     std::string phi_string=filterName+"Phi";
57     std::string eta_string=filterName+"Eta";
58     std::string mass_string=filterName+"Mass";
59     std::string id_string=filterName+"ID";
60     produces<std::vector<float> >(pT_string).setBranchAlias(pT_string);
61     produces<std::vector<float> >(phi_string).setBranchAlias(phi_string);
62     produces<std::vector<float> >(eta_string).setBranchAlias(eta_string);
63     produces<std::vector<float> >(mass_string).setBranchAlias(mass_string);
64     produces<std::vector<int> >(id_string).setBranchAlias(id_string);
65     }
66     }
67    
68     void HLTInfoDumperGeneral::produce(Event &evt, const EventSetup &iSetup)
69     {
70     auto_ptr<unsigned int> passed(new unsigned int); *passed=0;
71     std::vector<float> objPT;
72     std::vector<float> objPhi;
73     std::vector<float> objEta;
74     std::vector<float> objMass;
75     std::vector<int> objID;
76     std::map<std::string, std::vector<float> > objPT_filter;
77     std::map<std::string, std::vector<float> > objPhi_filter;
78     std::map<std::string, std::vector<float> > objEta_filter;
79     std::map<std::string, std::vector<float> > objMass_filter;
80     std::map<std::string, std::vector<int> > objID_filter;
81    
82     // Extract trigger result information
83     Handle<TriggerResults> triggerResults;
84     if (!evt.getByLabel(trigTag_, triggerResults))
85     {
86     std::cout<<"TriggerResults does not exist"<<std::endl;
87     return;
88     }
89     evt.getByLabel(trigTag_, triggerResults);
90     const edm::TriggerNames &trigNames = evt.triggerNames(*triggerResults);
91     for (size_t i=0; i<triggerResults->size(); ++i)
92     {
93     std::string trigName = trigNames.triggerName(i);
94     if (trigName.find(hltPath_)==0)
95     {
96     if (triggerResults->accept(i))
97     {
98     *passed=1;
99     }
100     break;
101     }
102     }
103    
104     objPT_filter.clear();
105     objPhi_filter.clear();
106     objEta_filter.clear();
107     objMass_filter.clear();
108     objID_filter.clear();
109     if (*passed==1)
110     {
111     // Now extract triggger object information
112     Handle<trigger::TriggerEvent> triggerEvent;
113     if (!evt.getByLabel(trigEvt_, triggerEvent))
114     {
115     std::cout<<"TriggerEvent does not exist"<<std::endl;
116     return;
117     }
118     evt.getByLabel(trigEvt_, triggerEvent);
119     const trigger::TriggerObjectCollection &toc=triggerEvent->getObjects();
120    
121     // Now extract 4-vectors for all the filter objects
122     for (int i=0; i<nFilters_; ++i)
123     {
124     edm::InputTag filterTag=filterNames_.at(i);
125     int filter_index=triggerEvent->filterIndex(filterTag);
126     if (filter_index!=triggerEvent->sizeFilters())
127     {
128     const trigger::Keys& keys=triggerEvent->filterKeys(filter_index);
129     objPT.clear();
130     objPhi.clear();
131     objEta.clear();
132     objMass.clear();
133     objID.clear();
134     for (unsigned int k=0; k<keys.size(); k++)
135     {
136     objPT.push_back(toc[keys[k]].pt());
137     objPhi.push_back(toc[keys[k]].phi());
138     objEta.push_back(toc[keys[k]].eta());
139     objMass.push_back(toc[keys[k]].mass());
140     objID.push_back(toc[keys[k]].id());
141     }
142     std::string filterName=filterTag.label();
143     std::string pT_string=filterName+"PT";
144     std::string phi_string=filterName+"Phi";
145     std::string eta_string=filterName+"Eta";
146     std::string mass_string=filterName+"Mass";
147     std::string id_string=filterName+"ID";
148     objPT_filter.insert(std::pair<std::string, std::vector<float> >(pT_string, objPT));
149     objPhi_filter.insert(std::pair<std::string, std::vector<float> >(phi_string, objPhi));
150     objEta_filter.insert(std::pair<std::string, std::vector<float> >(eta_string, objEta));
151     objMass_filter.insert(std::pair<std::string, std::vector<float> >(mass_string, objMass));
152     objID_filter.insert(std::pair<std::string, std::vector<int> >(id_string, objID));
153     }
154     else
155     {
156     std::cout<<"Filter name "<<filterTag<<" not found"<<std::endl;
157     }
158     }
159     }
160    
161     evt.put(passed, "Passed");
162     for (std::map<std::string, std::vector<float> >::iterator i=objPT_filter.begin(); i!=objPT_filter.end(); ++i)
163     {
164     auto_ptr<std::vector<float> >objPT_ptr(new std::vector<float>);
165     for (unsigned j=0; j<(i->second).size(); ++j) objPT_ptr->push_back((i->second).at(j));
166     evt.put(objPT_ptr, i->first);
167     }
168     for (std::map<std::string, std::vector<float> >::iterator i=objPhi_filter.begin(); i!=objPhi_filter.end(); ++i)
169     {
170     auto_ptr<std::vector<float> >objPhi_ptr(new std::vector<float>);
171     for (unsigned j=0; j<(i->second).size(); ++j) objPhi_ptr->push_back((i->second).at(j));
172     evt.put(objPhi_ptr, i->first);
173     }
174     for (std::map<std::string, std::vector<float> >::iterator i=objEta_filter.begin(); i!=objEta_filter.end(); ++i)
175     {
176     auto_ptr<std::vector<float> >objEta_ptr(new std::vector<float>);
177     for (unsigned j=0; j<(i->second).size(); ++j) objEta_ptr->push_back((i->second).at(j));
178     evt.put(objEta_ptr, i->first);
179     }
180     for (std::map<std::string, std::vector<float> >::iterator i=objMass_filter.begin(); i!=objMass_filter.end(); ++i)
181     {
182     auto_ptr<std::vector<float> >objMass_ptr(new std::vector<float>);
183     for (unsigned j=0; j<(i->second).size(); ++j) objMass_ptr->push_back((i->second).at(j));
184     evt.put(objMass_ptr, i->first);
185     }
186     for (std::map<std::string, std::vector<int> >::iterator i=objID_filter.begin(); i!=objID_filter.end(); ++i)
187     {
188     auto_ptr<std::vector<int> >objID_ptr(new std::vector<int>);
189     for (unsigned j=0; j<(i->second).size(); ++j) objID_ptr->push_back((i->second).at(j));
190     evt.put(objID_ptr, i->first);
191     }
192    
193     }
194    
195     #include "FWCore/Framework/interface/MakerMacros.h"
196    
197     DEFINE_FWK_MODULE( HLTInfoDumperGeneral );
198    
199     #endif