CMS 3D CMS Logo

Public Member Functions | Private Member Functions | Private Attributes

PFLinker Class Reference

#include <PFLinker.h>

Inheritance diagram for PFLinker:
edm::EDProducer edm::ProducerBase edm::ProductRegistryHelper

List of all members.

Public Member Functions

virtual void beginRun (edm::Run &run, const edm::EventSetup &es)
 PFLinker (const edm::ParameterSet &)
virtual void produce (edm::Event &, const edm::EventSetup &)
 ~PFLinker ()

Private Member Functions

template<typename T >
bool fetchCollection (edm::Handle< T > &c, const edm::InputTag &tag, const edm::Event &iEvent) const
template<typename TYPE >
edm::ValueMap
< reco::PFCandidatePtr
fillValueMap (edm::Event &event, std::string label, edm::Handle< TYPE > &inputObjCollection, const std::map< edm::Ref< TYPE >, reco::PFCandidatePtr > &mapToTheCandidate, const edm::OrphanHandle< reco::PFCandidateCollection > &newPFCandColl) const

Private Attributes

bool fillMuonRefs_
 Set muon refs and produce the value map?
edm::InputTag inputTagGsfElectrons_
 Input GsfElectrons.
edm::InputTag inputTagMuons_
 Input Muons.
std::vector< edm::InputTaginputTagPFCandidates_
 Input PFCandidates.
edm::InputTag inputTagPhotons_
 Input Photons.
std::string nameOutputElectronsPF_
 name of output ValueMap electrons
std::string nameOutputMergedPF_
 name of output merged ValueMap
std::string nameOutputPF_
 name of output collection of PFCandidate
std::string nameOutputPhotonsPF_
 name of output ValueMap photons
bool producePFCandidates_
 Flags - if true: References will be towards new collection ; if false to the original one.

Detailed Description

Producer meant for the Post PF reconstruction.

Fills the GsfElectron, Photon and Muon Ref into the PFCandidate Produces the ValueMap between GsfElectronRef/Photon/Mupns with PFCandidateRef

Date:
2011/07/20 09:57:25
Revision:
1.6
Author:
R. Bellan - UCSB <riccardo.bellan@cern.ch>, F. Beaudette - CERN <Florian.Beaudette@cern.ch>

Definition at line 31 of file PFLinker.h.


Constructor & Destructor Documentation

PFLinker::PFLinker ( const edm::ParameterSet iConfig) [explicit]

Definition at line 14 of file PFLinker.cc.

References fillMuonRefs_, edm::ParameterSet::getParameter(), inputTagGsfElectrons_, inputTagMuons_, inputTagPFCandidates_, inputTagPhotons_, edm::InputTag::label(), nameOutputElectronsPF_, nameOutputMergedPF_, nameOutputPF_, nameOutputPhotonsPF_, and producePFCandidates_.

                                                  {
  // vector of InputTag; more than 1 is not for RECO, it is for analysis
  inputTagPFCandidates_ 
    = iConfig.getParameter<std::vector<edm::InputTag> >("PFCandidate");
  inputTagGsfElectrons_
    = iConfig.getParameter<edm::InputTag>("GsfElectrons");
  inputTagPhotons_
    = iConfig.getParameter<edm::InputTag>("Photons");
  inputTagMuons_
    = iConfig.getParameter<edm::InputTag>("Muons");
  
  nameOutputPF_ 
    = iConfig.getParameter<std::string>("OutputPF");
  
  nameOutputElectronsPF_ 
    = iConfig.getParameter<std::string>("ValueMapElectrons");

  nameOutputPhotonsPF_ 
    = iConfig.getParameter<std::string>("ValueMapPhotons");

  producePFCandidates_  
    = iConfig.getParameter<bool>("ProducePFCandidates");

  nameOutputMergedPF_ 
    = iConfig.getParameter<std::string>("ValueMapMerged"); 
  
  fillMuonRefs_
    = iConfig.getParameter<bool>("FillMuonRefs");

  // should not produce PFCandidates and read seve
  if(producePFCandidates_ && inputTagPFCandidates_.size()>1) {
    edm::LogError("PFLinker") << " cannot read several collections of PFCandidates and produce a new collection at the same time. " << std::endl;
    assert(0);
  }
  if(producePFCandidates_) {
    produces<reco::PFCandidateCollection>(nameOutputPF_);
  }
  produces<edm::ValueMap<reco::PFCandidatePtr> > (nameOutputElectronsPF_);
  produces<edm::ValueMap<reco::PFCandidatePtr> > (nameOutputPhotonsPF_);
  produces<edm::ValueMap<reco::PFCandidatePtr> > (nameOutputMergedPF_);
  if(fillMuonRefs_)  produces<edm::ValueMap<reco::PFCandidatePtr> > (inputTagMuons_.label());

}
PFLinker::~PFLinker ( )

