CMS 3D CMS Logo

List of all members | Public Member Functions | Protected Member Functions | Private Member Functions | Private Attributes
BPHDecayVertex Class Reference

#include <BPHDecayVertex.h>

Inheritance diagram for BPHDecayVertex:
BPHDecayMomentum BPHKinematicFit BPHPlusMinusVertex BPHRecoCandidate BPHPlusMinusCandidate BPHPlusMinusCandidate

Public Member Functions

 BPHDecayVertex (const BPHDecayVertex &x)=delete
 
const BPHEventSetupWrappergetEventSetup () const
 retrieve EventSetup More...
 
char getTMode (const reco::Candidate *cand) const
 get Track mode for a daughter More...
 
const reco::TrackgetTrack (const reco::Candidate *cand) const
 get Track for a daughter More...
 
const std::string & getTrackSearchList (const reco::Candidate *cand) const
 retrieve track search list More...
 
reco::TransientTrackgetTransientTrack (const reco::Candidate *cand) const
 get TransientTrack for a daughter More...
 
BPHDecayVertexoperator= (const BPHDecayVertex &x)=delete
 
const std::vector< const reco::Track * > & tracks () const
 get list of Tracks More...
 
const std::vector< reco::TransientTrack > & transientTracks () const
 get list of TransientTracks More...
 
virtual bool validTracks () const
 check for valid reconstructed vertex More...
 
virtual bool validVertex () const
 
virtual const reco::Vertexvertex (VertexFitter< 5 > *fitter=nullptr, const reco::BeamSpot *bs=nullptr, const GlobalPoint *priorPos=nullptr, const GlobalError *priorError=nullptr) const
 get reconstructed vertex More...
 
 ~BPHDecayVertex () override
 
- Public Member Functions inherited from BPHDecayMomentum
 BPHDecayMomentum (const BPHDecayMomentum &x)=delete
 
const std::map< std::string, BPHRecoConstCandPtr > & compMap () const
 
virtual const std::vector< std::string > & compNames () const
 
virtual const pat::CompositeCandidatecomposite () const
 get a composite by the simple sum of simple particles More...
 
virtual const std::vector< BPHRecoConstCandPtr > & daughComp () const
 
virtual const std::vector< const reco::Candidate * > & daughFull () const
 
virtual const std::vector< const reco::Candidate * > & daughters () const
 
const std::map< std::string, const reco::Candidate * > & daugMap () const
 
virtual const std::vector< std::string > & daugNames () const
 
virtual BPHRecoConstCandPtr getComp (const std::string &name) const
 
virtual const reco::CandidategetDaug (const std::string &name) const
 
BPHDecayMomentumoperator= (const BPHDecayMomentum &x)=delete
 
virtual const reco::CandidateoriginalReco (const reco::Candidate *daug) const
 get the original particle from the clone More...
 
virtual ~BPHDecayMomentum ()
 

Protected Member Functions

virtual void addV (const std::string &name, const reco::Candidate *daug, const std::string &searchList, double mass)
 
virtual void addV (const std::string &name, const BPHRecoConstCandPtr &comp)
 
 BPHDecayVertex (const BPHEventSetupWrapper *es, int daugNum=2, int compNum=2)
 
 BPHDecayVertex (const BPHDecayVertex *ptr, const BPHEventSetupWrapper *es)
 
void setNotUpdated () const override
 
- Protected Member Functions inherited from BPHDecayMomentum
virtual void addP (const std::string &name, const reco::Candidate *daug, double mass=-1.0)
 
virtual void addP (const std::string &name, const BPHRecoConstCandPtr &comp)
 
 BPHDecayMomentum (int daugNum=2, int compNum=2)
 
 BPHDecayMomentum (const std::map< std::string, Component > &daugMap, int compNum=2)
 
 BPHDecayMomentum (const std::map< std::string, Component > &daugMap, const std::map< std::string, BPHRecoConstCandPtr > compMap)
 
