CMS 3D CMS Logo

TransientVertex.cc
Go to the documentation of this file.
6 #include <algorithm>
7 
8 using namespace std;
9 using namespace reco;
10 
11 TransientVertex::TransientVertex() : theVertexState(), theOriginalTracks(),
12  theChi2(0), theNDF(0), vertexValid(false), withPrior(false),
13  theWeightMapIsAvailable(false), theCovMapAvailable(false),
14  withRefittedTracks(false)
15 {}
16 
17 
19  const std::vector<TransientTrack> & tracks, float chi2) :
20  theVertexState(pos, posError), theOriginalTracks(tracks),
24 {
25  theNDF = 2.*theOriginalTracks.size() - 3.;
26 }
27 
29  const GlobalError & posTimeError,
30  const std::vector<TransientTrack> & tracks, float chi2) :
31  theVertexState(pos, time, posTimeError), theOriginalTracks(tracks),
35 {
36  theNDF = 2.*theOriginalTracks.size() - 3.;
37 }
38 
39 
41  const std::vector<TransientTrack> & tracks, float chi2, float ndf) :
42  theVertexState(pos, posError), theOriginalTracks(tracks),
43  theChi2(chi2), theNDF(ndf), vertexValid(true), withPrior(false),
46 {
47 }
48 
50  const GlobalError & posTimeError,
51  const std::vector<TransientTrack> & tracks, float chi2, float ndf) :
52  theVertexState(pos, time, posTimeError), theOriginalTracks(tracks),
53  theChi2(chi2), theNDF(ndf), vertexValid(true), withPrior(false),
56 {
57 }
58 
59 
60 TransientVertex::TransientVertex(const GlobalPoint & priorPos, const GlobalError & priorErr,
61  const GlobalPoint & pos, const GlobalError & posError,
62  const std::vector<TransientTrack> & tracks, float chi2) :
63  thePriorVertexState(priorPos, priorErr), theVertexState(pos, posError),
64  theOriginalTracks(tracks), theChi2(chi2), theNDF(0), vertexValid(true),
67 {
68  theNDF = 2.*theOriginalTracks.size();
69 }
70 
71 TransientVertex::TransientVertex(const GlobalPoint & priorPos, const double priorTime, const GlobalError & priorErr,
72  const GlobalPoint & pos, const double time, const GlobalError & posError,
73  const std::vector<TransientTrack> & tracks, float chi2) :
74  thePriorVertexState(priorPos, priorTime, priorErr), theVertexState(pos, time, posError),
75  theOriginalTracks(tracks), theChi2(chi2), theNDF(0), vertexValid(true),
78 {
79  theNDF = 2.*theOriginalTracks.size();
80 }
81 
82 
83 TransientVertex::TransientVertex(const GlobalPoint & priorPos, const GlobalError & priorErr,
84  const GlobalPoint & pos, const GlobalError & posError,
85  const std::vector<TransientTrack> & tracks, float chi2, float ndf) :
86  thePriorVertexState(priorPos, priorErr), theVertexState(pos, posError),
87  theOriginalTracks(tracks), theChi2(chi2), theNDF(ndf), vertexValid(true),
90 {
91 }
92 
93 TransientVertex::TransientVertex(const GlobalPoint & priorPos, const double priorTime, const GlobalError & priorErr,
94  const GlobalPoint & pos, const double time, const GlobalError & posError,
95  const std::vector<TransientTrack> & tracks, float chi2, float ndf) :
96  thePriorVertexState(priorPos, priorTime, priorErr), theVertexState(pos, time, posError),
97  theOriginalTracks(tracks), theChi2(chi2), theNDF(ndf), vertexValid(true),
100 {
101 }
102 
103 
105  const std::vector<TransientTrack> & tracks, float chi2) :
106  theVertexState(state), theOriginalTracks(tracks),
110 {
111  theNDF = 2.*theOriginalTracks.size() - 3.;
112 }
113 
114 
116  const std::vector<TransientTrack> & tracks, float chi2, float ndf) :
117  theVertexState(state), theOriginalTracks(tracks),
118  theChi2(chi2), theNDF(ndf), vertexValid(true), withPrior(false),
121 {
122 }
123 
124 
126  const VertexState & state,
127  const std::vector<TransientTrack> & tracks,
128  float chi2) :
129  thePriorVertexState(prior), theVertexState(state),
130  theOriginalTracks(tracks), theChi2(chi2), theNDF(0), vertexValid(true),
133 {
134  theNDF = 2.*theOriginalTracks.size();
135 }
136 
137 
139  const VertexState & state,
140  const std::vector<TransientTrack> & tracks,
141  float chi2, float ndf) :
142  thePriorVertexState(prior), theVertexState(state),
143  theOriginalTracks(tracks), theChi2(chi2), theNDF(ndf), vertexValid(true),
146 {
147 }
148 
150 {
151  theWeightMap = theMap;
153 }
154 
156  const std::vector<reco::TransientTrack> & refittedTracks)
157 {
158  if (refittedTracks.empty())
159  throw VertexException("TransientVertex::refittedTracks: No refitted tracks stored in input container");
161  withRefittedTracks = true;
162 }
163 
164 
166 {
167  theCovMap = covMap;
168  theCovMapAvailable = true;
169 }
170 
173  std::vector<TransientTrack>::const_iterator foundTrack = find(theOriginalTracks.begin(),
174  theOriginalTracks.end(), track);
175  return ((foundTrack != theOriginalTracks.end()) ? 1. : 0.);
176  }
177  TransientTrackToFloatMap::const_iterator it = theWeightMap.find(track);
178  if (it != theWeightMap.end()) {
179  return(it->second);
180  }
181  return 0.;
182 
183 }
184 
187 {
188  if (!theCovMapAvailable) {
189  throw VertexException("TransientVertex::Track-to-track covariance matrices not available");
190  }
191  const TransientTrack* tr1;
192  const TransientTrack* tr2;
193  if (t1<t2) {
194  tr1 = &t1;
195  tr2 = &t2;
196  } else {
197  tr1 = &t2;
198  tr2 = &t1;
199  }
200  TTtoTTmap::const_iterator it = theCovMap.find(*tr1);
201  if (it != theCovMap.end()) {
202  const TTmap & tm = it->second;
203  TTmap::const_iterator nit = tm.find(*tr2);
204  if (nit != tm.end()) {
205  return( nit->second);
206  }
207  else {
208  throw VertexException("TransientVertex::requested Track-to-track covariance matrix does not exist");
209  }
210  }
211  else {
212  throw VertexException("TransientVertex::requested Track-to-track covariance matrix does not exist");
213  }
214 }
215 
216 
218 {
219  if (theRefittedTracks.empty())
220  throw VertexException("No refitted tracks stored in vertex");
221  std::vector<TransientTrack>::const_iterator it =
222  find(theRefittedTracks.begin(), theRefittedTracks.end(), refTrack);
223  if (it==theRefittedTracks.end())
224  throw VertexException("Refitted track not found in list");
225  size_t pos = it - theRefittedTracks.begin();
226  return theOriginalTracks[pos];
227 }
228 
230 {
231  if (theRefittedTracks.empty())
232  throw VertexException("No refitted tracks stored in vertex");
233  std::vector<TransientTrack>::const_iterator it =
235  if (it==theOriginalTracks.end())
236  throw VertexException("Track not found in list");
237  size_t pos = it - theOriginalTracks.begin();
238  return theRefittedTracks[pos];
239 }
240 
242 {
243  //If the vertex is invalid, return an invalid TV !
244  if (!isValid()) return Vertex();
245 
247  // RecoVertex::convertError(theVertexState.error()),
251  for (std::vector<TransientTrack>::const_iterator i = theOriginalTracks.begin();
252  i != theOriginalTracks.end(); ++i) {
253 // const TrackTransientTrack* ttt = dynamic_cast<const TrackTransientTrack*>((*i).basicTransientTrack());
254 // if ((ttt!=0) && (ttt->persistentTrackRef().isNonnull()))
255 // {
256 // TrackRef tr = ttt->persistentTrackRef();
257 // TrackBaseRef tbr(tr);
258  if (withRefittedTracks) {
259 
260  vertex.add((*i).trackBaseRef(), refittedTrack(*i).track(), trackWeight ( *i ) );
261  } else {
262  vertex.add((*i).trackBaseRef(), trackWeight ( *i ) );
263  }
264  //}
265  }
266  return vertex;
267 }
268 
270 {
271  using namespace reco;
272  if (!isValid()) return VertexCompositePtrCandidate();
273 
274  VertexCompositePtrCandidate vtxCompPtrCand;
275 
276  vtxCompPtrCand.setTime(vertexState().time());
277  vtxCompPtrCand.setCovariance(vertexState().error4D().matrix4D());
278  vtxCompPtrCand.setChi2AndNdof(totalChiSquared(), degreesOfFreedom());
279  vtxCompPtrCand.setVertex(Candidate::Point(position().x(),position().y(),position().z()));
280 
282  for(std::vector<reco::TransientTrack>::const_iterator tt = theOriginalTracks.begin(); tt != theOriginalTracks.end(); ++tt)
283  {
284  if (trackWeight(*tt) < 0.5)
285  continue;
286 
287  const CandidatePtrTransientTrack* cptt = dynamic_cast<const CandidatePtrTransientTrack*>(tt->basicTransientTrack());
288  if ( cptt==nullptr )
289  edm::LogError("DynamicCastingFailed") << "Casting of TransientTrack to CandidatePtrTransientTrack failed!";
290  else
291  {
292  p4 += cptt->candidate()->p4();
293  vtxCompPtrCand.addDaughter(cptt->candidate());
294  }
295  }
296 
297  //TODO: if has refitted tracks we should scale the candidate p4 to the refitted one
298  vtxCompPtrCand.setP4(p4);
299  return vtxCompPtrCand;
300 }
301 
double priorTime() const
Common base class.
float totalChiSquared() const
std::map< reco::TransientTrack, float > TransientTrackToFloatMap
std::map< reco::TransientTrack, TTmap > TTtoTTmap
Definition: TTtoTTmap.h:12
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:20
GlobalPoint position() const
Definition: VertexState.h:69
std::map< reco::TransientTrack, AlgebraicMatrix33 > TTmap
Definition: TTtoTTmap.h:11
const AlgebraicSymMatrix44 & matrix4D() const
void setVertex(const Point &vertex) override
set vertex
reco::TransientTrack refittedTrack(const reco::TransientTrack &track) const
VertexState thePriorVertexState
virtual const LorentzVector & p4() const =0
four-momentum Lorentz vector
float degreesOfFreedom() const
TransientTrackToFloatMap theWeightMap
double p4[4]
Definition: TauolaWrapper.h:92
GlobalPoint position() const
double time() const
Definition: VertexState.h:95
math::XYZPoint Point
point in the space
Definition: Vertex.h:39
std::vector< reco::TransientTrack > theOriginalTracks
float trackWeight(const reco::TransientTrack &track) const
CandidatePtr candidate() const override
void add(const TrackBaseRef &r, float w=1.0)
add a reference to a Track
GlobalError error4D() const
Definition: VertexState.h:80
const Track & track() const
reco::TransientTrack originalTrack(const reco::TransientTrack &refTrack) const
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:37
double time() const
fixed size matrix
VertexState theVertexState
void setCovariance(const CovarianceMatrix &m)
set covariance matrix
void addDaughter(const CandidatePtr &)
add a daughter via a reference
math::XYZPoint Point
point in the space
Definition: Candidate.h:41
std::vector< reco::TransientTrack > const & refittedTracks() const
std::vector< reco::TransientTrack > theRefittedTracks
void setChi2AndNdof(double chi2, double ndof)
set chi2 and ndof
VertexState const & vertexState() const
AlgebraicMatrix33 tkToTkCovariance(const reco::TransientTrack &t1, const reco::TransientTrack &t2) const
bool isValid() const
void setP4(const LorentzVector &p4) final
set 4-momentum
TransientTrackToFloatMap weightMap() const
ROOT::Math::SMatrix< double, 3, 3, ROOT::Math::MatRepStd< double, 3, 3 > > AlgebraicMatrix33