CMS 3D CMS Logo

List of all members | Public Types | Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes
TemplatedInclusiveVertexFinder< InputContainer, VTX > Class Template Reference

#include <InclusiveVertexFinder.h>

Inheritance diagram for TemplatedInclusiveVertexFinder< InputContainer, VTX >:
edm::stream::EDProducer<>

Public Types

typedef std::vector< VTX > Product
 
typedef InputContainer::value_type TRK
 
- Public Types inherited from edm::stream::EDProducer<>
typedef CacheContexts< T... > CacheTypes
 
typedef CacheTypes::GlobalCache GlobalCache
 
typedef AbilityChecker< T... > HasAbility
 
typedef CacheTypes::LuminosityBlockCache LuminosityBlockCache
 
typedef LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCacheLuminosityBlockContext
 
typedef CacheTypes::LuminosityBlockSummaryCache LuminosityBlockSummaryCache
 
typedef CacheTypes::RunCache RunCache
 
typedef RunContextT< RunCache, GlobalCacheRunContext
 
typedef CacheTypes::RunSummaryCache RunSummaryCache
 

Public Member Functions

virtual void produce (edm::Event &event, const edm::EventSetup &es) override
 
 TemplatedInclusiveVertexFinder (const edm::ParameterSet &params)
 
- Public Member Functions inherited from edm::stream::EDProducer<>
 EDProducer ()=default
 
bool hasAbilityToProduceInLumis () const final
 
bool hasAbilityToProduceInRuns () const final
 

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &cdesc)
 

Private Member Functions

std::pair< std::vector< reco::TransientTrack >, GlobalPointnearTracks (const reco::TransientTrack &seed, const std::vector< reco::TransientTrack > &tracks, const reco::Vertex &primaryVertex) const
 
bool trackFilter (const reco::Track &track) const
 

Private Attributes

std::auto_ptr< TracksClusteringFromDisplacedSeedclusterizer
 
double fitterRatio
 
double fitterSigmacut
 
double fitterTini
 
double maxLIP
 
unsigned int maxNTracks
 
double maxTimeSig
 
unsigned int minHits
 
double minPt
 
edm::EDGetTokenT< reco::BeamSpottoken_beamSpot
 
edm::EDGetTokenT< reco::VertexCollectiontoken_primaryVertex
 
edm::EDGetTokenT< InputContainer > token_tracks
 
bool useVertexFitter
 
bool useVertexReco
 
double vertexMinAngleCosine
 
double vertexMinDLen2DSig
 
double vertexMinDLenSig
 
std::auto_ptr< VertexReconstructorvtxReco
 

Detailed Description

template<class InputContainer, class VTX>
class TemplatedInclusiveVertexFinder< InputContainer, VTX >

Definition at line 41 of file InclusiveVertexFinder.h.

Member Typedef Documentation

template<class InputContainer , class VTX >
typedef std::vector<VTX> TemplatedInclusiveVertexFinder< InputContainer, VTX >::Product

Definition at line 43 of file InclusiveVertexFinder.h.

template<class InputContainer , class VTX >
typedef InputContainer::value_type TemplatedInclusiveVertexFinder< InputContainer, VTX >::TRK

Definition at line 44 of file InclusiveVertexFinder.h.

Constructor & Destructor Documentation

template<class InputContainer , class VTX >
TemplatedInclusiveVertexFinder< InputContainer, VTX >::TemplatedInclusiveVertexFinder ( const edm::ParameterSet params)

Definition at line 129 of file InclusiveVertexFinder.h.

References edm::ParameterSet::getParameter(), TemplatedInclusiveVertexFinder< InputContainer, VTX >::token_beamSpot, TemplatedInclusiveVertexFinder< InputContainer, VTX >::token_primaryVertex, and TemplatedInclusiveVertexFinder< InputContainer, VTX >::token_tracks.