Definition at line 58 of file PFLinker.cc.

{;}

Member Function Documentation

void PFLinker::beginRun ( edm::Run run,
const edm::EventSetup es 
) [virtual]

Reimplemented from edm::EDProducer.

Definition at line 60 of file PFLinker.cc.

{;}
template<typename T >
bool PFLinker::fetchCollection ( edm::Handle< T > &  c,
const edm::InputTag tag,
const edm::Event iEvent 
) const [private]

Definition at line 200 of file PFLinker.cc.

References newFWLiteAna::found, and edm::Event::getByLabel().

                                                             {  

  bool found = iEvent.getByLabel(tag, c);
  
  if(!found )
    {
      std::ostringstream  err;
      err<<" cannot get " <<tag<<std::endl;
      edm::LogError("PFLinker")<<err.str();
    }
  return found;
}
template<typename TYPE >
edm::ValueMap< reco::PFCandidatePtr > PFLinker::fillValueMap ( edm::Event event,
std::string  label,
edm::Handle< TYPE > &  inputObjCollection,
const std::map< edm::Ref< TYPE >, reco::PFCandidatePtr > &  mapToTheCandidate,
const edm::OrphanHandle< reco::PFCandidateCollection > &  newPFCandColl 
) const [private]

Definition at line 218 of file PFLinker.cc.

References edm::helper::Filler< Map >::fill(), edm::helper::Filler< Map >::insert(), producePFCandidates_, and makeHLTPrescaleTable::values.

                                                                                                                                    {
  
  std::auto_ptr<edm::ValueMap<reco::PFCandidatePtr> > pfMap_p(new edm::ValueMap<reco::PFCandidatePtr>());
  edm::ValueMap<reco::PFCandidatePtr>::Filler filler(*pfMap_p);
  
  typedef typename std::map<edm::Ref<TYPE>, reco::PFCandidatePtr>::const_iterator MapTYPE_it; 
  
  unsigned nObj=inputObjCollection->size();
  std::vector<reco::PFCandidatePtr> values(nObj);

  for(unsigned iobj=0; iobj < nObj; ++iobj) {

    edm::Ref<TYPE> objRef(inputObjCollection, iobj);
    MapTYPE_it  itcheck = mapToTheCandidate.find(objRef);

    reco::PFCandidatePtr candPtr;

    if(itcheck != mapToTheCandidate.end())
      candPtr = producePFCandidates_ ? reco::PFCandidatePtr(newPFCandColl,itcheck->second.key()) : itcheck->second;
    
    values[iobj] = candPtr;    
  }

  filler.insert(inputObjCollection,values.begin(),values.end());
  filler.fill();
  edm::ValueMap<reco::PFCandidatePtr> returnValue = *pfMap_p;
  event.put(pfMap_p,label);
  return returnValue;
}
void PFLinker::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
) [virtual]

Implements edm::EDProducer.

Definition at line 62 of file PFLinker.cc.

