CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Attributes
V0Fitter Class Reference

#include <RecoVertex/V0Producer/interface/V0Fitter.h>

Public Member Functions

void fitAll (const edm::Event &iEvent, const edm::EventSetup &iSetup, reco::VertexCompositeCandidateCollection &k, reco::VertexCompositeCandidateCollection &l)
 
 V0Fitter (const edm::ParameterSet &theParams, edm::ConsumesCollector &&iC)
 

Private Attributes

bool doKShorts_
 
bool doLambdas_
 
double innerHitPosCut_
 
double kShortMassCut_
 
double lambdaMassCut_
 
double mPiPiCut_
 
double tkChi2Cut_
 
double tkDCACut_
 
double tkIPSigCut_
 
int tkNHitsCut_
 
double tkPtCut_
 
edm::EDGetTokenT< reco::BeamSpottoken_beamSpot
 
edm::EDGetTokenT
< reco::TrackCollection
token_tracks
 
bool useRefTracks_
 
double v0CosThetaCut_
 
bool vertexFitter_
 
double vtxChi2Cut_
 
double vtxDecayRSigCut_
 

Detailed Description

Description: <one line="" class="" summary>="">

Implementation: <Notes on="" implementation>="">

Definition at line 57 of file V0Fitter.h.

Constructor & Destructor Documentation

V0Fitter::V0Fitter ( const edm::ParameterSet theParams,
edm::ConsumesCollector &&  iC 
)

Definition at line 49 of file V0Fitter.cc.

References doKShorts_, doLambdas_, edm::ParameterSet::getParameter(), innerHitPosCut_, kShortMassCut_, lambdaMassCut_, mPiPiCut_, AlCaHLTBitMon_QueryRunRegistry::string, tkChi2Cut_, tkDCACut_, tkIPSigCut_, tkNHitsCut_, tkPtCut_, token_beamSpot, token_tracks, useRefTracks_, v0CosThetaCut_, vertexFitter_, vtxChi2Cut_, and vtxDecayRSigCut_.

50 {
51  using std::string;
52 
53  token_beamSpot = iC.consumes<reco::BeamSpot>(edm::InputTag("offlineBeamSpot"));
54  token_tracks = iC.consumes<reco::TrackCollection>(theParameters.getParameter<edm::InputTag>("trackRecoAlgorithm"));
55 
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  tkIPSigCut_ = theParameters.getParameter<double>("tkIPSigCut");
69  // cuts on vertex
70  vtxChi2Cut_ = theParameters.getParameter<double>("vtxChi2Cut");
71  vtxDecayRSigCut_ = theParameters.getParameter<double>("vtxDecayRSigCut");
72  // miscellaneous cuts
73  tkDCACut_ = theParameters.getParameter<double>("tkDCACut");
74  mPiPiCut_ = theParameters.getParameter<double>("mPiPiCut");
75  innerHitPosCut_ = theParameters.getParameter<double>("innerHitPosCut");
76  v0CosThetaCut_ = theParameters.getParameter<double>("v0CosThetaCut");
77  // cuts on the V0 candidate mass
78  kShortMassCut_ = theParameters.getParameter<double>("kShortMassCut");
79  lambdaMassCut_ = theParameters.getParameter<double>("lambdaMassCut");
80 }
bool doLambdas_
Definition: V0Fitter.h:70
double tkChi2Cut_
Definition: V0Fitter.h:73
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
double lambdaMassCut_
Definition: V0Fitter.h:87
double tkIPSigCut_
Definition: V0Fitter.h:76
bool vertexFitter_
Definition: V0Fitter.h:67
double v0CosThetaCut_
Definition: V0Fitter.h:84
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:13
double tkDCACut_
Definition: V0Fitter.h:81
int tkNHitsCut_
Definition: V0Fitter.h:74
bool useRefTracks_
Definition: V0Fitter.h:68
edm::EDGetTokenT< reco::TrackCollection > token_tracks
Definition: V0Fitter.h:89
double mPiPiCut_
Definition: V0Fitter.h:82
double kShortMassCut_
Definition: V0Fitter.h:86
double innerHitPosCut_
Definition: V0Fitter.h:83
double vtxDecayRSigCut_
Definition: V0Fitter.h:79
double tkPtCut_
Definition: V0Fitter.h:75
double vtxChi2Cut_
Definition: V0Fitter.h:78
edm::EDGetTokenT< reco::BeamSpot > token_beamSpot
Definition: V0Fitter.h:90
bool doKShorts_
Definition: V0Fitter.h:69

