CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Classes | Public Member Functions | Private Member Functions | Private Attributes
CombinedSVComputer Class Reference

#include <CombinedSVComputer.h>

Classes

struct  IterationRange
 

Public Member Functions

 CombinedSVComputer (const edm::ParameterSet &params)
 
reco::TaggingVariableList operator() (const reco::TrackIPTagInfo &ipInfo, const reco::SecondaryVertexTagInfo &svInfo) const
 

Private Member Functions

IterationRange flipIterate (int size, bool vertex) const
 
double flipValue (double value, bool vertex) const
 
const
reco::TrackIPTagInfo::TrackIPData
threshTrack (const reco::TrackIPTagInfo &trackIPTagInfo, const reco::TrackIPTagInfo::SortCriteria sort, const reco::Jet &jet, const GlobalPoint &pv) const
 

Private Attributes

double charmCut
 
double minTrackWeight
 
unsigned int pseudoMultiplicityMin
 
reco::V0Filter pseudoVertexV0Filter
 
reco::TrackIPTagInfo::SortCriteria sortCriterium
 
bool trackFlip
 
unsigned int trackMultiplicityMin
 
reco::TrackSelector trackNoDeltaRSelector
 
reco::V0Filter trackPairV0Filter
 
reco::TrackSelector trackPseudoSelector
 
reco::TrackSelector trackSelector
 
bool useTrackWeights
 
bool vertexFlip
 
bool vertexMassCorrection
 

Detailed Description

Definition at line 14 of file CombinedSVComputer.h.

Constructor & Destructor Documentation

CombinedSVComputer::CombinedSVComputer ( const edm::ParameterSet params)

Definition at line 51 of file CombinedSVComputer.cc.

51  :
52  trackFlip(params.getParameter<bool>("trackFlip")),
53  vertexFlip(params.getParameter<bool>("vertexFlip")),
54  charmCut(params.getParameter<double>("charmCut")),
56  trackSelector(params.getParameter<edm::ParameterSet>("trackSelection")),
58  trackPseudoSelector(params.getParameter<edm::ParameterSet>("trackPseudoSelection")),
59  pseudoMultiplicityMin(params.getParameter<unsigned int>("pseudoMultiplicityMin")),
60  trackMultiplicityMin(params.getParameter<unsigned int>("trackMultiplicityMin")),
61  minTrackWeight(params.getParameter<double>("minimumTrackWeight")),
62  useTrackWeights(params.getParameter<bool>("useTrackWeights")),
63  vertexMassCorrection(params.getParameter<bool>("correctVertexMass")),
64  pseudoVertexV0Filter(params.getParameter<edm::ParameterSet>("pseudoVertexV0Filter")),
65  trackPairV0Filter(params.getParameter<edm::ParameterSet>("trackPairV0Filter"))
66 {
67 }
reco::TrackSelector trackSelector
T getParameter(std::string const &) const
static edm::ParameterSet dropDeltaR(const edm::ParameterSet &pset)
reco::TrackIPTagInfo::SortCriteria getCriterium(const std::string &name)
Definition: TrackSorting.cc:11
unsigned int pseudoMultiplicityMin
reco::TrackSelector trackNoDeltaRSelector
reco::V0Filter pseudoVertexV0Filter
reco::TrackSelector trackPseudoSelector
reco::V0Filter trackPairV0Filter
reco::TrackIPTagInfo::SortCriteria sortCriterium
unsigned int trackMultiplicityMin

Member Function Documentation

CombinedSVComputer::IterationRange CombinedSVComputer::flipIterate ( int  size,
bool  vertex 
) const
inlineprivate

Definition at line 74 of file CombinedSVComputer.cc.

References CombinedSVComputer::IterationRange::begin, CombinedSVComputer::IterationRange::end, CombinedSVComputer::IterationRange::increment, findQualityFiles::size, trackFlip, and vertexFlip.

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

76 {
77  IterationRange range;
78  if (vertex ? vertexFlip : trackFlip) {
79  range.begin = size - 1;
80  range.end = -1;
81  range.increment = -1;
82  } else {
83  range.begin = 0;
84  range.end = size;
85  range.increment = +1;
86  }
87 
88  return range;
89 }
tuple size
Write out results.
double CombinedSVComputer::flipValue ( double  value,
bool  vertex 
) const
inlineprivate