129  :
130  minHits(params.getParameter<unsigned int>("minHits")),
131  maxNTracks(params.getParameter<unsigned int>("maxNTracks")),
132  maxLIP(params.getParameter<double>("maximumLongitudinalImpactParameter")),
133  maxTimeSig(params.getParameter<double>("maximumTimeSignificance")),
134  minPt(params.getParameter<double>("minPt")), //0.8
135  vertexMinAngleCosine(params.getParameter<double>("vertexMinAngleCosine")), //0.98
136  vertexMinDLen2DSig(params.getParameter<double>("vertexMinDLen2DSig")), //2.5
137  vertexMinDLenSig(params.getParameter<double>("vertexMinDLenSig")), //0.5
138  fitterSigmacut(params.getParameter<double>("fitterSigmacut")),
139  fitterTini(params.getParameter<double>("fitterTini")),
140  fitterRatio(params.getParameter<double>("fitterRatio")),
141  useVertexFitter(params.getParameter<bool>("useDirectVertexFitter")),
142  useVertexReco(params.getParameter<bool>("useVertexReco")),
145 
146 {
147  token_beamSpot = consumes<reco::BeamSpot>(params.getParameter<edm::InputTag>("beamSpot"));
148  token_primaryVertex = consumes<reco::VertexCollection>(params.getParameter<edm::InputTag>("primaryVertices"));
149  token_tracks = consumes<InputContainer>(params.getParameter<edm::InputTag>("tracks"));
150  produces<Product>();
151  //produces<reco::VertexCollection>("multi");
152 }
T getParameter(std::string const &) const
std::auto_ptr< TracksClusteringFromDisplacedSeed > clusterizer
edm::EDGetTokenT< InputContainer > token_tracks
edm::EDGetTokenT< reco::BeamSpot > token_beamSpot
edm::EDGetTokenT< reco::VertexCollection > token_primaryVertex
std::auto_ptr< VertexReconstructor > vtxReco

Member Function Documentation

template<class InputContainer , class VTX >
static void TemplatedInclusiveVertexFinder< InputContainer, VTX >::fillDescriptions ( edm::ConfigurationDescriptions cdesc)
inlinestatic

Definition at line 47 of file InclusiveVertexFinder.h.

References edm::ConfigurationDescriptions::add(), edm::ParameterSetDescription::add(), edm::ConfigurationDescriptions::addDefault(), TemplatedInclusiveVertexFinder< InputContainer, VTX >::clusterizer, TemplatedInclusiveVertexFinder< InputContainer, VTX >::nearTracks(), impactParameterTagInfos_cfi::primaryVertex, TemplatedInclusiveVertexFinder< InputContainer, VTX >::produce(), SurveyInfoScenario_cff::seed, AlCaHLTBitMon_QueryRunRegistry::string, HiIsolationCommonParameters_cff::track, TemplatedInclusiveVertexFinder< InputContainer, VTX >::trackFilter(), l1t::tracks, relativeConstraints::value, and ghostTrackVertexReco_cff::vertexReco.

