CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_4/src/TauAnalysis/MCEmbeddingTools/plugins/ZmumuPFEmbedder.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:    ZmumuPFEmbedder
00004 // Class:      ZmumuPFEmbedder
00005 // 
00013 //
00014 // Original Author:  Tomasz Maciej Frueboes
00015 //         Created:  Wed Dec  9 16:14:56 CET 2009
00016 // $Id: ZmumuPFEmbedder.cc,v 1.12 2012/01/27 13:17:12 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/View.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/VertexReco/interface/Vertex.h"
00041 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00042 #include "DataFormats/TrackCandidate/interface/TrackCandidateCollection.h"
00043 
00044 #include "DataFormats/Candidate/interface/CompositeCandidate.h"
00045 #include <DataFormats/Math/interface/deltaR.h>
00046 //
00047 // class decleration
00048 //
00049 
00050 class ZmumuPFEmbedder : public edm::EDProducer {
00051    public:
00052       explicit ZmumuPFEmbedder(const edm::ParameterSet&);
00053       ~ZmumuPFEmbedder();
00054 
00055    private:
00056       virtual void beginJob() ;
00057       virtual void produce(edm::Event&, const edm::EventSetup&);
00058       void producePFCandColl(edm::Event&, const std::vector< reco::Particle::LorentzVector > * toBeAdded );
00059       void produceTrackColl(edm::Event&, const std::vector< reco::Particle::LorentzVector > * toBeAdded );
00060       virtual void endJob() ;
00061       
00062       edm::InputTag _tracks;
00063       edm::InputTag _selectedMuons;
00064       bool _useCombinedCandidate;
00065 
00066       // ----------member data ---------------------------
00067 };
00068 
00069 //
00070 // constants, enums and typedefs
00071 //
00072 
00073 
00074 //
00075 // static data member definitions
00076 //
00077 
00078 //
00079 // constructors and destructor
00080 //
00081 ZmumuPFEmbedder::ZmumuPFEmbedder(const edm::ParameterSet& iConfig)
00082   : _tracks(iConfig.getParameter<edm::InputTag>("tracks")),
00083     _selectedMuons(iConfig.getParameter<edm::InputTag>("selectedMuons")),
00084     _useCombinedCandidate(iConfig.getUntrackedParameter<bool>("useCombinedCandidate", false))
00085 {
00086 
00087    //register your products
00088    // produces< std::vector< reco::Muon >  >("zMusExtracted"); // 
00089    produces<reco::TrackCollection>("tracks");
00090    produces< std::vector< reco::PFCandidate >  >("pfCands");
00091 
00092    
00093 }
00094 
00095 ZmumuPFEmbedder::~ZmumuPFEmbedder()
00096 {
00097  
00098    // do anything here that needs to be done at desctruction time
00099    // (e.g. close files, deallocate resources etc.)
00100 
00101 }
00102 
00103 
00104 //
00105 // member functions
00106 //
00107 
00108 // ------------ method called to produce the data  ------------
00109 void
00110 ZmumuPFEmbedder::produce(edm::Event& iEvent, const edm::EventSetup& iSetup)
00111 {
00112    std::vector< reco::Particle::LorentzVector > toBeAdded;
00113    
00114    if (_useCombinedCandidate)
00115    {
00116       edm::Handle< std::vector< reco::CompositeCandidate > > combCandidatesHandle;
00117       if (iEvent.getByLabel(_selectedMuons, combCandidatesHandle) && combCandidatesHandle->size()>0)
00118          for (size_t idx = 0; idx < combCandidatesHandle->at(0).numberOfDaughters(); ++idx)                     // use only the first combined candidate
00119             toBeAdded.push_back(combCandidatesHandle->at(0).daughter(idx)->p4());
00120    }
00121    else
00122    {
00123       edm::Handle< edm::View< reco::Muon > > selectedZMuonsHandle;
00124       if (iEvent.getByLabel(_selectedMuons, selectedZMuonsHandle))
00125         for (size_t idx = 0; idx < selectedZMuonsHandle->size(); ++idx)
00126           toBeAdded.push_back(selectedZMuonsHandle->at(idx).p4());
00127    }
00128 
00129    if (toBeAdded.size() == 0)
00130       return;
00131 
00132    producePFCandColl(iEvent, &toBeAdded);
00133    produceTrackColl(iEvent, &toBeAdded);
00134 }
00135 
00136 // produces clean PFCandidate collection w/o muon candidates.
00137 void ZmumuPFEmbedder::producePFCandColl(edm::Event & iEvent, const std::vector< reco::Particle::LorentzVector > * toBeAdded )
00138 {
00139    edm::Handle<reco::PFCandidateCollection> pfIn;
00140    iEvent.getByLabel("particleFlow",pfIn);
00141 
00142    //std::vector< reco::Muon > toBeAdded;
00143    // TODO - check col size
00144    //reco::Muon l1 = *(toBeAdded->begin());
00145    //reco::Muon l2 = *(toBeAdded->rbegin());
00146    
00147    std::auto_ptr<std::vector< reco::PFCandidate > > newCol(new std::vector< reco::PFCandidate  > );   
00148    
00149    //get selected muons
00150    // iterate over pfcandidates, make copy if its not a selected muon
00151    reco::PFCandidateConstIterator it = pfIn->begin();
00152    reco::PFCandidateConstIterator itE = pfIn->end();
00153 
00154    for (;it!=itE;++it) {
00155      int pdg = std::abs( it->pdgId() );
00156      double minDR = 10;
00157      /* 
00158      double dr1 = reco::deltaR( *it, l1); 
00159      double dr2 = reco::deltaR( *it, l2); 
00160      */
00161      std::vector< reco::Particle::LorentzVector >::const_iterator itSelectedMu = toBeAdded->begin();
00162      std::vector< reco::Particle::LorentzVector >::const_iterator itSelectedMuE = toBeAdded->end();
00163      for (; itSelectedMu != itSelectedMuE; ++itSelectedMu ){
00164        double dr = reco::deltaR( *it, *itSelectedMu);
00165        if (dr < minDR)  minDR = dr;
00166      }
00167 
00168      //if ( pdg == 13 && (dr1 < 0.001 || dr2 < 0.002 ) ) { // it is a selected muon, do not copy
00169      if ( pdg == 13 && (minDR < 0.001 ) ) { // it is a selected muon, do not copy
00170        
00171        
00172      } else {
00173        newCol->push_back(*it);
00174      }
00175    }
00176 
00177    iEvent.put(newCol, "pfCands");
00178 }
00179 
00180 // produces clean track collection w/o muon tracks.
00181 void ZmumuPFEmbedder::produceTrackColl(edm::Event & iEvent, const std::vector< reco::Particle::LorentzVector > * toBeAdded )
00182 {
00183    edm::Handle<reco::TrackCollection> tks;
00184    iEvent.getByLabel( _tracks, tks);
00185 
00186    std::auto_ptr< reco::TrackCollection  > newCol(new reco::TrackCollection );
00187 
00188    double epsilon = 0.00001;
00189    unsigned int nMatched = 0;
00190    
00191    for ( reco::TrackCollection::const_iterator it = tks->begin() ; it != tks->end(); ++it) 
00192    {
00193      bool ok = true;
00194      for ( std::vector< reco::Particle::LorentzVector >::const_iterator itTBA = toBeAdded->begin();
00195                                                             itTBA != toBeAdded->end();
00196                                                             ++itTBA)
00197      {
00198        double dr = reco::deltaR( *it, *itTBA);
00199        //std::cout << "TTTT " << dr << std::endl;
00200        if (dr < epsilon) {
00201          ++ nMatched;
00202          ok = false;
00203        }
00204      }
00205      if (ok)  newCol->push_back(*it);  
00206    }
00207    if (nMatched!=toBeAdded->size() ) std::cout << "TTT ARGGGHGH " << nMatched << std::endl;
00208 
00209    iEvent.put(newCol, "tracks");
00210 }
00211 // ------------ method called once each job just before starting event loop  ------------
00212 void 
00213 ZmumuPFEmbedder::beginJob()
00214 {
00215 }
00216 
00217 // ------------ method called once each job just after ending the event loop  ------------
00218 void 
00219 ZmumuPFEmbedder::endJob() {
00220 }
00221 
00222 //define this as a plug-in
00223 DEFINE_FWK_MODULE(ZmumuPFEmbedder);