CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes | Static Private Attributes
reco::tau::PFRecoTauChargedHadronFromGenericTrackPlugin< TrackClass > Class Template Reference

#include <PFRecoTauChargedHadron.h>

Inheritance diagram for reco::tau::PFRecoTauChargedHadronFromGenericTrackPlugin< TrackClass >:
reco::tau::PFRecoTauChargedHadronBuilderPlugin 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 chargedHadrons from objects in the input jet. More...
 
 PFRecoTauChargedHadronFromGenericTrackPlugin (const edm::ParameterSet &, edm::ConsumesCollector &&iC)
 
 ~PFRecoTauChargedHadronFromGenericTrackPlugin () override
 
- Public Member Functions inherited from reco::tau::PFRecoTauChargedHadronBuilderPlugin
 PFRecoTauChargedHadronBuilderPlugin (const edm::ParameterSet &pset, edm::ConsumesCollector &&iC)
 
 ~PFRecoTauChargedHadronBuilderPlugin () 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 Member Functions

bool filterTrack (const edm::Handle< std::vector< TrackClass > > &, size_t iTrack) const
 
template<>
bool filterTrack (const edm::Handle< std::vector< reco::Track > > &tracks, size_t iTrack) const
 
template<>
bool filterTrack (const edm::Handle< std::vector< pat::PackedCandidate > > &tracks, size_t iTrack) const
 
XYZTLorentzVector getTrackPos (const TrackClass &track) const
 
template<>
XYZTLorentzVector getTrackPos (const reco::Track &track) const
 
template<>
XYZTLorentzVector getTrackPos (const pat::PackedCandidate &track) const
 
double getTrackPtError (const TrackClass &track) const
 
template<>
double getTrackPtError (const reco::Track &track) const
 
template<>
double getTrackPtError (const pat::PackedCandidate &cand) const
 
void setChargedHadronTrack (PFRecoTauChargedHadron &chargedHadron, const edm::Ptr< TrackClass > &track) const
 
template<>
void setChargedHadronTrack (PFRecoTauChargedHadron &chargedHadron, const edm::Ptr< reco::Track > &track) const
 
template<>
void setChargedHadronTrack (PFRecoTauChargedHadron &chargedHadron, const edm::Ptr< pat::PackedCandidate > &track) const
 

Private Attributes

double dRcone_
 
bool dRconeLimitedToJetArea_
 
double dRmergeNeutralHadron_
 
double dRmergePhoton_
 
math::XYZVector magneticFieldStrength_
 
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecordmagneticFieldToken_
 
RecoTauQualityCutsqcuts_
 
edm::InputTag srcTracks_
 
edm::EDGetTokenT< std::vector< TrackClass > > Tracks_token
 
int verbosity_
 
RecoTauVertexAssociator vertexAssociator_
 

Static Private Attributes

static constexpr unsigned int maxWarnings_ = 3
 
static std::atomic< unsigned int > numWarnings_ {0}
 

Additional Inherited Members

- Public Types inherited from reco::tau::PFRecoTauChargedHadronBuilderPlugin
typedef std::vector< std::unique_ptr< PFRecoTauChargedHadron > > ChargedHadronVector
 
typedef ChargedHadronVector return_type
 

Detailed Description

template<class TrackClass>
class reco::tau::PFRecoTauChargedHadronFromGenericTrackPlugin< TrackClass >

Definition at line 13 of file PFRecoTauChargedHadron.h.

Constructor & Destructor Documentation

◆ PFRecoTauChargedHadronFromGenericTrackPlugin()

Definition at line 95 of file PFRecoTauChargedHadronFromGenericTrackPlugin.cc.

