CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_5/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::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
00061   
00062   std::auto_ptr<reco::PFCandidateCollection>
00063     pfCandidates_p(new reco::PFCandidateCollection);
00064   
00065   edm::Handle<reco::GsfElectronCollection> gsfElectrons;
00066   bool status=fetchCollection<reco::GsfElectronCollection>(gsfElectrons,
00067                                                       inputTagGsfElectrons_,
00068                                                       iEvent );
00069   std::map<reco::GsfElectronRef,reco::PFCandidatePtr> electronCandidateMap;  
00070 
00071   if(!status) {
00072     std::ostringstream err;
00073     err << " Problem in PFLinker: no electron collection called " << inputTagGsfElectrons_ << std::endl;
00074     edm::LogError("PFLinker") << err.str();
00075   }
00076 
00077 
00078   edm::Handle<reco::PhotonCollection> photons;
00079   status=fetchCollection<reco::PhotonCollection>(photons,
00080                                                  inputTagPhotons_,
00081                                                  iEvent );
00082   if(!status) {
00083     std::ostringstream err;
00084     err << " Problem in PFLinker: no photon collection called " << inputTagPhotons_ << std::endl;
00085     edm::LogError("PFLinker") << err.str();
00086   }
00087 
00088 
00089   std::map<reco::PhotonRef,reco::PFCandidatePtr> photonCandidateMap;
00090 
00091 
00092   edm::Handle<reco::MuonToMuonMap> muonMap;
00093   if(fillMuonRefs_)
00094     status=fetchCollection<reco::MuonToMuonMap>(muonMap,
00095                                                 inputTagMuons_,
00096                                                 iEvent);
00097   std::map<reco::MuonRef,reco::PFCandidatePtr> muonCandidateMap;
00098   
00099   unsigned nColPF=inputTagPFCandidates_.size();
00100   
00101   if(!status) {
00102     std::ostringstream err;
00103     err << " Problem in PFLinker: no muon collection called " << inputTagMuons_ << std::endl;
00104     edm::LogError("PFLinker") << err.str();
00105   }
00106 
00107   
00108   edm::Handle<reco::PFCandidateCollection> pfCandidates;
00109   for(unsigned icol=0;icol<nColPF;++icol) {
00110     
00111     bool status=fetchCollection<reco::PFCandidateCollection>(pfCandidates, 
00112                                                              inputTagPFCandidates_[icol], 
00113                                                              iEvent );
00114     
00115     unsigned ncand=(status)?pfCandidates->size():0;
00116     
00117     for( unsigned i=0; i<ncand; ++i) {
00118       edm::Ptr<reco::PFCandidate> candPtr(pfCandidates,i);
00119       reco::PFCandidate cand(candPtr);
00120       
00121       bool isphoton   = cand.particleId() == reco::PFCandidate::gamma && cand.mva_nothing_gamma()>0.;
00122       bool iselectron = cand.particleId() == reco::PFCandidate::e;
00123       // PFCandidates may have a valid MuonRef though they are not muons.
00124       bool hasNonNullMuonRef  = cand.muonRef().isNonnull() && fillMuonRefs_;
00125       
00126       // if not an electron or a photon or a muon just fill the PFCandidate collection
00127       if ( !(isphoton || iselectron || hasNonNullMuonRef)){pfCandidates_p->push_back(cand); continue;}
00128       
00129       
00130       if (hasNonNullMuonRef) {
00131         reco::MuonRef muRef = (*muonMap)[cand.muonRef()];
00132         cand.setMuonRef(muRef);
00133         muonCandidateMap[muRef] = candPtr;
00134       }
00135                  
00136       
00137       // if it is an electron. Find the GsfElectron with the same GsfTrack
00138       if (iselectron) {
00139         const reco::GsfTrackRef & gsfTrackRef(cand.gsfTrackRef());
00140         GsfElectronEqual myEqual(gsfTrackRef);
00141         std::vector<reco::GsfElectron>::const_iterator itcheck=find_if(gsfElectrons->begin(),gsfElectrons->end(),myEqual);
00142         if(itcheck==gsfElectrons->end()) {
00143           std::ostringstream err;
00144           err << " Problem in PFLinker: no GsfElectron " << std::endl;
00145           edm::LogError("PFLinker") << err.str();
00146           continue; // Watch out ! Continue
00147         } 
00148         reco::GsfElectronRef electronRef(gsfElectrons,itcheck-gsfElectrons->begin());
00149         cand.setGsfElectronRef(electronRef);
00150         cand.setSuperClusterRef(electronRef->superCluster());
00151         electronCandidateMap[electronRef] = candPtr;
00152       }  
00153       
00154       // if it is a photon, find the one with the same PF super-cluster
00155       if (isphoton) {
00156         const reco::SuperClusterRef & scRef(cand.superClusterRef());
00157         PhotonEqual myEqual(scRef);
00158         std::vector<reco::Photon>::const_iterator itcheck=find_if(photons->begin(),photons->end(),myEqual);
00159         if(itcheck==photons->end()) {
00160           std::ostringstream err;
00161           err << " Problem in PFLinker: no Photon " << std::endl;
00162           edm::LogError("PFLinker") << err.str();
00163           continue; // Watch out ! Continue
00164         } 
00165         reco::PhotonRef photonRef(photons,itcheck-photons->begin());
00166         cand.setPhotonRef(photonRef);
00167         cand.setSuperClusterRef(photonRef->superCluster());
00168         photonCandidateMap[photonRef] = candPtr;
00169       }      
00170 
00171       pfCandidates_p->push_back(cand);
00172     }      
00173     // save the PFCandidates and get a valid handle
00174   }
00175   const edm::OrphanHandle<reco::PFCandidateCollection> pfCandidateRefProd = (producePFCandidates_) ? iEvent.put(pfCandidates_p,nameOutputPF_) :
00176     edm::OrphanHandle<reco::PFCandidateCollection>();
00177   
00178   
00179   // now make the valuemaps
00180 
00181   edm::ValueMap<reco::PFCandidatePtr> pfMapGsfElectrons = fillValueMap<reco::GsfElectronCollection>(iEvent, 
00182                                                                 nameOutputElectronsPF_, 
00183                                                                 gsfElectrons, 
00184                                                                 electronCandidateMap,
00185                                                                 pfCandidateRefProd);
00186   
00187   edm::ValueMap<reco::PFCandidatePtr> pfMapPhotons = fillValueMap<reco::PhotonCollection>(iEvent, 
00188                                                       nameOutputPhotonsPF_, 
00189                                                       photons, 
00190                                                       photonCandidateMap,
00191                                                       pfCandidateRefProd);
00192   
00193 
00194   edm::ValueMap<reco::PFCandidatePtr> pfMapMuons;
00195 
00196   if(fillMuonRefs_){
00197     edm::Handle<reco::MuonCollection> muons; 
00198     iEvent.getByLabel(inputTagMuons_.label(), muons);
00199     
00200     pfMapMuons = fillValueMap<reco::MuonCollection>(iEvent, 
00201                                                     inputTagMuons_.label(), 
00202                                                     muons, 
00203                                                     muonCandidateMap,
00204                                                     pfCandidateRefProd);
00205   }
00206   
00207   std::auto_ptr<edm::ValueMap<reco::PFCandidatePtr> > 
00208     pfMapMerged(new edm::ValueMap<reco::PFCandidatePtr>());
00209   edm::ValueMap<reco::PFCandidatePtr>::Filler pfMapMergedFiller(*pfMapMerged);
00210   
00211   *pfMapMerged                   += pfMapGsfElectrons;
00212   *pfMapMerged                   += pfMapPhotons;
00213   if(fillMuonRefs_) *pfMapMerged += pfMapMuons;
00214   
00215   iEvent.put(pfMapMerged,nameOutputMergedPF_);
00216 }
00217 
00218 template<typename T>
00219 bool PFLinker::fetchCollection(edm::Handle<T>& c, 
00220                                const edm::InputTag& tag, 
00221                                const edm::Event& iEvent) const {  
00222 
00223   bool found = iEvent.getByLabel(tag, c);
00224   
00225   if(!found )
00226     {
00227       std::ostringstream  err;
00228       err<<" cannot get " <<tag<<std::endl;
00229       edm::LogError("PFLinker")<<err.str();
00230     }
00231   return found;
00232 }
00233 
00234 
00235 
00236 template<typename TYPE>
00237 edm::ValueMap<reco::PFCandidatePtr>  PFLinker::fillValueMap(edm::Event & event,
00238                                                             std::string label,
00239                                                             edm::Handle<TYPE>& inputObjCollection,
00240                                                             const std::map<edm::Ref<TYPE>, reco::PFCandidatePtr> & mapToTheCandidate,
00241                                                             const edm::OrphanHandle<reco::PFCandidateCollection> & newPFCandColl) const {
00242   
00243   std::auto_ptr<edm::ValueMap<reco::PFCandidatePtr> > pfMap_p(new edm::ValueMap<reco::PFCandidatePtr>());
00244   edm::ValueMap<reco::PFCandidatePtr>::Filler filler(*pfMap_p);
00245   
00246   typedef typename std::map<edm::Ref<TYPE>, reco::PFCandidatePtr>::const_iterator MapTYPE_it; 
00247   
00248   unsigned nObj=inputObjCollection->size();
00249   std::vector<reco::PFCandidatePtr> values(nObj);
00250 
00251   for(unsigned iobj=0; iobj < nObj; ++iobj) {
00252 
00253     edm::Ref<TYPE> objRef(inputObjCollection, iobj);
00254     MapTYPE_it  itcheck = mapToTheCandidate.find(objRef);
00255 
00256     reco::PFCandidatePtr candPtr;
00257 
00258     if(itcheck != mapToTheCandidate.end())
00259       candPtr = producePFCandidates_ ? reco::PFCandidatePtr(newPFCandColl,itcheck->second.key()) : itcheck->second;
00260     
00261     values[iobj] = candPtr;    
00262   }
00263 
00264   filler.insert(inputObjCollection,values.begin(),values.end());
00265   filler.fill();
00266   edm::ValueMap<reco::PFCandidatePtr> returnValue = *pfMap_p;
00267   event.put(pfMap_p,label);
00268   return returnValue;
00269 }