CMS 3D CMS Logo

List of all members | Public Member Functions | Private Types | Private Member Functions | Private Attributes
reco::tau::RecoTauPiZeroStripPlugin3 Class Reference
Inheritance diagram for reco::tau::RecoTauPiZeroStripPlugin3:
reco::tau::RecoTauPiZeroBuilderPlugin reco::tau::RecoTauEventHolderPlugin reco::tau::RecoTauNamedPlugin

Public Member Functions

void beginEvent () override
 Hook called at the beginning of the event. More...
 
return_type operator() (const reco::Jet &) const override
 Build a collection of piZeros from objects in the input jet. More...
 
 RecoTauPiZeroStripPlugin3 (const edm::ParameterSet &, edm::ConsumesCollector &&iC)
 
 ~RecoTauPiZeroStripPlugin3 () override
 
- Public Member Functions inherited from reco::tau::RecoTauPiZeroBuilderPlugin
 RecoTauPiZeroBuilderPlugin (const edm::ParameterSet &pset, edm::ConsumesCollector &&iC)
 
 ~RecoTauPiZeroBuilderPlugin () override
 
- Public Member Functions inherited from reco::tau::RecoTauEventHolderPlugin
const edm::Eventevt () const
 
edm::Eventevt ()
 
const edm::EventSetupevtSetup () const
 
 RecoTauEventHolderPlugin (const edm::ParameterSet &pset)
 
void setup (edm::Event &, const edm::EventSetup &)
 
 ~RecoTauEventHolderPlugin () override
 
- Public Member Functions inherited from reco::tau::RecoTauNamedPlugin
const std::string & name () const
 
 RecoTauNamedPlugin (const edm::ParameterSet &pset)
 
virtual ~RecoTauNamedPlugin ()
 

Private Types

typedef std::vector< reco::CandidatePtrCandPtrs
 

Private Member Functions

void addCandsToStrip (RecoTauPiZero &, CandPtrs &, const std::vector< bool > &, std::set< size_t > &, bool &) const
 

Private Attributes

bool applyElecTrackQcuts_
 
double combinatoricStripMassHypo_
 
bool combineStrips_
 
std::unique_ptr< const TFormula > etaAssociationDistance_
 
std::vector< int > inputParticleIds_
 
int maxStripBuildIterations_
 
int maxStrips_
 
double minGammaEtStripAdd_
 
double minGammaEtStripSeed_
 
double minStripEt_
 
AddFourMomenta p4Builder_
 
std::unique_ptr< const TFormula > phiAssociationDistance_
 
std::unique_ptr< RecoTauQualityCutsqcuts_
 
bool updateStripAfterEachDaughter_
 
int verbosity_
 
RecoTauVertexAssociator vertexAssociator_
 

Additional Inherited Members

- Public Types inherited from reco::tau::RecoTauPiZeroBuilderPlugin
typedef std::vector< std::unique_ptr< RecoTauPiZero > > PiZeroVector
 
typedef PiZeroVector return_type
 

Detailed Description

Definition at line 55 of file RecoTauPiZeroStripPlugin3.cc.

Member Typedef Documentation

◆ CandPtrs

Definition at line 65 of file RecoTauPiZeroStripPlugin3.cc.

Constructor & Destructor Documentation

◆ RecoTauPiZeroStripPlugin3()

reco::tau::RecoTauPiZeroStripPlugin3::RecoTauPiZeroStripPlugin3 ( const edm::ParameterSet pset,
edm::ConsumesCollector &&  iC 
)
explicit

Definition at line 109 of file RecoTauPiZeroStripPlugin3.cc.

References applyElecTrackQcuts_, combinatoricStripMassHypo_, combineStrips_, etaAssociationDistance_, inputParticleIds_, maxStripBuildIterations_, maxStrips_, SiStripPI::min, minGammaEtStripAdd_, minGammaEtStripSeed_, minStripEt_, phiAssociationDistance_, muonDTDigis_cfi::pset, qcuts_, updateStripAfterEachDaughter_, and verbosity_.

