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 boost::ptr_vector< RecoTauPiZeroPiZeroVector
 
typedef std::auto_ptr< PiZeroVectorreturn_type
 

Detailed Description

Definition at line 56 of file RecoTauPiZeroStripPlugin3.cc.

Member Typedef Documentation

Definition at line 67 of file RecoTauPiZeroStripPlugin3.cc.

Constructor & Destructor Documentation

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

Definition at line 113 of file RecoTauPiZeroStripPlugin3.cc.

References applyElecTrackQcuts_, combinatoricStripMassHypo_, combineStrips_, etaAssociationDistance_, edm::ParameterSet::getParameter(), edm::ParameterSet::getParameterSet(), inputParticleIds_, maxStripBuildIterations_, maxStrips_, min(), minGammaEtStripAdd_, minGammaEtStripSeed_, minStripEt_, phiAssociationDistance_, qcuts_, updateStripAfterEachDaughter_, and verbosity_.

116 {
117  minGammaEtStripSeed_ = pset.getParameter<double>("minGammaEtStripSeed");
118  minGammaEtStripAdd_ = pset.getParameter<double>("minGammaEtStripAdd");
119 
120  minStripEt_ = pset.getParameter<double>("minStripEt");
121 
122  edm::ParameterSet qcuts_pset = pset.getParameterSet("qualityCuts").getParameterSet("signalQualityCuts");
123 //-------------------------------------------------------------------------------
124 // CV: disable track quality cuts for PFElectronsPFElectron
125 // (treat PFElectrons like PFGammas for the purpose of building eta-phi strips)
126  applyElecTrackQcuts_ = pset.getParameter<bool>("applyElecTrackQcuts");
127  if ( !applyElecTrackQcuts_ ) {
128  qcuts_pset.addParameter<double>("minTrackPt", std::min(minGammaEtStripSeed_, minGammaEtStripAdd_));
129  qcuts_pset.addParameter<double>("maxTrackChi2", 1.e+9);
130  qcuts_pset.addParameter<double>("maxTransverseImpactParameter", 1.e+9);
131  qcuts_pset.addParameter<double>("maxDeltaZ", 1.e+9);
132  qcuts_pset.addParameter<double>("minTrackVertexWeight", -1.);
133  qcuts_pset.addParameter<unsigned>("minTrackPixelHits", 0);
134  qcuts_pset.addParameter<unsigned>("minTrackHits", 0);
135  }
136 //-------------------------------------------------------------------------------
137  qcuts_pset.addParameter<double>("minGammaEt", std::min(minGammaEtStripSeed_, minGammaEtStripAdd_));
138  //qcuts_ = new RecoTauQualityCuts(qcuts_pset);
139  //std::unique_ptr<RecoTauQualityCuts> qcuts_(new RecoTauQualityCuts(qcuts_pset));
140 
141  qcuts_.reset(new RecoTauQualityCuts(qcuts_pset));
142 
143  inputParticleIds_ = pset.getParameter<std::vector<int> >("stripCandidatesParticleIds");
144  const edm::ParameterSet& stripSize_eta_pset = pset.getParameterSet("stripEtaAssociationDistanceFunc");
145  etaAssociationDistance_ = makeFunction("etaAssociationDistance", stripSize_eta_pset);
146  const edm::ParameterSet& stripSize_phi_pset = pset.getParameterSet("stripPhiAssociationDistanceFunc");
147  phiAssociationDistance_ = makeFunction("phiAssociationDistance", stripSize_phi_pset);
148 
149  updateStripAfterEachDaughter_ = pset.getParameter<bool>("updateStripAfterEachDaughter");
150  maxStripBuildIterations_ = pset.getParameter<int>("maxStripBuildIterations");
151 
152  combineStrips_ = pset.getParameter<bool>("makeCombinatoricStrips");
153  if ( combineStrips_ ) {
154  maxStrips_ = pset.getParameter<int>("maxInputStrips");
155  combinatoricStripMassHypo_ = pset.getParameter<double>("stripMassWhenCombining");
156  }
157 
158  verbosity_ = pset.getParameter<int>("verbosity");
159 }
T getParameter(std::string const &) const
RecoTauPiZeroBuilderPlugin(const edm::ParameterSet &pset, edm::ConsumesCollector &&iC)
T min(T a, T b)
Definition: MathUtil.h:58
std::unique_ptr< const TFormula > etaAssociationDistance_
std::unique_ptr< RecoTauQualityCuts > qcuts_
ParameterSet const & getParameterSet(std::string const &) const
std::unique_ptr< const TFormula > phiAssociationDistance_
def move(src, dest)
Definition: eostools.py:511
reco::tau::RecoTauPiZeroStripPlugin3::~RecoTauPiZeroStripPlugin3 ( )
override