47  {
49  pdesc.add<edm::InputTag>("beamSpot",edm::InputTag("offlineBeamSpot"));
50  pdesc.add<edm::InputTag>("primaryVertices",edm::InputTag("offlinePrimaryVertices"));
52  pdesc.add<edm::InputTag>("tracks",edm::InputTag("generalTracks"));
53  pdesc.add<unsigned int>("minHits",8);
55  pdesc.add<edm::InputTag>("tracks",edm::InputTag("particleFlow"));
56  pdesc.add<unsigned int>("minHits",0);
57  } else {
58  pdesc.add<edm::InputTag>("tracks",edm::InputTag("generalTracks"));
59  }
60 
61  pdesc.add<double>("maximumLongitudinalImpactParameter",0.3);
62  pdesc.add<double>("maximumTimeSignificance",3.0);
63  pdesc.add<double>("minPt",0.8);
64  pdesc.add<unsigned int>("maxNTracks",30);
65  //clusterizer pset
67  clusterizer.add<double>("seedMax3DIPSignificance",9999.0);
68  clusterizer.add<double>("seedMax3DIPValue",9999.0);
69  clusterizer.add<double>("seedMin3DIPSignificance",1.2);
70  clusterizer.add<double>("seedMin3DIPValue",0.005);
71  clusterizer.add<double>("clusterMaxDistance",0.05);
72  clusterizer.add<double>("clusterMaxSignificance",4.5);
73  clusterizer.add<double>("distanceRatio",20.0);
74  clusterizer.add<double>("clusterMinAngleCosine",0.5);
75  clusterizer.add<double>("maxTimeSignificance",3.5);
76  pdesc.add<edm::ParameterSetDescription>("clusterizer", clusterizer);
77  // vertex and fitter config
78  pdesc.add<double>("vertexMinAngleCosine",0.95);
79  pdesc.add<double>("vertexMinDLen2DSig",2.5);
80  pdesc.add<double>("vertexMinDLenSig",0.5);
81  pdesc.add<double>("fitterSigmacut",3.0);
82  pdesc.add<double>("fitterTini",256.0);
83  pdesc.add<double>("fitterRatio",0.25);
84  pdesc.add<bool>("useDirectVertexFitter",true);
85  pdesc.add<bool>("useVertexReco",true);
86  // vertexReco pset
88  vertexReco.add<std::string>("finder", std::string("avr"));
89  vertexReco.add<double>("primcut",1.0);
90  vertexReco.add<double>("seccut",3.0);
91  vertexReco.add<bool>("smoothing",true);
92  pdesc.add<edm::ParameterSetDescription>("vertexReco", vertexReco);
94  cdesc.add("inclusiveVertexFinderDefault",pdesc);
96  cdesc.add("inclusiveCandidateVertexFinderDefault",pdesc);
97  } else {
98  cdesc.addDefault(pdesc);
99  }
100  }
std::auto_ptr< TracksClusteringFromDisplacedSeed > clusterizer
void addDefault(ParameterSetDescription const &psetDescription)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
template<class InputContainer , class VTX >
std::pair<std::vector<reco::TransientTrack>,GlobalPoint> TemplatedInclusiveVertexFinder< InputContainer, VTX >::nearTracks ( const reco::TransientTrack seed,
const std::vector< reco::TransientTrack > &  tracks,
const reco::Vertex primaryVertex 
) const
private
template<class InputContainer , class VTX >
void TemplatedInclusiveVertexFinder< InputContainer, VTX >::produce ( edm::Event event,
const edm::EventSetup es 
)
overridevirtual

Definition at line 165 of file InclusiveVertexFinder.h.

References funct::abs(), ecalDrivenElectronSeedsParameters_cff::beamSpot, tthelpers::buildTT(), TemplatedInclusiveVertexFinder< InputContainer, VTX >::clusterizer, fastPrimaryVertexProducer_cfi::clusters, gather_cfg::cout, reco::Vertex::covariance(), dir, VertexDistance3D::distance(), VertexDistanceXY::distance(), Vector3DBase< T, FrameTag >::dot(), reco::TransientTrack::dtErrorExt(), reco::TrackBase::dz(), Measurement1D::error(), TemplatedInclusiveVertexFinder< InputContainer, VTX >::fitterRatio, TemplatedInclusiveVertexFinder< InputContainer, VTX >::fitterSigmacut, TemplatedInclusiveVertexFinder< InputContainer, VTX >::fitterTini, edm::EventSetup::get(), mps_fire::i, edm::isFinite(), reco::TransientTrack::isValid(), TransientVertex::isValid(), TemplatedInclusiveVertexFinder< InputContainer, VTX >::maxLIP, TemplatedInclusiveVertexFinder< InputContainer, VTX >::maxNTracks, TemplatedInclusiveVertexFinder< InputContainer, VTX >::maxTimeSig, eostools::move(), reco::Vertex::position(), funct::pow(), jets_cff::primaryVertices, MetAnalyzer::pv(), reco::TransientTrack::setBeamSpot(), Measurement1D::significance(), mathSSE::sqrt(), pfDeepBoostedJetPreprocessParams_cfi::sv, reco::Vertex::t(), reco::TransientTrack::timeExt(), TemplatedInclusiveVertexFinder< InputContainer, VTX >::token_beamSpot, TemplatedInclusiveVertexFinder< InputContainer, VTX >::token_primaryVertex, TemplatedInclusiveVertexFinder< InputContainer, VTX >::token_tracks, HiIsolationCommonParameters_cff::track, reco::TransientTrack::track(), TemplatedInclusiveVertexFinder< InputContainer, VTX >::trackFilter(), l1t::tracks, groupFilesInBlocks::tt, csvLumiCalc::unit, Vector3DBase< T, FrameTag >::unit(), reco::Unknown, svhelper::updateVertexTime(), TemplatedInclusiveVertexFinder< InputContainer, VTX >::useVertexFitter, TemplatedInclusiveVertexFinder< InputContainer, VTX >::useVertexReco, findQualityFiles::v, Measurement1D::value(), AdaptiveVertexFitter::vertex(), TemplatedInclusiveVertexFinder< InputContainer, VTX >::vertexMinAngleCosine, TemplatedInclusiveVertexFinder< InputContainer, VTX >::vertexMinDLen2DSig, TemplatedInclusiveVertexFinder< InputContainer, VTX >::vertexMinDLenSig, electrons_cff::vertices, badGlobalMuonTaggersAOD_cff::vtx, TemplatedInclusiveVertexFinder< InputContainer, VTX >::vtxReco, and w.