const std::vector< Component > & componentList () const
 
virtual void fill (BPHRecoCandidate *ptr, int level) const =0
 

Private Member Functions

virtual void fitVertex (VertexFitter< 5 > *fitter, const reco::BeamSpot *bs, const GlobalPoint *priorPos, const GlobalError *priorError) const
 
virtual void fTracks () const
 
virtual void fTTracks () const
 

Private Attributes

const BPHEventSetupWrapperevSetup
 
reco::Vertex fittedVertex
 
bool oldTracks
 
bool oldTTracks
 
bool oldVertex
 
std::vector< const reco::Track * > rTracks
 
const reco::BeamSpotsavedBS
 
VertexFitter< 5 > * savedFitter
 
const GlobalErrorsavedPE
 
const GlobalPointsavedPP
 
std::map< const reco::Candidate *, std::string > searchMap
 
std::map< const reco::Candidate *, const reco::Track * > tkMap
 
std::map< const reco::Candidate *, char > tmMap
 
std::vector< reco::TransientTracktrTracks
 
std::map< const reco::Candidate *, reco::TransientTrack * > ttMap
 
bool validTks
 

Detailed Description

Description: Mid-level base class to reconstruct decay vertex

Author
Paolo Ronchese INFN Padova

Definition at line 45 of file BPHDecayVertex.h.

Constructor & Destructor Documentation

◆ BPHDecayVertex() [1/3]

BPHDecayVertex::BPHDecayVertex ( const BPHDecayVertex x)
delete

Constructors are protected this object can exist only as part of a derived class

◆ ~BPHDecayVertex()

BPHDecayVertex::~BPHDecayVertex ( )
override

Destructor

Definition at line 90 of file BPHDecayVertex.cc.

References evSetup.

90 { delete evSetup; }
const BPHEventSetupWrapper * evSetup

◆ BPHDecayVertex() [2/3]

BPHDecayVertex::BPHDecayVertex ( const BPHEventSetupWrapper es,
int  daugNum = 2,
int  compNum = 2 
)
protected

Definition at line 43 of file BPHDecayVertex.cc.

44  : BPHDecayMomentum(daugNum, compNum),
46  oldTracks(true),
47  oldTTracks(true),
48  oldVertex(true),
49  validTks(false),
50  savedFitter(nullptr),
51  savedBS(nullptr),
52  savedPP(nullptr),
53  savedPE(nullptr) {}
const BPHEventSetupWrapper * evSetup
VertexFitter< 5 > * savedFitter
BPHDecayMomentum(const BPHDecayMomentum &x)=delete
const GlobalError * savedPE
const GlobalPoint * savedPP
const reco::BeamSpot * savedBS

◆ BPHDecayVertex() [3/3]

BPHDecayVertex::BPHDecayVertex ( const BPHDecayVertex ptr,
const BPHEventSetupWrapper es 
)
protected

Definition at line 55 of file BPHDecayVertex.cc.

References DummyCfis::c, BPHDecayMomentum::componentList(), BPHDecayMomentum::daughComp(), BPHDecayMomentum::daughters(), BPHDecayMomentum::dMap, mps_fire::i, dqmiolumiharvest::j, visualization-live-secondInstance_cfg::m, dqmiodumpmetadata::n, BPHDecayMomentum::originalReco(), and searchMap.