111  vertexAssociator_(pset.getParameter<edm::ParameterSet>("qualityCuts"), std::move(iC)),
112  qcuts_(nullptr),
113  etaAssociationDistance_(nullptr),
114  phiAssociationDistance_(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_ = std::make_unique<RecoTauQualityCuts>(qcuts_pset);
137 
138  inputParticleIds_ = pset.getParameter<std::vector<int> >("stripCandidatesParticleIds");
139  const edm::ParameterSet& stripSize_eta_pset = pset.getParameterSet("stripEtaAssociationDistanceFunc");
140  etaAssociationDistance_ = makeFunction("etaAssociationDistance", stripSize_eta_pset);
141  const edm::ParameterSet& stripSize_phi_pset = pset.getParameterSet("stripPhiAssociationDistanceFunc");
142  phiAssociationDistance_ = makeFunction("phiAssociationDistance", stripSize_phi_pset);
143 
144  updateStripAfterEachDaughter_ = pset.getParameter<bool>("updateStripAfterEachDaughter");
145  maxStripBuildIterations_ = pset.getParameter<int>("maxStripBuildIterations");
146 
147  combineStrips_ = pset.getParameter<bool>("makeCombinatoricStrips");
148  if (combineStrips_) {
149  maxStrips_ = pset.getParameter<int>("maxInputStrips");
150  combinatoricStripMassHypo_ = pset.getParameter<double>("stripMassWhenCombining");
151  }
152 
153  verbosity_ = pset.getParameter<int>("verbosity");
154  }
RecoTauPiZeroBuilderPlugin(const edm::ParameterSet &pset, edm::ConsumesCollector &&iC)
std::unique_ptr< const TFormula > etaAssociationDistance_
std::unique_ptr< RecoTauQualityCuts > qcuts_
std::unique_ptr< const TFormula > phiAssociationDistance_
def move(src, dest)
Definition: eostools.py:511

◆ ~RecoTauPiZeroStripPlugin3()

reco::tau::RecoTauPiZeroStripPlugin3::~RecoTauPiZeroStripPlugin3 ( )
override

Definition at line 155 of file RecoTauPiZeroStripPlugin3.cc.

155 {}

Member Function Documentation

◆ addCandsToStrip()

void reco::tau::RecoTauPiZeroStripPlugin3::addCandsToStrip ( RecoTauPiZero strip,
CandPtrs cands,
const std::vector< bool > &  candFlags,
std::set< size_t > &  candIdsCurrentStrip,
bool &  isCandAdded 
) const
private

Definition at line 160 of file RecoTauPiZeroStripPlugin3.cc.

References HLT_2024v12_cff::cands, reco::deltaPhi(), etaAssociationDistance_, p4Builder_, phiAssociationDistance_, AddFourMomenta::set(), nano_mu_digi_cff::strip, updateStripAfterEachDaughter_, and verbosity_.

Referenced by operator()().

164  {
165  if (verbosity_ >= 1) {
166  edm::LogPrint("RecoTauPiZeroStripPlugin3") << "<RecoTauPiZeroStripPlugin3::addCandsToStrip>:";
167  }
168  size_t numCands = cands.size();
169  for (size_t candId = 0; candId < numCands; ++candId) {
170  if ((!candFlags[candId]) &&
171  candIdsCurrentStrip.find(candId) == candIdsCurrentStrip.end()) { // do not include same cand twice
172  reco::CandidatePtr cand = cands[candId];
173  double etaAssociationDistance_value =
174  etaAssociationDistance_->Eval(strip.pt()) + etaAssociationDistance_->Eval(cand->pt());
175  double phiAssociationDistance_value =
176  phiAssociationDistance_->Eval(strip.pt()) + phiAssociationDistance_->Eval(cand->pt());
177  if (fabs(strip.eta() - cand->eta()) <
178  etaAssociationDistance_value && // check if cand is within eta-phi window centered on strip
179  reco::deltaPhi(strip.phi(), cand->phi()) < phiAssociationDistance_value) {
180  if (verbosity_ >= 2) {
181  edm::LogPrint("RecoTauPiZeroStripPlugin3")
182  << "--> adding PFCand #" << candId << " (" << cand.id() << ":" << cand.key()
183  << "): Et = " << cand->et() << ", eta = " << cand->eta() << ", phi = " << cand->phi();
184  }
185  strip.addDaughter(cand);
188  isCandAdded = true;
189  candIdsCurrentStrip.insert(candId);
190  }
191  }
192  }
193  }
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:26
void set(reco::Candidate &c) const
set up a candidate
std::unique_ptr< const TFormula > etaAssociationDistance_
Log< level::Warning, true > LogPrint
std::unique_ptr< const TFormula > phiAssociationDistance_

◆ beginEvent()

void reco::tau::RecoTauPiZeroStripPlugin3::beginEvent ( )
overridevirtual

Hook called at the beginning of the event.

