CMS 3D CMS Logo

PFRecoTauDiscriminationAgainstElectron2.cc
Go to the documentation of this file.
1 /* class PFRecoTauDiscriminationAgainstElectron2
2  * created : Apr 3, 2014
3  * revised : ,
4  * Authorss : Aruna Nayak (DESY)
5  */
6 
8 
11 
14 
20 
21 #include <TMath.h>
22 #include "TPRegexp.h"
23 
24 #include <iostream>
25 #include <sstream>
26 #include <fstream>
27 #include <vector>
28 
29 using namespace reco;
30 typedef std::pair<double, double> pdouble;
31 
33 public:
36  LeadPFChargedHadrEoP_barrel_min_ = iConfig.getParameter<double>("LeadPFChargedHadrEoP_barrel_min");
37  LeadPFChargedHadrEoP_barrel_max_ = iConfig.getParameter<double>("LeadPFChargedHadrEoP_barrel_max");
38  Hcal3x3OverPLead_barrel_max_ = iConfig.getParameter<double>("Hcal3x3OverPLead_barrel_max");
39  GammaEtaMom_barrel_max_ = iConfig.getParameter<double>("GammaEtaMom_barrel_max");
40  GammaPhiMom_barrel_max_ = iConfig.getParameter<double>("GammaPhiMom_barrel_max");
41  GammaEnFrac_barrel_max_ = iConfig.getParameter<double>("GammaEnFrac_barrel_max");
42  LeadPFChargedHadrEoP_endcap_min1_ = iConfig.getParameter<double>("LeadPFChargedHadrEoP_endcap_min1");
43  LeadPFChargedHadrEoP_endcap_max1_ = iConfig.getParameter<double>("LeadPFChargedHadrEoP_endcap_max1");
44  LeadPFChargedHadrEoP_endcap_min2_ = iConfig.getParameter<double>("LeadPFChargedHadrEoP_endcap_min2");
45  LeadPFChargedHadrEoP_endcap_max2_ = iConfig.getParameter<double>("LeadPFChargedHadrEoP_endcap_max2");
46  Hcal3x3OverPLead_endcap_max_ = iConfig.getParameter<double>("Hcal3x3OverPLead_endcap_max");
47  GammaEtaMom_endcap_max_ = iConfig.getParameter<double>("GammaEtaMom_endcap_max");
48  GammaPhiMom_endcap_max_ = iConfig.getParameter<double>("GammaPhiMom_endcap_max");
49  GammaEnFrac_endcap_max_ = iConfig.getParameter<double>("GammaEnFrac_endcap_max");
50  keepTausInEcalCrack_ = iConfig.getParameter<bool>("keepTausInEcalCrack");
51  rejectTausInEcalCrack_ = iConfig.getParameter<bool>("rejectTausInEcalCrack");
52 
53  applyCut_hcal3x3OverPLead_ = iConfig.getParameter<bool>("applyCut_hcal3x3OverPLead");
54  applyCut_leadPFChargedHadrEoP_ = iConfig.getParameter<bool>("applyCut_leadPFChargedHadrEoP");
55  applyCut_GammaEtaMom_ = iConfig.getParameter<bool>("applyCut_GammaEtaMom");
56  applyCut_GammaPhiMom_ = iConfig.getParameter<bool>("applyCut_GammaPhiMom");
57  applyCut_GammaEnFrac_ = iConfig.getParameter<bool>("applyCut_GammaEnFrac");
58  applyCut_HLTSpecific_ = iConfig.getParameter<bool>("applyCut_HLTSpecific");
59 
60  etaCracks_string_ = iConfig.getParameter<std::vector<std::string>>("etaCracks");
61 
62  verbosity_ = iConfig.getParameter<int>("verbosity");
63 
64  cuts2_ = new AntiElectronIDCut2();
65  }
66 
67  void beginEvent(const edm::Event&, const edm::EventSetup&) override;
68 
69  double discriminate(const PFTauRef&) const override;
70 
72 
73  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
74 
75 private:
76  bool isInEcalCrack(double) const;
77  std::vector<pdouble> etaCracks_;
78  std::vector<std::string> etaCracks_string_;
79 
81 
96 
99 
106 
108 };
109 
111  cuts2_->SetBarrelCutValues(LeadPFChargedHadrEoP_barrel_min_,
112  LeadPFChargedHadrEoP_barrel_max_,
113  Hcal3x3OverPLead_barrel_max_,
114  GammaEtaMom_barrel_max_,
115  GammaPhiMom_barrel_max_,
116  GammaEnFrac_barrel_max_);
117 
118  cuts2_->SetEndcapCutValues(LeadPFChargedHadrEoP_endcap_min1_,
119  LeadPFChargedHadrEoP_endcap_max1_,
120  LeadPFChargedHadrEoP_endcap_min2_,
121  LeadPFChargedHadrEoP_endcap_max2_,
122  Hcal3x3OverPLead_endcap_max_,
123  GammaEtaMom_endcap_max_,
124  GammaPhiMom_endcap_max_,
125  GammaEnFrac_endcap_max_);
126 
127  cuts2_->ApplyCut_EcalCrack(keepTausInEcalCrack_, rejectTausInEcalCrack_);
128 
129  cuts2_->ApplyCuts(applyCut_hcal3x3OverPLead_,
130  applyCut_leadPFChargedHadrEoP_,
131  applyCut_GammaEtaMom_,
132  applyCut_GammaPhiMom_,
133  applyCut_GammaEnFrac_,
134  applyCut_HLTSpecific_);
135  //Ecal cracks in eta
136  etaCracks_.clear();
137  TPRegexp regexpParser_range("([0-9.e+/-]+):([0-9.e+/-]+)");
138  for (std::vector<std::string>::const_iterator etaCrack = etaCracks_string_.begin();
139  etaCrack != etaCracks_string_.end();
140  ++etaCrack) {
141  TObjArray* subStrings = regexpParser_range.MatchS(etaCrack->data());
142  if (subStrings->GetEntries() == 3) {
143  //std::cout << "substrings(1) = " << ((TObjString*)subStrings->At(1))->GetString() << std::endl;
144  double range_begin = ((TObjString*)subStrings->At(1))->GetString().Atof();
145  //std::cout << "substrings(2) = " << ((TObjString*)subStrings->At(2))->GetString() << std::endl;
146  double range_end = ((TObjString*)subStrings->At(2))->GetString().Atof();
147  etaCracks_.push_back(pdouble(range_begin, range_end));
148  }
149  }
150 
151  cuts2_->SetEcalCracks(etaCracks_);
152 }
153 
155  double discriminator = 0.;
156 
157  // ensure tau has at least one charged object
158 
159  if ((*thePFTauRef).leadChargedHadrCand().isNull()) {
160  return 0.;
161  } else {
162  discriminator = cuts2_->Discriminator(*thePFTauRef);
163  }
164 
165  if (verbosity_) {
166  std::cout << " Taus : " << TauProducer_ << std::endl;
167  std::cout << "<PFRecoTauDiscriminationAgainstElectron2::discriminate>:" << std::endl;
168  std::cout << " tau: Pt = " << thePFTauRef->pt() << ", eta = " << thePFTauRef->eta()
169  << ", phi = " << thePFTauRef->phi() << std::endl;
170  std::cout << " discriminator value = " << discriminator << std::endl;
171  std::cout << " Prongs in tau: " << thePFTauRef->signalChargedHadrCands().size() << std::endl;
172  }
173 
174  return discriminator;
175 }
176 
178  eta = fabs(eta);
179  return (eta < 0.018 || (eta > 0.423 && eta < 0.461) || (eta > 0.770 && eta < 0.806) || (eta > 1.127 && eta < 1.163) ||
180  (eta > 1.460 && eta < 1.558));
181 }
182 
184  // pfRecoTauDiscriminationAgainstElectron2
186  desc.add<bool>("rejectTausInEcalCrack", false);
187  desc.add<edm::InputTag>("PFTauProducer", edm::InputTag("pfRecoTauProducer"));
188  desc.add<bool>("applyCut_GammaEnFrac", true);
189  desc.add<bool>("applyCut_HLTSpecific", true);
190  desc.add<double>("GammaEnFrac_barrel_max", 0.15);
191  desc.add<bool>("keepTausInEcalCrack", true);
192  {
194  psd0.add<std::string>("BooleanOperator", "and");
195  {
197  psd1.add<double>("cut");
198  psd1.add<edm::InputTag>("Producer");
199  psd0.addOptional<edm::ParameterSetDescription>("leadTrack", psd1);
200  }
201  desc.add<edm::ParameterSetDescription>("Prediscriminants", psd0);
202  }
203  desc.add<bool>("applyCut_GammaPhiMom", false);
204  desc.add<double>("GammaPhiMom_endcap_max", 1.5);
205  desc.add<double>("GammaPhiMom_barrel_max", 1.5);
206  desc.add<bool>("applyCut_leadPFChargedHadrEoP", true);
207  desc.add<double>("LeadPFChargedHadrEoP_barrel_max", 1.01);
208  desc.add<double>("GammaEtaMom_endcap_max", 1.5);
209  desc.add<double>("GammaEtaMom_barrel_max", 1.5);
210  desc.add<double>("Hcal3x3OverPLead_endcap_max", 0.1);
211  desc.add<double>("LeadPFChargedHadrEoP_barrel_min", 0.99);
212  desc.add<double>("LeadPFChargedHadrEoP_endcap_max2", 1.01);
213  desc.add<double>("LeadPFChargedHadrEoP_endcap_min1", 0.7);
214  desc.add<double>("LeadPFChargedHadrEoP_endcap_min2", 0.99);
215  desc.add<double>("LeadPFChargedHadrEoP_endcap_max1", 1.3);
216  desc.add<int>("verbosity", 0);
217  desc.add<double>("GammaEnFrac_endcap_max", 0.2);
218  desc.add<bool>("applyCut_hcal3x3OverPLead", true);
219  desc.add<bool>("applyCut_GammaEtaMom", false);
220  desc.add<std::vector<std::string>>("etaCracks",
221  {
222  "0.0:0.018",
223  "0.423:0.461",
224  "0.770:0.806",
225  "1.127:1.163",
226  "1.460:1.558",
227  });
228  desc.add<double>("Hcal3x3OverPLead_barrel_max", 0.2);
229  descriptions.add("pfRecoTauDiscriminationAgainstElectron2", desc);
230 }
231 
T getParameter(std::string const &) const
ParameterDescriptionBase * addOptional(U const &iLabel, T const &value)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
std::pair< double, double > pdouble
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void beginEvent(const edm::Event &, const edm::EventSetup &) override
ParameterDescriptionBase * add(U const &iLabel, T const &value)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
fixed size matrix
PFRecoTauDiscriminationAgainstElectron2(const edm::ParameterSet &iConfig)