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

double cosThetaXYCut_
 
double cosThetaXYZCut_
 
bool doKShorts_
 
bool doLambdas_
 
double innerHitPosCut_
 
double kShortMassCut_
 
double lambdaMassCut_
 
double mPiPiCut_
 
double tkChi2Cut_
 
double tkDCACut_
 
double tkIPSigXYCut_
 
double tkIPSigZCut_
 
int tkNHitsCut_
 
double tkPtCut_
 
edm::EDGetTokenT< reco::BeamSpottoken_beamSpot
 
edm::EDGetTokenT
< reco::TrackCollection
token_tracks
 
edm::EDGetTokenT< std::vector
< reco::Vertex > > 
token_vertices
 
bool useRefTracks_
 
bool useVertex_
 
bool vertexFitter_
 
double vtxChi2Cut_
 
double vtxDecaySigXYCut_
 
double vtxDecaySigXYZCut_
 

Detailed Description

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

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

Definition at line 46 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 cosThetaXYCut_, cosThetaXYZCut_, doKShorts_, doLambdas_, edm::ParameterSet::getParameter(), innerHitPosCut_, kShortMassCut_, lambdaMassCut_, mPiPiCut_, tkChi2Cut_, tkDCACut_, tkIPSigXYCut_, tkIPSigZCut_, tkNHitsCut_, tkPtCut_, token_beamSpot, token_tracks, token_vertices, useRefTracks_, useVertex_, vertexFitter_, vtxChi2Cut_, vtxDecaySigXYCut_, and vtxDecaySigXYZCut_.

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 }
bool doLambdas_
Definition: V0Fitter.h:58
double tkChi2Cut_
Definition: V0Fitter.h:61
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
double cosThetaXYZCut_
Definition: V0Fitter.h:75
double lambdaMassCut_
Definition: V0Fitter.h:78
double tkIPSigXYCut_
Definition: V0Fitter.h:64
bool vertexFitter_
Definition: V0Fitter.h:55
edm::EDGetTokenT< std::vector< reco::Vertex > > token_vertices
Definition: V0Fitter.h:83
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:13
double vtxDecaySigXYZCut_
Definition: V0Fitter.h:69
double tkDCACut_
Definition: V0Fitter.h:71
int tkNHitsCut_
Definition: V0Fitter.h:62
bool useRefTracks_
Definition: V0Fitter.h:56
edm::EDGetTokenT< reco::TrackCollection > token_tracks
Definition: V0Fitter.h:80
bool useVertex_
Definition: V0Fitter.h:82
double mPiPiCut_
Definition: V0Fitter.h:72
double kShortMassCut_
Definition: V0Fitter.h:77
double tkIPSigZCut_
Definition: V0Fitter.h:65
double cosThetaXYCut_
Definition: V0Fitter.h:74
double innerHitPosCut_
Definition: V0Fitter.h:73
double vtxDecaySigXYCut_
Definition: V0Fitter.h:68
double tkPtCut_
Definition: V0Fitter.h:63
double vtxChi2Cut_
Definition: V0Fitter.h:67
edm::EDGetTokenT< reco::BeamSpot > token_beamSpot
Definition: V0Fitter.h:81
bool doKShorts_
Definition: V0Fitter.h:57

Member Function Documentation

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

Definition at line 87 of file V0Fitter.cc.

