CMS 3D CMS Logo

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  */
13 
14 #include <algorithm>
15 #include <memory>
16 
17 #include "boost/bind.hpp"
18 
20 
29 
34 
35 //-------------------------------------------------------------------------------
36 // CV: the following headers are needed only for debug print-out
39 //-------------------------------------------------------------------------------
40 
41 namespace reco { namespace tau {
42 
43 namespace {
44  // Apply a hypothesis on the mass of the strips.
45  math::XYZTLorentzVector applyMassConstraint(
46  const math::XYZTLorentzVector& vec,double mass) {
47  double factor = sqrt(vec.energy()*vec.energy()-mass*mass)/vec.P();
49  vec.px()*factor,vec.py()*factor,vec.pz()*factor,vec.energy());
50  }
51 }
52 
54 {
55  public:
57  ~RecoTauPiZeroStripPlugin2() override;
58  // Return type is auto_ptr<PiZeroVector>
59  return_type operator()(const reco::Jet&) const override;
60  // Hook to update PV information
61  void beginEvent() override;
62 
63  private:
64  typedef std::vector<reco::CandidatePtr> CandPtrs;
65  void addCandsToStrip(RecoTauPiZero&, CandPtrs&, const std::vector<bool>&, std::set<size_t>&, bool&) const;
66 
68 
73 
74  double minStripEt_;
75 
76  std::vector<int> inputParticleIds_; // type of candidates to clusterize
77  double etaAssociationDistance_; // size of strip clustering window in eta direction
78  double phiAssociationDistance_; // size of strip clustering window in phi direction
79 
82 
83  // Parameters for build strip combinations
87 
89 
91 };
92 
94  : RecoTauPiZeroBuilderPlugin(pset, std::move(iC)),
95  vertexAssociator_(pset.getParameter<edm::ParameterSet>("qualityCuts"), std::move(iC)),
97 {
98  minGammaEtStripSeed_ = pset.getParameter<double>("minGammaEtStripSeed");
99  minGammaEtStripAdd_ = pset.getParameter<double>("minGammaEtStripAdd");
100 
101  minStripEt_ = pset.getParameter<double>("minStripEt");
102 
103  edm::ParameterSet qcuts_pset = pset.getParameterSet("qualityCuts").getParameterSet("signalQualityCuts");
104 //-------------------------------------------------------------------------------
105 // CV: disable track quality cuts for PFElectronsPFElectron
106 // (treat PFElectrons like PFGammas for the purpose of building eta-phi strips)
107  applyElecTrackQcuts_ = pset.getParameter<bool>("applyElecTrackQcuts");
108  if ( !applyElecTrackQcuts_ ) {
109  qcuts_pset.addParameter<double>("minTrackPt", std::min(minGammaEtStripSeed_, minGammaEtStripAdd_));
110  qcuts_pset.addParameter<double>("maxTrackChi2", 1.e+9);
111  qcuts_pset.addParameter<double>("maxTransverseImpactParameter", 1.e+9);
112  qcuts_pset.addParameter<double>("maxDeltaZ", 1.e+9);
113  qcuts_pset.addParameter<double>("minTrackVertexWeight", -1.);
114  qcuts_pset.addParameter<unsigned>("minTrackPixelHits", 0);
115  qcuts_pset.addParameter<unsigned>("minTrackHits", 0);
116  }
117 //-------------------------------------------------------------------------------
118  qcuts_pset.addParameter<double>("minGammaEt", std::min(minGammaEtStripSeed_, minGammaEtStripAdd_));
119  qcuts_ = new RecoTauQualityCuts(qcuts_pset);
120 
121  inputParticleIds_ = pset.getParameter<std::vector<int> >("stripCandidatesParticleIds");
122  etaAssociationDistance_ = pset.getParameter<double>("stripEtaAssociationDistance");
123  phiAssociationDistance_ = pset.getParameter<double>("stripPhiAssociationDistance");
124 
125  updateStripAfterEachDaughter_ = pset.getParameter<bool>("updateStripAfterEachDaughter");
126  maxStripBuildIterations_ = pset.getParameter<int>("maxStripBuildIterations");
127 
128  combineStrips_ = pset.getParameter<bool>("makeCombinatoricStrips");
129  if ( combineStrips_ ) {
130  maxStrips_ = pset.getParameter<int>("maxInputStrips");
131  combinatoricStripMassHypo_ = pset.getParameter<double>("stripMassWhenCombining");
132  }
133 
134  verbosity_ = pset.getParameter<int>("verbosity");
135 }
136 
138 {
139  delete qcuts_;
140 }
141 
142 // Update the primary vertex
144 {
146 }
147 
148 void RecoTauPiZeroStripPlugin2::addCandsToStrip(RecoTauPiZero& strip, CandPtrs& cands, const std::vector<bool>& candFlags,
149  std::set<size_t>& candIdsCurrentStrip, bool& isCandAdded) const
150 {
151  if ( verbosity_ >= 1 ) {
152  edm::LogPrint("RecoTauPiZeroStripPlugin2") << "<RecoTauPiZeroStripPlugin2::addCandsToStrip>:" ;
153  }
154  size_t numCands = cands.size();
155  for ( size_t candId = 0; candId < numCands; ++candId ) {
156  if ( (!candFlags[candId]) && candIdsCurrentStrip.find(candId) == candIdsCurrentStrip.end() ) { // do not include same cand twice
157  reco::CandidatePtr cand = cands[candId];
158  if ( fabs(strip.eta() - cand->eta()) < etaAssociationDistance_ && // check if cand is within eta-phi window centered on strip
159  fabs(strip.phi() - cand->phi()) < phiAssociationDistance_ ) {
160  if ( verbosity_ >= 2 ) {
161  edm::LogPrint("RecoTauPiZeroStripPlugin2") << "--> adding PFCand #" << candId << " (" << cand.id() << ":" << cand.key() << "): Et = " << cand->et() << ", eta = " << cand->eta() << ", phi = " << cand->phi() ;
162  }
163  strip.addDaughter(cand);
165  isCandAdded = true;
166  candIdsCurrentStrip.insert(candId);
167  }
168  }
169  }
170 }
171 
172 namespace
173 {
174  void markCandsInStrip(std::vector<bool>& candFlags, const std::set<size_t>& candIds)
175  {
176  for ( std::set<size_t>::const_iterator candId = candIds.begin();
177  candId != candIds.end(); ++candId ) {
178  candFlags[*candId] = true;
179  }
180  }
181 
182  inline const reco::TrackBaseRef getTrack(const Candidate& cand)
183  {
184  const PFCandidate* pfCandPtr = dynamic_cast<const PFCandidate*>(&cand);
185  if (pfCandPtr) {
186  if ( pfCandPtr->trackRef().isNonnull() ) return reco::TrackBaseRef(pfCandPtr->trackRef());
187  else if ( pfCandPtr->gsfTrackRef().isNonnull() ) return reco::TrackBaseRef(pfCandPtr->gsfTrackRef());
188  else return reco::TrackBaseRef();
189  }
190 
191  return reco::TrackBaseRef();
192  }
193 }
194 
196 {
197  if ( verbosity_ >= 1 ) {
198  edm::LogPrint("RecoTauPiZeroStripPlugin2") << "<RecoTauPiZeroStripPlugin2::operator()>:" ;
199  edm::LogPrint("RecoTauPiZeroStripPlugin2") << " minGammaEtStripSeed = " << minGammaEtStripSeed_ ;
200  edm::LogPrint("RecoTauPiZeroStripPlugin2") << " minGammaEtStripAdd = " << minGammaEtStripAdd_ ;
201  edm::LogPrint("RecoTauPiZeroStripPlugin2") << " minStripEt = " << minStripEt_ ;
202  }
203 
205 
206  // Get the candidates passing our quality cuts
209 
210  // Convert to stl::list to allow fast deletions
211  CandPtrs seedCands;
212  CandPtrs addCands;
213  int idx = 0;
214  for ( CandPtrs::iterator cand = candsVector.begin();
215  cand != candsVector.end(); ++cand ) {
216  if ( verbosity_ >= 1 ) {
217  edm::LogPrint("RecoTauPiZeroStripPlugin2") << "PFGamma #" << idx << " (" << cand->id() << ":" << cand->key() << "): Et = " << (*cand)->et() << ", eta = " << (*cand)->eta() << ", phi = " << (*cand)->phi() ;
218  }
219  if ( (*cand)->et() > minGammaEtStripSeed_ ) {
220  if ( verbosity_ >= 2 ) {
221  edm::LogPrint("RecoTauPiZeroStripPlugin2") << "--> assigning seedCandId = " << seedCands.size() ;
222  const reco::TrackBaseRef candTrack = getTrack(**cand);
223  if ( candTrack.isNonnull() ) {
224  edm::LogPrint("RecoTauPiZeroStripPlugin2") << "track: Pt = " << candTrack->pt() << " eta = " << candTrack->eta() << ", phi = " << candTrack->phi() << ", charge = " << candTrack->charge() ;
225  edm::LogPrint("RecoTauPiZeroStripPlugin2") << " (dZ = " << candTrack->dz(vertexAssociator_.associatedVertex(jet)->position()) << ", dXY = " << candTrack->dxy(vertexAssociator_.associatedVertex(jet)->position()) << ","
226  << " numHits = " << candTrack->hitPattern().numberOfValidTrackerHits() << ", numPxlHits = " << candTrack->hitPattern().numberOfValidPixelHits() << ","
227  << " chi2 = " << candTrack->normalizedChi2() << ", dPt/Pt = " << (candTrack->ptError()/candTrack->pt()) << ")" ;
228  }
229  }
230  seedCands.push_back(*cand);
231  } else if ( (*cand)->et() > minGammaEtStripAdd_ ) {
232  if ( verbosity_ >= 2 ) {
233  edm::LogPrint("RecoTauPiZeroStripPlugin2") << "--> assigning addCandId = " << addCands.size() ;
234  }
235  addCands.push_back(*cand);
236  }
237  ++idx;
238  }
239 
240  std::vector<bool> seedCandFlags(seedCands.size()); // true/false: seedCand is already/not yet included in strip
241  std::vector<bool> addCandFlags(addCands.size()); // true/false: addCand is already/not yet included in strip
242 
243  std::set<size_t> seedCandIdsCurrentStrip;
244  std::set<size_t> addCandIdsCurrentStrip;
245 
246  size_t idxSeed = 0;
247  while ( idxSeed < seedCands.size() ) {
248  if ( verbosity_ >= 2 ) edm::LogPrint("RecoTauPiZeroStripPlugin2") << "processing seed #" << idxSeed ;
249 
250  seedCandIdsCurrentStrip.clear();
251  addCandIdsCurrentStrip.clear();
252 
253  std::auto_ptr<RecoTauPiZero> strip(new RecoTauPiZero(*seedCands[idxSeed], RecoTauPiZero::kStrips));
254  strip->addDaughter(seedCands[idxSeed]);
255  seedCandIdsCurrentStrip.insert(idxSeed);
256 
257  bool isCandAdded;
258  int stripBuildIteration = 0;
259  do {
260  isCandAdded = false;
261 
262  //if ( verbosity_ >= 2 ) edm::LogPrint("RecoTauPiZeroStripPlugin2") << " adding seedCands to strip..." ;
263  addCandsToStrip(*strip, seedCands, seedCandFlags, seedCandIdsCurrentStrip, isCandAdded);
264  //if ( verbosity_ >= 2 ) edm::LogPrint("RecoTauPiZeroStripPlugin2") << " adding addCands to strip..." ;
265  addCandsToStrip(*strip, addCands, addCandFlags, addCandIdsCurrentStrip, isCandAdded);
266 
268 
269  ++stripBuildIteration;
270  } while ( isCandAdded && (stripBuildIteration < maxStripBuildIterations_ || maxStripBuildIterations_ == -1) );
271 
272  if ( strip->et() > minStripEt_ ) { // strip passed Et cuts, add it to the event
273  if ( verbosity_ >= 2 ) edm::LogPrint("RecoTauPiZeroStripPlugin2") << "Building strip: Et = " << strip->et() << ", eta = " << strip->eta() << ", phi = " << strip->phi() ;
274 
275  // Update the vertex
276  if ( strip->daughterPtr(0).isNonnull() ) strip->setVertex(strip->daughterPtr(0)->vertex());
277  output.push_back(strip);
278 
279  // Mark daughters as being part of this strip
280  markCandsInStrip(seedCandFlags, seedCandIdsCurrentStrip);
281  markCandsInStrip(addCandFlags, addCandIdsCurrentStrip);
282  } else { // strip failed Et cuts, just skip it
283  if ( verbosity_ >= 2 ) edm::LogPrint("RecoTauPiZeroStripPlugin2") << "Discarding strip: Et = " << strip->et() << ", eta = " << strip->eta() << ", phi = " << strip->phi() ;
284  }
285 
286  ++idxSeed;
287  while ( idxSeed < seedCands.size() && seedCandFlags[idxSeed] ) {
288  ++idxSeed; // fast-forward to next seed cand not yet included in any strip
289  }
290  }
291 
292  // Check if we want to combine our strips
293  if ( combineStrips_ && output.size() > 1 ) {
294  PiZeroVector stripCombinations;
295  // Sort the output by descending pt
296  output.sort(output.begin(), output.end(),
297  boost::bind(&RecoTauPiZero::pt, _1) >
298  boost::bind(&RecoTauPiZero::pt, _2));
299  // Get the end of interesting set of strips to try and combine
300  PiZeroVector::const_iterator end_iter = takeNElements(
301  output.begin(), output.end(), maxStrips_);
302 
303  // Look at all the combinations
304  for ( PiZeroVector::const_iterator first = output.begin();
305  first != end_iter-1; ++first ) {
306  for ( PiZeroVector::const_iterator second = first+1;
307  second != end_iter; ++second ) {
308  Candidate::LorentzVector firstP4 = first->p4();
309  Candidate::LorentzVector secondP4 = second->p4();
310  // If we assume a certain mass for each strip apply it here.
311  firstP4 = applyMassConstraint(firstP4, combinatoricStripMassHypo_);
312  secondP4 = applyMassConstraint(secondP4, combinatoricStripMassHypo_);
313  Candidate::LorentzVector totalP4 = firstP4 + secondP4;
314  // Make our new combined strip
315  std::auto_ptr<RecoTauPiZero> combinedStrips(
316  new RecoTauPiZero(0, totalP4,
317  Candidate::Point(0, 0, 0),
318  //111, 10001, true, RecoTauPiZero::kCombinatoricStrips));
319  111, 10001, true, RecoTauPiZero::kUndefined));
320 
321  // Now loop over the strip members
322  for(auto const& gamma : first->daughterPtrVector()) {
323  combinedStrips->addDaughter(gamma);
324  }
325  for(auto const& gamma : second->daughterPtrVector()) {
326  combinedStrips->addDaughter(gamma);
327  }
328  // Update the vertex
329  if ( combinedStrips->daughterPtr(0).isNonnull() )
330  combinedStrips->setVertex(combinedStrips->daughterPtr(0)->vertex());
331  // Add to our collection of combined strips
332  stripCombinations.push_back(combinedStrips);
333  }
334  }
335  // When done doing all the combinations, add the combined strips to the
336  // output.
337  output.transfer(output.end(), stripCombinations);
338  }
339 
340  return output.release();
341 }
342 }} // end namespace reco::tau
343 
346  reco::tau::RecoTauPiZeroStripPlugin2, "RecoTauPiZeroStripPlugin2");
RecoTauPiZeroStripPlugin2(const edm::ParameterSet &, edm::ConsumesCollector &&iC)
std::vector< CandidatePtr > pfCandidates(const Jet &jet, int particleId, bool sort=true)
T getParameter(std::string const &) const
InputIterator takeNElements(const InputIterator &begin, const InputIterator &end, size_t N)
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:251
double eta() const final
momentum pseudorapidity
key_type key() const
Definition: Ptr.h:185
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:600
Coll filterCandRefs(const Coll &refcoll, bool invert=false) const
Filter a ref vector of Candidates.
Base class for all types of Jets.
Definition: Jet.h:20
#define nullptr
double phi() const
azimuthal angle of momentum vector
Definition: TrackBase.h:684
double pt() const final
transverse momentum
void setEvent(const edm::Event &evt)
Load the vertices from the event.
bool isNonnull() const
Checks for non-null.
Definition: RefToBase.h:340
reco::TrackRef trackRef() const
Definition: PFCandidate.cc:442
U second(std::pair< T, U > const &p)
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:690
edm::RefToBase< reco::Track > TrackBaseRef
persistent reference to a Track, using views
Definition: TrackFwd.h:36
void setPV(const reco::VertexRef &vtx) const
Update the primary vertex.
T sqrt(T t)
Definition: SSEVec.h:18
double pt() const
track transverse momentum
Definition: TrackBase.h:660
double ptError() const
error on Pt (set to 1000 TeV if charge==0 for safety)
Definition: TrackBase.h:814
T min(T a, T b)
Definition: MathUtil.h:58
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:648
static const TrackGhostTrackState * getTrack(const BasicGhostTrackState *basic)
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:485
ParameterSet const & getParameterSet(std::string const &) const
ProductID id() const
Accessor for product ID.
Definition: Ptr.h:180
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:37
int numberOfValidTrackerHits() const
Definition: HitPattern.h:901
std::auto_ptr< PiZeroVector > return_type
return_type operator()(const reco::Jet &) const override
Build a collection of piZeros from objects in the input jet.
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:40
fixed size matrix
HLT enums.
reco::VertexRef associatedVertex(const Jet &jet) const
boost::ptr_vector< RecoTauPiZero > PiZeroVector
reco::GsfTrackRef gsfTrackRef() const
Definition: PFCandidate.cc:480
void addDaughter(const CandidatePtr &)
add a daughter via a reference
int numberOfValidPixelHits() const
Definition: HitPattern.h:916
std::vector< reco::CandidatePtr > CandPtrs
math::XYZPoint Point
point in the space
Definition: Candidate.h:41
int charge() const
track electric charge
Definition: TrackBase.h:606
#define DEFINE_EDM_PLUGIN(factory, type, name)
void set(reco::Candidate &c) const
set up a candidate
void addCandsToStrip(RecoTauPiZero &, CandPtrs &, const std::vector< bool > &, std::set< size_t > &, bool &) const
void beginEvent() override
Hook called at the beginning of the event.
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:630
double phi() const final
momentum azimuthal angle
def move(src, dest)
Definition: eostools.py:511