Member Function Documentation

void V0Fitter::fitAll ( const edm::Event iEvent,
const edm::EventSetup iSetup,
reco::VertexCompositeCandidateCollection k,
reco::VertexCompositeCandidateCollection l 
)

Definition at line 83 of file V0Fitter.cc.

References funct::abs(), reco::CompositeCandidate::addDaughter(), ClosestApproachInRPhi::calculate(), reco::Vertex::chi2(), reco::Vertex::covariance(), ClosestApproachInRPhi::crossingPoint(), ClosestApproachInRPhi::distance(), doKShorts_, doLambdas_, Vector3DBase< T, FrameTag >::dot(), reco::TrackBase::dxy(), reco::TrackBase::dxyError(), edm::EventSetup::get(), edm::Event::getByToken(), TransientVertex::hasRefittedTracks(), reco::TransientTrack::impactPointTSCP(), innerHitPosCut_, TrajectoryStateClosestToPoint::isValid(), TransientVertex::isValid(), kShortMassCut_, lambdaMassCut_, PV3DBase< T, PVType, FrameType >::mag2(), mag2(), reco::LeafCandidate::mass(), TrajectoryStateClosestToPoint::momentum(), mPiPiCut_, reco::Vertex::ndof(), reco::Vertex::normalizedChi2(), reco::TrackBase::normalizedChi2(), reco::TrackBase::numberOfValidHits(), reco::TrackBase::pt(), dt_dqm_sourceclient_common_cff::reco, TransientVertex::refittedTracks(), AddFourMomenta::set(), reco::LeafCandidate::setPdgId(), reco::RecoChargedCandidate::setTrack(), mathSSE::sqrt(), ClosestApproachInRPhi::status(), TrajectoryStateClosestToPoint::theState(), tkChi2Cut_, tkDCACut_, tkIPSigCut_, tkNHitsCut_, tkPtCut_, token_beamSpot, token_tracks, reco::TransientTrack::trajectoryStateClosestToPoint(), useRefTracks_, v0CosThetaCut_, KalmanVertexFitter::vertex(), AdaptiveVertexFitter::vertex(), vertexFitter_, vtxChi2Cut_, vtxDecayRSigCut_, PV3DBase< T, PVType, FrameType >::x(), reco::Vertex::x(), PV3DBase< T, PVType, FrameType >::y(), reco::Vertex::y(), PV3DBase< T, PVType, FrameType >::z(), and reco::Vertex::z().

Referenced by V0Producer::produce().