Reimplemented from reco::tau::RecoTauPiZeroBuilderPlugin.

Definition at line 158 of file RecoTauPiZeroStripPlugin3.cc.

References reco::tau::RecoTauEventHolderPlugin::evt(), reco::tau::RecoTauVertexAssociator::setEvent(), and vertexAssociator_.

void setEvent(const edm::Event &evt)
Load the vertices from the event.

◆ operator()()

RecoTauPiZeroStripPlugin3::return_type reco::tau::RecoTauPiZeroStripPlugin3::operator() ( const reco::Jet ) const
overridevirtual

Build a collection of piZeros from objects in the input jet.

Implements reco::tau::RecoTauPiZeroBuilderPlugin.

Definition at line 220 of file RecoTauPiZeroStripPlugin3.cc.

References addCandsToStrip(), reco::tau::RecoTauVertexAssociator::associatedVertex(), reco::TrackBase::charge(), combinatoricStripMassHypo_, combineStrips_, reco::TrackBase::dxy(), reco::TrackBase::dz(), MillePedeFileConverter_cfg::e, CastorDataFrameFilter_impl::energySum(), reco::TrackBase::eta(), etaAssociationDistance_, dqmdumpme::first, CustomPhysics_cfi::gamma, getTrack(), reco::TrackBase::hitPattern(), heavyIonCSV_trainingSettings::idx, inputParticleIds_, edm::RefToBase< T >::isNonnull(), metsig::jet, reco::RecoTauPiZero::kStrips, reco::RecoTauPiZero::kUndefined, maxStripBuildIterations_, maxStrips_, minGammaEtStripAdd_, minGammaEtStripSeed_, minStripEt_, eostools::move(), reco::TrackBase::normalizedChi2(), reco::HitPattern::numberOfValidPixelHits(), reco::HitPattern::numberOfValidTrackerHits(), convertSQLitetoXML_cfg::output, p4Builder_, reco::tau::pfCandidates(), reco::TrackBase::phi(), phiAssociationDistance_, reco::TrackBase::pt(), reco::TrackBase::ptError(), qcuts_, edm::second(), AddFourMomenta::set(), jetUpdater_cfi::sort, nano_mu_digi_cff::strip, reco::tau::takeNElements(), updateStripAfterEachDaughter_, verbosity_, and vertexAssociator_.

