CMS 3D CMS Logo

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