References funct::abs(), reco::CompositeCandidate::addDaughter(), ClosestApproachInRPhi::calculate(), RecoTauCleanerPlugins::charge, reco::Vertex::chi2(), cosThetaXYCut_, cosThetaXYZCut_, reco::Vertex::covariance(), ClosestApproachInRPhi::crossingPoint(), ClosestApproachInRPhi::distance(), HLT_25ns14e33_v1_cff::distance, doKShorts_, doLambdas_, Vector3DBase< T, FrameTag >::dot(), reco::TrackBase::dxy(), reco::TrackBase::dxyError(), reco::TrackBase::dz(), reco::TrackBase::dzError(), 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::BeamSpot::position(), edm::Handle< T >::product(), edm::ESHandle< class >::product(), reco::TrackBase::pt(), TransientVertex::refittedTracks(), reco::BeamSpot::rotatedCovariance3D(), AddFourMomenta::set(), reco::LeafCandidate::setPdgId(), reco::RecoChargedCandidate::setTrack(), mathSSE::sqrt(), ClosestApproachInRPhi::status(), TrajectoryStateClosestToPoint::theState(), tkChi2Cut_, tkDCACut_, tkIPSigXYCut_, tkIPSigZCut_, tkNHitsCut_, tkPtCut_, token_beamSpot, token_tracks, token_vertices, reco::TransientTrack::trajectoryStateClosestToPoint(), useRefTracks_, useVertex_, KalmanVertexFitter::vertex(), AdaptiveVertexFitter::vertex(), vertexFitter_, HLT_25ns14e33_v1_cff::vertices, vtxChi2Cut_, vtxDecaySigXYCut_, vtxDecaySigXYZCut_, 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().

