CMS 3D CMS Logo

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