CMS 3D CMS Logo

PFRecoTauDiscriminationByHPSSelection.cc
Go to the documentation of this file.
3 
6 
9 
10 namespace {
11  // Apply a hypothesis on the mass of the strips.
12  math::XYZTLorentzVector applyMassConstraint(const math::XYZTLorentzVector& vec, double mass)
13  {
14  double factor = sqrt(vec.energy()*vec.energy() - mass*mass)/vec.P();
15  return math::XYZTLorentzVector(vec.px()*factor, vec.py()*factor, vec.pz()*factor, vec.energy());
16  }
17 }
18 
20 {
21  public:
24  double discriminate(const reco::PFTauRef&) const override;
25 
26  private:
28 
29  struct DecayModeCuts
30  {
32  : maxMass_(nullptr)
33  {}
34  ~DecayModeCuts() {} // CV: maxMass object gets deleted by PFRecoTauDiscriminationByHPSSelection destructor
35  unsigned nTracksMin_;
37  double minMass_;
38  TauFunc* maxMass_;
42  double minPi0Mass_;
43  double maxPi0Mass_;
45  };
46 
47  typedef std::pair<unsigned int, unsigned int> IntPair;
48  typedef std::pair<double, double> DoublePair;
49  typedef std::map<IntPair, DecayModeCuts> DecayModeCutMap;
50 
51  DecayModeCutMap decayModeCuts_;
52  double matchingCone_;
53  double minPt_;
54 
56 
58 
60 };
61 
64 {
65  // Get the matchign cut
66  matchingCone_ = pset.getParameter<double>("matchingCone");
67  minPt_ = pset.getParameter<double>("minTauPt");
68  // Get the mass cuts for each decay mode
69  typedef std::vector<edm::ParameterSet> VPSet;
70  const VPSet& decayModes = pset.getParameter<VPSet>("decayModes");
71  for(auto const& decayMode : decayModes ) {
72  // The mass window(s)
74  if ( decayMode.exists("nTracksMin") ) {
75  cuts.nTracksMin_ = decayMode.getParameter<unsigned>("nTracksMin");
76  } else {
77  cuts.nTracksMin_ = 0;
78  }
79  if ( decayMode.exists("nChargedPFCandsMin") ) {
80  cuts.nChargedPFCandsMin_ = decayMode.getParameter<unsigned>("nChargedPFCandsMin");
81  } else {
82  cuts.nChargedPFCandsMin_ = 0;
83  }
84  cuts.minMass_ = decayMode.getParameter<double>("minMass");
85  cuts.maxMass_ = new TauFunc(decayMode.getParameter<std::string>("maxMass"));
86  edm::ParameterSet applyBendCorrection = decayMode.getParameter<edm::ParameterSet>("applyBendCorrection");
87  cuts.applyBendCorrection_eta_ = applyBendCorrection.getParameter<bool>("eta");
88  cuts.applyBendCorrection_phi_ = applyBendCorrection.getParameter<bool>("phi");
89  cuts.applyBendCorrection_mass_ = applyBendCorrection.getParameter<bool>("mass");
90  if ( decayMode.exists("minPi0Mass") ) {
91  cuts.minPi0Mass_ = decayMode.getParameter<double>("minPi0Mass");
92  cuts.maxPi0Mass_ = decayMode.getParameter<double>("maxPi0Mass");
93  } else {
94  cuts.minPi0Mass_ = -1.e3;
95  cuts.maxPi0Mass_ = 1.e9;
96  }
97  if ( decayMode.exists("assumeStripMass") ) {
98  cuts.assumeStripMass_ = decayMode.getParameter<double>("assumeStripMass");
99  } else {
100  cuts.assumeStripMass_ = -1.0;
101  }
102  decayModeCuts_.insert(std::make_pair(
103  // The decay mode as a key
104  std::make_pair(
105  decayMode.getParameter<uint32_t>("nCharged"),
106  decayMode.getParameter<uint32_t>("nPiZeros")),
107  cuts
108  ));
109  }
110  requireTauChargedHadronsToBeChargedPFCands_ = pset.getParameter<bool>("requireTauChargedHadronsToBeChargedPFCands");
111  minPixelHits_ = pset.getParameter<int>("minPixelHits");
112  verbosity_ = pset.exists("verbosity") ?
113  pset.getParameter<int>("verbosity") : 0;
114 
115 
116 }
117 
119 {
120  for ( DecayModeCutMap::iterator it = decayModeCuts_.begin();
121  it != decayModeCuts_.end(); ++it ) {
122  delete it->second.maxMass_;
123  }
124 }
125 
126 double
128 {
129  if ( verbosity_ ) {
130  edm::LogPrint("PFTauByHPSSelect") << "<PFRecoTauDiscriminationByHPSSelection::discriminate>:" ;
131  edm::LogPrint("PFTauByHPSSelect") << " nCharged = " << tau->signalTauChargedHadronCandidates().size() ;
132  edm::LogPrint("PFTauByHPSSelect") << " nPiZeros = " << tau->signalPiZeroCandidates().size() ;
133  }
134 
135  // Check if we pass the min pt
136  if ( tau->pt() < minPt_ ) {
137  if ( verbosity_ ) {
138  edm::LogPrint("PFTauByHPSSelect") << " fails minPt cut." ;
139  }
140  return 0.0;
141  }
142 
143  // See if we select this decay mode
144  DecayModeCutMap::const_iterator massWindowIter =
145  decayModeCuts_.find(std::make_pair(tau->signalTauChargedHadronCandidates().size(),
146  tau->signalPiZeroCandidates().size()));
147  // Check if decay mode is supported
148  if ( massWindowIter == decayModeCuts_.end() ) {
149  if ( verbosity_ ) {
150  edm::LogPrint("PFTauByHPSSelect") << " fails mass-window definition requirement." ;
151  }
152  return 0.0;
153  }
154  const DecayModeCuts& massWindow = massWindowIter->second;
155 
156  if ( massWindow.nTracksMin_ > 0 ) {
157  unsigned int nTracks = 0;
158  const std::vector<reco::PFRecoTauChargedHadron>& chargedHadrons = tau->signalTauChargedHadronCandidates();
159  for ( std::vector<reco::PFRecoTauChargedHadron>::const_iterator chargedHadron = chargedHadrons.begin();
160  chargedHadron != chargedHadrons.end(); ++chargedHadron ) {
162  ++nTracks;
163  }
164  }
165  if ( verbosity_ ) {
166  edm::LogPrint("PFTauByHPSSelect") << "nTracks = " << nTracks << " (min = " << massWindow.nTracksMin_ << ")" ;
167  }
168  if ( nTracks < massWindow.nTracksMin_ ) {
169  if ( verbosity_ ) {
170  edm::LogPrint("PFTauByHPSSelect") << " fails nTracks requirement for mass-window." ;
171  }
172  return 0.0;
173  }
174  }
175  if ( massWindow.nChargedPFCandsMin_ > 0 ) {
176  unsigned int nChargedPFCands = 0;
177  const std::vector<reco::PFRecoTauChargedHadron>& chargedHadrons = tau->signalTauChargedHadronCandidates();
178  for ( std::vector<reco::PFRecoTauChargedHadron>::const_iterator chargedHadron = chargedHadrons.begin();
179  chargedHadron != chargedHadrons.end(); ++chargedHadron ) {
181  ++nChargedPFCands;
182  }
183  }
184  if ( verbosity_ ) {
185  edm::LogPrint("PFTauByHPSSelect") << "nChargedPFCands = " << nChargedPFCands << " (min = " << massWindow.nChargedPFCandsMin_ << ")" ;
186  }
187  if ( nChargedPFCands < massWindow.nChargedPFCandsMin_ ) {
188  if ( verbosity_ ) {
189  edm::LogPrint("PFTauByHPSSelect") << " fails nChargedPFCands requirement for mass-window." ;
190  }
191  return 0.0;
192  }
193  }
194 
195  math::XYZTLorentzVector tauP4 = tau->p4();
196  if ( verbosity_ ) {
197  edm::LogPrint("PFTauByHPSSelect") << "tau: Pt = " << tauP4.pt() << ", eta = " << tauP4.eta() << ", phi = " << tauP4.phi() << ", mass = " << tauP4.mass() ;
198  }
199  // Find the total pizero p4
201  for(auto const& cand : tau->signalPiZeroCandidates()){
202  const math::XYZTLorentzVector& candP4 = cand.p4();
203  stripsP4 += candP4;
204  }
205 
206  // Apply strip mass assumption corrections
207  if (massWindow.assumeStripMass_ >= 0) {
208  for(auto const& cand : tau->signalPiZeroCandidates()){
209  const math::XYZTLorentzVector& uncorrected = cand.p4();
210  math::XYZTLorentzVector corrected =
211  applyMassConstraint(uncorrected, massWindow.assumeStripMass_);
212  math::XYZTLorentzVector correction = corrected - uncorrected;
213  tauP4 += correction;
214  stripsP4 += correction;
215  }
216  }
217  if ( verbosity_ ) {
218  edm::LogPrint("PFTauByHPSSelect") << "strips: Pt = " << stripsP4.pt() << ", eta = " << stripsP4.eta() << ", phi = " << stripsP4.phi() << ", mass = " << stripsP4.mass() ;
219  }
220 
221  // Check if tau fails mass cut
222  double maxMass_value = (*massWindow.maxMass_)(*tau);
223  double bendCorrection_mass = ( massWindow.applyBendCorrection_mass_ ) ? tau->bendCorrMass() : 0.;
224  if ( !((tauP4.M() - bendCorrection_mass) < maxMass_value && (tauP4.M() + bendCorrection_mass) > massWindow.minMass_) ) {
225  if ( verbosity_ ) {
226  edm::LogPrint("PFTauByHPSSelect") << " fails tau mass-window cut." ;
227  }
228  return 0.0;
229  }
230  // CV: require that mass of charged tau decay products is always within specified mass window,
231  // irrespective of bendCorrection_mass
232  reco::Candidate::LorentzVector tauP4_charged;
233  const std::vector<reco::PFRecoTauChargedHadron>& signalChargedHadrons = tau->signalTauChargedHadronCandidates();
234  for ( std::vector<reco::PFRecoTauChargedHadron>::const_iterator signalChargedHadron = signalChargedHadrons.begin();
235  signalChargedHadron != signalChargedHadrons.end(); ++signalChargedHadron ) {
236  tauP4_charged += signalChargedHadron->p4();
237  }
238  if ( !(tauP4_charged.mass() < maxMass_value) ) {
239  if ( verbosity_ ) {
240  edm::LogPrint("PFTauByHPSSelect") << " fails tau mass-window cut." ;
241  }
242  return 0.0;
243  }
244 
245  // Check if it fails the pi0 IM cut
246  if ( stripsP4.M() > massWindow.maxPi0Mass_ ||
247  stripsP4.M() < massWindow.minPi0Mass_ ) {
248  if ( verbosity_ ) {
249  edm::LogPrint("PFTauByHPSSelect") << " fails strip mass-window cut." ;
250  }
251  return 0.0;
252  }
253 
254  // Check if tau passes matching cone cut
255  //edm::LogPrint("PFTauByHPSSelect") << "dR(tau, jet) = " << deltaR(tauP4, tau->jetRef()->p4()) ;
256  if ( deltaR(tauP4, tau->jetRef()->p4()) > matchingCone_ ) {
257  if ( verbosity_ ) {
258  edm::LogPrint("PFTauByHPSSelect") << " fails matching-cone cut." ;
259  }
260  return 0.0;
261  }
262 
263  // Check if tau passes cone cut
264  double cone_size = tau->signalConeSize();
265  // Check if any charged objects fail the signal cone cut
266  for (auto const& cand : tau->signalTauChargedHadronCandidates()) {
267  if ( verbosity_ ) {
268  edm::LogPrint("PFTauByHPSSelect") << "dR(tau, signalPFChargedHadr) = " << deltaR(cand.p4(), tauP4) ;
269  }
270  if ( deltaR(cand.p4(), tauP4) > cone_size ) {
271  if ( verbosity_ ) {
272  edm::LogPrint("PFTauByHPSSelect") << " fails signal-cone cut for charged hadron(s)." ;
273  }
274  return 0.0;
275  }
276  }
277  // Now check the pizeros
278  for (auto const& cand : tau->signalPiZeroCandidates()) {
279  double bendCorrection_eta = ( massWindow.applyBendCorrection_eta_ ) ? cand.bendCorrEta() : 0.;
280  double dEta = std::max(0., fabs(cand.eta() - tauP4.eta()) - bendCorrection_eta);
281  double bendCorrection_phi = ( massWindow.applyBendCorrection_phi_ ) ? cand.bendCorrPhi() : 0.;
282  double dPhi = std::max(0., std::abs(reco::deltaPhi(cand.phi(), tauP4.phi())) - bendCorrection_phi);
283  double dR2 = dEta*dEta + dPhi*dPhi;
284  if ( verbosity_ ) {
285  edm::LogPrint("PFTauByHPSSelect") << "dR2(tau, signalPiZero) = " << dR2 ;
286  }
287  if ( dR2 > cone_size*cone_size ) {
288  if ( verbosity_ ) {
289  edm::LogPrint("PFTauByHPSSelect") << " fails signal-cone cut for strip(s)." ;
290  }
291  return 0.0;
292  }
293  }
294 
296  for(auto const& cand : tau->signalTauChargedHadronCandidates()) {
297  if ( verbosity_ ) {
298  std::string algo_string;
299  if ( cand.algo() == reco::PFRecoTauChargedHadron::kChargedPFCandidate ) algo_string = "ChargedPFCandidate";
300  else if ( cand.algo() == reco::PFRecoTauChargedHadron::kTrack ) algo_string = "Track";
301  else if ( cand.algo() == reco::PFRecoTauChargedHadron::kPFNeutralHadron ) algo_string = "PFNeutralHadron";
302  else algo_string = "Undefined";
303  edm::LogPrint("PFTauByHPSSelect") << "algo(signalPFChargedHadr) = " << algo_string ;
304  }
306  if ( verbosity_ ) {
307  edm::LogPrint("PFTauByHPSSelect") << " fails cut on PFRecoTauChargedHadron algo." ;
308  }
309  return 0.0;
310  }
311  }
312  }
313 
314  if ( minPixelHits_ > 0 ) {
315  int numPixelHits = 0;
316  const std::vector<reco::PFCandidatePtr>& chargedHadrCands = tau->signalPFChargedHadrCands();
317  for ( std::vector<reco::PFCandidatePtr>::const_iterator chargedHadrCand = chargedHadrCands.begin();
318  chargedHadrCand != chargedHadrCands.end(); ++chargedHadrCand ) {
319  const reco::Track* track = nullptr;
320  if ( (*chargedHadrCand)->trackRef().isNonnull() ) track = (*chargedHadrCand)->trackRef().get();
321  else if ( (*chargedHadrCand)->gsfTrackRef().isNonnull() ) track = (*chargedHadrCand)->gsfTrackRef().get();
322  if ( track ) {
323  numPixelHits += track->hitPattern().numberOfValidPixelHits();
324  }
325  }
326  if ( !(numPixelHits >= minPixelHits_) ) {
327  if ( verbosity_ ) {
328  edm::LogPrint("PFTauByHPSSelect") << " fails cut on sum of pixel hits." ;
329  }
330  return 0.0;
331  }
332  }
333 
334  // Otherwise, we pass!
335  if ( verbosity_ ) {
336  edm::LogPrint("PFTauByHPSSelect") << " passes all cuts." ;
337  }
338  return 1.0;
339 }
340 
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:22
T getParameter(std::string const &) const
const unsigned int nTracks(const reco::Vertex &sv)
double discriminate(const reco::PFTauRef &) const override
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
bool exists(std::string const &parameterName) const
checks if a parameter exists
#define nullptr
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
T sqrt(T t)
Definition: SSEVec.h:18
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:479
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:37
int numberOfValidPixelHits() const
Definition: HitPattern.h:908