CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_2_9_HLT1_bphpatch4/src/RecoJets/JetProducers/plugins/FastjetJetProducer.cc

Go to the documentation of this file.
00001 
00002 //
00003 // FastjetJetProducer
00004 // ------------------
00005 //
00006 //            04/21/2009 Philipp Schieferdecker <philipp.schieferdecker@cern.ch>
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 
00039 
00040 #include <iostream>
00041 #include <memory>
00042 #include <algorithm>
00043 #include <limits>
00044 #include <cmath>
00045 
00046 using namespace std;
00047 
00048 
00049 
00051 // construction / destruction
00053 
00054 //______________________________________________________________________________
00055 FastjetJetProducer::FastjetJetProducer(const edm::ParameterSet& iConfig)
00056   : VirtualJetProducer( iConfig )
00057 {
00058   
00059   if ( iConfig.exists("UseOnlyVertexTracks") )
00060     useOnlyVertexTracks_ = iConfig.getParameter<bool>("UseOnlyVertexTracks");
00061   else 
00062     useOnlyVertexTracks_ = false;
00063   
00064   if ( iConfig.exists("UseOnlyOnePV") )
00065     useOnlyOnePV_        = iConfig.getParameter<bool>("UseOnlyOnePV");
00066   else
00067     useOnlyOnePV_ = false;
00068 
00069   if ( iConfig.exists("DzTrVtxMax") )
00070     dzTrVtxMax_          = iConfig.getParameter<double>("DzTrVtxMax");
00071   else
00072     dzTrVtxMax_ = 999999.;
00073   if ( iConfig.exists("DxyTrVtxMax") )
00074     dxyTrVtxMax_          = iConfig.getParameter<double>("DxyTrVtxMax");
00075   else
00076     dxyTrVtxMax_ = 999999.;
00077   if ( iConfig.exists("MinVtxNdof") )
00078     minVtxNdof_ = iConfig.getParameter<int>("MinVtxNdof");
00079   else
00080     minVtxNdof_ = 5;
00081   if ( iConfig.exists("MaxVtxZ") )
00082     maxVtxZ_ = iConfig.getParameter<double>("MaxVtxZ");
00083   else
00084     maxVtxZ_ = 15;
00085 
00086 }
00087 
00088 
00089 //______________________________________________________________________________
00090 FastjetJetProducer::~FastjetJetProducer()
00091 {
00092 } 
00093 
00094 
00096 // implementation of member functions
00098 
00099 void FastjetJetProducer::produce( edm::Event & iEvent, const edm::EventSetup & iSetup )
00100 {
00101 
00102   // for everything but track jets
00103   if (!makeTrackJet(jetTypeE)) {
00104  
00105      // use the default production from one collection
00106      VirtualJetProducer::produce( iEvent, iSetup );
00107 
00108   } else { // produce trackjets from tracks grouped per primary vertex
00109 
00110     produceTrackJets(iEvent, iSetup);
00111   
00112   }
00113 
00114 }
00115 
00116 
00117 void FastjetJetProducer::produceTrackJets( edm::Event & iEvent, const edm::EventSetup & iSetup )
00118 {
00119 
00120     // read in the track candidates
00121     edm::Handle<edm::View<reco::RecoChargedRefCandidate> > inputsHandle;
00122     iEvent.getByLabel(src_, inputsHandle);
00123     // make collection with pointers so we can play around with it
00124     std::vector<edm::Ptr<reco::RecoChargedRefCandidate> > allInputs;
00125     std::vector<edm::Ptr<reco::Candidate> > origInputs;
00126     for (size_t i = 0; i < inputsHandle->size(); ++i) {
00127       allInputs.push_back(inputsHandle->ptrAt(i));
00128       origInputs.push_back(inputsHandle->ptrAt(i));
00129     }
00130 
00131     // read in the PV collection
00132     edm::Handle<reco::VertexCollection> pvCollection;
00133     iEvent.getByLabel(srcPVs_, pvCollection);
00134     // define the overall output jet container
00135     std::auto_ptr<std::vector<reco::TrackJet> > jets(new std::vector<reco::TrackJet>() );
00136 
00137     // loop over the good vertices, clustering for each vertex separately
00138     for (reco::VertexCollection::const_iterator itVtx = pvCollection->begin(); itVtx != pvCollection->end(); ++itVtx) {
00139       if (itVtx->isFake() || itVtx->ndof() < minVtxNdof_ || fabs(itVtx->z()) > maxVtxZ_) continue;
00140 
00141       // clear the intermediate containers
00142       inputs_.clear();
00143       fjInputs_.clear();
00144       fjJets_.clear();
00145 
00146       // if only vertex-associated tracks should be used
00147       if (useOnlyVertexTracks_) {
00148         // loop over the tracks associated to the vertex
00149         for (reco::Vertex::trackRef_iterator itTr = itVtx->tracks_begin(); itTr != itVtx->tracks_end(); ++itTr) {
00150           // whether a match was found in the track candidate input
00151           bool found = false;
00152           // loop over input track candidates
00153           for (std::vector<edm::Ptr<reco::RecoChargedRefCandidate> >::iterator itIn = allInputs.begin(); itIn != allInputs.end(); ++itIn) {
00154             // match the input track candidate to the track from the vertex
00155             reco::TrackRef trref(itTr->castTo<reco::TrackRef>());
00156             // check if the tracks match
00157             if ((*itIn)->track() == trref) {
00158               found = true;
00159               // add this track candidate to the input for clustering
00160               inputs_.push_back(*itIn);
00161               // erase the track candidate from the total list of input, so we don't reuse it later
00162               allInputs.erase(itIn);
00163               // found the candidate track corresponding to the vertex track, so stop the loop
00164               break;
00165             } // end if match found
00166           } // end loop over input tracks
00167           // give an info message in case no match is found (can happen if candidates are subset of tracks used for clustering)
00168           if (!found) edm::LogInfo("FastjetTrackJetProducer") << "Ignoring a track at vertex which is not in input track collection!";
00169         } // end loop over tracks associated to vertex
00170       // if all inpt track candidates should be used
00171       } else {
00172         // loop over input track candidates
00173         for (std::vector<edm::Ptr<reco::RecoChargedRefCandidate> >::iterator itIn = allInputs.begin(); itIn != allInputs.end(); ++itIn) {
00174           // check if the track is close enough to the vertex
00175           float dz = (*itIn)->track()->dz(itVtx->position());
00176           float dxy = (*itIn)->track()->dxy(itVtx->position());
00177           if (fabs(dz) > dzTrVtxMax_) continue;
00178           if (fabs(dxy) > dxyTrVtxMax_) continue;
00179           bool closervtx = false;
00180           // now loop over the good vertices a second time
00181           for (reco::VertexCollection::const_iterator itVtx2 = pvCollection->begin(); itVtx2 != pvCollection->end(); ++itVtx2) {
00182             if (itVtx->isFake() || itVtx->ndof() < minVtxNdof_ || fabs(itVtx->z()) > maxVtxZ_) continue;
00183             // and check this track is closer to any other vertex (if more than 1 vertex considered)
00184             if (!useOnlyOnePV_ &&
00185                 itVtx != itVtx2 &&
00186                 fabs((*itIn)->track()->dz(itVtx2->position())) < fabs(dz)) {
00187               closervtx = true;
00188               break; // 1 closer vertex makes the track already not matched, so break
00189             }
00190           }
00191           // don't add this track if another vertex is found closer
00192           if (closervtx) continue;
00193           // add this track candidate to the input for clustering
00194           inputs_.push_back(*itIn);
00195           // erase the track candidate from the total list of input, so we don't reuse it later
00196           allInputs.erase(itIn);
00197           // take a step back in the loop since we just erased
00198           --itIn;
00199         }
00200       }
00201 
00202       // convert candidates in inputs_ to fastjet::PseudoJets in fjInputs_
00203       fjInputs_.reserve(inputs_.size());
00204       inputTowers();
00205       LogDebug("FastjetTrackJetProducer") << "Inputted towers\n";
00206 
00207       // run algorithm, using fjInputs_, modifying fjJets_ and allocating fjClusterSeq_
00208       runAlgorithm(iEvent, iSetup);
00209       LogDebug("FastjetTrackJetProducer") << "Ran algorithm\n";
00210 
00211       // convert our jets and add to the overall jet vector
00212       for (unsigned int ijet=0;ijet<fjJets_.size();++ijet) {
00213         // get the constituents from fastjet
00214         std::vector<fastjet::PseudoJet> fjConstituents = sorted_by_pt(fjClusterSeq_->constituents(fjJets_[ijet]));
00215         // convert them to CandidatePtr vector
00216         std::vector<reco::CandidatePtr> constituents = getConstituents(fjConstituents);
00217         // fill the trackjet
00218         reco::TrackJet jet;
00219         // write the specifics to the jet (simultaneously sets 4-vector, vertex).
00220         writeSpecific( jet,
00221                        reco::Particle::LorentzVector(fjJets_[ijet].px(), fjJets_[ijet].py(), fjJets_[ijet].pz(), fjJets_[ijet].E()),
00222                        vertex_, constituents, iSetup);
00223         jet.setJetArea(0);
00224         jet.setPileup(0);
00225         jet.setPrimaryVertex(edm::Ref<reco::VertexCollection>(pvCollection, (int) (itVtx-pvCollection->begin())));
00226         jet.setVertex(itVtx->position());
00227         jets->push_back(jet);
00228       }
00229 
00230       if (useOnlyOnePV_) break; // stop vertex loop if only one vertex asked for
00231     } // end loop over vertices
00232 
00233     // put the jets in the collection
00234     LogDebug("FastjetTrackJetProducer") << "Put " << jets->size() << " jets in the event.\n";
00235     iEvent.put(jets);
00236 
00237 }
00238 
00239 
00240 //______________________________________________________________________________
00241 void FastjetJetProducer::runAlgorithm( edm::Event & iEvent, edm::EventSetup const& iSetup)
00242 {
00243   // run algorithm
00244   if ( !doAreaFastjet_ && !doRhoFastjet_) {
00245     fjClusterSeq_ = ClusterSequencePtr( new fastjet::ClusterSequence( fjInputs_, *fjJetDefinition_ ) );
00246   } else if (voronoiRfact_ <= 0) {
00247     fjClusterSeq_ = ClusterSequencePtr( new fastjet::ClusterSequenceArea( fjInputs_, *fjJetDefinition_ , *fjActiveArea_ ) );
00248   } else {
00249     fjClusterSeq_ = ClusterSequencePtr( new fastjet::ClusterSequenceVoronoiArea( fjInputs_, *fjJetDefinition_ , fastjet::VoronoiAreaSpec(voronoiRfact_) ) );
00250   }
00251 
00252   fjJets_ = fastjet::sorted_by_pt(fjClusterSeq_->inclusive_jets(jetPtMin_));
00253 
00254 }
00255 
00256 
00257 
00259 // define as cmssw plugin
00261 
00262 DEFINE_FWK_MODULE(FastjetJetProducer);
00263