CMS 3D CMS Logo

V0Fitter.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: V0Producer
4 // Class: V0Fitter
5 //
13 //
14 // Original Author: Brian Drell
15 // Created: Fri May 18 22:57:40 CEST 2007
16 //
17 //
18 
19 #include "V0Fitter.h"
20 
28 #include <Math/Functions.h>
29 #include <Math/SMatrix.h>
30 #include <Math/SVector.h>
33 #include <memory>
34 #include <typeinfo>
35 
36 // pdg mass constants
37 namespace {
38  const double piMass = 0.13957018;
39  const double piMassSquared = piMass * piMass;
40  const double protonMass = 0.938272046;
41  const double protonMassSquared = protonMass * protonMass;
42  const double kShortMass = 0.497614;
43  const double lambdaMass = 1.115683;
44 } // namespace
45 
46 typedef ROOT::Math::SMatrix<double, 3, 3, ROOT::Math::MatRepSym<double, 3>> SMatrixSym3D;
47 typedef ROOT::Math::SVector<double, 3> SVector3;
48 
49 V0Fitter::V0Fitter(const edm::ParameterSet& theParameters, edm::ConsumesCollector&& iC) : esTokenMF_(iC.esConsumes()) {
50  token_beamSpot = iC.consumes<reco::BeamSpot>(theParameters.getParameter<edm::InputTag>("beamSpot"));
51  useVertex_ = theParameters.getParameter<bool>("useVertex");
52  if (useVertex_)
53  token_vertices = iC.consumes<std::vector<reco::Vertex>>(theParameters.getParameter<edm::InputTag>("vertices"));
54 
55  token_tracks = iC.consumes<reco::TrackCollection>(theParameters.getParameter<edm::InputTag>("trackRecoAlgorithm"));
56  doFit_ = theParameters.getParameter<bool>("doFit");
57  vertexFitter_ = theParameters.getParameter<bool>("vertexFitter");
58  useRefTracks_ = theParameters.getParameter<bool>("useRefTracks");
59 
60  // whether to reconstruct KShorts
61  doKShorts_ = theParameters.getParameter<bool>("doKShorts");
62  // whether to reconstruct Lambdas
63  doLambdas_ = theParameters.getParameter<bool>("doLambdas");
64 
65  // cuts on initial track selection
66  tkChi2Cut_ = theParameters.getParameter<double>("tkChi2Cut");
67  tkNHitsCut_ = theParameters.getParameter<int>("tkNHitsCut");
68  tkPtCut_ = theParameters.getParameter<double>("tkPtCut");
69  tkIPSigXYCut_ = theParameters.getParameter<double>("tkIPSigXYCut");
70  tkIPSigZCut_ = theParameters.getParameter<double>("tkIPSigZCut");
71 
72  // cuts on vertex
73  vtxChi2Cut_ = theParameters.getParameter<double>("vtxChi2Cut");
74  vtxDecaySigXYZCut_ = theParameters.getParameter<double>("vtxDecaySigXYZCut");
75  vtxDecaySigXYCut_ = theParameters.getParameter<double>("vtxDecaySigXYCut");
76  vtxDecayXYCut_ = theParameters.getParameter<double>("vtxDecayXYCut");
77  ssVtxDecayXYCut_ = theParameters.getParameter<double>("ssVtxDecayXYCut");
78  // miscellaneous cuts
79  allowSS_ = theParameters.getParameter<bool>("allowSS");
80  innerOuterTkDCAThreshold_ = theParameters.getParameter<double>("innerOuterTkDCAThreshold");
81  innerTkDCACut_ = theParameters.getParameter<double>("innerTkDCACut");
82  outerTkDCACut_ = theParameters.getParameter<double>("outerTkDCACut");
83  allowWideAngleVtx_ = theParameters.getParameter<bool>("allowWideAngleVtx");
84  mPiPiCut_ = theParameters.getParameter<double>("mPiPiCut");
85  innerHitPosCut_ = theParameters.getParameter<double>("innerHitPosCut");
86  cosThetaXYCut_ = theParameters.getParameter<double>("cosThetaXYCut");
87  cosThetaXYZCut_ = theParameters.getParameter<double>("cosThetaXYZCut");
88  // cuts on the V0 candidate mass
89  kShortMassCut_ = theParameters.getParameter<double>("kShortMassCut");
90  lambdaMassCut_ = theParameters.getParameter<double>("lambdaMassCut");
91 }
92 
93 // method containing the algorithm for vertex reconstruction
95  const edm::EventSetup& iSetup,
98  using std::vector;
99 
100  edm::Handle<reco::TrackCollection> theTrackHandle;
101  iEvent.getByToken(token_tracks, theTrackHandle);
102  if (theTrackHandle->empty())
103  return;
104  const reco::TrackCollection* theTrackCollection = theTrackHandle.product();
105 
106  edm::Handle<reco::BeamSpot> theBeamSpotHandle;
107  iEvent.getByToken(token_beamSpot, theBeamSpotHandle);
108  const reco::BeamSpot* theBeamSpot = theBeamSpotHandle.product();
109  math::XYZPoint referencePos(theBeamSpot->position());
110 
111  reco::Vertex referenceVtx;
112  if (useVertex_) {
114  iEvent.getByToken(token_vertices, vertices);
115  referenceVtx = vertices->at(0);
116  referencePos = referenceVtx.position();
117  }
118 
119  const MagneticField* theMagneticField = &iSetup.getData(esTokenMF_);
120 
121  std::vector<reco::TrackRef> theTrackRefs;
122  std::vector<reco::TransientTrack> theTransTracks;
123 
124  // fill vectors of TransientTracks and TrackRefs after applying preselection cuts
125  for (reco::TrackCollection::const_iterator iTk = theTrackCollection->begin(); iTk != theTrackCollection->end();
126  ++iTk) {
127  const reco::Track* tmpTrack = &(*iTk);
128  double ipsigXY = std::abs(tmpTrack->dxy(*theBeamSpot) / tmpTrack->dxyError());
129  if (useVertex_)
130  ipsigXY = std::abs(tmpTrack->dxy(referencePos) / tmpTrack->dxyError());
131  double ipsigZ = std::abs(tmpTrack->dz(referencePos) / tmpTrack->dzError());
132  if (tmpTrack->normalizedChi2() < tkChi2Cut_ && tmpTrack->numberOfValidHits() >= tkNHitsCut_ &&
133  tmpTrack->pt() > tkPtCut_ && ipsigXY > tkIPSigXYCut_ && ipsigZ > tkIPSigZCut_) {
134  reco::TrackRef tmpRef(theTrackHandle, std::distance(theTrackCollection->begin(), iTk));
135  theTrackRefs.push_back(std::move(tmpRef));
136  reco::TransientTrack tmpTransient(*tmpRef, theMagneticField);
137  theTransTracks.push_back(std::move(tmpTransient));
138  }
139  }
140  // good tracks have now been selected for vertexing
141 
142  // loop over tracks and vertex good charged track pairs
143  for (unsigned int trdx1 = 0; trdx1 < theTrackRefs.size(); ++trdx1) {
144  for (unsigned int trdx2 = trdx1 + 1; trdx2 < theTrackRefs.size(); ++trdx2) {
145  reco::TrackRef positiveTrackRef;
146  reco::TrackRef negativeTrackRef;
147  reco::TransientTrack* posTransTkPtr = nullptr;
148  reco::TransientTrack* negTransTkPtr = nullptr;
149 
150  if (theTrackRefs[trdx1]->charge() < 0. && theTrackRefs[trdx2]->charge() > 0.) {
151  negativeTrackRef = theTrackRefs[trdx1];
152  positiveTrackRef = theTrackRefs[trdx2];
153  negTransTkPtr = &theTransTracks[trdx1];
154  posTransTkPtr = &theTransTracks[trdx2];
155  } else if (theTrackRefs[trdx1]->charge() > 0. && theTrackRefs[trdx2]->charge() < 0.) {
156  negativeTrackRef = theTrackRefs[trdx2];
157  positiveTrackRef = theTrackRefs[trdx1];
158  negTransTkPtr = &theTransTracks[trdx2];
159  posTransTkPtr = &theTransTracks[trdx1];
160  } else {
161  if (!allowSS_)
162  continue;
163 
164  // if same-sign pairs are allowed, assign the negative and positive tracks arbitrarily
165  negativeTrackRef = theTrackRefs[trdx1];
166  positiveTrackRef = theTrackRefs[trdx2];
167  negTransTkPtr = &theTransTracks[trdx1];
168  posTransTkPtr = &theTransTracks[trdx2];
169  }
170 
171  // measure distance between tracks at their closest approach
172 
173  //these two variables are needed to 'pin' the temporary value returned to the stack
174  // in order to keep posState and negState from pointing to destructed objects
175  auto const& posImpact = posTransTkPtr->impactPointTSCP();
176  auto const& negImpact = negTransTkPtr->impactPointTSCP();
177  if (!posImpact.isValid() || !negImpact.isValid())
178  continue;
179  FreeTrajectoryState const& posState = posImpact.theState();
180  FreeTrajectoryState const& negState = negImpact.theState();
182  cApp.calculate(posState, negState);
183  if (!cApp.status())
184  continue;
185  float dca = std::abs(cApp.distance());
186 
187  // the POCA should at least be in the sensitive volume
188  GlobalPoint cxPt = cApp.crossingPoint();
189  const double cxPtR2 = cxPt.x() * cxPt.x() + cxPt.y() * cxPt.y();
190  if (cxPtR2 > 120. * 120. || std::abs(cxPt.z()) > 300.)
191  continue;
192 
193  // allow for different DCA cuts depending on position of POCA
195  if (dca > innerTkDCACut_)
196  continue;
197  } else {
198  if (dca > outerTkDCACut_)
199  continue;
200  }
201 
202  // the tracks should at least point in the same quadrant
203  TrajectoryStateClosestToPoint const& posTSCP = posTransTkPtr->trajectoryStateClosestToPoint(cxPt);
204  TrajectoryStateClosestToPoint const& negTSCP = negTransTkPtr->trajectoryStateClosestToPoint(cxPt);
205  if (!posTSCP.isValid() || !negTSCP.isValid())
206  continue;
207  if (!allowWideAngleVtx_ && posTSCP.momentum().dot(negTSCP.momentum()) < 0)
208  continue;
209 
210  // calculate mPiPi
211  double totalE = sqrt(posTSCP.momentum().mag2() + piMassSquared) + sqrt(negTSCP.momentum().mag2() + piMassSquared);
212  double totalESq = totalE * totalE;
213  double totalPSq = (posTSCP.momentum() + negTSCP.momentum()).mag2();
214  double massSquared = totalESq - totalPSq;
215  if (massSquared > mPiPiCut_ * mPiPiCut_)
216  continue;
217 
218  // Fill the vector of TransientTracks to send to KVF
219  std::vector<reco::TransientTrack> transTracks;
220  transTracks.reserve(2);
221  transTracks.push_back(*posTransTkPtr);
222  transTracks.push_back(*negTransTkPtr);
223 
224  // create the vertex fitter object and vertex the tracks
225  const GlobalError dummyError(1.0e-3, 0.0, 1.0e-3, 0.0, 0.0, 1.0e-3);
226  TransientVertex theRecoVertex(cxPt, dummyError, transTracks, 1.0e-3);
227  if (doFit_) {
228  if (vertexFitter_) {
229  KalmanVertexFitter theKalmanFitter(useRefTracks_ == 0 ? false : true);
230  theRecoVertex = theKalmanFitter.vertex(transTracks);
231  } else {
232  useRefTracks_ = false;
233  AdaptiveVertexFitter theAdaptiveFitter;
234  theRecoVertex = theAdaptiveFitter.vertex(transTracks);
235  }
236  }
237  if (!theRecoVertex.isValid())
238  continue;
239 
240  reco::Vertex theVtx = theRecoVertex;
241  if (theVtx.normalizedChi2() > vtxChi2Cut_)
242  continue;
243  GlobalPoint vtxPos(theVtx.x(), theVtx.y(), theVtx.z());
244 
245  // 2D decay significance
246  SMatrixSym3D totalCov = theBeamSpot->rotatedCovariance3D() + theVtx.covariance();
247  if (useVertex_)
248  totalCov = referenceVtx.covariance() + theVtx.covariance();
249  SVector3 distVecXY(vtxPos.x() - referencePos.x(), vtxPos.y() - referencePos.y(), 0.);
250  double distMagXY = ROOT::Math::Mag(distVecXY);
251  if (distMagXY < vtxDecayXYCut_)
252  continue;
253  if (posTransTkPtr->charge() * negTransTkPtr->charge() > 0 && distMagXY < ssVtxDecayXYCut_)
254  continue;
255  double sigmaDistMagXY = sqrt(ROOT::Math::Similarity(totalCov, distVecXY)) / distMagXY;
256  if (distMagXY < vtxDecaySigXYCut_ * sigmaDistMagXY)
257  continue;
258 
259  // 3D decay significance
260  if (vtxDecaySigXYZCut_ > 0.) {
261  SVector3 distVecXYZ(
262  vtxPos.x() - referencePos.x(), vtxPos.y() - referencePos.y(), vtxPos.z() - referencePos.z());
263  double distMagXYZ = ROOT::Math::Mag(distVecXYZ);
264  double sigmaDistMagXYZ = sqrt(ROOT::Math::Similarity(totalCov, distVecXYZ)) / distMagXYZ;
265  if (distMagXYZ / sigmaDistMagXYZ < vtxDecaySigXYZCut_)
266  continue;
267  }
268 
269  // make sure the vertex radius is within the inner track hit radius
270  double tkHitPosLimitSquared =
271  (distMagXY - sigmaDistMagXY * innerHitPosCut_) * (distMagXY - sigmaDistMagXY * innerHitPosCut_);
272  if (innerHitPosCut_ > 0. && positiveTrackRef->innerOk()) {
273  reco::Vertex::Point posTkHitPos = positiveTrackRef->innerPosition();
274  double posTkHitPosD2 = (posTkHitPos.x() - referencePos.x()) * (posTkHitPos.x() - referencePos.x()) +
275  (posTkHitPos.y() - referencePos.y()) * (posTkHitPos.y() - referencePos.y());
276  if (posTkHitPosD2 < tkHitPosLimitSquared)
277  continue;
278  }
279  if (innerHitPosCut_ > 0. && negativeTrackRef->innerOk()) {
280  reco::Vertex::Point negTkHitPos = negativeTrackRef->innerPosition();
281  double negTkHitPosD2 = (negTkHitPos.x() - referencePos.x()) * (negTkHitPos.x() - referencePos.x()) +
282  (negTkHitPos.y() - referencePos.y()) * (negTkHitPos.y() - referencePos.y());
283  if (negTkHitPosD2 < tkHitPosLimitSquared)
284  continue;
285  }
286 
287  std::unique_ptr<TrajectoryStateClosestToPoint> trajPlus;
288  std::unique_ptr<TrajectoryStateClosestToPoint> trajMins;
289  std::vector<reco::TransientTrack> theRefTracks;
290  if (theRecoVertex.hasRefittedTracks()) {
291  theRefTracks = theRecoVertex.refittedTracks();
292  }
293 
294  if (useRefTracks_ && theRefTracks.size() > 1) {
295  reco::TransientTrack* thePositiveRefTrack = nullptr;
296  reco::TransientTrack* theNegativeRefTrack = nullptr;
297  for (std::vector<reco::TransientTrack>::iterator iTrack = theRefTracks.begin(); iTrack != theRefTracks.end();
298  ++iTrack) {
299  if (iTrack->track().charge() > 0.) {
300  thePositiveRefTrack = &*iTrack;
301  } else if (iTrack->track().charge() < 0.) {
302  theNegativeRefTrack = &*iTrack;
303  }
304  }
305  if (thePositiveRefTrack == nullptr || theNegativeRefTrack == nullptr)
306  continue;
307  trajPlus =
308  std::make_unique<TrajectoryStateClosestToPoint>(thePositiveRefTrack->trajectoryStateClosestToPoint(vtxPos));
309  trajMins =
310  std::make_unique<TrajectoryStateClosestToPoint>(theNegativeRefTrack->trajectoryStateClosestToPoint(vtxPos));
311  } else {
312  trajPlus =
313  std::make_unique<TrajectoryStateClosestToPoint>(posTransTkPtr->trajectoryStateClosestToPoint(vtxPos));
314  trajMins =
315  std::make_unique<TrajectoryStateClosestToPoint>(negTransTkPtr->trajectoryStateClosestToPoint(vtxPos));
316  }
317 
318  if (trajPlus.get() == nullptr || trajMins.get() == nullptr || !trajPlus->isValid() || !trajMins->isValid())
319  continue;
320 
321  GlobalVector positiveP(trajPlus->momentum());
322  GlobalVector negativeP(trajMins->momentum());
323  GlobalVector totalP(positiveP + negativeP);
324 
325  // 2D pointing angle
326  double dx = theVtx.x() - referencePos.x();
327  double dy = theVtx.y() - referencePos.y();
328  double px = totalP.x();
329  double py = totalP.y();
330  double angleXY = (dx * px + dy * py) / (sqrt(dx * dx + dy * dy) * sqrt(px * px + py * py));
331  if (angleXY < cosThetaXYCut_)
332  continue;
333 
334  // 3D pointing angle
335  if (cosThetaXYZCut_ > -1.) {
336  double dz = theVtx.z() - referencePos.z();
337  double pz = totalP.z();
338  double angleXYZ =
339  (dx * px + dy * py + dz * pz) / (sqrt(dx * dx + dy * dy + dz * dz) * sqrt(px * px + py * py + pz * pz));
340  if (angleXYZ < cosThetaXYZCut_)
341  continue;
342  }
343 
344  // calculate total energy of V0 3 ways: assume it's a kShort, a Lambda, or a LambdaBar.
345  double piPlusE = sqrt(positiveP.mag2() + piMassSquared);
346  double piMinusE = sqrt(negativeP.mag2() + piMassSquared);
347  double protonE = sqrt(positiveP.mag2() + protonMassSquared);
348  double antiProtonE = sqrt(negativeP.mag2() + protonMassSquared);
349  double kShortETot = piPlusE + piMinusE;
350  double lambdaEtot = protonE + piMinusE;
351  double lambdaBarEtot = antiProtonE + piPlusE;
352 
353  // Create momentum 4-vectors for the 3 candidate types
354  const reco::Particle::LorentzVector kShortP4(totalP.x(), totalP.y(), totalP.z(), kShortETot);
355  const reco::Particle::LorentzVector lambdaP4(totalP.x(), totalP.y(), totalP.z(), lambdaEtot);
356  const reco::Particle::LorentzVector lambdaBarP4(totalP.x(), totalP.y(), totalP.z(), lambdaBarEtot);
357 
358  reco::Particle::Point vtx(theVtx.x(), theVtx.y(), theVtx.z());
359  const reco::Vertex::CovarianceMatrix vtxCov(theVtx.covariance());
360  double vtxChi2(theVtx.chi2());
361  double vtxNdof(theVtx.ndof());
362 
363  // Create the VertexCompositeCandidate object that will be stored in the Event
364  reco::VertexCompositeCandidate* theKshort = nullptr;
365  reco::VertexCompositeCandidate* theLambda = nullptr;
366  reco::VertexCompositeCandidate* theLambdaBar = nullptr;
367 
368  if (doKShorts_) {
369  theKshort = new reco::VertexCompositeCandidate(0, kShortP4, vtx, vtxCov, vtxChi2, vtxNdof);
370  }
371  if (doLambdas_) {
372  if (positiveP.mag2() > negativeP.mag2()) {
373  theLambda = new reco::VertexCompositeCandidate(0, lambdaP4, vtx, vtxCov, vtxChi2, vtxNdof);
374  } else {
375  theLambdaBar = new reco::VertexCompositeCandidate(0, lambdaBarP4, vtx, vtxCov, vtxChi2, vtxNdof);
376  }
377  }
378 
379  // Create daughter candidates for the VertexCompositeCandidates
380  reco::RecoChargedCandidate thePiPlusCand(
381  1, reco::Particle::LorentzVector(positiveP.x(), positiveP.y(), positiveP.z(), piPlusE), vtx);
382  thePiPlusCand.setTrack(positiveTrackRef);
383 
384  reco::RecoChargedCandidate thePiMinusCand(
385  -1, reco::Particle::LorentzVector(negativeP.x(), negativeP.y(), negativeP.z(), piMinusE), vtx);
386  thePiMinusCand.setTrack(negativeTrackRef);
387 
388  reco::RecoChargedCandidate theProtonCand(
389  1, reco::Particle::LorentzVector(positiveP.x(), positiveP.y(), positiveP.z(), protonE), vtx);
390  theProtonCand.setTrack(positiveTrackRef);
391 
392  reco::RecoChargedCandidate theAntiProtonCand(
393  -1, reco::Particle::LorentzVector(negativeP.x(), negativeP.y(), negativeP.z(), antiProtonE), vtx);
394  theAntiProtonCand.setTrack(negativeTrackRef);
395 
396  AddFourMomenta addp4;
397  // Store the daughter Candidates in the VertexCompositeCandidates if they pass mass cuts
398  if (doKShorts_) {
399  theKshort->addDaughter(thePiPlusCand);
400  theKshort->addDaughter(thePiMinusCand);
401  theKshort->setPdgId(310);
402  addp4.set(*theKshort);
403  if (theKshort->mass() < kShortMass + kShortMassCut_ && theKshort->mass() > kShortMass - kShortMassCut_) {
404  theKshorts.push_back(std::move(*theKshort));
405  }
406  }
407  if (doLambdas_ && theLambda) {
408  theLambda->addDaughter(theProtonCand);
409  theLambda->addDaughter(thePiMinusCand);
410  theLambda->setPdgId(3122);
411  addp4.set(*theLambda);
412  if (theLambda->mass() < lambdaMass + lambdaMassCut_ && theLambda->mass() > lambdaMass - lambdaMassCut_) {
413  theLambdas.push_back(std::move(*theLambda));
414  }
415  } else if (doLambdas_ && theLambdaBar) {
416  theLambdaBar->addDaughter(theAntiProtonCand);
417  theLambdaBar->addDaughter(thePiPlusCand);
418  theLambdaBar->setPdgId(-3122);
419  addp4.set(*theLambdaBar);
420  if (theLambdaBar->mass() < lambdaMass + lambdaMassCut_ && theLambdaBar->mass() > lambdaMass - lambdaMassCut_) {
421  theLambdas.push_back(std::move(*theLambdaBar));
422  }
423  }
424 
425  delete theKshort;
426  delete theLambda;
427  delete theLambdaBar;
428  theKshort = theLambda = theLambdaBar = nullptr;
429  }
430  }
431 }
bool doLambdas_
Definition: V0Fitter.h:55
double tkChi2Cut_
Definition: V0Fitter.h:58
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
std::vector< VertexCompositeCandidate > VertexCompositeCandidateCollection
collection of Candidate objects
ROOT::Math::SMatrix< double, 3, 3, ROOT::Math::MatRepSym< double, 3 > > SMatrixSym3D
Definition: V0Fitter.cc:46
double cosThetaXYZCut_
Definition: V0Fitter.h:78
double outerTkDCACut_
Definition: V0Fitter.h:73
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
double lambdaMassCut_
Definition: V0Fitter.h:81
void set(reco::Candidate &c) const
set up a candidate
double z() const
z coordinate
Definition: Vertex.h:133
float distance() const override
double tkIPSigXYCut_
Definition: V0Fitter.h:61
unsigned short numberOfValidHits() const
number of valid hits found
Definition: TrackBase.h:798
const Point & position() const
position
Definition: BeamSpot.h:59
bool vertexFitter_
Definition: V0Fitter.h:52
edm::EDGetTokenT< std::vector< reco::Vertex > > token_vertices
Definition: V0Fitter.h:86
bool allowWideAngleVtx_
Definition: V0Fitter.h:74
PreciseFloatType< T, U >::Type dot(const Vector3DBase< U, FrameTag > &v) const
Definition: Vector3DBase.h:99
T const * product() const
Definition: Handle.h:70
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
bool allowSS_
Definition: V0Fitter.h:70
bool doFit_
Definition: V0Fitter.h:51
T mag2() const
Definition: PV3DBase.h:63
CachingVertex< 5 > vertex(const std::vector< reco::TransientTrack > &) const override
double ndof() const
Definition: Vertex.h:123
double vtxDecaySigXYZCut_
Definition: V0Fitter.h:66
math::Error< dimension >::type CovarianceMatrix
covariance error matrix (3x3)
Definition: Vertex.h:46
CachingVertex< 5 > vertex(const std::vector< reco::TransientTrack > &tracks) const override
double pt() const
track transverse momentum
Definition: TrackBase.h:637
TrajectoryStateClosestToPoint trajectoryStateClosestToPoint(const GlobalPoint &point) const
T x() const
Definition: PV3DBase.h:59
double covariance(int i, int j) const
(i, j)-th element of error matrix, i, j = 0, ... 2
Definition: Vertex.h:148
int tkNHitsCut_
Definition: V0Fitter.h:59
double dz() const
dz parameter (= dsz/cos(lambda)). This is the track z0 w.r.t (0,0,0) only if the refPoint is close to...
Definition: TrackBase.h:622
GlobalPoint crossingPoint() const override
int iEvent
Definition: GenABIO.cc:224
bool useRefTracks_
Definition: V0Fitter.h:53
double dxyError() const
error on dxy
Definition: TrackBase.h:769
edm::EDGetTokenT< reco::TrackCollection > token_tracks
Definition: V0Fitter.h:83
bool isValid() const
double dzError() const
error on dz
Definition: TrackBase.h:778
TrajectoryStateClosestToPoint impactPointTSCP() const
bool calculate(const TrajectoryStateOnSurface &sta, const TrajectoryStateOnSurface &stb) override
T sqrt(T t)
Definition: SSEVec.h:19
math::XYZPoint Point
point in the space
Definition: Particle.h:25
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
void fitAll(const edm::Event &iEvent, const edm::EventSetup &iSetup, reco::VertexCompositeCandidateCollection &k, reco::VertexCompositeCandidateCollection &l)
Definition: V0Fitter.cc:94
bool useVertex_
Definition: V0Fitter.h:85
math::XYZPoint Point
point in the space
Definition: Vertex.h:40
Covariance3DMatrix rotatedCovariance3D() const
Definition: BeamSpot.cc:73
double mPiPiCut_
Definition: V0Fitter.h:75
double x() const
x coordinate
Definition: Vertex.h:129
double vtxDecayXYCut_
Definition: V0Fitter.h:67
double kShortMassCut_
Definition: V0Fitter.h:80
V0Fitter(const edm::ParameterSet &theParams, edm::ConsumesCollector &&iC)
Definition: V0Fitter.cc:49
double tkIPSigZCut_
Definition: V0Fitter.h:62
double y() const
y coordinate
Definition: Vertex.h:131
TrackCharge charge() const
void addDaughter(const Candidate &, const std::string &s="")
add a clone of the passed candidate as daughter
double normalizedChi2() const
chi-squared divided by n.d.o.f.
Definition: Vertex.h:125
double cosThetaXYCut_
Definition: V0Fitter.h:77
double innerHitPosCut_
Definition: V0Fitter.h:76
double vtxDecaySigXYCut_
Definition: V0Fitter.h:65
bool hasRefittedTracks() const
T mag2() const
The vector magnitude squared. Equivalent to vec.dot(vec)
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
double ssVtxDecayXYCut_
Definition: V0Fitter.h:68
double normalizedChi2() const
chi-squared divided by n.d.o.f. (or chi-squared * 1e6 if n.d.o.f. is zero)
Definition: TrackBase.h:593
double innerOuterTkDCAThreshold_
Definition: V0Fitter.h:71
double chi2() const
chi-squares
Definition: Vertex.h:116
double tkPtCut_
Definition: V0Fitter.h:60
double mass() const final
mass
double vtxChi2Cut_
Definition: V0Fitter.h:64
edm::EDGetTokenT< reco::BeamSpot > token_beamSpot
Definition: V0Fitter.h:84
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > esTokenMF_
Definition: V0Fitter.h:49
ROOT::Math::SVector< double, 3 > SVector3
Definition: V0Fitter.cc:47
void setTrack(const reco::TrackRef &r)
set reference to track
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Particle.h:21
void setPdgId(int pdgId) final
bool status() const override
double innerTkDCACut_
Definition: V0Fitter.h:72
def move(src, dest)
Definition: eostools.py:511
std::vector< reco::TransientTrack > const & refittedTracks() const
double dxy() const
dxy parameter. (This is the transverse impact parameter w.r.t. to (0,0,0) ONLY if refPoint is close t...
Definition: TrackBase.h:608
bool doKShorts_
Definition: V0Fitter.h:54