1 |
// -*- C++ -*-
|
2 |
//
|
3 |
// Package: FinalPlots
|
4 |
// Class: FinalPlots
|
5 |
//
|
6 |
/**\class FinalPlots FinalPlots.cc RA2/FinalPlots/src/FinalPlots.cc
|
7 |
|
8 |
Description: <one line class summary>
|
9 |
|
10 |
Implementation:
|
11 |
<Notes on implementation>
|
12 |
*/
|
13 |
//
|
14 |
// Original Author: Christian Autermann,68/112,2115,
|
15 |
// Created: Sun Oct 18 20:00:45 CEST 2009
|
16 |
// $Id: FinalPlots.cc,v 1.7 2010/01/06 13:42:37 auterman Exp $
|
17 |
//
|
18 |
//
|
19 |
|
20 |
|
21 |
// system include files
|
22 |
#include <memory>
|
23 |
#include <string>
|
24 |
#include <cmath>
|
25 |
#include <numeric>
|
26 |
|
27 |
// user include files
|
28 |
#include "FWCore/Framework/interface/Frameworkfwd.h"
|
29 |
#include "FWCore/Framework/interface/EDAnalyzer.h"
|
30 |
#include "FWCore/Framework/interface/Event.h"
|
31 |
#include "FWCore/Framework/interface/MakerMacros.h"
|
32 |
#include "FWCore/ParameterSet/interface/ParameterSet.h"
|
33 |
#include "DataFormats/Math/interface/LorentzVector.h"
|
34 |
#include "DataFormats/Math/interface/deltaPhi.h"
|
35 |
#include "DataFormats/METReco/interface/CaloMET.h"
|
36 |
#include "DataFormats/METReco/interface/MET.h"
|
37 |
#include "DataFormats/JetReco/interface/CaloJet.h"
|
38 |
#include "DataFormats/JetReco/interface/PFJet.h"
|
39 |
#include "DataFormats/Candidate/interface/Candidate.h"
|
40 |
#include "FWCore/ServiceRegistry/interface/Service.h"
|
41 |
#include "PhysicsTools/UtilAlgos/interface/TFileService.h"
|
42 |
|
43 |
#include "TH1F.h"
|
44 |
//
|
45 |
// class decleration
|
46 |
//
|
47 |
|
48 |
class FinalPlots : public edm::EDAnalyzer {
|
49 |
public:
|
50 |
explicit FinalPlots(const edm::ParameterSet&);
|
51 |
~FinalPlots();
|
52 |
|
53 |
|
54 |
private:
|
55 |
virtual void beginJob(const edm::EventSetup&);
|
56 |
virtual void analyze(const edm::Event&, const edm::EventSetup&);
|
57 |
virtual void endJob() ;
|
58 |
double dPhiBiased(const edm::View<reco::Candidate>*jets,
|
59 |
const edm::View<reco::Candidate>::const_iterator jet);
|
60 |
|
61 |
|
62 |
typedef ROOT::Math::LorentzVector<ROOT::Math::PxPyPzE4D<double> > LorentzV;
|
63 |
std::vector<double> deltaSumPt_permutations(const std::vector<LorentzV>& p4s);
|
64 |
double alphaT(const std::vector<LorentzV>& p4s);
|
65 |
|
66 |
// ----------member data ---------------------------
|
67 |
edm::InputTag weightName_;
|
68 |
edm::InputTag Jet_;
|
69 |
edm::InputTag Met_;
|
70 |
std::string name_;
|
71 |
double weight_, jetptmin_, jetetamax_;
|
72 |
|
73 |
TH1F * HT_;
|
74 |
TH1F * MHT_;
|
75 |
TH1F * MET_;
|
76 |
TH1F * Jet1Pt_, *Jet2Pt_, *Jet3Pt_, *Jet4Pt_;
|
77 |
TH1F * JetMult_;
|
78 |
TH1F * alphaT_;
|
79 |
TH1F * dPhiMin_;
|
80 |
TH1F * dPhiBiased_;
|
81 |
|
82 |
};
|
83 |
|
84 |
FinalPlots::FinalPlots(const edm::ParameterSet& iConfig):
|
85 |
weightName_(iConfig.getParameter<edm::InputTag>("weightName") ),
|
86 |
Jet_( iConfig.getParameter<edm::InputTag>("Jet") ),
|
87 |
Met_( iConfig.getParameter<edm::InputTag>("MET") ),
|
88 |
name_( iConfig.getParameter<std::string>("uncertainty_name") ),
|
89 |
jetptmin_ ( iConfig.getParameter<double>("JetPtMin") ),
|
90 |
jetetamax_( iConfig.getParameter<double>("JetEtaMax") )
|
91 |
{
|
92 |
}
|
93 |
|
94 |
// ------------ method called once each job just before starting event loop ------------
|
95 |
void
|
96 |
FinalPlots::beginJob(const edm::EventSetup&)
|
97 |
{
|
98 |
edm::Service<TFileService> fs;
|
99 |
if( !fs ){
|
100 |
throw edm::Exception( edm::errors::Configuration,
|
101 |
"TFile Service is not registered in cfg file" );
|
102 |
}
|
103 |
|
104 |
std::string histname = "HT"+name_;
|
105 |
HT_ = fs->make<TH1F>(histname.c_str(),";HT [GeV];events", 200, 0.0, 2000.0);
|
106 |
HT_->Sumw2();
|
107 |
|
108 |
histname = "MHT"+name_;
|
109 |
MHT_ = fs->make<TH1F>(histname.c_str(),";MHT [GeV];events", 200, 0.0, 1000.0);
|
110 |
MHT_->Sumw2();
|
111 |
|
112 |
histname = "MET"+name_;
|
113 |
MET_ = fs->make<TH1F>(histname.c_str(),";MET [GeV];events", 200, 0.0, 1000.0);
|
114 |
MET_->Sumw2();
|
115 |
|
116 |
histname = "Jet1Pt"+name_;
|
117 |
Jet1Pt_ = fs->make<TH1F>(histname.c_str(),";1. Jet Pt [GeV];events", 200, 0.0, 1000.0);
|
118 |
Jet1Pt_->Sumw2();
|
119 |
|
120 |
histname = "Jet2Pt"+name_;
|
121 |
Jet2Pt_ = fs->make<TH1F>(histname.c_str(),";2. Jet Pt [GeV];events", 200, 0.0, 1000.0);
|
122 |
Jet2Pt_->Sumw2();
|
123 |
|
124 |
histname = "Jet3Pt"+name_;
|
125 |
Jet3Pt_ = fs->make<TH1F>(histname.c_str(),";3. Jet Pt [GeV];events", 200, 0.0, 1000.0);
|
126 |
Jet3Pt_->Sumw2();
|
127 |
|
128 |
histname = "Jet4Pt"+name_;
|
129 |
Jet4Pt_ = fs->make<TH1F>(histname.c_str(),";4. Jet Pt [GeV];events", 200, 0.0, 1000.0);
|
130 |
Jet4Pt_->Sumw2();
|
131 |
|
132 |
histname = "JetMult"+name_;
|
133 |
JetMult_ = fs->make<TH1F>(histname.c_str(),";Jet multiplicity;events", 21, -0.5, 20.5);
|
134 |
JetMult_->Sumw2();
|
135 |
|
136 |
histname = "alphaT"+name_;
|
137 |
alphaT_ = fs->make<TH1F>(histname.c_str(),";alphaT;events", 200, 0.0, 2.0);
|
138 |
alphaT_->Sumw2();
|
139 |
|
140 |
histname = "dPhiMin"+name_;
|
141 |
dPhiMin_ = fs->make<TH1F>(histname.c_str(),";#Delta#phi_{min};events", 200, 0.0, 3.2);
|
142 |
dPhiMin_->Sumw2();
|
143 |
|
144 |
histname = "dPhiBiased"+name_;
|
145 |
dPhiBiased_ = fs->make<TH1F>(histname.c_str(),";#Delta#phi_{biased};events", 200, 0.0, 3.2);
|
146 |
dPhiBiased_->Sumw2();
|
147 |
}
|
148 |
|
149 |
|
150 |
FinalPlots::~FinalPlots()
|
151 |
{
|
152 |
|
153 |
// do anything here that needs to be done at desctruction time
|
154 |
// (e.g. close files, deallocate resources etc.)
|
155 |
|
156 |
}
|
157 |
|
158 |
|
159 |
//
|
160 |
// member functions
|
161 |
//
|
162 |
// ------------ biased delta phi ----------------------------
|
163 |
double FinalPlots::dPhiBiased(const edm::View<reco::Candidate> * jets, const edm::View<reco::Candidate>::const_iterator excljet)
|
164 |
{
|
165 |
math::PtEtaPhiMLorentzVector htvec(0.0, 0.0, 0.0, 0.0);
|
166 |
for (edm::View<reco::Candidate>::const_iterator jet=jets->begin();
|
167 |
jet!=jets->end(); ++jet){
|
168 |
if ( fabs(jet->eta())>jetetamax_ ) continue;
|
169 |
if ( jet->pt()<jetptmin_ ) continue;
|
170 |
if ( jet==excljet ) continue;
|
171 |
htvec+=jet->p4();
|
172 |
}
|
173 |
return deltaPhi(htvec, *excljet);
|
174 |
}
|
175 |
|
176 |
std::vector<double> FinalPlots::deltaSumPt_permutations(const std::vector<LorentzV>& p4s)
|
177 |
{
|
178 |
std::vector<std::vector<double> > ht( 1<<(p4s.size()-1) , std::vector<double>( 2, 0.) );
|
179 |
for(unsigned i=0; i < ht.size(); i++) {
|
180 |
for(unsigned j=0; j < p4s.size(); j++) {
|
181 |
ht [i] [(i/(1<<j))%2] += p4s[j].pt();
|
182 |
}
|
183 |
}
|
184 |
std::vector<double> deltaHT; for(unsigned i=0; i<ht.size(); i++) deltaHT.push_back(fabs(ht[i][0]-ht[i][1]));
|
185 |
return deltaHT;
|
186 |
}
|
187 |
|
188 |
double FinalPlots::alphaT(const std::vector<LorentzV>& p4s)
|
189 |
{
|
190 |
if(p4s.size()<2) return 0;
|
191 |
|
192 |
std::vector<double> pTs;
|
193 |
for(unsigned i=0; i<p4s.size(); i++) pTs.push_back(p4s[i].pt());
|
194 |
for(unsigned i=0; i<p4s.size(); i++) pTs.push_back(p4s[i].pt());
|
195 |
const std::vector<double> DHTper( deltaSumPt_permutations(p4s) );
|
196 |
|
197 |
const double mDHT = *(std::min_element( DHTper.begin(), DHTper.end() ));
|
198 |
const double sumPT = accumulate( pTs.begin(), pTs.end(), double(0) );
|
199 |
const LorentzV sumP4 = accumulate( p4s.begin(), p4s.end(), LorentzV() );
|
200 |
|
201 |
return 0.5 * ( sumPT - mDHT ) / sqrt( sumPT*sumPT - sumP4.perp2() );
|
202 |
}
|
203 |
|
204 |
|
205 |
// ------------ method called to for each event ------------
|
206 |
void
|
207 |
FinalPlots::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup)
|
208 |
{
|
209 |
using namespace edm;
|
210 |
using namespace std;
|
211 |
|
212 |
edm::Handle< double > event_weight;
|
213 |
iEvent.getByLabel(weightName_, event_weight);
|
214 |
weight_ = (event_weight.isValid() ? (*event_weight) : 1.0);
|
215 |
|
216 |
edm::Handle<edm::View<reco::Candidate> > jet_hnd;
|
217 |
iEvent.getByLabel(Jet_, jet_hnd);
|
218 |
|
219 |
edm::Handle<edm::View<reco::MET> > met_hnd;
|
220 |
iEvent.getByLabel(Met_, met_hnd);
|
221 |
|
222 |
MET_->Fill( met_hnd->front().et(), weight_ );
|
223 |
|
224 |
double ht=0.0, dphimin=999., dphibiased=999.;
|
225 |
int jetmult = 0;
|
226 |
math::PtEtaPhiMLorentzVector htvec(0.0, 0.0, 0.0, 0.0);
|
227 |
std::vector <LorentzV> forAlphaT;
|
228 |
for (edm::View<reco::Candidate>::const_iterator jet=jet_hnd->begin();
|
229 |
jet!=jet_hnd->end(); ++jet){
|
230 |
if ( fabs(jet->eta())>jetetamax_ ) continue;
|
231 |
if ( jet->pt()<jetptmin_ ) continue;
|
232 |
forAlphaT.push_back(jet->p4());
|
233 |
htvec+=jet->p4();
|
234 |
ht+=jet->pt();
|
235 |
if (jetmult==0) Jet1Pt_->Fill( jet->pt(), weight_ );
|
236 |
else if (jetmult==1) Jet2Pt_->Fill( jet->pt(), weight_ );
|
237 |
else if (jetmult==2) Jet3Pt_->Fill( jet->pt(), weight_ );
|
238 |
else if (jetmult==3) Jet4Pt_->Fill( jet->pt(), weight_ );
|
239 |
++jetmult;
|
240 |
|
241 |
if ( deltaPhi(*jet, met_hnd->front()) < dphimin )
|
242 |
dphimin = deltaPhi(*jet, met_hnd->front());
|
243 |
|
244 |
double thisdphi = dPhiBiased(&(*jet_hnd), jet);
|
245 |
if (thisdphi < dphibiased)
|
246 |
dphibiased = thisdphi;
|
247 |
|
248 |
}
|
249 |
HT_ ->Fill( ht, weight_ );
|
250 |
MHT_->Fill( htvec.pt(), weight_ );
|
251 |
JetMult_->Fill( jetmult, weight_ );
|
252 |
alphaT_->Fill( alphaT(forAlphaT), weight_ );
|
253 |
dPhiMin_->Fill( dphimin, weight_ );
|
254 |
dPhiBiased_->Fill( dphibiased, weight_ );
|
255 |
|
256 |
//other variables
|
257 |
//MPT...
|
258 |
}
|
259 |
|
260 |
|
261 |
// ------------ method called once each job just after ending the event loop ------------
|
262 |
void
|
263 |
FinalPlots::endJob() {
|
264 |
}
|
265 |
|
266 |
//define this as a plug-in
|
267 |
DEFINE_FWK_MODULE(FinalPlots);
|