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<>
using CacheTypes = CacheContexts< T... >
 
using GlobalCache = typename CacheTypes::GlobalCache
 
using HasAbility = AbilityChecker< T... >
 
using InputProcessBlockCache = typename CacheTypes::InputProcessBlockCache
 
using LuminosityBlockCache = typename CacheTypes::LuminosityBlockCache
 
using LuminosityBlockContext = LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCache >
 
using LuminosityBlockSummaryCache = typename CacheTypes::LuminosityBlockSummaryCache
 
using RunCache = typename CacheTypes::RunCache
 
using RunContext = RunContextT< RunCache, GlobalCache >
 
using RunSummaryCache = typename CacheTypes::RunSummaryCache
 

Public Member Functions

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
 
 EDProducer (const EDProducer &)=delete
 
bool hasAbilityToProduceInBeginLumis () const final
 
bool hasAbilityToProduceInBeginProcessBlocks () const final
 
bool hasAbilityToProduceInBeginRuns () const final
 
bool hasAbilityToProduceInEndLumis () const final
 
bool hasAbilityToProduceInEndProcessBlocks () const final
 
bool hasAbilityToProduceInEndRuns () const final
 
const EDProduceroperator= (const EDProducer &)=delete
 

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::unique_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::ESGetToken< TransientTrackBuilder, TransientTrackRecordtoken_trackBuilder
 
edm::EDGetTokenT< InputContainer > token_tracks
 
bool useVertexFitter
 
bool useVertexReco
 
double vertexMinAngleCosine
 
double vertexMinDLen2DSig
 
double vertexMinDLenSig
 
std::unique_ptr< VertexReconstructorvtxReco
 

Detailed Description

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

Definition at line 42 of file InclusiveVertexFinder.h.

Member Typedef Documentation

◆ Product

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

Definition at line 44 of file InclusiveVertexFinder.h.

◆ TRK

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

Definition at line 45 of file InclusiveVertexFinder.h.

Constructor & Destructor Documentation

◆ TemplatedInclusiveVertexFinder()

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

Definition at line 132 of file InclusiveVertexFinder.h.

References submitPVValidationJobs::params, TemplatedInclusiveVertexFinder< InputContainer, VTX >::token_beamSpot, TemplatedInclusiveVertexFinder< InputContainer, VTX >::token_primaryVertex, TemplatedInclusiveVertexFinder< InputContainer, VTX >::token_trackBuilder, and TemplatedInclusiveVertexFinder< InputContainer, VTX >::token_tracks.

133  : minHits(params.getParameter<unsigned int>("minHits")),
134  maxNTracks(params.getParameter<unsigned int>("maxNTracks")),
135  maxLIP(params.getParameter<double>("maximumLongitudinalImpactParameter")),
136  maxTimeSig(params.getParameter<double>("maximumTimeSignificance")),
137  minPt(params.getParameter<double>("minPt")), //0.8
138  vertexMinAngleCosine(params.getParameter<double>("vertexMinAngleCosine")), //0.98
139  vertexMinDLen2DSig(params.getParameter<double>("vertexMinDLen2DSig")), //2.5
140  vertexMinDLenSig(params.getParameter<double>("vertexMinDLenSig")), //0.5
141  fitterSigmacut(params.getParameter<double>("fitterSigmacut")),
142  fitterTini(params.getParameter<double>("fitterTini")),
143  fitterRatio(params.getParameter<double>("fitterRatio")),
144  useVertexFitter(params.getParameter<bool>("useDirectVertexFitter")),
145  useVertexReco(params.getParameter<bool>("useVertexReco")),
146  vtxReco(new ConfigurableVertexReconstructor(params.getParameter<edm::ParameterSet>("vertexReco"))),
147  clusterizer(new TracksClusteringFromDisplacedSeed(params.getParameter<edm::ParameterSet>("clusterizer")))
148 
149 {
150  token_beamSpot = consumes<reco::BeamSpot>(params.getParameter<edm::InputTag>("beamSpot"));
151  token_primaryVertex = consumes<reco::VertexCollection>(params.getParameter<edm::InputTag>("primaryVertices"));
152  token_tracks = consumes<InputContainer>(params.getParameter<edm::InputTag>("tracks"));
154  esConsumes<TransientTrackBuilder, TransientTrackRecord>(edm::ESInputTag("", "TransientTrackBuilder"));
155  produces<Product>();
156  //produces<reco::VertexCollection>("multi");
157 }
std::unique_ptr< VertexReconstructor > vtxReco
std::unique_ptr< TracksClusteringFromDisplacedSeed > clusterizer
edm::ESGetToken< TransientTrackBuilder, TransientTrackRecord > token_trackBuilder
edm::EDGetTokenT< InputContainer > token_tracks
edm::EDGetTokenT< reco::BeamSpot > token_beamSpot
edm::EDGetTokenT< reco::VertexCollection > token_primaryVertex

Member Function Documentation

◆ fillDescriptions()

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

Definition at line 48 of file InclusiveVertexFinder.h.

