Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:25:06

0001 #include "PhysicsTools/SelectorUtils/interface/CutApplicatorWithEventContentBase.h"
0002 #include "DataFormats/EgammaCandidates/interface/GsfElectron.h"
0003 #include "DataFormats/EgammaCandidates/interface/ConversionFwd.h"
0004 #include "DataFormats/EgammaCandidates/interface/Conversion.h"
0005 #include "CommonTools/Egamma/interface/ConversionTools.h"
0006 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0007 
0008 class GsfEleConversionVetoCut : public CutApplicatorWithEventContentBase {
0009 public:
0010   GsfEleConversionVetoCut(const edm::ParameterSet& c);
0011 
0012   result_type operator()(const reco::GsfElectronPtr&) const final;
0013 
0014   void setConsumes(edm::ConsumesCollector&) final;
0015   void getEventContent(const edm::EventBase&) final;
0016 
0017   double value(const reco::CandidatePtr& cand) const final;
0018 
0019   CandidateType candidateType() const final { return ELECTRON; }
0020 
0021 private:
0022   edm::Handle<reco::ConversionCollection> _convs;
0023   edm::Handle<reco::BeamSpot> _thebs;
0024 };
0025 
0026 DEFINE_EDM_PLUGIN(CutApplicatorFactory, GsfEleConversionVetoCut, "GsfEleConversionVetoCut");
0027 
0028 GsfEleConversionVetoCut::GsfEleConversionVetoCut(const edm::ParameterSet& c) : CutApplicatorWithEventContentBase(c) {
0029   edm::InputTag conversiontag = c.getParameter<edm::InputTag>("conversionSrc");
0030   contentTags_.emplace("conversions", conversiontag);
0031 
0032   edm::InputTag conversiontagMiniAOD = c.getParameter<edm::InputTag>("conversionSrcMiniAOD");
0033   contentTags_.emplace("conversionsMiniAOD", conversiontagMiniAOD);
0034 
0035   edm::InputTag beamspottag = c.getParameter<edm::InputTag>("beamspotSrc");
0036   contentTags_.emplace("beamspot", beamspottag);
0037 }
0038 
0039 void GsfEleConversionVetoCut::setConsumes(edm::ConsumesCollector& cc) {
0040   auto convs = cc.mayConsume<reco::ConversionCollection>(contentTags_["conversions"]);
0041   auto convsMiniAOD = cc.mayConsume<reco::ConversionCollection>(contentTags_["conversionsMiniAOD"]);
0042   auto thebs = cc.consumes<reco::BeamSpot>(contentTags_["beamspot"]);
0043   contentTokens_.emplace("conversions", convs);
0044   contentTokens_.emplace("conversionsMiniAOD", convsMiniAOD);
0045   contentTokens_.emplace("beamspot", thebs);
0046 }
0047 
0048 void GsfEleConversionVetoCut::getEventContent(const edm::EventBase& ev) {
0049   // First try AOD, then go to miniAOD. Use the same Handle since collection class is the same.
0050   ev.getByLabel(contentTags_["conversions"], _convs);
0051   if (!_convs.isValid())
0052     ev.getByLabel(contentTags_["conversionsMiniAOD"], _convs);
0053 
0054   ev.getByLabel(contentTags_["beamspot"], _thebs);
0055 }
0056 
0057 CutApplicatorBase::result_type GsfEleConversionVetoCut::operator()(const reco::GsfElectronPtr& cand) const {
0058   if (_thebs.isValid() && _convs.isValid()) {
0059     return !ConversionTools::hasMatchedConversion(*cand, *_convs, _thebs->position());
0060   } else {
0061     edm::LogWarning("GsfEleConversionVetoCut") << "Couldn't find a necessary collection, returning true!";
0062   }
0063   return true;
0064 }
0065 
0066 double GsfEleConversionVetoCut::value(const reco::CandidatePtr& cand) const {
0067   reco::GsfElectronPtr ele(cand);
0068   if (_thebs.isValid() && _convs.isValid()) {
0069     return !ConversionTools::hasMatchedConversion(*ele, *_convs, _thebs->position());
0070   } else {
0071     edm::LogWarning("GsfEleConversionVetoCut") << "Couldn't find a necessary collection, returning true!";
0072     return true;
0073   }
0074 }