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&) 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 
52  double matchingCone_;
53  double minPt_;
54 
56 
58 };
59 
62  signalConeFun_(pset.getParameter<std::string>("coneSizeFormula"))
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  verbosity_ = pset.exists("verbosity") ?
107  pset.getParameter<int>("verbosity") : 0;
108 }
109 
111 {
112  for ( DecayModeCutMap::iterator it = decayModeCuts_.begin();
113  it != decayModeCuts_.end(); ++it ) {
114  delete it->second.maxMass_;
115  }
116 }
117 
118 double
120 {
121  if ( verbosity_ ) {
122  edm::LogPrint("PFTauByHPSSelect") << "<PFRecoTauDiscriminationByHPSSelection::discriminate>:" ;
123  edm::LogPrint("PFTauByHPSSelect") << " nCharged = " << tau->signalTauChargedHadronCandidates().size() ;
124  edm::LogPrint("PFTauByHPSSelect") << " nPiZeros = " << tau->signalPiZeroCandidates().size() ;
125  }
126 
127  // Check if we pass the min pt
128  if ( tau->pt() < minPt_ ) {
129  if ( verbosity_ ) {
130  edm::LogPrint("PFTauByHPSSelect") << " fails minPt cut." ;
131  }
132  return 0.0;
133  }
134 
135  // See if we select this decay mode
136  DecayModeCutMap::const_iterator massWindowIter =
137  decayModeCuts_.find(std::make_pair(tau->signalTauChargedHadronCandidates().size(),
138  tau->signalPiZeroCandidates().size()));
139  // Check if decay mode is supported
140  if ( massWindowIter == decayModeCuts_.end() ) {
141  if ( verbosity_ ) {
142  edm::LogPrint("PFTauByHPSSelect") << " fails mass-window definition requirement." ;
143  }
144  return 0.0;
145  }
146  const DecayModeCuts& massWindow = massWindowIter->second;
147 
148  if ( massWindow.nTracksMin_ > 0 ) {
149  unsigned int nTracks = 0;
150  const std::vector<reco::PFRecoTauChargedHadron>& chargedHadrons = tau->signalTauChargedHadronCandidates();
151  for ( std::vector<reco::PFRecoTauChargedHadron>::const_iterator chargedHadron = chargedHadrons.begin();
152  chargedHadron != chargedHadrons.end(); ++chargedHadron ) {
153  if ( chargedHadron->algoIs(reco::PFRecoTauChargedHadron::kChargedPFCandidate) || chargedHadron->algoIs(reco::PFRecoTauChargedHadron::kTrack) ) {
154  ++nTracks;
155  }
156  }
157  if ( verbosity_ ) {
158  edm::LogPrint("PFTauByHPSSelect") << "nTracks = " << nTracks << " (min = " << massWindow.nTracksMin_ << ")" ;
159  }
160  if ( nTracks < massWindow.nTracksMin_ ) {
161  if ( verbosity_ ) {
162  edm::LogPrint("PFTauByHPSSelect") << " fails nTracks requirement for mass-window." ;
163  }
164  return 0.0;
165  }
166  }
167  if ( massWindow.nChargedPFCandsMin_ > 0 ) {
168  unsigned int nChargedPFCands = 0;
169  const std::vector<reco::PFRecoTauChargedHadron>& chargedHadrons = tau->signalTauChargedHadronCandidates();
170  for ( std::vector<reco::PFRecoTauChargedHadron>::const_iterator chargedHadron = chargedHadrons.begin();
171  chargedHadron != chargedHadrons.end(); ++chargedHadron ) {
172  if ( chargedHadron->algoIs(reco::PFRecoTauChargedHadron::kChargedPFCandidate) ) {
173  ++nChargedPFCands;
174  }
175  }
176  if ( verbosity_ ) {
177  edm::LogPrint("PFTauByHPSSelect") << "nChargedPFCands = " << nChargedPFCands << " (min = " << massWindow.nChargedPFCandsMin_ << ")" ;
178  }
179  if ( nChargedPFCands < massWindow.nChargedPFCandsMin_ ) {
180  if ( verbosity_ ) {
181  edm::LogPrint("PFTauByHPSSelect") << " fails nChargedPFCands requirement for mass-window." ;
182  }
183  return 0.0;
184  }
185  }
186 
187  math::XYZTLorentzVector tauP4 = tau->p4();
188  if ( verbosity_ ) {
189  edm::LogPrint("PFTauByHPSSelect") << "tau: Pt = " << tauP4.pt() << ", eta = " << tauP4.eta() << ", phi = " << tauP4.phi() << ", mass = " << tauP4.mass() ;
190  }
191  // Find the total pizero p4
193  BOOST_FOREACH(const reco::RecoTauPiZero& cand,
194  tau->signalPiZeroCandidates()){
195  math::XYZTLorentzVector candP4 = cand.p4();
196  stripsP4 += candP4;
197  }
198 
199  // Apply strip mass assumption corrections
200  if (massWindow.assumeStripMass_ >= 0) {
201  BOOST_FOREACH(const reco::RecoTauPiZero& cand,
202  tau->signalPiZeroCandidates()){
203  math::XYZTLorentzVector uncorrected = cand.p4();
204  math::XYZTLorentzVector corrected =
205  applyMassConstraint(uncorrected, massWindow.assumeStripMass_);
206  math::XYZTLorentzVector correction = corrected - uncorrected;
207  tauP4 += correction;
208  stripsP4 += correction;
209  }
210  }
211  if ( verbosity_ ) {
212  edm::LogPrint("PFTauByHPSSelect") << "strips: Pt = " << stripsP4.pt() << ", eta = " << stripsP4.eta() << ", phi = " << stripsP4.phi() << ", mass = " << stripsP4.mass() ;
213  }
214 
215  // Check if tau fails mass cut
216  double maxMass_value = (*massWindow.maxMass_)(*tau);
217  if ( tauP4.M() > maxMass_value || tauP4.M() < massWindow.minMass_ ) {
218  if ( verbosity_ ) {
219  edm::LogPrint("PFTauByHPSSelect") << " fails tau mass-window cut." ;
220  }
221  return 0.0;
222  }
223 
224  // Check if it fails the pi0 IM cut
225  if ( stripsP4.M() > massWindow.maxPi0Mass_ ||
226  stripsP4.M() < massWindow.minPi0Mass_ ) {
227  if ( verbosity_ ) {
228  edm::LogPrint("PFTauByHPSSelect") << " fails strip mass-window cut." ;
229  }
230  return 0.0;
231  }
232 
233  // Check if tau passes matching cone cut
234  //edm::LogPrint("PFTauByHPSSelect") << "dR(tau, jet) = " << deltaR(tauP4, tau->jetRef()->p4()) ;
235  if ( deltaR(tauP4, tau->jetRef()->p4()) > matchingCone_ ) {
236  if ( verbosity_ ) {
237  edm::LogPrint("PFTauByHPSSelect") << " fails matching-cone cut." ;
238  }
239  return 0.0;
240  }
241 
242  // Check if tau passes cone cut
243  double cone_size = signalConeFun_(*tau);
244  // Check if any charged objects fail the signal cone cut
245  BOOST_FOREACH(const reco::PFRecoTauChargedHadron& cand, tau->signalTauChargedHadronCandidates()) {
246  if ( verbosity_ ) {
247  edm::LogPrint("PFTauByHPSSelect") << "dR(tau, signalPFChargedHadr) = " << deltaR(cand.p4(), tauP4) ;
248  }
249  if ( deltaR(cand.p4(), tauP4) > cone_size ) {
250  if ( verbosity_ ) {
251  edm::LogPrint("PFTauByHPSSelect") << " fails signal-cone cut for charged hadron(s)." ;
252  }
253  return 0.0;
254  }
255  }
256  // Now check the pizeros
257  BOOST_FOREACH(const reco::RecoTauPiZero& cand, tau->signalPiZeroCandidates()) {
258  if ( verbosity_ ) {
259  edm::LogPrint("PFTauByHPSSelect") << "dR(tau, signalPiZero) = " << deltaR(cand.p4(), tauP4) ;
260  }
261  if ( deltaR(cand.p4(), tauP4) > cone_size ) {
262  if ( verbosity_ ) {
263  edm::LogPrint("PFTauByHPSSelect") << " fails signal-cone cut for strip(s)." ;
264  }
265  return 0.0;
266  }
267  }
268 
270  BOOST_FOREACH(const reco::PFRecoTauChargedHadron& cand, tau->signalTauChargedHadronCandidates()) {
271  if ( verbosity_ ) {
272  std::string algo_string;
273  if ( cand.algo() == reco::PFRecoTauChargedHadron::kChargedPFCandidate ) algo_string = "ChargedPFCandidate";
274  else if ( cand.algo() == reco::PFRecoTauChargedHadron::kTrack ) algo_string = "Track";
275  else if ( cand.algo() == reco::PFRecoTauChargedHadron::kPFNeutralHadron ) algo_string = "PFNeutralHadron";
276  else algo_string = "Undefined";
277  edm::LogPrint("PFTauByHPSSelect") << "algo(signalPFChargedHadr) = " << algo_string ;
278  }
280  if ( verbosity_ ) {
281  edm::LogPrint("PFTauByHPSSelect") << " fails cut on PFRecoTauChargedHadron algo." ;
282  }
283  return 0.0;
284  }
285  }
286  }
287 
288  // Otherwise, we pass!
289  if ( verbosity_ ) {
290  edm::LogPrint("PFTauByHPSSelect") << " passes all cuts." ;
291  }
292  return 1.0;
293 }
294 
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
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
double discriminate(const reco::PFTauRef &) override
T sqrt(T t)
Definition: SSEVec.h:48
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:41
virtual const LorentzVector & p4() const
four-momentum Lorentz vector