Definition at line 69 of file CombinedSVComputer.cc.

References trackFlip, relativeConstraints::value, and vertexFlip.

Referenced by operator()().

70 {
71  return (vertex ? vertexFlip : trackFlip) ? -value : value;
72 }
TaggingVariableList CombinedSVComputer::operator() ( const reco::TrackIPTagInfo ipInfo,
const reco::SecondaryVertexTagInfo svInfo 
) const

Definition at line 140 of file CombinedSVComputer.cc.

References edm::RefVector< C, T, F >::begin(), reco::TrackIPTagInfo::TrackIPData::closestToJetAxis, data, deltaR(), dir, reco::TrackIPTagInfo::TrackIPData::distanceToJetAxis, edm::RefVector< C, T, F >::end(), reco::LeafCandidate::et(), etaRel(), reco::SecondaryVertexTagInfo::flightDirection(), reco::SecondaryVertexTagInfo::flightDistance(), reco::btau::flightDistance2dSig, reco::btau::flightDistance2dVal, reco::btau::flightDistance3dSig, reco::btau::flightDistance3dVal, flipIterate(), flipValue(), reco::Vertex::hasRefittedTracks(), i, customizeTrackingMonitorSeedNumber::idx, reco::TrackIPTagInfo::impactParameterData(), reco::TaggingVariableList::insert(), reco::TrackIPTagInfo::TrackIPData::ip2d, reco::TrackIPTagInfo::TrackIPData::ip3d, edm::Ref< C, T, F >::isNonnull(), j, reco::JTATagInfo::jet(), metsig::jet, reco::btau::jetEta, reco::btau::jetNSecondaryVertices, reco::btau::jetPt, PV3DBase< T, PVType, FrameType >::mag2(), minTrackWeight, reco::TrackBase::momentum(), reco::btag::Vertices::NoVertex, reco::SecondaryVertexTagInfo::nVertices(), convertSQLiteXML::ok, reco::TrackIPTagInfo::primaryVertex(), pseudoMultiplicityMin, reco::btag::Vertices::PseudoVertex, pseudoVertexV0Filter, range_for, reco::btag::Vertices::RecoVertex, reco::Vertex::refittedTrack(), reco::SecondaryVertexTagInfo::secondaryVertex(), reco::TrackIPTagInfo::selectedTracks(), Measurement1D::significance(), edm::RefVector< C, T, F >::size(), sortCriterium, reco::TrackIPTagInfo::sortedIndexes(), mathSSE::sqrt(), threshTrack(), reco::btau::trackDecayLenVal, reco::btau::trackDeltaR, reco::btau::trackEta, reco::btau::trackEtaRel, reco::btau::trackJetDistVal, reco::btau::trackMomentum, trackMultiplicityMin, trackPairV0Filter, reco::btau::trackPPar, reco::btau::trackPParRatio, trackPseudoSelector, reco::btau::trackPtRatio, reco::btau::trackPtRel, testEve_cfg::tracks, trackSelector, reco::btau::trackSip2dSig, reco::btau::trackSip2dSigAboveCharm, reco::btau::trackSip2dVal, reco::btau::trackSip3dSig, reco::btau::trackSip3dSigAboveCharm, reco::btau::trackSip3dVal, reco::btau::trackSumJetDeltaR, reco::btau::trackSumJetEtRatio, reco::SecondaryVertexTagInfo::trackWeight(), useTrackWeights, Measurement1D::value(), reco::btau::vertexCategory, reco::btau::vertexEnergyRatio, reco::btau::vertexJetDeltaR, reco::btau::vertexMass, vertexMassCorrection, reco::btau::vertexNTracks, reco::SecondaryVertexTagInfo::vertexTracks(), w(), PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

