CMS 3D CMS Logo

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