CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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"
21 
29 
30 #include <Math/Functions.h>
31 #include <Math/SVector.h>
32 #include <Math/SMatrix.h>
33 #include <typeinfo>
34 #include <memory>
35 
36 // Constants
37 
38 namespace {
39 const double piMass = 0.13957018;
40 const double piMassSquared = piMass*piMass;
41 const double protonMass = 0.938272013;
42 const double protonMassSquared = protonMass*protonMass;
43 const double kShortMass = 0.497614;
44 const double lambdaMass = 1.115683;
45 }
46 
47 
48 // Constructor and (empty) destructor
50  edm::ConsumesCollector && iC) {
51  using std::string;
52 
53  // Get the track reco algorithm from the ParameterSet
54  token_beamSpot = iC.consumes<reco::BeamSpot>(edm::InputTag("offlineBeamSpot"));
55  token_tracks = iC.consumes<reco::TrackCollection>(theParameters.getParameter<edm::InputTag>("trackRecoAlgorithm"));
56 
57  // ------> Initialize parameters from PSet. ALL TRACKED, so no defaults.
58  // First set bits to do various things:
59  // -decide whether to use the KVF track smoother, and whether to store those
60  // tracks in the reco::Vertex
61  useRefTrax = theParameters.getParameter<bool>(string("useSmoothing"));
62 
63  // -whether to reconstruct K0s
64  doKshorts = theParameters.getParameter<bool>(string("selectKshorts"));
65  // -whether to reconstruct Lambdas
66  doLambdas = theParameters.getParameter<bool>(string("selectLambdas"));
67 
68  // Second, initialize post-fit cuts
69  chi2Cut = theParameters.getParameter<double>(string("vtxChi2Cut"));
70  tkChi2Cut = theParameters.getParameter<double>(string("tkChi2Cut"));
71  tkNhitsCut = theParameters.getParameter<int>(string("tkNhitsCut"));
72  rVtxCut = theParameters.getParameter<double>(string("rVtxCut"));
73  vtxSigCut = theParameters.getParameter<double>(string("vtxSignificance2DCut"));
74  collinCut = theParameters.getParameter<double>(string("collinearityCut"));
75  kShortMassCut = theParameters.getParameter<double>(string("kShortMassCut"));
76  lambdaMassCut = theParameters.getParameter<double>(string("lambdaMassCut"));
77  impactParameterSigCut = theParameters.getParameter<double>(string("impactParameterSigCut"));
78  mPiPiCut = theParameters.getParameter<double>(string("mPiPiCut"));
79  tkDCACut = theParameters.getParameter<double>(string("tkDCACut"));
80  vtxFitter = theParameters.getParameter<edm::InputTag>("vertexFitter");
81  innerHitPosCut = theParameters.getParameter<double>(string("innerHitPosCut"));
82  std::vector<std::string> qual = theParameters.getParameter<std::vector<std::string> >("trackQualities");
83  for (unsigned int ndx = 0; ndx < qual.size(); ndx++) {
84  qualities.push_back(reco::TrackBase::qualityByName(qual[ndx]));
85  }
86 
87  //edm::LogInfo("V0Producer") << "Using " << vtxFitter << " to fit V0 vertices.\n";
88  //std::cout << "Using " << vtxFitter << " to fit V0 vertices." << std::endl;
89  // FOR DEBUG:
90  //initFileOutput();
91  //--------------------
92 
93  //std::cout << "Entering V0Producer" << std::endl;
94 
95 
96  // FOR DEBUG:
97  //cleanupFileOutput();
98  //--------------------
99 
100 }
101 
103 }
104 
105 // Method containing the algorithm for vertex reconstruction
109 
110 
111  using std::vector;
112  using std::cout;
113  using std::endl;
114  using namespace reco;
115  using namespace edm;
116 
117 
118  // Create std::vectors for Tracks and TrackRefs (required for
119  // passing to the KalmanVertexFitter)
120  std::vector<TrackRef> theTrackRefs;
121  std::vector<TransientTrack> theTransTracks;
122 
123  // Handles for tracks, B-field, and tracker geometry
124  Handle<reco::TrackCollection> theTrackHandle;
125  Handle<reco::BeamSpot> theBeamSpotHandle;
126  ESHandle<MagneticField> bFieldHandle;
127  ESHandle<TrackerGeometry> trackerGeomHandle;
128  ESHandle<GlobalTrackingGeometry> globTkGeomHandle;
129  //cout << "Check 0" << endl;
130 
131  // Get the tracks from the event, and get the B-field record
132  // from the EventSetup
133  iEvent.getByToken(token_tracks, theTrackHandle);
134  iEvent.getByToken(token_beamSpot,theBeamSpotHandle);
135  if( !theTrackHandle->size() ) return;
136  iSetup.get<IdealMagneticFieldRecord>().get(bFieldHandle);
137  iSetup.get<TrackerDigiGeometryRecord>().get(trackerGeomHandle);
138  iSetup.get<GlobalTrackingGeometryRecord>().get(globTkGeomHandle);
139 
140  trackerGeom = trackerGeomHandle.product();
141  magField = bFieldHandle.product();
142 
143  // Fill vectors of TransientTracks and TrackRefs after applying preselection cuts.
144  for(unsigned int indx = 0; indx < theTrackHandle->size(); indx++) {
145  TrackRef tmpRef( theTrackHandle, indx );
146  bool quality_ok = true;
147  if (qualities.size()!=0) {
148  quality_ok = false;
149  for (unsigned int ndx_ = 0; ndx_ < qualities.size(); ndx_++) {
150  if (tmpRef->quality(qualities[ndx_])){
151  quality_ok = true;
152  break;
153  }
154  }
155  }
156  if( !quality_ok ) continue;
157 
158 
159  if( tmpRef->normalizedChi2() < tkChi2Cut &&
160  tmpRef->numberOfValidHits() >= tkNhitsCut ) {
161  TransientTrack tmpTk( *tmpRef, &(*bFieldHandle), globTkGeomHandle );
162 
164  TSCBLBuilderNoMaterial blsBuilder;
165  TrajectoryStateClosestToBeamLine tscb( blsBuilder(initialFTS, *theBeamSpotHandle) );
166 
167  if( tscb.isValid() ) {
169  theTrackRefs.push_back( std::move(tmpRef) );
170  theTransTracks.push_back( std::move(tmpTk) );
171  }
172  }
173  }
174  }
175 
176  // Good tracks have now been selected for vertexing. Move on to vertex fitting.
177 
178 
179  // Loop over tracks and vertex good charged track pairs
180  for(unsigned int trdx1 = 0; trdx1 < theTrackRefs.size(); trdx1++) {
181 
182  for(unsigned int trdx2 = trdx1 + 1; trdx2 < theTrackRefs.size(); trdx2++) {
183 
184 
185  TrackRef positiveTrackRef;
186  TrackRef negativeTrackRef;
187  TransientTrack* posTransTkPtr = nullptr;
188  TransientTrack* negTransTkPtr = nullptr;
189 
190  // Look at the two tracks we're looping over. If they're oppositely
191  // charged, load them into the hypothesized positive and negative tracks
192  // and references to be sent to the KalmanVertexFitter
193  if(theTrackRefs[trdx1]->charge() < 0. &&
194  theTrackRefs[trdx2]->charge() > 0.) {
195  negativeTrackRef = theTrackRefs[trdx1];
196  positiveTrackRef = theTrackRefs[trdx2];
197  negTransTkPtr = &theTransTracks[trdx1];
198  posTransTkPtr = &theTransTracks[trdx2];
199  }
200  else if(theTrackRefs[trdx1]->charge() > 0. &&
201  theTrackRefs[trdx2]->charge() < 0.) {
202  negativeTrackRef = theTrackRefs[trdx2];
203  positiveTrackRef = theTrackRefs[trdx1];
204  negTransTkPtr = &theTransTracks[trdx2];
205  posTransTkPtr = &theTransTracks[trdx1];
206  }
207  // If they're not 2 oppositely charged tracks, loop back to the
208  // beginning and try the next pair.
209  else continue;
210 
211  // Fill the vector of TransientTracks to send to KVF
212  std::vector<TransientTrack> transTracks; transTracks.reserve(2);
213  transTracks.push_back(*posTransTkPtr);
214  transTracks.push_back(*negTransTkPtr);
215 
216  // Trajectory states to calculate DCA for the 2 tracks
217  FreeTrajectoryState const & posState = posTransTkPtr->impactPointTSCP().theState();
218  FreeTrajectoryState const & negState = negTransTkPtr->impactPointTSCP().theState();
219 
220  if( !posTransTkPtr->impactPointTSCP().isValid() || !negTransTkPtr->impactPointTSCP().isValid() ) continue;
221 
222  // Measure distance between tracks at their closest approach
224  cApp.calculate(posState, negState);
225  if( !cApp.status() ) continue;
226  float dca = fabs( cApp.distance() );
227  GlobalPoint cxPt = cApp.crossingPoint();
228 
229  if (dca < 0. || dca > tkDCACut) continue;
230  if (sqrt( cxPt.x()*cxPt.x() + cxPt.y()*cxPt.y() ) > 120.
231  || std::abs(cxPt.z()) > 300.) continue;
232 
233  // Get trajectory states for the tracks at POCA for later cuts
234  TrajectoryStateClosestToPoint const & posTSCP =
235  posTransTkPtr->trajectoryStateClosestToPoint( cxPt );
236  TrajectoryStateClosestToPoint const & negTSCP =
237  negTransTkPtr->trajectoryStateClosestToPoint( cxPt );
238 
239  if( !posTSCP.isValid() || !negTSCP.isValid() ) continue;
240 
241 
242  /*double posESq = posTSCP.momentum().mag2() + piMassSquared;
243  double negESq = negTSCP.momentum().mag2() + piMassSquared;
244  double posE = sqrt(posESq);
245  double negE = sqrt(negESq);
246  double totalE = posE + negE;*/
247  double totalE = sqrt( posTSCP.momentum().mag2() + piMassSquared ) +
248  sqrt( negTSCP.momentum().mag2() + piMassSquared );
249  double totalESq = totalE*totalE;
250  double totalPSq =
251  ( posTSCP.momentum() + negTSCP.momentum() ).mag2();
252  double mass = sqrt( totalESq - totalPSq);
253 
254  //mPiPiMassOut << mass << std::endl;
255 
256  if( mass > mPiPiCut ) continue;
257 
258  // Create the vertex fitter object and vertex the tracks
259  TransientVertex theRecoVertex;
260  if(vtxFitter == std::string("KalmanVertexFitter")) {
261  KalmanVertexFitter theKalmanFitter(useRefTrax == 0 ? false : true);
262  theRecoVertex = theKalmanFitter.vertex(transTracks);
263  }
264  else if (vtxFitter == std::string("AdaptiveVertexFitter")) {
265  useRefTrax = false;
266  AdaptiveVertexFitter theAdaptiveFitter;
267  theRecoVertex = theAdaptiveFitter.vertex(transTracks);
268  }
269 
270  // If the vertex is valid, make a VertexCompositeCandidate with it
271 
272  if( !theRecoVertex.isValid() || theRecoVertex.totalChiSquared() < 0. ) {
273  continue;
274  }
275 
276  // Create reco::Vertex object for use in creating the Candidate
277  reco::Vertex theVtx = theRecoVertex;
278  // Create and fill vector of refitted TransientTracks
279  // (iff they've been created by the KVF)
280  std::vector<TransientTrack> refittedTrax;
281  if( theRecoVertex.hasRefittedTracks() ) {
282  refittedTrax = theRecoVertex.refittedTracks();
283  }
284 
285  // Do post-fit cuts if specified in config file.
286 
287  // Find the vertex d0 and its error
288 
289  typedef ROOT::Math::SMatrix<double, 3, 3, ROOT::Math::MatRepSym<double, 3> > SMatrixSym3D;
290  typedef ROOT::Math::SVector<double, 3> SVector3;
291 
292  GlobalPoint vtxPos(theVtx.x(), theVtx.y(), theVtx.z());
293 
294  GlobalPoint beamSpotPos(theBeamSpotHandle->position().x(),
295  theBeamSpotHandle->position().y(),
296  theBeamSpotHandle->position().z());
297 
298  SMatrixSym3D totalCov = theBeamSpotHandle->rotatedCovariance3D() + theVtx.covariance();
299  SVector3 distanceVector(vtxPos.x() - beamSpotPos.x(),
300  vtxPos.y() - beamSpotPos.y(),
301  0.);//so that we get radial values only,
302  //since z beamSpot uncertainty is huge
303 
304  double rVtxMag = ROOT::Math::Mag(distanceVector);
305  double sigmaRvtxMag = sqrt(ROOT::Math::Similarity(totalCov, distanceVector)) / rVtxMag;
306 
307  // The methods innerOk() and innerPosition() require TrackExtra, which
308  // is only available in the RECO data tier, not AOD. Setting innerHitPosCut
309  // to -1 avoids this problem and allows to run on AOD.
310  if( innerHitPosCut > 0. && positiveTrackRef->innerOk() ) {
311  reco::Vertex::Point posTkHitPos = positiveTrackRef->innerPosition();
312  double posTkHitPosD2 =
313  (posTkHitPos.x()-beamSpotPos.x())*(posTkHitPos.x()-beamSpotPos.x()) +
314  (posTkHitPos.y()-beamSpotPos.y())*(posTkHitPos.y()-beamSpotPos.y());
315  if( sqrt( posTkHitPosD2 ) < ( rVtxMag - sigmaRvtxMag*innerHitPosCut )
316  ) {
317  continue;
318  }
319  }
320  if( innerHitPosCut > 0. && negativeTrackRef->innerOk() ) {
321  reco::Vertex::Point negTkHitPos = negativeTrackRef->innerPosition();
322  double negTkHitPosD2 =
323  (negTkHitPos.x()-beamSpotPos.x())*(negTkHitPos.x()-beamSpotPos.x()) +
324  (negTkHitPos.y()-beamSpotPos.y())*(negTkHitPos.y()-beamSpotPos.y());
325  if( sqrt( negTkHitPosD2 ) < ( rVtxMag - sigmaRvtxMag*innerHitPosCut )
326  ) {
327  continue;
328  }
329  }
330 
331  if( theVtx.normalizedChi2() > chi2Cut ||
332  rVtxMag < rVtxCut ||
333  rVtxMag / sigmaRvtxMag < vtxSigCut ) {
334  continue;
335  }
336 
337  // Cuts finished, now we create the candidates and push them back into the collections.
338 
339  std::auto_ptr<TrajectoryStateClosestToPoint> trajPlus;
340  std::auto_ptr<TrajectoryStateClosestToPoint> trajMins;
341 
342  if( useRefTrax && refittedTrax.size() > 1 ) {
343  // Need an iterator over the refitted tracks for below
344  std::vector<TransientTrack>::iterator traxIter = refittedTrax.begin(),
345  traxEnd = refittedTrax.end();
346 
347  // TransientTrack objects to hold the positive and negative
348  // refitted tracks
349  TransientTrack* thePositiveRefTrack = 0;
350  TransientTrack* theNegativeRefTrack = 0;
351 
352  for( ; traxIter != traxEnd; ++traxIter) {
353  if( traxIter->track().charge() > 0. ) {
354  thePositiveRefTrack = &*traxIter;
355  }
356  else if (traxIter->track().charge() < 0.) {
357  theNegativeRefTrack = &*traxIter;
358  }
359  }
360  if (thePositiveRefTrack == 0 || theNegativeRefTrack == 0) continue;
361  trajPlus.reset(new TrajectoryStateClosestToPoint(thePositiveRefTrack->trajectoryStateClosestToPoint(vtxPos)));
362  trajMins.reset(new TrajectoryStateClosestToPoint(theNegativeRefTrack->trajectoryStateClosestToPoint(vtxPos)));
363  }
364  else {
365  trajPlus.reset(new TrajectoryStateClosestToPoint(posTransTkPtr->trajectoryStateClosestToPoint(vtxPos)));
366  trajMins.reset(new TrajectoryStateClosestToPoint(negTransTkPtr->trajectoryStateClosestToPoint(vtxPos)));
367 
368  }
369 
370  if( trajPlus.get() == 0 || trajMins.get() == 0 || !trajPlus->isValid() || !trajMins->isValid() ) continue;
371 
372  posTransTkPtr = negTransTkPtr = 0;
373 
374  GlobalVector positiveP(trajPlus->momentum());
375  GlobalVector negativeP(trajMins->momentum());
376  GlobalVector totalP(positiveP + negativeP);
377 
378  //cleanup stuff we don't need anymore
379  trajPlus.reset();
380  trajMins.reset();
381 
382  // calculate total energy of V0 3 ways:
383  // Assume it's a kShort, a Lambda, or a LambdaBar.
384  double piPlusE = sqrt( positiveP.mag2() + piMassSquared );
385  double piMinusE = sqrt( negativeP.mag2() + piMassSquared );
386  double protonE = sqrt( positiveP.mag2() + protonMassSquared );
387  double antiProtonE = sqrt( negativeP.mag2() + protonMassSquared );
388  double kShortETot = piPlusE + piMinusE;
389  double lambdaEtot = protonE + piMinusE;
390  double lambdaBarEtot = antiProtonE + piPlusE;
391 
392  using namespace reco;
393 
394  // Create momentum 4-vectors for the 3 candidate types
395  const Particle::LorentzVector kShortP4(totalP.x(),
396  totalP.y(), totalP.z(),
397  kShortETot);
398  const Particle::LorentzVector lambdaP4(totalP.x(),
399  totalP.y(), totalP.z(),
400  lambdaEtot);
401  const Particle::LorentzVector lambdaBarP4(totalP.x(),
402  totalP.y(), totalP.z(),
403  lambdaBarEtot);
404 
405  Particle::Point vtx(theVtx.x(), theVtx.y(), theVtx.z());
406  const Vertex::CovarianceMatrix vtxCov(theVtx.covariance());
407  double vtxChi2(theVtx.chi2());
408  double vtxNdof(theVtx.ndof());
409 
410  // Create the VertexCompositeCandidate object that will be stored in the Event
411  VertexCompositeCandidate* theKshort = nullptr;
412  VertexCompositeCandidate* theLambda = nullptr;
413  VertexCompositeCandidate* theLambdaBar = nullptr;
414 
415  if( doKshorts ) {
416  theKshort = new VertexCompositeCandidate(0, kShortP4, vtx, vtxCov, vtxChi2, vtxNdof);
417  }
418  if( doLambdas ) {
419  if( positiveP.mag2() > negativeP.mag2() ) {
420  theLambda =
421  new VertexCompositeCandidate(0, lambdaP4, vtx, vtxCov, vtxChi2, vtxNdof);
422  }
423  else {
424  theLambdaBar =
425  new VertexCompositeCandidate(0, lambdaBarP4, vtx, vtxCov, vtxChi2, vtxNdof);
426  }
427  }
428 
429  // Create daughter candidates for the VertexCompositeCandidates
431  thePiPlusCand(1, Particle::LorentzVector(positiveP.x(),
432  positiveP.y(), positiveP.z(),
433  piPlusE), vtx);
434  thePiPlusCand.setTrack(positiveTrackRef);
435 
437  thePiMinusCand(-1, Particle::LorentzVector(negativeP.x(),
438  negativeP.y(), negativeP.z(),
439  piMinusE), vtx);
440  thePiMinusCand.setTrack(negativeTrackRef);
441 
442 
444  theProtonCand(1, Particle::LorentzVector(positiveP.x(),
445  positiveP.y(), positiveP.z(),
446  protonE), vtx);
447  theProtonCand.setTrack(positiveTrackRef);
448 
450  theAntiProtonCand(-1, Particle::LorentzVector(negativeP.x(),
451  negativeP.y(), negativeP.z(),
452  antiProtonE), vtx);
453  theAntiProtonCand.setTrack(negativeTrackRef);
454 
455 
456  AddFourMomenta addp4;
457  // Store the daughter Candidates in the VertexCompositeCandidates
458  // if they pass mass cuts
459  if( doKshorts ) {
460  theKshort->addDaughter(thePiPlusCand);
461  theKshort->addDaughter(thePiMinusCand);
462  theKshort->setPdgId(310);
463  addp4.set( *theKshort );
464  if( theKshort->mass() < kShortMass + kShortMassCut &&
465  theKshort->mass() > kShortMass - kShortMassCut ) {
466  theKshorts.push_back( std::move(*theKshort) );
467  }
468  }
469 
470  if( doLambdas && theLambda ) {
471  theLambda->addDaughter(theProtonCand);
472  theLambda->addDaughter(thePiMinusCand);
473  theLambda->setPdgId(3122);
474  addp4.set( *theLambda );
475  if( theLambda->mass() < lambdaMass + lambdaMassCut &&
476  theLambda->mass() > lambdaMass - lambdaMassCut ) {
477  theLambdas.push_back( std::move(*theLambda) );
478  }
479  }
480  else if ( doLambdas && theLambdaBar ) {
481  theLambdaBar->addDaughter(theAntiProtonCand);
482  theLambdaBar->addDaughter(thePiPlusCand);
483  theLambdaBar->setPdgId(-3122);
484  addp4.set( *theLambdaBar );
485  if( theLambdaBar->mass() < lambdaMass + lambdaMassCut &&
486  theLambdaBar->mass() > lambdaMass - lambdaMassCut ) {
487  theLambdas.push_back( std::move(*theLambdaBar) );
488  }
489  }
490 
491  delete theKshort;
492  delete theLambda;
493  delete theLambdaBar;
494  theKshort = theLambda = theLambdaBar = nullptr;
495 
496  }
497  }
498 }
499 
500 
501 
502 // Experimental
503 double V0Fitter::findV0MassError(const GlobalPoint &vtxPos, const std::vector<reco::TransientTrack> &dauTracks) {
504  return -1.;
505 }
506 
507 /*
508 double V0Fitter::findV0MassError(const GlobalPoint &vtxPos, std::vector<reco::TransientTrack> const & dauTracks) {
509  // Returns -99999. if trajectory states fail at vertex position
510 
511  // Load positive track trajectory at vertex into vector, then negative track
512  std::vector<TrajectoryStateClosestToPoint> sortedTrajStatesAtVtx;
513  for( unsigned int ndx = 0; ndx < dauTracks.size(); ndx++ ) {
514  if( dauTracks[ndx].trajectoryStateClosestToPoint(vtxPos).isValid() ) {
515  std::cout << "From TSCP: "
516  << dauTracks[ndx].trajectoryStateClosestToPoint(vtxPos).perigeeParameters().transverseCurvature()
517  << "; From Track: " << dauTracks[ndx].track().qoverp() << std::endl;
518  }
519  if( sortedTrajStatesAtVtx.size() == 0 ) {
520  if( dauTracks[ndx].charge() > 0 ) {
521  sortedTrajStatesAtVtx.push_back( dauTracks[ndx].trajectoryStateClosestToPoint(vtxPos) );
522  }
523  else {
524  sortedTrajStatesAtVtx.push_back( dauTracks[ndx].trajectoryStateClosestToPoint(vtxPos) );
525  }
526  }
527  }
528  std::vector<PerigeeTrajectoryParameters> param;
529  std::vector<PerigeeTrajectoryError> paramError;
530  std::vector<GlobalVector> momenta;
531 
532  for( unsigned int ndx2 = 0; ndx2 < sortedTrajStatesAtVtx.size(); ndx2++ ) {
533  if( sortedTrajStatesAtVtx[ndx2].isValid() ) {
534  param.push_back( sortedTrajStatesAtVtx[ndx2].perigeeParameters() );
535  paramError.push_back( sortedTrajStatesAtVtx[ndx2].perigeeError() );
536  momenta.push_back( sortedTrajStatesAtVtx[ndx2].momentum() );
537  }
538  else return -99999.;
539  }
540  return 0;
541 }
542 */
543 
544 
545 
T getParameter(std::string const &) const
bool useRefTrax
Definition: V0Fitter.h:76
T mag2() const
Definition: PV3DBase.h:66
std::vector< VertexCompositeCandidate > VertexCompositeCandidateCollection
collection of Candidate objects
const double protonMass
int tkNhitsCut
Definition: V0Fitter.h:87
~V0Fitter()
Definition: V0Fitter.cc:102
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:446
TrajectoryStateClosestToPoint impactPointTSCP() const
const FreeTrajectoryState & theState() const
double y() const
y coordinate
Definition: Vertex.h:110
virtual CachingVertex< 5 > vertex(const std::vector< reco::TransientTrack > &tracks) const
float totalChiSquared() const
const MagneticField * magField
Definition: V0Fitter.h:74
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:13
double chi2Cut
Definition: V0Fitter.h:85
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
edm::InputTag vtxFitter
Definition: V0Fitter.h:103
double rVtxCut
Definition: V0Fitter.h:88
double collinCut
Definition: V0Fitter.h:91
math::Error< dimension >::type CovarianceMatrix
covariance error matrix (3x3)
Definition: Vertex.h:45
virtual GlobalPoint crossingPoint() const
bool doKshorts
Definition: V0Fitter.h:78
double findV0MassError(const GlobalPoint &vtxPos, const std::vector< reco::TransientTrack > &dauTracks)
Definition: V0Fitter.cc:503
double tkChi2Cut
Definition: V0Fitter.h:86
int iEvent
Definition: GenABIO.cc:230
double impactParameterSigCut
Definition: V0Fitter.h:94
edm::EDGetTokenT< reco::TrackCollection > token_tracks
Definition: V0Fitter.h:101
T sqrt(T t)
Definition: SSEVec.h:48
bool doLambdas
Definition: V0Fitter.h:79
T mag2() const
The vector magnitude squared. Equivalent to vec.dot(vec)
const TrackerGeometry * trackerGeom
Definition: V0Fitter.h:72
virtual CachingVertex< 5 > vertex(const std::vector< reco::TransientTrack > &) const
math::XYZPoint Point
point in the space
Definition: Particle.h:31
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
void fitAll(const edm::Event &iEvent, const edm::EventSetup &iSetup, reco::VertexCompositeCandidateCollection &k, reco::VertexCompositeCandidateCollection &l)
Definition: V0Fitter.cc:106
math::XYZPoint Point
point in the space
Definition: Vertex.h:39
V0Fitter(const edm::ParameterSet &theParams, edm::ConsumesCollector &&iC)
Definition: V0Fitter.cc:49
std::vector< reco::TrackBase::TrackQuality > qualities
Definition: V0Fitter.h:99
double ndof() const
Definition: Vertex.h:102
double x() const
x coordinate
Definition: Vertex.h:108
double significance() const
Definition: Measurement1D.h:32
static TrackQuality qualityByName(const std::string &name)
Definition: TrackBase.cc:102
virtual bool calculate(const TrajectoryStateOnSurface &sta, const TrajectoryStateOnSurface &stb)
double tkDCACut
Definition: V0Fitter.h:96
FreeTrajectoryState initialFreeState(const reco::Track &tk, const MagneticField *field)
TrajectoryStateClosestToPoint trajectoryStateClosestToPoint(const GlobalPoint &point) const
const T & get() const
Definition: EventSetup.h:55
double kShortMassCut
Definition: V0Fitter.h:92
double mPiPiCut
Definition: V0Fitter.h:95
tuple cout
Definition: gather_cfg.py:121
edm::EDGetTokenT< reco::BeamSpot > token_beamSpot
Definition: V0Fitter.h:102
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
double vtxSigCut
Definition: V0Fitter.h:89
void setTrack(const reco::TrackRef &r)
set reference to track
virtual float distance() const
double innerHitPosCut
Definition: V0Fitter.h:97
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Particle.h:27
bool isValid() const
virtual bool status() const
double lambdaMassCut
Definition: V0Fitter.h:93