References reco::PFCandidate::e, fillMuonRefs_, reco::PFCandidate::gamma, gsfElectrons_cfi::gsfElectrons, reco::PFCandidate::gsfTrackRef(), i, iEvent, inputTagGsfElectrons_, inputTagMuons_, inputTagPFCandidates_, inputTagPhotons_, edm::Ref< C, T, F >::isNonnull(), reco::PFCandidate::muonRef(), patZpeak::muons, reco::PFCandidate::mva_nothing_gamma(), nameOutputElectronsPF_, nameOutputMergedPF_, nameOutputPF_, nameOutputPhotonsPF_, reco::PFCandidate::particleId(), reco::tau::pfCandidates(), interactiveExample::photons, producePFCandidates_, edm::Event::put(), reco::PFCandidate::setGsfElectronRef(), reco::PFCandidate::setMuonRef(), reco::PFCandidate::setPhotonRef(), reco::PFCandidate::setSuperClusterRef(), ntuplemaker::status, and reco::PFCandidate::superClusterRef().

                                                                    {
  
  std::auto_ptr<reco::PFCandidateCollection>
    pfCandidates_p(new reco::PFCandidateCollection);
  
  edm::Handle<reco::GsfElectronCollection> gsfElectrons;
  bool status=fetchCollection<reco::GsfElectronCollection>(gsfElectrons,
                                                      inputTagGsfElectrons_,
                                                      iEvent );
  std::map<reco::GsfElectronRef,reco::PFCandidatePtr> electronCandidateMap;  


  edm::Handle<reco::PhotonCollection> photons;
  status=fetchCollection<reco::PhotonCollection>(photons,
                                                 inputTagPhotons_,
                                                 iEvent );
  std::map<reco::PhotonRef,reco::PFCandidatePtr> photonCandidateMap;


  edm::Handle<reco::MuonToMuonMap> muonMap;
  if(fillMuonRefs_)
    status=fetchCollection<reco::MuonToMuonMap>(muonMap,
                                                inputTagMuons_,
                                                iEvent);
  std::map<reco::MuonRef,reco::PFCandidatePtr> muonCandidateMap;
  
  unsigned nColPF=inputTagPFCandidates_.size();
  edm::Handle<reco::PFCandidateCollection> pfCandidates;
  for(unsigned icol=0;icol<nColPF;++icol) {
    
    bool status=fetchCollection<reco::PFCandidateCollection>(pfCandidates, 
                                                             inputTagPFCandidates_[icol], 
                                                             iEvent );
    
    unsigned ncand=(status)?pfCandidates->size():0;
    
    for( unsigned i=0; i<ncand; ++i) {
      edm::Ptr<reco::PFCandidate> candPtr(pfCandidates,i);
      reco::PFCandidate cand(candPtr);
      
      bool isphoton   = cand.particleId() == reco::PFCandidate::gamma && cand.mva_nothing_gamma()>0.;
      bool iselectron = cand.particleId() == reco::PFCandidate::e;
      // PFCandidates may have a valid MuonRef though they are not muons.
      bool hasNonNullMuonRef  = cand.muonRef().isNonnull() && fillMuonRefs_;
      
      // if not an electron or a photon or a muon just fill the PFCandidate collection
      if ( !(isphoton || iselectron || hasNonNullMuonRef)){pfCandidates_p->push_back(cand); continue;}
      
      
      if (hasNonNullMuonRef) {
        reco::MuonRef muRef = (*muonMap)[cand.muonRef()];
        cand.setMuonRef(muRef);
        muonCandidateMap[muRef] = candPtr;
      }
                 
      
      // if it is an electron. Find the GsfElectron with the same GsfTrack
      if (iselectron) {
        const reco::GsfTrackRef & gsfTrackRef(cand.gsfTrackRef());
        GsfElectronEqual myEqual(gsfTrackRef);
        std::vector<reco::GsfElectron>::const_iterator itcheck=find_if(gsfElectrons->begin(),gsfElectrons->end(),myEqual);
        if(itcheck==gsfElectrons->end()) {
          std::ostringstream err;
          err << " Problem in PFLinker: no GsfElectron " << std::endl;
          edm::LogError("PFLinker") << err.str();
          continue; // Watch out ! Continue
        } 
        reco::GsfElectronRef electronRef(gsfElectrons,itcheck-gsfElectrons->begin());
        cand.setGsfElectronRef(electronRef);
        cand.setSuperClusterRef(electronRef->superCluster());
        electronCandidateMap[electronRef] = candPtr;
      }  
      
      // if it is a photon, find the one with the same PF super-cluster
      if (isphoton) {
        const reco::SuperClusterRef & scRef(cand.superClusterRef());
        PhotonEqual myEqual(scRef);
        std::vector<reco::Photon>::const_iterator itcheck=find_if(photons->begin(),photons->end(),myEqual);
        if(itcheck==photons->end()) {
          std::ostringstream err;
          err << " Problem in PFLinker: no Photon " << std::endl;
          edm::LogError("PFLinker") << err.str();
          continue; // Watch out ! Continue
        } 
        reco::PhotonRef photonRef(photons,itcheck-photons->begin());
        cand.setPhotonRef(photonRef);
        cand.setSuperClusterRef(photonRef->superCluster());
        photonCandidateMap[photonRef] = candPtr;
      }      

      pfCandidates_p->push_back(cand);
    }      
    // save the PFCandidates and get a valid handle
  }
  const edm::OrphanHandle<reco::PFCandidateCollection> pfCandidateRefProd = (producePFCandidates_) ? iEvent.put(pfCandidates_p,nameOutputPF_) :
    edm::OrphanHandle<reco::PFCandidateCollection>();
  
  
  // now make the valuemaps

  edm::ValueMap<reco::PFCandidatePtr> pfMapGsfElectrons = fillValueMap<reco::GsfElectronCollection>(iEvent, 
                                                                nameOutputElectronsPF_, 
                                                                gsfElectrons, 
                                                                electronCandidateMap,
                                                                pfCandidateRefProd);
  
  edm::ValueMap<reco::PFCandidatePtr> pfMapPhotons = fillValueMap<reco::PhotonCollection>(iEvent, 
                                                      nameOutputPhotonsPF_, 
                                                      photons, 
                                                      photonCandidateMap,
                                                      pfCandidateRefProd);
  

  edm::ValueMap<reco::PFCandidatePtr> pfMapMuons;

  if(fillMuonRefs_){
    edm::Handle<reco::MuonCollection> muons; 
    iEvent.getByLabel(inputTagMuons_.label(), muons);
    
    pfMapMuons = fillValueMap<reco::MuonCollection>(iEvent, 
                                                    inputTagMuons_.label(), 
                                                    muons, 
                                                    muonCandidateMap,
                                                    pfCandidateRefProd);
  }
  
  std::auto_ptr<edm::ValueMap<reco::PFCandidatePtr> > 
    pfMapMerged(new edm::ValueMap<reco::PFCandidatePtr>());
  edm::ValueMap<reco::PFCandidatePtr>::Filler pfMapMergedFiller(*pfMapMerged);
  
  *pfMapMerged                   += pfMapGsfElectrons;
  *pfMapMerged                   += pfMapPhotons;
  if(fillMuonRefs_) *pfMapMerged += pfMapMuons;
  
  iEvent.put(pfMapMerged,nameOutputMergedPF_);
}