89 {
90  using std::vector;
91 
93  iEvent.getByToken(token_tracks, theTrackHandle);
94  if (!theTrackHandle->size()) 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  if (!posTransTkPtr->impactPointTSCP().isValid() || !negTransTkPtr->impactPointTSCP().isValid()) continue;
158  FreeTrajectoryState const & posState = posTransTkPtr->impactPointTSCP().theState();
159  FreeTrajectoryState const & negState = negTransTkPtr->impactPointTSCP().theState();
161  cApp.calculate(posState, negState);
162  if (!cApp.status()) continue;
163  float dca = std::abs(cApp.distance());
164  if (dca > tkDCACut_) continue;
165 
166  // the POCA should at least be in the sensitive volume
167  GlobalPoint cxPt = cApp.crossingPoint();
168  if (sqrt(cxPt.x()*cxPt.x() + cxPt.y()*cxPt.y()) > 120. || std::abs(cxPt.z()) > 300.) continue;
169 
170  // the tracks should at least point in the same quadrant
171  TrajectoryStateClosestToPoint const & posTSCP = posTransTkPtr->trajectoryStateClosestToPoint(cxPt);
172  TrajectoryStateClosestToPoint const & negTSCP = negTransTkPtr->trajectoryStateClosestToPoint(cxPt);
173  if (!posTSCP.isValid() || !negTSCP.isValid()) continue;
174  if (posTSCP.momentum().dot(negTSCP.momentum()) < 0) continue;
175 
176  // calculate mPiPi
177  double totalE = sqrt(posTSCP.momentum().mag2() + piMassSquared) + sqrt(negTSCP.momentum().mag2() + piMassSquared);
178  double totalESq = totalE*totalE;
179  double totalPSq = (posTSCP.momentum() + negTSCP.momentum()).mag2();
180  double mass = sqrt(totalESq - totalPSq);
181  if (mass > mPiPiCut_) continue;
182 
183  // Fill the vector of TransientTracks to send to KVF
184  std::vector<reco::TransientTrack> transTracks;
185  transTracks.reserve(2);
186  transTracks.push_back(*posTransTkPtr);
187  transTracks.push_back(*negTransTkPtr);
188 
189  // create the vertex fitter object and vertex the tracks
190  TransientVertex theRecoVertex;
191  if (vertexFitter_) {
192  KalmanVertexFitter theKalmanFitter(useRefTracks_ == 0 ? false : true);
193  theRecoVertex = theKalmanFitter.vertex(transTracks);
194  } else if (!vertexFitter_) {
195  useRefTracks_ = false;
196  AdaptiveVertexFitter theAdaptiveFitter;
197  theRecoVertex = theAdaptiveFitter.vertex(transTracks);
198  }
199  if (!theRecoVertex.isValid()) continue;
200 
201  reco::Vertex theVtx = theRecoVertex;
202  if (theVtx.normalizedChi2() > vtxChi2Cut_) continue;
203  GlobalPoint vtxPos(theVtx.x(), theVtx.y(), theVtx.z());
204 
205  // 2D decay significance
206  SMatrixSym3D totalCov = theBeamSpot->rotatedCovariance3D() + theVtx.covariance();
207  if (useVertex_) totalCov = referenceVtx.covariance() + theVtx.covariance();
208  SVector3 distVecXY(vtxPos.x()-referencePos.x(), vtxPos.y()-referencePos.y(), 0.);
209  double distMagXY = ROOT::Math::Mag(distVecXY);
210  double sigmaDistMagXY = sqrt(ROOT::Math::Similarity(totalCov, distVecXY)) / distMagXY;
211  if (distMagXY/sigmaDistMagXY < vtxDecaySigXYCut_) continue;
212 
213  // 3D decay significance
214  SVector3 distVecXYZ(vtxPos.x()-referencePos.x(), vtxPos.y()-referencePos.y(), vtxPos.z()-referencePos.z());
215  double distMagXYZ = ROOT::Math::Mag(distVecXYZ);
216  double sigmaDistMagXYZ = sqrt(ROOT::Math::Similarity(totalCov, distVecXYZ)) / distMagXYZ;
217  if (distMagXYZ/sigmaDistMagXYZ < vtxDecaySigXYZCut_) continue;
218 
219  // make sure the vertex radius is within the inner track hit radius
220  if (innerHitPosCut_ > 0. && positiveTrackRef->innerOk()) {
221  reco::Vertex::Point posTkHitPos = positiveTrackRef->innerPosition();
222  double posTkHitPosD2 = (posTkHitPos.x()-referencePos.x())*(posTkHitPos.x()-referencePos.x()) +
223  (posTkHitPos.y()-referencePos.y())*(posTkHitPos.y()-referencePos.y());
224  if (sqrt(posTkHitPosD2) < (distMagXY - sigmaDistMagXY*innerHitPosCut_)) continue;
225  }
226  if (innerHitPosCut_ > 0. && negativeTrackRef->innerOk()) {
227  reco::Vertex::Point negTkHitPos = negativeTrackRef->innerPosition();
228  double negTkHitPosD2 = (negTkHitPos.x()-referencePos.x())*(negTkHitPos.x()-referencePos.x()) +
229  (negTkHitPos.y()-referencePos.y())*(negTkHitPos.y()-referencePos.y());
230  if (sqrt(negTkHitPosD2) < (distMagXY - sigmaDistMagXY*innerHitPosCut_)) continue;
231  }
232 
233  std::auto_ptr<TrajectoryStateClosestToPoint> trajPlus;
234  std::auto_ptr<TrajectoryStateClosestToPoint> trajMins;
235  std::vector<reco::TransientTrack> theRefTracks;
236  if (theRecoVertex.hasRefittedTracks()) {
237  theRefTracks = theRecoVertex.refittedTracks();
238  }
239 
240  if (useRefTracks_ && theRefTracks.size() > 1) {
241  reco::TransientTrack* thePositiveRefTrack = 0;
242  reco::TransientTrack* theNegativeRefTrack = 0;
243  for (std::vector<reco::TransientTrack>::iterator iTrack = theRefTracks.begin(); iTrack != theRefTracks.end(); ++iTrack) {
244  if (iTrack->track().charge() > 0.) {
245  thePositiveRefTrack = &*iTrack;
246  } else if (iTrack->track().charge() < 0.) {
247  theNegativeRefTrack = &*iTrack;
248  }
249  }
250  if (thePositiveRefTrack == 0 || theNegativeRefTrack == 0) continue;
251  trajPlus.reset(new TrajectoryStateClosestToPoint(thePositiveRefTrack->trajectoryStateClosestToPoint(vtxPos)));
252  trajMins.reset(new TrajectoryStateClosestToPoint(theNegativeRefTrack->trajectoryStateClosestToPoint(vtxPos)));
253  } else {
254  trajPlus.reset(new TrajectoryStateClosestToPoint(posTransTkPtr->trajectoryStateClosestToPoint(vtxPos)));
255  trajMins.reset(new TrajectoryStateClosestToPoint(negTransTkPtr->trajectoryStateClosestToPoint(vtxPos)));
256  }
257 
258  if (trajPlus.get() == 0 || trajMins.get() == 0 || !trajPlus->isValid() || !trajMins->isValid()) continue;
259 
260  GlobalVector positiveP(trajPlus->momentum());
261  GlobalVector negativeP(trajMins->momentum());
262  GlobalVector totalP(positiveP + negativeP);
263 
264  // 2D pointing angle
265  double dx = theVtx.x()-referencePos.x();
266  double dy = theVtx.y()-referencePos.y();
267  double px = totalP.x();
268  double py = totalP.y();
269  double angleXY = (dx*px+dy*py)/(sqrt(dx*dx+dy*dy)*sqrt(px*px+py*py));
270  if (angleXY < cosThetaXYCut_) continue;
271 
272  // 3D pointing angle
273  double dz = theVtx.z()-referencePos.z();
274  double pz = totalP.z();
275  double angleXYZ = (dx*px+dy*py+dz*pz)/(sqrt(dx*dx+dy*dy+dz*dz)*sqrt(px*px+py*py+pz*pz));
276  if (angleXYZ < cosThetaXYZCut_) continue;
277 
278  // calculate total energy of V0 3 ways: assume it's a kShort, a Lambda, or a LambdaBar.
279  double piPlusE = sqrt(positiveP.mag2() + piMassSquared);
280  double piMinusE = sqrt(negativeP.mag2() + piMassSquared);
281  double protonE = sqrt(positiveP.mag2() + protonMassSquared);
282  double antiProtonE = sqrt(negativeP.mag2() + protonMassSquared);
283  double kShortETot = piPlusE + piMinusE;
284  double lambdaEtot = protonE + piMinusE;
285  double lambdaBarEtot = antiProtonE + piPlusE;
286 
287  // Create momentum 4-vectors for the 3 candidate types
288  const reco::Particle::LorentzVector kShortP4(totalP.x(), totalP.y(), totalP.z(), kShortETot);
289  const reco::Particle::LorentzVector lambdaP4(totalP.x(), totalP.y(), totalP.z(), lambdaEtot);
290  const reco::Particle::LorentzVector lambdaBarP4(totalP.x(), totalP.y(), totalP.z(), lambdaBarEtot);
291 
292  reco::Particle::Point vtx(theVtx.x(), theVtx.y(), theVtx.z());
293  const reco::Vertex::CovarianceMatrix vtxCov(theVtx.covariance());
294  double vtxChi2(theVtx.chi2());
295  double vtxNdof(theVtx.ndof());
296 
297  // Create the VertexCompositeCandidate object that will be stored in the Event
298  reco::VertexCompositeCandidate* theKshort = nullptr;
299  reco::VertexCompositeCandidate* theLambda = nullptr;
300  reco::VertexCompositeCandidate* theLambdaBar = nullptr;
301 
302  if (doKShorts_) {
303  theKshort = new reco::VertexCompositeCandidate(0, kShortP4, vtx, vtxCov, vtxChi2, vtxNdof);
304  }
305  if (doLambdas_) {
306  if (positiveP.mag2() > negativeP.mag2()) {
307  theLambda = new reco::VertexCompositeCandidate(0, lambdaP4, vtx, vtxCov, vtxChi2, vtxNdof);
308  } else {
309  theLambdaBar = new reco::VertexCompositeCandidate(0, lambdaBarP4, vtx, vtxCov, vtxChi2, vtxNdof);
310  }
311  }
312 
313  // Create daughter candidates for the VertexCompositeCandidates
314  reco::RecoChargedCandidate thePiPlusCand(
315  1, reco::Particle::LorentzVector(positiveP.x(), positiveP.y(), positiveP.z(), piPlusE), vtx);
316  thePiPlusCand.setTrack(positiveTrackRef);
317 
318  reco::RecoChargedCandidate thePiMinusCand(
319  -1, reco::Particle::LorentzVector(negativeP.x(), negativeP.y(), negativeP.z(), piMinusE), vtx);
320  thePiMinusCand.setTrack(negativeTrackRef);
321 
322  reco::RecoChargedCandidate theProtonCand(
323  1, reco::Particle::LorentzVector(positiveP.x(), positiveP.y(), positiveP.z(), protonE), vtx);
324  theProtonCand.setTrack(positiveTrackRef);
325 
326  reco::RecoChargedCandidate theAntiProtonCand(
327  -1, reco::Particle::LorentzVector(negativeP.x(), negativeP.y(), negativeP.z(), antiProtonE), vtx);
328  theAntiProtonCand.setTrack(negativeTrackRef);
329 
330  AddFourMomenta addp4;
331  // Store the daughter Candidates in the VertexCompositeCandidates if they pass mass cuts
332  if (doKShorts_) {
333  theKshort->addDaughter(thePiPlusCand);
334  theKshort->addDaughter(thePiMinusCand);
335  theKshort->setPdgId(310);
336  addp4.set(*theKshort);
337  if (theKshort->mass() < kShortMass + kShortMassCut_ && theKshort->mass() > kShortMass - kShortMassCut_) {
338  theKshorts.push_back(std::move(*theKshort));
339  }
340  }
341  if (doLambdas_ && theLambda) {
342  theLambda->addDaughter(theProtonCand);
343  theLambda->addDaughter(thePiMinusCand);
344  theLambda->setPdgId(3122);
345  addp4.set( *theLambda );
346  if (theLambda->mass() < lambdaMass + lambdaMassCut_ && theLambda->mass() > lambdaMass - lambdaMassCut_) {
347  theLambdas.push_back(std::move(*theLambda));
348  }
349  } else if (doLambdas_ && theLambdaBar) {
350  theLambdaBar->addDaughter(theAntiProtonCand);
351  theLambdaBar->addDaughter(thePiPlusCand);
352  theLambdaBar->setPdgId(-3122);
353  addp4.set(*theLambdaBar);
354  if (theLambdaBar->mass() < lambdaMass + lambdaMassCut_ && theLambdaBar->mass() > lambdaMass - lambdaMassCut_) {
355  theLambdas.push_back(std::move(*theLambdaBar));
356  }
357  }
358 
359  delete theKshort;
360  delete theLambda;
361  delete theLambdaBar;
362  theKshort = theLambda = theLambdaBar = nullptr;
363 
364  }
365  }
366 }
bool doLambdas_
Definition: V0Fitter.h:58
double tkChi2Cut_
Definition: V0Fitter.h:61
T mag2() const
Definition: PV3DBase.h:66
ROOT::Math::SMatrix< double, 3, 3, ROOT::Math::MatRepSym< double, 3 > > SMatrixSym3D
Definition: V0Fitter.cc:46
double cosThetaXYZCut_
Definition: V0Fitter.h:75
double lambdaMassCut_
Definition: V0Fitter.h:78
double tkIPSigXYCut_
Definition: V0Fitter.h:64
bool vertexFitter_
Definition: V0Fitter.h:55
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:524
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
TrajectoryStateClosestToPoint impactPointTSCP() const
const FreeTrajectoryState & theState() const
double dxyError() const
error on dxy
Definition: TrackBase.h:759
edm::EDGetTokenT< std::vector< reco::Vertex > > token_vertices
Definition: V0Fitter.h:83
double y() const
y coordinate
Definition: Vertex.h:110
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
double vtxDecaySigXYZCut_
Definition: V0Fitter.h:69
math::Error< dimension >::type CovarianceMatrix
covariance error matrix (3x3)
Definition: Vertex.h:45
virtual GlobalPoint crossingPoint() const
double tkDCACut_
Definition: V0Fitter.h:71
int tkNHitsCut_
Definition: V0Fitter.h:62
virtual double mass() const
mass
bool useRefTracks_
Definition: V0Fitter.h:56
edm::EDGetTokenT< reco::TrackCollection > token_tracks
Definition: V0Fitter.h:80
T sqrt(T t)
Definition: SSEVec.h:48
double pt() const
track transverse momentum
Definition: TrackBase.h:584
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
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:95
double z() const
y coordinate
Definition: Vertex.h:112
bool useVertex_
Definition: V0Fitter.h:82
math::XYZPoint Point
point in the space
Definition: Vertex.h:39
unsigned short numberOfValidHits() const
number of valid hits found
Definition: TrackBase.h:783
double mPiPiCut_
Definition: V0Fitter.h:72
double kShortMassCut_
Definition: V0Fitter.h:77
double tkIPSigZCut_
Definition: V0Fitter.h:65
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 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:572
double dzError() const
error on dz
Definition: TrackBase.h:777
double cosThetaXYCut_
Definition: V0Fitter.h:74
double x() const
x coordinate
Definition: Vertex.h:108
double innerHitPosCut_
Definition: V0Fitter.h:73
double vtxDecaySigXYCut_
Definition: V0Fitter.h:68
virtual bool calculate(const TrajectoryStateOnSurface &sta, const TrajectoryStateOnSurface &stb)
T const * product() const
Definition: Handle.h:81
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
TrajectoryStateClosestToPoint trajectoryStateClosestToPoint(const GlobalPoint &point) const
const T & get() const
Definition: EventSetup.h:55
T const * product() const
Definition: ESHandle.h:86
double tkPtCut_
Definition: V0Fitter.h:63
double vtxChi2Cut_
Definition: V0Fitter.h:67
virtual void setPdgId(int pdgId)
const Point & position() const
position
Definition: BeamSpot.h:62
edm::EDGetTokenT< reco::BeamSpot > token_beamSpot
Definition: V0Fitter.h:81
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
Covariance3DMatrix rotatedCovariance3D() const
Definition: BeamSpot.cc:78
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:554
T x() const
Definition: PV3DBase.h:62
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Particle.h:21
bool isValid() const
virtual bool status() const
bool doKShorts_
Definition: V0Fitter.h:57