142 {
143  using namespace ROOT::Math;
144 
145  edm::RefToBase<Jet> jet = ipInfo.jet();
146  math::XYZVector jetDir = jet->momentum().Unit();
147  bool havePv = ipInfo.primaryVertex().isNonnull();
148  GlobalPoint pv;
149  if (havePv)
150  pv = GlobalPoint(ipInfo.primaryVertex()->x(),
151  ipInfo.primaryVertex()->y(),
152  ipInfo.primaryVertex()->z());
153 
155 
156  TaggingVariableList vars; // = ipInfo.taggingVariables();
157 
158  vars.insert(btau::jetPt, jet->pt(), true);
159  vars.insert(btau::jetEta, jet->eta(), true);
160 
161  if (ipInfo.selectedTracks().size() < trackMultiplicityMin)
162  return vars;
163 
164  TrackKinematics allKinematics;
165  TrackKinematics vertexKinematics;
166 
167  int vtx = -1;
168  unsigned int numberofvertextracks = 0;
169 
170  IterationRange range = flipIterate(svInfo.nVertices(), true);
171  range_for(i, range) {
172  if (vtx < 0)
173  vtx = i;
174 
175  numberofvertextracks = numberofvertextracks + (svInfo.secondaryVertex(i)).nTracks();
176 
177  const Vertex &vertex = svInfo.secondaryVertex(i);
178  bool hasRefittedTracks = vertex.hasRefittedTracks();
180  for(TrackRefVector::const_iterator track = tracks.begin();
181  track != tracks.end(); track++) {
182  double w = svInfo.trackWeight(i, *track);
183  if (w < minTrackWeight)
184  continue;
185  if (hasRefittedTracks) {
186  Track actualTrack =
187  vertex.refittedTrack(*track);
188  vertexKinematics.add(actualTrack, w);
189  vars.insert(btau::trackEtaRel, etaRel(jetDir,
190  actualTrack.momentum()), true);
191  } else {
192  vertexKinematics.add(**track, w);
193  vars.insert(btau::trackEtaRel, etaRel(jetDir,
194  (*track)->momentum()), true);
195  }
196  }
197  }
198 
199  if (vtx >= 0) {
200  vtxType = btag::Vertices::RecoVertex;
201 
202  vars.insert(btau::flightDistance2dVal,
203  flipValue(
204  svInfo.flightDistance(vtx, true).value(),
205  true),
206  true);
207  vars.insert(btau::flightDistance2dSig,
208  flipValue(
209  svInfo.flightDistance(vtx, true).significance(),
210  true),
211  true);
212  vars.insert(btau::flightDistance3dVal,
213  flipValue(
214  svInfo.flightDistance(vtx, false).value(),
215  true),
216  true);
217  vars.insert(btau::flightDistance3dSig,
218  flipValue(
219  svInfo.flightDistance(vtx, false).significance(),
220  true),
221  true);
222  vars.insert(btau::vertexJetDeltaR,
223  Geom::deltaR(svInfo.flightDirection(vtx), jetDir),true);
224  vars.insert(btau::jetNSecondaryVertices, svInfo.nVertices(), true);
225  vars.insert(btau::vertexNTracks, numberofvertextracks, true);
226  }
227 
228  std::vector<std::size_t> indices = ipInfo.sortedIndexes(sortCriterium);
229  const std::vector<TrackIPTagInfo::TrackIPData> &ipData =
230  ipInfo.impactParameterData();
231  const edm::RefVector<TrackCollection> &tracks =
232  ipInfo.selectedTracks();
233  std::vector<TrackRef> pseudoVertexTracks;
234 
235  TrackRef trackPairV0Test[2];
236  range = flipIterate(indices.size(), false);
237  range_for(i, range) {
238  std::size_t idx = indices[i];
239  const TrackIPTagInfo::TrackIPData &data = ipData[idx];
240  const TrackRef &trackRef = tracks[idx];
241  const Track &track = *trackRef;
242 
243  // filter track
244 
245  if (!trackSelector(track, data, *jet, pv))
246  continue;
247 
248  // add track to kinematics for all tracks in jet
249 
250  allKinematics.add(track);
251 
252  // if no vertex was reconstructed, attempt pseudo vertex
253 
254  if (vtxType == btag::Vertices::NoVertex &&
255  trackPseudoSelector(track, data, *jet, pv)) {
256  pseudoVertexTracks.push_back(trackRef);
257  vertexKinematics.add(track);
258  }
259 
260  // check against all other tracks for V0 track pairs
261 
262  trackPairV0Test[0] = tracks[idx];
263  bool ok = true;
264  range_for(j, range) {
265  if (i == j)
266  continue;
267 
268  std::size_t pairIdx = indices[j];
269  const TrackIPTagInfo::TrackIPData &pairTrackData =
270  ipData[pairIdx];
271  const TrackRef &pairTrackRef = tracks[pairIdx];
272  const Track &pairTrack = *pairTrackRef;
273 
274  if (!trackSelector(pairTrack, pairTrackData, *jet, pv))
275  continue;
276 
277  trackPairV0Test[1] = pairTrackRef;
278  if (!trackPairV0Filter(trackPairV0Test, 2)) {
279  ok = false;
280  break;
281  }
282  }
283  if (!ok)
284  continue;
285 
286  // add track variables
287 
288  math::XYZVector trackMom = track.momentum();
289  double trackMag = std::sqrt(trackMom.Mag2());
290 
291  vars.insert(btau::trackSip3dVal,
292  flipValue(data.ip3d.value(), false), true);
293  vars.insert(btau::trackSip3dSig,
294  flipValue(data.ip3d.significance(), false), true);
295  vars.insert(btau::trackSip2dVal,
296  flipValue(data.ip2d.value(), false), true);
297  vars.insert(btau::trackSip2dSig,
298  flipValue(data.ip2d.significance(), false), true);
299  vars.insert(btau::trackJetDistVal, data.distanceToJetAxis.value(), true);
300 // vars.insert(btau::trackJetDistSig, data.distanceToJetAxis.significance(), true);
301 // vars.insert(btau::trackFirstTrackDist, data.distanceToFirstTrack, true);
302 // vars.insert(btau::trackGhostTrackVal, data.distanceToGhostTrack.value(), true);
303 // vars.insert(btau::trackGhostTrackSig, data.distanceToGhostTrack.significance(), true);
304  vars.insert(btau::trackDecayLenVal, havePv ? (data.closestToJetAxis - pv).mag() : -1.0, true);
305 
306  vars.insert(btau::trackMomentum, trackMag, true);
307  vars.insert(btau::trackEta, trackMom.Eta(), true);
308 
309  vars.insert(btau::trackPtRel, VectorUtil::Perp(trackMom, jetDir), true);
310  vars.insert(btau::trackPPar, jetDir.Dot(trackMom), true);
311  vars.insert(btau::trackDeltaR, VectorUtil::DeltaR(trackMom, jetDir), true);
312  vars.insert(btau::trackPtRatio, VectorUtil::Perp(trackMom, jetDir) / trackMag, true);
313  vars.insert(btau::trackPParRatio, jetDir.Dot(trackMom) / trackMag, true);
314  }
315 
316  if (vtxType == btag::Vertices::NoVertex &&
317  vertexKinematics.numberOfTracks() >= pseudoMultiplicityMin &&
318  pseudoVertexV0Filter(pseudoVertexTracks)) {
320  for(std::vector<TrackRef>::const_iterator track =
321  pseudoVertexTracks.begin();
322  track != pseudoVertexTracks.end(); ++track)
323  vars.insert(btau::trackEtaRel, etaRel(jetDir,
324  (*track)->momentum()), true);
325  }
326 
327  vars.insert(btau::vertexCategory, vtxType, true);
328  vars.insert(btau::trackSumJetDeltaR,
329  VectorUtil::DeltaR(allKinematics.vectorSum(), jetDir), true);
330  vars.insert(btau::trackSumJetEtRatio,
331  allKinematics.vectorSum().Et() / ipInfo.jet()->et(), true);
332  vars.insert(btau::trackSip3dSigAboveCharm,
333  flipValue(
334  threshTrack(ipInfo, TrackIPTagInfo::IP3DSig, *jet, pv)
335  .ip3d.significance(),
336  false),
337  true);
338  vars.insert(btau::trackSip2dSigAboveCharm,
339  flipValue(
340  threshTrack(ipInfo, TrackIPTagInfo::IP2DSig, *jet, pv)
341  .ip2d.significance(),
342  false),
343  true);
344 
345  if (vtxType != btag::Vertices::NoVertex) {
347  ? allKinematics.weightedVectorSum()
348  : allKinematics.vectorSum();
350  ? vertexKinematics.weightedVectorSum()
351  : vertexKinematics.vectorSum();
352 
353  if (vtxType != btag::Vertices::RecoVertex) {
354  vars.insert(btau::vertexNTracks,
355  vertexKinematics.numberOfTracks(), true);
356  vars.insert(btau::vertexJetDeltaR,
357  VectorUtil::DeltaR(vertexSum, jetDir), true);
358  }
359 
360  double vertexMass = vertexSum.M();
361  if (vtxType == btag::Vertices::RecoVertex &&
363  GlobalVector dir = svInfo.flightDirection(vtx);
364  double vertexPt2 =
365  math::XYZVector(dir.x(), dir.y(), dir.z()).
366  Cross(vertexSum).Mag2() / dir.mag2();
367  vertexMass = std::sqrt(vertexMass * vertexMass +
368  vertexPt2) + std::sqrt(vertexPt2);
369  }
370  vars.insert(btau::vertexMass, vertexMass, true);
371  if (allKinematics.numberOfTracks())
372  vars.insert(btau::vertexEnergyRatio,
373  vertexSum.E() / allSum.E(), true);
374  else
375  vars.insert(btau::vertexEnergyRatio, 1, true);
376  }
377 
378  vars.finalize();
379 
380  return vars;
381 }
reco::TrackSelector trackSelector
virtual double et() const GCC11_FINAL
transverse energy
int i
Definition: DBlmapReader.cc:9
T mag2() const
Definition: PV3DBase.h:66
const Vector & momentum() const
track momentum vector
Definition: TrackBase.h:148
Measurement1D flightDistance(unsigned int index, bool in2d=false) const
unsigned int pseudoMultiplicityMin
const edm::Ref< VertexCollection > & primaryVertex() const
Track refittedTrack(const TrackBaseRef &track) const
bool hasRefittedTracks() const
Checks whether refitted tracks are stored.
Definition: Vertex.h:134
reco::V0Filter pseudoVertexV0Filter
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
#define range_for(i, x)
T y() const
Definition: PV3DBase.h:63
const_iterator end() const
Termination of iteration.
Definition: RefVector.h:249
const_iterator begin() const
Initialize an iterator over the RefVector.
Definition: RefVector.h:244
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:250
const GlobalVector & flightDirection(unsigned int index) const
const Vertex & secondaryVertex(unsigned int index) const
reco::TrackSelector trackPseudoSelector
T sqrt(T t)
Definition: SSEVec.h:48
T z() const
Definition: PV3DBase.h:64
int j
Definition: DBlmapReader.cc:9
const reco::TrackIPTagInfo::TrackIPData & threshTrack(const reco::TrackIPTagInfo &trackIPTagInfo, const reco::TrackIPTagInfo::SortCriteria sort, const reco::Jet &jet, const GlobalPoint &pv) const
float trackWeight(unsigned int svIndex, unsigned int trackindex) const
IterationRange flipIterate(int size, bool vertex) const
virtual edm::RefToBase< Jet > jet(void) const
returns a polymorphic reference to the tagged jet
Definition: JTATagInfo.h:20
TrackRefVector vertexTracks() const
reco::V0Filter trackPairV0Filter
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
std::vector< size_t > sortedIndexes(SortCriteria mode=IP3DSig) const
static double etaRel(const math::XYZVector &dir, const math::XYZVector &track)
double significance() const
Definition: Measurement1D.h:32
tuple tracks
Definition: testEve_cfg.py:39
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
tuple idx
DEBUGGING if hasattr(process,&quot;trackMonIterativeTracking2012&quot;): print &quot;trackMonIterativeTracking2012 D...
double flipValue(double value, bool vertex) const
double value() const
Definition: Measurement1D.h:28
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:82
const std::vector< TrackIPData > & impactParameterData() const
size_type size() const
Size of the RefVector.
Definition: RefVector.h:89
T w() const
dbl *** dir
Definition: mlp_gen.cc:35
const edm::RefVector< TrackCollection > & selectedTracks() const
T x() const
Definition: PV3DBase.h:62
reco::TrackIPTagInfo::SortCriteria sortCriterium
void insert(const TaggingVariable &variable, bool delayed=false)
unsigned int trackMultiplicityMin
const TrackIPTagInfo::TrackIPData & CombinedSVComputer::threshTrack ( const reco::TrackIPTagInfo trackIPTagInfo,
const reco::TrackIPTagInfo::SortCriteria  sort,
const reco::Jet jet,
const GlobalPoint pv 
) const
private

