Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:27:49

0001 #include "RecoTauTag/RecoTau/interface/TauDiscriminationProducerBase.h"
0002 #include <FWCore/ParameterSet/interface/ConfigurationDescriptions.h>
0003 #include <FWCore/ParameterSet/interface/ParameterSetDescription.h>
0004 
0005 /* 
0006  * class PFRecoTauDiscriminationByLeadingObjectPtCut
0007  * created : October 08 2008,
0008  * revised : Wed Aug 19 17:13:04 PDT 2009
0009  * Authors : Simone Gennai (SNS), Evan Friis (UC Davis)
0010  */
0011 
0012 using namespace reco;
0013 
0014 class PFRecoTauDiscriminationByLeadingObjectPtCut : public PFTauDiscriminationProducerBase {
0015 public:
0016   explicit PFRecoTauDiscriminationByLeadingObjectPtCut(const edm::ParameterSet& iConfig)
0017       : PFTauDiscriminationProducerBase(iConfig) {
0018     chargedOnly_ = iConfig.getParameter<bool>("UseOnlyChargedHadrons");
0019     minPtLeadObject_ = iConfig.getParameter<double>("MinPtLeadingObject");
0020   }
0021   ~PFRecoTauDiscriminationByLeadingObjectPtCut() override {}
0022   double discriminate(const PFTauRef& pfTau) const override;
0023   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0024 
0025 private:
0026   bool chargedOnly_;
0027   double minPtLeadObject_;
0028 };
0029 
0030 double PFRecoTauDiscriminationByLeadingObjectPtCut::discriminate(const PFTauRef& thePFTauRef) const {
0031   double leadObjectPt = -1.;
0032   if (chargedOnly_) {
0033     // consider only charged hadrons.  note that the leadChargedHadrCand is the highest pt
0034     // charged signal cone object above the quality cut level (typically 0.5 GeV).
0035     if (thePFTauRef->leadChargedHadrCand().isNonnull()) {
0036       leadObjectPt = thePFTauRef->leadChargedHadrCand()->pt();
0037     }
0038   } else {
0039     // If using the 'leading pion' option, require that:
0040     //   1) at least one charged hadron exists above threshold (thePFTauRef->leadChargedHadrCand().isNonnull())
0041     //   2) the lead PFCand exists.  In the case that the highest pt charged hadron is above the PFRecoTauProducer threshold
0042     //      (typically 5 GeV), the leadCand and the leadChargedHadrCand are the same object.  If the leadChargedHadrCand
0043     //      is below 5GeV, but there exists a neutral PF particle > 5 GeV, it is set to be the leadCand
0044     if (thePFTauRef->leadCand().isNonnull() && thePFTauRef->leadChargedHadrCand().isNonnull()) {
0045       leadObjectPt = thePFTauRef->leadCand()->pt();
0046     }
0047   }
0048 
0049   return (leadObjectPt > minPtLeadObject_ ? 1. : 0.);
0050 }
0051 
0052 void PFRecoTauDiscriminationByLeadingObjectPtCut::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0053   // pfRecoTauDiscriminationByLeadingObjectPtCut
0054   edm::ParameterSetDescription desc;
0055   desc.add<double>("MinPtLeadingObject", 5.0);
0056   {
0057     edm::ParameterSetDescription psd0;
0058     psd0.add<std::string>("BooleanOperator", "and");
0059     desc.add<edm::ParameterSetDescription>("Prediscriminants", psd0);
0060   }
0061   desc.add<bool>("UseOnlyChargedHadrons", false);
0062   desc.add<edm::InputTag>("PFTauProducer", edm::InputTag("pfRecoTauProducer"));
0063   descriptions.add("pfRecoTauDiscriminationByLeadingObjectPtCut", desc);
0064 }
0065 
0066 DEFINE_FWK_MODULE(PFRecoTauDiscriminationByLeadingObjectPtCut);