CMS 3D CMS Logo

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