1 |
|
2 |
#include "FWCore/Framework/interface/EDAnalyzer.h"
|
3 |
#include "FWCore/ParameterSet/interface/ParameterSet.h"
|
4 |
#include "FWCore/ServiceRegistry/interface/Service.h"
|
5 |
|
6 |
#include <vector>
|
7 |
|
8 |
class TFileService;
|
9 |
class TH1F;
|
10 |
//class TH2F;
|
11 |
|
12 |
class SimpleTrackAnalyzer: public edm::EDAnalyzer
|
13 |
{
|
14 |
public:
|
15 |
typedef std::map<std::string,TH1F*> HistoMapTH1F;
|
16 |
|
17 |
explicit SimpleTrackAnalyzer(edm::ParameterSet const&);
|
18 |
~SimpleTrackAnalyzer();
|
19 |
|
20 |
void beginJob();
|
21 |
void analyze(edm::Event const&, edm::EventSetup const&);
|
22 |
private:
|
23 |
void bookHistos(HistoMapTH1F&, edm::Service<TFileService> const&);
|
24 |
|
25 |
static const int nVarBinsMax_;
|
26 |
|
27 |
edm::InputTag trackTag_;
|
28 |
|
29 |
double minPt_;
|
30 |
double maxPt_;
|
31 |
unsigned int nBinsPt_;
|
32 |
double minEta_;
|
33 |
double maxEta_;
|
34 |
unsigned int nBinsEta_;
|
35 |
double minPtSum_;
|
36 |
double maxPtSum_;
|
37 |
unsigned int nBinsPtSum_;
|
38 |
unsigned int minNTracks_;
|
39 |
unsigned int maxNTracks_;
|
40 |
unsigned int nBinsNTracks_;
|
41 |
|
42 |
std::vector<double> boundariesPt;
|
43 |
std::vector<double> boundariesEta;
|
44 |
std::vector<double> boundariesPtSum;
|
45 |
std::vector<double> boundariesNTracks;
|
46 |
|
47 |
HistoMapTH1F histosTH1F_;
|
48 |
//HistoMapTH2F histosTH2F_;
|
49 |
};
|
50 |
|
51 |
#include "FWCore/Framework/interface/Frameworkfwd.h"
|
52 |
#include "FWCore/Framework/interface/Event.h"
|
53 |
#include "FWCore/Framework/interface/EDAnalyzer.h"
|
54 |
#include "FWCore/Framework/interface/MakerMacros.h"
|
55 |
#include "FWCore/ParameterSet/interface/ParameterSet.h"
|
56 |
#include "FWCore/MessageLogger/interface/MessageLogger.h"
|
57 |
|
58 |
#include "DataFormats/TrackReco/interface/TrackFwd.h"
|
59 |
#include "DataFormats/TrackReco/interface/Track.h"
|
60 |
|
61 |
#include "FWCore/ServiceRegistry/interface/Service.h"
|
62 |
#include "CommonTools/UtilAlgos/interface/TFileService.h"
|
63 |
|
64 |
#include "TTree.h"
|
65 |
#include "TH1F.h"
|
66 |
#include "TH2F.h"
|
67 |
|
68 |
const int SimpleTrackAnalyzer::nVarBinsMax_ = 1000;
|
69 |
|
70 |
SimpleTrackAnalyzer::SimpleTrackAnalyzer(edm::ParameterSet const& pset):
|
71 |
trackTag_(pset.getParameter<edm::InputTag>("TrackTag")),
|
72 |
minPt_(pset.getParameter<double>("MinPt")),
|
73 |
maxPt_(pset.getParameter<double>("MaxPt")),
|
74 |
nBinsPt_(pset.getParameter<unsigned int>("NBinsPt")),
|
75 |
minEta_(pset.getParameter<double>("MinEta")),
|
76 |
maxEta_(pset.getParameter<double>("MaxEta")),
|
77 |
nBinsEta_(pset.getParameter<unsigned int>("NBinsEta")),
|
78 |
minPtSum_(pset.getParameter<double>("MinPtSum")),
|
79 |
maxPtSum_(pset.getParameter<double>("MaxPtSum")),
|
80 |
nBinsPtSum_(pset.getParameter<unsigned int>("NBinsPtSum")),
|
81 |
minNTracks_(pset.getParameter<unsigned int>("MinNTracks")),
|
82 |
maxNTracks_(pset.getParameter<unsigned int>("MaxNTracks")),
|
83 |
nBinsNTracks_(pset.getParameter<unsigned int>("NBinsNTracks")) {
|
84 |
if( pset.exists("VarBin") ){
|
85 |
edm::ParameterSet const& psetVarBin = pset.getParameter<edm::ParameterSet>("VarBin");
|
86 |
if( psetVarBin.exists("pt") )
|
87 |
boundariesPt = psetVarBin.getParameter<std::vector<double> >("pt");
|
88 |
if( psetVarBin.exists("eta") )
|
89 |
boundariesEta = psetVarBin.getParameter<std::vector<double> >("eta");
|
90 |
if( psetVarBin.exists("ptSum") )
|
91 |
boundariesPtSum = psetVarBin.getParameter<std::vector<double> >("ptSum");
|
92 |
if( psetVarBin.exists("nTracks") )
|
93 |
boundariesNTracks = psetVarBin.getParameter<std::vector<double> >("nTracks");
|
94 |
}
|
95 |
}
|
96 |
|
97 |
SimpleTrackAnalyzer::~SimpleTrackAnalyzer(){}
|
98 |
|
99 |
void SimpleTrackAnalyzer::beginJob(){
|
100 |
edm::Service<TFileService> fs;
|
101 |
|
102 |
bool sumw2 = TH1::GetDefaultSumw2();
|
103 |
TH1::SetDefaultSumw2(false);
|
104 |
|
105 |
bookHistos(histosTH1F_,fs);
|
106 |
//bookHistos(histosTH2F_,fs);
|
107 |
for(HistoMapTH1F::const_iterator it = histosTH1F_.begin(); it != histosTH1F_.end(); ++it) it->second->Sumw2();
|
108 |
|
109 |
TH1::SetDefaultSumw2(sumw2);
|
110 |
}
|
111 |
|
112 |
void SimpleTrackAnalyzer::analyze(const edm::Event& event, const edm::EventSetup& setup){
|
113 |
edm::Handle<edm::View<reco::Track> > trackCollectionH;
|
114 |
event.getByLabel(trackTag_,trackCollectionH);
|
115 |
const edm::View<reco::Track>& trackColl = *(trackCollectionH.product());
|
116 |
|
117 |
int nTracks = 0;
|
118 |
double ptSum = 0.;
|
119 |
edm::View<reco::Track>::const_iterator track = trackColl.begin();
|
120 |
edm::View<reco::Track>::const_iterator tracks_end = trackColl.end();
|
121 |
for(; track != tracks_end; ++track){
|
122 |
// Any selection here
|
123 |
|
124 |
++nTracks;
|
125 |
ptSum += track->pt();
|
126 |
histosTH1F_["TrackPt"]->Fill( track->pt() );
|
127 |
histosTH1F_["TrackEta"]->Fill( track->eta() );
|
128 |
if(histosTH1F_.find("TrackPtVarBin") != histosTH1F_.end())
|
129 |
histosTH1F_["TrackPtVarBin"]->Fill( track->pt() );
|
130 |
if(histosTH1F_.find("TrackEtaVarBin") != histosTH1F_.end())
|
131 |
histosTH1F_["TrackEtaVarBin"]->Fill( track->eta() );
|
132 |
|
133 |
}
|
134 |
histosTH1F_["TrackPtSum"]->Fill( ptSum );
|
135 |
histosTH1F_["NTracks"]->Fill( nTracks );
|
136 |
if(histosTH1F_.find("TrackPtSumVarBin") != histosTH1F_.end())
|
137 |
histosTH1F_["TrackPtSumVarBin"]->Fill( ptSum );
|
138 |
if(histosTH1F_.find("NTracksVarBin") != histosTH1F_.end())
|
139 |
histosTH1F_["NTracksVarBin"]->Fill( nTracks );
|
140 |
}
|
141 |
|
142 |
void SimpleTrackAnalyzer::bookHistos(HistoMapTH1F& histos, edm::Service<TFileService> const& fs){
|
143 |
histos["TrackPt"] = fs->make<TH1F>("TrackPt","TrackPt",nBinsPt_,minPt_,maxPt_);
|
144 |
histos["TrackEta"] = fs->make<TH1F>("TrackEta","TrackEta",nBinsEta_,minEta_,maxEta_);
|
145 |
histos["TrackPtSum"] = fs->make<TH1F>("TrackPtSum","TrackPtSum",nBinsPtSum_,minPtSum_,maxPtSum_);
|
146 |
histos["NTracks"] = fs->make<TH1F>("NTracks","NTracks",nBinsNTracks_,minNTracks_,maxNTracks_);
|
147 |
|
148 |
if( boundariesPt.size() ){
|
149 |
float ptBins[nVarBinsMax_];
|
150 |
std::copy(boundariesPt.begin(),boundariesPt.end(),ptBins);
|
151 |
histos["TrackPtVarBin"] = fs->make<TH1F>("TrackPtVarBin","TrackPtVarBin",(boundariesPt.size() - 1),ptBins);
|
152 |
}
|
153 |
if( boundariesEta.size() ){
|
154 |
float etaBins[nVarBinsMax_];
|
155 |
std::copy(boundariesEta.begin(),boundariesEta.end(),etaBins);
|
156 |
histos["TrackEtaVarBin"] = fs->make<TH1F>("TrackEtaVarBin","TrackEtaVarBin",(boundariesEta.size() - 1),etaBins);
|
157 |
}
|
158 |
if( boundariesPtSum.size() ){
|
159 |
float ptSumBins[nVarBinsMax_];
|
160 |
std::copy(boundariesPtSum.begin(),boundariesPtSum.end(),ptSumBins);
|
161 |
histos["TrackPtSumVarBin"] = fs->make<TH1F>("TrackPtSumVarBin","TrackPtSumVarBin",(boundariesPtSum.size() - 1),ptSumBins);
|
162 |
}
|
163 |
if( boundariesNTracks.size() ){
|
164 |
float nTracksBins[nVarBinsMax_];
|
165 |
std::copy(boundariesNTracks.begin(),boundariesNTracks.end(),nTracksBins);
|
166 |
histos["NTracksVarBin"] = fs->make<TH1F>("NTracksVarBin","NTracksVarBin",(boundariesNTracks.size() - 1),nTracksBins);
|
167 |
}
|
168 |
}
|
169 |
|
170 |
DEFINE_FWK_MODULE(SimpleTrackAnalyzer);
|