1 |
yilmaz |
1.1 |
#ifndef RecoJets_JetProducers_plugins_JetAlgorithmAnalyzer_h
|
2 |
|
|
#define RecoJets_JetProducers_plugins_JetAlgorithmAnalyzer_h
|
3 |
|
|
|
4 |
|
|
#include "RecoJets/JetProducers/plugins/VirtualJetProducer.h"
|
5 |
|
|
|
6 |
|
|
#include "FWCore/ServiceRegistry/interface/Service.h"
|
7 |
|
|
#include "CommonTools/UtilAlgos/interface/TFileService.h"
|
8 |
|
|
#include "TNtuple.h"
|
9 |
|
|
|
10 |
|
|
class JetAlgorithmAnalyzer : public VirtualJetProducer
|
11 |
|
|
{
|
12 |
|
|
|
13 |
|
|
public:
|
14 |
|
|
//
|
15 |
|
|
// construction/destruction
|
16 |
|
|
//
|
17 |
|
|
explicit JetAlgorithmAnalyzer(const edm::ParameterSet& iConfig);
|
18 |
|
|
virtual ~JetAlgorithmAnalyzer();
|
19 |
|
|
|
20 |
|
|
protected:
|
21 |
|
|
|
22 |
|
|
//
|
23 |
|
|
// member functions
|
24 |
|
|
//
|
25 |
|
|
|
26 |
|
|
virtual void produce( edm::Event& iEvent, const edm::EventSetup& iSetup );
|
27 |
|
|
virtual void runAlgorithm( edm::Event& iEvent, const edm::EventSetup& iSetup );
|
28 |
|
|
virtual void output( edm::Event & iEvent, edm::EventSetup const& iSetup );
|
29 |
|
|
template< typename T >
|
30 |
|
|
void writeBkgJets( edm::Event & iEvent, edm::EventSetup const& iSetup );
|
31 |
|
|
|
32 |
|
|
void fillNtuple(bool output, const std::vector<fastjet::PseudoJet>& jets, int step);
|
33 |
|
|
void fillTowerNtuple(const std::vector<fastjet::PseudoJet>& jets, int step);
|
34 |
|
|
void fillJetNtuple(const std::vector<fastjet::PseudoJet>& jets, int step);
|
35 |
|
|
void fillBkgNtuple(const PileUpSubtractor* subtractor, int step);
|
36 |
|
|
|
37 |
|
|
private:
|
38 |
|
|
|
39 |
|
|
// trackjet clustering parameters
|
40 |
|
|
bool useOnlyVertexTracks_;
|
41 |
|
|
bool useOnlyOnePV_;
|
42 |
|
|
float dzTrVtxMax_;
|
43 |
|
|
|
44 |
|
|
int nFill_;
|
45 |
|
|
float etaMax_;
|
46 |
|
|
int iev_;
|
47 |
|
|
bool avoidNegative_;
|
48 |
|
|
|
49 |
|
|
bool doAnalysis_;
|
50 |
|
|
|
51 |
|
|
TNtuple* ntTowers;
|
52 |
|
|
TNtuple* ntJets;
|
53 |
|
|
TNtuple* ntPU;
|
54 |
|
|
TNtuple* ntRandom;
|
55 |
|
|
|
56 |
|
|
const CaloGeometry *geo;
|
57 |
|
|
edm::Service<TFileService> f;
|
58 |
|
|
};
|
59 |
|
|
|
60 |
|
|
|
61 |
|
|
#endif
|
62 |
|
|
////////////////////////////////////////////////////////////////////////////////
|
63 |
|
|
//
|
64 |
|
|
// JetAlgorithmAnalyzer
|
65 |
|
|
// ------------------
|
66 |
|
|
//
|
67 |
|
|
// 04/21/2009 Philipp Schieferdecker <philipp.schieferdecker@cern.ch>
|
68 |
|
|
////////////////////////////////////////////////////////////////////////////////
|
69 |
|
|
|
70 |
|
|
//#include "RecoJets/JetProducers/plugins/JetAlgorithmAnalyzer.h"
|
71 |
|
|
//#include "JetAlgorithmAnalyzer.h"
|
72 |
|
|
|
73 |
|
|
#include "RecoJets/JetProducers/interface/JetSpecific.h"
|
74 |
|
|
|
75 |
|
|
#include "FWCore/Framework/interface/Event.h"
|
76 |
|
|
#include "FWCore/Framework/interface/EventSetup.h"
|
77 |
|
|
#include "FWCore/Framework/interface/ESHandle.h"
|
78 |
|
|
#include "FWCore/Utilities/interface/CodedException.h"
|
79 |
|
|
#include "FWCore/MessageLogger/interface/MessageLogger.h"
|
80 |
|
|
#include "FWCore/Framework/interface/MakerMacros.h"
|
81 |
|
|
#include "FWCore/Utilities/interface/RandomNumberGenerator.h"
|
82 |
|
|
#include "FWCore/ServiceRegistry/interface/Service.h"
|
83 |
|
|
|
84 |
|
|
#include "DataFormats/Common/interface/Handle.h"
|
85 |
|
|
#include "DataFormats/VertexReco/interface/Vertex.h"
|
86 |
|
|
#include "DataFormats/VertexReco/interface/VertexFwd.h"
|
87 |
|
|
#include "DataFormats/JetReco/interface/CaloJetCollection.h"
|
88 |
|
|
#include "DataFormats/JetReco/interface/GenJetCollection.h"
|
89 |
|
|
#include "DataFormats/JetReco/interface/PFJetCollection.h"
|
90 |
|
|
#include "DataFormats/JetReco/interface/BasicJetCollection.h"
|
91 |
|
|
#include "DataFormats/Candidate/interface/CandidateFwd.h"
|
92 |
|
|
#include "DataFormats/RecoCandidate/interface/RecoChargedCandidate.h"
|
93 |
|
|
#include "DataFormats/Candidate/interface/LeafCandidate.h"
|
94 |
|
|
#include "DataFormats/Math/interface/deltaR.h"
|
95 |
|
|
|
96 |
|
|
#include "DataFormats/CaloTowers/interface/CaloTower.h"
|
97 |
|
|
#include "Geometry/CaloGeometry/interface/CaloGeometry.h"
|
98 |
|
|
#include "Geometry/Records/interface/CaloGeometryRecord.h"
|
99 |
|
|
|
100 |
|
|
#include "fastjet/SISConePlugin.hh"
|
101 |
|
|
#include "fastjet/CMSIterativeConePlugin.hh"
|
102 |
|
|
#include "fastjet/ATLASConePlugin.hh"
|
103 |
|
|
#include "fastjet/CDFMidPointPlugin.hh"
|
104 |
|
|
|
105 |
|
|
#include "CLHEP/Random/RandomEngine.h"
|
106 |
|
|
|
107 |
|
|
#include <iostream>
|
108 |
|
|
#include <memory>
|
109 |
|
|
#include <algorithm>
|
110 |
|
|
#include <limits>
|
111 |
|
|
#include <cmath>
|
112 |
|
|
|
113 |
|
|
using namespace std;
|
114 |
|
|
using namespace edm;
|
115 |
|
|
|
116 |
|
|
static const double pi = 3.14159265358979;
|
117 |
|
|
|
118 |
|
|
CLHEP::HepRandomEngine* randomEngine;
|
119 |
|
|
|
120 |
|
|
|
121 |
|
|
////////////////////////////////////////////////////////////////////////////////
|
122 |
|
|
// construction / destruction
|
123 |
|
|
////////////////////////////////////////////////////////////////////////////////
|
124 |
|
|
|
125 |
|
|
//______________________________________________________________________________
|
126 |
|
|
JetAlgorithmAnalyzer::JetAlgorithmAnalyzer(const edm::ParameterSet& iConfig)
|
127 |
|
|
: VirtualJetProducer( iConfig ),
|
128 |
|
|
nFill_(5),
|
129 |
|
|
etaMax_(3),
|
130 |
|
|
iev_(0),
|
131 |
|
|
geo(0)
|
132 |
|
|
{
|
133 |
|
|
edm::Service<RandomNumberGenerator> rng;
|
134 |
|
|
randomEngine = &(rng->getEngine());
|
135 |
|
|
|
136 |
|
|
doAnalysis_ = iConfig.getUntrackedParameter<bool>("doAnalysis",true);
|
137 |
|
|
avoidNegative_ = iConfig.getParameter<bool>("avoidNegative");
|
138 |
|
|
|
139 |
|
|
if ( iConfig.exists("UseOnlyVertexTracks") )
|
140 |
|
|
useOnlyVertexTracks_ = iConfig.getParameter<bool>("UseOnlyVertexTracks");
|
141 |
|
|
else
|
142 |
|
|
useOnlyVertexTracks_ = false;
|
143 |
|
|
|
144 |
|
|
if ( iConfig.exists("UseOnlyOnePV") )
|
145 |
|
|
useOnlyOnePV_ = iConfig.getParameter<bool>("UseOnlyOnePV");
|
146 |
|
|
else
|
147 |
|
|
useOnlyOnePV_ = false;
|
148 |
|
|
|
149 |
|
|
if ( iConfig.exists("DrTrVtxMax") )
|
150 |
|
|
dzTrVtxMax_ = iConfig.getParameter<double>("DzTrVtxMax");
|
151 |
|
|
else
|
152 |
|
|
dzTrVtxMax_ = false;
|
153 |
|
|
|
154 |
|
|
produces<std::vector<bool> >("directions");
|
155 |
|
|
|
156 |
|
|
if(doAnalysis_){
|
157 |
|
|
ntTowers = f->make<TNtuple>("ntTowers","Algorithm Analysis Towers","eta:phi:et:step:event");
|
158 |
|
|
ntJets = f->make<TNtuple>("ntJets","Algorithm Analysis Jets","eta:phi:et:step:event");
|
159 |
|
|
ntPU = f->make<TNtuple>("ntPU","Algorithm Analysis Background","eta:mean:sigma:step:event");
|
160 |
|
|
ntRandom = f->make<TNtuple>("ntRandom","Algorithm Analysis Background","eta:phi:et:pu:event");
|
161 |
|
|
}
|
162 |
|
|
|
163 |
|
|
}
|
164 |
|
|
|
165 |
|
|
|
166 |
|
|
//______________________________________________________________________________
|
167 |
|
|
JetAlgorithmAnalyzer::~JetAlgorithmAnalyzer()
|
168 |
|
|
{
|
169 |
|
|
}
|
170 |
|
|
|
171 |
|
|
|
172 |
|
|
void JetAlgorithmAnalyzer::fillNtuple(bool output, const std::vector<fastjet::PseudoJet>& jets, int step){
|
173 |
|
|
if(!doAnalysis_) return;
|
174 |
|
|
|
175 |
|
|
TNtuple* nt;
|
176 |
|
|
if(output) nt = ntJets;
|
177 |
|
|
else nt = ntTowers;
|
178 |
|
|
|
179 |
|
|
for(unsigned int i = 0; i < jets.size(); ++i){
|
180 |
|
|
const fastjet::PseudoJet& jet = jets[i];
|
181 |
|
|
nt->Fill(jet.eta(),jet.phi(),jet.perp(),step,iev_);
|
182 |
|
|
}
|
183 |
|
|
|
184 |
|
|
}
|
185 |
|
|
|
186 |
|
|
|
187 |
|
|
void JetAlgorithmAnalyzer::fillTowerNtuple(const std::vector<fastjet::PseudoJet>& jets, int step){
|
188 |
|
|
fillNtuple(0,jets,step);
|
189 |
|
|
}
|
190 |
|
|
|
191 |
|
|
void JetAlgorithmAnalyzer::fillJetNtuple(const std::vector<fastjet::PseudoJet>& jets, int step){
|
192 |
|
|
fillNtuple(1,jets,step);
|
193 |
|
|
}
|
194 |
|
|
|
195 |
|
|
void JetAlgorithmAnalyzer::fillBkgNtuple(const PileUpSubtractor* subtractor, int step){
|
196 |
|
|
if(!doAnalysis_) return;
|
197 |
|
|
|
198 |
|
|
CaloTowerCollection col;
|
199 |
|
|
|
200 |
|
|
for(int ieta = 1; ieta < 29; ++ieta){
|
201 |
|
|
|
202 |
|
|
CaloTowerDetId id(ieta,0);
|
203 |
|
|
const GlobalPoint& hitpoint = geo->getPosition(id);
|
204 |
|
|
double eta = hitpoint.eta();
|
205 |
|
|
|
206 |
|
|
// CaloTowerDetId id(int ieta = 0, int iphi = 0);
|
207 |
|
|
math::PtEtaPhiMLorentzVector p4(1,eta,0,0);
|
208 |
|
|
GlobalPoint pos(1,1,1);
|
209 |
|
|
CaloTower c(id,1.,1.,1.,1,1, p4, pos,pos);
|
210 |
|
|
col.push_back(c);
|
211 |
|
|
reco::CandidatePtr ptr(&col,col.size() - 1);
|
212 |
|
|
double mean = subtractor->getMeanAtTower(ptr);
|
213 |
|
|
double sigma = subtractor->getMeanAtTower(ptr);
|
214 |
|
|
|
215 |
|
|
ntPU->Fill(eta,mean,sigma,step,iev_);
|
216 |
|
|
}
|
217 |
|
|
}
|
218 |
|
|
|
219 |
|
|
void JetAlgorithmAnalyzer::produce(edm::Event& iEvent,const edm::EventSetup& iSetup)
|
220 |
|
|
{
|
221 |
|
|
if(!geo){
|
222 |
|
|
edm::ESHandle<CaloGeometry> pGeo;
|
223 |
|
|
iSetup.get<CaloGeometryRecord>().get(pGeo);
|
224 |
|
|
geo = pGeo.product();
|
225 |
|
|
}
|
226 |
|
|
|
227 |
|
|
LogDebug("VirtualJetProducer") << "Entered produce\n";
|
228 |
|
|
//determine signal vertex
|
229 |
|
|
vertex_=reco::Jet::Point(0,0,0);
|
230 |
|
|
if (makeCaloJet(jetTypeE)&&doPVCorrection_) {
|
231 |
|
|
LogDebug("VirtualJetProducer") << "Adding PV info\n";
|
232 |
|
|
edm::Handle<reco::VertexCollection> pvCollection;
|
233 |
|
|
iEvent.getByLabel(srcPVs_,pvCollection);
|
234 |
|
|
if (pvCollection->size()>0) vertex_=pvCollection->begin()->position();
|
235 |
|
|
}
|
236 |
|
|
|
237 |
|
|
// For Pileup subtraction using offset correction:
|
238 |
|
|
// set up geometry map
|
239 |
|
|
if ( doPUOffsetCorr_ ) {
|
240 |
|
|
subtractor_->setupGeometryMap(iEvent, iSetup);
|
241 |
|
|
}
|
242 |
|
|
|
243 |
|
|
// clear data
|
244 |
|
|
LogDebug("VirtualJetProducer") << "Clear data\n";
|
245 |
|
|
fjInputs_.clear();
|
246 |
|
|
fjJets_.clear();
|
247 |
|
|
inputs_.clear();
|
248 |
|
|
|
249 |
|
|
// get inputs and convert them to the fastjet format (fastjet::PeudoJet)
|
250 |
|
|
edm::Handle<reco::CandidateView> inputsHandle;
|
251 |
|
|
iEvent.getByLabel(src_,inputsHandle);
|
252 |
|
|
for (size_t i = 0; i < inputsHandle->size(); ++i) {
|
253 |
|
|
inputs_.push_back(inputsHandle->ptrAt(i));
|
254 |
|
|
}
|
255 |
|
|
LogDebug("VirtualJetProducer") << "Got inputs\n";
|
256 |
|
|
|
257 |
|
|
// Convert candidates to fastjet::PseudoJets.
|
258 |
|
|
// Also correct to Primary Vertex. Will modify fjInputs_
|
259 |
|
|
// and use inputs_
|
260 |
|
|
fjInputs_.reserve(inputs_.size());
|
261 |
|
|
inputTowers();
|
262 |
|
|
LogDebug("VirtualJetProducer") << "Inputted towers\n";
|
263 |
|
|
|
264 |
|
|
fillTowerNtuple(fjInputs_,0);
|
265 |
|
|
fillBkgNtuple(subtractor_.get(),0);
|
266 |
|
|
|
267 |
|
|
// For Pileup subtraction using offset correction:
|
268 |
|
|
// Subtract pedestal.
|
269 |
|
|
if ( doPUOffsetCorr_ ) {
|
270 |
|
|
subtractor_->calculatePedestal(fjInputs_);
|
271 |
|
|
|
272 |
|
|
fillTowerNtuple(fjInputs_,1);
|
273 |
|
|
fillBkgNtuple(subtractor_.get(),1);
|
274 |
|
|
|
275 |
|
|
subtractor_->subtractPedestal(fjInputs_);
|
276 |
|
|
|
277 |
|
|
fillTowerNtuple(fjInputs_,2);
|
278 |
|
|
fillBkgNtuple(subtractor_.get(),2);
|
279 |
|
|
|
280 |
|
|
LogDebug("VirtualJetProducer") << "Subtracted pedestal\n";
|
281 |
|
|
}
|
282 |
|
|
|
283 |
|
|
// Run algorithm. Will modify fjJets_ and allocate fjClusterSeq_.
|
284 |
|
|
// This will use fjInputs_
|
285 |
|
|
runAlgorithm( iEvent, iSetup );
|
286 |
|
|
|
287 |
|
|
fillTowerNtuple(fjInputs_,3);
|
288 |
|
|
fillBkgNtuple(subtractor_.get(),3);
|
289 |
|
|
fillJetNtuple(fjJets_,3);
|
290 |
|
|
|
291 |
|
|
if ( doPUOffsetCorr_ ) {
|
292 |
|
|
subtractor_->setAlgorithm(fjClusterSeq_);
|
293 |
|
|
}
|
294 |
|
|
|
295 |
|
|
LogDebug("VirtualJetProducer") << "Ran algorithm\n";
|
296 |
|
|
|
297 |
|
|
// For Pileup subtraction using offset correction:
|
298 |
|
|
// Now we find jets and need to recalculate their energy,
|
299 |
|
|
// mark towers participated in jet,
|
300 |
|
|
// remove occupied towers from the list and recalculate mean and sigma
|
301 |
|
|
// put the initial towers collection to the jet,
|
302 |
|
|
// and subtract from initial towers in jet recalculated mean and sigma of towers
|
303 |
|
|
if ( doPUOffsetCorr_ ) {
|
304 |
|
|
vector<fastjet::PseudoJet> orphanInput;
|
305 |
|
|
subtractor_->calculateOrphanInput(orphanInput);
|
306 |
|
|
fillTowerNtuple(orphanInput,4);
|
307 |
|
|
fillBkgNtuple(subtractor_.get(),4);
|
308 |
|
|
fillJetNtuple(fjJets_,4);
|
309 |
|
|
|
310 |
|
|
subtractor_->calculatePedestal(orphanInput);
|
311 |
|
|
fillTowerNtuple(orphanInput,5);
|
312 |
|
|
fillBkgNtuple(subtractor_.get(),5);
|
313 |
|
|
fillJetNtuple(fjJets_,5);
|
314 |
|
|
|
315 |
|
|
subtractor_->offsetCorrectJets();
|
316 |
|
|
fillTowerNtuple(orphanInput,6);
|
317 |
|
|
fillBkgNtuple(subtractor_.get(),6);
|
318 |
|
|
fillJetNtuple(fjJets_,6);
|
319 |
|
|
|
320 |
|
|
}
|
321 |
|
|
|
322 |
|
|
// Write the output jets.
|
323 |
|
|
// This will (by default) call the member function template
|
324 |
|
|
// "writeJets", but can be overridden.
|
325 |
|
|
// this will use inputs_
|
326 |
|
|
output( iEvent, iSetup );
|
327 |
|
|
fillBkgNtuple(subtractor_.get(),7);
|
328 |
|
|
fillJetNtuple(fjJets_,7);
|
329 |
|
|
|
330 |
|
|
LogDebug("VirtualJetProducer") << "Wrote jets\n";
|
331 |
|
|
|
332 |
|
|
++iev_;
|
333 |
|
|
|
334 |
|
|
return;
|
335 |
|
|
}
|
336 |
|
|
|
337 |
|
|
void JetAlgorithmAnalyzer::output(edm::Event & iEvent, edm::EventSetup const& iSetup)
|
338 |
|
|
{
|
339 |
|
|
// Write jets and constitutents. Will use fjJets_, inputs_
|
340 |
|
|
// and fjClusterSeq_
|
341 |
|
|
switch( jetTypeE ) {
|
342 |
|
|
case JetType::CaloJet :
|
343 |
|
|
writeBkgJets<reco::CaloJet>( iEvent, iSetup);
|
344 |
|
|
break;
|
345 |
|
|
case JetType::PFJet :
|
346 |
|
|
writeBkgJets<reco::PFJet>( iEvent, iSetup);
|
347 |
|
|
break;
|
348 |
|
|
case JetType::GenJet :
|
349 |
|
|
writeBkgJets<reco::GenJet>( iEvent, iSetup);
|
350 |
|
|
break;
|
351 |
|
|
case JetType::TrackJet :
|
352 |
|
|
writeBkgJets<reco::TrackJet>( iEvent, iSetup);
|
353 |
|
|
break;
|
354 |
|
|
case JetType::BasicJet :
|
355 |
|
|
writeBkgJets<reco::BasicJet>( iEvent, iSetup);
|
356 |
|
|
break;
|
357 |
|
|
default:
|
358 |
|
|
edm::LogError("InvalidInput") << " invalid jet type in VirtualJetProducer\n";
|
359 |
|
|
break;
|
360 |
|
|
};
|
361 |
|
|
|
362 |
|
|
}
|
363 |
|
|
|
364 |
|
|
template< typename T >
|
365 |
|
|
void JetAlgorithmAnalyzer::writeBkgJets( edm::Event & iEvent, edm::EventSetup const& iSetup )
|
366 |
|
|
{
|
367 |
|
|
// produce output jet collection
|
368 |
|
|
|
369 |
|
|
using namespace reco;
|
370 |
|
|
|
371 |
|
|
std::vector<fastjet::PseudoJet> fjFakeJets_;
|
372 |
|
|
std::vector<std::vector<reco::CandidatePtr> > constituents_;
|
373 |
|
|
vector<double> phiRandom;
|
374 |
|
|
vector<double> etaRandom;
|
375 |
|
|
vector<double> et;
|
376 |
|
|
vector<double> pileUp;
|
377 |
|
|
std::auto_ptr<std::vector<bool> > directions(new std::vector<bool>());
|
378 |
|
|
directions->reserve(nFill_);
|
379 |
|
|
|
380 |
|
|
phiRandom.reserve(nFill_);
|
381 |
|
|
etaRandom.reserve(nFill_);
|
382 |
|
|
et.reserve(nFill_);
|
383 |
|
|
pileUp.reserve(nFill_);
|
384 |
|
|
|
385 |
|
|
fjFakeJets_.reserve(nFill_);
|
386 |
|
|
constituents_.reserve(nFill_);
|
387 |
|
|
|
388 |
|
|
constituents_.reserve(nFill_);
|
389 |
|
|
for(int ijet = 0; ijet < nFill_; ++ijet){
|
390 |
|
|
vector<reco::CandidatePtr> vec;
|
391 |
|
|
constituents_.push_back(vec);
|
392 |
|
|
directions->push_back(true);
|
393 |
|
|
}
|
394 |
|
|
|
395 |
|
|
for(int ijet = 0; ijet < nFill_; ++ijet){
|
396 |
|
|
phiRandom[ijet] = 2*pi*randomEngine->flat() - pi;
|
397 |
|
|
etaRandom[ijet] = 2*etaMax_*randomEngine->flat() - etaMax_;
|
398 |
|
|
et[ijet] = 0;
|
399 |
|
|
pileUp[ijet] = 0;
|
400 |
|
|
}
|
401 |
|
|
|
402 |
|
|
for (vector<fastjet::PseudoJet>::const_iterator input_object = fjInputs_.begin (),
|
403 |
|
|
fjInputsEnd = fjInputs_.end();
|
404 |
|
|
input_object != fjInputsEnd; ++input_object) {
|
405 |
|
|
|
406 |
|
|
const reco::CandidatePtr & tower=inputs_[input_object->user_index()];
|
407 |
|
|
for(int ir = 0; ir < nFill_; ++ir){
|
408 |
|
|
if(reco::deltaR(tower->eta(),tower->phi(),etaRandom[ir],phiRandom[ir]) > rParam_) continue;
|
409 |
|
|
|
410 |
|
|
constituents_[ir].push_back(tower);
|
411 |
|
|
|
412 |
|
|
double towet = tower->et();
|
413 |
|
|
double putow = subtractor_->getPileUpAtTower(tower);
|
414 |
|
|
double etadd = towet - putow;
|
415 |
|
|
if(avoidNegative_ && etadd < 0.) etadd = 0;
|
416 |
|
|
et[ir] += etadd;
|
417 |
|
|
pileUp[ir] += towet - etadd;
|
418 |
|
|
}
|
419 |
|
|
}
|
420 |
|
|
cout<<"Start filling jets"<<endl;
|
421 |
|
|
|
422 |
|
|
for(int ir = 0; ir < nFill_; ++ir){
|
423 |
|
|
|
424 |
|
|
if(doAnalysis_) ntRandom->Fill(etaRandom[ir],phiRandom[ir],et[ir],pileUp[ir],iev_);
|
425 |
|
|
|
426 |
|
|
if(et[ir] < 0){
|
427 |
|
|
cout<<"Flipping vector"<<endl;
|
428 |
|
|
(*directions)[ir] = false;
|
429 |
|
|
et[ir] = -et[ir];
|
430 |
|
|
}else{
|
431 |
|
|
cout<<"Keep vector same"<<endl;
|
432 |
|
|
(*directions)[ir] = true;
|
433 |
|
|
}
|
434 |
|
|
cout<<"Lorentz"<<endl;
|
435 |
|
|
|
436 |
|
|
math::PtEtaPhiMLorentzVector p(et[ir],etaRandom[ir],phiRandom[ir],0);
|
437 |
|
|
fastjet::PseudoJet jet(p.px(),p.py(),p.pz(),p.energy());
|
438 |
|
|
fjFakeJets_.push_back(jet);
|
439 |
|
|
}
|
440 |
|
|
|
441 |
|
|
std::auto_ptr<std::vector<T> > jets(new std::vector<T>() );
|
442 |
|
|
jets->reserve(fjFakeJets_.size());
|
443 |
|
|
|
444 |
|
|
for (unsigned int ijet=0;ijet<fjFakeJets_.size();++ijet) {
|
445 |
|
|
// allocate this jet
|
446 |
|
|
T jet;
|
447 |
|
|
// get the fastjet jet
|
448 |
|
|
const fastjet::PseudoJet& fjJet = fjFakeJets_[ijet];
|
449 |
|
|
|
450 |
|
|
// convert them to CandidatePtr vector
|
451 |
|
|
std::vector<CandidatePtr> constituents =
|
452 |
|
|
constituents_[ijet];
|
453 |
|
|
|
454 |
|
|
writeSpecific(jet,
|
455 |
|
|
Particle::LorentzVector(fjJet.px(),
|
456 |
|
|
fjJet.py(),
|
457 |
|
|
fjJet.pz(),
|
458 |
|
|
fjJet.E()),
|
459 |
|
|
vertex_,
|
460 |
|
|
constituents, iSetup);
|
461 |
|
|
|
462 |
|
|
// calcuate the jet area
|
463 |
|
|
double jetArea=0.0;
|
464 |
|
|
jet.setJetArea (jetArea);
|
465 |
|
|
if(doPUOffsetCorr_){
|
466 |
|
|
jet.setPileup(pileUp[ijet]);
|
467 |
|
|
}else{
|
468 |
|
|
jet.setPileup (0.0);
|
469 |
|
|
}
|
470 |
|
|
|
471 |
|
|
// add to the list
|
472 |
|
|
jets->push_back(jet);
|
473 |
|
|
}
|
474 |
|
|
|
475 |
|
|
// put the jets in the collection
|
476 |
|
|
iEvent.put(jets);
|
477 |
|
|
iEvent.put(directions,"directions");
|
478 |
|
|
// calculate rho (median pT per unit area, for PU&UE subtraction down the line
|
479 |
|
|
std::auto_ptr<double> rho(new double(0.0));
|
480 |
|
|
std::auto_ptr<double> sigma(new double(0.0));
|
481 |
|
|
|
482 |
|
|
if (doRhoFastjet_) {
|
483 |
|
|
fastjet::ClusterSequenceArea const * clusterSequenceWithArea =
|
484 |
|
|
dynamic_cast<fastjet::ClusterSequenceArea const *> ( &*fjClusterSeq_ );
|
485 |
|
|
double mean_area = 0;
|
486 |
|
|
clusterSequenceWithArea->get_median_rho_and_sigma(*fjRangeDef_,false,*rho,*sigma,mean_area);
|
487 |
|
|
iEvent.put(rho,"rho");
|
488 |
|
|
iEvent.put(sigma,"sigma");
|
489 |
|
|
}
|
490 |
|
|
}
|
491 |
|
|
|
492 |
|
|
void JetAlgorithmAnalyzer::runAlgorithm( edm::Event & iEvent, edm::EventSetup const& iSetup)
|
493 |
|
|
{
|
494 |
|
|
if ( !doAreaFastjet_ && !doRhoFastjet_) {
|
495 |
|
|
fjClusterSeq_ = ClusterSequencePtr( new fastjet::ClusterSequence( fjInputs_, *fjJetDefinition_ ) );
|
496 |
|
|
} else {
|
497 |
|
|
fjClusterSeq_ = ClusterSequencePtr( new fastjet::ClusterSequenceArea( fjInputs_, *fjJetDefinition_ , *fjActiveArea_ ) );
|
498 |
|
|
}
|
499 |
|
|
fjJets_ = fastjet::sorted_by_pt(fjClusterSeq_->inclusive_jets(jetPtMin_));
|
500 |
|
|
|
501 |
|
|
}
|
502 |
|
|
|
503 |
|
|
|
504 |
|
|
|
505 |
|
|
|
506 |
|
|
|
507 |
|
|
////////////////////////////////////////////////////////////////////////////////
|
508 |
|
|
// define as cmssw plugin
|
509 |
|
|
////////////////////////////////////////////////////////////////////////////////
|
510 |
|
|
|
511 |
|
|
DEFINE_FWK_MODULE(JetAlgorithmAnalyzer);
|
512 |
|
|
|