Member Data Documentation

double V0Fitter::cosThetaXYCut_
private

Definition at line 74 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

double V0Fitter::cosThetaXYZCut_
private

Definition at line 75 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

bool V0Fitter::doKShorts_
private

Definition at line 57 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

bool V0Fitter::doLambdas_
private

Definition at line 58 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

double V0Fitter::innerHitPosCut_
private

Definition at line 73 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

double V0Fitter::kShortMassCut_
private

Definition at line 77 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

double V0Fitter::lambdaMassCut_
private

Definition at line 78 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

double V0Fitter::mPiPiCut_
private

Definition at line 72 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

double V0Fitter::tkChi2Cut_
private

Definition at line 61 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

double V0Fitter::tkDCACut_
private

Definition at line 71 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

double V0Fitter::tkIPSigXYCut_
private

Definition at line 64 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

double V0Fitter::tkIPSigZCut_
private

Definition at line 65 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

int V0Fitter::tkNHitsCut_
private

Definition at line 62 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

double V0Fitter::tkPtCut_
private

Definition at line 63 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

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

Definition at line 81 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

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

Definition at line 80 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

edm::EDGetTokenT<std::vector<reco::Vertex> > V0Fitter::token_vertices
private

Definition at line 83 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

bool V0Fitter::useRefTracks_
private

Definition at line 56 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

bool V0Fitter::useVertex_
private

Definition at line 82 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

bool V0Fitter::vertexFitter_
private

Definition at line 55 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

double V0Fitter::vtxChi2Cut_
private

Definition at line 67 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

double V0Fitter::vtxDecaySigXYCut_
private

Definition at line 68 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().

double V0Fitter::vtxDecaySigXYZCut_
private

Definition at line 69 of file V0Fitter.h.

Referenced by fitAll(), and V0Fitter().