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/Utilities/interface/InputTag.h"
0003 #include <FWCore/ParameterSet/interface/ConfigurationDescriptions.h>
0004 #include <FWCore/ParameterSet/interface/ParameterSetDescription.h>
0005 
0006 #include <memory>
0007 
0008 #include "DataFormats/VertexReco/interface/Vertex.h"
0009 #include "RecoTauTag/RecoTau/interface/RecoTauQualityCuts.h"
0010 #include "RecoTauTag/RecoTau/interface/RecoTauVertexAssociator.h"
0011 
0012 /* class PFRecoTauDiscriminationByNProngs
0013  * created : August 30 2010,
0014  * contributors : Sami Lehti (sami.lehti@cern.ch ; HIP, Helsinki)
0015  * based on H+ tau ID by Lauri Wendland
0016  * Modified April 16 2014 by S.Lehti
0017  */
0018 
0019 using namespace reco;
0020 using namespace std;
0021 using namespace edm;
0022 
0023 class PFRecoTauDiscriminationByNProngs : public PFTauDiscriminationProducerBase {
0024 public:
0025   explicit PFRecoTauDiscriminationByNProngs(const ParameterSet&);
0026   ~PFRecoTauDiscriminationByNProngs() override {}
0027 
0028   void beginEvent(const edm::Event&, const edm::EventSetup&) override;
0029   double discriminate(const reco::PFTauRef&) const override;
0030 
0031   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0032 
0033 private:
0034   std::unique_ptr<tau::RecoTauQualityCuts> qcuts_;
0035   std::unique_ptr<tau::RecoTauVertexAssociator> vertexAssociator_;
0036 
0037   uint32_t minN, maxN;
0038   bool booleanOutput;
0039   edm::ParameterSet qualityCuts;
0040 };
0041 
0042 PFRecoTauDiscriminationByNProngs::PFRecoTauDiscriminationByNProngs(const ParameterSet& iConfig)
0043     : PFTauDiscriminationProducerBase(iConfig), qualityCuts(iConfig.getParameterSet("qualityCuts")) {
0044   minN = iConfig.getParameter<uint32_t>("MinN");
0045   maxN = iConfig.getParameter<uint32_t>("MaxN");
0046   booleanOutput = iConfig.getParameter<bool>("BooleanOutput");
0047 
0048   qcuts_ = std::make_unique<tau::RecoTauQualityCuts>(qualityCuts.getParameterSet("signalQualityCuts"));
0049   vertexAssociator_ = std::make_unique<tau::RecoTauVertexAssociator>(qualityCuts, consumesCollector());
0050 }
0051 
0052 void PFRecoTauDiscriminationByNProngs::beginEvent(const Event& iEvent, const EventSetup& iSetup) {
0053   vertexAssociator_->setEvent(iEvent);
0054 }
0055 
0056 double PFRecoTauDiscriminationByNProngs::discriminate(const PFTauRef& tau) const {
0057   reco::VertexRef pv = vertexAssociator_->associatedVertex(*tau);
0058   const CandidatePtr leadingTrack = tau->leadChargedHadrCand();
0059 
0060   uint np = 0;
0061   if (leadingTrack.isNonnull() && pv.isNonnull()) {
0062     qcuts_->setPV(pv);
0063     qcuts_->setLeadTrack(*tau->leadChargedHadrCand());
0064 
0065     for (auto const& cand : tau->signalChargedHadrCands()) {
0066       if (qcuts_->filterCandRef(cand))
0067         np++;
0068     }
0069   }
0070 
0071   bool accepted = false;
0072   if (maxN == 0) {
0073     if (np == 1 || np == 3)
0074       accepted = true;
0075   } else {
0076     if (np >= minN && np <= maxN)
0077       accepted = true;
0078   }
0079 
0080   if (!accepted)
0081     np = 0;
0082   if (booleanOutput)
0083     return accepted;
0084   return np;
0085 }
0086 
0087 void PFRecoTauDiscriminationByNProngs::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0088   // pfRecoTauDiscriminationByNProngs
0089   edm::ParameterSetDescription desc;
0090 
0091   edm::ParameterSetDescription desc_qualityCuts;
0092   reco::tau::RecoTauQualityCuts::fillDescriptions(desc_qualityCuts);
0093   desc.add<edm::ParameterSetDescription>("qualityCuts", desc_qualityCuts);
0094 
0095   {
0096     edm::ParameterSetDescription psd0;
0097     psd0.add<std::string>("BooleanOperator", "and");
0098     desc.add<edm::ParameterSetDescription>("Prediscriminants", psd0);
0099   }
0100 
0101   desc.add<bool>("BooleanOutput", true);
0102   desc.add<edm::InputTag>("PFTauProducer", edm::InputTag("combinatoricRecoTaus"));
0103   desc.add<unsigned int>("MinN", 1);
0104   desc.add<unsigned int>("MaxN", 0);
0105   descriptions.add("pfRecoTauDiscriminationByNProngs", desc);
0106 }
0107 
0108 DEFINE_FWK_MODULE(PFRecoTauDiscriminationByNProngs);