CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_3/src/PhysicsTools/TagAndProbe/plugins/DeltaRNearestObjectComputer.cc

Go to the documentation of this file.
00001 /*****************************************************************************
00002  * Project: CMS detector at the CERN
00003  *
00004  * Package: PhysicsTools/TagAndProbe
00005  *
00006  *
00007  * Authors:
00008  *   Giovanni Petrucciani, UCSD - Giovanni.Petrucciani@cern.ch
00009  *
00010  * Description:
00011  *   - Matches a given object with other objects using deltaR-matching.
00012  *   - For example: can match a photon with track within a given deltaR.
00013  *   - Saves collection of the reference vectors of matched objects.
00014  * History:
00015  *   
00016  * Kalanand Mishra, Fermilab - kalanand@fnal.gov
00017  * Extended the class to compute deltaR with respect to any object 
00018  * (i.e., Candidate, Jet, Muon, Electron, or Photon). The previous 
00019  * version of this class could deltaR only with respect to reco::Candidates.
00020  * This didn't work if one wanted to apply selection cuts on the Candidate's 
00021  * RefToBase object.
00022  *****************************************************************************/
00023 
00024 #include "FWCore/Framework/interface/EDProducer.h"
00025 #include "FWCore/Framework/interface/Event.h"
00026 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00027 #include "FWCore/Utilities/interface/InputTag.h"
00028 
00029 #include "DataFormats/Math/interface/deltaR.h"
00030 #include "DataFormats/Common/interface/ValueMap.h"
00031 #include "DataFormats/Common/interface/View.h"
00032 
00033 #include "DataFormats/Candidate/interface/CandidateFwd.h"
00034 #include "DataFormats/Candidate/interface/Candidate.h"
00035 #include "DataFormats/MuonReco/interface/Muon.h"
00036 #include "DataFormats/EgammaCandidates/interface/GsfElectron.h"
00037 #include "DataFormats/EgammaCandidates/interface/Electron.h"
00038 #include "DataFormats/JetReco/interface/Jet.h"
00039 #include "DataFormats/EgammaCandidates/interface/Photon.h"
00040 
00041 #include "CommonTools/Utils/interface/StringCutObjectSelector.h"
00042 
00043 template<typename T>
00044 class DeltaRNearestObjectComputer : public edm::EDProducer {
00045     public:
00046         explicit DeltaRNearestObjectComputer(const edm::ParameterSet & iConfig);
00047         virtual ~DeltaRNearestObjectComputer() ;
00048 
00049         virtual void produce(edm::Event & iEvent, const edm::EventSetup & iSetup);
00050 
00051     private:
00052         edm::InputTag probes_;            
00053         edm::InputTag objects_; 
00054         StringCutObjectSelector<T,true> objCut_; // lazy parsing, to allow cutting on variables not in reco::Candidate class
00055 };
00056 
00057 template<typename T>
00058 DeltaRNearestObjectComputer<T>::DeltaRNearestObjectComputer(const edm::ParameterSet & iConfig) :
00059     probes_(iConfig.getParameter<edm::InputTag>("probes")),
00060     objects_(iConfig.getParameter<edm::InputTag>("objects")),
00061     objCut_(iConfig.existsAs<std::string>("objectSelection") ? iConfig.getParameter<std::string>("objectSelection") : "", true)
00062 {
00063     produces<edm::ValueMap<float> >();
00064 }
00065 
00066 
00067 template<typename T>
00068 DeltaRNearestObjectComputer<T>::~DeltaRNearestObjectComputer()
00069 {
00070 }
00071 
00072 template<typename T>
00073 void 
00074 DeltaRNearestObjectComputer<T>::produce(edm::Event & iEvent, const edm::EventSetup & iSetup) {
00075     using namespace edm;
00076 
00077     // read input
00078     Handle<View<reco::Candidate> > probes;
00079     iEvent.getByLabel(probes_,  probes);
00080 
00081     Handle<View<T> > objects;
00082     iEvent.getByLabel(objects_, objects);
00083 
00084     // prepare vector for output    
00085     std::vector<float> values;
00086     
00087     // fill
00088     View<reco::Candidate>::const_iterator probe, endprobes = probes->end();
00089     for (probe = probes->begin(); probe != endprobes; ++probe) {
00090         double dr2min = 10000;
00091         for (unsigned int iObj=0; iObj<objects->size(); iObj++) {
00092           const T& obj = objects->at(iObj);
00093           if (!objCut_(obj)) continue;
00094             double dr2 = deltaR2(*probe, obj);
00095             if (dr2 < dr2min) { dr2min = dr2; }
00096         }
00097         values.push_back(sqrt(dr2min));
00098     }
00099 
00100     // convert into ValueMap and store
00101     std::auto_ptr<ValueMap<float> > valMap(new ValueMap<float>());
00102     ValueMap<float>::Filler filler(*valMap);
00103     filler.insert(probes, values.begin(), values.end());
00104     filler.fill();
00105     iEvent.put(valMap);
00106 }
00107 
00108 
00110 // plugin definition
00112 
00113 typedef DeltaRNearestObjectComputer<reco::Candidate>     DeltaRNearestCandidateComputer;
00114 typedef DeltaRNearestObjectComputer<reco::Muon>          DeltaRNearestMuonComputer;
00115 typedef DeltaRNearestObjectComputer<reco::Electron>      DeltaRNearestElectronComputer;
00116 typedef DeltaRNearestObjectComputer<reco::GsfElectron>   DeltaRNearestGsfElectronComputer;
00117 typedef DeltaRNearestObjectComputer<reco::Photon>        DeltaRNearestPhotonComputer;
00118 typedef DeltaRNearestObjectComputer<reco::Jet>           DeltaRNearestJetComputer;
00119 
00120 #include "FWCore/Framework/interface/MakerMacros.h"
00121 DEFINE_FWK_MODULE(DeltaRNearestCandidateComputer);          
00122 DEFINE_FWK_MODULE(DeltaRNearestMuonComputer);          
00123 DEFINE_FWK_MODULE(DeltaRNearestElectronComputer);          
00124 DEFINE_FWK_MODULE(DeltaRNearestGsfElectronComputer);          
00125 DEFINE_FWK_MODULE(DeltaRNearestPhotonComputer);          
00126 DEFINE_FWK_MODULE(DeltaRNearestJetComputer);