85 {
86  using std::vector;
87  using namespace reco;
88  using namespace edm;
89 
90  Handle<reco::TrackCollection> theTrackHandle;
91  iEvent.getByToken(token_tracks, theTrackHandle);
92  if (!theTrackHandle->size()) return;
93  const reco::TrackCollection* theTrackCollection = theTrackHandle.product();
94 
95  Handle<reco::BeamSpot> theBeamSpotHandle;
96  iEvent.getByToken(token_beamSpot, theBeamSpotHandle);
97  const reco::BeamSpot* theBeamSpot = theBeamSpotHandle.product();
98  const GlobalPoint beamSpotPos(theBeamSpot->position().x(), theBeamSpot->position().y(), theBeamSpot->position().z());
99 
100  ESHandle<MagneticField> theMagneticFieldHandle;
101  iSetup.get<IdealMagneticFieldRecord>().get(theMagneticFieldHandle);
102  const MagneticField* theMagneticField = theMagneticFieldHandle.product();
103 
104  std::vector<TrackRef> theTrackRefs;
105  std::vector<TransientTrack> theTransTracks;
106 
107  // fill vectors of TransientTracks and TrackRefs after applying preselection cuts
108  for (reco::TrackCollection::const_iterator iTk = theTrackCollection->begin(); iTk != theTrackCollection->end(); ++iTk) {
109  const reco::Track* tmpTrack = &(*iTk);
110  double ipsig = std::abs(tmpTrack->dxy(*theBeamSpot)/tmpTrack->dxyError());
111  if (tmpTrack->normalizedChi2() < tkChi2Cut_ && tmpTrack->numberOfValidHits() >= tkNHitsCut_ &&
112  tmpTrack->pt() > tkPtCut_ && ipsig > tkIPSigCut_) {
113  TrackRef tmpRef(theTrackHandle, std::distance(theTrackCollection->begin(), iTk));
114  theTrackRefs.push_back(std::move(tmpRef));
115  TransientTrack tmpTransient(*tmpRef, theMagneticField);
116  theTransTracks.push_back(std::move(tmpTransient));
117  }
118  }
119  // good tracks have now been selected for vertexing
120 
121  // loop over tracks and vertex good charged track pairs
122  for (unsigned int trdx1 = 0; trdx1 < theTrackRefs.size(); ++trdx1) {
123  for (unsigned int trdx2 = trdx1 + 1; trdx2 < theTrackRefs.size(); ++trdx2) {
124 
125  TrackRef positiveTrackRef;
126  TrackRef negativeTrackRef;
127  TransientTrack* posTransTkPtr = nullptr;
128  TransientTrack* negTransTkPtr = nullptr;
129 
130  if (theTrackRefs[trdx1]->charge() < 0. && theTrackRefs[trdx2]->charge() > 0.) {
131  negativeTrackRef = theTrackRefs[trdx1];
132  positiveTrackRef = theTrackRefs[trdx2];
133  negTransTkPtr = &theTransTracks[trdx1];
134  posTransTkPtr = &theTransTracks[trdx2];
135  } else if (theTrackRefs[trdx1]->charge() > 0. && theTrackRefs[trdx2]->charge() < 0.) {
136  negativeTrackRef = theTrackRefs[trdx2];
137  positiveTrackRef = theTrackRefs[trdx1];
138  negTransTkPtr = &theTransTracks[trdx2];
139  posTransTkPtr = &theTransTracks[trdx1];
140  } else {
141  continue;
142  }
143 
144  // measure distance between tracks at their closest approach
145  if (!posTransTkPtr->impactPointTSCP().isValid() || !negTransTkPtr->impactPointTSCP().isValid()) continue;
146  FreeTrajectoryState const & posState = posTransTkPtr->impactPointTSCP().theState();
147  FreeTrajectoryState const & negState = negTransTkPtr->impactPointTSCP().theState();
149  cApp.calculate(posState, negState);
150  if (!cApp.status()) continue;
151  float dca = std::abs(cApp.distance());
152  if (dca > tkDCACut_) continue;
153 
154  // the POCA should at least be in the sensitive volume
155  GlobalPoint cxPt = cApp.crossingPoint();
156  if (sqrt(cxPt.x()*cxPt.x() + cxPt.y()*cxPt.y()) > 120. || std::abs(cxPt.z()) > 300.) continue;
157 
158  // the tracks should at least point in the same quadrant
159  TrajectoryStateClosestToPoint const & posTSCP = posTransTkPtr->trajectoryStateClosestToPoint(cxPt);
160  TrajectoryStateClosestToPoint const & negTSCP = negTransTkPtr->trajectoryStateClosestToPoint(cxPt);
161  if (!posTSCP.isValid() || !negTSCP.isValid()) continue;
162  if (posTSCP.momentum().dot(negTSCP.momentum()) < 0) continue;
163 
164  // calculate mPiPi
165  double totalE = sqrt(posTSCP.momentum().mag2() + piMassSquared) + sqrt(negTSCP.momentum().mag2() + piMassSquared);
166  double totalESq = totalE*totalE;
167  double totalPSq = (posTSCP.momentum() + negTSCP.momentum()).mag2();
168  double mass = sqrt(totalESq - totalPSq);
169  if (mass > mPiPiCut_) continue;
170 
171  // Fill the vector of TransientTracks to send to KVF
172  std::vector<TransientTrack> transTracks;
173  transTracks.reserve(2);
174  transTracks.push_back(*posTransTkPtr);
175  transTracks.push_back(*negTransTkPtr);
176 
177  // Create the vertex fitter object and vertex the tracks
178  TransientVertex theRecoVertex;
179  if (vertexFitter_) {
180  KalmanVertexFitter theKalmanFitter(useRefTracks_ == 0 ? false : true);
181  theRecoVertex = theKalmanFitter.vertex(transTracks);
182  } else if (!vertexFitter_) {
183  useRefTracks_ = false;
184  AdaptiveVertexFitter theAdaptiveFitter;
185  theRecoVertex = theAdaptiveFitter.vertex(transTracks);
186  }
187  if (!theRecoVertex.isValid()) continue;
188 
189  reco::Vertex theVtx = theRecoVertex;
190  if (theVtx.normalizedChi2() > vtxChi2Cut_) continue;
191 
192  GlobalPoint vtxPos(theVtx.x(), theVtx.y(), theVtx.z());
193  SMatrixSym3D totalCov = theBeamSpot->rotatedCovariance3D() + theVtx.covariance();
194  SVector3 distanceVector(vtxPos.x() - beamSpotPos.x(), vtxPos.y() - beamSpotPos.y(), 0.);
195  double rVtxMag = ROOT::Math::Mag(distanceVector);
196  double sigmaRvtxMag = sqrt(ROOT::Math::Similarity(totalCov, distanceVector)) / rVtxMag;
197  if (rVtxMag/sigmaRvtxMag < vtxDecayRSigCut_) continue;
198 
199  if (innerHitPosCut_ > 0. && positiveTrackRef->innerOk()) {
200  reco::Vertex::Point posTkHitPos = positiveTrackRef->innerPosition();
201  double posTkHitPosD2 = (posTkHitPos.x()-beamSpotPos.x())*(posTkHitPos.x()-beamSpotPos.x()) +
202  (posTkHitPos.y()-beamSpotPos.y())*(posTkHitPos.y()-beamSpotPos.y());
203  if (sqrt(posTkHitPosD2) < (rVtxMag - sigmaRvtxMag*innerHitPosCut_)) continue;
204  }
205  if (innerHitPosCut_ > 0. && negativeTrackRef->innerOk()) {
206  reco::Vertex::Point negTkHitPos = negativeTrackRef->innerPosition();
207  double negTkHitPosD2 = (negTkHitPos.x()-beamSpotPos.x())*(negTkHitPos.x()-beamSpotPos.x()) +
208  (negTkHitPos.y()-beamSpotPos.y())*(negTkHitPos.y()-beamSpotPos.y());
209  if (sqrt(negTkHitPosD2) < (rVtxMag - sigmaRvtxMag*innerHitPosCut_)) continue;
210  }
211 
212  std::auto_ptr<TrajectoryStateClosestToPoint> trajPlus;
213  std::auto_ptr<TrajectoryStateClosestToPoint> trajMins;
214  std::vector<TransientTrack> theRefTracks;
215  if (theRecoVertex.hasRefittedTracks()) {
216  theRefTracks = theRecoVertex.refittedTracks();
217  }
218 
219  if (useRefTracks_ && theRefTracks.size() > 1) {
220  TransientTrack* thePositiveRefTrack = 0;
221  TransientTrack* theNegativeRefTrack = 0;
222  for (std::vector<TransientTrack>::iterator iTrack = theRefTracks.begin(); iTrack != theRefTracks.end(); ++iTrack) {
223  if (iTrack->track().charge() > 0.) {
224  thePositiveRefTrack = &*iTrack;
225  } else if (iTrack->track().charge() < 0.) {
226  theNegativeRefTrack = &*iTrack;
227  }
228  }
229  if (thePositiveRefTrack == 0 || theNegativeRefTrack == 0) continue;
230  trajPlus.reset(new TrajectoryStateClosestToPoint(thePositiveRefTrack->trajectoryStateClosestToPoint(vtxPos)));
231  trajMins.reset(new TrajectoryStateClosestToPoint(theNegativeRefTrack->trajectoryStateClosestToPoint(vtxPos)));
232  } else {
233  trajPlus.reset(new TrajectoryStateClosestToPoint(posTransTkPtr->trajectoryStateClosestToPoint(vtxPos)));
234  trajMins.reset(new TrajectoryStateClosestToPoint(negTransTkPtr->trajectoryStateClosestToPoint(vtxPos)));
235  }
236 
237  if (trajPlus.get() == 0 || trajMins.get() == 0 || !trajPlus->isValid() || !trajMins->isValid()) continue;
238 
239  GlobalVector positiveP(trajPlus->momentum());
240  GlobalVector negativeP(trajMins->momentum());
241  GlobalVector totalP(positiveP + negativeP);
242 
243  // calculate the pointing angle
244  double posx = theVtx.x() - beamSpotPos.x();
245  double posy = theVtx.y() - beamSpotPos.y();
246  double momx = totalP.x();
247  double momy = totalP.y();
248  double pointangle = (posx*momx+posy*momy)/(sqrt(posx*posx+posy*posy)*sqrt(momx*momx+momy*momy));
249  if (pointangle < v0CosThetaCut_) continue;
250 
251  // calculate total energy of V0 3 ways: assume it's a kShort, a Lambda, or a LambdaBar.
252  double piPlusE = sqrt(positiveP.mag2() + piMassSquared);
253  double piMinusE = sqrt(negativeP.mag2() + piMassSquared);
254  double protonE = sqrt(positiveP.mag2() + protonMassSquared);
255  double antiProtonE = sqrt(negativeP.mag2() + protonMassSquared);
256  double kShortETot = piPlusE + piMinusE;
257  double lambdaEtot = protonE + piMinusE;
258  double lambdaBarEtot = antiProtonE + piPlusE;
259 
260  // Create momentum 4-vectors for the 3 candidate types
261  const Particle::LorentzVector kShortP4(totalP.x(), totalP.y(), totalP.z(), kShortETot);
262  const Particle::LorentzVector lambdaP4(totalP.x(), totalP.y(), totalP.z(), lambdaEtot);
263  const Particle::LorentzVector lambdaBarP4(totalP.x(), totalP.y(), totalP.z(), lambdaBarEtot);
264 
265  Particle::Point vtx(theVtx.x(), theVtx.y(), theVtx.z());
266  const Vertex::CovarianceMatrix vtxCov(theVtx.covariance());
267  double vtxChi2(theVtx.chi2());
268  double vtxNdof(theVtx.ndof());
269 
270  // Create the VertexCompositeCandidate object that will be stored in the Event
271  VertexCompositeCandidate* theKshort = nullptr;
272  VertexCompositeCandidate* theLambda = nullptr;
273  VertexCompositeCandidate* theLambdaBar = nullptr;
274 
275  if (doKShorts_) {
276  theKshort = new VertexCompositeCandidate(0, kShortP4, vtx, vtxCov, vtxChi2, vtxNdof);
277  }
278  if (doLambdas_) {
279  if (positiveP.mag2() > negativeP.mag2()) {
280  theLambda = new VertexCompositeCandidate(0, lambdaP4, vtx, vtxCov, vtxChi2, vtxNdof);
281  } else {
282  theLambdaBar = new VertexCompositeCandidate(0, lambdaBarP4, vtx, vtxCov, vtxChi2, vtxNdof);
283  }
284  }
285 
286  // Create daughter candidates for the VertexCompositeCandidates
287  RecoChargedCandidate thePiPlusCand(
288  1, Particle::LorentzVector(positiveP.x(), positiveP.y(), positiveP.z(), piPlusE), vtx);
289  thePiPlusCand.setTrack(positiveTrackRef);
290 
291  RecoChargedCandidate thePiMinusCand(
292  -1, Particle::LorentzVector(negativeP.x(), negativeP.y(), negativeP.z(), piMinusE), vtx);
293  thePiMinusCand.setTrack(negativeTrackRef);
294 
295  RecoChargedCandidate theProtonCand(
296  1, Particle::LorentzVector(positiveP.x(), positiveP.y(), positiveP.z(), protonE), vtx);
297  theProtonCand.setTrack(positiveTrackRef);
298 
299  RecoChargedCandidate theAntiProtonCand(
300  -1, Particle::LorentzVector(negativeP.x(), negativeP.y(), negativeP.z(), antiProtonE), vtx);
301  theAntiProtonCand.setTrack(negativeTrackRef);
302 
303  AddFourMomenta addp4;
304  // Store the daughter Candidates in the VertexCompositeCandidates if they pass mass cuts
305  if (doKShorts_) {
306  theKshort->addDaughter(thePiPlusCand);
307  theKshort->addDaughter(thePiMinusCand);
308  theKshort->setPdgId(310);
309  addp4.set(*theKshort);
310  if (theKshort->mass() < kShortMass + kShortMassCut_ && theKshort->mass() > kShortMass - kShortMassCut_) {
311  theKshorts.push_back(std::move(*theKshort));
312  }
313  }
314  if (doLambdas_ && theLambda) {
315  theLambda->addDaughter(theProtonCand);
316  theLambda->addDaughter(thePiMinusCand);
317  theLambda->setPdgId(3122);
318  addp4.set( *theLambda );
319  if (theLambda->mass() < lambdaMass + lambdaMassCut_ && theLambda->mass() > lambdaMass - lambdaMassCut_) {
320  theLambdas.push_back(std::move(*theLambda));
321  }
322  } else if (doLambdas_ && theLambdaBar) {
323  theLambdaBar->addDaughter(theAntiProtonCand);
324  theLambdaBar->addDaughter(thePiPlusCand);
325  theLambdaBar->setPdgId(-3122);
326  addp4.set(*theLambdaBar);
327  if (theLambdaBar->mass() < lambdaMass + lambdaMassCut_ && theLambdaBar->mass() > lambdaMass - lambdaMassCut_) {
328  theLambdas.push_back(std::move(*theLambdaBar));
329  }
330  }
331 
332  delete theKshort;
333  delete theLambda;
334  delete theLambdaBar;
335  theKshort = theLambda = theLambdaBar = nullptr;
336 
337  }
338  }
339 }
bool doLambdas_
Definition: V0Fitter.h:70
double tkChi2Cut_
Definition: V0Fitter.h:73
T mag2() const
Definition: PV3DBase.h:66
ROOT::Math::SMatrix< double, 3, 3, ROOT::Math::MatRepSym< double, 3 > > SMatrixSym3D
Definition: V0Fitter.cc:46
double lambdaMassCut_
Definition: V0Fitter.h:87
double tkIPSigCut_
Definition: V0Fitter.h:76
bool vertexFitter_
Definition: V0Fitter.h:67
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:609
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:446
TrajectoryStateClosestToPoint impactPointTSCP() const
const FreeTrajectoryState & theState() const
double dxyError() const
error on dxy
Definition: TrackBase.h:844
double y() const
y coordinate
Definition: Vertex.h:110
double v0CosThetaCut_
Definition: V0Fitter.h:84
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:13
T y() const
Definition: PV3DBase.h:63
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:123
bool hasRefittedTracks() const
std::pair< double, double > Point
Definition: CaloEllipse.h:18
math::Error< dimension >::type CovarianceMatrix
covariance error matrix (3x3)
Definition: Vertex.h:45
virtual GlobalPoint crossingPoint() const
double tkDCACut_
Definition: V0Fitter.h:81
int tkNHitsCut_
Definition: V0Fitter.h:74
bool useRefTracks_
Definition: V0Fitter.h:68
edm::EDGetTokenT< reco::TrackCollection > token_tracks
Definition: V0Fitter.h:89
T sqrt(T t)
Definition: SSEVec.h:48
double pt() const
track transverse momentum
Definition: TrackBase.h:669
T z() const
Definition: PV3DBase.h:64
T mag2() const
The vector magnitude squared. Equivalent to vec.dot(vec)
virtual CachingVertex< 5 > vertex(const std::vector< reco::TransientTrack > &) const
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double chi2() const
chi-squares
Definition: Vertex.h:95
double z() const
y coordinate
Definition: Vertex.h:112
math::XYZPoint Point
point in the space
Definition: Vertex.h:39
unsigned short numberOfValidHits() const
number of valid hits found
Definition: TrackBase.h:868
double mPiPiCut_
Definition: V0Fitter.h:82
double kShortMassCut_
Definition: V0Fitter.h:86
void addDaughter(const Candidate &, const std::string &s="")
add a clone of the passed candidate as daughter
double ndof() const
Definition: Vertex.h:102
double x() const
x coordinate
Definition: Vertex.h:108
double innerHitPosCut_
Definition: V0Fitter.h:83
virtual bool calculate(const TrajectoryStateOnSurface &sta, const TrajectoryStateOnSurface &stb)
TrajectoryStateClosestToPoint trajectoryStateClosestToPoint(const GlobalPoint &point) const
const T & get() const
Definition: EventSetup.h:55
double vtxDecayRSigCut_
Definition: V0Fitter.h:79
double tkPtCut_
Definition: V0Fitter.h:75
double vtxChi2Cut_
Definition: V0Fitter.h:78
virtual void setPdgId(int pdgId)
virtual float mass() const
mass
edm::EDGetTokenT< reco::BeamSpot > token_beamSpot
Definition: V0Fitter.h:90
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:104
ROOT::Math::SVector< double, 3 > SVector3
Definition: V0Fitter.cc:47
virtual float distance() 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:639
T x() const
Definition: PV3DBase.h:62
bool isValid() const
virtual bool status() const
math::PtEtaPhiELorentzVectorF LorentzVector
bool doKShorts_
Definition: V0Fitter.h:69

