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/SVector.h>
30 #include <Math/SMatrix.h>
31 #include <typeinfo>
32 #include <memory>
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 }
45 
46 typedef ROOT::Math::SMatrix<double, 3, 3, ROOT::Math::MatRepSym<double, 3> > SMatrixSym3D;
47 typedef ROOT::Math::SVector<double, 3> SVector3;
48 
50 {
51  token_beamSpot = iC.consumes<reco::BeamSpot>(theParameters.getParameter<edm::InputTag>("beamSpot"));
52  useVertex_ = theParameters.getParameter<bool>("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  vertexFitter_ = theParameters.getParameter<bool>("vertexFitter");
57  useRefTracks_ = theParameters.getParameter<bool>("useRefTracks");
58 
59  // whether to reconstruct KShorts
60  doKShorts_ = theParameters.getParameter<bool>("doKShorts");
61  // whether to reconstruct Lambdas
62  doLambdas_ = theParameters.getParameter<bool>("doLambdas");
63 
64  // cuts on initial track selection
65  tkChi2Cut_ = theParameters.getParameter<double>("tkChi2Cut");
66  tkNHitsCut_ = theParameters.getParameter<int>("tkNHitsCut");
67  tkPtCut_ = theParameters.getParameter<double>("tkPtCut");
68  tkIPSigXYCut_ = theParameters.getParameter<double>("tkIPSigXYCut");
69  tkIPSigZCut_ = theParameters.getParameter<double>("tkIPSigZCut");
70 
71  // cuts on vertex
72  vtxChi2Cut_ = theParameters.getParameter<double>("vtxChi2Cut");
73  vtxDecaySigXYZCut_ = theParameters.getParameter<double>("vtxDecaySigXYZCut");
74  vtxDecaySigXYCut_ = theParameters.getParameter<double>("vtxDecaySigXYCut");
75  // miscellaneous cuts
76  tkDCACut_ = theParameters.getParameter<double>("tkDCACut");
77  mPiPiCut_ = theParameters.getParameter<double>("mPiPiCut");
78  innerHitPosCut_ = theParameters.getParameter<double>("innerHitPosCut");
79  cosThetaXYCut_ = theParameters.getParameter<double>("cosThetaXYCut");
80  cosThetaXYZCut_ = theParameters.getParameter<double>("cosThetaXYZCut");
81  // cuts on the V0 candidate mass
82  kShortMassCut_ = theParameters.getParameter<double>("kShortMassCut");
83  lambdaMassCut_ = theParameters.getParameter<double>("lambdaMassCut");
84 }
85 
86 // method containing the algorithm for vertex reconstruction
89 {
90  using std::vector;
91 
93  iEvent.getByToken(token_tracks, theTrackHandle);
94  if (theTrackHandle->empty()) return;
95  const reco::TrackCollection* theTrackCollection = theTrackHandle.product();
96 
97  edm::Handle<reco::BeamSpot> theBeamSpotHandle;
98  iEvent.getByToken(token_beamSpot, theBeamSpotHandle);
99  const reco::BeamSpot* theBeamSpot = theBeamSpotHandle.product();
100  math::XYZPoint referencePos(theBeamSpot->position());
101 
102  reco::Vertex referenceVtx;
103  if (useVertex_) {
105  iEvent.getByToken(token_vertices, vertices);
106  referenceVtx = vertices->at(0);
107  referencePos = referenceVtx.position();
108  }
109 
110  edm::ESHandle<MagneticField> theMagneticFieldHandle;
111  iSetup.get<IdealMagneticFieldRecord>().get(theMagneticFieldHandle);
112  const MagneticField* theMagneticField = theMagneticFieldHandle.product();
113 
114  std::vector<reco::TrackRef> theTrackRefs;
115  std::vector<reco::TransientTrack> theTransTracks;
116 
117  // fill vectors of TransientTracks and TrackRefs after applying preselection cuts
118  for (reco::TrackCollection::const_iterator iTk = theTrackCollection->begin(); iTk != theTrackCollection->end(); ++iTk) {
119  const reco::Track* tmpTrack = &(*iTk);
120  double ipsigXY = std::abs(tmpTrack->dxy(*theBeamSpot)/tmpTrack->dxyError());
121  if (useVertex_) ipsigXY = std::abs(tmpTrack->dxy(referencePos)/tmpTrack->dxyError());
122  double ipsigZ = std::abs(tmpTrack->dz(referencePos)/tmpTrack->dzError());
123  if (tmpTrack->normalizedChi2() < tkChi2Cut_ && tmpTrack->numberOfValidHits() >= tkNHitsCut_ &&
124  tmpTrack->pt() > tkPtCut_ && ipsigXY > tkIPSigXYCut_ && ipsigZ > tkIPSigZCut_) {
125  reco::TrackRef tmpRef(theTrackHandle, std::distance(theTrackCollection->begin(), iTk));
126  theTrackRefs.push_back(std::move(tmpRef));
127  reco::TransientTrack tmpTransient(*tmpRef, theMagneticField);
128  theTransTracks.push_back(std::move(tmpTransient));
129  }
130  }
131  // good tracks have now been selected for vertexing
132 
133  // loop over tracks and vertex good charged track pairs
134  for (unsigned int trdx1 = 0; trdx1 < theTrackRefs.size(); ++trdx1) {
135  for (unsigned int trdx2 = trdx1 + 1; trdx2 < theTrackRefs.size(); ++trdx2) {
136 
137  reco::TrackRef positiveTrackRef;
138  reco::TrackRef negativeTrackRef;
139  reco::TransientTrack* posTransTkPtr = nullptr;
140  reco::TransientTrack* negTransTkPtr = nullptr;
141 
142  if (theTrackRefs[trdx1]->charge() < 0. && theTrackRefs[trdx2]->charge() > 0.) {
143  negativeTrackRef = theTrackRefs[trdx1];
144  positiveTrackRef = theTrackRefs[trdx2];
145  negTransTkPtr = &theTransTracks[trdx1];
146  posTransTkPtr = &theTransTracks[trdx2];
147  } else if (theTrackRefs[trdx1]->charge() > 0. && theTrackRefs[trdx2]->charge() < 0.) {
148  negativeTrackRef = theTrackRefs[trdx2];
149  positiveTrackRef = theTrackRefs[trdx1];
150  negTransTkPtr = &theTransTracks[trdx2];
151  posTransTkPtr = &theTransTracks[trdx1];
152  } else {
153  continue;
154  }
155 
156  // measure distance between tracks at their closest approach
157 
158  //these two variables are needed to 'pin' the temporary value returned to the stack
159  // in order to keep posState and negState from pointing to destructed objects
160  auto const& posImpact = posTransTkPtr->impactPointTSCP();
161  auto const& negImpact = negTransTkPtr->impactPointTSCP();
162  if (!posImpact.isValid() || !negImpact.isValid()) continue;
163  FreeTrajectoryState const & posState = posImpact.theState();
164  FreeTrajectoryState const & negState = negImpact.theState();
166  cApp.calculate(posState, negState);
167  if (!cApp.status()) continue;
168  float dca = std::abs(cApp.distance());
169  if (dca > tkDCACut_) continue;
170 
171  // the POCA should at least be in the sensitive volume
172  GlobalPoint cxPt = cApp.crossingPoint();
173  if (sqrt(cxPt.x()*cxPt.x() + cxPt.y()*cxPt.y()) > 120. || std::abs(cxPt.z()) > 300.) continue;
174 
175  // the tracks should at least point in the same quadrant
176  TrajectoryStateClosestToPoint const & posTSCP = posTransTkPtr->trajectoryStateClosestToPoint(cxPt);
177  TrajectoryStateClosestToPoint const & negTSCP = negTransTkPtr->trajectoryStateClosestToPoint(cxPt);
178  if (!posTSCP.isValid() || !negTSCP.isValid()) continue;
179  if (posTSCP.momentum().dot(negTSCP.momentum()) < 0) continue;
180 
181  // calculate mPiPi
182  double totalE = sqrt(posTSCP.momentum().mag2() + piMassSquared) + sqrt(negTSCP.momentum().mag2() + piMassSquared);
183  double totalESq = totalE*totalE;
184  double totalPSq = (posTSCP.momentum() + negTSCP.momentum()).mag2();
185  double mass = sqrt(totalESq - totalPSq);
186  if (mass > mPiPiCut_) continue;
187 
188  // Fill the vector of TransientTracks to send to KVF
189  std::vector<reco::TransientTrack> transTracks;
190  transTracks.reserve(2);
191  transTracks.push_back(*posTransTkPtr);
192  transTracks.push_back(*negTransTkPtr);
193 
194  // create the vertex fitter object and vertex the tracks
195  TransientVertex theRecoVertex;
196  if (vertexFitter_) {
197  KalmanVertexFitter theKalmanFitter(useRefTracks_ == 0 ? false : true);
198  theRecoVertex = theKalmanFitter.vertex(transTracks);
199  } else if (!vertexFitter_) {
200  useRefTracks_ = false;
201  AdaptiveVertexFitter theAdaptiveFitter;
202  theRecoVertex = theAdaptiveFitter.vertex(transTracks);
203  }
204  if (!theRecoVertex.isValid()) continue;
205 
206  reco::Vertex theVtx = theRecoVertex;
207  if (theVtx.normalizedChi2() > vtxChi2Cut_) continue;
208  GlobalPoint vtxPos(theVtx.x(), theVtx.y(), theVtx.z());
209 
210  // 2D decay significance
211  SMatrixSym3D totalCov = theBeamSpot->rotatedCovariance3D() + theVtx.covariance();
212  if (useVertex_) totalCov = referenceVtx.covariance() + theVtx.covariance();
213  SVector3 distVecXY(vtxPos.x()-referencePos.x(), vtxPos.y()-referencePos.y(), 0.);
214  double distMagXY = ROOT::Math::Mag(distVecXY);
215  double sigmaDistMagXY = sqrt(ROOT::Math::Similarity(totalCov, distVecXY)) / distMagXY;
216  if (distMagXY/sigmaDistMagXY < vtxDecaySigXYCut_) continue;
217 
218  // 3D decay significance
219  if (vtxDecaySigXYZCut_ > 0.) {
220  SVector3 distVecXYZ(vtxPos.x()-referencePos.x(), vtxPos.y()-referencePos.y(), vtxPos.z()-referencePos.z());
221  double distMagXYZ = ROOT::Math::Mag(distVecXYZ);
222  double sigmaDistMagXYZ = sqrt(ROOT::Math::Similarity(totalCov, distVecXYZ)) / distMagXYZ;
223  if (distMagXYZ/sigmaDistMagXYZ < vtxDecaySigXYZCut_) continue;
224  }
225 
226  // make sure the vertex radius is within the inner track hit radius
227  if (innerHitPosCut_ > 0. && positiveTrackRef->innerOk()) {
228  reco::Vertex::Point posTkHitPos = positiveTrackRef->innerPosition();
229  double posTkHitPosD2 = (posTkHitPos.x()-referencePos.x())*(posTkHitPos.x()-referencePos.x()) +
230  (posTkHitPos.y()-referencePos.y())*(posTkHitPos.y()-referencePos.y());
231  if (sqrt(posTkHitPosD2) < (distMagXY - sigmaDistMagXY*innerHitPosCut_)) continue;
232  }
233  if (innerHitPosCut_ > 0. && negativeTrackRef->innerOk()) {
234  reco::Vertex::Point negTkHitPos = negativeTrackRef->innerPosition();
235  double negTkHitPosD2 = (negTkHitPos.x()-referencePos.x())*(negTkHitPos.x()-referencePos.x()) +
236  (negTkHitPos.y()-referencePos.y())*(negTkHitPos.y()-referencePos.y());
237  if (sqrt(negTkHitPosD2) < (distMagXY - sigmaDistMagXY*innerHitPosCut_)) continue;
238  }
239 
240  std::unique_ptr<TrajectoryStateClosestToPoint> trajPlus;
241  std::unique_ptr<TrajectoryStateClosestToPoint> trajMins;
242  std::vector<reco::TransientTrack> theRefTracks;
243  if (theRecoVertex.hasRefittedTracks()) {
244  theRefTracks = theRecoVertex.refittedTracks();
245  }
246 
247  if (useRefTracks_ && theRefTracks.size() > 1) {
248  reco::TransientTrack* thePositiveRefTrack = nullptr;
249  reco::TransientTrack* theNegativeRefTrack = nullptr;
250  for (std::vector<reco::TransientTrack>::iterator iTrack = theRefTracks.begin(); iTrack != theRefTracks.end(); ++iTrack) {
251  if (iTrack->track().charge() > 0.) {
252  thePositiveRefTrack = &*iTrack;
253  } else if (iTrack->track().charge() < 0.) {
254  theNegativeRefTrack = &*iTrack;
255  }
256  }
257  if (thePositiveRefTrack == nullptr || theNegativeRefTrack == nullptr) continue;
258  trajPlus.reset(new TrajectoryStateClosestToPoint(thePositiveRefTrack->trajectoryStateClosestToPoint(vtxPos)));
259  trajMins.reset(new TrajectoryStateClosestToPoint(theNegativeRefTrack->trajectoryStateClosestToPoint(vtxPos)));
260  } else {
261  trajPlus.reset(new TrajectoryStateClosestToPoint(posTransTkPtr->trajectoryStateClosestToPoint(vtxPos)));
262  trajMins.reset(new TrajectoryStateClosestToPoint(negTransTkPtr->trajectoryStateClosestToPoint(vtxPos)));
263  }
264 
265  if (trajPlus.get() == nullptr || trajMins.get() == nullptr || !trajPlus->isValid() || !trajMins->isValid()) continue;
266 
267  GlobalVector positiveP(trajPlus->momentum());
268  GlobalVector negativeP(trajMins->momentum());
269  GlobalVector totalP(positiveP + negativeP);
270 
271  // 2D pointing angle
272  double dx = theVtx.x()-referencePos.x();
273  double dy = theVtx.y()-referencePos.y();
274  double px = totalP.x();
275  double py = totalP.y();
276  double angleXY = (dx*px+dy*py)/(sqrt(dx*dx+dy*dy)*sqrt(px*px+py*py));
277  if (angleXY < cosThetaXYCut_) continue;
278 
279  // 3D pointing angle
280  if (cosThetaXYZCut_ > -1.) {
281  double dz = theVtx.z()-referencePos.z();
282  double pz = totalP.z();
283  double angleXYZ = (dx*px+dy*py+dz*pz)/(sqrt(dx*dx+dy*dy+dz*dz)*sqrt(px*px+py*py+pz*pz));
284  if (angleXYZ < cosThetaXYZCut_) continue;
285  }
286 
287  // calculate total energy of V0 3 ways: assume it's a kShort, a Lambda, or a LambdaBar.
288  double piPlusE = sqrt(positiveP.mag2() + piMassSquared);
289  double piMinusE = sqrt(negativeP.mag2() + piMassSquared);
290  double protonE = sqrt(positiveP.mag2() + protonMassSquared);
291  double antiProtonE = sqrt(negativeP.mag2() + protonMassSquared);
292  double kShortETot = piPlusE + piMinusE;
293  double lambdaEtot = protonE + piMinusE;
294  double lambdaBarEtot = antiProtonE + piPlusE;
295 
296  // Create momentum 4-vectors for the 3 candidate types
297  const reco::Particle::LorentzVector kShortP4(totalP.x(), totalP.y(), totalP.z(), kShortETot);
298  const reco::Particle::LorentzVector lambdaP4(totalP.x(), totalP.y(), totalP.z(), lambdaEtot);
299  const reco::Particle::LorentzVector lambdaBarP4(totalP.x(), totalP.y(), totalP.z(), lambdaBarEtot);
300 
301  reco::Particle::Point vtx(theVtx.x(), theVtx.y(), theVtx.z());
302  const reco::Vertex::CovarianceMatrix vtxCov(theVtx.covariance());
303  double vtxChi2(theVtx.chi2());
304  double vtxNdof(theVtx.ndof());
305 
306  // Create the VertexCompositeCandidate object that will be stored in the Event
307  reco::VertexCompositeCandidate* theKshort = nullptr;
308  reco::VertexCompositeCandidate* theLambda = nullptr;
309  reco::VertexCompositeCandidate* theLambdaBar = nullptr;
310 
311  if (doKShorts_) {
312  theKshort = new reco::VertexCompositeCandidate(0, kShortP4, vtx, vtxCov, vtxChi2, vtxNdof);
313  }
314  if (doLambdas_) {
315  if (positiveP.mag2() > negativeP.mag2()) {
316  theLambda = new reco::VertexCompositeCandidate(0, lambdaP4, vtx, vtxCov, vtxChi2, vtxNdof);
317  } else {
318  theLambdaBar = new reco::VertexCompositeCandidate(0, lambdaBarP4, vtx, vtxCov, vtxChi2, vtxNdof);
319  }
320  }
321 
322  // Create daughter candidates for the VertexCompositeCandidates
323  reco::RecoChargedCandidate thePiPlusCand(
324  1, reco::Particle::LorentzVector(positiveP.x(), positiveP.y(), positiveP.z(), piPlusE), vtx);
325  thePiPlusCand.setTrack(positiveTrackRef);
326 
327  reco::RecoChargedCandidate thePiMinusCand(
328  -1, reco::Particle::LorentzVector(negativeP.x(), negativeP.y(), negativeP.z(), piMinusE), vtx);
329  thePiMinusCand.setTrack(negativeTrackRef);
330 
331  reco::RecoChargedCandidate theProtonCand(
332  1, reco::Particle::LorentzVector(positiveP.x(), positiveP.y(), positiveP.z(), protonE), vtx);
333  theProtonCand.setTrack(positiveTrackRef);
334 
335  reco::RecoChargedCandidate theAntiProtonCand(
336  -1, reco::Particle::LorentzVector(negativeP.x(), negativeP.y(), negativeP.z(), antiProtonE), vtx);
337  theAntiProtonCand.setTrack(negativeTrackRef);
338 
339  AddFourMomenta addp4;
340  // Store the daughter Candidates in the VertexCompositeCandidates if they pass mass cuts
341  if (doKShorts_) {
342  theKshort->addDaughter(thePiPlusCand);
343  theKshort->addDaughter(thePiMinusCand);
344  theKshort->setPdgId(310);
345  addp4.set(*theKshort);
346  if (theKshort->mass() < kShortMass + kShortMassCut_ && theKshort->mass() > kShortMass - kShortMassCut_) {
347  theKshorts.push_back(std::move(*theKshort));
348  }
349  }
350  if (doLambdas_ && theLambda) {
351  theLambda->addDaughter(theProtonCand);
352  theLambda->addDaughter(thePiMinusCand);
353  theLambda->setPdgId(3122);
354  addp4.set( *theLambda );
355  if (theLambda->mass() < lambdaMass + lambdaMassCut_ && theLambda->mass() > lambdaMass - lambdaMassCut_) {
356  theLambdas.push_back(std::move(*theLambda));
357  }
358  } else if (doLambdas_ && theLambdaBar) {
359  theLambdaBar->addDaughter(theAntiProtonCand);
360  theLambdaBar->addDaughter(thePiPlusCand);
361  theLambdaBar->setPdgId(-3122);
362  addp4.set(*theLambdaBar);
363  if (theLambdaBar->mass() < lambdaMass + lambdaMassCut_ && theLambdaBar->mass() > lambdaMass - lambdaMassCut_) {
364  theLambdas.push_back(std::move(*theLambdaBar));
365  }
366  }
367 
368  delete theKshort;
369  delete theLambda;
370  delete theLambdaBar;
371  theKshort = theLambda = theLambdaBar = nullptr;
372 
373  }
374  }
375 }
376 
bool doLambdas_
Definition: V0Fitter.h:52
double tkChi2Cut_
Definition: V0Fitter.h:55
T getParameter(std::string const &) const
T mag2() const
Definition: PV3DBase.h:66
void reset()
Definition: ProxyBase11.h:52
float distance() const override
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:69
double lambdaMassCut_
Definition: V0Fitter.h:72
double tkIPSigXYCut_
Definition: V0Fitter.h:58
bool vertexFitter_
Definition: V0Fitter.h:49
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:600
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
TrajectoryStateClosestToPoint impactPointTSCP() const
double dxyError() const
error on dxy
Definition: TrackBase.h:847
edm::EDGetTokenT< std::vector< reco::Vertex > > token_vertices
Definition: V0Fitter.h:77
double y() const
y coordinate
Definition: Vertex.h:113
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:15
T y() const
Definition: PV3DBase.h:63
CachingVertex< 5 > vertex(const std::vector< reco::TransientTrack > &) const override
PreciseFloatType< T, U >::Type dot(const Vector3DBase< U, FrameTag > &v) const
Definition: Vector3DBase.h:107
double covariance(int i, int j) const
(i, j)-th element of error matrix, i, j = 0, ... 2
Definition: Vertex.h:130
bool hasRefittedTracks() const
double vtxDecaySigXYZCut_
Definition: V0Fitter.h:63
math::Error< dimension >::type CovarianceMatrix
covariance error matrix (3x3)
Definition: Vertex.h:45
double tkDCACut_
Definition: V0Fitter.h:65
int tkNHitsCut_
Definition: V0Fitter.h:56
int iEvent
Definition: GenABIO.cc:224
bool status() const override
bool useRefTracks_
Definition: V0Fitter.h:50
edm::EDGetTokenT< reco::TrackCollection > token_tracks
Definition: V0Fitter.h:74
bool calculate(const TrajectoryStateOnSurface &sta, const TrajectoryStateOnSurface &stb) override
T sqrt(T t)
Definition: SSEVec.h:18
double pt() const
track transverse momentum
Definition: TrackBase.h:660
T z() const
Definition: PV3DBase.h:64
T mag2() const
The vector magnitude squared. Equivalent to vec.dot(vec)
math::XYZPoint Point
point in the space
Definition: Particle.h:25
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double chi2() const
chi-squares
Definition: Vertex.h:98
double z() const
z coordinate
Definition: Vertex.h:115
void fitAll(const edm::Event &iEvent, const edm::EventSetup &iSetup, reco::VertexCompositeCandidateCollection &k, reco::VertexCompositeCandidateCollection &l)
Definition: V0Fitter.cc:87
CachingVertex< 5 > vertex(const std::vector< reco::TransientTrack > &tracks) const override
bool useVertex_
Definition: V0Fitter.h:76
math::XYZPoint Point
point in the space
Definition: Vertex.h:39
unsigned short numberOfValidHits() const
number of valid hits found
Definition: TrackBase.h:901
double mPiPiCut_
Definition: V0Fitter.h:66
double kShortMassCut_
Definition: V0Fitter.h:71
V0Fitter(const edm::ParameterSet &theParams, edm::ConsumesCollector &&iC)
Definition: V0Fitter.cc:49
double tkIPSigZCut_
Definition: V0Fitter.h:59
void addDaughter(const Candidate &, const std::string &s="")
add a clone of the passed candidate as daughter
double ndof() const
Definition: Vertex.h:105
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:648
double dzError() const
error on dz
Definition: TrackBase.h:865
double cosThetaXYCut_
Definition: V0Fitter.h:68
double x() const
x coordinate
Definition: Vertex.h:111
double innerHitPosCut_
Definition: V0Fitter.h:67
double vtxDecaySigXYCut_
Definition: V0Fitter.h:62
T const * product() const
Definition: Handle.h:74
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
TrajectoryStateClosestToPoint trajectoryStateClosestToPoint(const GlobalPoint &point) const
double tkPtCut_
Definition: V0Fitter.h:57
T get() const
Definition: EventSetup.h:71
double vtxChi2Cut_
Definition: V0Fitter.h:61
const Point & position() const
position
Definition: BeamSpot.h:62
edm::EDGetTokenT< reco::BeamSpot > token_beamSpot
Definition: V0Fitter.h:75
std::vector< reco::TransientTrack > const & refittedTracks() const
void set(reco::Candidate &c) const
set up a candidate
double normalizedChi2() const
chi-squared divided by n.d.o.f.
Definition: Vertex.h:107
Covariance3DMatrix rotatedCovariance3D() const
Definition: BeamSpot.cc:78
ROOT::Math::SVector< double, 3 > SVector3
Definition: V0Fitter.cc:47
GlobalPoint crossingPoint() const override
void setTrack(const reco::TrackRef &r)
set reference to track
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:630
T x() const
Definition: PV3DBase.h:62
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Particle.h:21
void setPdgId(int pdgId) final
T const * product() const
Definition: ESHandle.h:86
bool isValid() const
def move(src, dest)
Definition: eostools.py:511
double mass() const final
mass
bool doKShorts_
Definition: V0Fitter.h:51