Referenced by TemplatedInclusiveVertexFinder< InputContainer, VTX >::fillDescriptions().

166 {
167  using namespace reco;
168 
169  VertexDistance3D vdist;
170  VertexDistanceXY vdist2d;
171  MultiVertexFitter theMultiVertexFitter;
172  AdaptiveVertexFitter theAdaptiveFitter(
178 
179 
181  event.getByToken(token_beamSpot,beamSpot);
182 
184  event.getByToken(token_primaryVertex, primaryVertices);
185 
187  event.getByToken(token_tracks, tracks);
188 
190  es.get<TransientTrackRecord>().get("TransientTrackBuilder",
191  trackBuilder);
192 
193 
194  auto recoVertices = std::make_unique<Product>();
195  if(primaryVertices->size()!=0) {
196 
197  const reco::Vertex &pv = (*primaryVertices)[0];
198  GlobalPoint ppv(pv.position().x(),pv.position().y(),pv.position().z());
199 
200  std::vector<TransientTrack> tts;
201  //Fill transient track vector
202  for(typename InputContainer::const_iterator track = tracks->begin();
203  track != tracks->end(); ++track) {
204 //TransientTrack tt = trackBuilder->build(ref);
205 //TrackRef ref(tracks, track - tracks->begin());
206  TransientTrack tt(tthelpers::buildTT(tracks,trackBuilder, track - tracks->begin()));
207  if(!tt.isValid()) continue;
208  if (!trackFilter(tt.track()))
209  continue;
210  if( std::abs(tt.track().dz(pv.position())) > maxLIP)
211  continue;
212  if( edm::isFinite(tt.timeExt()) && pv.covariance(3,3) > 0. ) { // only apply if time available
213  auto tError = std::sqrt( std::pow(tt.dtErrorExt(),2) + pv.covariance(3,3) );
214  auto dtSig = std::abs(tt.timeExt() - pv.t())/tError;
215  if( dtSig > maxTimeSig ) continue;
216  }
217  tt.setBeamSpot(*beamSpot);
218  tts.push_back(tt);
219  }
220  std::vector<TracksClusteringFromDisplacedSeed::Cluster> clusters = clusterizer->clusters(pv,tts);
221 
222  //Create BS object from PV to feed in the AVR
224  for(unsigned int i = 0; i < 7; i++) {
225  for(unsigned int j = 0; j < 7; j++) {
226  if (i < 3 && j < 3)
227  cov(i, j) = pv.covariance(i, j);
228  else
229  cov(i, j) = 0.0;
230  }
231  }
232  BeamSpot bs(pv.position(), 0.0, 0.0, 0.0, 0.0, cov, BeamSpot::Unknown);
233 
234 
235  int i=0;
236 #ifdef VTXDEBUG
237 
238  std::cout << "CLUSTERS " << clusters.size() << std::endl;
239 #endif
240 
241  for(std::vector<TracksClusteringFromDisplacedSeed::Cluster>::iterator cluster = clusters.begin();
242  cluster != clusters.end(); ++cluster,++i)
243  {
244  if(cluster->tracks.size() < 2 || cluster->tracks.size() > maxNTracks )
245  continue;
246  std::vector<TransientVertex> vertices;
247  if(useVertexReco) {
248  vertices = vtxReco->vertices(cluster->tracks, bs); // attempt with config given reconstructor
249  }
250  TransientVertex singleFitVertex;
251  if(useVertexFitter) {
252  singleFitVertex = theAdaptiveFitter.vertex(cluster->tracks,cluster->seedPoint); //attempt with direct fitting
253  if(singleFitVertex.isValid())
254  vertices.push_back(singleFitVertex);
255  }
256 
257  // for each transient vertex state determine if a time can be measured and fill covariance
258  if( pv.covariance(3,3) > 0. ) {
259  for(auto& vtx : vertices) {
261  }
262  }
263 
264  for(std::vector<TransientVertex>::const_iterator v = vertices.begin();
265  v != vertices.end(); ++v) {
266  Measurement1D dlen= vdist.distance(pv,*v);
267  Measurement1D dlen2= vdist2d.distance(pv,*v);
268 #ifdef VTXDEBUG
269  VTX vv(*v);
270  std::cout << "V chi2/n: " << v->normalisedChiSquared() << " ndof: " <<v->degreesOfFreedom() ;
271  std::cout << " dlen: " << dlen.value() << " error: " << dlen.error() << " signif: " << dlen.significance();
272  std::cout << " dlen2: " << dlen2.value() << " error2: " << dlen2.error() << " signif2: " << dlen2.significance();
273  std::cout << " pos: " << vv.position() << " error: " <<vv.xError() << " " << vv.yError() << " " << vv.zError() << std::endl;
274  std::cout << " time: " << vv.time() << " error: " << vv.tError() << std::endl;
275 #endif
276  GlobalVector dir;
277  std::vector<reco::TransientTrack> ts = v->originalTracks();
278  for(std::vector<reco::TransientTrack>::const_iterator i = ts.begin();
279  i != ts.end(); ++i) {
280  float w = v->trackWeight(*i);
281  if (w > 0.5) dir+=i->impactPointState().globalDirection();
282 #ifdef VTXDEBUG
283  std::cout << "\t[" << (*i).track().pt() << ": "
284  << (*i).track().eta() << ", "
285  << (*i).track().phi() << "], "
286  << w << std::endl;
287 #endif
288  }
289  GlobalPoint sv((*v).position().x(),(*v).position().y(),(*v).position().z());
290  float vscal = dir.unit().dot((sv-ppv).unit());
291  if(dlen.significance() > vertexMinDLenSig &&
292  ( (vertexMinAngleCosine > 0) ? (vscal > vertexMinAngleCosine) : (vscal < vertexMinAngleCosine) )
293  && v->normalisedChiSquared() < 10 && dlen2.significance() > vertexMinDLen2DSig)
294  {
295  recoVertices->push_back(*v);
296 
297 #ifdef VTXDEBUG
298  std::cout << "ADDED" << std::endl;
299 #endif
300  }
301 
302  }
303  }
304 #ifdef VTXDEBUG
305 
306  std::cout << "Final put " << recoVertices->size() << std::endl;
307 #endif
308  }
309 
310  event.put(std::move(recoVertices));
311 
312 }
math::Error< dimension >::type CovarianceMatrix
Definition: BeamSpot.h:31
std::auto_ptr< TracksClusteringFromDisplacedSeed > clusterizer
const double w
Definition: UKUtility.cc:23
double error() const
Definition: Measurement1D.h:30
PreciseFloatType< T, U >::Type dot(const Vector3DBase< U, FrameTag > &v) const
Definition: Vector3DBase.h:107
double covariance(int i, int j) const
(i, j)-th element of error matrix, i, j = 0, ... 2
Definition: Vertex.h:130
const Point & position() const
position
Definition: Vertex.h:109
Measurement1D distance(const GlobalPoint &vtx1Position, const GlobalError &vtx1PositionError, const GlobalPoint &vtx2Position, const GlobalError &vtx2PositionError) const override
bool isFinite(T x)
Measurement1D distance(const GlobalPoint &vtx1Position, const GlobalError &vtx1PositionError, const GlobalPoint &vtx2Position, const GlobalError &vtx2PositionError) const override
T sqrt(T t)
Definition: SSEVec.h:18
def pv(vc)
Definition: MetAnalyzer.py:6
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
edm::EDGetTokenT< InputContainer > token_tracks
bool trackFilter(const reco::Track &track) const
Vector3DBase unit() const
Definition: Vector3DBase.h:57
double significance() const
Definition: Measurement1D.h:32
primaryVertices
Definition: jets_cff.py:130
double value() const
Definition: Measurement1D.h:28
void updateVertexTime(TransientVertex &vtx)
Definition: SVTimeHelpers.h:11
significance
Definition: met_cff.py:30
fixed size matrix
T get() const
Definition: EventSetup.h:63
edm::EDGetTokenT< reco::BeamSpot > token_beamSpot
dbl *** dir
Definition: mlp_gen.cc:35
edm::EDGetTokenT< reco::VertexCollection > token_primaryVertex
std::auto_ptr< VertexReconstructor > vtxReco
reco::TransientTrack buildTT(edm::Handle< reco::TrackCollection > &tracks, edm::ESHandle< TransientTrackBuilder > &trackbuilder, unsigned int k)
Definition: TTHelpers.h:10
bool isValid() const
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
def move(src, dest)
Definition: eostools.py:510
double t() const
t coordinate
Definition: Vertex.h:117
template<class InputContainer , class VTX >
bool TemplatedInclusiveVertexFinder< InputContainer, VTX >::trackFilter ( const reco::Track track) const
private

Definition at line 154 of file InclusiveVertexFinder.h.

References reco::TrackBase::hitPattern(), createfilelist::int, TemplatedInclusiveVertexFinder< InputContainer, VTX >::minHits, TemplatedInclusiveVertexFinder< InputContainer, VTX >::minPt, reco::HitPattern::numberOfValidHits(), and reco::TrackBase::pt().

Referenced by TemplatedInclusiveVertexFinder< InputContainer, VTX >::fillDescriptions(), and TemplatedInclusiveVertexFinder< InputContainer, VTX >::produce().

155 {
156  if (track.hitPattern().numberOfValidHits() < (int)minHits)
157  return false;
158  if (track.pt() < minPt )
159  return false;
160 
161  return true;
162 }
int numberOfValidHits() const
Definition: HitPattern.h:824
double pt() const
track transverse momentum
Definition: TrackBase.h:621
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:446

Member Data Documentation

template<class InputContainer , class VTX >
std::auto_ptr<TracksClusteringFromDisplacedSeed> TemplatedInclusiveVertexFinder< InputContainer, VTX >::clusterizer
private
template<class InputContainer , class VTX >
double TemplatedInclusiveVertexFinder< InputContainer, VTX >::fitterRatio
private
template<class InputContainer , class VTX >
double TemplatedInclusiveVertexFinder< InputContainer, VTX >::fitterSigmacut
private
template<class InputContainer , class VTX >
double TemplatedInclusiveVertexFinder< InputContainer, VTX >::fitterTini
private
template<class InputContainer , class VTX >
double TemplatedInclusiveVertexFinder< InputContainer, VTX >::maxLIP
private
template<class InputContainer , class VTX >
unsigned int TemplatedInclusiveVertexFinder< InputContainer, VTX >::maxNTracks
private
template<class InputContainer , class VTX >
double TemplatedInclusiveVertexFinder< InputContainer, VTX >::maxTimeSig
private
template<class InputContainer , class VTX >
unsigned int TemplatedInclusiveVertexFinder< InputContainer, VTX >::minHits
private
template<class InputContainer , class VTX >
double TemplatedInclusiveVertexFinder< InputContainer, VTX >::minPt
private
template<class InputContainer , class VTX >
edm::EDGetTokenT<reco::BeamSpot> TemplatedInclusiveVertexFinder< InputContainer, VTX >::token_beamSpot
private
template<class InputContainer , class VTX >
edm::EDGetTokenT<reco::VertexCollection> TemplatedInclusiveVertexFinder< InputContainer, VTX >::token_primaryVertex
private
template<class InputContainer , class VTX >
edm::EDGetTokenT<InputContainer> TemplatedInclusiveVertexFinder< InputContainer, VTX >::token_tracks
private
template<class InputContainer , class VTX >
bool TemplatedInclusiveVertexFinder< InputContainer, VTX >::useVertexFitter
private
template<class InputContainer , class VTX >
bool TemplatedInclusiveVertexFinder< InputContainer, VTX >::useVertexReco
private
template<class InputContainer , class VTX >
double TemplatedInclusiveVertexFinder< InputContainer, VTX >::vertexMinAngleCosine
private
template<class InputContainer , class VTX >
double TemplatedInclusiveVertexFinder< InputContainer, VTX >::vertexMinDLen2DSig
private
template<class InputContainer , class VTX >
double TemplatedInclusiveVertexFinder< InputContainer, VTX >::vertexMinDLenSig
private
template<class InputContainer , class VTX >
std::auto_ptr<VertexReconstructor> TemplatedInclusiveVertexFinder< InputContainer, VTX >::vtxReco
private