CMS 3D CMS Logo

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