56  : evSetup(new BPHEventSetupWrapper(es)),
57  oldTracks(true),
58  oldTTracks(true),
59  oldVertex(true),
60  validTks(false),
61  savedFitter(nullptr),
62  savedBS(nullptr),
63  savedPP(nullptr),
64  savedPE(nullptr) {
65  map<const reco::Candidate*, const reco::Candidate*> iMap;
66  const vector<const reco::Candidate*>& daug = daughters();
67  const vector<Component>& list = ptr->componentList();
68  int i;
69  int n = daug.size();
70  for (i = 0; i < n; ++i) {
71  const reco::Candidate* cand = daug[i];
72  iMap[originalReco(cand)] = cand;
73  }
74  for (i = 0; i < n; ++i) {
75  const Component& c = list[i];
76  searchMap[iMap[c.cand]] = c.searchList;
77  }
78  const vector<BPHRecoConstCandPtr>& dComp = daughComp();
79  int j;
80  int m = dComp.size();
81  for (j = 0; j < m; ++j) {
82  const map<const reco::Candidate*, string>& dMap = dComp[j]->searchMap;
83  searchMap.insert(dMap.begin(), dMap.end());
84  }
85 }
virtual const reco::Candidate * originalReco(const reco::Candidate *daug) const
get the original particle from the clone
const BPHEventSetupWrapper * evSetup
VertexFitter< 5 > * savedFitter
const std::vector< Component > & componentList() const
std::map< std::string, const reco::Candidate * > dMap
std::map< const reco::Candidate *, std::string > searchMap
const GlobalError * savedPE
virtual const std::vector< BPHRecoConstCandPtr > & daughComp() const
const GlobalPoint * savedPP
virtual const std::vector< const reco::Candidate * > & daughters() const
const reco::BeamSpot * savedBS

Member Function Documentation

◆ addV() [1/2]

virtual void BPHDecayVertex::addV ( const std::string &  name,
const reco::Candidate daug,
const std::string &  searchList,
double  mass 
)
protectedvirtual

◆ addV() [2/2]

virtual void BPHDecayVertex::addV ( const std::string &  name,
const BPHRecoConstCandPtr comp 
)
protectedvirtual

◆ fitVertex()

void BPHDecayVertex::fitVertex ( VertexFitter< 5 > *  fitter,
const reco::BeamSpot bs,
const GlobalPoint priorPos,
const GlobalError priorError 
) const
privatevirtual

Definition at line 254 of file BPHDecayVertex.cc.

References cms::cuda::bs, cppFunctionSkipper::exception, fittedVertex, fTTracks(), oldTTracks, oldVertex, savedBS, savedFitter, savedPE, savedPP, trTracks, and VertexFitter< N >::vertex().

Referenced by vertex().

257  {
258  oldVertex = false;
259  savedFitter = fitter;
260  savedBS = bs;
261  savedPP = priorPos;
262  savedPE = priorError;
263  if (oldTTracks)
264  fTTracks();
265  if (trTracks.size() < 2)
266  return;
267  try {
268  if (bs == nullptr) {
269  if (priorPos == nullptr) {
270  TransientVertex tv = fitter->vertex(trTracks);
271  fittedVertex = tv;
272  } else {
273  if (priorError == nullptr) {
274  TransientVertex tv = fitter->vertex(trTracks, *priorPos);
275  fittedVertex = tv;
276  } else {
277  TransientVertex tv = fitter->vertex(trTracks, *priorPos, *priorError);
278  fittedVertex = tv;
279  }
280  }
281  } else {
282  TransientVertex tv = fitter->vertex(trTracks, *bs);
283  fittedVertex = tv;
284  }
285  } catch (std::exception const&) {
286  reco::Vertex tv;
287  fittedVertex = tv;
288  edm::LogPrint("FitFailed") << "BPHDecayVertex::fitVertex: "
289  << "vertex fit failed";
290  }
291  return;
292 }
VertexFitter< 5 > * savedFitter
std::vector< reco::TransientTrack > trTracks
virtual CachingVertex< N > vertex(const std::vector< reco::TransientTrack > &tracks) const =0
Log< level::Warning, true > LogPrint
reco::Vertex fittedVertex
const GlobalError * savedPE
const GlobalPoint * savedPP
virtual void fTTracks() const
const reco::BeamSpot * savedBS

◆ fTracks()

void BPHDecayVertex::fTracks ( ) const
privatevirtual

