CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_7/src/TauAnalysis/MCEmbeddingTools/plugins/PFCandidateMixer.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:    PFCandidateMixer
00004 // Class:      PFCandidateMixer
00005 // 
00013 //
00014 // Original Author:  Tomasz Maciej Frueboes
00015 //         Created:  Wed Dec  9 16:14:56 CET 2009
00016 // $Id: PFCandidateMixer.cc,v 1.5.2.1 2012/06/04 12:19:00 aburgmei Exp $
00017 //
00018 //
00019 
00020 
00021 // system include files
00022 #include <memory>
00023 
00024 // user include files
00025 #include "FWCore/Framework/interface/Frameworkfwd.h"
00026 #include "FWCore/Framework/interface/EDProducer.h"
00027 
00028 #include "FWCore/Framework/interface/Event.h"
00029 #include "FWCore/Framework/interface/MakerMacros.h"
00030 
00031 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00032 #include <DataFormats/Common/interface/ValueMap.h>
00033 #include <DataFormats/RecoCandidate/interface/RecoCandidate.h>
00034 #include <DataFormats/Candidate/interface/CompositeRefCandidate.h>
00035 #include <DataFormats/MuonReco/interface/Muon.h>
00036 
00037 #include <DataFormats/ParticleFlowCandidate/interface/PFCandidate.h>
00038 #include <DataFormats/ParticleFlowCandidate/interface/PFCandidateFwd.h>
00039 
00040 #include "DataFormats/EgammaCandidates/interface/GsfElectronFwd.h"
00041 #include "DataFormats/EgammaCandidates/interface/GsfElectron.h"
00042 
00043 #include "DataFormats/TrackReco/interface/Track.h"
00044 #include "DataFormats/Math/interface/deltaR.h"
00045 
00046 #include "RecoParticleFlow/PFProducer/interface/GsfElectronEqual.h"
00047 
00048 //
00049 // class decleration
00050 //
00051 
00052 class PFCandidateMixer : public edm::EDProducer {
00053    public:
00054       explicit PFCandidateMixer(const edm::ParameterSet&);
00055       ~PFCandidateMixer();
00056 
00057    private:
00058       virtual void beginJob() ;
00059       virtual void produce(edm::Event&, const edm::EventSetup&);
00060       virtual void endJob() ;
00061 
00062       void mix(edm::Event& iEvent, const edm::Handle<reco::TrackCollection>& trackCol, const edm::Handle<reco::MuonCollection>& muonCol, const edm::Handle<reco::GsfElectronCollection>& electronCol, const reco::PFCandidateCollection& pfIn1, const reco::PFCandidateCollection& pfIn2);
00063      
00064       edm::InputTag _col1;
00065       edm::InputTag _col2;
00066 
00067       edm::InputTag _trackCol;
00068       edm::InputTag _muonCol;
00069       edm::InputTag _electronCol;
00070       // ----------member data ---------------------------
00071 };
00072 
00073 //
00074 // constants, enums and typedefs
00075 //
00076 
00077 
00078 //
00079 // static data member definitions
00080 //
00081 
00082 //
00083 // constructors and destructor
00084 //
00085 PFCandidateMixer::PFCandidateMixer(const edm::ParameterSet& iConfig):
00086    _col1(iConfig.getUntrackedParameter<edm::InputTag>("col1") ),
00087    _col2(iConfig.getUntrackedParameter<edm::InputTag>("col2") ),
00088    _trackCol(iConfig.getUntrackedParameter<edm::InputTag>("trackCol") ),
00089    _muonCol(iConfig.getUntrackedParameter<edm::InputTag>("muons") ),
00090    _electronCol(iConfig.getUntrackedParameter<edm::InputTag>("gsfElectrons"))
00091 {
00092    produces< std::vector< reco::PFCandidate >  >(); 
00093 
00094    if(!_electronCol.label().empty())
00095       produces< edm::ValueMap<edm::Ptr<reco::PFCandidate> > >("electrons"); 
00096 #if 0
00097    if(!_muonCol.label().empty())
00098       produces< edm::ValueMap<edm::Ptr<reco::PFCandidate> > >("muons"); 
00099 #endif
00100 }
00101 
00102 PFCandidateMixer::~PFCandidateMixer()
00103 {
00104  
00105    // do anything here that needs to be done at desctruction time
00106    // (e.g. close files, deallocate resources etc.)
00107 
00108 }
00109 
00110 
00111 //
00112 // member functions
00113 //
00114 
00115 // ------------ method called to produce the data  ------------
00116 void
00117 PFCandidateMixer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup)
00118 {
00119    const bool muonColGiven = !_muonCol.label().empty();
00120    const bool electronColGiven = !_electronCol.label().empty();
00121 
00122    edm::Handle< std::vector<reco::Track> > trackCol;
00123    edm::Handle<reco::MuonCollection> muonCol;
00124    edm::Handle<reco::GsfElectronCollection> electronCol;
00125    iEvent.getByLabel( _trackCol, trackCol);
00126 
00127    if(muonColGiven)
00128       if(!iEvent.getByLabel(_muonCol, muonCol))
00129          throw cms::Exception("PFCandidateMixer") << "Muon Collection not found!";
00130    if(electronColGiven)
00131       if(!iEvent.getByLabel(_electronCol, electronCol))
00132          throw cms::Exception("PFCandidateMixer") << "GsfElectron Collection not found!";
00133 
00134    edm::Handle<reco::PFCandidateCollection> pfIn1;
00135    edm::Handle<reco::PFCandidateCollection> pfIn2;
00136    iEvent.getByLabel(_col1, pfIn1);
00137    iEvent.getByLabel(_col2, pfIn2);
00138 
00139    mix(iEvent, trackCol, muonCol, electronCol, *pfIn1, *pfIn2);
00140 }
00141 
00142 void PFCandidateMixer::mix(edm::Event& iEvent, const edm::Handle<reco::TrackCollection>& trackCol, const edm::Handle<reco::MuonCollection>& muonCol, const edm::Handle<reco::GsfElectronCollection>& electronCol, const reco::PFCandidateCollection& pfIn1, const reco::PFCandidateCollection& pfIn2)
00143 {
00144    using namespace edm;
00145    using namespace reco;
00146 
00147    std::vector<const reco::PFCandidateCollection*> colVec;
00148    
00149    colVec.push_back(&pfIn1);
00150    colVec.push_back(&pfIn2);
00151 
00152    std::auto_ptr<std::vector< reco::PFCandidate > > pOut(new std::vector< reco::PFCandidate  > );
00153    
00154    std::vector<const reco::PFCandidateCollection*>::iterator itCol= colVec.begin();
00155    std::vector<const reco::PFCandidateCollection*>::iterator itColE= colVec.end();
00156 
00157    std::map<reco::GsfElectronRef, reco::PFCandidatePtr> electronCandidateMap; 
00158 
00159    int iCol = 0;
00160    for (;itCol!=itColE; ++itCol){
00161      PFCandidateConstIterator it = (*itCol)->begin();
00162      PFCandidateConstIterator itE = (*itCol)->end();
00163      for (;it!=itE;++it) {
00164       edm::Ptr<reco::PFCandidate> candPtr(*itCol, it - (*itCol)->begin());
00165       reco::PFCandidate cand(*it);
00166       //reco::PFCandidate cand(candPtr); // This breaks the output module, I'm not sure why
00167        
00168        //const bool isphoton   = cand.particleId() == reco::PFCandidate::gamma && cand.mva_nothing_gamma()>0.;
00169        const bool iselectron = cand.particleId() == reco::PFCandidate::e;
00170        //const bool hasNonNullMuonRef  = cand.muonRef().isNonnull();
00171 
00172        // if it is an electron. Find the GsfElectron with the same GsfTrack
00173        if (electronCol.isValid() && iselectron) {
00174          const reco::GsfTrackRef& gsfTrackRef(cand.gsfTrackRef());
00175          GsfElectronEqual myEqual(gsfTrackRef);
00176          std::vector<reco::GsfElectron>::const_iterator itcheck=find_if(electronCol->begin(), electronCol->end(), myEqual);
00177          if(itcheck==electronCol->end())
00178            throw cms::Exception("PFCandidateMixer") << "GsfElectron for candidate not found!";
00179 
00180          reco::GsfElectronRef electronRef(electronCol, itcheck - electronCol->begin());
00181          cand.setGsfElectronRef(electronRef);
00182          cand.setSuperClusterRef(electronRef->superCluster());
00183          electronCandidateMap[electronRef] = candPtr;
00184        }
00185 
00186        // TODO: Do the same for muons -- see PFLinker.cc, we just don't have the MuonToMuonMap available.
00187        // But it might work with the pure muons collections as well.
00188 
00189        size_t i = 0;
00190        bool found = false;
00191        double minDR = 9999.;
00192        int iMinDr = -1;
00193        if (it->trackRef().isNonnull()) {
00194          for ( i = 0 ; i < trackCol->size(); ++i){
00195            if ( reco::deltaR( *(it->trackRef()), (*trackCol)[i] )<0.001 ) {
00196                 found = true;
00197                 break; 
00198            }
00199            double dr = reco::deltaR( *(it->trackRef()), (*trackCol)[i] );
00200            if ( dr < minDR) {
00201               iMinDr = i;
00202               minDR = dr;
00203            } 
00204          } 
00205        } 
00206        if ( found ){ // ref was found, overwrite in PFCand
00207          reco::TrackRef trref(trackCol,i);
00208          cand.setTrackRef(trref);
00209          //std::cout << " YY track ok"<<std::endl;
00210 
00211        } else { // keep orginall ref
00212          if (it->trackRef().isNonnull()) {
00213            std::cout << " XXXXXXXXXXX track not found " 
00214                  << " col " << iCol
00215                  << " ch " << it->charge()
00216                  << " id " << it->pdgId() 
00217                  << " pt " << it->pt() 
00218                  << " track: eta " << it->trackRef()->eta()
00219                  << " pt:  " << it->trackRef()->pt()
00220                  << " charge:  " << it->trackRef()->charge()
00221                  <<  std::endl;
00222            std::cout << " minDR=" << minDR << std::endl; 
00223            if ( iMinDr >= 0 ) {
00224                 std::cout 
00225                      << " closest track pt=" << (*trackCol)[iMinDr].pt()
00226                      << " ch=" << (*trackCol)[iMinDr].charge()
00227                      <<  std::endl; 
00228            } 
00229            edm::Provenance prov=iEvent.getProvenance(it->trackRef().id());
00230            edm::InputTag tag(prov.moduleLabel(),  prov.productInstanceName(),   prov.processName());
00231            std::cout << " trackref in PFCand came from: "   << tag.encode() << std::endl;
00232          }
00233        }
00234        pOut->push_back(cand);
00235      }
00236      ++iCol;
00237    }
00238 
00239    edm::OrphanHandle<reco::PFCandidateCollection> newColl = iEvent.put(pOut);
00240 
00241    // Now fixup the references and write the valuemap
00242    if(electronCol.isValid())
00243    {
00244       std::vector<reco::PFCandidatePtr> values(electronCol->size());
00245       for(unsigned int i = 0; i < electronCol->size(); ++i)
00246       {
00247          edm::Ref<reco::GsfElectronCollection> objRef(electronCol, i);
00248          std::map<reco::GsfElectronRef, reco::PFCandidatePtr>::const_iterator iter = electronCandidateMap.find(objRef);
00249 
00250          reco::PFCandidatePtr candPtr;
00251          if(iter != electronCandidateMap.end())
00252             candPtr = reco::PFCandidatePtr(newColl, iter->second.key());
00253          values[i] = candPtr;
00254       }
00255 
00256       std::auto_ptr<edm::ValueMap<reco::PFCandidatePtr> > pfMap_p(new edm::ValueMap<reco::PFCandidatePtr>());
00257       edm::ValueMap<reco::PFCandidatePtr>::Filler filler(*pfMap_p);
00258       filler.insert(electronCol, values.begin(), values.end());
00259       filler.fill();
00260       iEvent.put(pfMap_p, "electrons");
00261    }
00262 
00263    // TODO: Do the same for muons
00264 }
00265 
00266 // ------------ method called once each job just before starting event loop  ------------
00267 void 
00268 PFCandidateMixer::beginJob()
00269 {
00270 }
00271 
00272 // ------------ method called once each job just after ending the event loop  ------------
00273 void 
00274 PFCandidateMixer::endJob() {
00275 }
00276 
00277 //define this as a plug-in
00278 DEFINE_FWK_MODULE(PFCandidateMixer);