References reco::tau::PFRecoTauChargedHadronFromGenericTrackPlugin< TrackClass >::dRcone_, reco::tau::PFRecoTauChargedHadronFromGenericTrackPlugin< TrackClass >::dRconeLimitedToJetArea_, reco::tau::PFRecoTauChargedHadronFromGenericTrackPlugin< TrackClass >::dRmergeNeutralHadron_, reco::tau::PFRecoTauChargedHadronFromGenericTrackPlugin< TrackClass >::dRmergePhoton_, muonDTDigis_cfi::pset, reco::tau::PFRecoTauChargedHadronFromGenericTrackPlugin< TrackClass >::qcuts_, reco::tau::PFRecoTauChargedHadronFromGenericTrackPlugin< TrackClass >::srcTracks_, reco::tau::PFRecoTauChargedHadronFromGenericTrackPlugin< TrackClass >::Tracks_token, and reco::tau::PFRecoTauChargedHadronFromGenericTrackPlugin< TrackClass >::verbosity_.

98  vertexAssociator_(pset.getParameter<edm::ParameterSet>("qualityCuts"), std::move(iC)),
99  qcuts_(nullptr),
101  edm::ParameterSet qcuts_pset = pset.getParameterSet("qualityCuts").getParameterSet("signalQualityCuts");
102  qcuts_ = new RecoTauQualityCuts(qcuts_pset);
103 
104  srcTracks_ = pset.getParameter<edm::InputTag>("srcTracks");
105  Tracks_token = iC.consumes<std::vector<TrackClass> >(srcTracks_);
106  dRcone_ = pset.getParameter<double>("dRcone");
107  dRconeLimitedToJetArea_ = pset.getParameter<bool>("dRconeLimitedToJetArea");
108 
109  dRmergeNeutralHadron_ = pset.getParameter<double>("dRmergeNeutralHadron");
110  dRmergePhoton_ = pset.getParameter<double>("dRmergePhoton");
111 
112  verbosity_ = pset.getParameter<int>("verbosity");
113  }
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > magneticFieldToken_
def move(src, dest)
Definition: eostools.py:511
PFRecoTauChargedHadronBuilderPlugin(const edm::ParameterSet &pset, edm::ConsumesCollector &&iC)

◆ ~PFRecoTauChargedHadronFromGenericTrackPlugin()

Member Function Documentation

◆ beginEvent()

template<class TrackClass >
void reco::tau::PFRecoTauChargedHadronFromGenericTrackPlugin< TrackClass >::beginEvent ( )
overridevirtual

Hook called at the beginning of the event.

Reimplemented from reco::tau::PFRecoTauChargedHadronBuilderPlugin.

Definition at line 122 of file PFRecoTauChargedHadronFromGenericTrackPlugin.cc.

122  {
123  vertexAssociator_.setEvent(*this->evt());
124 
126  }
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
const edm::EventSetup * evtSetup() const
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
void setEvent(const edm::Event &evt)
Load the vertices from the event.
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > magneticFieldToken_

◆ filterTrack() [1/3]

template<class TrackClass >
bool reco::tau::PFRecoTauChargedHadronFromGenericTrackPlugin< TrackClass >::filterTrack ( const edm::Handle< std::vector< TrackClass > > &  ,
size_t  iTrack 
) const
private

◆ filterTrack() [2/3]

template<>
bool reco::tau::PFRecoTauChargedHadronFromGenericTrackPlugin< reco::Track >::filterTrack ( const edm::Handle< std::vector< reco::Track > > &  tracks,
size_t  iTrack 
) const
private

Definition at line 129 of file PFRecoTauChargedHadronFromGenericTrackPlugin.cc.

References pwdgSkimBPark_cfi::tracks.

130  {
131  // ignore tracks which fail quality cuts
132  reco::TrackRef trackRef(tracks, iTrack);
133  return qcuts_->filterTrack(trackRef);
134  }
bool filterTrack(const reco::TrackBaseRef &track) const
Filter a single Track.

◆ filterTrack() [3/3]

template<>
bool reco::tau::PFRecoTauChargedHadronFromGenericTrackPlugin< pat::PackedCandidate >::filterTrack ( const edm::Handle< std::vector< pat::PackedCandidate > > &  tracks,
size_t  iTrack 
) const
private

Definition at line 137 of file PFRecoTauChargedHadronFromGenericTrackPlugin.cc.