Definition at line 194 of file BPHDecayVertex.cc.

References reco::Candidate::charge(), BPHDecayMomentum::daughFull(), BPHTrackReference::getTrack(), dqmiodumpmetadata::n, oldTracks, oldTTracks, BPHDecayMomentum::originalReco(), rTracks, newFWLiteAna::searchList, searchMap, tkMap, tmMap, cmsswSequenceInfo::tp, trTracks, and validTks.

Referenced by fTTracks(), getTMode(), getTrack(), tracks(), and validTracks().

194  {
195  oldTTracks = true;
196  rTracks.clear();
197  tkMap.clear();
198  tmMap.clear();
199  const vector<const reco::Candidate*>& dL = daughFull();
200  int n = dL.size();
201  trTracks.reserve(n);
202  validTks = true;
203  while (n--) {
204  const reco::Candidate* rp = dL[n];
205  tkMap[rp] = nullptr;
206  tmMap[rp] = '.';
207  if (!rp->charge())
208  continue;
209  const char* searchList = "cfhp";
210  char usedMode;
211  map<const reco::Candidate*, string>::const_iterator iter = searchMap.find(rp);
212  if (iter != searchMap.end())
213  searchList = iter->second.c_str();
214  const reco::Track* tp = tkMap[rp] = BPHTrackReference::getTrack(*originalReco(rp), searchList, &usedMode);
215  if (tp == nullptr) {
216  edm::LogPrint("DataNotFound") << "BPHDecayVertex::tTracks: "
217  << "no track for reco::(PF)Candidate";
218  validTks = false;
219  continue;
220  }
221  rTracks.push_back(tp);
222  tmMap[rp] = usedMode;
223  }
224  oldTracks = false;
225  return;
226 }
virtual const reco::Candidate * originalReco(const reco::Candidate *daug) const
get the original particle from the clone
std::vector< reco::TransientTrack > trTracks
std::map< const reco::Candidate *, std::string > searchMap
virtual const std::vector< const reco::Candidate * > & daughFull() const
static const reco::Track * getTrack(const reco::Candidate &rc, const char *modeList="cfhbpmnigset", char *modeFlag=nullptr)
virtual int charge() const =0
electric charge
Log< level::Warning, true > LogPrint
std::map< const reco::Candidate *, const reco::Track * > tkMap
std::map< const reco::Candidate *, char > tmMap
std::vector< const reco::Track * > rTracks

◆ fTTracks()

void BPHDecayVertex::fTTracks ( ) const
privatevirtual

Definition at line 228 of file BPHDecayVertex.cc.

References TransientTrackBuilder::build(), BPHDecayMomentum::daughFull(), SiStripBadComponentsDQMServiceTemplate_cfg::ep, evSetup, fTracks(), BPHEventSetupWrapper::get(), dqmiodumpmetadata::n, oldTracks, oldTTracks, tkMap, unpackBuffers-CaloStage2::token, cmsswSequenceInfo::tp, BPHRecoCandidate::transientTrackBuilder, trTracks, and ttMap.

Referenced by fitVertex(), getTransientTrack(), and transientTracks().

228  {
229  if (oldTracks)
230  fTracks();
231  trTracks.clear();
234  const edm::EventSetup* ep = evSetup->get();
236  token->get(*ep, ttB);
237  ttMap.clear();
238  const vector<const reco::Candidate*>& dL = daughFull();
239  int n = dL.size();
240  trTracks.reserve(n);
241  while (n--) {
242  const reco::Candidate* rp = dL[n];
243  ttMap[rp] = nullptr;
244  map<const reco::Candidate*, const reco::Track*>::const_iterator iter = tkMap.find(rp);
245  const reco::Track* tp = iter->second;
246  if (tp == nullptr)
247  continue;
248  trTracks.push_back(ttB->build(tp));
249  ttMap[rp] = &trTracks.back();
250  }
251  oldTTracks = false;
252 }
const BPHEventSetupWrapper * evSetup
std::vector< reco::TransientTrack > trTracks
virtual const std::vector< const reco::Candidate * > & daughFull() const
reco::TransientTrack build(const reco::Track *p) const
std::map< const reco::Candidate *, reco::TransientTrack * > ttMap
virtual void fTracks() const
std::map< const reco::Candidate *, const reco::Track * > tkMap
const edm::EventSetup * get() const

