CMS 3D CMS Logo

PFRecoTauDiscriminationByNProngs.cc
Go to the documentation of this file.
8 
9 /* class PFRecoTauDiscriminationByNProngs
10  * created : August 30 2010,
11  * contributors : Sami Lehti (sami.lehti@cern.ch ; HIP, Helsinki)
12  * based on H+ tau ID by Lauri Wendland
13  * Modified April 16 2014 by S.Lehti
14  */
15 
16 using namespace reco;
17 using namespace std;
18 using namespace edm;
19 
21 public:
24 
25  void beginEvent(const edm::Event&, const edm::EventSetup&) override;
26  double discriminate(const reco::PFTauRef&) const override;
27 
28  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
29 
30 private:
31  std::unique_ptr<tau::RecoTauQualityCuts> qcuts_;
32  std::unique_ptr<tau::RecoTauVertexAssociator> vertexAssociator_;
33 
34  uint32_t minN, maxN;
37 };
38 
40  : PFTauDiscriminationProducerBase(iConfig), qualityCuts(iConfig.getParameterSet("qualityCuts")) {
41  minN = iConfig.getParameter<uint32_t>("MinN");
42  maxN = iConfig.getParameter<uint32_t>("MaxN");
43  booleanOutput = iConfig.getParameter<bool>("BooleanOutput");
44 
45  qcuts_.reset(new tau::RecoTauQualityCuts(qualityCuts.getParameterSet("signalQualityCuts")));
46  vertexAssociator_.reset(new tau::RecoTauVertexAssociator(qualityCuts, consumesCollector()));
47 }
48 
50  vertexAssociator_->setEvent(iEvent);
51 }
52 
54  reco::VertexRef pv = vertexAssociator_->associatedVertex(*tau);
55  const CandidatePtr leadingTrack = tau->leadChargedHadrCand();
56 
57  uint np = 0;
58  if (leadingTrack.isNonnull() && pv.isNonnull()) {
59  qcuts_->setPV(pv);
60  qcuts_->setLeadTrack(*tau->leadChargedHadrCand());
61 
62  for (auto const& cand : tau->signalChargedHadrCands()) {
63  if (qcuts_->filterCandRef(cand))
64  np++;
65  }
66  }
67 
68  bool accepted = false;
69  if (maxN == 0) {
70  if (np == 1 || np == 3)
71  accepted = true;
72  } else {
73  if (np >= minN && np <= maxN)
74  accepted = true;
75  }
76 
77  if (!accepted)
78  np = 0;
79  if (booleanOutput)
80  return accepted;
81  return np;
82 }
83 
85  // pfRecoTauDiscriminationByNProngs
87 
88  {
89  edm::ParameterSetDescription pset_signalQualityCuts;
90  pset_signalQualityCuts.add<double>("maxDeltaZ", 0.4);
91  pset_signalQualityCuts.add<double>("minTrackPt", 0.5);
92  pset_signalQualityCuts.add<double>("minTrackVertexWeight", -1.0);
93  pset_signalQualityCuts.add<double>("maxTrackChi2", 100.0);
94  pset_signalQualityCuts.add<unsigned int>("minTrackPixelHits", 0);
95  pset_signalQualityCuts.add<double>("minGammaEt", 1.0);
96  pset_signalQualityCuts.add<unsigned int>("minTrackHits", 3);
97  pset_signalQualityCuts.add<double>("minNeutralHadronEt", 30.0);
98  pset_signalQualityCuts.add<double>("maxTransverseImpactParameter", 0.1);
99  pset_signalQualityCuts.addOptional<bool>("useTracksInsteadOfPFHadrons");
100 
101  edm::ParameterSetDescription pset_vxAssocQualityCuts;
102  pset_vxAssocQualityCuts.add<double>("minTrackPt", 0.5);
103  pset_vxAssocQualityCuts.add<double>("minTrackVertexWeight", -1.0);
104  pset_vxAssocQualityCuts.add<double>("maxTrackChi2", 100.0);
105  pset_vxAssocQualityCuts.add<unsigned int>("minTrackPixelHits", 0);
106  pset_vxAssocQualityCuts.add<double>("minGammaEt", 1.0);
107  pset_vxAssocQualityCuts.add<unsigned int>("minTrackHits", 3);
108  pset_vxAssocQualityCuts.add<double>("maxTransverseImpactParameter", 0.1);
109  pset_vxAssocQualityCuts.addOptional<bool>("useTracksInsteadOfPFHadrons");
110 
111  edm::ParameterSetDescription pset_isolationQualityCuts;
112  pset_isolationQualityCuts.add<double>("maxDeltaZ", 0.2);
113  pset_isolationQualityCuts.add<double>("minTrackPt", 1.0);
114  pset_isolationQualityCuts.add<double>("minTrackVertexWeight", -1.0);
115  pset_isolationQualityCuts.add<double>("maxTrackChi2", 100.0);
116  pset_isolationQualityCuts.add<unsigned int>("minTrackPixelHits", 0);
117  pset_isolationQualityCuts.add<double>("minGammaEt", 1.5);
118  pset_isolationQualityCuts.add<unsigned int>("minTrackHits", 8);
119  pset_isolationQualityCuts.add<double>("maxTransverseImpactParameter", 0.03);
120  pset_isolationQualityCuts.addOptional<bool>("useTracksInsteadOfPFHadrons");
121 
122  edm::ParameterSetDescription pset_qualityCuts;
123  pset_qualityCuts.add<edm::ParameterSetDescription>("signalQualityCuts", pset_signalQualityCuts);
124  pset_qualityCuts.add<edm::ParameterSetDescription>("vxAssocQualityCuts", pset_vxAssocQualityCuts);
125  pset_qualityCuts.add<edm::ParameterSetDescription>("isolationQualityCuts", pset_isolationQualityCuts);
126  pset_qualityCuts.add<std::string>("leadingTrkOrPFCandOption", "leadPFCand");
127  pset_qualityCuts.add<std::string>("pvFindingAlgo", "closestInDeltaZ");
128  pset_qualityCuts.add<edm::InputTag>("primaryVertexSrc", edm::InputTag("offlinePrimaryVertices"));
129  pset_qualityCuts.add<bool>("vertexTrackFiltering", false);
130  pset_qualityCuts.add<bool>("recoverLeadingTrk", false);
131 
132  desc.add<edm::ParameterSetDescription>("qualityCuts", pset_qualityCuts);
133  }
134 
135  {
137  psd0.add<std::string>("BooleanOperator", "and");
138  desc.add<edm::ParameterSetDescription>("Prediscriminants", psd0);
139  }
140 
141  desc.add<bool>("BooleanOutput", true);
142  desc.add<edm::InputTag>("PFTauProducer", edm::InputTag("combinatoricRecoTaus"));
143  desc.add<unsigned int>("MinN", 1);
144  desc.add<unsigned int>("MaxN", 0);
145  descriptions.add("pfRecoTauDiscriminationByNProngs", desc);
146 }
147 
ConfigurationDescriptions.h
PFRecoTauDiscriminationByNProngs::discriminate
double discriminate(const reco::PFTauRef &) const override
Definition: PFRecoTauDiscriminationByNProngs.cc:53
edm::ParameterSetDescription::add
ParameterDescriptionBase * add(U const &iLabel, T const &value)
Definition: ParameterSetDescription.h:95
TauDiscriminationProducerBase.h
reco::tau::RecoTauQualityCuts
Definition: RecoTauQualityCuts.h:34
metsig::tau
Definition: SignAlgoResolutions.h:49
edm
HLT enums.
Definition: AlignableModifier.h:19
PFRecoTauDiscriminationByNProngs
Definition: PFRecoTauDiscriminationByNProngs.cc:20
np
int np
Definition: AMPTWrapper.h:43
edm::ParameterSetDescription
Definition: ParameterSetDescription.h:52
PFRecoTauDiscriminationByNProngs::qualityCuts
edm::ParameterSet qualityCuts
Definition: PFRecoTauDiscriminationByNProngs.cc:36
reco
fixed size matrix
Definition: AlignmentAlgorithmBase.h:45
edm::ParameterSetDescription::addOptional
ParameterDescriptionBase * addOptional(U const &iLabel, T const &value)
Definition: ParameterSetDescription.h:105
parallelization.uint
uint
Definition: parallelization.py:124
RecoTauVertexAssociator.h
edm::Ref< PFTauCollection >
PFRecoTauDiscriminationByNProngs::maxN
uint32_t maxN
Definition: PFRecoTauDiscriminationByNProngs.cc:34
reco::tau::RecoTauVertexAssociator
Definition: RecoTauVertexAssociator.h:50
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::ConfigurationDescriptions::add
void add(std::string const &label, ParameterSetDescription const &psetDescription)
Definition: ConfigurationDescriptions.cc:57
fillDescriptions
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
cms::dd::accepted
bool accepted(std::vector< std::regex > const &, std::string_view)
numberPerLSFilter_cff.maxN
maxN
Definition: numberPerLSFilter_cff.py:5
TauDiscriminationProducerBase
Definition: TauDiscriminationProducerBase.h:55
PFRecoTauDiscriminationByNProngs::~PFRecoTauDiscriminationByNProngs
~PFRecoTauDiscriminationByNProngs() override
Definition: PFRecoTauDiscriminationByNProngs.cc:23
ParameterSetDescription.h
edm::ConfigurationDescriptions
Definition: ConfigurationDescriptions.h:28
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
Vertex.h
RecoTauQualityCuts.h
HLT_2018_cff.InputTag
InputTag
Definition: HLT_2018_cff.py:79016
edm::ParameterSet
Definition: ParameterSet.h:36
PFRecoTauDiscriminationByNProngs::qcuts_
std::unique_ptr< tau::RecoTauQualityCuts > qcuts_
Definition: PFRecoTauDiscriminationByNProngs.cc:31
PFRecoTauDiscriminationByNProngs::booleanOutput
bool booleanOutput
Definition: PFRecoTauDiscriminationByNProngs.cc:35
cand
Definition: decayParser.h:34
MetAnalyzer.pv
def pv(vc)
Definition: MetAnalyzer.py:7
iEvent
int iEvent
Definition: GenABIO.cc:224
edm::EventSetup
Definition: EventSetup.h:57
PFRecoTauDiscriminationByNProngs::vertexAssociator_
std::unique_ptr< tau::RecoTauVertexAssociator > vertexAssociator_
Definition: PFRecoTauDiscriminationByNProngs.cc:32
PFRecoTauDiscriminationByNProngs::PFRecoTauDiscriminationByNProngs
PFRecoTauDiscriminationByNProngs(const ParameterSet &)
Definition: PFRecoTauDiscriminationByNProngs.cc:39
InputTag.h
edm::Ptr< Candidate >
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
PFRecoTauDiscriminationByNProngs::minN
uint32_t minN
Definition: PFRecoTauDiscriminationByNProngs.cc:34
edm::getParameterSet
ParameterSet const & getParameterSet(ParameterSetID const &id)
Definition: ParameterSet.cc:855
std
Definition: JetResolutionObject.h:76
PFRecoTauDiscriminationByNProngs::beginEvent
void beginEvent(const edm::Event &, const edm::EventSetup &) override
Definition: PFRecoTauDiscriminationByNProngs.cc:49
edm::Ptr::isNonnull
bool isNonnull() const
Checks for non-null.
Definition: Ptr.h:146
edm::Event
Definition: Event.h:73
edm::InputTag
Definition: InputTag.h:15
PFRecoTauDiscriminationByNProngs::fillDescriptions
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: PFRecoTauDiscriminationByNProngs.cc:84
edm::ParameterSet::getParameterSet
ParameterSet const & getParameterSet(std::string const &) const
Definition: ParameterSet.cc:2121
beam_dqm_sourceclient-live_cfg.qualityCuts
qualityCuts
Definition: beam_dqm_sourceclient-live_cfg.py:115