1 |
auterman |
1.1 |
// -*- C++ -*-
|
2 |
|
|
//
|
3 |
|
|
// Package: FinalPlots
|
4 |
|
|
// Class: FinalPlots
|
5 |
|
|
//
|
6 |
|
|
/**\class FinalPlots FinalPlots.cc Analysis/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.2 2010/01/27 13:34:08 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/PatCandidates/interface/Jet.h"
|
39 |
|
|
#include "DataFormats/JetReco/interface/PFJet.h"
|
40 |
|
|
#include "DataFormats/Candidate/interface/Candidate.h"
|
41 |
|
|
#include "FWCore/ServiceRegistry/interface/Service.h"
|
42 |
|
|
#include "PhysicsTools/UtilAlgos/interface/TFileService.h"
|
43 |
|
|
|
44 |
|
|
#include "FWCore/Framework/interface/EventSetup.h"
|
45 |
|
|
#include "Geometry/CaloGeometry/interface/CaloGeometry.h"
|
46 |
|
|
#include "Geometry/Records/interface/IdealGeometryRecord.h"
|
47 |
|
|
#include "Geometry/CaloGeometry/interface/CaloSubdetectorGeometry.h"
|
48 |
|
|
#include "FWCore/Framework/interface/ESHandle.h"
|
49 |
|
|
#include "MagneticField/Engine/interface/MagneticField.h"
|
50 |
|
|
#include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
|
51 |
|
|
#include "PhysicsTools/IsolationAlgos/interface/PropagateToCal.h"
|
52 |
|
|
|
53 |
|
|
#include "TH1F.h"
|
54 |
|
|
|
55 |
|
|
//#include "CommonTools/Utils/interface/PtComparator.h"
|
56 |
|
|
template<typename T>
|
57 |
|
|
struct PtrGreaterByPt {
|
58 |
|
|
typedef T first_argument_type;
|
59 |
|
|
typedef T second_argument_type;
|
60 |
|
|
bool operator()( const T & t1, const T & t2 ) const {
|
61 |
|
|
return t1->pt() > t2->pt();
|
62 |
|
|
}
|
63 |
|
|
};
|
64 |
|
|
|
65 |
|
|
//
|
66 |
|
|
// class decleration
|
67 |
|
|
//
|
68 |
|
|
|
69 |
|
|
class FinalPlots : public edm::EDAnalyzer {
|
70 |
|
|
public:
|
71 |
|
|
explicit FinalPlots(const edm::ParameterSet&);
|
72 |
|
|
~FinalPlots();
|
73 |
|
|
|
74 |
|
|
|
75 |
|
|
private:
|
76 |
|
|
virtual void beginJob(const edm::EventSetup&);
|
77 |
|
|
virtual void analyze(const edm::Event&, const edm::EventSetup&);
|
78 |
|
|
virtual void endJob() ;
|
79 |
|
|
double dPhiBiased(const std::vector<const reco::Candidate*>& jets, const reco::Candidate& excljet);
|
80 |
|
|
double summedTowerEnergy(CaloTowerDetId& towerId, CaloTowerCollection towers, int ntw);
|
81 |
|
|
|
82 |
|
|
typedef ROOT::Math::LorentzVector<ROOT::Math::PxPyPzE4D<double> > LorentzV;
|
83 |
|
|
std::vector<double> deltaSumPt_permutations(const std::vector<LorentzV>& p4s);
|
84 |
|
|
double alphaT(const std::vector<LorentzV>& p4s);
|
85 |
|
|
|
86 |
|
|
// ----------member data ---------------------------
|
87 |
|
|
edm::InputTag weightName_;
|
88 |
|
|
edm::InputTag Jet_;
|
89 |
|
|
edm::InputTag Met_;
|
90 |
|
|
std::string name_;
|
91 |
|
|
double weight_, jetptmin_, jetetamax_;
|
92 |
|
|
|
93 |
|
|
//Jet & MET plots
|
94 |
|
|
TH1F * HT_;
|
95 |
|
|
TH1F * MHT_;
|
96 |
|
|
TH1F * MET_;
|
97 |
|
|
TH1F * AllJetsPt_;
|
98 |
|
|
TH1F * Jet1Pt_, *Jet2Pt_, *Jet3Pt_, *Jet4Pt_;
|
99 |
|
|
TH1F * Jet1E_, *Jet2E_, *Jet3E_, *Jet4E_;
|
100 |
|
|
TH1F * Jet1Phi_,*Jet2Phi_,*Jet3Phi_,*Jet4Phi_;
|
101 |
|
|
TH1F * Jet1Eta_,*Jet2Eta_,*Jet3Eta_,*Jet4Eta_;
|
102 |
|
|
TH1F * JetMult_;
|
103 |
|
|
TH1F * dPhiMin_;
|
104 |
|
|
TH1F * dPhiBiased_;
|
105 |
|
|
|
106 |
|
|
//Track & Tower plots
|
107 |
|
|
edm::InputTag tracks_;
|
108 |
|
|
edm::InputTag towers_;
|
109 |
|
|
PropagateToCal* toEcal_;
|
110 |
|
|
PropagateToCal* toHcal_;
|
111 |
|
|
TH1F * TrackTower_dEta_, *TrackTower_dPhi_, *TrackTower_dPt_;
|
112 |
|
|
int nTowerGroup_;
|
113 |
|
|
|
114 |
|
|
};
|
115 |
|
|
|
116 |
|
|
FinalPlots::FinalPlots(const edm::ParameterSet& iConfig):
|
117 |
|
|
weightName_(iConfig.getParameter<edm::InputTag>("weightName") ),
|
118 |
|
|
Jet_( iConfig.getParameter<edm::InputTag>("Jet") ),
|
119 |
|
|
Met_( iConfig.getParameter<edm::InputTag>("MET") ),
|
120 |
|
|
name_( iConfig.getParameter<std::string>("uncertainty_name") ),
|
121 |
|
|
jetptmin_ ( iConfig.getParameter<double>("JetPtMin") ),
|
122 |
|
|
jetetamax_( iConfig.getParameter<double>("JetEtaMax") ),
|
123 |
|
|
//tracks & Towers:
|
124 |
|
|
tracks_( iConfig.getParameter<edm::InputTag>( "Tracks" ) ),
|
125 |
|
|
towers_( iConfig.getParameter<edm::InputTag>( "Towers" ) ),
|
126 |
|
|
nTowerGroup_( iConfig.getParameter<int>( "GroupNTowers" ) )
|
127 |
|
|
{
|
128 |
|
|
}
|
129 |
|
|
|
130 |
|
|
// ------------ method called once each job just before starting event loop ------------
|
131 |
|
|
void
|
132 |
|
|
FinalPlots::beginJob(const edm::EventSetup&)
|
133 |
|
|
{
|
134 |
|
|
edm::Service<TFileService> fs;
|
135 |
|
|
if( !fs ){
|
136 |
|
|
throw edm::Exception( edm::errors::Configuration,
|
137 |
|
|
"TFile Service is not registered in cfg file" );
|
138 |
|
|
}
|
139 |
|
|
|
140 |
|
|
std::string histname = "HT"+name_;
|
141 |
|
|
HT_ = fs->make<TH1F>(histname.c_str(),";HT [GeV];events", 40, 0.0, 100.0);
|
142 |
|
|
HT_->Sumw2();
|
143 |
|
|
|
144 |
|
|
histname = "MHT"+name_;
|
145 |
|
|
MHT_ = fs->make<TH1F>(histname.c_str(),";MHT [GeV];events", 40, 0.0, 50.0);
|
146 |
|
|
MHT_->Sumw2();
|
147 |
|
|
histname = "MET"+name_;
|
148 |
|
|
MET_ = fs->make<TH1F>(histname.c_str(),";MET [GeV];events", 40, 0.0, 50.0);
|
149 |
|
|
MET_->Sumw2();
|
150 |
|
|
|
151 |
|
|
histname = "AllJetsPt"+name_;
|
152 |
|
|
AllJetsPt_ = fs->make<TH1F>(histname.c_str(),";All Jets Pt [GeV];events", 40, 0.0, 50.0);
|
153 |
|
|
AllJetsPt_->Sumw2();
|
154 |
|
|
|
155 |
|
|
histname = "Jet1Pt"+name_;
|
156 |
|
|
Jet1Pt_ = fs->make<TH1F>(histname.c_str(),";1. Jet Pt [GeV];events", 40, 0.0, 20.0);
|
157 |
|
|
Jet1Pt_->Sumw2();
|
158 |
|
|
histname = "Jet2Pt"+name_;
|
159 |
|
|
Jet2Pt_ = fs->make<TH1F>(histname.c_str(),";2. Jet Pt [GeV];events", 40, 0.0, 20.0);
|
160 |
|
|
Jet2Pt_->Sumw2();
|
161 |
|
|
histname = "Jet3Pt"+name_;
|
162 |
|
|
Jet3Pt_ = fs->make<TH1F>(histname.c_str(),";3. Jet Pt [GeV];events", 40, 0.0, 20.0);
|
163 |
|
|
Jet3Pt_->Sumw2();
|
164 |
|
|
histname = "Jet4Pt"+name_;
|
165 |
|
|
Jet4Pt_ = fs->make<TH1F>(histname.c_str(),";4. Jet Pt [GeV];events", 40, 0.0, 20.0);
|
166 |
|
|
Jet4Pt_->Sumw2();
|
167 |
|
|
|
168 |
|
|
histname = "Jet1E"+name_;
|
169 |
|
|
Jet1E_ = fs->make<TH1F>(histname.c_str(),";1. Jet E [GeV];events", 40, 0.0, 20.0);
|
170 |
|
|
Jet1E_->Sumw2();
|
171 |
|
|
histname = "Jet2E"+name_;
|
172 |
|
|
Jet2E_ = fs->make<TH1F>(histname.c_str(),";2. Jet E [GeV];events", 40, 0.0, 20.0);
|
173 |
|
|
Jet2E_->Sumw2();
|
174 |
|
|
histname = "Jet3E"+name_;
|
175 |
|
|
Jet3E_ = fs->make<TH1F>(histname.c_str(),";3. Jet E [GeV];events", 40, 0.0, 20.0);
|
176 |
|
|
Jet3E_->Sumw2();
|
177 |
|
|
histname = "Jet4E"+name_;
|
178 |
|
|
Jet4E_ = fs->make<TH1F>(histname.c_str(),";4. Jet E [GeV];events", 40, 0.0, 20.0);
|
179 |
|
|
Jet4E_->Sumw2();
|
180 |
|
|
|
181 |
|
|
histname = "Jet1Phi"+name_;
|
182 |
|
|
Jet1Phi_ = fs->make<TH1F>(histname.c_str(),";1. Jet #Phi;events", 40, -3.2, 3.2);
|
183 |
|
|
Jet1Phi_->Sumw2();
|
184 |
|
|
histname = "Jet2Phi"+name_;
|
185 |
|
|
Jet2Phi_ = fs->make<TH1F>(histname.c_str(),";2. Jet #Phi;events", 40, -3.2, 3.2);
|
186 |
|
|
Jet2Phi_->Sumw2();
|
187 |
|
|
histname = "Jet3Phi"+name_;
|
188 |
|
|
Jet3Phi_ = fs->make<TH1F>(histname.c_str(),";3. Jet #Phi;events", 40, -3.2, 3.2);
|
189 |
|
|
Jet3Phi_->Sumw2();
|
190 |
|
|
histname = "Jet4Phi"+name_;
|
191 |
|
|
Jet4Phi_ = fs->make<TH1F>(histname.c_str(),";4. Jet #Phi;events", 40, -3.2, 3.2);
|
192 |
|
|
Jet4Phi_->Sumw2();
|
193 |
|
|
|
194 |
|
|
histname = "Jet1Eta"+name_;
|
195 |
|
|
Jet1Eta_ = fs->make<TH1F>(histname.c_str(),";1. Jet #Eta;events", 40, -4.0, 4.0);
|
196 |
|
|
Jet1Eta_->Sumw2();
|
197 |
|
|
histname = "Jet2Eta"+name_;
|
198 |
|
|
Jet2Eta_ = fs->make<TH1F>(histname.c_str(),";2. Jet #Eta;events", 40, -4.0, 4.0);
|
199 |
|
|
Jet2Eta_->Sumw2();
|
200 |
|
|
histname = "Jet3Eta"+name_;
|
201 |
|
|
Jet3Eta_ = fs->make<TH1F>(histname.c_str(),";3. Jet #Eta;events", 40, -4.0, 4.0);
|
202 |
|
|
Jet3Eta_->Sumw2();
|
203 |
|
|
histname = "Jet4Eta"+name_;
|
204 |
|
|
Jet4Eta_ = fs->make<TH1F>(histname.c_str(),";4. Jet #Eta;events", 40, -4.0, 4.0);
|
205 |
|
|
Jet4Eta_->Sumw2();
|
206 |
|
|
|
207 |
|
|
histname = "JetMult"+name_;
|
208 |
|
|
JetMult_ = fs->make<TH1F>(histname.c_str(),";Jet multiplicity;events", 21, -0.5, 20.5);
|
209 |
|
|
JetMult_->Sumw2();
|
210 |
|
|
histname = "dPhiMin"+name_;
|
211 |
|
|
dPhiMin_ = fs->make<TH1F>(histname.c_str(),";#Delta#phi_{min};events", 20, 0.0, 3.2);
|
212 |
|
|
dPhiMin_->Sumw2();
|
213 |
|
|
histname = "dPhiBiased"+name_;
|
214 |
|
|
dPhiBiased_ = fs->make<TH1F>(histname.c_str(),";#Delta#phi_{biased};events", 20, 0.0, 3.2);
|
215 |
|
|
dPhiBiased_->Sumw2();
|
216 |
|
|
|
217 |
|
|
|
218 |
|
|
// Track & Tower
|
219 |
|
|
histname = "TrackTower_dEta"+name_;
|
220 |
|
|
TrackTower_dEta_ = fs->make<TH1F>(histname.c_str(),";#Delta#Eta(track, tower);events", 40, -5.0, 5.0);
|
221 |
|
|
TrackTower_dEta_->Sumw2();
|
222 |
|
|
|
223 |
|
|
histname = "TrackTower_dPhi"+name_;
|
224 |
|
|
TrackTower_dPhi_ = fs->make<TH1F>(histname.c_str(),";#Delta#Phi(track, tower);events", 40, 0.0, 4.0);
|
225 |
|
|
TrackTower_dPhi_->Sumw2();
|
226 |
|
|
|
227 |
|
|
histname = "TrackTower_dPt"+name_;
|
228 |
|
|
TrackTower_dPt_ = fs->make<TH1F>(histname.c_str(),";#Delta#Pt(track, tower);events", 40, -20.0, 20.0);
|
229 |
|
|
TrackTower_dPt_->Sumw2();
|
230 |
|
|
}
|
231 |
|
|
|
232 |
|
|
|
233 |
|
|
FinalPlots::~FinalPlots()
|
234 |
|
|
{
|
235 |
|
|
|
236 |
|
|
// do anything here that needs to be done at desctruction time
|
237 |
|
|
// (e.g. close files, deallocate resources etc.)
|
238 |
|
|
|
239 |
|
|
}
|
240 |
|
|
|
241 |
|
|
|
242 |
|
|
//
|
243 |
|
|
// member functions
|
244 |
|
|
//
|
245 |
|
|
// ------------ biased delta phi ----------------------------
|
246 |
|
|
double FinalPlots::dPhiBiased(const std::vector<const reco::Candidate*>& jets, const reco::Candidate& excljet)
|
247 |
|
|
{
|
248 |
|
|
math::PtEtaPhiMLorentzVector htvec(0.0, 0.0, 0.0, 0.0);
|
249 |
|
|
for (std::vector<const reco::Candidate*>::const_iterator jet=jets.begin();
|
250 |
|
|
jet!=jets.end(); ++jet){
|
251 |
|
|
if ( fabs( (*jet)->eta())>jetetamax_ ) continue;
|
252 |
|
|
if ( (*jet)->pt()<jetptmin_ ) continue;
|
253 |
|
|
if ( (*jet)->phi()==excljet.phi()&&(*jet)->eta()==excljet.eta() ) continue;
|
254 |
|
|
htvec+=(*jet)->p4();
|
255 |
|
|
}
|
256 |
|
|
return fabs(deltaPhi(htvec, excljet));
|
257 |
|
|
}
|
258 |
|
|
|
259 |
|
|
std::vector<double> FinalPlots::deltaSumPt_permutations(const std::vector<LorentzV>& p4s)
|
260 |
|
|
{
|
261 |
|
|
std::vector<std::vector<double> > ht( 1<<(p4s.size()-1) , std::vector<double>( 2, 0.) );
|
262 |
|
|
for(unsigned i=0; i < ht.size(); i++) {
|
263 |
|
|
for(unsigned j=0; j < p4s.size(); j++) {
|
264 |
|
|
ht [i] [(i/(1<<j))%2] += p4s[j].pt();
|
265 |
|
|
}
|
266 |
|
|
}
|
267 |
|
|
std::vector<double> deltaHT; for(unsigned i=0; i<ht.size(); i++) deltaHT.push_back(fabs(ht[i][0]-ht[i][1]));
|
268 |
|
|
return deltaHT;
|
269 |
|
|
}
|
270 |
|
|
|
271 |
|
|
double FinalPlots::summedTowerEnergy(CaloTowerDetId& towerId, CaloTowerCollection towers, int ntw)
|
272 |
|
|
{
|
273 |
|
|
double summedTowerEnergy=0;
|
274 |
|
|
for(int iphi = towerId.iphi()-ntw; iphi<=towerId.iphi()+ntw; ++iphi){
|
275 |
|
|
for(int ieta = towerId.ieta()-ntw; ieta<=towerId.ieta()+ntw; ++ieta){
|
276 |
|
|
CaloTowerDetId id(ieta,iphi);
|
277 |
|
|
CaloTowerCollection::const_iterator tower = towers.find(id);
|
278 |
|
|
if(tower == towers.end() ) continue;
|
279 |
|
|
summedTowerEnergy+=tower->et();
|
280 |
|
|
}
|
281 |
|
|
}
|
282 |
|
|
return summedTowerEnergy;
|
283 |
|
|
}
|
284 |
|
|
|
285 |
|
|
|
286 |
|
|
// ------------ method called to for each event ------------
|
287 |
|
|
void
|
288 |
|
|
FinalPlots::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup)
|
289 |
|
|
{
|
290 |
|
|
using namespace edm;
|
291 |
|
|
using namespace std;
|
292 |
|
|
|
293 |
|
|
//Make Jet & MET plots: -------------------------------------------------
|
294 |
|
|
edm::Handle< double > event_weight;
|
295 |
|
|
iEvent.getByLabel(weightName_, event_weight);
|
296 |
|
|
weight_ = (event_weight.isValid() ? (*event_weight) : 1.0);
|
297 |
|
|
|
298 |
|
|
edm::Handle<edm::View<reco::Candidate> > jet_hnd;
|
299 |
|
|
iEvent.getByLabel(Jet_, jet_hnd);
|
300 |
|
|
edm::View<reco::Candidate> Jets = *jet_hnd;
|
301 |
|
|
|
302 |
|
|
|
303 |
|
|
edm::Handle<edm::View<reco::MET> > met_hnd;
|
304 |
|
|
iEvent.getByLabel(Met_, met_hnd);
|
305 |
|
|
|
306 |
|
|
MET_->Fill( met_hnd->front().et(), weight_ );
|
307 |
|
|
|
308 |
|
|
//apply JEC for PAT jets:
|
309 |
|
|
std::vector<const reco::Candidate*> cJets;
|
310 |
|
|
for (edm::View<reco::Candidate>::const_iterator jet=Jets.begin();
|
311 |
|
|
jet!=Jets.end(); ++jet){
|
312 |
|
|
const pat::Jet * patJet = 0;
|
313 |
|
|
patJet = dynamic_cast <const pat::Jet*> ( &(*jet) );
|
314 |
|
|
const reco::Candidate * cJet;
|
315 |
|
|
if (patJet)
|
316 |
|
|
cJet = new pat::Jet( patJet->correctedJet("part", "b") );
|
317 |
|
|
else
|
318 |
|
|
cJet = &(*jet);
|
319 |
|
|
cJets.push_back( cJet );
|
320 |
|
|
}
|
321 |
|
|
// sort jets in Et
|
322 |
|
|
PtrGreaterByPt<const reco::Candidate*> pTComparator_;
|
323 |
|
|
std::sort(cJets.begin(), cJets.end(), pTComparator_);
|
324 |
|
|
|
325 |
|
|
double ht=0.0, dphimin=999., dphibiased=999.;
|
326 |
|
|
int jetmult = 0;
|
327 |
|
|
math::PtEtaPhiMLorentzVector htvec(0.0, 0.0, 0.0, 0.0);
|
328 |
|
|
std::vector <LorentzV> forAlphaT;
|
329 |
|
|
|
330 |
|
|
for (std::vector<const reco::Candidate*>::const_iterator cJet=cJets.begin();
|
331 |
|
|
cJet!=cJets.end(); ++cJet) {
|
332 |
|
|
|
333 |
|
|
if ( fabs((*cJet)->eta())>jetetamax_ ) continue;
|
334 |
|
|
if ( (*cJet)->pt()<jetptmin_ ) continue;
|
335 |
|
|
forAlphaT.push_back((*cJet)->p4());
|
336 |
|
|
htvec+=(*cJet)->p4();
|
337 |
|
|
ht+=(*cJet)->pt();
|
338 |
|
|
AllJetsPt_->Fill( (*cJet)->pt() );
|
339 |
|
|
if (jetmult==0) {
|
340 |
|
|
Jet1E_->Fill( (*cJet)->energy(), weight_ );
|
341 |
|
|
Jet1Pt_->Fill( (*cJet)->pt(), weight_ );
|
342 |
|
|
Jet1Phi_->Fill( (*cJet)->phi(), weight_ );
|
343 |
|
|
Jet1Eta_->Fill( (*cJet)->eta(), weight_ );
|
344 |
|
|
}
|
345 |
|
|
else if (jetmult==1) {
|
346 |
|
|
Jet2E_->Fill( (*cJet)->energy(), weight_ );
|
347 |
|
|
Jet2Pt_->Fill( (*cJet)->pt(), weight_ );
|
348 |
|
|
Jet2Phi_->Fill( (*cJet)->phi(), weight_ );
|
349 |
|
|
Jet2Eta_->Fill( (*cJet)->eta(), weight_ );
|
350 |
|
|
}
|
351 |
|
|
else if (jetmult==2) {
|
352 |
|
|
Jet3E_->Fill( (*cJet)->energy(), weight_ );
|
353 |
|
|
Jet3Pt_->Fill( (*cJet)->pt(), weight_ );
|
354 |
|
|
Jet3Phi_->Fill( (*cJet)->phi(), weight_ );
|
355 |
|
|
Jet3Eta_->Fill( (*cJet)->eta(), weight_ );
|
356 |
|
|
}
|
357 |
|
|
else if (jetmult==3) {
|
358 |
|
|
Jet4E_->Fill( (*cJet)->energy(), weight_ );
|
359 |
|
|
Jet4Pt_->Fill( (*cJet)->pt(), weight_ );
|
360 |
|
|
Jet4Phi_->Fill( (*cJet)->phi(), weight_ );
|
361 |
|
|
Jet4Eta_->Fill( (*cJet)->eta(), weight_ );
|
362 |
|
|
}
|
363 |
|
|
++jetmult;
|
364 |
|
|
|
365 |
|
|
if ( fabs(deltaPhi(**cJet, met_hnd->front())) < dphimin )
|
366 |
|
|
dphimin = fabs(deltaPhi( **cJet, met_hnd->front()));
|
367 |
|
|
|
368 |
|
|
double thisdphi = dPhiBiased(cJets, **cJet);
|
369 |
|
|
if (thisdphi < dphibiased)
|
370 |
|
|
dphibiased = thisdphi;
|
371 |
|
|
|
372 |
|
|
}
|
373 |
|
|
HT_ ->Fill( ht, weight_ );
|
374 |
|
|
MHT_->Fill( htvec.pt(), weight_ );
|
375 |
|
|
JetMult_->Fill( jetmult, weight_ );
|
376 |
|
|
dPhiMin_->Fill( dphimin, weight_ );
|
377 |
|
|
dPhiBiased_->Fill( dphibiased, weight_ );
|
378 |
|
|
|
379 |
|
|
//other variables
|
380 |
|
|
//MPT...
|
381 |
|
|
|
382 |
|
|
|
383 |
|
|
//Analyze Tracks <-> Towers ----------------------------------------------------------------------
|
384 |
|
|
|
385 |
|
|
edm::Handle<CaloTowerCollection> towers;
|
386 |
|
|
iEvent.getByLabel(towers_,towers);
|
387 |
|
|
|
388 |
|
|
edm::Handle<reco::TrackCollection> tracks;
|
389 |
|
|
iEvent.getByLabel(tracks_,tracks);
|
390 |
|
|
|
391 |
|
|
//Propagate a track to the calorimeter surface:
|
392 |
|
|
/*
|
393 |
|
|
edm::ESHandle<MagneticField> field;
|
394 |
|
|
iSetup.get<IdealMagneticFieldRecord>().get(field);
|
395 |
|
|
|
396 |
|
|
edm::ESHandle<CaloGeometry> geometry;
|
397 |
|
|
iSetup.get<IdealGeometryRecord>().get(geometry);
|
398 |
|
|
const CaloSubdetectorGeometry* towerGeometry = geometry->getSubdetectorGeometry(DetId::Calo, CaloTowerDetId::SubdetId);
|
399 |
|
|
|
400 |
|
|
|
401 |
|
|
double trkIso=-1.;
|
402 |
|
|
reco::Track selTrack;
|
403 |
|
|
//----------------------------------------------------------------
|
404 |
|
|
//search for most isolated track of good quality
|
405 |
|
|
//----------------------------------------------------------------
|
406 |
|
|
//Isolation::const_iterator iso = isolation->begin();
|
407 |
|
|
|
408 |
|
|
//these parameters should be defined in the cfg file:
|
409 |
|
|
double maxEta_ = 2.4;
|
410 |
|
|
double minPt_ = 5.0;
|
411 |
|
|
double maxChiSquare_ = 40.0; //just a guess!
|
412 |
|
|
int minNumOfHits_ = 10;
|
413 |
|
|
|
414 |
|
|
reco::TrackCollection::const_iterator track = tracks->begin();
|
415 |
|
|
for( ; track!=tracks->end(); //, iso!=isolation->end();
|
416 |
|
|
++track//,++iso
|
417 |
|
|
){
|
418 |
|
|
if( !(fabs(track->eta())<maxEta_) )
|
419 |
|
|
continue;
|
420 |
|
|
|
421 |
|
|
if( !(minPt_<track->pt()) )
|
422 |
|
|
continue;
|
423 |
|
|
|
424 |
|
|
if( !(track->normalizedChi2()<maxChiSquare_) )
|
425 |
|
|
continue;
|
426 |
|
|
|
427 |
|
|
if( !(track->found()>minNumOfHits_) )
|
428 |
|
|
continue;
|
429 |
|
|
|
430 |
|
|
//double isolation = *iso-track->pt();
|
431 |
|
|
//if(trkIso<0 || isolation<trkIso){
|
432 |
|
|
// selTrack =*track;
|
433 |
|
|
// trkIso=isolation;
|
434 |
|
|
//}
|
435 |
|
|
}
|
436 |
|
|
|
437 |
|
|
//these parameters should be defined in the cfg file:
|
438 |
|
|
double maxDeltaEta_ = 0.3;
|
439 |
|
|
double maxDeltaPhi_ = 0.3;
|
440 |
|
|
|
441 |
|
|
double en1x1=0, enSum=0, dPhi=0, dEta=0;
|
442 |
|
|
if( 0<=trkIso ){ // && trkIso<maxIsolation_ ){
|
443 |
|
|
//----------------------------------------------------------------
|
444 |
|
|
//propagate track to hcal and
|
445 |
|
|
//search for closest calo tower
|
446 |
|
|
//----------------------------------------------------------------
|
447 |
|
|
const GlobalPoint vtx(selTrack.vx(), selTrack.vy(), selTrack.vz());
|
448 |
|
|
GlobalVector atHcal(selTrack.px(), selTrack.py(), selTrack.pz());//track at hcal
|
449 |
|
|
toHcal_->propagate( vtx, atHcal, selTrack.charge(), &*field);
|
450 |
|
|
const GlobalPoint poI(atHcal.x(), atHcal.y(), atHcal.z()); //point of incidence
|
451 |
|
|
CaloTowerDetId towerId = towerGeometry->getClosestCell(poI);
|
452 |
|
|
|
453 |
|
|
CaloTowerDetId id(towerId.ieta(),towerId.iphi());
|
454 |
|
|
CaloTowerCollection::const_iterator selTower = towers->find(id);
|
455 |
|
|
if(selTower!=towers->end()){
|
456 |
|
|
dEta = selTower->eta()-selTrack.eta();
|
457 |
|
|
dPhi = deltaPhi(selTower->phi(),selTrack.phi());
|
458 |
|
|
|
459 |
|
|
if( fabs(dEta)<maxDeltaEta_ && fabs(dPhi<maxDeltaPhi_) ){
|
460 |
|
|
for(int iphi = towerId.iphi()-nTowerGroup_; iphi<=towerId.iphi()+nTowerGroup_; ++iphi){
|
461 |
|
|
for(int ieta = towerId.ieta()-nTowerGroup_; ieta<=towerId.ieta()+nTowerGroup_; ++ieta){
|
462 |
|
|
CaloTowerDetId id(ieta,iphi);
|
463 |
|
|
CaloTowerCollection::const_iterator twr = towers->find(id);
|
464 |
|
|
if(twr == towers->end() ) continue;
|
465 |
|
|
}
|
466 |
|
|
}
|
467 |
|
|
en1x1 = summedTowerEnergy(towerId, *towers, 0);
|
468 |
|
|
enSum = summedTowerEnergy(towerId, *towers, nTowerGroup_);
|
469 |
|
|
}
|
470 |
|
|
}
|
471 |
|
|
}
|
472 |
|
|
*/
|
473 |
|
|
|
474 |
|
|
//these parameters should be defined in the cfg file:
|
475 |
|
|
double maxEta_ = 2.4;
|
476 |
|
|
double minPt_ = 5.0;
|
477 |
|
|
double maxChiSquare_ = 40.0; //just a guess!
|
478 |
|
|
int minNumOfHits_ = 10;
|
479 |
|
|
reco::TrackCollection::const_iterator track = tracks->begin();
|
480 |
|
|
for( ; track!=tracks->end(); ++track ){ //loop over all tracks
|
481 |
|
|
if( !(fabs(track->eta())<maxEta_) ) continue;
|
482 |
|
|
if( !(minPt_<track->pt()) ) continue;
|
483 |
|
|
if( !(track->normalizedChi2()<maxChiSquare_) ) continue;
|
484 |
|
|
if( !(track->found()>minNumOfHits_) ) continue;
|
485 |
|
|
//check also for isolation
|
486 |
|
|
//e.g. check also if track is a Pion (if running on MC)
|
487 |
|
|
|
488 |
|
|
double dRmin2 = 99999.;
|
489 |
|
|
double dEta=0, dPhi=0, dPt=0, dEtaMin, dPhiMin;
|
490 |
|
|
|
491 |
|
|
//loop over all towers
|
492 |
|
|
for (CaloTowerCollection::const_iterator tower = towers->begin(); tower!=towers->end(); ++tower) {
|
493 |
|
|
|
494 |
|
|
//dEta = track->eta() - tower->p4( track.vertex() ).eta();
|
495 |
|
|
dEta = track->outerEta() - tower->eta();
|
496 |
|
|
dPhi = deltaPhi(track->outerPhi(),tower->phi());
|
497 |
|
|
if (dEta*dEta+dPhi*dPhi<dRmin2) { //find closest tower to track
|
498 |
|
|
dRmin2 = dEta*dEta+dPhi*dPhi;
|
499 |
|
|
dPt = track->pt() - tower->pt();
|
500 |
|
|
dEtaMin = dEta;
|
501 |
|
|
dPhiMin = dPhi;
|
502 |
|
|
}
|
503 |
|
|
}
|
504 |
|
|
|
505 |
|
|
//fill histograms
|
506 |
|
|
TrackTower_dEta_->Fill( dEtaMin );
|
507 |
|
|
TrackTower_dPhi_->Fill( dPhiMin );
|
508 |
|
|
TrackTower_dPt_ ->Fill( dPt );
|
509 |
|
|
}
|
510 |
|
|
}
|
511 |
|
|
|
512 |
|
|
|
513 |
|
|
// ------------ method called once each job just after ending the event loop ------------
|
514 |
|
|
void
|
515 |
|
|
FinalPlots::endJob() {
|
516 |
|
|
}
|
517 |
|
|
|
518 |
|
|
//define this as a plug-in
|
519 |
|
|
DEFINE_FWK_MODULE(FinalPlots);
|