◆ getEventSetup()

const BPHEventSetupWrapper * BPHDecayVertex::getEventSetup ( ) const

retrieve EventSetup

Definition at line 164 of file BPHDecayVertex.cc.

References evSetup.

Referenced by BPHPlusMinusCandidate::clone(), and BPHRecoCandidate::clone().

164 { return evSetup; }
const BPHEventSetupWrapper * evSetup

◆ getTMode()

char BPHDecayVertex::getTMode ( const reco::Candidate cand) const

get Track mode for a daughter

Definition at line 141 of file BPHDecayVertex.cc.

References fTracks(), oldTracks, and tmMap.

141  {
142  if (oldTracks)
143  fTracks();
144  map<const reco::Candidate*, char>::const_iterator iter = tmMap.find(cand);
145  map<const reco::Candidate*, char>::const_iterator iend = tmMap.end();
146  return (iter != iend ? iter->second : '.');
147 }
virtual void fTracks() const
std::map< const reco::Candidate *, char > tmMap

◆ getTrack()

const reco::Track * BPHDecayVertex::getTrack ( const reco::Candidate cand) const

get Track for a daughter

Definition at line 133 of file BPHDecayVertex.cc.

References fTracks(), oldTracks, and tkMap.

133  {
134  if (oldTracks)
135  fTracks();
136  map<const reco::Candidate*, const reco::Track*>::const_iterator iter = tkMap.find(cand);
137  map<const reco::Candidate*, const reco::Track*>::const_iterator iend = tkMap.end();
138  return (iter != iend ? iter->second : nullptr);
139 }
virtual void fTracks() const
std::map< const reco::Candidate *, const reco::Track * > tkMap

◆ getTrackSearchList()

const string & BPHDecayVertex::getTrackSearchList ( const reco::Candidate cand) const

retrieve track search list

Definition at line 166 of file BPHDecayVertex.cc.

References searchMap.

Referenced by BPHRecoCandidate::fill().

166  {
167  static const string dum = "";
168  map<const reco::Candidate*, string>::const_iterator iter = searchMap.find(cand);
169  if (iter != searchMap.end())
170  return iter->second;
171  return dum;
172 }
std::map< const reco::Candidate *, std::string > searchMap

◆ getTransientTrack()

reco::TransientTrack * BPHDecayVertex::getTransientTrack ( const reco::Candidate cand) const

get TransientTrack for a daughter

Definition at line 155 of file BPHDecayVertex.cc.

References fTTracks(), oldTTracks, and ttMap.

Referenced by BPHKinematicFit::addParticles().

155  {
156  if (oldTTracks)
157  fTTracks();
158  map<const reco::Candidate*, reco::TransientTrack*>::const_iterator iter = ttMap.find(cand);
159  map<const reco::Candidate*, reco::TransientTrack*>::const_iterator iend = ttMap.end();
160  return (iter != iend ? iter->second : nullptr);
161 }
std::map< const reco::Candidate *, reco::TransientTrack * > ttMap
virtual void fTTracks() const

◆ operator=()

BPHDecayVertex& BPHDecayVertex::operator= ( const BPHDecayVertex x)
delete

◆ setNotUpdated()

void BPHDecayVertex::setNotUpdated ( ) const
overrideprotectedvirtual

Reimplemented from BPHDecayMomentum.

Reimplemented in BPHKinematicFit, BPHPlusMinusCandidate, and BPHPlusMinusVertex.

