Go to the documentation of this file.00001
00002
00003
00004
00005
00006
00008
00009 #include "RecoJets/JetProducers/plugins/FastjetJetProducer.h"
00010
00011 #include "RecoJets/JetProducers/interface/JetSpecific.h"
00012
00013 #include "FWCore/Framework/interface/Event.h"
00014 #include "FWCore/Framework/interface/EventSetup.h"
00015 #include "FWCore/Framework/interface/ESHandle.h"
00016 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00017 #include "FWCore/Framework/interface/MakerMacros.h"
00018
00019 #include "DataFormats/Common/interface/Handle.h"
00020 #include "DataFormats/VertexReco/interface/Vertex.h"
00021 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00022 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
00023 #include "DataFormats/JetReco/interface/GenJetCollection.h"
00024 #include "DataFormats/JetReco/interface/PFJetCollection.h"
00025 #include "DataFormats/JetReco/interface/BasicJetCollection.h"
00026 #include "DataFormats/Candidate/interface/CandidateFwd.h"
00027 #include "DataFormats/RecoCandidate/interface/RecoChargedCandidate.h"
00028 #include "DataFormats/Candidate/interface/LeafCandidate.h"
00029
00030
00031 #include "Geometry/CaloGeometry/interface/CaloGeometry.h"
00032 #include "Geometry/Records/interface/CaloGeometryRecord.h"
00033
00034 #include "fastjet/SISConePlugin.hh"
00035 #include "fastjet/CMSIterativeConePlugin.hh"
00036 #include "fastjet/ATLASConePlugin.hh"
00037 #include "fastjet/CDFMidPointPlugin.hh"
00038 #include "fastjet/tools/Filter.hh"
00039 #include "fastjet/tools/Pruner.hh"
00040 #include "fastjet/tools/MassDropTagger.hh"
00041
00042 #include <iostream>
00043 #include <memory>
00044 #include <algorithm>
00045 #include <limits>
00046 #include <cmath>
00047
00048
00049 using namespace std;
00050
00051
00052
00054
00056
00057
00058 FastjetJetProducer::FastjetJetProducer(const edm::ParameterSet& iConfig)
00059 : VirtualJetProducer( iConfig ),
00060 useMassDropTagger_(false),
00061 useFiltering_(false),
00062 useTrimming_(false),
00063 usePruning_(false),
00064 muCut_(-1.0),
00065 yCut_(-1.0),
00066 rFilt_(-1.0),
00067 nFilt_(-1),
00068 trimPtFracMin_(-1.0),
00069 zCut_(-1.0),
00070 RcutFactor_(-1.0)
00071 {
00072
00073 if ( iConfig.exists("UseOnlyVertexTracks") )
00074 useOnlyVertexTracks_ = iConfig.getParameter<bool>("UseOnlyVertexTracks");
00075 else
00076 useOnlyVertexTracks_ = false;
00077
00078 if ( iConfig.exists("UseOnlyOnePV") )
00079 useOnlyOnePV_ = iConfig.getParameter<bool>("UseOnlyOnePV");
00080 else
00081 useOnlyOnePV_ = false;
00082
00083 if ( iConfig.exists("DzTrVtxMax") )
00084 dzTrVtxMax_ = iConfig.getParameter<double>("DzTrVtxMax");
00085 else
00086 dzTrVtxMax_ = 999999.;
00087 if ( iConfig.exists("DxyTrVtxMax") )
00088 dxyTrVtxMax_ = iConfig.getParameter<double>("DxyTrVtxMax");
00089 else
00090 dxyTrVtxMax_ = 999999.;
00091 if ( iConfig.exists("MinVtxNdof") )
00092 minVtxNdof_ = iConfig.getParameter<int>("MinVtxNdof");
00093 else
00094 minVtxNdof_ = 5;
00095 if ( iConfig.exists("MaxVtxZ") )
00096 maxVtxZ_ = iConfig.getParameter<double>("MaxVtxZ");
00097 else
00098 maxVtxZ_ = 15;
00099
00100
00101 if ( iConfig.exists("useFiltering") ||
00102 iConfig.exists("useTrimming") ||
00103 iConfig.exists("usePruning") ||
00104 iConfig.exists("useMassDropTagger") ) {
00105 useMassDropTagger_=false;
00106 useFiltering_=false;
00107 useTrimming_=false;
00108 usePruning_=false;
00109 rFilt_=-1.0;
00110 nFilt_=-1;
00111 trimPtFracMin_=-1.0;
00112 zCut_=-1.0;
00113 RcutFactor_=-1.0;
00114 muCut_=-1.0;
00115 yCut_=-1.0;
00116 useExplicitGhosts_ = true;
00117
00118
00119 if ( iConfig.exists("useMassDropTagger") ) {
00120 useMassDropTagger_ = true;
00121 muCut_ = iConfig.getParameter<double>("muCut");
00122 yCut_ = iConfig.getParameter<double>("yCut");
00123 }
00124
00125 if ( iConfig.exists("useFiltering") ) {
00126 useFiltering_ = true;
00127 rFilt_ = iConfig.getParameter<double>("rFilt");
00128 nFilt_ = iConfig.getParameter<int>("nFilt");
00129 }
00130
00131 if ( iConfig.exists("useTrimming") ) {
00132 useTrimming_ = true;
00133 rFilt_ = iConfig.getParameter<double>("rFilt");
00134 trimPtFracMin_ = iConfig.getParameter<double>("trimPtFracMin");
00135 }
00136
00137 if ( iConfig.exists("usePruning") ) {
00138 usePruning_ = true;
00139 zCut_ = iConfig.getParameter<double>("zcut");
00140 RcutFactor_ = iConfig.getParameter<double>("rcut_factor");
00141 nFilt_ = iConfig.getParameter<int>("nFilt");
00142 }
00143
00144 }
00145
00146 }
00147
00148
00149
00150 FastjetJetProducer::~FastjetJetProducer()
00151 {
00152 }
00153
00154
00156
00158
00159 void FastjetJetProducer::produce( edm::Event & iEvent, const edm::EventSetup & iSetup )
00160 {
00161
00162
00163 if (!makeTrackJet(jetTypeE)) {
00164
00165
00166 VirtualJetProducer::produce( iEvent, iSetup );
00167
00168 } else {
00169
00170 produceTrackJets(iEvent, iSetup);
00171
00172 }
00173
00174 }
00175
00176
00177 void FastjetJetProducer::produceTrackJets( edm::Event & iEvent, const edm::EventSetup & iSetup )
00178 {
00179
00180
00181 edm::Handle<edm::View<reco::RecoChargedRefCandidate> > inputsHandle;
00182 iEvent.getByLabel(src_, inputsHandle);
00183
00184 std::vector<edm::Ptr<reco::RecoChargedRefCandidate> > allInputs;
00185 std::vector<edm::Ptr<reco::Candidate> > origInputs;
00186 for (size_t i = 0; i < inputsHandle->size(); ++i) {
00187 allInputs.push_back(inputsHandle->ptrAt(i));
00188 origInputs.push_back(inputsHandle->ptrAt(i));
00189 }
00190
00191
00192 edm::Handle<reco::VertexCollection> pvCollection;
00193 iEvent.getByLabel(srcPVs_, pvCollection);
00194
00195 std::auto_ptr<std::vector<reco::TrackJet> > jets(new std::vector<reco::TrackJet>() );
00196
00197
00198 for (reco::VertexCollection::const_iterator itVtx = pvCollection->begin(); itVtx != pvCollection->end(); ++itVtx) {
00199 if (itVtx->isFake() || itVtx->ndof() < minVtxNdof_ || fabs(itVtx->z()) > maxVtxZ_) continue;
00200
00201
00202 inputs_.clear();
00203 fjInputs_.clear();
00204 fjJets_.clear();
00205
00206
00207 if (useOnlyVertexTracks_) {
00208
00209 for (reco::Vertex::trackRef_iterator itTr = itVtx->tracks_begin(); itTr != itVtx->tracks_end(); ++itTr) {
00210
00211 bool found = false;
00212
00213 for (std::vector<edm::Ptr<reco::RecoChargedRefCandidate> >::iterator itIn = allInputs.begin(); itIn != allInputs.end(); ++itIn) {
00214
00215 reco::TrackRef trref(itTr->castTo<reco::TrackRef>());
00216
00217 if ((*itIn)->track() == trref) {
00218 found = true;
00219
00220 inputs_.push_back(*itIn);
00221
00222 allInputs.erase(itIn);
00223
00224 break;
00225 }
00226 }
00227
00228 if (!found) edm::LogInfo("FastjetTrackJetProducer") << "Ignoring a track at vertex which is not in input track collection!";
00229 }
00230
00231 } else {
00232
00233 for (std::vector<edm::Ptr<reco::RecoChargedRefCandidate> >::iterator itIn = allInputs.begin(); itIn != allInputs.end(); ++itIn) {
00234
00235 float dz = (*itIn)->track()->dz(itVtx->position());
00236 float dxy = (*itIn)->track()->dxy(itVtx->position());
00237 if (fabs(dz) > dzTrVtxMax_) continue;
00238 if (fabs(dxy) > dxyTrVtxMax_) continue;
00239 bool closervtx = false;
00240
00241 for (reco::VertexCollection::const_iterator itVtx2 = pvCollection->begin(); itVtx2 != pvCollection->end(); ++itVtx2) {
00242 if (itVtx->isFake() || itVtx->ndof() < minVtxNdof_ || fabs(itVtx->z()) > maxVtxZ_) continue;
00243
00244 if (!useOnlyOnePV_ &&
00245 itVtx != itVtx2 &&
00246 fabs((*itIn)->track()->dz(itVtx2->position())) < fabs(dz)) {
00247 closervtx = true;
00248 break;
00249 }
00250 }
00251
00252 if (closervtx) continue;
00253
00254 inputs_.push_back(*itIn);
00255
00256 allInputs.erase(itIn);
00257
00258 --itIn;
00259 }
00260 }
00261
00262
00263 fjInputs_.reserve(inputs_.size());
00264 inputTowers();
00265 LogDebug("FastjetTrackJetProducer") << "Inputted towers\n";
00266
00267
00268 runAlgorithm(iEvent, iSetup);
00269 LogDebug("FastjetTrackJetProducer") << "Ran algorithm\n";
00270
00271
00272 for (unsigned int ijet=0;ijet<fjJets_.size();++ijet) {
00273
00274 std::vector<fastjet::PseudoJet> fjConstituents = sorted_by_pt(fjClusterSeq_->constituents(fjJets_[ijet]));
00275
00276 std::vector<reco::CandidatePtr> constituents = getConstituents(fjConstituents);
00277
00278 reco::TrackJet jet;
00279
00280 writeSpecific( jet,
00281 reco::Particle::LorentzVector(fjJets_[ijet].px(), fjJets_[ijet].py(), fjJets_[ijet].pz(), fjJets_[ijet].E()),
00282 vertex_, constituents, iSetup);
00283 jet.setJetArea(0);
00284 jet.setPileup(0);
00285 jet.setPrimaryVertex(edm::Ref<reco::VertexCollection>(pvCollection, (int) (itVtx-pvCollection->begin())));
00286 jet.setVertex(itVtx->position());
00287 jets->push_back(jet);
00288 }
00289
00290 if (useOnlyOnePV_) break;
00291 }
00292
00293
00294 LogDebug("FastjetTrackJetProducer") << "Put " << jets->size() << " jets in the event.\n";
00295 iEvent.put(jets);
00296
00297 }
00298
00299
00300
00301 void FastjetJetProducer::runAlgorithm( edm::Event & iEvent, edm::EventSetup const& iSetup)
00302 {
00303
00304
00305
00306
00307
00308
00309
00310
00311
00312
00313
00314
00315
00316
00317
00318
00319
00320
00321
00322
00323 if ( !doAreaFastjet_ && !doRhoFastjet_) {
00324 fjClusterSeq_ = ClusterSequencePtr( new fastjet::ClusterSequence( fjInputs_, *fjJetDefinition_ ) );
00325 } else if (voronoiRfact_ <= 0) {
00326 fjClusterSeq_ = ClusterSequencePtr( new fastjet::ClusterSequenceArea( fjInputs_, *fjJetDefinition_ , *fjAreaDefinition_ ) );
00327 } else {
00328 fjClusterSeq_ = ClusterSequencePtr( new fastjet::ClusterSequenceVoronoiArea( fjInputs_, *fjJetDefinition_ , fastjet::VoronoiAreaSpec(voronoiRfact_) ) );
00329 }
00330
00331 if ( !useTrimming_ && !useFiltering_ && !usePruning_ ) {
00332 fjJets_ = fastjet::sorted_by_pt(fjClusterSeq_->inclusive_jets(jetPtMin_));
00333 }
00334 else {
00335 fjJets_.clear();
00336 std::vector<fastjet::PseudoJet> tempJets = fastjet::sorted_by_pt(fjClusterSeq_->inclusive_jets(jetPtMin_));
00337
00338 fastjet::MassDropTagger md_tagger( muCut_, yCut_ );
00339 fastjet::Filter trimmer( fastjet::Filter(fastjet::JetDefinition(fastjet::kt_algorithm, rFilt_), fastjet::SelectorPtFractionMin(trimPtFracMin_)));
00340 fastjet::Filter filter( fastjet::Filter(fastjet::JetDefinition(fastjet::cambridge_algorithm, rFilt_), fastjet::SelectorNHardest(nFilt_)));
00341 fastjet::Pruner pruner(fastjet::cambridge_algorithm, zCut_, RcutFactor_);
00342
00343 std::vector<fastjet::Transformer const *> transformers;
00344
00345 if ( useMassDropTagger_ ) {
00346 transformers.push_back(&md_tagger);
00347 }
00348 if ( useTrimming_ ) {
00349 transformers.push_back(&trimmer);
00350 }
00351 if ( useFiltering_ ) {
00352 transformers.push_back(&filter);
00353 }
00354 if ( usePruning_ ) {
00355 transformers.push_back(&pruner);
00356 }
00357
00358
00359
00360 for ( std::vector<fastjet::PseudoJet>::const_iterator ijet = tempJets.begin(),
00361 ijetEnd = tempJets.end(); ijet != ijetEnd; ++ijet ) {
00362
00363 fastjet::PseudoJet transformedJet = *ijet;
00364 bool passed = true;
00365 for ( std::vector<fastjet::Transformer const *>::const_iterator itransf = transformers.begin(),
00366 itransfEnd = transformers.end(); itransf != itransfEnd; ++itransf ) {
00367 if ( transformedJet != 0 ) {
00368 transformedJet = (**itransf)(transformedJet);
00369 }
00370 else {
00371 passed=false;
00372 }
00373 }
00374 if ( passed )
00375 fjJets_.push_back( transformedJet );
00376 }
00377 }
00378
00379 }
00380
00381
00382
00384
00386
00387 DEFINE_FWK_MODULE(FastjetJetProducer);
00388