Definition at line 160 of file RecoTauPiZeroStripPlugin3.cc.

161 {
162 }

Member Function Documentation

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 170 of file RecoTauPiZeroStripPlugin3.cc.

References reco::CompositePtrCandidate::addDaughter(), reco::deltaPhi(), reco::LeafCandidate::eta(), etaAssociationDistance_, getTrack(), reco::PFCandidate::gsfTrackRef(), edm::Ptr< T >::id(), edm::Ref< C, T, F >::isNonnull(), edm::Ptr< T >::key(), p4Builder_, reco::LeafCandidate::phi(), phiAssociationDistance_, reco::LeafCandidate::pt(), AddFourMomenta::set(), reco::PFCandidate::trackRef(), updateStripAfterEachDaughter_, and verbosity_.

Referenced by operator()().

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

Hook called at the beginning of the event.

Reimplemented from reco::tau::RecoTauPiZeroBuilderPlugin.

Definition at line 165 of file RecoTauPiZeroStripPlugin3.cc.

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

166 {
168 }
void setEvent(const edm::Event &evt)
Load the vertices from the event.
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 222 of file RecoTauPiZeroStripPlugin3.cc.

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

223 {
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();
242  cand != candsVector.end(); ++cand ) {
243  if ( verbosity_ >= 1 ) {
244  edm::LogPrint("RecoTauPiZeroStripPlugin3") << "PFGamma #" << idx << " (" << cand->id() << ":" << cand->key() << "): Et = " << (*cand)->et() << ", eta = " << (*cand)->eta() << ", phi = " << (*cand)->phi() ;
245  }
246  if ( (*cand)->et() > minGammaEtStripSeed_ ) {
247  if ( verbosity_ >= 2 ) {
248  edm::LogPrint("RecoTauPiZeroStripPlugin3") << "--> assigning seedCandId = " << seedCands.size() ;
249  const reco::TrackBaseRef candTrack = getTrack(**cand);
250  if ( candTrack.isNonnull() ) {
251  edm::LogPrint("RecoTauPiZeroStripPlugin3") << "track: Pt = " << candTrack->pt() << " eta = " << candTrack->eta() << ", phi = " << candTrack->phi() << ", charge = " << candTrack->charge() ;
252  edm::LogPrint("RecoTauPiZeroStripPlugin3") << " (dZ = " << candTrack->dz(vertexAssociator_.associatedVertex(jet)->position()) << ", dXY = " << candTrack->dxy(vertexAssociator_.associatedVertex(jet)->position()) << ","
253  << " numHits = " << candTrack->hitPattern().numberOfValidTrackerHits() << ", numPxlHits = " << candTrack->hitPattern().numberOfValidPixelHits() << ","
254  << " chi2 = " << candTrack->normalizedChi2() << ", dPt/Pt = " << (candTrack->ptError()/candTrack->pt()) << ")" ;
255  }
256  }
257  seedCands.push_back(*cand);
258  } else if ( (*cand)->et() > minGammaEtStripAdd_ ) {
259  if ( verbosity_ >= 2 ) {
260  edm::LogPrint("RecoTauPiZeroStripPlugin3") << "--> assigning addCandId = " << addCands.size() ;
261  }
262  addCands.push_back(*cand);
263  }
264  ++idx;
265  }
266 
267  std::vector<bool> seedCandFlags(seedCands.size()); // true/false: seedCand is already/not yet included in strip
268  std::vector<bool> addCandFlags(addCands.size()); // true/false: addCand is already/not yet included in strip
269 
270  std::set<size_t> seedCandIdsCurrentStrip;
271  std::set<size_t> addCandIdsCurrentStrip;
272 
273  size_t idxSeed = 0;
274  while ( idxSeed < seedCands.size() ) {
275  if ( verbosity_ >= 2 ) edm::LogPrint("RecoTauPiZeroStripPlugin3") << "processing seed #" << idxSeed ;
276 
277  seedCandIdsCurrentStrip.clear();
278  addCandIdsCurrentStrip.clear();
279 
280  std::auto_ptr<RecoTauPiZero> strip(new RecoTauPiZero(*seedCands[idxSeed], RecoTauPiZero::kStrips));
281  strip->addDaughter(seedCands[idxSeed]);
282  seedCandIdsCurrentStrip.insert(idxSeed);
283 
284  bool isCandAdded;
285  int stripBuildIteration = 0;
286  do {
287  isCandAdded = false;
288 
289  //if ( verbosity_ >= 2 ) edm::LogPrint("RecoTauPiZeroStripPlugin3") << " adding seedCands to strip..." ;
290  addCandsToStrip(*strip, seedCands, seedCandFlags, seedCandIdsCurrentStrip, isCandAdded);
291  //if ( verbosity_ >= 2 ) edm::LogPrint("RecoTauPiZeroStripPlugin3") << " adding addCands to strip..." ;
292  addCandsToStrip(*strip, addCands, addCandFlags, addCandIdsCurrentStrip, isCandAdded);
293 
295 
296  ++stripBuildIteration;
297  } while ( isCandAdded && (stripBuildIteration < maxStripBuildIterations_ || maxStripBuildIterations_ == -1) );
298 
299  if ( strip->et() > minStripEt_ ) { // strip passed Et cuts, add it to the event
300  if ( verbosity_ >= 2 ) edm::LogPrint("RecoTauPiZeroStripPlugin3") << "Building strip: Et = " << strip->et() << ", eta = " << strip->eta() << ", phi = " << strip->phi() ;
301 
302  // Update the vertex
303  if ( strip->daughterPtr(0).isNonnull() ) strip->setVertex(strip->daughterPtr(0)->vertex());
304  output.push_back(strip);
305 
306  // Mark daughters as being part of this strip
307  markCandsInStrip(seedCandFlags, seedCandIdsCurrentStrip);
308  markCandsInStrip(addCandFlags, addCandIdsCurrentStrip);
309  } else { // strip failed Et cuts, just skip it
310  if ( verbosity_ >= 2 ) edm::LogPrint("RecoTauPiZeroStripPlugin3") << "Discarding strip: Et = " << strip->et() << ", eta = " << strip->eta() << ", phi = " << strip->phi() ;
311  }
312 
313  ++idxSeed;
314  while ( idxSeed < seedCands.size() && seedCandFlags[idxSeed] ) {
315  ++idxSeed; // fast-forward to next seed cand not yet included in any strip
316  }
317  }
318 
319  // Check if we want to combine our strips
320  if ( combineStrips_ && output.size() > 1 ) {
321  PiZeroVector stripCombinations;
322  // Sort the output by descending pt
323  output.sort(output.begin(), output.end(),
324  boost::bind(&RecoTauPiZero::pt, _1) >
325  boost::bind(&RecoTauPiZero::pt, _2));
326  // Get the end of interesting set of strips to try and combine
327  PiZeroVector::const_iterator end_iter = takeNElements(
328  output.begin(), output.end(), maxStrips_);
329 
330  // Look at all the combinations
331  for ( PiZeroVector::const_iterator first = output.begin();
332  first != end_iter-1; ++first ) {
333  for ( PiZeroVector::const_iterator second = first+1;
334  second != end_iter; ++second ) {
335  Candidate::LorentzVector firstP4 = first->p4();
336  Candidate::LorentzVector secondP4 = second->p4();
337  // If we assume a certain mass for each strip apply it here.
338  firstP4 = applyMassConstraint(firstP4, combinatoricStripMassHypo_);
339  secondP4 = applyMassConstraint(secondP4, combinatoricStripMassHypo_);
340  Candidate::LorentzVector totalP4 = firstP4 + secondP4;
341  // Make our new combined strip
342  std::auto_ptr<RecoTauPiZero> combinedStrips(
343  new RecoTauPiZero(0, totalP4,
344  Candidate::Point(0, 0, 0),
345  //111, 10001, true, RecoTauPiZero::kCombinatoricStrips));
346  111, 10001, true, RecoTauPiZero::kUndefined));
347 
348  // Now loop over the strip members
349  for ( auto const& gamma : first->daughterPtrVector()) {
350  combinedStrips->addDaughter(gamma);
351  }
352  for ( auto const& gamma : second->daughterPtrVector()) {
353  combinedStrips->addDaughter(gamma);
354  }
355  // Update the vertex
356  if ( combinedStrips->daughterPtr(0).isNonnull() ) {
357  combinedStrips->setVertex(combinedStrips->daughterPtr(0)->vertex());
358  }
359 
360  // Add to our collection of combined strips
361  stripCombinations.push_back(combinedStrips);
362  }
363  }
364  // When done doing all the combinations, add the combined strips to the
365  // output.
366  output.transfer(output.end(), stripCombinations);
367  }
368 
369  // Compute correction to account for spread of photon energy in eta and phi
370  // in case charged pions make nuclear interactions or photons convert within the tracking detector
371  for ( PiZeroVector::iterator strip = output.begin();
372  strip != output.end(); ++strip ) {
373  double bendCorrEta = 0.;
374  double bendCorrPhi = 0.;
375  double energySum = 0.;
376  for (auto const& gamma : strip->daughterPtrVector()) {
377  bendCorrEta += (gamma->energy()*etaAssociationDistance_->Eval(gamma->pt()));
378  bendCorrPhi += (gamma->energy()*phiAssociationDistance_->Eval(gamma->pt()));
379  energySum += gamma->energy();
380  }
381  if ( energySum > 1.e-2 ) {
382  bendCorrEta /= energySum;
383  bendCorrPhi /= energySum;
384  }
385  //std::cout << "stripPt = " << strip->pt() << ": bendCorrEta = " << bendCorrEta << ", bendCorrPhi = " << bendCorrPhi << std::endl;
386  strip->setBendCorrEta(bendCorrEta);
387  strip->setBendCorrPhi(bendCorrPhi);
388  }
389 
390  return output.release();
391 }
std::vector< CandidatePtr > pfCandidates(const Jet &jet, int particleId, bool sort=true)
InputIterator takeNElements(const InputIterator &begin, const InputIterator &end, size_t N)
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:594
double phi() const
azimuthal angle of momentum vector
Definition: TrackBase.h:678
double pt() const final
transverse momentum
bool isNonnull() const
Checks for non-null.
Definition: RefToBase.h:340
U second(std::pair< T, U > const &p)
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:684
std::vector< reco::CandidatePtr > CandPtrs
double pt() const
track transverse momentum
Definition: TrackBase.h:654
double ptError() const
error on Pt (set to 1000 TeV if charge==0 for safety)
Definition: TrackBase.h:808
void addCandsToStrip(RecoTauPiZero &, CandPtrs &, const std::vector< bool > &, std::set< size_t > &, bool &) const
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:642
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:479
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:37
int numberOfValidTrackerHits() const
Definition: HitPattern.h:901
reco::VertexRef associatedVertex(const Jet &jet) const
boost::ptr_vector< RecoTauPiZero > PiZeroVector
std::unique_ptr< const TFormula > phiAssociationDistance_
int numberOfValidPixelHits() const
Definition: HitPattern.h:916
math::XYZPoint Point
point in the space
Definition: Candidate.h:41
int charge() const
track electric charge
Definition: TrackBase.h:600
void set(reco::Candidate &c) const
set up a candidate
double energySum(const DataFrame &df, int fs, int ls)
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:624