Definition at line 187 of file BPHDecayVertex.cc.

References oldTracks, oldVertex, BPHDecayMomentum::setNotUpdated(), and validTks.

Referenced by BPHPlusMinusVertex::setNotUpdated(), and BPHKinematicFit::setNotUpdated().

187  {
189  oldTracks = oldVertex = true;
190  validTks = false;
191  return;
192 }
virtual void setNotUpdated() const

◆ tracks()

const vector< const reco::Track * > & BPHDecayVertex::tracks ( void  ) const

get list of Tracks

Definition at line 127 of file BPHDecayVertex.cc.

References fTracks(), oldTracks, and rTracks.

127  {
128  if (oldTracks)
129  fTracks();
130  return rTracks;
131 }
virtual void fTracks() const
std::vector< const reco::Track * > rTracks

◆ transientTracks()

const vector< reco::TransientTrack > & BPHDecayVertex::transientTracks ( ) const

get list of TransientTracks

Definition at line 149 of file BPHDecayVertex.cc.

References fTTracks(), oldTTracks, and trTracks.

Referenced by BPHPlusMinusVertex::computeApp().

149  {
150  if (oldTTracks)
151  fTTracks();
152  return trTracks;
153 }
std::vector< reco::TransientTrack > trTracks
virtual void fTTracks() const

◆ validTracks()

bool BPHDecayVertex::validTracks ( ) const
virtual

check for valid reconstructed vertex

Operations

Definition at line 95 of file BPHDecayVertex.cc.

References fTracks(), oldTracks, and validTks.

95  {
96  if (oldTracks)
97  fTracks();
98  return validTks;
99 }
virtual void fTracks() const

◆ validVertex()

bool BPHDecayVertex::validVertex ( ) const
virtual

Definition at line 101 of file BPHDecayVertex.cc.

References fittedVertex, reco::Vertex::isValid(), validTks, and vertex().

101  {
102  vertex();
103  return validTks && fittedVertex.isValid();
104 }
reco::Vertex fittedVertex
virtual const reco::Vertex & vertex(VertexFitter< 5 > *fitter=nullptr, const reco::BeamSpot *bs=nullptr, const GlobalPoint *priorPos=nullptr, const GlobalError *priorError=nullptr) const
get reconstructed vertex
bool isValid() const
Tells whether the vertex is valid.
Definition: Vertex.h:73

◆ vertex()

const reco::Vertex & BPHDecayVertex::vertex ( VertexFitter< 5 > *  fitter = nullptr,
const reco::BeamSpot bs = nullptr,
const GlobalPoint priorPos = nullptr,
const GlobalError priorError = nullptr 
) const
virtual

get reconstructed vertex

Definition at line 106 of file BPHDecayVertex.cc.

References cms::cuda::bs, fittedVertex, fitVertex(), oldVertex, savedBS, savedFitter, savedPE, and savedPP.

Referenced by Tau.Tau::dxy(), BPHWriteSpecificDecay::fill(), and validVertex().

109  {
110  if ((fitter == nullptr) && (bs == nullptr) && (priorPos == nullptr) && (priorError == nullptr)) {
111  fitter = savedFitter;
112  bs = savedBS;
113  priorPos = savedPP;
114  priorError = savedPE;
115  }
116  if (oldVertex || (fitter != savedFitter) || (bs != savedBS) || (priorPos != savedPP) || (priorError != savedPE)) {
117  if (fitter != nullptr) {
118  fitVertex(fitter, bs, priorPos, priorError);
119  } else {
120  KalmanVertexFitter kvf(true);
121  fitVertex(&kvf, bs, priorPos, priorError);
122  }
123  }
124  return fittedVertex;
125 }
VertexFitter< 5 > * savedFitter
virtual void fitVertex(VertexFitter< 5 > *fitter, const reco::BeamSpot *bs, const GlobalPoint *priorPos, const GlobalError *priorError) const
reco::Vertex fittedVertex
const GlobalError * savedPE
const GlobalPoint * savedPP
const reco::BeamSpot * savedBS

