CMS 3D CMS Logo

RecoTauPiZeroStripPlugin3.cc
Go to the documentation of this file.
1 /*
2  * RecoTauPiZeroStripPlugin3
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 #include "TString.h"
42 #include "TFormula.h"
43 
44 namespace reco {
45  namespace tau {
46 
47  namespace {
48  // Apply a hypothesis on the mass of the strips.
49  math::XYZTLorentzVector applyMassConstraint(const math::XYZTLorentzVector& vec, double mass) {
50  double factor = sqrt(vec.energy() * vec.energy() - mass * mass) / vec.P();
51  return math::XYZTLorentzVector(vec.px() * factor, vec.py() * factor, vec.pz() * factor, vec.energy());
52  }
53  } // namespace
54 
56  public:
58  ~RecoTauPiZeroStripPlugin3() override;
59  // Return type is unique_ptr<PiZeroVector>
60  return_type operator()(const reco::Jet&) const override;
61  // Hook to update PV information
62  void beginEvent() override;
63 
64  private:
65  typedef std::vector<reco::CandidatePtr> CandPtrs;
66  void addCandsToStrip(RecoTauPiZero&, CandPtrs&, const std::vector<bool>&, std::set<size_t>&, bool&) const;
67 
69 
70  std::unique_ptr<RecoTauQualityCuts> qcuts_;
74 
75  double minStripEt_;
76 
77  std::vector<int> inputParticleIds_; // type of candidates to clusterize
78  std::unique_ptr<const TFormula> etaAssociationDistance_; // size of strip clustering window in eta direction
79  std::unique_ptr<const TFormula> phiAssociationDistance_; // size of strip clustering window in phi direction
80 
83 
84  // Parameters for build strip combinations
88 
90 
92  };
93 
94  namespace {
95  std::unique_ptr<TFormula> makeFunction(const std::string& functionName, const edm::ParameterSet& pset) {
96  TString formula = pset.getParameter<std::string>("function");
97  formula = formula.ReplaceAll("pT", "x");
98  std::unique_ptr<TFormula> function(new TFormula(functionName.data(), formula.Data()));
99  int numParameter = function->GetNpar();
100  for (int idxParameter = 0; idxParameter < numParameter; ++idxParameter) {
101  std::string parameterName = Form("par%i", idxParameter);
102  double parameter = pset.getParameter<double>(parameterName);
103  function->SetParameter(idxParameter, parameter);
104  }
105  return function;
106  }
107  } // namespace
108 
110  : RecoTauPiZeroBuilderPlugin(pset, std::move(iC)),
111  vertexAssociator_(pset.getParameter<edm::ParameterSet>("qualityCuts"), std::move(iC)),
112  qcuts_(nullptr),
115  minGammaEtStripSeed_ = pset.getParameter<double>("minGammaEtStripSeed");
116  minGammaEtStripAdd_ = pset.getParameter<double>("minGammaEtStripAdd");
117 
118  minStripEt_ = pset.getParameter<double>("minStripEt");
119 
120  edm::ParameterSet qcuts_pset = pset.getParameterSet("qualityCuts").getParameterSet("signalQualityCuts");
121  //-------------------------------------------------------------------------------
122  // CV: disable track quality cuts for PFElectronsPFElectron
123  // (treat PFElectrons like PFGammas for the purpose of building eta-phi strips)
124  applyElecTrackQcuts_ = pset.getParameter<bool>("applyElecTrackQcuts");
125  if (!applyElecTrackQcuts_) {
126  qcuts_pset.addParameter<double>("minTrackPt", std::min(minGammaEtStripSeed_, minGammaEtStripAdd_));
127  qcuts_pset.addParameter<double>("maxTrackChi2", 1.e+9);
128  qcuts_pset.addParameter<double>("maxTransverseImpactParameter", 1.e+9);
129  qcuts_pset.addParameter<double>("maxDeltaZ", 1.e+9);
130  qcuts_pset.addParameter<double>("minTrackVertexWeight", -1.);
131  qcuts_pset.addParameter<unsigned>("minTrackPixelHits", 0);
132  qcuts_pset.addParameter<unsigned>("minTrackHits", 0);
133  }
134  //-------------------------------------------------------------------------------
135  qcuts_pset.addParameter<double>("minGammaEt", std::min(minGammaEtStripSeed_, minGammaEtStripAdd_));
136  //qcuts_ = new RecoTauQualityCuts(qcuts_pset);
137  //std::unique_ptr<RecoTauQualityCuts> qcuts_(new RecoTauQualityCuts(qcuts_pset));
138 
139  qcuts_.reset(new RecoTauQualityCuts(qcuts_pset));
140 
141  inputParticleIds_ = pset.getParameter<std::vector<int> >("stripCandidatesParticleIds");
142  const edm::ParameterSet& stripSize_eta_pset = pset.getParameterSet("stripEtaAssociationDistanceFunc");
143  etaAssociationDistance_ = makeFunction("etaAssociationDistance", stripSize_eta_pset);
144  const edm::ParameterSet& stripSize_phi_pset = pset.getParameterSet("stripPhiAssociationDistanceFunc");
145  phiAssociationDistance_ = makeFunction("phiAssociationDistance", stripSize_phi_pset);
146 
147  updateStripAfterEachDaughter_ = pset.getParameter<bool>("updateStripAfterEachDaughter");
148  maxStripBuildIterations_ = pset.getParameter<int>("maxStripBuildIterations");
149 
150  combineStrips_ = pset.getParameter<bool>("makeCombinatoricStrips");
151  if (combineStrips_) {
152  maxStrips_ = pset.getParameter<int>("maxInputStrips");
153  combinatoricStripMassHypo_ = pset.getParameter<double>("stripMassWhenCombining");
154  }
155 
156  verbosity_ = pset.getParameter<int>("verbosity");
157  }
159 
160  // Update the primary vertex
162 
164  CandPtrs& cands,
165  const std::vector<bool>& candFlags,
166  std::set<size_t>& candIdsCurrentStrip,
167  bool& isCandAdded) const {
168  if (verbosity_ >= 1) {
169  edm::LogPrint("RecoTauPiZeroStripPlugin3") << "<RecoTauPiZeroStripPlugin3::addCandsToStrip>:";
170  }
171  size_t numCands = cands.size();
172  for (size_t candId = 0; candId < numCands; ++candId) {
173  if ((!candFlags[candId]) &&
174  candIdsCurrentStrip.find(candId) == candIdsCurrentStrip.end()) { // do not include same cand twice
175  reco::CandidatePtr cand = cands[candId];
176  double etaAssociationDistance_value =
177  etaAssociationDistance_->Eval(strip.pt()) + etaAssociationDistance_->Eval(cand->pt());
178  double phiAssociationDistance_value =
179  phiAssociationDistance_->Eval(strip.pt()) + phiAssociationDistance_->Eval(cand->pt());
180  if (fabs(strip.eta() - cand->eta()) <
181  etaAssociationDistance_value && // check if cand is within eta-phi window centered on strip
182  reco::deltaPhi(strip.phi(), cand->phi()) < phiAssociationDistance_value) {
183  if (verbosity_ >= 2) {
184  edm::LogPrint("RecoTauPiZeroStripPlugin3")
185  << "--> adding PFCand #" << candId << " (" << cand.id() << ":" << cand.key()
186  << "): Et = " << cand->et() << ", eta = " << cand->eta() << ", phi = " << cand->phi();
187  }
188  strip.addDaughter(cand);
190  p4Builder_.set(strip);
191  isCandAdded = true;
192  candIdsCurrentStrip.insert(candId);
193  }
194  }
195  }
196  }
197 
198  namespace {
199  void markCandsInStrip(std::vector<bool>& candFlags, const std::set<size_t>& candIds) {
200  for (std::set<size_t>::const_iterator candId = candIds.begin(); candId != candIds.end(); ++candId) {
201  candFlags[*candId] = true;
202  }
203  }
204 
205  // The following method has not been adapted to work with PackedCandidates,
206  // however, it is only used for printing/debugging and can be adapted when
207  // needed.
208  inline const reco::TrackBaseRef getTrack(const Candidate& cand) {
209  const PFCandidate* pfCandPtr = dynamic_cast<const PFCandidate*>(&cand);
210  if (pfCandPtr) {
211  if (pfCandPtr->trackRef().isNonnull())
212  return reco::TrackBaseRef(pfCandPtr->trackRef());
213  else if (pfCandPtr->gsfTrackRef().isNonnull())
214  return reco::TrackBaseRef(pfCandPtr->gsfTrackRef());
215  else
216  return reco::TrackBaseRef();
217  }
218 
219  return reco::TrackBaseRef();
220  }
221  } // namespace
222 
224  if (verbosity_ >= 1) {
225  edm::LogPrint("RecoTauPiZeroStripPlugin3") << "<RecoTauPiZeroStripPlugin3::operator()>:";
226  edm::LogPrint("RecoTauPiZeroStripPlugin3") << " minGammaEtStripSeed = " << minGammaEtStripSeed_;
227  edm::LogPrint("RecoTauPiZeroStripPlugin3") << " minGammaEtStripAdd = " << minGammaEtStripAdd_;
228  edm::LogPrint("RecoTauPiZeroStripPlugin3") << " minStripEt = " << minStripEt_;
229  }
230 
232 
233  // Get the candidates passing our quality cuts
235  CandPtrs candsVector = qcuts_->filterCandRefs(pfCandidates(jet, inputParticleIds_));
236 
237  // Convert to stl::list to allow fast deletions
238  CandPtrs seedCands;
239  CandPtrs addCands;
240  int idx = 0;
241  for (CandPtrs::iterator cand = candsVector.begin(); cand != candsVector.end(); ++cand) {
242  if (verbosity_ >= 1) {
243  edm::LogPrint("RecoTauPiZeroStripPlugin3")
244  << "PFGamma #" << idx << " (" << cand->id() << ":" << cand->key() << "): Et = " << (*cand)->et()
245  << ", eta = " << (*cand)->eta() << ", phi = " << (*cand)->phi();
246  }
247  if ((*cand)->et() > minGammaEtStripSeed_) {
248  if (verbosity_ >= 2) {
249  edm::LogPrint("RecoTauPiZeroStripPlugin3") << "--> assigning seedCandId = " << seedCands.size();
250  const reco::TrackBaseRef candTrack = getTrack(**cand);
251  if (candTrack.isNonnull()) {
252  edm::LogPrint("RecoTauPiZeroStripPlugin3")
253  << "track: Pt = " << candTrack->pt() << " eta = " << candTrack->eta()
254  << ", phi = " << candTrack->phi() << ", charge = " << candTrack->charge();
255  edm::LogPrint("RecoTauPiZeroStripPlugin3")
256  << " (dZ = " << candTrack->dz(vertexAssociator_.associatedVertex(jet)->position())
257  << ", dXY = " << candTrack->dxy(vertexAssociator_.associatedVertex(jet)->position()) << ","
258  << " numHits = " << candTrack->hitPattern().numberOfValidTrackerHits()
259  << ", numPxlHits = " << candTrack->hitPattern().numberOfValidPixelHits() << ","
260  << " chi2 = " << candTrack->normalizedChi2()
261  << ", dPt/Pt = " << (candTrack->ptError() / candTrack->pt()) << ")";
262  }
263  }
264  seedCands.push_back(*cand);
265  } else if ((*cand)->et() > minGammaEtStripAdd_) {
266  if (verbosity_ >= 2) {
267  edm::LogPrint("RecoTauPiZeroStripPlugin3") << "--> assigning addCandId = " << addCands.size();
268  }
269  addCands.push_back(*cand);
270  }
271  ++idx;
272  }
273 
274  std::vector<bool> seedCandFlags(seedCands.size()); // true/false: seedCand is already/not yet included in strip
275  std::vector<bool> addCandFlags(addCands.size()); // true/false: addCand is already/not yet included in strip
276 
277  std::set<size_t> seedCandIdsCurrentStrip;
278  std::set<size_t> addCandIdsCurrentStrip;
279 
280  size_t idxSeed = 0;
281  while (idxSeed < seedCands.size()) {
282  if (verbosity_ >= 2)
283  edm::LogPrint("RecoTauPiZeroStripPlugin3") << "processing seed #" << idxSeed;
284 
285  seedCandIdsCurrentStrip.clear();
286  addCandIdsCurrentStrip.clear();
287 
288  std::unique_ptr<RecoTauPiZero> strip(new RecoTauPiZero(*seedCands[idxSeed], RecoTauPiZero::kStrips));
289  strip->addDaughter(seedCands[idxSeed]);
290  seedCandIdsCurrentStrip.insert(idxSeed);
291 
292  bool isCandAdded;
293  int stripBuildIteration = 0;
294  do {
295  isCandAdded = false;
296 
297  //if ( verbosity_ >= 2 ) edm::LogPrint("RecoTauPiZeroStripPlugin3") << " adding seedCands to strip..." ;
298  addCandsToStrip(*strip, seedCands, seedCandFlags, seedCandIdsCurrentStrip, isCandAdded);
299  //if ( verbosity_ >= 2 ) edm::LogPrint("RecoTauPiZeroStripPlugin3") << " adding addCands to strip..." ;
300  addCandsToStrip(*strip, addCands, addCandFlags, addCandIdsCurrentStrip, isCandAdded);
301 
303  p4Builder_.set(*strip);
304 
305  ++stripBuildIteration;
306  } while (isCandAdded && (stripBuildIteration < maxStripBuildIterations_ || maxStripBuildIterations_ == -1));
307 
308  if (strip->et() > minStripEt_) { // strip passed Et cuts, add it to the event
309  if (verbosity_ >= 2)
310  edm::LogPrint("RecoTauPiZeroStripPlugin3")
311  << "Building strip: Et = " << strip->et() << ", eta = " << strip->eta() << ", phi = " << strip->phi();
312 
313  // Update the vertex
314  if (strip->daughterPtr(0).isNonnull())
315  strip->setVertex(strip->daughterPtr(0)->vertex());
316  output.push_back(std::move(strip));
317 
318  // Mark daughters as being part of this strip
319  markCandsInStrip(seedCandFlags, seedCandIdsCurrentStrip);
320  markCandsInStrip(addCandFlags, addCandIdsCurrentStrip);
321  } else { // strip failed Et cuts, just skip it
322  if (verbosity_ >= 2)
323  edm::LogPrint("RecoTauPiZeroStripPlugin3")
324  << "Discarding strip: Et = " << strip->et() << ", eta = " << strip->eta() << ", phi = " << strip->phi();
325  }
326 
327  ++idxSeed;
328  while (idxSeed < seedCands.size() && seedCandFlags[idxSeed]) {
329  ++idxSeed; // fast-forward to next seed cand not yet included in any strip
330  }
331  }
332 
333  // Check if we want to combine our strips
334  if (combineStrips_ && output.size() > 1) {
335  PiZeroVector stripCombinations;
336  // Sort the output by descending pt
337  output.sort(
338  output.begin(), output.end(), boost::bind(&RecoTauPiZero::pt, _1) > boost::bind(&RecoTauPiZero::pt, _2));
339  // Get the end of interesting set of strips to try and combine
340  PiZeroVector::const_iterator end_iter = takeNElements(output.begin(), output.end(), maxStrips_);
341 
342  // Look at all the combinations
343  for (PiZeroVector::const_iterator first = output.begin(); first != end_iter - 1; ++first) {
344  for (PiZeroVector::const_iterator second = first + 1; second != end_iter; ++second) {
345  Candidate::LorentzVector firstP4 = first->p4();
346  Candidate::LorentzVector secondP4 = second->p4();
347  // If we assume a certain mass for each strip apply it here.
348  firstP4 = applyMassConstraint(firstP4, combinatoricStripMassHypo_);
349  secondP4 = applyMassConstraint(secondP4, combinatoricStripMassHypo_);
350  Candidate::LorentzVector totalP4 = firstP4 + secondP4;
351  // Make our new combined strip
352  std::unique_ptr<RecoTauPiZero> combinedStrips(
353  new RecoTauPiZero(0,
354  totalP4,
355  Candidate::Point(0, 0, 0),
356  //111, 10001, true, RecoTauPiZero::kCombinatoricStrips));
357  111,
358  10001,
359  true,
361 
362  // Now loop over the strip members
363  for (auto const& gamma : first->daughterPtrVector()) {
364  combinedStrips->addDaughter(gamma);
365  }
366  for (auto const& gamma : second->daughterPtrVector()) {
367  combinedStrips->addDaughter(gamma);
368  }
369  // Update the vertex
370  if (combinedStrips->daughterPtr(0).isNonnull()) {
371  combinedStrips->setVertex(combinedStrips->daughterPtr(0)->vertex());
372  }
373 
374  // Add to our collection of combined strips
375  stripCombinations.push_back(std::move(combinedStrips));
376  }
377  }
378  // When done doing all the combinations, add the combined strips to the
379  // output.
380  output.transfer(output.end(), stripCombinations);
381  }
382 
383  // Compute correction to account for spread of photon energy in eta and phi
384  // in case charged pions make nuclear interactions or photons convert within the tracking detector
385  for (PiZeroVector::iterator strip = output.begin(); strip != output.end(); ++strip) {
386  double bendCorrEta = 0.;
387  double bendCorrPhi = 0.;
388  double energySum = 0.;
389  for (auto const& gamma : strip->daughterPtrVector()) {
390  bendCorrEta += (gamma->energy() * etaAssociationDistance_->Eval(gamma->pt()));
391  bendCorrPhi += (gamma->energy() * phiAssociationDistance_->Eval(gamma->pt()));
392  energySum += gamma->energy();
393  }
394  if (energySum > 1.e-2) {
395  bendCorrEta /= energySum;
396  bendCorrPhi /= energySum;
397  }
398  //std::cout << "stripPt = " << strip->pt() << ": bendCorrEta = " << bendCorrEta << ", bendCorrPhi = " << bendCorrPhi << std::endl;
399  strip->setBendCorrEta(bendCorrEta);
400  strip->setBendCorrPhi(bendCorrPhi);
401  }
402 
403  return output.release();
404  }
405  } // namespace tau
406 } // namespace reco
407 
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:26
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:238
double eta() const final
momentum pseudorapidity
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:572
Base class for all types of Jets.
Definition: Jet.h:20
#define nullptr
double phi() const
azimuthal angle of momentum vector
Definition: TrackBase.h:614
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:301
reco::TrackRef trackRef() const
Definition: PFCandidate.cc:408
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:617
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:602
double ptError() const
error on Pt (set to 1000 TeV if charge==0 for safety)
Definition: TrackBase.h:696
void addCandsToStrip(RecoTauPiZero &, CandPtrs &, const std::vector< bool > &, std::set< size_t > &, bool &) const
std::vector< reco::CandidatePtr > CandPtrs
T min(T a, T b)
Definition: MathUtil.h:58
void beginEvent() override
Hook called at the beginning of the event.
std::unique_ptr< const TFormula > etaAssociationDistance_
std::unique_ptr< RecoTauQualityCuts > qcuts_
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:596
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:483
ParameterSet const & getParameterSet(std::string const &) const
ProductID id() const
Accessor for product ID.
Definition: Ptr.h:158
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:37
int numberOfValidTrackerHits() const
Definition: HitPattern.h:789
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:40
fixed size matrix
RecoTauPiZeroStripPlugin3(const edm::ParameterSet &, edm::ConsumesCollector &&iC)
HLT enums.
reco::VertexRef associatedVertex(const Jet &jet) const
boost::ptr_vector< RecoTauPiZero > PiZeroVector
reco::GsfTrackRef gsfTrackRef() const
Definition: PFCandidate.cc:440
std::unique_ptr< const TFormula > phiAssociationDistance_
void addDaughter(const CandidatePtr &)
add a daughter via a reference
int numberOfValidPixelHits() const
Definition: HitPattern.h:801
math::XYZPoint Point
point in the space
Definition: Candidate.h:41
int charge() const
track electric charge
Definition: TrackBase.h:575
#define DEFINE_EDM_PLUGIN(factory, type, name)
void set(reco::Candidate &c) const
set up a candidate
double energySum(const DataFrame &df, int fs, int ls)
return_type operator()(const reco::Jet &) const override
Build a collection of piZeros from objects in the input jet.
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:587
double phi() const final
momentum azimuthal angle
def move(src, dest)
Definition: eostools.py:511
std::unique_ptr< PiZeroVector > return_type