138  {
139  // ignore tracks which fail quality cuts
140  const pat::PackedCandidate& cand = (*tracks)[iTrack];
141  if (cand.charge() == 0)
142  return false;
143 
144  return qcuts_->filterChargedCand(cand);
145  }
bool filterChargedCand(const reco::Candidate &cand) const
or a single charged candidate

◆ getTrackPos() [1/3]

template<class TrackClass >
XYZTLorentzVector reco::tau::PFRecoTauChargedHadronFromGenericTrackPlugin< TrackClass >::getTrackPos ( const TrackClass track) const
private

◆ getTrackPos() [2/3]

template<>
XYZTLorentzVector reco::tau::PFRecoTauChargedHadronFromGenericTrackPlugin< reco::Track >::getTrackPos ( const reco::Track track) const
private

Definition at line 182 of file PFRecoTauChargedHadronFromGenericTrackPlugin.cc.

References HLT_2023v12_cff::track.

183  {
184  return XYZTLorentzVector(track.referencePoint().x(), track.referencePoint().y(), track.referencePoint().z(), 0.);
185  }
math::XYZTLorentzVector XYZTLorentzVector
Definition: RawParticle.h:25

◆ getTrackPos() [3/3]

Definition at line 188 of file PFRecoTauChargedHadronFromGenericTrackPlugin.cc.

References HLT_2023v12_cff::track.

189  {
190  return XYZTLorentzVector(track.vertex().x(), track.vertex().y(), track.vertex().z(), 0.);
191  }
math::XYZTLorentzVector XYZTLorentzVector
Definition: RawParticle.h:25

◆ getTrackPtError() [1/3]

template<class TrackClass >
double reco::tau::PFRecoTauChargedHadronFromGenericTrackPlugin< TrackClass >::getTrackPtError ( const TrackClass track) const
private

◆ getTrackPtError() [2/3]

template<>
double reco::tau::PFRecoTauChargedHadronFromGenericTrackPlugin< reco::Track >::getTrackPtError ( const reco::Track track) const
private

Definition at line 160 of file PFRecoTauChargedHadronFromGenericTrackPlugin.cc.

References HLT_2023v12_cff::track.

160  {
161  return track.ptError();
162  }

◆ getTrackPtError() [3/3]

template<>
double reco::tau::PFRecoTauChargedHadronFromGenericTrackPlugin< pat::PackedCandidate >::getTrackPtError ( const pat::PackedCandidate cand) const
private

Definition at line 165 of file PFRecoTauChargedHadronFromGenericTrackPlugin.cc.

References funct::abs(), and HLT_2023v12_cff::track.

166  {
167  double trackPtError =
168  0.06; // MB: Approximate avarage track PtError by 2.5% (barrel), 4% (transition), 6% (endcaps) lostTracks w/o detailed track information available (after TRK-11-001)
169  const reco::Track* track(cand.bestTrack());
170  if (track != nullptr) {
171  trackPtError = track->ptError();
172  } else {
173  if (std::abs(cand.eta()) < 0.9)
174  trackPtError = 0.025;
175  else if (std::abs(cand.eta()) < 1.4)
176  trackPtError = 0.04;
177  }
178  return trackPtError;
179  }
Abs< T >::type abs(const T &t)
Definition: Abs.h:22

◆ operator()()

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

Implements reco::tau::PFRecoTauChargedHadronBuilderPlugin.

Definition at line 206 of file PFRecoTauChargedHadronFromGenericTrackPlugin.cc.

References funct::abs(), reco::tau::atECALEntrance(), HPSPFTauProducerPuppi_cfi::chargedHadron, reco::deltaR(), HGC3DClusterGenMatchSelector_cfi::dR, edm::Event::getByToken(), metsig::jet, reco::PFRecoTauChargedHadron::kTrack, M_PI, SiStripPI::min, eostools::move(), Skims_PA_cff::name, convertSQLitetoXML_cfg::output, AlCaHLTBitMon_ParallelJobs::p, EgammaValidation_cff::pdgId, BaseParticlePropagator::propagateToEcalEntrance(), L1TObjectsTimingClient_cff::resolution, reco::tau::setChargedHadronP4(), jetUpdater_cfi::sort, mathSSE::sqrt(), HLT_2023v12_cff::track, pwdgSkimBPark_cfi::tracks, and extraflags_cff::vtx.