220  {
221  if (verbosity_ >= 1) {
222  edm::LogPrint("RecoTauPiZeroStripPlugin3") << "<RecoTauPiZeroStripPlugin3::operator()>:";
223  edm::LogPrint("RecoTauPiZeroStripPlugin3") << " minGammaEtStripSeed = " << minGammaEtStripSeed_;
224  edm::LogPrint("RecoTauPiZeroStripPlugin3") << " minGammaEtStripAdd = " << minGammaEtStripAdd_;
225  edm::LogPrint("RecoTauPiZeroStripPlugin3") << " minStripEt = " << minStripEt_;
226  }
227 
229 
230  // Get the candidates passing our quality cuts
232  CandPtrs candsVector = qcuts_->filterCandRefs(pfCandidates(jet, inputParticleIds_));
233 
234  // Convert to stl::list to allow fast deletions
235  CandPtrs seedCands;
236  CandPtrs addCands;
237  int idx = 0;
238  for (CandPtrs::iterator cand = candsVector.begin(); cand != candsVector.end(); ++cand) {
239  if (verbosity_ >= 1) {
240  edm::LogPrint("RecoTauPiZeroStripPlugin3")
241  << "PFGamma #" << idx << " (" << cand->id() << ":" << cand->key() << "): Et = " << (*cand)->et()
242  << ", eta = " << (*cand)->eta() << ", phi = " << (*cand)->phi();
243  }
244  if ((*cand)->et() > minGammaEtStripSeed_) {
245  if (verbosity_ >= 2) {
246  edm::LogPrint("RecoTauPiZeroStripPlugin3") << "--> assigning seedCandId = " << seedCands.size();
247  const reco::TrackBaseRef candTrack = getTrack(**cand);
248  if (candTrack.isNonnull()) {
249  edm::LogPrint("RecoTauPiZeroStripPlugin3")
250  << "track: Pt = " << candTrack->pt() << " eta = " << candTrack->eta()
251  << ", phi = " << candTrack->phi() << ", charge = " << candTrack->charge();
252  edm::LogPrint("RecoTauPiZeroStripPlugin3")
253  << " (dZ = " << candTrack->dz(vertexAssociator_.associatedVertex(jet)->position())
254  << ", dXY = " << candTrack->dxy(vertexAssociator_.associatedVertex(jet)->position()) << ","
255  << " numHits = " << candTrack->hitPattern().numberOfValidTrackerHits()
256  << ", numPxlHits = " << candTrack->hitPattern().numberOfValidPixelHits() << ","
257  << " chi2 = " << candTrack->normalizedChi2()
258  << ", dPt/Pt = " << (candTrack->ptError() / candTrack->pt()) << ")";
259  }
260  }
261  seedCands.push_back(*cand);
262  } else if ((*cand)->et() > minGammaEtStripAdd_) {
263  if (verbosity_ >= 2) {
264  edm::LogPrint("RecoTauPiZeroStripPlugin3") << "--> assigning addCandId = " << addCands.size();
265  }
266  addCands.push_back(*cand);
267  }
268  ++idx;
269  }
270 
271  std::vector<bool> seedCandFlags(seedCands.size()); // true/false: seedCand is already/not yet included in strip
272  std::vector<bool> addCandFlags(addCands.size()); // true/false: addCand is already/not yet included in strip
273 
274  std::set<size_t> seedCandIdsCurrentStrip;
275  std::set<size_t> addCandIdsCurrentStrip;
276 
277  size_t idxSeed = 0;
278  while (idxSeed < seedCands.size()) {
279  if (verbosity_ >= 2)
280  edm::LogPrint("RecoTauPiZeroStripPlugin3") << "processing seed #" << idxSeed;
281 
282  seedCandIdsCurrentStrip.clear();
283  addCandIdsCurrentStrip.clear();
284 
285  auto strip = std::make_unique<RecoTauPiZero>(*seedCands[idxSeed], RecoTauPiZero::kStrips);
286  strip->addDaughter(seedCands[idxSeed]);
287  seedCandIdsCurrentStrip.insert(idxSeed);
288 
289  bool isCandAdded;
290  int stripBuildIteration = 0;
291  do {
292  isCandAdded = false;
293 
294  //if ( verbosity_ >= 2 ) edm::LogPrint("RecoTauPiZeroStripPlugin3") << " adding seedCands to strip..." ;
295  addCandsToStrip(*strip, seedCands, seedCandFlags, seedCandIdsCurrentStrip, isCandAdded);
296  //if ( verbosity_ >= 2 ) edm::LogPrint("RecoTauPiZeroStripPlugin3") << " adding addCands to strip..." ;
297  addCandsToStrip(*strip, addCands, addCandFlags, addCandIdsCurrentStrip, isCandAdded);
298 
300  p4Builder_.set(*strip);
301 
302  ++stripBuildIteration;
303  } while (isCandAdded && (stripBuildIteration < maxStripBuildIterations_ || maxStripBuildIterations_ == -1));
304 
305  if (strip->et() > minStripEt_) { // strip passed Et cuts, add it to the event
306  if (verbosity_ >= 2)
307  edm::LogPrint("RecoTauPiZeroStripPlugin3")
308  << "Building strip: Et = " << strip->et() << ", eta = " << strip->eta() << ", phi = " << strip->phi();
309 
310  // Update the vertex
311  if (strip->daughterPtr(0).isNonnull())
312  strip->setVertex(strip->daughterPtr(0)->vertex());
313  output.push_back(std::move(strip));
314 
315  // Mark daughters as being part of this strip
316  markCandsInStrip(seedCandFlags, seedCandIdsCurrentStrip);
317  markCandsInStrip(addCandFlags, addCandIdsCurrentStrip);
318  } else { // strip failed Et cuts, just skip it
319  if (verbosity_ >= 2)
320  edm::LogPrint("RecoTauPiZeroStripPlugin3")
321  << "Discarding strip: Et = " << strip->et() << ", eta = " << strip->eta() << ", phi = " << strip->phi();
322  }
323 
324  ++idxSeed;
325  while (idxSeed < seedCands.size() && seedCandFlags[idxSeed]) {
326  ++idxSeed; // fast-forward to next seed cand not yet included in any strip
327  }
328  }
329 
330  // Check if we want to combine our strips
331  if (combineStrips_ && output.size() > 1) {
332  PiZeroVector stripCombinations;
333  // Sort the output by descending pt
334  std::sort(output.begin(), output.end(), [&](auto& arg1, auto& arg2) { return arg1->pt() > arg2->pt(); });
335  // Get the end of interesting set of strips to try and combine
336  PiZeroVector::const_iterator end_iter = takeNElements(output.begin(), output.end(), maxStrips_);
337 
338  // Look at all the combinations
339  for (PiZeroVector::const_iterator firstIter = output.begin(); firstIter != end_iter - 1; ++firstIter) {
340  for (PiZeroVector::const_iterator secondIter = firstIter + 1; secondIter != end_iter; ++secondIter) {
341  auto const& first = *firstIter;
342  auto const& second = *secondIter;
343  Candidate::LorentzVector firstP4 = first->p4();
344  Candidate::LorentzVector secondP4 = second->p4();
345  // If we assume a certain mass for each strip apply it here.
346  firstP4 = applyMassConstraint(firstP4, combinatoricStripMassHypo_);
347  secondP4 = applyMassConstraint(secondP4, combinatoricStripMassHypo_);
348  Candidate::LorentzVector totalP4 = firstP4 + secondP4;
349  // Make our new combined strip
350  auto combinedStrips =
351  std::make_unique<RecoTauPiZero>(0,
352  totalP4,
353  Candidate::Point(0, 0, 0),
354  //111, 10001, true, RecoTauPiZero::kCombinatoricStrips));
355  111,
356  10001,
357  true,
359 
360  // Now loop over the strip members
361  for (auto const& gamma : first->daughterPtrVector()) {
362  combinedStrips->addDaughter(gamma);
363  }
364  for (auto const& gamma : second->daughterPtrVector()) {
365  combinedStrips->addDaughter(gamma);
366  }
367  // Update the vertex
368  if (combinedStrips->daughterPtr(0).isNonnull()) {
369  combinedStrips->setVertex(combinedStrips->daughterPtr(0)->vertex());
370  }
371 
372  // Add to our collection of combined strips
373  stripCombinations.push_back(std::move(combinedStrips));
374  }
375  }
376  // When done doing all the combinations, add the combined strips to the
377  // output.
378  std::move(stripCombinations.begin(), stripCombinations.end(), std::back_inserter(output));
379  }
380 
381  // Compute correction to account for spread of photon energy in eta and phi
382  // in case charged pions make nuclear interactions or photons convert within the tracking detector
383  for (auto const& strip : output) {
384  double bendCorrEta = 0.;
385  double bendCorrPhi = 0.;
386  double energySum = 0.;
387  for (auto const& gamma : strip->daughterPtrVector()) {
388  bendCorrEta += (gamma->energy() * etaAssociationDistance_->Eval(gamma->pt()));
389  bendCorrPhi += (gamma->energy() * phiAssociationDistance_->Eval(gamma->pt()));
390  energySum += gamma->energy();
391  }
392  if (energySum > 1.e-2) {
393  bendCorrEta /= energySum;
394  bendCorrPhi /= energySum;
395  }
396  //std::cout << "stripPt = " << strip->pt() << ": bendCorrEta = " << bendCorrEta << ", bendCorrPhi = " << bendCorrPhi << std::endl;
397  strip->setBendCorrEta(bendCorrEta);
398  strip->setBendCorrPhi(bendCorrPhi);
399  }
400 
401  return output;
402  }
int numberOfValidPixelHits() const
Definition: HitPattern.h:831
std::vector< CandidatePtr > pfCandidates(const Jet &jet, int particleId, bool sort=true)
InputIterator takeNElements(const InputIterator &begin, const InputIterator &end, size_t N)
void set(reco::Candidate &c) const
set up a candidate
double ptError() const
error on Pt (set to 1000 TeV if charge==0 for safety)
Definition: TrackBase.h:754
bool isNonnull() const
Checks for non-null.
Definition: RefToBase.h:303
void addCandsToStrip(RecoTauPiZero &, CandPtrs &, const std::vector< bool > &, std::set< size_t > &, bool &) const
reco::VertexRef associatedVertex(const Jet &jet) const
U second(std::pair< T, U > const &p)
double pt() const
track transverse momentum
Definition: TrackBase.h:637
int charge() const
track electric charge
Definition: TrackBase.h:596
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
std::vector< reco::CandidatePtr > CandPtrs
double phi() const
azimuthal angle of momentum vector
Definition: TrackBase.h:649
std::unique_ptr< const TFormula > etaAssociationDistance_
Log< level::Warning, true > LogPrint
std::unique_ptr< RecoTauQualityCuts > qcuts_
static const TrackGhostTrackState * getTrack(const BasicGhostTrackState *basic)
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:652
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:504
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
std::vector< std::unique_ptr< RecoTauPiZero > > PiZeroVector
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:36
std::unique_ptr< const TFormula > phiAssociationDistance_
math::XYZPoint Point
point in the space
Definition: Candidate.h:40
Definition: output.py:1
double energySum(const DataFrame &df, int fs, int ls)
int numberOfValidTrackerHits() const
Definition: HitPattern.h:819
def move(src, dest)
Definition: eostools.py:511
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

Member Data Documentation

◆ applyElecTrackQcuts_

bool reco::tau::RecoTauPiZeroStripPlugin3::applyElecTrackQcuts_
private

Definition at line 71 of file RecoTauPiZeroStripPlugin3.cc.

Referenced by RecoTauPiZeroStripPlugin3().

◆ combinatoricStripMassHypo_

double reco::tau::RecoTauPiZeroStripPlugin3::combinatoricStripMassHypo_
private

Definition at line 87 of file RecoTauPiZeroStripPlugin3.cc.

Referenced by operator()(), and RecoTauPiZeroStripPlugin3().

◆ combineStrips_

bool reco::tau::RecoTauPiZeroStripPlugin3::combineStrips_
private

Definition at line 85 of file RecoTauPiZeroStripPlugin3.cc.

Referenced by operator()(), and RecoTauPiZeroStripPlugin3().

◆ etaAssociationDistance_

std::unique_ptr<const TFormula> reco::tau::RecoTauPiZeroStripPlugin3::etaAssociationDistance_
private

◆ inputParticleIds_

std::vector<int> reco::tau::RecoTauPiZeroStripPlugin3::inputParticleIds_
private

Definition at line 77 of file RecoTauPiZeroStripPlugin3.cc.

Referenced by operator()(), and RecoTauPiZeroStripPlugin3().

◆ maxStripBuildIterations_

int reco::tau::RecoTauPiZeroStripPlugin3::maxStripBuildIterations_
private

Definition at line 82 of file RecoTauPiZeroStripPlugin3.cc.

Referenced by operator()(), and RecoTauPiZeroStripPlugin3().

◆ maxStrips_

int reco::tau::RecoTauPiZeroStripPlugin3::maxStrips_
private

Definition at line 86 of file RecoTauPiZeroStripPlugin3.cc.

Referenced by operator()(), and RecoTauPiZeroStripPlugin3().

◆ minGammaEtStripAdd_

double reco::tau::RecoTauPiZeroStripPlugin3::minGammaEtStripAdd_
private

Definition at line 73 of file RecoTauPiZeroStripPlugin3.cc.

Referenced by operator()(), and RecoTauPiZeroStripPlugin3().

◆ minGammaEtStripSeed_

double reco::tau::RecoTauPiZeroStripPlugin3::minGammaEtStripSeed_
private

Definition at line 72 of file RecoTauPiZeroStripPlugin3.cc.

Referenced by operator()(), and RecoTauPiZeroStripPlugin3().

◆ minStripEt_

double reco::tau::RecoTauPiZeroStripPlugin3::minStripEt_
private

Definition at line 75 of file RecoTauPiZeroStripPlugin3.cc.

Referenced by operator()(), and RecoTauPiZeroStripPlugin3().

◆ p4Builder_

AddFourMomenta reco::tau::RecoTauPiZeroStripPlugin3::p4Builder_
private

Definition at line 89 of file RecoTauPiZeroStripPlugin3.cc.

Referenced by addCandsToStrip(), and operator()().

◆ phiAssociationDistance_

std::unique_ptr<const TFormula> reco::tau::RecoTauPiZeroStripPlugin3::phiAssociationDistance_
private

◆ qcuts_

std::unique_ptr<RecoTauQualityCuts> reco::tau::RecoTauPiZeroStripPlugin3::qcuts_
private

Definition at line 70 of file RecoTauPiZeroStripPlugin3.cc.

Referenced by operator()(), and RecoTauPiZeroStripPlugin3().

◆ updateStripAfterEachDaughter_

bool reco::tau::RecoTauPiZeroStripPlugin3::updateStripAfterEachDaughter_
private

◆ verbosity_

int reco::tau::RecoTauPiZeroStripPlugin3::verbosity_
private

◆ vertexAssociator_

RecoTauVertexAssociator reco::tau::RecoTauPiZeroStripPlugin3::vertexAssociator_
private

Definition at line 68 of file RecoTauPiZeroStripPlugin3.cc.

Referenced by beginEvent(), and operator()().