CMS 3D CMS Logo

PFRecoTauDiscriminationAgainstMuon2.cc
Go to the documentation of this file.
1 
12 
14 
17 
27 
29 
30 #include <vector>
31 #include <string>
32 #include <iostream>
33 #include <atomic>
34 
35 namespace {
36 
37  class PFRecoTauDiscriminationAgainstMuon2 final : public PFTauDiscriminationProducerBase {
38  public:
39  explicit PFRecoTauDiscriminationAgainstMuon2(const edm::ParameterSet& cfg)
40  : PFTauDiscriminationProducerBase(cfg), moduleLabel_(cfg.getParameter<std::string>("@module_label")) {
41  std::string discriminatorOption_string = cfg.getParameter<std::string>("discriminatorOption");
42  int discOption;
43  if (discriminatorOption_string == "loose")
45  else if (discriminatorOption_string == "medium")
47  else if (discriminatorOption_string == "tight")
49  else if (discriminatorOption_string == "custom")
51  else
53  << " Invalid Configuration parameter 'discriminatorOption' = " << discriminatorOption_string << " !!\n";
54  wpDef_ = std::make_unique<PFRecoTauDiscriminationAgainstMuonConfigSet>(
55  discOption,
56  cfg.getParameter<double>("HoPMin"),
57  cfg.getParameter<int>("maxNumberOfMatches"),
58  cfg.getParameter<bool>("doCaloMuonVeto"),
59  cfg.getParameter<int>("maxNumberOfHitsLast2Stations"));
60  srcMuons_ = cfg.getParameter<edm::InputTag>("srcMuons");
61  Muons_token = consumes<reco::MuonCollection>(srcMuons_);
62  dRmuonMatch_ = cfg.getParameter<double>("dRmuonMatch");
63  dRmuonMatchLimitedToJetArea_ = cfg.getParameter<bool>("dRmuonMatchLimitedToJetArea");
64  minPtMatchedMuon_ = cfg.getParameter<double>("minPtMatchedMuon");
65  typedef std::vector<int> vint;
66  maskMatchesDT_ = cfg.getParameter<vint>("maskMatchesDT");
67  maskMatchesCSC_ = cfg.getParameter<vint>("maskMatchesCSC");
68  maskMatchesRPC_ = cfg.getParameter<vint>("maskMatchesRPC");
69  maskHitsDT_ = cfg.getParameter<vint>("maskHitsDT");
70  maskHitsCSC_ = cfg.getParameter<vint>("maskHitsCSC");
71  maskHitsRPC_ = cfg.getParameter<vint>("maskHitsRPC");
72  verbosity_ = cfg.getParameter<int>("verbosity");
73  }
74  ~PFRecoTauDiscriminationAgainstMuon2() override {}
75 
76  void beginEvent(const edm::Event&, const edm::EventSetup&) override;
77 
78  double discriminate(const reco::PFTauRef&) const override;
79 
80  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
81 
82  private:
84  std::unique_ptr<PFRecoTauDiscriminationAgainstMuonConfigSet> wpDef_;
85  edm::InputTag srcMuons_;
88  double dRmuonMatch_;
89  bool dRmuonMatchLimitedToJetArea_;
90  double minPtMatchedMuon_;
91  std::vector<int> maskMatchesDT_;
92  std::vector<int> maskMatchesCSC_;
93  std::vector<int> maskMatchesRPC_;
94  std::vector<int> maskHitsDT_;
95  std::vector<int> maskHitsCSC_;
96  std::vector<int> maskHitsRPC_;
97  static std::atomic<unsigned int> numWarnings_;
98  static constexpr unsigned int maxWarnings_ = 3;
99  int verbosity_;
100  };
101 
102  std::atomic<unsigned int> PFRecoTauDiscriminationAgainstMuon2::numWarnings_{0};
103 
105  if (!srcMuons_.label().empty()) {
106  evt.getByToken(Muons_token, muons_);
107  }
108  }
109 
111  const reco::PFCandidatePtr& pfCand = pfTau->leadPFChargedHadrCand();
113  moduleLabel_,
114  srcMuons_.label().empty(),
115  minPtMatchedMuon_,
116  dRmuonMatch_,
117  dRmuonMatchLimitedToJetArea_,
118  numWarnings_,
119  maxWarnings_,
120  maskMatchesDT_,
121  maskMatchesCSC_,
122  maskMatchesRPC_,
123  maskHitsDT_,
124  maskHitsCSC_,
125  maskHitsRPC_,
126  muons_,
127  pfTau,
128  pfCand);
129  double discriminatorValue = helper.eval(*wpDef_, pfTau);
130  if (verbosity_)
131  edm::LogPrint("PFTauAgainstMuon2") << "--> returning discriminatorValue = " << discriminatorValue;
132  return discriminatorValue;
133  }
134 
135 } // namespace
136 
138  // pfRecoTauDiscriminationAgainstMuon2
140  desc.add<std::vector<int>>("maskHitsRPC",
141  {
142  0,
143  0,
144  0,
145  0,
146  });
147  desc.add<int>("maxNumberOfHitsLast2Stations", 0);
148  desc.add<std::vector<int>>("maskMatchesRPC",
149  {
150  0,
151  0,
152  0,
153  0,
154  });
155  desc.add<std::vector<int>>("maskMatchesCSC",
156  {
157  1,
158  0,
159  0,
160  0,
161  });
162  desc.add<std::vector<int>>("maskHitsCSC",
163  {
164  0,
165  0,
166  0,
167  0,
168  });
169  desc.add<edm::InputTag>("PFTauProducer", edm::InputTag("pfRecoTauProducer"));
170  desc.add<int>("verbosity", 0);
171  desc.add<std::vector<int>>("maskMatchesDT",
172  {
173  0,
174  0,
175  0,
176  0,
177  });
178  desc.add<double>("minPtMatchedMuon", 5.0);
179  desc.add<bool>("dRmuonMatchLimitedToJetArea", false);
180  {
182  psd0.add<std::string>("BooleanOperator", "and");
183  {
185  psd1.add<double>("cut"); //, 0.5);
186  psd1.add<edm::InputTag>("Producer"); //, edm::InputTag("pfRecoTauDiscriminationByLeadingTrackFinding"));
187  psd0.addOptional<edm::ParameterSetDescription>("leadTrack", psd1); // optional with default?
188  }
189  // Prediscriminants can be
190  // Prediscriminants = noPrediscriminants,
191  // as in RecoTauTag/Configuration/python/HPSPFTaus_cff.py
192  //
193  // and the definition is:
194  // RecoTauTag/RecoTau/python/TauDiscriminatorTools.py
195  // # Require no prediscriminants
196  // noPrediscriminants = cms.PSet(
197  // BooleanOperator = cms.string("and"),
198  // )
199  // -- so this is the minimum required definition
200  // otherwise it inserts the leadTrack with Producer = InpuTag(...) and does not find the corresponding output during the run
201  desc.add<edm::ParameterSetDescription>("Prediscriminants", psd0);
202  }
203  desc.add<std::vector<int>>("maskHitsDT",
204  {
205  0,
206  0,
207  0,
208  0,
209  });
210  desc.add<double>("HoPMin", 0.2);
211  desc.add<int>("maxNumberOfMatches", 0);
212  desc.add<std::string>("discriminatorOption", "loose");
213  desc.add<double>("dRmuonMatch", 0.3);
214  desc.add<edm::InputTag>("srcMuons", edm::InputTag("muons"));
215  desc.add<bool>("doCaloMuonVeto", false);
216  descriptions.add("pfRecoTauDiscriminationAgainstMuon2", desc);
217 }
218 
219 DEFINE_FWK_MODULE(PFRecoTauDiscriminationAgainstMuon2);
ParameterDescriptionBase * addOptional(U const &iLabel, T const &value)
Definition: helper.py:1
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
virtual TauDiscriminatorDataType discriminate(const TauRef &tau) const =0
double vint[400]
virtual void beginEvent(const edm::Event &, const edm::EventSetup &)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
Log< level::Warning, true > LogPrint
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)