References edm::ConfigurationDescriptions::add(), edm::ParameterSetDescription::add(), edm::ConfigurationDescriptions::addDefault(), TemplatedInclusiveVertexFinder< InputContainer, VTX >::clusterizer, ProducerED_cfi::InputTag, AlCaHLTBitMon_QueryRunRegistry::string, relativeConstraints::value, and HLT_2024v14_cff::vertexReco.

48  {
50  pdesc.add<edm::InputTag>("beamSpot", edm::InputTag("offlineBeamSpot"));
51  pdesc.add<edm::InputTag>("primaryVertices", edm::InputTag("offlinePrimaryVertices"));
53  pdesc.add<edm::InputTag>("tracks", edm::InputTag("generalTracks"));
54  pdesc.add<unsigned int>("minHits", 8);
56  pdesc.add<edm::InputTag>("tracks", edm::InputTag("particleFlow"));
57  pdesc.add<unsigned int>("minHits", 0);
58  } else {
59  pdesc.add<edm::InputTag>("tracks", edm::InputTag("generalTracks"));
60  }
61 
62  pdesc.add<double>("maximumLongitudinalImpactParameter", 0.3);
63  pdesc.add<double>("maximumTimeSignificance", 3.0);
64  pdesc.add<double>("minPt", 0.8);
65  pdesc.add<unsigned int>("maxNTracks", 30);
66  //clusterizer pset
68  clusterizer.add<double>("seedMax3DIPSignificance", 9999.0);
69  clusterizer.add<double>("seedMax3DIPValue", 9999.0);
70  clusterizer.add<double>("seedMin3DIPSignificance", 1.2);
71  clusterizer.add<double>("seedMin3DIPValue", 0.005);
72  clusterizer.add<double>("clusterMaxDistance", 0.05);
73  clusterizer.add<double>("clusterMaxSignificance", 4.5);
74  clusterizer.add<double>("distanceRatio", 20.0);
75  clusterizer.add<double>("clusterMinAngleCosine", 0.5);
76  clusterizer.add<double>("maxTimeSignificance", 3.5);
77  pdesc.add<edm::ParameterSetDescription>("clusterizer", clusterizer);
78  // vertex and fitter config
79  pdesc.add<double>("vertexMinAngleCosine", 0.95);
80  pdesc.add<double>("vertexMinDLen2DSig", 2.5);
81  pdesc.add<double>("vertexMinDLenSig", 0.5);
82  pdesc.add<double>("fitterSigmacut", 3.0);
83  pdesc.add<double>("fitterTini", 256.0);
84  pdesc.add<double>("fitterRatio", 0.25);
85  pdesc.add<bool>("useDirectVertexFitter", true);
86  pdesc.add<bool>("useVertexReco", true);
87  // vertexReco pset
89  vertexReco.add<std::string>("finder", std::string("avr"));
90  vertexReco.add<double>("primcut", 1.0);
91  vertexReco.add<double>("seccut", 3.0);
92  vertexReco.add<bool>("smoothing", true);
93  pdesc.add<edm::ParameterSetDescription>("vertexReco", vertexReco);
95  cdesc.add("inclusiveVertexFinderDefault", pdesc);
97  cdesc.add("inclusiveCandidateVertexFinderDefault", pdesc);
98  } else {
99  cdesc.addDefault(pdesc);
100  }
101  }
std::unique_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)

◆ nearTracks()

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

◆ produce()

template<class InputContainer , class VTX >
void TemplatedInclusiveVertexFinder< InputContainer, VTX >::produce ( edm::Event event,
const edm::EventSetup es 
)
override

Definition at line 169 of file InclusiveVertexFinder.h.

References funct::abs(), pwdgSkimBPark_cfi::beamSpot, cms::cuda::bs, tthelpers::buildTT(), HLT_2024v14_cff::clusterizer, bsc_activity_cfg::clusters, gather_cfg::cout, DeadROC_duringRun::dir, VertexDistanceXY::distance(), VertexDistance3D::distance(), Measurement1D::error(), HLT_2024v14_cff::fitterRatio, HLT_2024v14_cff::fitterSigmacut, HLT_2024v14_cff::fitterTini, edm::EventSetup::getHandle(), mps_fire::i, edm::isFinite(), TransientVertex::isValid(), dqmiolumiharvest::j, HLT_2024v14_cff::maxNTracks, eostools::move(), funct::pow(), HLT_2024v14_cff::primaryVertices, Measurement1D::significance(), mathSSE::sqrt(), pfDeepBoostedJetPreprocessParams_cfi::sv, HLT_2024v14_cff::track, MinBiasPDSkim_cfg::trackFilter, DiMuonV_cfg::tracks, groupFilesInBlocks::tt, unit(), reco::BeamSpot::Unknown, svhelper::updateVertexTime(), HLT_2024v14_cff::useVertexReco, findQualityFiles::v, Measurement1D::value(), AdaptiveVertexFitter::vertex(), HLT_2024v14_cff::vertexMinAngleCosine, HLT_2024v14_cff::vertexMinDLen2DSig, HLT_2024v14_cff::vertexMinDLenSig, AlignmentTracksFromVertexSelector_cfi::vertices, L1BJetProducer_cff::vtx, and w().

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