206  {
207  if (verbosity_) {
208  edm::LogPrint("TauChHFromTrack") << "<PFRecoTauChargedHadronFromGenericTrackPlugin::operator()>:";
209  edm::LogPrint("TauChHFromTrack") << " pluginName = " << name();
210  }
211 
213 
214  const edm::Event& evt = (*this->evt());
215 
218 
220  float jEta = jet.eta();
221  float jPhi = jet.phi();
222  size_t numTracks = tracks->size();
223  for (size_t iTrack = 0; iTrack < numTracks; ++iTrack) {
224  const TrackClass& track = (*tracks)[iTrack];
225 
226  // consider tracks in vicinity of tau-jet candidate only
227  double dR = deltaR(track.eta(), track.phi(), jEta, jPhi);
228  double dRmatch = dRcone_;
230  double jetArea = jet.jetArea();
231  if (jetArea > 0.) {
232  dRmatch = std::min(dRmatch, sqrt(jetArea / M_PI));
233  } else {
234  if (numWarnings_ < maxWarnings_) {
235  edm::LogInfo("PFRecoTauChargedHadronFromGenericTrackPlugin::operator()")
236  << "Jet: Pt = " << jet.pt() << ", eta = " << jet.eta() << ", phi = " << jet.phi()
237  << " has area = " << jetArea << " !!" << std::endl;
238  ++numWarnings_;
239  }
240  dRmatch = 0.1;
241  }
242  }
243  if (dR > dRmatch)
244  continue;
245 
246  if (!this->filterTrack(tracks, iTrack))
247  continue;
248 
249  reco::Candidate::Charge trackCharge_int = 0;
250  if (track.charge() > 0.)
251  trackCharge_int = +1;
252  else if (track.charge() < 0.)
253  trackCharge_int = -1;
254 
255  const double chargedPionMass = 0.13957; // GeV
256  double chargedPionP = track.p();
257  double chargedPionEn = TMath::Sqrt(chargedPionP * chargedPionP + chargedPionMass * chargedPionMass);
258  reco::Candidate::LorentzVector chargedPionP4(track.px(), track.py(), track.pz(), chargedPionEn);
259 
260  reco::Vertex::Point vtx(0., 0., 0.);
262  vtx = vertexAssociator_.associatedVertex(jet)->position();
263 
264  std::unique_ptr<PFRecoTauChargedHadron> chargedHadron(
265  new PFRecoTauChargedHadron(trackCharge_int, chargedPionP4, vtx, 0, true, PFRecoTauChargedHadron::kTrack));
266 
268 
269  // CV: Take code for propagating track to ECAL entrance
270  // from RecoParticleFlow/PFTracking/src/PFTrackTransformer.cc
271  // to make sure propagation is done in the same way as for charged PFCandidates.
272  //
273  // The following replacements need to be made
274  // outerMomentum -> momentum
275  // outerPosition -> referencePoint
276  // in order to run on AOD input
277  // (outerMomentum and outerPosition require access to reco::TrackExtra objects, which are available in RECO only)
278  //
279  XYZTLorentzVector chargedPionPos(getTrackPos(track));
280  RawParticle p(chargedPionP4, chargedPionPos);
281  p.setCharge(track.charge());
282  BaseParticlePropagator trackPropagator(p, 0., 0., magneticFieldStrength_.z());
283  trackPropagator.propagateToEcalEntrance(false);
284  if (trackPropagator.getSuccess() != 0) {
285  chargedHadron->positionAtECALEntrance_ = trackPropagator.particle().vertex();
286  } else {
287  if (chargedPionP4.pt() > 2. and std::abs(chargedPionP4.eta()) < 3.) {
288  edm::LogWarning("PFRecoTauChargedHadronFromGenericTrackPlugin::operator()")
289  << "Failed to propagate track: Pt = " << track.pt() << ", eta = " << track.eta()
290  << ", phi = " << track.phi() << " to ECAL entrance !!" << std::endl;
291  }
292  chargedHadron->positionAtECALEntrance_ = math::XYZPointF(0., 0., 0.);
293  }
294 
295  std::vector<Candidate_withDistance> neutralJetConstituents_withDistance;
296  for (const auto& jetConstituent : jet.daughterPtrVector()) {
297  int pdgId = jetConstituent->pdgId();
298  if (!(pdgId == 130 || pdgId == 22))
299  continue;
300  double dR = deltaR(atECALEntrance(&*jetConstituent, magneticFieldStrength_.z()),
301  chargedHadron->positionAtECALEntrance_);
302  double dRmerge = -1.;
303  if (pdgId == 130)
304  dRmerge = dRmergeNeutralHadron_;
305  else if (pdgId == 22)
306  dRmerge = dRmergePhoton_;
307  if (dR < dRmerge) {
308  Candidate_withDistance jetConstituent_withDistance;
309  jetConstituent_withDistance.pfCandidate_ = jetConstituent;
310  jetConstituent_withDistance.distance_ = dR;
311  neutralJetConstituents_withDistance.push_back(jetConstituent_withDistance);
312  chargedHadron->addDaughter(jetConstituent);
313  }
314  }
315  std::sort(
316  neutralJetConstituents_withDistance.begin(), neutralJetConstituents_withDistance.end(), isSmallerDistance);
317 
318  const double caloResolutionCoeff =
319  1.0; // CV: approximate ECAL + HCAL calorimeter resolution for hadrons by 100%*sqrt(E)
320  double resolutionTrackP = track.p() * (getTrackPtError(track) / track.pt());
321  double neutralEnSum = 0.;
322  for (std::vector<Candidate_withDistance>::const_iterator nextNeutral =
323  neutralJetConstituents_withDistance.begin();
324  nextNeutral != neutralJetConstituents_withDistance.end();
325  ++nextNeutral) {
326  double nextNeutralEn = nextNeutral->pfCandidate_->energy();
327  double resolutionCaloEn = caloResolutionCoeff * sqrt(neutralEnSum + nextNeutralEn);
328  double resolution = sqrt(resolutionTrackP * resolutionTrackP + resolutionCaloEn * resolutionCaloEn);
329  if ((neutralEnSum + nextNeutralEn) < (track.p() + 2. * resolution)) {
330  chargedHadron->neutralPFCandidates_.push_back(nextNeutral->pfCandidate_);
331  neutralEnSum += nextNeutralEn;
332  } else {
333  break;
334  }
335  }
336 
338 
339  if (verbosity_) {
340  edm::LogPrint("TauChHFromTrack") << *chargedHadron;
341  }
342 
343  output.push_back(std::move(chargedHadron));
344  }
345 
346  return output;
347  }
int Charge
electric charge type
Definition: Candidate.h:34
std::vector< std::unique_ptr< PFRecoTauChargedHadron > > ChargedHadronVector
math::XYZPointF atECALEntrance(const reco::Candidate *part, double bField)
bool filterTrack(const edm::Handle< std::vector< TrackClass > > &, size_t iTrack) const
XYZTLorentzVector getTrackPos(const TrackClass &track) const
void setChargedHadronTrack(PFRecoTauChargedHadron &chargedHadron, const edm::Ptr< TrackClass > &track) const
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:238
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:540
ROOT::Math::PositionVector3D< ROOT::Math::Cartesian3D< float > > XYZPointF
point in space with cartesian internal representation
Definition: Point3D.h:10
double getTrackPtError(const TrackClass &track) const
reco::VertexRef associatedVertex(const Jet &jet) const
T sqrt(T t)
Definition: SSEVec.h:19
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
math::XYZPoint Point
point in the space
Definition: Vertex.h:40
Log< level::Warning, true > LogPrint
bool propagateToEcalEntrance(bool first=true)
#define M_PI
constexpr auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:30
Log< level::Info, false > LogInfo
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:36
const std::string & name() const
Log< level::Warning, false > LogWarning
void setChargedHadronP4(reco::PFRecoTauChargedHadron &chargedHadron, double scaleFactor_neutralPFCands=1.0)
void setPV(const reco::VertexRef &vtx)
Update the primary vertex.
def move(src, dest)
Definition: eostools.py:511
math::XYZTLorentzVector XYZTLorentzVector
Definition: RawParticle.h:25