Member Data Documentation

bool V0Fitter::doKShorts_
private

Definition at line 69 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

bool V0Fitter::doLambdas_
private

Definition at line 70 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

double V0Fitter::innerHitPosCut_
private

Definition at line 83 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

double V0Fitter::kShortMassCut_
private

Definition at line 86 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

double V0Fitter::lambdaMassCut_
private

Definition at line 87 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

double V0Fitter::mPiPiCut_
private

Definition at line 82 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

double V0Fitter::tkChi2Cut_
private

Definition at line 73 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

double V0Fitter::tkDCACut_
private

Definition at line 81 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

double V0Fitter::tkIPSigCut_
private

Definition at line 76 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

int V0Fitter::tkNHitsCut_
private

Definition at line 74 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

double V0Fitter::tkPtCut_
private

Definition at line 75 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

edm::EDGetTokenT<reco::BeamSpot> V0Fitter::token_beamSpot
private

Definition at line 90 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

edm::EDGetTokenT<reco::TrackCollection> V0Fitter::token_tracks
private

Definition at line 89 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

bool V0Fitter::useRefTracks_
private

Definition at line 68 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

double V0Fitter::v0CosThetaCut_
private

Definition at line 84 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

bool V0Fitter::vertexFitter_
private

Definition at line 67 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

double V0Fitter::vtxChi2Cut_
private

Definition at line 78 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

double V0Fitter::vtxDecayRSigCut_
private

Definition at line 79 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().