Member Data Documentation

◆ evSetup

const BPHEventSetupWrapper* BPHDecayVertex::evSetup
private

Definition at line 109 of file BPHDecayVertex.h.

Referenced by fTTracks(), getEventSetup(), and ~BPHDecayVertex().

◆ fittedVertex

reco::Vertex BPHDecayVertex::fittedVertex
mutableprivate

Definition at line 124 of file BPHDecayVertex.h.

Referenced by fitVertex(), validVertex(), and vertex().

◆ oldTracks

bool BPHDecayVertex::oldTracks
mutableprivate

Definition at line 115 of file BPHDecayVertex.h.

Referenced by fTracks(), fTTracks(), getTMode(), getTrack(), setNotUpdated(), tracks(), and validTracks().

◆ oldTTracks

bool BPHDecayVertex::oldTTracks
mutableprivate

Definition at line 116 of file BPHDecayVertex.h.

Referenced by fitVertex(), fTracks(), fTTracks(), getTransientTrack(), and transientTracks().

◆ oldVertex

bool BPHDecayVertex::oldVertex
mutableprivate

Definition at line 117 of file BPHDecayVertex.h.

Referenced by fitVertex(), setNotUpdated(), and vertex().

◆ rTracks

std::vector<const reco::Track*> BPHDecayVertex::rTracks
mutableprivate

Definition at line 119 of file BPHDecayVertex.h.

Referenced by fTracks(), and tracks().

◆ savedBS

const reco::BeamSpot* BPHDecayVertex::savedBS
mutableprivate

Definition at line 126 of file BPHDecayVertex.h.

Referenced by fitVertex(), and vertex().

◆ savedFitter

VertexFitter<5>* BPHDecayVertex::savedFitter
mutableprivate

Definition at line 125 of file BPHDecayVertex.h.

Referenced by fitVertex(), and vertex().

◆ savedPE

const GlobalError* BPHDecayVertex::savedPE
mutableprivate

Definition at line 128 of file BPHDecayVertex.h.

Referenced by fitVertex(), and vertex().

◆ savedPP

const GlobalPoint* BPHDecayVertex::savedPP
mutableprivate

Definition at line 127 of file BPHDecayVertex.h.

Referenced by fitVertex(), and vertex().

◆ searchMap

std::map<const reco::Candidate*, std::string> BPHDecayVertex::searchMap
private

Definition at line 112 of file BPHDecayVertex.h.

Referenced by BPHDecayVertex(), fTracks(), and getTrackSearchList().

◆ tkMap

std::map<const reco::Candidate*, const reco::Track*> BPHDecayVertex::tkMap
mutableprivate

Definition at line 121 of file BPHDecayVertex.h.

Referenced by fTracks(), fTTracks(), and getTrack().

◆ tmMap

std::map<const reco::Candidate*, char> BPHDecayVertex::tmMap
mutableprivate

Definition at line 122 of file BPHDecayVertex.h.

Referenced by fTracks(), and getTMode().

◆ trTracks

std::vector<reco::TransientTrack> BPHDecayVertex::trTracks
mutableprivate

Definition at line 120 of file BPHDecayVertex.h.

Referenced by fitVertex(), fTracks(), fTTracks(), and transientTracks().

◆ ttMap

std::map<const reco::Candidate*, reco::TransientTrack*> BPHDecayVertex::ttMap
mutableprivate

Definition at line 123 of file BPHDecayVertex.h.

Referenced by fTTracks(), and getTransientTrack().

◆ validTks

bool BPHDecayVertex::validTks
mutableprivate

Definition at line 118 of file BPHDecayVertex.h.

Referenced by fTracks(), setNotUpdated(), validTracks(), and validVertex().