◆ setChargedHadronTrack() [1/3]

template<class TrackClass >
void reco::tau::PFRecoTauChargedHadronFromGenericTrackPlugin< TrackClass >::setChargedHadronTrack ( PFRecoTauChargedHadron chargedHadron,
const edm::Ptr< TrackClass > &  track 
) const
private

◆ setChargedHadronTrack() [2/3]

template<>
void reco::tau::PFRecoTauChargedHadronFromGenericTrackPlugin< reco::Track >::setChargedHadronTrack ( PFRecoTauChargedHadron chargedHadron,
const edm::Ptr< reco::Track > &  track 
) const
private

◆ setChargedHadronTrack() [3/3]

template<>
void reco::tau::PFRecoTauChargedHadronFromGenericTrackPlugin< pat::PackedCandidate >::setChargedHadronTrack ( PFRecoTauChargedHadron chargedHadron,
const edm::Ptr< pat::PackedCandidate > &  track 
) const
private

Member Data Documentation

◆ dRcone_

template<class TrackClass >
double reco::tau::PFRecoTauChargedHadronFromGenericTrackPlugin< TrackClass >::dRcone_
private

◆ dRconeLimitedToJetArea_

template<class TrackClass >
bool reco::tau::PFRecoTauChargedHadronFromGenericTrackPlugin< TrackClass >::dRconeLimitedToJetArea_
private

