00001
00002
00003
00004
00005
00006
00008
00009 #include "RecoJets/JetProducers/plugins/VirtualJetProducer.h"
00010 #include "RecoJets/JetProducers/interface/JetSpecific.h"
00011
00012 #include "FWCore/Framework/interface/Event.h"
00013 #include "FWCore/Framework/interface/EventSetup.h"
00014 #include "FWCore/Framework/interface/ESHandle.h"
00015 #include "FWCore/Utilities/interface/CodedException.h"
00016 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00017 #include "FWCore/Framework/interface/MakerMacros.h"
00018
00019 #include "DataFormats/Common/interface/View.h"
00020 #include "DataFormats/Common/interface/Handle.h"
00021 #include "DataFormats/VertexReco/interface/Vertex.h"
00022 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00023 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
00024 #include "DataFormats/JetReco/interface/GenJetCollection.h"
00025 #include "DataFormats/JetReco/interface/PFJetCollection.h"
00026 #include "DataFormats/JetReco/interface/BasicJetCollection.h"
00027 #include "DataFormats/JetReco/interface/TrackJetCollection.h"
00028 #include "DataFormats/Candidate/interface/CandidateFwd.h"
00029 #include "DataFormats/Candidate/interface/LeafCandidate.h"
00030
00031 #include "fastjet/SISConePlugin.hh"
00032 #include "fastjet/CMSIterativeConePlugin.hh"
00033 #include "fastjet/ATLASConePlugin.hh"
00034 #include "fastjet/CDFMidPointPlugin.hh"
00035
00036 #include <iostream>
00037 #include <memory>
00038 #include <algorithm>
00039 #include <limits>
00040 #include <cmath>
00041
00042 using namespace std;
00043
00044
00045 namespace reco {
00046 namespace helper {
00047 struct GreaterByPtPseudoJet {
00048 bool operator()( const fastjet::PseudoJet & t1, const fastjet::PseudoJet & t2 ) const {
00049 return t1.perp() > t2.perp();
00050 }
00051 };
00052
00053 }
00054 }
00055
00056
00057 const char *VirtualJetProducer::JetType::names[] = {
00058 "BasicJet","GenJet","CaloJet","PFJet","TrackJet"
00059 };
00060
00061
00062
00063 VirtualJetProducer::JetType::Type
00064 VirtualJetProducer::JetType::byName(const string &name)
00065 {
00066 const char **pos = std::find(names, names + LastJetType, name);
00067 if (pos == names + LastJetType) {
00068 std::string errorMessage="Requested jetType not supported: "+name+"\n";
00069 throw cms::Exception("Configuration",errorMessage);
00070 }
00071 return (Type)(pos-names);
00072 }
00073
00074
00075 void VirtualJetProducer::makeProduces( std::string alias, std::string tag )
00076 {
00077 if (makeCaloJet(jetTypeE)) {
00078 produces<reco::CaloJetCollection>(tag).setBranchAlias(alias);
00079 }
00080 else if (makePFJet(jetTypeE)) {
00081 produces<reco::PFJetCollection>(tag).setBranchAlias(alias);
00082 }
00083 else if (makeGenJet(jetTypeE)) {
00084 produces<reco::GenJetCollection>(tag).setBranchAlias(alias);
00085 }
00086 else if (makeTrackJet(jetTypeE)) {
00087 produces<reco::TrackJetCollection>(tag).setBranchAlias(alias);
00088 }
00089 else if (makeBasicJet(jetTypeE)) {
00090 produces<reco::BasicJetCollection>(tag).setBranchAlias(alias);
00091 }
00092 }
00093
00095
00097
00098
00099 VirtualJetProducer::VirtualJetProducer(const edm::ParameterSet& iConfig)
00100 : moduleLabel_ (iConfig.getParameter<string> ("@module_label"))
00101 , src_ (iConfig.getParameter<edm::InputTag>("src"))
00102 , srcPVs_ (iConfig.getParameter<edm::InputTag>("srcPVs"))
00103 , jetType_ (iConfig.getParameter<string> ("jetType"))
00104 , jetAlgorithm_ (iConfig.getParameter<string> ("jetAlgorithm"))
00105 , rParam_ (iConfig.getParameter<double> ("rParam"))
00106 , inputEtMin_ (iConfig.getParameter<double> ("inputEtMin"))
00107 , inputEMin_ (iConfig.getParameter<double> ("inputEMin"))
00108 , jetPtMin_ (iConfig.getParameter<double> ("jetPtMin"))
00109 , doPVCorrection_(iConfig.getParameter<bool> ("doPVCorrection"))
00110 , restrictInputs_(false)
00111 , maxInputs_(99999999)
00112 , doAreaFastjet_ (iConfig.getParameter<bool> ("doAreaFastjet"))
00113 , doRhoFastjet_ (iConfig.getParameter<bool> ("doRhoFastjet"))
00114 , doPUOffsetCorr_(iConfig.getParameter<bool> ("doPUOffsetCorr"))
00115 , maxBadEcalCells_ (iConfig.getParameter<unsigned>("maxBadEcalCells"))
00116 , maxRecoveredEcalCells_ (iConfig.getParameter<unsigned>("maxRecoveredEcalCells"))
00117 , maxProblematicEcalCells_(iConfig.getParameter<unsigned>("maxProblematicEcalCells"))
00118 , maxBadHcalCells_ (iConfig.getParameter<unsigned>("maxBadHcalCells"))
00119 , maxRecoveredHcalCells_ (iConfig.getParameter<unsigned>("maxRecoveredHcalCells"))
00120 , maxProblematicHcalCells_(iConfig.getParameter<unsigned>("maxProblematicHcalCells"))
00121 , jetCollInstanceName_ ("")
00122 {
00123
00124
00125
00126
00127
00128
00129
00130 if (jetAlgorithm_=="SISCone") {
00131 fjPlugin_ = PluginPtr( new fastjet::SISConePlugin(rParam_,0.75,0,0.0,false,
00132 fastjet::SISConePlugin::SM_pttilde) );
00133 fjJetDefinition_= JetDefPtr( new fastjet::JetDefinition(&*fjPlugin_) );
00134 }
00135 else if (jetAlgorithm_=="IterativeCone") {
00136 fjPlugin_ = PluginPtr(new fastjet::CMSIterativeConePlugin(rParam_,1.0));
00137 fjJetDefinition_= JetDefPtr(new fastjet::JetDefinition(&*fjPlugin_));
00138 }
00139 else if (jetAlgorithm_=="CDFMidPoint") {
00140 fjPlugin_ = PluginPtr(new fastjet::CDFMidPointPlugin(rParam_,0.75));
00141 fjJetDefinition_= JetDefPtr(new fastjet::JetDefinition(&*fjPlugin_));
00142 }
00143 else if (jetAlgorithm_=="ATLASCone") {
00144 fjPlugin_ = PluginPtr(new fastjet::ATLASConePlugin(rParam_));
00145 fjJetDefinition_= JetDefPtr(new fastjet::JetDefinition(&*fjPlugin_));
00146 }
00147 else if (jetAlgorithm_=="Kt")
00148 fjJetDefinition_= JetDefPtr(new fastjet::JetDefinition(fastjet::kt_algorithm,rParam_));
00149 else if (jetAlgorithm_=="CambridgeAachen")
00150 fjJetDefinition_= JetDefPtr(new fastjet::JetDefinition(fastjet::cambridge_algorithm,
00151 rParam_) );
00152 else if (jetAlgorithm_=="AntiKt")
00153 fjJetDefinition_= JetDefPtr( new fastjet::JetDefinition(fastjet::antikt_algorithm,rParam_) );
00154 else if (jetAlgorithm_=="GeneralizedKt")
00155 fjJetDefinition_= JetDefPtr( new fastjet::JetDefinition(fastjet::genkt_algorithm,
00156 rParam_,-2) );
00157 else
00158 throw cms::Exception("Invalid jetAlgorithm")
00159 <<"Jet algorithm for VirtualJetProducer is invalid, Abort!\n";
00160
00161 jetTypeE=JetType::byName(jetType_);
00162
00163 if ( iConfig.exists("jetCollInstanceName") ) {
00164 jetCollInstanceName_ = iConfig.getParameter<string>("jetCollInstanceName");
00165 }
00166
00167 if ( doPUOffsetCorr_ ) {
00168 if ( jetTypeE != JetType::CaloJet ) {
00169 throw cms::Exception("InvalidInput") << "Can only offset correct jets of type CaloJet";
00170 }
00171
00172 if(iConfig.exists("subtractorName")) puSubtractorName_ = iConfig.getParameter<string> ("subtractorName");
00173 else puSubtractorName_ = std::string();
00174
00175 if(puSubtractorName_.empty()){
00176 edm::LogWarning("VirtualJetProducer") << "Pile Up correction on; however, pile up type is not specified. Using default... \n";
00177 subtractor_ = boost::shared_ptr<PileUpSubtractor>(new PileUpSubtractor(iConfig));
00178 }else{
00179 subtractor_ = boost::shared_ptr<PileUpSubtractor>(PileUpSubtractorFactory::get()->create( puSubtractorName_, iConfig));
00180 }
00181 }
00182
00183
00184 if ( doAreaFastjet_ || doRhoFastjet_ ) {
00185
00186
00187 double rhoEtaMax=iConfig.getParameter<double>("Rho_EtaMax");
00188
00189 double ghostEtaMax = iConfig.getParameter<double>("Ghost_EtaMax");
00190
00191 int activeAreaRepeats = iConfig.getParameter<int> ("Active_Area_Repeats");
00192
00193 double ghostArea = iConfig.getParameter<double> ("GhostArea");
00194 fjActiveArea_ = ActiveAreaSpecPtr(new fastjet::ActiveAreaSpec(ghostEtaMax,
00195 activeAreaRepeats,
00196 ghostArea));
00197 fjRangeDef_ = RangeDefPtr( new fastjet::RangeDefinition(rhoEtaMax) );
00198 }
00199
00200
00201 if ( iConfig.exists("restrictInputs") ) {
00202 restrictInputs_ = iConfig.getParameter<bool>("restrictInputs");
00203 maxInputs_ = iConfig.getParameter<unsigned int>("maxInputs");
00204 }
00205
00206
00207 string alias=iConfig.getUntrackedParameter<string>("alias",moduleLabel_);
00208
00209
00210 makeProduces( alias, jetCollInstanceName_ );
00211
00212 doFastJetNonUniform_ = false;
00213 if(iConfig.exists("doFastJetNonUniform")) doFastJetNonUniform_ = iConfig.getParameter<bool> ("doFastJetNonUniform");
00214 if(doFastJetNonUniform_){
00215 puCenters_ = iConfig.getParameter<std::vector<double> >("puCenters");
00216 puWidth_ = iConfig.getParameter<double>("puWidth");
00217 }
00218 produces<std::vector<double> >("rhos");
00219 produces<std::vector<double> >("sigmas");
00220 produces<double>("rho");
00221 produces<double>("sigma");
00222
00223 }
00224
00225
00226 VirtualJetProducer::~VirtualJetProducer()
00227 {
00228 }
00229
00230
00232
00234
00235
00236 void VirtualJetProducer::produce(edm::Event& iEvent,const edm::EventSetup& iSetup)
00237 {
00238 LogDebug("VirtualJetProducer") << "Entered produce\n";
00239
00240 vertex_=reco::Jet::Point(0,0,0);
00241 if (makeCaloJet(jetTypeE)&&doPVCorrection_) {
00242 LogDebug("VirtualJetProducer") << "Adding PV info\n";
00243 edm::Handle<reco::VertexCollection> pvCollection;
00244 iEvent.getByLabel(srcPVs_,pvCollection);
00245 if (pvCollection->size()>0) vertex_=pvCollection->begin()->position();
00246 }
00247
00248
00249
00250 if ( doPUOffsetCorr_ ) {
00251 subtractor_->setupGeometryMap(iEvent, iSetup);
00252 }
00253
00254
00255 LogDebug("VirtualJetProducer") << "Clear data\n";
00256 fjInputs_.clear();
00257 fjJets_.clear();
00258 inputs_.clear();
00259
00260
00261 edm::Handle<reco::CandidateView> inputsHandle;
00262 iEvent.getByLabel(src_,inputsHandle);
00263 for (size_t i = 0; i < inputsHandle->size(); ++i) {
00264 inputs_.push_back(inputsHandle->ptrAt(i));
00265 }
00266 LogDebug("VirtualJetProducer") << "Got inputs\n";
00267
00268
00269
00270
00271 fjInputs_.reserve(inputs_.size());
00272 inputTowers();
00273 LogDebug("VirtualJetProducer") << "Inputted towers\n";
00274
00275
00276
00277 if ( doPUOffsetCorr_ ) {
00278 subtractor_->reset(inputs_,fjInputs_,fjJets_);
00279 subtractor_->calculatePedestal(fjInputs_);
00280 subtractor_->subtractPedestal(fjInputs_);
00281 LogDebug("VirtualJetProducer") << "Subtracted pedestal\n";
00282 }
00283
00284
00285
00286 runAlgorithm( iEvent, iSetup );
00287 if ( doPUOffsetCorr_ ) {
00288 subtractor_->setAlgorithm(fjClusterSeq_);
00289 }
00290
00291 LogDebug("VirtualJetProducer") << "Ran algorithm\n";
00292
00293
00294
00295
00296
00297
00298
00299 if ( doPUOffsetCorr_ ) {
00300 vector<fastjet::PseudoJet> orphanInput;
00301 subtractor_->calculateOrphanInput(orphanInput);
00302 subtractor_->calculatePedestal(orphanInput);
00303 subtractor_->offsetCorrectJets();
00304 }
00305
00306
00307
00308
00309
00310 output( iEvent, iSetup );
00311 LogDebug("VirtualJetProducer") << "Wrote jets\n";
00312
00313 return;
00314 }
00315
00316
00317
00318 void VirtualJetProducer::inputTowers( )
00319 {
00320 std::vector<edm::Ptr<reco::Candidate> >::const_iterator inBegin = inputs_.begin(),
00321 inEnd = inputs_.end(), i = inBegin;
00322 for (; i != inEnd; ++i ) {
00323 reco::CandidatePtr input = *i;
00324 if (isnan(input->pt())) continue;
00325 if (input->et() <inputEtMin_) continue;
00326 if (input->energy()<inputEMin_) continue;
00327 if (isAnomalousTower(input)) continue;
00328 if (makeCaloJet(jetTypeE)&&doPVCorrection_) {
00329 const CaloTower* tower=dynamic_cast<const CaloTower*>(input.get());
00330 math::PtEtaPhiMLorentzVector ct(tower->p4(vertex_));
00331 fjInputs_.push_back(fastjet::PseudoJet(ct.px(),ct.py(),ct.pz(),ct.energy()));
00332 }
00333 else {
00334 fjInputs_.push_back(fastjet::PseudoJet(input->px(),input->py(),input->pz(),
00335 input->energy()));
00336 }
00337 fjInputs_.back().set_user_index(i - inBegin);
00338 }
00339
00340 if ( restrictInputs_ && fjInputs_.size() > maxInputs_ ) {
00341 reco::helper::GreaterByPtPseudoJet pTComparator;
00342 std::sort(fjInputs_.begin(), fjInputs_.end(), pTComparator);
00343 fjInputs_.resize(maxInputs_);
00344 edm::LogWarning("JetRecoTooManyEntries") << "Too many inputs in the event, limiting to first " << maxInputs_ << ". Output is suspect.";
00345 }
00346 }
00347
00348
00349 bool VirtualJetProducer::isAnomalousTower(reco::CandidatePtr input)
00350 {
00351 if (!makeCaloJet(jetTypeE)) return false;
00352 const CaloTower* tower=dynamic_cast<const CaloTower*>(input.get());
00353 if (0==tower) return false;
00354 if (tower->numBadEcalCells() >maxBadEcalCells_ ||
00355 tower->numRecoveredEcalCells() >maxRecoveredEcalCells_ ||
00356 tower->numProblematicEcalCells()>maxProblematicEcalCells_||
00357 tower->numBadHcalCells() >maxBadHcalCells_ ||
00358 tower->numRecoveredHcalCells() >maxRecoveredHcalCells_ ||
00359 tower->numProblematicHcalCells()>maxProblematicHcalCells_) return true;
00360 return false;
00361 }
00362
00363
00364
00365
00366
00367
00368
00369
00370
00371 void VirtualJetProducer::copyConstituents(const vector<fastjet::PseudoJet>& fjConstituents,
00372 reco::Jet* jet)
00373 {
00374 for (unsigned int i=0;i<fjConstituents.size();++i)
00375 jet->addDaughter(inputs_[fjConstituents[i].user_index()]);
00376 }
00377
00378
00379
00380 vector<reco::CandidatePtr>
00381 VirtualJetProducer::getConstituents(const vector<fastjet::PseudoJet>&fjConstituents)
00382 {
00383 vector<reco::CandidatePtr> result;
00384 for (unsigned int i=0;i<fjConstituents.size();i++) {
00385 int index = fjConstituents[i].user_index();
00386 reco::CandidatePtr candidate = inputs_[index];
00387 result.push_back(candidate);
00388 }
00389 return result;
00390 }
00391
00392 void VirtualJetProducer::output(edm::Event & iEvent, edm::EventSetup const& iSetup)
00393 {
00394
00395
00396 switch( jetTypeE ) {
00397 case JetType::CaloJet :
00398 writeJets<reco::CaloJet>( iEvent, iSetup);
00399 break;
00400 case JetType::PFJet :
00401 writeJets<reco::PFJet>( iEvent, iSetup);
00402 break;
00403 case JetType::GenJet :
00404 writeJets<reco::GenJet>( iEvent, iSetup);
00405 break;
00406 case JetType::TrackJet :
00407 writeJets<reco::TrackJet>( iEvent, iSetup);
00408 break;
00409 case JetType::BasicJet :
00410 writeJets<reco::BasicJet>( iEvent, iSetup);
00411 break;
00412 default:
00413 throw cms::Exception("InvalidInput") << "invalid jet type in VirtualJetProducer\n";
00414 break;
00415 };
00416
00417 }
00418
00419 template< typename T >
00420 void VirtualJetProducer::writeJets( edm::Event & iEvent, edm::EventSetup const& iSetup )
00421 {
00422
00423
00424 using namespace reco;
00425
00426 std::auto_ptr<std::vector<T> > jets(new std::vector<T>() );
00427 jets->reserve(fjJets_.size());
00428
00429 for (unsigned int ijet=0;ijet<fjJets_.size();++ijet) {
00430
00431 T jet;
00432
00433 const fastjet::PseudoJet& fjJet = fjJets_[ijet];
00434
00435 std::vector<fastjet::PseudoJet> fjConstituents =
00436 sorted_by_pt(fjClusterSeq_->constituents(fjJet));
00437
00438 std::vector<CandidatePtr> constituents =
00439 getConstituents(fjConstituents);
00440
00441
00442 double jetArea=0.0;
00443 if ( doAreaFastjet_ ) {
00444 fastjet::ClusterSequenceArea const * clusterSequenceWithArea =
00445 dynamic_cast<fastjet::ClusterSequenceArea const *>(&*fjClusterSeq_);
00446 jetArea = clusterSequenceWithArea->area(fjJet);
00447 }
00448
00449
00450
00451
00452 writeSpecific(jet,
00453 Particle::LorentzVector(fjJet.px(),
00454 fjJet.py(),
00455 fjJet.pz(),
00456 fjJet.E()),
00457 vertex_,
00458 constituents, iSetup);
00459
00460 jet.setJetArea (jetArea);
00461
00462 if(doPUOffsetCorr_){
00463 jet.setPileup(subtractor_->getPileUpEnergy(ijet));
00464 }else{
00465 jet.setPileup (0.0);
00466 }
00467
00468
00469 jets->push_back(jet);
00470 }
00471
00472
00473 iEvent.put(jets);
00474
00475 if (doRhoFastjet_) {
00476 if(doFastJetNonUniform_){
00477 std::auto_ptr<std::vector<double> > rhos(new std::vector<double>);
00478 std::auto_ptr<std::vector<double> > sigmas(new std::vector<double>);
00479 int nEta = puCenters_.size();
00480 rhos->reserve(nEta);
00481 sigmas->reserve(nEta);
00482 fastjet::ClusterSequenceArea const * clusterSequenceWithArea =
00483 dynamic_cast<fastjet::ClusterSequenceArea const *> ( &*fjClusterSeq_ );
00484
00485 for(int ie = 0; ie < nEta; ++ie){
00486 double eta = puCenters_[ie];
00487 double rho = 0;
00488 double sigma = 0;
00489 double etamin=eta-puWidth_;
00490 double etamax=eta+puWidth_;
00491 fastjet::RangeDefinition range_rho(etamin,etamax);
00492 clusterSequenceWithArea->get_median_rho_and_sigma(range_rho, true, rho, sigma);
00493 rhos->push_back(rho);
00494 sigmas->push_back(sigma);
00495 }
00496 iEvent.put(rhos,"rhos");
00497 iEvent.put(sigmas,"sigmas");
00498 }else{
00499 std::auto_ptr<double> rho(new double(0.0));
00500 std::auto_ptr<double> sigma(new double(0.0));
00501 double mean_area = 0;
00502
00503 fastjet::ClusterSequenceArea const * clusterSequenceWithArea =
00504 dynamic_cast<fastjet::ClusterSequenceArea const *> ( &*fjClusterSeq_ );
00505 clusterSequenceWithArea->get_median_rho_and_sigma(*fjRangeDef_,false,*rho,*sigma,mean_area);
00506 iEvent.put(rho,"rho");
00507 iEvent.put(sigma,"sigma");
00508 }
00509
00510 }
00511 }
00512
00513
00514