CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_4_5_patch3/src/RecoParticleFlow/PFProducer/plugins/PFLinker.cc

Go to the documentation of this file.
00001 #include "RecoParticleFlow/PFProducer/plugins/PFLinker.h"
00002 
00003 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidate.h"
00004 #include "DataFormats/EgammaCandidates/interface/GsfElectron.h"
00005 #include "DataFormats/MuonReco/interface/MuonToMuonMap.h"
00006 #include "DataFormats/MuonReco/interface/Muon.h"
00007 #include "DataFormats/MuonReco/interface/MuonFwd.h"
00008 
00009 #include "RecoParticleFlow/PFProducer/interface/GsfElectronEqual.h"
00010 #include "RecoParticleFlow/PFProducer/interface/PhotonEqual.h"
00011 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00012 
00013 
00014 PFLinker::PFLinker(const edm::ParameterSet & iConfig) {
00015   // vector of InputTag; more than 1 is not for RECO, it is for analysis
00016   inputTagPFCandidates_ 
00017     = iConfig.getParameter<std::vector<edm::InputTag> >("PFCandidate");
00018   inputTagGsfElectrons_
00019     = iConfig.getParameter<edm::InputTag>("GsfElectrons");
00020   inputTagPhotons_
00021     = iConfig.getParameter<edm::InputTag>("Photons");
00022   inputTagMuons_
00023     = iConfig.getParameter<edm::InputTag>("Muons");
00024   
00025   nameOutputPF_ 
00026     = iConfig.getParameter<std::string>("OutputPF");
00027   
00028   nameOutputElectronsPF_ 
00029     = iConfig.getParameter<std::string>("ValueMapElectrons");
00030 
00031   nameOutputPhotonsPF_ 
00032     = iConfig.getParameter<std::string>("ValueMapPhotons");
00033 
00034   producePFCandidates_  
00035     = iConfig.getParameter<bool>("ProducePFCandidates");
00036 
00037   nameOutputMergedPF_ 
00038     = iConfig.getParameter<std::string>("ValueMapMerged"); 
00039   
00040   fillMuonRefs_
00041     = iConfig.getParameter<bool>("FillMuonRefs");
00042 
00043   // should not produce PFCandidates and read seve
00044   if(producePFCandidates_ && inputTagPFCandidates_.size()>1) {
00045     edm::LogError("PFLinker") << " cannot read several collections of PFCandidates and produce a new collection at the same time. " << std::endl;
00046     assert(0);
00047   }
00048   if(producePFCandidates_) {
00049     produces<reco::PFCandidateCollection>(nameOutputPF_);
00050   }
00051   produces<edm::ValueMap<reco::PFCandidatePtr> > (nameOutputElectronsPF_);
00052   produces<edm::ValueMap<reco::PFCandidatePtr> > (nameOutputPhotonsPF_);
00053   produces<edm::ValueMap<reco::PFCandidatePtr> > (nameOutputMergedPF_);
00054   if(fillMuonRefs_)  produces<edm::ValueMap<reco::PFCandidatePtr> > (inputTagMuons_.label());
00055 
00056 }
00057 
00058 PFLinker::~PFLinker() {;}
00059 
00060 void PFLinker::beginRun(edm::Run& run,const edm::EventSetup & es) {;}
00061 
00062 void PFLinker::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
00063   
00064   std::auto_ptr<reco::PFCandidateCollection>
00065     pfCandidates_p(new reco::PFCandidateCollection);
00066   
00067   edm::Handle<reco::GsfElectronCollection> gsfElectrons;
00068   bool status=fetchCollection<reco::GsfElectronCollection>(gsfElectrons,
00069                                                       inputTagGsfElectrons_,
00070                                                       iEvent );
00071   std::map<reco::GsfElectronRef,reco::PFCandidatePtr> electronCandidateMap;  
00072 
00073 
00074   edm::Handle<reco::PhotonCollection> photons;
00075   status=fetchCollection<reco::PhotonCollection>(photons,
00076                                                  inputTagPhotons_,
00077                                                  iEvent );
00078   std::map<reco::PhotonRef,reco::PFCandidatePtr> photonCandidateMap;
00079 
00080 
00081   edm::Handle<reco::MuonToMuonMap> muonMap;
00082   if(fillMuonRefs_)
00083     status=fetchCollection<reco::MuonToMuonMap>(muonMap,
00084                                                 inputTagMuons_,
00085                                                 iEvent);
00086   std::map<reco::MuonRef,reco::PFCandidatePtr> muonCandidateMap;
00087   
00088   unsigned nColPF=inputTagPFCandidates_.size();
00089   edm::Handle<reco::PFCandidateCollection> pfCandidates;
00090   for(unsigned icol=0;icol<nColPF;++icol) {
00091     
00092     bool status=fetchCollection<reco::PFCandidateCollection>(pfCandidates, 
00093                                                              inputTagPFCandidates_[icol], 
00094                                                              iEvent );
00095     
00096     unsigned ncand=(status)?pfCandidates->size():0;
00097     
00098     for( unsigned i=0; i<ncand; ++i) {
00099       edm::Ptr<reco::PFCandidate> candPtr(pfCandidates,i);
00100       reco::PFCandidate cand(candPtr);
00101       
00102       bool isphoton   = cand.particleId() == reco::PFCandidate::gamma && cand.mva_nothing_gamma()>0.;
00103       bool iselectron = cand.particleId() == reco::PFCandidate::e;
00104       // PFCandidates may have a valid MuonRef though they are not muons.
00105       bool hasNonNullMuonRef  = cand.muonRef().isNonnull() && fillMuonRefs_;
00106       
00107       // if not an electron or a photon or a muon just fill the PFCandidate collection
00108       if ( !(isphoton || iselectron || hasNonNullMuonRef)){pfCandidates_p->push_back(cand); continue;}
00109       
00110       
00111       if (hasNonNullMuonRef) {
00112         reco::MuonRef muRef = (*muonMap)[cand.muonRef()];
00113         cand.setMuonRef(muRef);
00114         muonCandidateMap[muRef] = candPtr;
00115       }
00116                  
00117       
00118       // if it is an electron. Find the GsfElectron with the same GsfTrack
00119       if (iselectron) {
00120         const reco::GsfTrackRef & gsfTrackRef(cand.gsfTrackRef());
00121         GsfElectronEqual myEqual(gsfTrackRef);
00122         std::vector<reco::GsfElectron>::const_iterator itcheck=find_if(gsfElectrons->begin(),gsfElectrons->end(),myEqual);
00123         if(itcheck==gsfElectrons->end()) {
00124           std::ostringstream err;
00125           err << " Problem in PFLinker: no GsfElectron " << std::endl;
00126           edm::LogError("PFLinker") << err.str();
00127           continue; // Watch out ! Continue
00128         } 
00129         reco::GsfElectronRef electronRef(gsfElectrons,itcheck-gsfElectrons->begin());
00130         cand.setGsfElectronRef(electronRef);
00131         cand.setSuperClusterRef(electronRef->superCluster());
00132         electronCandidateMap[electronRef] = candPtr;
00133       }  
00134       
00135       // if it is a photon, find the one with the same PF super-cluster
00136       if (isphoton) {
00137         const reco::SuperClusterRef & scRef(cand.superClusterRef());
00138         PhotonEqual myEqual(scRef);
00139         std::vector<reco::Photon>::const_iterator itcheck=find_if(photons->begin(),photons->end(),myEqual);
00140         if(itcheck==photons->end()) {
00141           std::ostringstream err;
00142           err << " Problem in PFLinker: no Photon " << std::endl;
00143           edm::LogError("PFLinker") << err.str();
00144           continue; // Watch out ! Continue
00145         } 
00146         reco::PhotonRef photonRef(photons,itcheck-photons->begin());
00147         cand.setPhotonRef(photonRef);
00148         cand.setSuperClusterRef(photonRef->superCluster());
00149         photonCandidateMap[photonRef] = candPtr;
00150       }      
00151 
00152       pfCandidates_p->push_back(cand);
00153     }      
00154     // save the PFCandidates and get a valid handle
00155   }
00156   const edm::OrphanHandle<reco::PFCandidateCollection> pfCandidateRefProd = (producePFCandidates_) ? iEvent.put(pfCandidates_p,nameOutputPF_) :
00157     edm::OrphanHandle<reco::PFCandidateCollection>();
00158   
00159   
00160   // now make the valuemaps
00161 
00162   edm::ValueMap<reco::PFCandidatePtr> pfMapGsfElectrons = fillValueMap<reco::GsfElectronCollection>(iEvent, 
00163                                                                 nameOutputElectronsPF_, 
00164                                                                 gsfElectrons, 
00165                                                                 electronCandidateMap,
00166                                                                 pfCandidateRefProd);
00167   
00168   edm::ValueMap<reco::PFCandidatePtr> pfMapPhotons = fillValueMap<reco::PhotonCollection>(iEvent, 
00169                                                       nameOutputPhotonsPF_, 
00170                                                       photons, 
00171                                                       photonCandidateMap,
00172                                                       pfCandidateRefProd);
00173   
00174 
00175   edm::ValueMap<reco::PFCandidatePtr> pfMapMuons;
00176 
00177   if(fillMuonRefs_){
00178     edm::Handle<reco::MuonCollection> muons; 
00179     iEvent.getByLabel(inputTagMuons_.label(), muons);
00180     
00181     pfMapMuons = fillValueMap<reco::MuonCollection>(iEvent, 
00182                                                     inputTagMuons_.label(), 
00183                                                     muons, 
00184                                                     muonCandidateMap,
00185                                                     pfCandidateRefProd);
00186   }
00187   
00188   std::auto_ptr<edm::ValueMap<reco::PFCandidatePtr> > 
00189     pfMapMerged(new edm::ValueMap<reco::PFCandidatePtr>());
00190   edm::ValueMap<reco::PFCandidatePtr>::Filler pfMapMergedFiller(*pfMapMerged);
00191   
00192   *pfMapMerged                   += pfMapGsfElectrons;
00193   *pfMapMerged                   += pfMapPhotons;
00194   if(fillMuonRefs_) *pfMapMerged += pfMapMuons;
00195   
00196   iEvent.put(pfMapMerged,nameOutputMergedPF_);
00197 }
00198 
00199 template<typename T>
00200 bool PFLinker::fetchCollection(edm::Handle<T>& c, 
00201                                const edm::InputTag& tag, 
00202                                const edm::Event& iEvent) const {  
00203 
00204   bool found = iEvent.getByLabel(tag, c);
00205   
00206   if(!found )
00207     {
00208       std::ostringstream  err;
00209       err<<" cannot get " <<tag<<std::endl;
00210       edm::LogError("PFLinker")<<err.str();
00211     }
00212   return found;
00213 }
00214 
00215 
00216 
00217 template<typename TYPE>
00218 edm::ValueMap<reco::PFCandidatePtr>  PFLinker::fillValueMap(edm::Event & event,
00219                                                             std::string label,
00220                                                             edm::Handle<TYPE>& inputObjCollection,
00221                                                             const std::map<edm::Ref<TYPE>, reco::PFCandidatePtr> & mapToTheCandidate,
00222                                                             const edm::OrphanHandle<reco::PFCandidateCollection> & newPFCandColl) const {
00223   
00224   std::auto_ptr<edm::ValueMap<reco::PFCandidatePtr> > pfMap_p(new edm::ValueMap<reco::PFCandidatePtr>());
00225   edm::ValueMap<reco::PFCandidatePtr>::Filler filler(*pfMap_p);
00226   
00227   typedef typename std::map<edm::Ref<TYPE>, reco::PFCandidatePtr>::const_iterator MapTYPE_it; 
00228   
00229   unsigned nObj=inputObjCollection->size();
00230   std::vector<reco::PFCandidatePtr> values(nObj);
00231 
00232   for(unsigned iobj=0; iobj < nObj; ++iobj) {
00233 
00234     edm::Ref<TYPE> objRef(inputObjCollection, iobj);
00235     MapTYPE_it  itcheck = mapToTheCandidate.find(objRef);
00236 
00237     reco::PFCandidatePtr candPtr;
00238 
00239     if(itcheck != mapToTheCandidate.end())
00240       candPtr = producePFCandidates_ ? reco::PFCandidatePtr(newPFCandColl,itcheck->second.key()) : itcheck->second;
00241     
00242     values[iobj] = candPtr;    
00243   }
00244 
00245   filler.insert(inputObjCollection,values.begin(),values.end());
00246   filler.fill();
00247   edm::ValueMap<reco::PFCandidatePtr> returnValue = *pfMap_p;
00248   event.put(pfMap_p,label);
00249   return returnValue;
00250 }