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  token_vertices = iC.consumes<std::vector<reco::Vertex>>(theParameters.getParameter<edm::InputTag>("vertices"));
53 
54  token_tracks = iC.consumes<reco::TrackCollection>(theParameters.getParameter<edm::InputTag>("trackRecoAlgorithm"));
55  vertexFitter_ = theParameters.getParameter<bool>("vertexFitter");
56  useRefTracks_ = theParameters.getParameter<bool>("useRefTracks");
57 
58  // whether to reconstruct KShorts
59  doKShorts_ = theParameters.getParameter<bool>("doKShorts");
60  // whether to reconstruct Lambdas
61  doLambdas_ = theParameters.getParameter<bool>("doLambdas");
62 
63  // cuts on initial track selection
64  tkChi2Cut_ = theParameters.getParameter<double>("tkChi2Cut");
65  tkNHitsCut_ = theParameters.getParameter<int>("tkNHitsCut");
66  tkPtCut_ = theParameters.getParameter<double>("tkPtCut");
67  tkIPSigXYCut_ = theParameters.getParameter<double>("tkIPSigXYCut");
68  tkIPSigZCut_ = theParameters.getParameter<double>("tkIPSigZCut");
69 
70  // cuts on vertex
71  vtxChi2Cut_ = theParameters.getParameter<double>("vtxChi2Cut");
72  vtxDecaySigXYZCut_ = theParameters.getParameter<double>("vtxDecaySigXYZCut");
73  vtxDecaySigXYCut_ = theParameters.getParameter<double>("vtxDecaySigXYCut");
74  // miscellaneous cuts
75  tkDCACut_ = theParameters.getParameter<double>("tkDCACut");
76  mPiPiCut_ = theParameters.getParameter<double>("mPiPiCut");
77  innerHitPosCut_ = theParameters.getParameter<double>("innerHitPosCut");
78  cosThetaXYCut_ = theParameters.getParameter<double>("cosThetaXYCut");
79  cosThetaXYZCut_ = theParameters.getParameter<double>("cosThetaXYZCut");
80  // cuts on the V0 candidate mass
81  kShortMassCut_ = theParameters.getParameter<double>("kShortMassCut");
82  lambdaMassCut_ = theParameters.getParameter<double>("lambdaMassCut");
83 }
84 
85 // method containing the algorithm for vertex reconstruction
87  const edm::EventSetup& iSetup,
90  using std::vector;
91 
93  iEvent.getByToken(token_tracks, theTrackHandle);
94  if (theTrackHandle->empty())
95  return;
96  const reco::TrackCollection* theTrackCollection = theTrackHandle.product();
97 
98  edm::Handle<reco::BeamSpot> theBeamSpotHandle;
99  iEvent.getByToken(token_beamSpot, theBeamSpotHandle);
100  const reco::BeamSpot* theBeamSpot = theBeamSpotHandle.product();
101  math::XYZPoint referencePos(theBeamSpot->position());
102 
103  reco::Vertex referenceVtx;
104  if (useVertex_) {
106  iEvent.getByToken(token_vertices, vertices);
107  referenceVtx = vertices->at(0);
108  referencePos = referenceVtx.position();
109  }
110 
111  const MagneticField* theMagneticField = &iSetup.getData(esTokenMF_);
112 
113  std::vector<reco::TrackRef> theTrackRefs;
114  std::vector<reco::TransientTrack> theTransTracks;
115 
116  // fill vectors of TransientTracks and TrackRefs after applying preselection cuts
117  for (reco::TrackCollection::const_iterator iTk = theTrackCollection->begin(); iTk != theTrackCollection->end();
118  ++iTk) {
119  const reco::Track* tmpTrack = &(*iTk);
120  double ipsigXY = std::abs(tmpTrack->dxy(*theBeamSpot) / tmpTrack->dxyError());
121  if (useVertex_)
122  ipsigXY = std::abs(tmpTrack->dxy(referencePos) / tmpTrack->dxyError());
123  double ipsigZ = std::abs(tmpTrack->dz(referencePos) / tmpTrack->dzError());
124  if (tmpTrack->normalizedChi2() < tkChi2Cut_ && tmpTrack->numberOfValidHits() >= tkNHitsCut_ &&
125  tmpTrack->pt() > tkPtCut_ && ipsigXY > tkIPSigXYCut_ && ipsigZ > tkIPSigZCut_) {
126  reco::TrackRef tmpRef(theTrackHandle, std::distance(theTrackCollection->begin(), iTk));
127  theTrackRefs.push_back(std::move(tmpRef));
128  reco::TransientTrack tmpTransient(*tmpRef, theMagneticField);
129  theTransTracks.push_back(std::move(tmpTransient));
130  }
131  }
132  // good tracks have now been selected for vertexing
133 
134  // loop over tracks and vertex good charged track pairs
135  for (unsigned int trdx1 = 0; trdx1 < theTrackRefs.size(); ++trdx1) {
136  for (unsigned int trdx2 = trdx1 + 1; trdx2 < theTrackRefs.size(); ++trdx2) {
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())
163  continue;
164  FreeTrajectoryState const& posState = posImpact.theState();
165  FreeTrajectoryState const& negState = negImpact.theState();
167  cApp.calculate(posState, negState);
168  if (!cApp.status())
169  continue;
170  float dca = std::abs(cApp.distance());
171  if (dca > tkDCACut_)
172  continue;
173 
174  // the POCA should at least be in the sensitive volume
175  GlobalPoint cxPt = cApp.crossingPoint();
176  if ((cxPt.x() * cxPt.x() + cxPt.y() * cxPt.y()) > 120. * 120. || std::abs(cxPt.z()) > 300.)
177  continue;
178 
179  // the tracks should at least point in the same quadrant
180  TrajectoryStateClosestToPoint const& posTSCP = posTransTkPtr->trajectoryStateClosestToPoint(cxPt);
181  TrajectoryStateClosestToPoint const& negTSCP = negTransTkPtr->trajectoryStateClosestToPoint(cxPt);
182  if (!posTSCP.isValid() || !negTSCP.isValid())
183  continue;
184  if (posTSCP.momentum().dot(negTSCP.momentum()) < 0)
185  continue;
186 
187  // calculate mPiPi
188  double totalE = sqrt(posTSCP.momentum().mag2() + piMassSquared) + sqrt(negTSCP.momentum().mag2() + piMassSquared);
189  double totalESq = totalE * totalE;
190  double totalPSq = (posTSCP.momentum() + negTSCP.momentum()).mag2();
191  double massSquared = totalESq - totalPSq;
192  if (massSquared > mPiPiCut_ * mPiPiCut_)
193  continue;
194 
195  // Fill the vector of TransientTracks to send to KVF
196  std::vector<reco::TransientTrack> transTracks;
197  transTracks.reserve(2);
198  transTracks.push_back(*posTransTkPtr);
199  transTracks.push_back(*negTransTkPtr);
200 
201  // create the vertex fitter object and vertex the tracks
202  TransientVertex theRecoVertex;
203  if (vertexFitter_) {
204  KalmanVertexFitter theKalmanFitter(useRefTracks_ == 0 ? false : true);
205  theRecoVertex = theKalmanFitter.vertex(transTracks);
206  } else if (!vertexFitter_) {
207  useRefTracks_ = false;
208  AdaptiveVertexFitter theAdaptiveFitter;
209  theRecoVertex = theAdaptiveFitter.vertex(transTracks);
210  }
211  if (!theRecoVertex.isValid())
212  continue;
213 
214  reco::Vertex theVtx = theRecoVertex;
215  if (theVtx.normalizedChi2() > vtxChi2Cut_)
216  continue;
217  GlobalPoint vtxPos(theVtx.x(), theVtx.y(), theVtx.z());
218 
219  // 2D decay significance
220  SMatrixSym3D totalCov = theBeamSpot->rotatedCovariance3D() + theVtx.covariance();
221  if (useVertex_)
222  totalCov = referenceVtx.covariance() + theVtx.covariance();
223  SVector3 distVecXY(vtxPos.x() - referencePos.x(), vtxPos.y() - referencePos.y(), 0.);
224  double distMagXY = ROOT::Math::Mag(distVecXY);
225  double sigmaDistMagXY = sqrt(ROOT::Math::Similarity(totalCov, distVecXY)) / distMagXY;
226  if (distMagXY / sigmaDistMagXY < vtxDecaySigXYCut_)
227  continue;
228 
229  // 3D decay significance
230  if (vtxDecaySigXYZCut_ > 0.) {
231  SVector3 distVecXYZ(
232  vtxPos.x() - referencePos.x(), vtxPos.y() - referencePos.y(), vtxPos.z() - referencePos.z());
233  double distMagXYZ = ROOT::Math::Mag(distVecXYZ);
234  double sigmaDistMagXYZ = sqrt(ROOT::Math::Similarity(totalCov, distVecXYZ)) / distMagXYZ;
235  if (distMagXYZ / sigmaDistMagXYZ < vtxDecaySigXYZCut_)
236  continue;
237  }
238 
239  // make sure the vertex radius is within the inner track hit radius
240  double tkHitPosLimitSquared =
241  (distMagXY - sigmaDistMagXY * innerHitPosCut_) * (distMagXY - sigmaDistMagXY * innerHitPosCut_);
242  if (innerHitPosCut_ > 0. && positiveTrackRef->innerOk()) {
243  reco::Vertex::Point posTkHitPos = positiveTrackRef->innerPosition();
244  double posTkHitPosD2 = (posTkHitPos.x() - referencePos.x()) * (posTkHitPos.x() - referencePos.x()) +
245  (posTkHitPos.y() - referencePos.y()) * (posTkHitPos.y() - referencePos.y());
246  if (posTkHitPosD2 < tkHitPosLimitSquared)
247  continue;
248  }
249  if (innerHitPosCut_ > 0. && negativeTrackRef->innerOk()) {
250  reco::Vertex::Point negTkHitPos = negativeTrackRef->innerPosition();
251  double negTkHitPosD2 = (negTkHitPos.x() - referencePos.x()) * (negTkHitPos.x() - referencePos.x()) +
252  (negTkHitPos.y() - referencePos.y()) * (negTkHitPos.y() - referencePos.y());
253  if (negTkHitPosD2 < tkHitPosLimitSquared)
254  continue;
255  }
256 
257  std::unique_ptr<TrajectoryStateClosestToPoint> trajPlus;
258  std::unique_ptr<TrajectoryStateClosestToPoint> trajMins;
259  std::vector<reco::TransientTrack> theRefTracks;
260  if (theRecoVertex.hasRefittedTracks()) {
261  theRefTracks = theRecoVertex.refittedTracks();
262  }
263 
264  if (useRefTracks_ && theRefTracks.size() > 1) {
265  reco::TransientTrack* thePositiveRefTrack = nullptr;
266  reco::TransientTrack* theNegativeRefTrack = nullptr;
267  for (std::vector<reco::TransientTrack>::iterator iTrack = theRefTracks.begin(); iTrack != theRefTracks.end();
268  ++iTrack) {
269  if (iTrack->track().charge() > 0.) {
270  thePositiveRefTrack = &*iTrack;
271  } else if (iTrack->track().charge() < 0.) {
272  theNegativeRefTrack = &*iTrack;
273  }
274  }
275  if (thePositiveRefTrack == nullptr || theNegativeRefTrack == nullptr)
276  continue;
277  trajPlus =
278  std::make_unique<TrajectoryStateClosestToPoint>(thePositiveRefTrack->trajectoryStateClosestToPoint(vtxPos));
279  trajMins =
280  std::make_unique<TrajectoryStateClosestToPoint>(theNegativeRefTrack->trajectoryStateClosestToPoint(vtxPos));
281  } else {
282  trajPlus =
283  std::make_unique<TrajectoryStateClosestToPoint>(posTransTkPtr->trajectoryStateClosestToPoint(vtxPos));
284  trajMins =
285  std::make_unique<TrajectoryStateClosestToPoint>(negTransTkPtr->trajectoryStateClosestToPoint(vtxPos));
286  }
287 
288  if (trajPlus.get() == nullptr || trajMins.get() == nullptr || !trajPlus->isValid() || !trajMins->isValid())
289  continue;
290 
291  GlobalVector positiveP(trajPlus->momentum());
292  GlobalVector negativeP(trajMins->momentum());
293  GlobalVector totalP(positiveP + negativeP);
294 
295  // 2D pointing angle
296  double dx = theVtx.x() - referencePos.x();
297  double dy = theVtx.y() - referencePos.y();
298  double px = totalP.x();
299  double py = totalP.y();
300  double angleXY = (dx * px + dy * py) / (sqrt(dx * dx + dy * dy) * sqrt(px * px + py * py));
301  if (angleXY < cosThetaXYCut_)
302  continue;
303 
304  // 3D pointing angle
305  if (cosThetaXYZCut_ > -1.) {
306  double dz = theVtx.z() - referencePos.z();
307  double pz = totalP.z();
308  double angleXYZ =
309  (dx * px + dy * py + dz * pz) / (sqrt(dx * dx + dy * dy + dz * dz) * sqrt(px * px + py * py + pz * pz));
310  if (angleXYZ < cosThetaXYZCut_)
311  continue;
312  }
313 
314  // calculate total energy of V0 3 ways: assume it's a kShort, a Lambda, or a LambdaBar.
315  double piPlusE = sqrt(positiveP.mag2() + piMassSquared);
316  double piMinusE = sqrt(negativeP.mag2() + piMassSquared);
317  double protonE = sqrt(positiveP.mag2() + protonMassSquared);
318  double antiProtonE = sqrt(negativeP.mag2() + protonMassSquared);
319  double kShortETot = piPlusE + piMinusE;
320  double lambdaEtot = protonE + piMinusE;
321  double lambdaBarEtot = antiProtonE + piPlusE;
322 
323  // Create momentum 4-vectors for the 3 candidate types
324  const reco::Particle::LorentzVector kShortP4(totalP.x(), totalP.y(), totalP.z(), kShortETot);
325  const reco::Particle::LorentzVector lambdaP4(totalP.x(), totalP.y(), totalP.z(), lambdaEtot);
326  const reco::Particle::LorentzVector lambdaBarP4(totalP.x(), totalP.y(), totalP.z(), lambdaBarEtot);
327 
328  reco::Particle::Point vtx(theVtx.x(), theVtx.y(), theVtx.z());
329  const reco::Vertex::CovarianceMatrix vtxCov(theVtx.covariance());
330  double vtxChi2(theVtx.chi2());
331  double vtxNdof(theVtx.ndof());
332 
333  // Create the VertexCompositeCandidate object that will be stored in the Event
334  reco::VertexCompositeCandidate* theKshort = nullptr;
335  reco::VertexCompositeCandidate* theLambda = nullptr;
336  reco::VertexCompositeCandidate* theLambdaBar = nullptr;
337 
338  if (doKShorts_) {
339  theKshort = new reco::VertexCompositeCandidate(0, kShortP4, vtx, vtxCov, vtxChi2, vtxNdof);
340  }
341  if (doLambdas_) {
342  if (positiveP.mag2() > negativeP.mag2()) {
343  theLambda = new reco::VertexCompositeCandidate(0, lambdaP4, vtx, vtxCov, vtxChi2, vtxNdof);
344  } else {
345  theLambdaBar = new reco::VertexCompositeCandidate(0, lambdaBarP4, vtx, vtxCov, vtxChi2, vtxNdof);
346  }
347  }
348 
349  // Create daughter candidates for the VertexCompositeCandidates
350  reco::RecoChargedCandidate thePiPlusCand(
351  1, reco::Particle::LorentzVector(positiveP.x(), positiveP.y(), positiveP.z(), piPlusE), vtx);
352  thePiPlusCand.setTrack(positiveTrackRef);
353 
354  reco::RecoChargedCandidate thePiMinusCand(
355  -1, reco::Particle::LorentzVector(negativeP.x(), negativeP.y(), negativeP.z(), piMinusE), vtx);
356  thePiMinusCand.setTrack(negativeTrackRef);
357 
358  reco::RecoChargedCandidate theProtonCand(
359  1, reco::Particle::LorentzVector(positiveP.x(), positiveP.y(), positiveP.z(), protonE), vtx);
360  theProtonCand.setTrack(positiveTrackRef);
361 
362  reco::RecoChargedCandidate theAntiProtonCand(
363  -1, reco::Particle::LorentzVector(negativeP.x(), negativeP.y(), negativeP.z(), antiProtonE), vtx);
364  theAntiProtonCand.setTrack(negativeTrackRef);
365 
366  AddFourMomenta addp4;
367  // Store the daughter Candidates in the VertexCompositeCandidates if they pass mass cuts
368  if (doKShorts_) {
369  theKshort->addDaughter(thePiPlusCand);
370  theKshort->addDaughter(thePiMinusCand);
371  theKshort->setPdgId(310);
372  addp4.set(*theKshort);
373  if (theKshort->mass() < kShortMass + kShortMassCut_ && theKshort->mass() > kShortMass - kShortMassCut_) {
374  theKshorts.push_back(std::move(*theKshort));
375  }
376  }
377  if (doLambdas_ && theLambda) {
378  theLambda->addDaughter(theProtonCand);
379  theLambda->addDaughter(thePiMinusCand);
380  theLambda->setPdgId(3122);
381  addp4.set(*theLambda);
382  if (theLambda->mass() < lambdaMass + lambdaMassCut_ && theLambda->mass() > lambdaMass - lambdaMassCut_) {
383  theLambdas.push_back(std::move(*theLambda));
384  }
385  } else if (doLambdas_ && theLambdaBar) {
386  theLambdaBar->addDaughter(theAntiProtonCand);
387  theLambdaBar->addDaughter(thePiPlusCand);
388  theLambdaBar->setPdgId(-3122);
389  addp4.set(*theLambdaBar);
390  if (theLambdaBar->mass() < lambdaMass + lambdaMassCut_ && theLambdaBar->mass() > lambdaMass - lambdaMassCut_) {
391  theLambdas.push_back(std::move(*theLambdaBar));
392  }
393  }
394 
395  delete theKshort;
396  delete theLambda;
397  delete theLambdaBar;
398  theKshort = theLambda = theLambdaBar = nullptr;
399  }
400  }
401 }
bool doLambdas_
Definition: V0Fitter.h:54
double tkChi2Cut_
Definition: V0Fitter.h:57
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:71
double lambdaMassCut_
Definition: V0Fitter.h:74
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:60
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:51
T z() const
Definition: PV3DBase.h:61
edm::EDGetTokenT< std::vector< reco::Vertex > > token_vertices
Definition: V0Fitter.h:79
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
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:65
math::Error< dimension >::type CovarianceMatrix
covariance error matrix (3x3)
Definition: Vertex.h:46
double tkDCACut_
Definition: V0Fitter.h:67
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
T y() const
Definition: PV3DBase.h:60
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:58
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:52
double dxyError() const
error on dxy
Definition: TrackBase.h:769
edm::EDGetTokenT< reco::TrackCollection > token_tracks
Definition: V0Fitter.h:76
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:86
bool useVertex_
Definition: V0Fitter.h:78
math::XYZPoint Point
point in the space
Definition: Vertex.h:40
Covariance3DMatrix rotatedCovariance3D() const
Definition: BeamSpot.cc:73
double mPiPiCut_
Definition: V0Fitter.h:68
bool getData(T &iHolder) const
Definition: EventSetup.h:122
double x() const
x coordinate
Definition: Vertex.h:129
double kShortMassCut_
Definition: V0Fitter.h:73
V0Fitter(const edm::ParameterSet &theParams, edm::ConsumesCollector &&iC)
Definition: V0Fitter.cc:49
double tkIPSigZCut_
Definition: V0Fitter.h:61
double y() const
y coordinate
Definition: Vertex.h:131
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:70
double innerHitPosCut_
Definition: V0Fitter.h:69
double vtxDecaySigXYCut_
Definition: V0Fitter.h:64
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 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 chi2() const
chi-squares
Definition: Vertex.h:116
double tkPtCut_
Definition: V0Fitter.h:59
double mass() const final
mass
double vtxChi2Cut_
Definition: V0Fitter.h:63
edm::EDGetTokenT< reco::BeamSpot > token_beamSpot
Definition: V0Fitter.h:77
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
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:53