Definition at line 92 of file CombinedSVComputer.cc.

References charmCut, data, flipIterate(), i, customizeTrackingMonitorSeedNumber::idx, reco::TrackIPTagInfo::impactParameterData(), range_for, reco::TrackIPTagInfo::selectedTracks(), reco::TrackIPTagInfo::sortedIndexes(), trackNoDeltaRSelector, and testEve_cfg::tracks.

Referenced by operator()().

96 {
98  trackIPTagInfo.selectedTracks();
99  const std::vector<TrackIPTagInfo::TrackIPData> &ipData =
100  trackIPTagInfo.impactParameterData();
101  std::vector<std::size_t> indices = trackIPTagInfo.sortedIndexes(sort);
102 
103  IterationRange range = flipIterate(indices.size(), false);
104  TrackKinematics kin;
105  range_for(i, range) {
106  std::size_t idx = indices[i];
107  const TrackIPTagInfo::TrackIPData &data = ipData[idx];
108  const Track &track = *tracks[idx];
109 
110  if (!trackNoDeltaRSelector(track, data, jet, pv))
111  continue;
112 
113  kin.add(track);
114  if (kin.vectorSum().M() > charmCut)
115  return data;
116  }
117 
118  static const TrackIPTagInfo::TrackIPData dummy = {
119  GlobalPoint(),
120  GlobalPoint(),
121  Measurement1D(-1.0, 1.0),
122  Measurement1D(-1.0, 1.0),
123  Measurement1D(-1.0, 1.0),
124  Measurement1D(-1.0, 1.0),
125  0.
126  };
127  return dummy;
128 }
int i
Definition: DBlmapReader.cc:9
reco::TrackSelector trackNoDeltaRSelector
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
#define range_for(i, x)
IterationRange flipIterate(int size, bool vertex) const
std::vector< size_t > sortedIndexes(SortCriteria mode=IP3DSig) const
tuple tracks
Definition: testEve_cfg.py:39
tuple idx
DEBUGGING if hasattr(process,&quot;trackMonIterativeTracking2012&quot;): print &quot;trackMonIterativeTracking2012 D...
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:82
const std::vector< TrackIPData > & impactParameterData() const
const edm::RefVector< TrackCollection > & selectedTracks() const