◆ trackFilter()

template<class InputContainer , class VTX >
bool TemplatedInclusiveVertexFinder< InputContainer, VTX >::trackFilter ( const reco::Track track) const
private

Definition at line 159 of file InclusiveVertexFinder.h.

References createfilelist::int, reco_skim_cfg_mod::minHits, PV_cfg::minPt, and HLT_2024v14_cff::track.

159  {
160  if (track.hitPattern().numberOfValidHits() < (int)minHits)
161  return false;
162  if (track.pt() < minPt)
163  return false;
164 
165  return true;
166 }

Member Data Documentation

◆ clusterizer

template<class InputContainer , class VTX >
std::unique_ptr<TracksClusteringFromDisplacedSeed> TemplatedInclusiveVertexFinder< InputContainer, VTX >::clusterizer
private

◆ fitterRatio

template<class InputContainer , class VTX >
double TemplatedInclusiveVertexFinder< InputContainer, VTX >::fitterRatio
private

Definition at line 125 of file InclusiveVertexFinder.h.

◆ fitterSigmacut

template<class InputContainer , class VTX >
double TemplatedInclusiveVertexFinder< InputContainer, VTX >::fitterSigmacut
private

Definition at line 123 of file InclusiveVertexFinder.h.

◆ fitterTini

template<class InputContainer , class VTX >
double TemplatedInclusiveVertexFinder< InputContainer, VTX >::fitterTini
private

Definition at line 124 of file InclusiveVertexFinder.h.

◆ maxLIP

template<class InputContainer , class VTX >
double TemplatedInclusiveVertexFinder< InputContainer, VTX >::maxLIP
private

Definition at line 117 of file InclusiveVertexFinder.h.

◆ maxNTracks

template<class InputContainer , class VTX >
unsigned int TemplatedInclusiveVertexFinder< InputContainer, VTX >::maxNTracks
private

Definition at line 116 of file InclusiveVertexFinder.h.

◆ maxTimeSig

template<class InputContainer , class VTX >
double TemplatedInclusiveVertexFinder< InputContainer, VTX >::maxTimeSig
private

Definition at line 118 of file InclusiveVertexFinder.h.

◆ minHits

template<class InputContainer , class VTX >
unsigned int TemplatedInclusiveVertexFinder< InputContainer, VTX >::minHits
private

Definition at line 115 of file InclusiveVertexFinder.h.

◆ minPt

template<class InputContainer , class VTX >
double TemplatedInclusiveVertexFinder< InputContainer, VTX >::minPt
private

Definition at line 119 of file InclusiveVertexFinder.h.

◆ token_beamSpot

template<class InputContainer , class VTX >
edm::EDGetTokenT<reco::BeamSpot> TemplatedInclusiveVertexFinder< InputContainer, VTX >::token_beamSpot
private

◆ token_primaryVertex

template<class InputContainer , class VTX >
edm::EDGetTokenT<reco::VertexCollection> TemplatedInclusiveVertexFinder< InputContainer, VTX >::token_primaryVertex
private

◆ token_trackBuilder

template<class InputContainer , class VTX >
edm::ESGetToken<TransientTrackBuilder, TransientTrackRecord> TemplatedInclusiveVertexFinder< InputContainer, VTX >::token_trackBuilder
private

◆ token_tracks

template<class InputContainer , class VTX >
edm::EDGetTokenT<InputContainer> TemplatedInclusiveVertexFinder< InputContainer, VTX >::token_tracks
private

◆ useVertexFitter

template<class InputContainer , class VTX >
bool TemplatedInclusiveVertexFinder< InputContainer, VTX >::useVertexFitter
private

Definition at line 126 of file InclusiveVertexFinder.h.

◆ useVertexReco

template<class InputContainer , class VTX >
bool TemplatedInclusiveVertexFinder< InputContainer, VTX >::useVertexReco
private

Definition at line 127 of file InclusiveVertexFinder.h.

◆ vertexMinAngleCosine

template<class InputContainer , class VTX >
double TemplatedInclusiveVertexFinder< InputContainer, VTX >::vertexMinAngleCosine
private

Definition at line 120 of file InclusiveVertexFinder.h.

◆ vertexMinDLen2DSig

template<class InputContainer , class VTX >
double TemplatedInclusiveVertexFinder< InputContainer, VTX >::vertexMinDLen2DSig
private

Definition at line 121 of file InclusiveVertexFinder.h.

◆ vertexMinDLenSig

template<class InputContainer , class VTX >
double TemplatedInclusiveVertexFinder< InputContainer, VTX >::vertexMinDLenSig
private

Definition at line 122 of file InclusiveVertexFinder.h.

◆ vtxReco

template<class InputContainer , class VTX >
std::unique_ptr<VertexReconstructor> TemplatedInclusiveVertexFinder< InputContainer, VTX >::vtxReco
private

Definition at line 128 of file InclusiveVertexFinder.h.