CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
RecoTauPiZeroStripPlugin2.cc
Go to the documentation of this file.
1 /*
2  * RecoTauPiZeroStripPlugin2
3  *
4  * Merges PFGammas in a PFJet into Candidate piZeros defined as
5  * strips in eta-phi.
6  *
7  * Author: Michail Bachtis (University of Wisconsin)
8  *
9  * Code modifications: Evan Friis (UC Davis),
10  * Christian Veelken (LLR)
11  *
12  * $Id $
13  */
14 #include <algorithm>
15 #include <memory>
16 
18 
27 
32 
33 //-------------------------------------------------------------------------------
34 // CV: the following headers are needed only for debug print-out
37 //-------------------------------------------------------------------------------
38 
39 namespace reco { namespace tau {
40 
41 namespace {
42  // Apply a hypothesis on the mass of the strips.
43  math::XYZTLorentzVector applyMassConstraint(
44  const math::XYZTLorentzVector& vec,double mass) {
45  double factor = sqrt(vec.energy()*vec.energy()-mass*mass)/vec.P();
47  vec.px()*factor,vec.py()*factor,vec.pz()*factor,vec.energy());
48  }
49 }
50 
52 {
53  public:
56  // Return type is auto_ptr<PiZeroVector>
57  return_type operator()(const reco::PFJet&) const;
58  // Hook to update PV information
59  virtual void beginEvent();
60 
61  private:
62  typedef std::vector<reco::PFCandidatePtr> PFCandPtrs;
63  void addCandsToStrip(RecoTauPiZero&, PFCandPtrs&, const std::vector<bool>&, std::set<size_t>&, bool&) const;
64 
66 
71 
72  double minStripEt_;
73 
74  std::vector<int> inputPdgIds_; // type of candidates to clusterize
75  double etaAssociationDistance_; // size of strip clustering window in eta direction
76  double phiAssociationDistance_; // size of strip clustering window in phi direction
77 
80 
81  // Parameters for build strip combinations
85 
87 };
88 
91  vertexAssociator_(pset.getParameter<edm::ParameterSet>("qualityCuts")),
92  qcuts_(0)
93 {
94  LogTrace("RecoTauPiZeroStripPlugin2") << "<RecoTauPiZeroStripPlugin2::RecoTauPiZeroStripPlugin2>:" ;
95 
96  minGammaEtStripSeed_ = pset.getParameter<double>("minGammaEtStripSeed");
97  LogTrace("RecoTauPiZeroStripPlugin2") << " minGammaEtStripSeed = " << minGammaEtStripSeed_ ;
98  minGammaEtStripAdd_ = pset.getParameter<double>("minGammaEtStripAdd");
99  LogTrace("RecoTauPiZeroStripPlugin2") << " minGammaEtStripAdd = " << minGammaEtStripAdd_ ;
100 
101  minStripEt_ = pset.getParameter<double>("minStripEt");
102  LogTrace("RecoTauPiZeroStripPlugin2") << " minStripEt = " << minStripEt_ ;
103 
104  edm::ParameterSet qcuts_pset = pset.getParameterSet("qualityCuts").getParameterSet("signalQualityCuts");
105 //-------------------------------------------------------------------------------
106 // CV: disable track quality cuts for PFElectronsPFElectron
107 // (treat PFElectrons like PFGammas for the purpose of building eta-phi strips)
108  applyElecTrackQcuts_ = pset.getParameter<bool>("applyElecTrackQcuts");
109  if ( !applyElecTrackQcuts_ ) {
110  qcuts_pset.addParameter<double>("minTrackPt", std::min(minGammaEtStripSeed_, minGammaEtStripAdd_));
111  qcuts_pset.addParameter<double>("maxTrackChi2", 1.e+9);
112  qcuts_pset.addParameter<double>("maxTransverseImpactParameter", 1.e+9);
113  qcuts_pset.addParameter<double>("maxDeltaZ", 1.);
114  qcuts_pset.addParameter<double>("minTrackVertexWeight", -1.);
115  qcuts_pset.addParameter<unsigned>("minTrackPixelHits", 0);
116  qcuts_pset.addParameter<unsigned>("minTrackHits", 0);
117  }
118 //-------------------------------------------------------------------------------
119  qcuts_pset.addParameter<double>("minGammaEt", std::min(minGammaEtStripSeed_, minGammaEtStripAdd_));
120  qcuts_ = new RecoTauQualityCuts(qcuts_pset);
121 
122  inputPdgIds_ = pset.getParameter<std::vector<int> >("stripCandidatesParticleIds");
123  etaAssociationDistance_ = pset.getParameter<double>("stripEtaAssociationDistance");
124  phiAssociationDistance_ = pset.getParameter<double>("stripPhiAssociationDistance");
125 
126  updateStripAfterEachDaughter_ = pset.getParameter<bool>("updateStripAfterEachDaughter");
127  maxStripBuildIterations_ = pset.getParameter<int>("maxStripBuildIterations");
128 
129  combineStrips_ = pset.getParameter<bool>("makeCombinatoricStrips");
130  if ( combineStrips_ ) {
131  maxStrips_ = pset.getParameter<int>("maxInputStrips");
132  combinatoricStripMassHypo_ = pset.getParameter<double>("stripMassWhenCombining");
133  }
134 }
135 
137 {
138  delete qcuts_;
139 }
140 
141 // Update the primary vertex
143 {
145 }
146 
147 void RecoTauPiZeroStripPlugin2::addCandsToStrip(RecoTauPiZero& strip, PFCandPtrs& cands, const std::vector<bool>& candFlags,
148  std::set<size_t>& candIdsCurrentStrip, bool& isCandAdded) const
149 {
150  size_t numCands = cands.size();
151  for ( size_t candId = 0; candId < numCands; ++candId ) {
152  if ( (!candFlags[candId]) && candIdsCurrentStrip.find(candId) == candIdsCurrentStrip.end() ) { // do not include same cand twice
153  reco::PFCandidatePtr cand = cands[candId];
154  if ( fabs(strip.eta() - cand->eta()) < etaAssociationDistance_ && // check if cand is within eta-phi window centered on strip
155  fabs(strip.phi() - cand->phi()) < phiAssociationDistance_ ) {
156  LogTrace("RecoTauPiZeroStripPlugin2") << "--> candId = " << candId << " has been added." ;
157  strip.addDaughter(cand);
159  isCandAdded = true;
160  candIdsCurrentStrip.insert(candId);
161  }
162  }
163  }
164 }
165 
166 void markCandsInStrip(std::vector<bool>& candFlags, const std::set<size_t>& candIds)
167 {
168  for ( std::set<size_t>::const_iterator candId = candIds.begin();
169  candId != candIds.end(); ++candId ) {
170  candFlags[*candId] = true;
171  }
172 }
173 
174 namespace {
175  const reco::TrackBaseRef getTrack(const PFCandidate& cand)
176  {
177  if ( cand.trackRef().isNonnull() ) return reco::TrackBaseRef(cand.trackRef());
178  else if ( cand.gsfTrackRef().isNonnull() ) return reco::TrackBaseRef(cand.gsfTrackRef());
179  else return reco::TrackBaseRef();
180  }
181 }
182 
184 {
186 
187  // Get the candidates passing our quality cuts
189  PFCandPtrs candsVector = qcuts_->filterRefs(pfCandidates(jet, inputPdgIds_));
190 
191  // Convert to stl::list to allow fast deletions
192  PFCandPtrs seedCands;
193  PFCandPtrs addCands;
194  int idx = 0;
195  for ( PFCandPtrs::iterator cand = candsVector.begin();
196  cand != candsVector.end(); ++cand ) {
197  LogTrace("RecoTauPiZeroStripPlugin2") << "PFGamma (" << idx << "): Et = " << (*cand)->et() << ","
198  << " eta = " << (*cand)->eta() << ", phi = " << (*cand)->phi();
199  if ( (*cand)->et() > minGammaEtStripSeed_ ) {
200  LogTrace("RecoTauPiZeroStripPlugin2") << " --> assigning seedCandId = " << seedCands.size() ;
201 #ifdef EDM_ML_DEBUG
202  const reco::TrackBaseRef candTrack = getTrack(*cand);
203  if ( candTrack.isNonnull() ) {
204  LogTrace("RecoTauPiZeroStripPlugin2") << "has Track: pt = " << candTrack->pt() << " +/- " << candTrack->ptError() << ","
205  << " eta = " << candTrack->eta() << ", phi = " << candTrack->phi() << ","
206  << " charge = " << candTrack->charge() ;
207  LogTrace("RecoTauPiZeroStripPlugin2") << " chi2 = " << candTrack->normalizedChi2() ;
208  LogTrace("RecoTauPiZeroStripPlugin2") << " dIP = " << std::abs(candTrack->dxy(vertexAssociator_.associatedVertex(jet)->position())) ;
209  LogTrace("RecoTauPiZeroStripPlugin2") << " dZ = " << std::abs(candTrack->dz(vertexAssociator_.associatedVertex(jet)->position())) ;
210  LogTrace("RecoTauPiZeroStripPlugin2") << " vtxAssocWeight = " << vertexAssociator_.associatedVertex(jet)->trackWeight(candTrack) ;
211  LogTrace("RecoTauPiZeroStripPlugin2") << " numPxlHits = " << candTrack->hitPattern().numberOfValidPixelHits() ;
212  LogTrace("RecoTauPiZeroStripPlugin2") << " numTrkHits = " << candTrack->hitPattern().numberOfValidHits() ;
213  }
214 #endif
215  LogTrace("RecoTauPiZeroStripPlugin2") << "ECAL Et: calibrated = " << (*cand)->ecalEnergy()*sin((*cand)->theta()) << ","
216  << " raw = " << (*cand)->rawEcalEnergy()*sin((*cand)->theta()) ;
217  LogTrace("RecoTauPiZeroStripPlugin2") << "HCAL Et: calibrated = " << (*cand)->hcalEnergy()*sin((*cand)->theta()) << ","
218  << " raw = " << (*cand)->rawHcalEnergy()*sin((*cand)->theta()) ;
219  seedCands.push_back(*cand);
220  } else if ( (*cand)->et() > minGammaEtStripAdd_ ) {
221  LogTrace("RecoTauPiZeroStripPlugin2") << " --> assigning addCandId = " << addCands.size();
222  addCands.push_back(*cand);
223  }
224  LogTrace("RecoTauPiZeroStripPlugin2") ;
225  ++idx;
226  }
227 
228  std::vector<bool> seedCandFlags(seedCands.size()); // true/false: seedCand is already/not yet included in strip
229  std::vector<bool> addCandFlags(addCands.size()); // true/false: addCand is already/not yet included in strip
230 
231  std::set<size_t> seedCandIdsCurrentStrip;
232  std::set<size_t> addCandIdsCurrentStrip;
233 
234  size_t idxSeed = 0;
235  while ( idxSeed < seedCands.size() ) {
236  LogTrace("RecoTauPiZeroStripPlugin2") << "idxSeed = " << idxSeed ;
237 
238  seedCandIdsCurrentStrip.clear();
239  addCandIdsCurrentStrip.clear();
240 
241  std::auto_ptr<RecoTauPiZero> strip(new RecoTauPiZero(*seedCands[idxSeed], RecoTauPiZero::kStrips));
242  strip->addDaughter(seedCands[idxSeed]);
243  seedCandIdsCurrentStrip.insert(idxSeed);
244 
245  bool isCandAdded;
246  int stripBuildIteration = 0;
247  do {
248  isCandAdded = false;
249 
250  LogTrace("RecoTauPiZeroStripPlugin2") << " adding seedCands to strip..." ;
251  addCandsToStrip(*strip, seedCands, seedCandFlags, seedCandIdsCurrentStrip, isCandAdded);
252  LogTrace("RecoTauPiZeroStripPlugin2") << " adding addCands to strip..." ;
253  addCandsToStrip(*strip, addCands, addCandFlags, addCandIdsCurrentStrip, isCandAdded);
254 
256 
257  ++stripBuildIteration;
258  } while ( isCandAdded && (stripBuildIteration < maxStripBuildIterations_ || maxStripBuildIterations_ == -1) );
259 
260  if ( strip->et() > minStripEt_ ) { // strip passed Et cuts, add it to the event
261  LogTrace("RecoTauPiZeroStripPlugin2") << "Strip: Et = " << strip->et() << ","
262  << " eta = " << strip->eta() << ", phi = " << strip->phi()
263  << " --> building it !!" ;
264 
265  // Update the vertex
266  if ( strip->daughterPtr(0).isNonnull() ) strip->setVertex(strip->daughterPtr(0)->vertex());
267  output.push_back(strip);
268 
269  // Mark daughters as being part of this strip
270  markCandsInStrip(seedCandFlags, seedCandIdsCurrentStrip);
271  markCandsInStrip(addCandFlags, addCandIdsCurrentStrip);
272  } else { // strip failed Et cuts, just skip it
273  LogTrace("RecoTauPiZeroStripPlugin2") << "Strip: Et = " << strip->et() << ","
274  << " eta = " << strip->eta() << ", phi = " << strip->phi()
275  << " --> discarding it !!" ;
276  }
277 
278  ++idxSeed;
279  while ( idxSeed < seedCands.size() && seedCandFlags[idxSeed] ) {
280  ++idxSeed; // fast-forward to next seed cand not yet included in any strip
281  }
282  }
283 
284  // Check if we want to combine our strips
285  if ( combineStrips_ && output.size() > 1 ) {
286  PiZeroVector stripCombinations;
287  // Sort the output by descending pt
288  output.sort(output.begin(), output.end(),
289  boost::bind(&RecoTauPiZero::pt, _1) >
290  boost::bind(&RecoTauPiZero::pt, _2));
291  // Get the end of interesting set of strips to try and combine
292  PiZeroVector::const_iterator end_iter = takeNElements(
293  output.begin(), output.end(), maxStrips_);
294 
295  // Look at all the combinations
296  for ( PiZeroVector::const_iterator first = output.begin();
297  first != end_iter-1; ++first ) {
298  for ( PiZeroVector::const_iterator second = first+1;
299  second != end_iter; ++second ) {
300  Candidate::LorentzVector firstP4 = first->p4();
301  Candidate::LorentzVector secondP4 = second->p4();
302  // If we assume a certain mass for each strip apply it here.
303  firstP4 = applyMassConstraint(firstP4, combinatoricStripMassHypo_);
304  secondP4 = applyMassConstraint(secondP4, combinatoricStripMassHypo_);
305  Candidate::LorentzVector totalP4 = firstP4 + secondP4;
306  // Make our new combined strip
307  std::auto_ptr<RecoTauPiZero> combinedStrips(
308  new RecoTauPiZero(0, totalP4,
309  Candidate::Point(0, 0, 0),
310  //111, 10001, true, RecoTauPiZero::kCombinatoricStrips));
311  111, 10001, true, RecoTauPiZero::kUndefined));
312 
313  // Now loop over the strip members
314  BOOST_FOREACH(const RecoTauPiZero::daughters::value_type& gamma,
315  first->daughterPtrVector()) {
316  combinedStrips->addDaughter(gamma);
317  }
318  BOOST_FOREACH(const RecoTauPiZero::daughters::value_type& gamma,
319  second->daughterPtrVector()) {
320  combinedStrips->addDaughter(gamma);
321  }
322  // Update the vertex
323  if ( combinedStrips->daughterPtr(0).isNonnull() )
324  combinedStrips->setVertex(combinedStrips->daughterPtr(0)->vertex());
325  // Add to our collection of combined strips
326  stripCombinations.push_back(combinedStrips);
327  }
328  }
329  // When done doing all the combinations, add the combined strips to the
330  // output.
331  output.transfer(output.end(), stripCombinations);
332  }
333 
334  return output.release();
335 }
336 }} // end namespace reco::tau
337 
340  reco::tau::RecoTauPiZeroStripPlugin2, "RecoTauPiZeroStripPlugin2");
T getParameter(std::string const &) const
reco::VertexRef associatedVertex(const PFJet &tau) const
InputIterator takeNElements(const InputIterator &begin, const InputIterator &end, size_t N)
Coll filterRefs(const Coll &refcoll, bool invert=false) const
Filter a ref vector of PFCandidates.
void strip(std::string &input, const std::string &blanks=" \n\t")
Definition: stringTools.cc:16
double normalizedChi2() const
chi-squared divided by n.d.o.f. (or chi-squared * 1e6 if n.d.o.f. is zero)
Definition: TrackBase.h:111
void markCandsInStrip(std::vector< bool > &candFlags, const std::set< size_t > &candIds)
int numberOfValidHits() const
Definition: HitPattern.h:554
virtual const_iterator end() const
last daughter const_iterator
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
#define abs(x)
Definition: mlp_lapack.h:159
double phi() const
azimuthal angle of momentum vector
Definition: TrackBase.h:139
void setEvent(const edm::Event &evt)
Load the vertices from the event.
#define min(a, b)
Definition: mlp_lapack.h:161
Jets made from PFObjects.
Definition: PFJet.h:22
virtual double eta() const
momentum pseudorapidity
std::vector< PFCandidatePtr > pfCandidates(const PFJet &jet, int particleId, bool sort=true)
reco::TrackRef trackRef() const
Definition: PFCandidate.cc:339
U second(std::pair< T, U > const &p)
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:30
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:250
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:141
edm::RefToBase< reco::Track > TrackBaseRef
persistent reference to a Track, using views
Definition: TrackFwd.h:22
void setPV(const reco::VertexRef &vtx) const
Update the primary vertex.
T sqrt(T t)
Definition: SSEVec.h:46
double pt() const
track transverse momentum
Definition: TrackBase.h:131
void addParameter(std::string const &name, T const &value)
Definition: ParameterSet.h:145
virtual const_iterator begin() const
first daughter const_iterator
double ptError() const
error on Pt (set to 1000 TeV if charge==0 for safety)
Definition: TrackBase.h:194
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:223
bool first
Definition: L1TdeRCT.cc:94
Container::value_type value_type
#define LogTrace(id)
double dz() const
dz parameter (= dsz/cos(lambda)). This is the track z0 w.r.t (0,0,0) only if the refPoint is close to...
Definition: TrackBase.h:127
RecoTauPiZeroStripPlugin2(const edm::ParameterSet &)
ParameterSet const & getParameterSet(std::string const &) const
virtual double pt() const
transverse momentum
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:38
std::auto_ptr< PiZeroVector > return_type
virtual void beginEvent()
Hook called at the beginning of the event.
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:33
boost::ptr_vector< RecoTauPiZero > PiZeroVector
reco::GsfTrackRef gsfTrackRef() const
Definition: PFCandidate.cc:376
tuple mass
Definition: scaleCards.py:27
const reco::TrackBaseRef getTrack(const reco::PFCandidate &cand)
std::vector< reco::PFCandidatePtr > PFCandPtrs
void addDaughter(const CandidatePtr &)
add a daughter via a reference
int numberOfValidPixelHits() const
Definition: HitPattern.h:566
math::XYZPoint Point
point in the space
Definition: Candidate.h:42
int charge() const
track electric charge
Definition: TrackBase.h:113
#define DEFINE_EDM_PLUGIN(factory, type, name)
void set(reco::Candidate &c) const
set up a candidate
double dxy() const
dxy parameter. (This is the transverse impact parameter w.r.t. to (0,0,0) ONLY if refPoint is close t...
Definition: TrackBase.h:121
virtual double phi() const
momentum azimuthal angle
bool isNonnull() const
Checks for non-null.
Definition: RefToBase.h:279
void addCandsToStrip(RecoTauPiZero &, PFCandPtrs &, const std::vector< bool > &, std::set< size_t > &, bool &) const
return_type operator()(const reco::PFJet &) const
Build a collection of piZeros from objects in the input jet.