◆ dRmergeNeutralHadron_

template<class TrackClass >
double reco::tau::PFRecoTauChargedHadronFromGenericTrackPlugin< TrackClass >::dRmergeNeutralHadron_
private

◆ dRmergePhoton_

template<class TrackClass >
double reco::tau::PFRecoTauChargedHadronFromGenericTrackPlugin< TrackClass >::dRmergePhoton_
private

◆ magneticFieldStrength_

template<class TrackClass >
math::XYZVector reco::tau::PFRecoTauChargedHadronFromGenericTrackPlugin< TrackClass >::magneticFieldStrength_
private

◆ magneticFieldToken_

◆ maxWarnings_

template<class TrackClass >
constexpr unsigned int reco::tau::PFRecoTauChargedHadronFromGenericTrackPlugin< TrackClass >::maxWarnings_ = 3
staticprivate

◆ numWarnings_

template<class TrackClass >
std::atomic< unsigned int > reco::tau::PFRecoTauChargedHadronFromGenericTrackPlugin< TrackClass >::numWarnings_ {0}
staticprivate

◆ qcuts_

◆ srcTracks_

template<class TrackClass >
edm::InputTag reco::tau::PFRecoTauChargedHadronFromGenericTrackPlugin< TrackClass >::srcTracks_
private

◆ Tracks_token

template<class TrackClass >
edm::EDGetTokenT<std::vector<TrackClass> > reco::tau::PFRecoTauChargedHadronFromGenericTrackPlugin< TrackClass >::Tracks_token
private

◆ verbosity_

template<class TrackClass >
int reco::tau::PFRecoTauChargedHadronFromGenericTrackPlugin< TrackClass >::verbosity_
private

◆ vertexAssociator_

template<class TrackClass >
RecoTauVertexAssociator reco::tau::PFRecoTauChargedHadronFromGenericTrackPlugin< TrackClass >::vertexAssociator_
private