ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/CmsHi/JetAnalysis/modifiedFiles/VirtualJetProducer.cc
Revision: 1.1
Committed: Tue Jul 12 14:13:07 2011 UTC (13 years, 10 months ago) by frankma
Content type: text/plain
Branch: MAIN
CVS Tags: tag_d20110915, cmssw39x_base, cmssw39X_base
Branch point for: cmssw39x_branch
Log Message:
pf and calo jet analysis in 39X

File Contents

# Content
1 ////////////////////////////////////////////////////////////////////////////////
2 //
3 // VirtualJetProducer
4 // ------------------
5 //
6 // 04/21/2009 Philipp Schieferdecker <philipp.schieferdecker@cern.ch>
7 ////////////////////////////////////////////////////////////////////////////////
8
9 #include "RecoJets/JetProducers/plugins/VirtualJetProducer.h"
10 #include "RecoJets/JetProducers/interface/JetSpecific.h"
11
12 #include "FWCore/Framework/interface/Event.h"
13 #include "FWCore/Framework/interface/EventSetup.h"
14 #include "FWCore/Framework/interface/ESHandle.h"
15 #include "FWCore/Utilities/interface/CodedException.h"
16 #include "FWCore/MessageLogger/interface/MessageLogger.h"
17 #include "FWCore/Framework/interface/MakerMacros.h"
18
19 #include "DataFormats/Common/interface/View.h"
20 #include "DataFormats/Common/interface/Handle.h"
21 #include "DataFormats/VertexReco/interface/Vertex.h"
22 #include "DataFormats/VertexReco/interface/VertexFwd.h"
23 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
24 #include "DataFormats/JetReco/interface/GenJetCollection.h"
25 #include "DataFormats/JetReco/interface/PFJetCollection.h"
26 #include "DataFormats/JetReco/interface/BasicJetCollection.h"
27 #include "DataFormats/JetReco/interface/TrackJetCollection.h"
28 #include "DataFormats/Candidate/interface/CandidateFwd.h"
29 #include "DataFormats/Candidate/interface/LeafCandidate.h"
30
31 #include "fastjet/SISConePlugin.hh"
32 #include "fastjet/CMSIterativeConePlugin.hh"
33 #include "fastjet/ATLASConePlugin.hh"
34 #include "fastjet/CDFMidPointPlugin.hh"
35
36 #include <iostream>
37 #include <memory>
38 #include <algorithm>
39 #include <limits>
40 #include <cmath>
41
42 using namespace std;
43
44
45 namespace reco {
46 namespace helper {
47 struct GreaterByPtPseudoJet {
48 bool operator()( const fastjet::PseudoJet & t1, const fastjet::PseudoJet & t2 ) const {
49 return t1.perp() > t2.perp();
50 }
51 };
52
53 }
54 }
55
56 //______________________________________________________________________________
57 const char *VirtualJetProducer::JetType::names[] = {
58 "BasicJet","GenJet","CaloJet","PFJet","TrackJet"
59 };
60
61
62 //______________________________________________________________________________
63 VirtualJetProducer::JetType::Type
64 VirtualJetProducer::JetType::byName(const string &name)
65 {
66 const char **pos = std::find(names, names + LastJetType, name);
67 if (pos == names + LastJetType) {
68 std::string errorMessage="Requested jetType not supported: "+name+"\n";
69 throw cms::Exception("Configuration",errorMessage);
70 }
71 return (Type)(pos-names);
72 }
73
74
75 void VirtualJetProducer::makeProduces( std::string alias, std::string tag )
76 {
77 if (makeCaloJet(jetTypeE)) {
78 produces<reco::CaloJetCollection>(tag).setBranchAlias(alias);
79 }
80 else if (makePFJet(jetTypeE)) {
81 produces<reco::PFJetCollection>(tag).setBranchAlias(alias);
82 }
83 else if (makeGenJet(jetTypeE)) {
84 produces<reco::GenJetCollection>(tag).setBranchAlias(alias);
85 }
86 else if (makeTrackJet(jetTypeE)) {
87 produces<reco::TrackJetCollection>(tag).setBranchAlias(alias);
88 }
89 else if (makeBasicJet(jetTypeE)) {
90 produces<reco::BasicJetCollection>(tag).setBranchAlias(alias);
91 }
92 }
93
94 ////////////////////////////////////////////////////////////////////////////////
95 // construction / destruction
96 ////////////////////////////////////////////////////////////////////////////////
97
98 //______________________________________________________________________________
99 VirtualJetProducer::VirtualJetProducer(const edm::ParameterSet& iConfig)
100 : moduleLabel_ (iConfig.getParameter<string> ("@module_label"))
101 , src_ (iConfig.getParameter<edm::InputTag>("src"))
102 , srcPVs_ (iConfig.getParameter<edm::InputTag>("srcPVs"))
103 , jetType_ (iConfig.getParameter<string> ("jetType"))
104 , jetAlgorithm_ (iConfig.getParameter<string> ("jetAlgorithm"))
105 , rParam_ (iConfig.getParameter<double> ("rParam"))
106 , inputEtMin_ (iConfig.getParameter<double> ("inputEtMin"))
107 , inputEMin_ (iConfig.getParameter<double> ("inputEMin"))
108 , jetPtMin_ (iConfig.getParameter<double> ("jetPtMin"))
109 , doPVCorrection_(iConfig.getParameter<bool> ("doPVCorrection"))
110 , restrictInputs_(false)
111 , maxInputs_(99999999)
112 , doAreaFastjet_ (iConfig.getParameter<bool> ("doAreaFastjet"))
113 , doRhoFastjet_ (iConfig.getParameter<bool> ("doRhoFastjet"))
114 , doPUOffsetCorr_(iConfig.getParameter<bool> ("doPUOffsetCorr"))
115 , maxBadEcalCells_ (iConfig.getParameter<unsigned>("maxBadEcalCells"))
116 , maxRecoveredEcalCells_ (iConfig.getParameter<unsigned>("maxRecoveredEcalCells"))
117 , maxProblematicEcalCells_(iConfig.getParameter<unsigned>("maxProblematicEcalCells"))
118 , maxBadHcalCells_ (iConfig.getParameter<unsigned>("maxBadHcalCells"))
119 , maxRecoveredHcalCells_ (iConfig.getParameter<unsigned>("maxRecoveredHcalCells"))
120 , maxProblematicHcalCells_(iConfig.getParameter<unsigned>("maxProblematicHcalCells"))
121 , jetCollInstanceName_ ("")
122 {
123
124 //
125 // additional parameters to think about:
126 // - overlap threshold (set to 0.75 for the time being)
127 // - p parameter for generalized kT (set to -2 for the time being)
128 // - fastjet PU subtraction parameters (not yet considered)
129 //
130 if (jetAlgorithm_=="SISCone") {
131 fjPlugin_ = PluginPtr( new fastjet::SISConePlugin(rParam_,0.75,0,0.0,false,
132 fastjet::SISConePlugin::SM_pttilde) );
133 fjJetDefinition_= JetDefPtr( new fastjet::JetDefinition(&*fjPlugin_) );
134 }
135 else if (jetAlgorithm_=="IterativeCone") {
136 fjPlugin_ = PluginPtr(new fastjet::CMSIterativeConePlugin(rParam_,1.0));
137 fjJetDefinition_= JetDefPtr(new fastjet::JetDefinition(&*fjPlugin_));
138 }
139 else if (jetAlgorithm_=="CDFMidPoint") {
140 fjPlugin_ = PluginPtr(new fastjet::CDFMidPointPlugin(rParam_,0.75));
141 fjJetDefinition_= JetDefPtr(new fastjet::JetDefinition(&*fjPlugin_));
142 }
143 else if (jetAlgorithm_=="ATLASCone") {
144 fjPlugin_ = PluginPtr(new fastjet::ATLASConePlugin(rParam_));
145 fjJetDefinition_= JetDefPtr(new fastjet::JetDefinition(&*fjPlugin_));
146 }
147 else if (jetAlgorithm_=="Kt")
148 fjJetDefinition_= JetDefPtr(new fastjet::JetDefinition(fastjet::kt_algorithm,rParam_));
149 else if (jetAlgorithm_=="CambridgeAachen")
150 fjJetDefinition_= JetDefPtr(new fastjet::JetDefinition(fastjet::cambridge_algorithm,
151 rParam_) );
152 else if (jetAlgorithm_=="AntiKt")
153 fjJetDefinition_= JetDefPtr( new fastjet::JetDefinition(fastjet::antikt_algorithm,rParam_) );
154 else if (jetAlgorithm_=="GeneralizedKt")
155 fjJetDefinition_= JetDefPtr( new fastjet::JetDefinition(fastjet::genkt_algorithm,
156 rParam_,-2) );
157 else
158 throw cms::Exception("Invalid jetAlgorithm")
159 <<"Jet algorithm for VirtualJetProducer is invalid, Abort!\n";
160
161 jetTypeE=JetType::byName(jetType_);
162
163 if ( iConfig.exists("jetCollInstanceName") ) {
164 jetCollInstanceName_ = iConfig.getParameter<string>("jetCollInstanceName");
165 }
166
167 if ( doPUOffsetCorr_ ) {
168 if ( jetTypeE != JetType::CaloJet ) {
169 //throw cms::Exception("InvalidInput") << "Can only offset correct jets of type CaloJet";
170 std::cout<<" Implementing background subtraction using a pseudo calo grid"<<std::endl;
171 }
172
173 if(iConfig.exists("subtractorName")) puSubtractorName_ = iConfig.getParameter<string> ("subtractorName");
174 else puSubtractorName_ = std::string();
175
176 if(puSubtractorName_.empty()){
177 edm::LogWarning("VirtualJetProducer") << "Pile Up correction on; however, pile up type is not specified. Using default... \n";
178 subtractor_ = boost::shared_ptr<PileUpSubtractor>(new PileUpSubtractor(iConfig));
179 }else{
180 subtractor_ = boost::shared_ptr<PileUpSubtractor>(PileUpSubtractorFactory::get()->create( puSubtractorName_, iConfig));
181 }
182 }
183
184 // do fasjet area / rho calcluation? => accept corresponding parameters
185 if ( doAreaFastjet_ || doRhoFastjet_ ) {
186 // Eta range of jets to be considered for Rho calculation
187 // Should be at most (jet acceptance - jet radius)
188 double rhoEtaMax=iConfig.getParameter<double>("Rho_EtaMax");
189 // default Ghost_EtaMax should be 5
190 double ghostEtaMax = iConfig.getParameter<double>("Ghost_EtaMax");
191 // default Active_Area_Repeats 1
192 int activeAreaRepeats = iConfig.getParameter<int> ("Active_Area_Repeats");
193 // default GhostArea 0.01
194 double ghostArea = iConfig.getParameter<double> ("GhostArea");
195 fjActiveArea_ = ActiveAreaSpecPtr(new fastjet::ActiveAreaSpec(ghostEtaMax,
196 activeAreaRepeats,
197 ghostArea));
198 fjRangeDef_ = RangeDefPtr( new fastjet::RangeDefinition(rhoEtaMax) );
199 }
200
201 // restrict inputs to first "maxInputs" towers?
202 if ( iConfig.exists("restrictInputs") ) {
203 restrictInputs_ = iConfig.getParameter<bool>("restrictInputs");
204 maxInputs_ = iConfig.getParameter<unsigned int>("maxInputs");
205 }
206
207
208 string alias=iConfig.getUntrackedParameter<string>("alias",moduleLabel_);
209
210 // make the "produces" statements
211 makeProduces( alias, jetCollInstanceName_ );
212
213 doFastJetNonUniform_ = false;
214 if(iConfig.exists("doFastJetNonUniform")) doFastJetNonUniform_ = iConfig.getParameter<bool> ("doFastJetNonUniform");
215 if(doFastJetNonUniform_){
216 puCenters_ = iConfig.getParameter<std::vector<double> >("puCenters");
217 puWidth_ = iConfig.getParameter<double>("puWidth");
218 }
219 produces<std::vector<double> >("rhos");
220 produces<std::vector<double> >("sigmas");
221 produces<double>("rho");
222 produces<double>("sigma");
223
224 }
225
226 //______________________________________________________________________________
227 VirtualJetProducer::~VirtualJetProducer()
228 {
229 }
230
231
232 ////////////////////////////////////////////////////////////////////////////////
233 // implementation of member functions
234 ////////////////////////////////////////////////////////////////////////////////
235
236 //______________________________________________________________________________
237 void VirtualJetProducer::produce(edm::Event& iEvent,const edm::EventSetup& iSetup)
238 {
239 LogDebug("VirtualJetProducer") << "Entered produce\n";
240 //determine signal vertex
241 vertex_=reco::Jet::Point(0,0,0);
242 if (makeCaloJet(jetTypeE)&&doPVCorrection_) {
243 LogDebug("VirtualJetProducer") << "Adding PV info\n";
244 edm::Handle<reco::VertexCollection> pvCollection;
245 iEvent.getByLabel(srcPVs_,pvCollection);
246 if (pvCollection->size()>0) vertex_=pvCollection->begin()->position();
247 }
248
249 // For Pileup subtraction using offset correction:
250 // set up geometry map
251 if ( doPUOffsetCorr_ ) {
252 subtractor_->setupGeometryMap(iEvent, iSetup);
253 }
254
255 // clear data
256 LogDebug("VirtualJetProducer") << "Clear data\n";
257 fjInputs_.clear();
258 fjJets_.clear();
259 inputs_.clear();
260
261 // get inputs and convert them to the fastjet format (fastjet::PeudoJet)
262 edm::Handle<reco::CandidateView> inputsHandle;
263 iEvent.getByLabel(src_,inputsHandle);
264 for (size_t i = 0; i < inputsHandle->size(); ++i) {
265 inputs_.push_back(inputsHandle->ptrAt(i));
266 }
267 LogDebug("VirtualJetProducer") << "Got inputs\n";
268
269 // Convert candidates to fastjet::PseudoJets.
270 // Also correct to Primary Vertex. Will modify fjInputs_
271 // and use inputs_
272 fjInputs_.reserve(inputs_.size());
273 inputTowers();
274 LogDebug("VirtualJetProducer") << "Inputted towers\n";
275
276 // For Pileup subtraction using offset correction:
277 // Subtract pedestal.
278 if ( doPUOffsetCorr_ ) {
279 subtractor_->reset(inputs_,fjInputs_,fjJets_);
280 subtractor_->calculatePedestal(fjInputs_);
281 subtractor_->subtractPedestal(fjInputs_);
282 LogDebug("VirtualJetProducer") << "Subtracted pedestal\n";
283 }
284
285 // Run algorithm. Will modify fjJets_ and allocate fjClusterSeq_.
286 // This will use fjInputs_
287 runAlgorithm( iEvent, iSetup );
288 if ( doPUOffsetCorr_ ) {
289 subtractor_->setAlgorithm(fjClusterSeq_);
290 }
291
292 LogDebug("VirtualJetProducer") << "Ran algorithm\n";
293
294 // For Pileup subtraction using offset correction:
295 // Now we find jets and need to recalculate their energy,
296 // mark towers participated in jet,
297 // remove occupied towers from the list and recalculate mean and sigma
298 // put the initial towers collection to the jet,
299 // and subtract from initial towers in jet recalculated mean and sigma of towers
300 if ( doPUOffsetCorr_ ) {
301 vector<fastjet::PseudoJet> orphanInput;
302 subtractor_->calculateOrphanInput(orphanInput);
303 subtractor_->calculatePedestal(orphanInput);
304 subtractor_->offsetCorrectJets();
305 }
306
307 // Write the output jets.
308 // This will (by default) call the member function template
309 // "writeJets", but can be overridden.
310 // this will use inputs_
311 output( iEvent, iSetup );
312 LogDebug("VirtualJetProducer") << "Wrote jets\n";
313
314 return;
315 }
316
317 //______________________________________________________________________________
318
319 void VirtualJetProducer::inputTowers( )
320 {
321 std::vector<edm::Ptr<reco::Candidate> >::const_iterator inBegin = inputs_.begin(),
322 inEnd = inputs_.end(), i = inBegin;
323 for (; i != inEnd; ++i ) {
324 reco::CandidatePtr input = *i;
325 if (isnan(input->pt())) continue;
326 if (input->et() <inputEtMin_) continue;
327 if (input->energy()<inputEMin_) continue;
328 if (isAnomalousTower(input)) continue;
329 if (makeCaloJet(jetTypeE)&&doPVCorrection_) {
330 const CaloTower* tower=dynamic_cast<const CaloTower*>(input.get());
331 math::PtEtaPhiMLorentzVector ct(tower->p4(vertex_));
332 fjInputs_.push_back(fastjet::PseudoJet(ct.px(),ct.py(),ct.pz(),ct.energy()));
333 }
334 else {
335 fjInputs_.push_back(fastjet::PseudoJet(input->px(),input->py(),input->pz(),
336 input->energy()));
337 }
338 fjInputs_.back().set_user_index(i - inBegin);
339 }
340
341 if ( restrictInputs_ && fjInputs_.size() > maxInputs_ ) {
342 reco::helper::GreaterByPtPseudoJet pTComparator;
343 std::sort(fjInputs_.begin(), fjInputs_.end(), pTComparator);
344 fjInputs_.resize(maxInputs_);
345 edm::LogWarning("JetRecoTooManyEntries") << "Too many inputs in the event, limiting to first " << maxInputs_ << ". Output is suspect.";
346 }
347 }
348
349 //______________________________________________________________________________
350 bool VirtualJetProducer::isAnomalousTower(reco::CandidatePtr input)
351 {
352 if (!makeCaloJet(jetTypeE)) return false;
353 const CaloTower* tower=dynamic_cast<const CaloTower*>(input.get());
354 if (0==tower) return false;
355 if (tower->numBadEcalCells() >maxBadEcalCells_ ||
356 tower->numRecoveredEcalCells() >maxRecoveredEcalCells_ ||
357 tower->numProblematicEcalCells()>maxProblematicEcalCells_||
358 tower->numBadHcalCells() >maxBadHcalCells_ ||
359 tower->numRecoveredHcalCells() >maxRecoveredHcalCells_ ||
360 tower->numProblematicHcalCells()>maxProblematicHcalCells_) return true;
361 return false;
362 }
363
364
365 //------------------------------------------------------------------------------
366 // This is pure virtual.
367 //______________________________________________________________________________
368 // void VirtualJetProducer::runAlgorithm( edm::Event & iEvent, edm::EventSetup const& iSetup,
369 // std::vector<edm::Ptr<reco::Candidate> > const & inputs_);
370
371 //______________________________________________________________________________
372 void VirtualJetProducer::copyConstituents(const vector<fastjet::PseudoJet>& fjConstituents,
373 reco::Jet* jet)
374 {
375 for (unsigned int i=0;i<fjConstituents.size();++i)
376 jet->addDaughter(inputs_[fjConstituents[i].user_index()]);
377 }
378
379
380 //______________________________________________________________________________
381 vector<reco::CandidatePtr>
382 VirtualJetProducer::getConstituents(const vector<fastjet::PseudoJet>&fjConstituents)
383 {
384 vector<reco::CandidatePtr> result;
385 for (unsigned int i=0;i<fjConstituents.size();i++) {
386 int index = fjConstituents[i].user_index();
387 reco::CandidatePtr candidate = inputs_[index];
388 result.push_back(candidate);
389 }
390 return result;
391 }
392
393 void VirtualJetProducer::output(edm::Event & iEvent, edm::EventSetup const& iSetup)
394 {
395 // Write jets and constitutents. Will use fjJets_, inputs_
396 // and fjClusterSeq_
397 switch( jetTypeE ) {
398 case JetType::CaloJet :
399 writeJets<reco::CaloJet>( iEvent, iSetup);
400 break;
401 case JetType::PFJet :
402 writeJets<reco::PFJet>( iEvent, iSetup);
403 break;
404 case JetType::GenJet :
405 writeJets<reco::GenJet>( iEvent, iSetup);
406 break;
407 case JetType::TrackJet :
408 writeJets<reco::TrackJet>( iEvent, iSetup);
409 break;
410 case JetType::BasicJet :
411 writeJets<reco::BasicJet>( iEvent, iSetup);
412 break;
413 default:
414 throw cms::Exception("InvalidInput") << "invalid jet type in VirtualJetProducer\n";
415 break;
416 };
417
418 }
419
420 template< typename T >
421 void VirtualJetProducer::writeJets( edm::Event & iEvent, edm::EventSetup const& iSetup )
422 {
423 // produce output jet collection
424
425 using namespace reco;
426
427 std::auto_ptr<std::vector<T> > jets(new std::vector<T>() );
428 jets->reserve(fjJets_.size());
429
430 for (unsigned int ijet=0;ijet<fjJets_.size();++ijet) {
431 // allocate this jet
432 T jet;
433 // get the fastjet jet
434 const fastjet::PseudoJet& fjJet = fjJets_[ijet];
435 // get the constituents from fastjet
436 std::vector<fastjet::PseudoJet> fjConstituents =
437 sorted_by_pt(fjClusterSeq_->constituents(fjJet));
438 // convert them to CandidatePtr vector
439 std::vector<CandidatePtr> constituents =
440 getConstituents(fjConstituents);
441
442 // calcuate the jet area
443 double jetArea=0.0;
444 if ( doAreaFastjet_ ) {
445 fastjet::ClusterSequenceArea const * clusterSequenceWithArea =
446 dynamic_cast<fastjet::ClusterSequenceArea const *>(&*fjClusterSeq_);
447 jetArea = clusterSequenceWithArea->area(fjJet);
448 }
449
450 // write the specifics to the jet (simultaneously sets 4-vector, vertex).
451 // These are overridden functions that will call the appropriate
452 // specific allocator.
453 writeSpecific(jet,
454 Particle::LorentzVector(fjJet.px(),
455 fjJet.py(),
456 fjJet.pz(),
457 fjJet.E()),
458 vertex_,
459 constituents, iSetup);
460
461 jet.setJetArea (jetArea);
462
463 if(doPUOffsetCorr_){
464 jet.setPileup(subtractor_->getPileUpEnergy(ijet));
465 }else{
466 jet.setPileup (0.0);
467 }
468
469 // add to the list
470 jets->push_back(jet);
471 }
472
473 // put the jets in the collection
474 iEvent.put(jets);
475
476 if (doRhoFastjet_) {
477 if(doFastJetNonUniform_){
478 std::auto_ptr<std::vector<double> > rhos(new std::vector<double>);
479 std::auto_ptr<std::vector<double> > sigmas(new std::vector<double>);
480 int nEta = puCenters_.size();
481 rhos->reserve(nEta);
482 sigmas->reserve(nEta);
483 fastjet::ClusterSequenceArea const * clusterSequenceWithArea =
484 dynamic_cast<fastjet::ClusterSequenceArea const *> ( &*fjClusterSeq_ );
485
486 for(int ie = 0; ie < nEta; ++ie){
487 double eta = puCenters_[ie];
488 double rho = 0;
489 double sigma = 0;
490 double etamin=eta-puWidth_;
491 double etamax=eta+puWidth_;
492 fastjet::RangeDefinition range_rho(etamin,etamax);
493 clusterSequenceWithArea->get_median_rho_and_sigma(range_rho, true, rho, sigma);
494 rhos->push_back(rho);
495 sigmas->push_back(sigma);
496 }
497 iEvent.put(rhos,"rhos");
498 iEvent.put(sigmas,"sigmas");
499 }else{
500 std::auto_ptr<double> rho(new double(0.0));
501 std::auto_ptr<double> sigma(new double(0.0));
502 double mean_area = 0;
503
504 fastjet::ClusterSequenceArea const * clusterSequenceWithArea =
505 dynamic_cast<fastjet::ClusterSequenceArea const *> ( &*fjClusterSeq_ );
506 clusterSequenceWithArea->get_median_rho_and_sigma(*fjRangeDef_,false,*rho,*sigma,mean_area);
507 iEvent.put(rho,"rho");
508 iEvent.put(sigma,"sigma");
509 }
510
511 }
512 }
513
514
515