Member Data Documentation

bool PFLinker::fillMuonRefs_ [private]

Set muon refs and produce the value map?

Definition at line 86 of file PFLinker.h.

Referenced by PFLinker(), and produce().

Input GsfElectrons.

Definition at line 62 of file PFLinker.h.

Referenced by PFLinker(), and produce().

Input Muons.

Definition at line 68 of file PFLinker.h.

Referenced by PFLinker(), and produce().

Input PFCandidates.

Definition at line 59 of file PFLinker.h.

Referenced by PFLinker(), and produce().

Input Photons.

Definition at line 65 of file PFLinker.h.

Referenced by PFLinker(), and produce().

std::string PFLinker::nameOutputElectronsPF_ [private]

name of output ValueMap electrons

Definition at line 74 of file PFLinker.h.

Referenced by PFLinker(), and produce().

std::string PFLinker::nameOutputMergedPF_ [private]

name of output merged ValueMap

Definition at line 80 of file PFLinker.h.

Referenced by PFLinker(), and produce().

std::string PFLinker::nameOutputPF_ [private]

name of output collection of PFCandidate

Definition at line 71 of file PFLinker.h.

Referenced by PFLinker(), and produce().

std::string PFLinker::nameOutputPhotonsPF_ [private]

name of output ValueMap photons

Definition at line 77 of file PFLinker.h.

Referenced by PFLinker(), and produce().

Flags - if true: References will be towards new collection ; if false to the original one.

Definition at line 83 of file PFLinker.h.

Referenced by fillValueMap(), PFLinker(), and produce().