Member Data Documentation

double CombinedSVComputer::charmCut
private

Definition at line 36 of file CombinedSVComputer.h.

Referenced by threshTrack().

double CombinedSVComputer::minTrackWeight
private

Definition at line 43 of file CombinedSVComputer.h.

Referenced by operator()().

unsigned int CombinedSVComputer::pseudoMultiplicityMin
private

Definition at line 41 of file CombinedSVComputer.h.

Referenced by operator()().

reco::V0Filter CombinedSVComputer::pseudoVertexV0Filter
private

Definition at line 46 of file CombinedSVComputer.h.

Referenced by operator()().

reco::TrackIPTagInfo::SortCriteria CombinedSVComputer::sortCriterium
private

Definition at line 37 of file CombinedSVComputer.h.

Referenced by operator()().

bool CombinedSVComputer::trackFlip
private

Definition at line 34 of file CombinedSVComputer.h.

Referenced by flipIterate(), and flipValue().

unsigned int CombinedSVComputer::trackMultiplicityMin
private

Definition at line 42 of file CombinedSVComputer.h.

Referenced by operator()().

reco::TrackSelector CombinedSVComputer::trackNoDeltaRSelector
private

Definition at line 39 of file CombinedSVComputer.h.

Referenced by threshTrack().

reco::V0Filter CombinedSVComputer::trackPairV0Filter
private

Definition at line 47 of file CombinedSVComputer.h.

Referenced by operator()().

reco::TrackSelector CombinedSVComputer::trackPseudoSelector
private

Definition at line 40 of file CombinedSVComputer.h.

Referenced by operator()().

reco::TrackSelector CombinedSVComputer::trackSelector
private

Definition at line 38 of file CombinedSVComputer.h.

Referenced by operator()().

bool CombinedSVComputer::useTrackWeights
private

Definition at line 44 of file CombinedSVComputer.h.

Referenced by operator()().

bool CombinedSVComputer::vertexFlip
private

Definition at line 35 of file CombinedSVComputer.h.

Referenced by flipIterate(), and flipValue().

bool CombinedSVComputer::vertexMassCorrection
private

Definition at line 45 of file CombinedSVComputer.h.

Referenced by operator()().