Member Data Documentation

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

Definition at line 73 of file RecoTauPiZeroStripPlugin3.cc.

Referenced by RecoTauPiZeroStripPlugin3().

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

Definition at line 89 of file RecoTauPiZeroStripPlugin3.cc.

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

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

Definition at line 87 of file RecoTauPiZeroStripPlugin3.cc.

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

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

Definition at line 79 of file RecoTauPiZeroStripPlugin3.cc.

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

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

Definition at line 84 of file RecoTauPiZeroStripPlugin3.cc.

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

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

Definition at line 88 of file RecoTauPiZeroStripPlugin3.cc.

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

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

Definition at line 75 of file RecoTauPiZeroStripPlugin3.cc.

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

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

Definition at line 74 of file RecoTauPiZeroStripPlugin3.cc.

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

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

Definition at line 77 of file RecoTauPiZeroStripPlugin3.cc.

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

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

Definition at line 91 of file RecoTauPiZeroStripPlugin3.cc.

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

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

Definition at line 72 of file RecoTauPiZeroStripPlugin3.cc.

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

bool reco::tau::RecoTauPiZeroStripPlugin3::updateStripAfterEachDaughter_
private
int reco::tau::RecoTauPiZeroStripPlugin3::verbosity_
private
RecoTauVertexAssociator reco::tau::RecoTauPiZeroStripPlugin3::vertexAssociator_
private

Definition at line 70 of file RecoTauPiZeroStripPlugin3.cc.

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