CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_7/src/RecoVertex/V0Producer/src/V0Fitter.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:    V0Producer
00004 // Class:      V0Fitter
00005 // 
00013 //
00014 // Original Author:  Brian Drell
00015 //         Created:  Fri May 18 22:57:40 CEST 2007
00016 // $Id: V0Fitter.cc,v 1.57 2013/04/11 23:17:44 wmtan Exp $
00017 //
00018 //
00019 
00020 #include "RecoVertex/V0Producer/interface/V0Fitter.h"
00021 #include "CommonTools/CandUtils/interface/AddFourMomenta.h"
00022 
00023 #include "TrackingTools/TransientTrack/interface/TransientTrackBuilder.h"
00024 #include "TrackingTools/Records/interface/TransientTrackRecord.h"
00025 #include "TrackingTools/PatternTools/interface/ClosestApproachInRPhi.h"
00026 #include "Geometry/CommonDetUnit/interface/GlobalTrackingGeometry.h"
00027 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00028 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateTransform.h"
00029 #include "TrackingTools/PatternTools/interface/TSCBLBuilderNoMaterial.h"
00030 
00031 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
00032 
00033 #include <Math/Functions.h>
00034 #include <Math/SVector.h>
00035 #include <Math/SMatrix.h>
00036 #include <typeinfo>
00037 #include <memory>
00038 
00039 // Constants
00040 
00041 const double piMass = 0.13957018;
00042 const double piMassSquared = piMass*piMass;
00043 const double protonMass = 0.938272013;
00044 const double protonMassSquared = protonMass*protonMass;
00045 const double kShortMass = 0.497614;
00046 const double lambdaMass = 1.115683;
00047 
00048 // Constructor and (empty) destructor
00049 V0Fitter::V0Fitter(const edm::ParameterSet& theParameters,
00050                    const edm::Event& iEvent, const edm::EventSetup& iSetup) {
00051   using std::string;
00052 
00053   // Get the track reco algorithm from the ParameterSet
00054   recoAlg = theParameters.getParameter<edm::InputTag>("trackRecoAlgorithm");
00055 
00056   // ------> Initialize parameters from PSet. ALL TRACKED, so no defaults.
00057   // First set bits to do various things:
00058   //  -decide whether to use the KVF track smoother, and whether to store those
00059   //     tracks in the reco::Vertex
00060   useRefTrax = theParameters.getParameter<bool>(string("useSmoothing"));
00061 
00062   //  -whether to reconstruct K0s
00063   doKshorts = theParameters.getParameter<bool>(string("selectKshorts"));
00064   //  -whether to reconstruct Lambdas
00065   doLambdas = theParameters.getParameter<bool>(string("selectLambdas"));
00066 
00067   // Second, initialize post-fit cuts
00068   chi2Cut = theParameters.getParameter<double>(string("vtxChi2Cut"));
00069   tkChi2Cut = theParameters.getParameter<double>(string("tkChi2Cut"));
00070   tkNhitsCut = theParameters.getParameter<int>(string("tkNhitsCut"));
00071   rVtxCut = theParameters.getParameter<double>(string("rVtxCut"));
00072   vtxSigCut = theParameters.getParameter<double>(string("vtxSignificance2DCut"));
00073   collinCut = theParameters.getParameter<double>(string("collinearityCut"));
00074   kShortMassCut = theParameters.getParameter<double>(string("kShortMassCut"));
00075   lambdaMassCut = theParameters.getParameter<double>(string("lambdaMassCut"));
00076   impactParameterSigCut = theParameters.getParameter<double>(string("impactParameterSigCut"));
00077   mPiPiCut = theParameters.getParameter<double>(string("mPiPiCut"));
00078   tkDCACut = theParameters.getParameter<double>(string("tkDCACut"));
00079   vtxFitter = theParameters.getParameter<edm::InputTag>("vertexFitter");
00080   innerHitPosCut = theParameters.getParameter<double>(string("innerHitPosCut"));
00081   std::vector<std::string> qual = theParameters.getParameter<std::vector<std::string> >("trackQualities");
00082   for (unsigned int ndx = 0; ndx < qual.size(); ndx++) {
00083     qualities.push_back(reco::TrackBase::qualityByName(qual[ndx]));
00084   }
00085 
00086   //edm::LogInfo("V0Producer") << "Using " << vtxFitter << " to fit V0 vertices.\n";
00087   //std::cout << "Using " << vtxFitter << " to fit V0 vertices." << std::endl;
00088   // FOR DEBUG:
00089   //initFileOutput();
00090   //--------------------
00091 
00092   //std::cout << "Entering V0Producer" << std::endl;
00093 
00094   fitAll(iEvent, iSetup);
00095 
00096   // FOR DEBUG:
00097   //cleanupFileOutput();
00098   //--------------------
00099 
00100 }
00101 
00102 V0Fitter::~V0Fitter() {
00103 }
00104 
00105 // Method containing the algorithm for vertex reconstruction
00106 void V0Fitter::fitAll(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
00107 
00108   using std::vector;
00109   using std::cout;
00110   using std::endl;
00111   using namespace reco;
00112   using namespace edm;
00113 
00114   // Create std::vectors for Tracks and TrackRefs (required for
00115   //  passing to the KalmanVertexFitter)
00116   std::vector<TrackRef> theTrackRefs;
00117   std::vector<TransientTrack> theTransTracks;
00118 
00119   // Handles for tracks, B-field, and tracker geometry
00120   Handle<reco::TrackCollection> theTrackHandle;
00121   Handle<reco::BeamSpot> theBeamSpotHandle;
00122   ESHandle<MagneticField> bFieldHandle;
00123   ESHandle<TrackerGeometry> trackerGeomHandle;
00124   ESHandle<GlobalTrackingGeometry> globTkGeomHandle;
00125   //cout << "Check 0" << endl;
00126 
00127   // Get the tracks from the event, and get the B-field record
00128   //  from the EventSetup
00129   iEvent.getByLabel(recoAlg, theTrackHandle);
00130   iEvent.getByLabel(std::string("offlineBeamSpot"), theBeamSpotHandle);
00131   if( !theTrackHandle->size() ) return;
00132   iSetup.get<IdealMagneticFieldRecord>().get(bFieldHandle);
00133   iSetup.get<TrackerDigiGeometryRecord>().get(trackerGeomHandle);
00134   iSetup.get<GlobalTrackingGeometryRecord>().get(globTkGeomHandle);
00135 
00136   trackerGeom = trackerGeomHandle.product();
00137   magField = bFieldHandle.product();
00138 
00139   // Fill vectors of TransientTracks and TrackRefs after applying preselection cuts.
00140   for(unsigned int indx = 0; indx < theTrackHandle->size(); indx++) {
00141     TrackRef tmpRef( theTrackHandle, indx );
00142     bool quality_ok = true;
00143     if (qualities.size()!=0) {
00144       quality_ok = false;
00145       for (unsigned int ndx_ = 0; ndx_ < qualities.size(); ndx_++) {
00146         if (tmpRef->quality(qualities[ndx_])){
00147           quality_ok = true;
00148           break;          
00149         }
00150       }
00151     }
00152     if( !quality_ok ) continue;
00153 
00154 
00155     if( tmpRef->normalizedChi2() < tkChi2Cut &&
00156         tmpRef->numberOfValidHits() >= tkNhitsCut ) {
00157       TransientTrack tmpTk( *tmpRef, &(*bFieldHandle), globTkGeomHandle );
00158       
00159       FreeTrajectoryState initialFTS = trajectoryStateTransform::initialFreeState(*tmpRef, magField);
00160       TSCBLBuilderNoMaterial blsBuilder;
00161       TrajectoryStateClosestToBeamLine tscb( blsBuilder(initialFTS, *theBeamSpotHandle) );
00162       
00163       if( tscb.isValid() ) {
00164         if( tscb.transverseImpactParameter().significance() > impactParameterSigCut ) {
00165           theTrackRefs.push_back( tmpRef );
00166           theTransTracks.push_back( tmpTk );
00167         }
00168       }
00169     }
00170   }
00171 
00172   // Good tracks have now been selected for vertexing.  Move on to vertex fitting.
00173 
00174 
00175   // Loop over tracks and vertex good charged track pairs
00176   for(unsigned int trdx1 = 0; trdx1 < theTrackRefs.size(); trdx1++) {
00177 
00178     for(unsigned int trdx2 = trdx1 + 1; trdx2 < theTrackRefs.size(); trdx2++) {
00179 
00180       //This vector holds the pair of oppositely-charged tracks to be vertexed
00181       std::vector<TransientTrack> transTracks;
00182 
00183       TrackRef positiveTrackRef;
00184       TrackRef negativeTrackRef;
00185       TransientTrack* posTransTkPtr = 0;
00186       TransientTrack* negTransTkPtr = 0;
00187 
00188       // Look at the two tracks we're looping over.  If they're oppositely
00189       //  charged, load them into the hypothesized positive and negative tracks
00190       //  and references to be sent to the KalmanVertexFitter
00191       if(theTrackRefs[trdx1]->charge() < 0. && 
00192          theTrackRefs[trdx2]->charge() > 0.) {
00193         negativeTrackRef = theTrackRefs[trdx1];
00194         positiveTrackRef = theTrackRefs[trdx2];
00195         negTransTkPtr = &theTransTracks[trdx1];
00196         posTransTkPtr = &theTransTracks[trdx2];
00197       }
00198       else if(theTrackRefs[trdx1]->charge() > 0. &&
00199               theTrackRefs[trdx2]->charge() < 0.) {
00200         negativeTrackRef = theTrackRefs[trdx2];
00201         positiveTrackRef = theTrackRefs[trdx1];
00202         negTransTkPtr = &theTransTracks[trdx2];
00203         posTransTkPtr = &theTransTracks[trdx1];
00204       }
00205       // If they're not 2 oppositely charged tracks, loop back to the
00206       //  beginning and try the next pair.
00207       else continue;
00208 
00209       // Fill the vector of TransientTracks to send to KVF
00210       transTracks.push_back(*posTransTkPtr);
00211       transTracks.push_back(*negTransTkPtr);
00212 
00213       // Trajectory states to calculate DCA for the 2 tracks
00214       FreeTrajectoryState posState = posTransTkPtr->impactPointTSCP().theState();
00215       FreeTrajectoryState negState = negTransTkPtr->impactPointTSCP().theState();
00216 
00217       if( !posTransTkPtr->impactPointTSCP().isValid() || !negTransTkPtr->impactPointTSCP().isValid() ) continue;
00218 
00219       // Measure distance between tracks at their closest approach
00220       ClosestApproachInRPhi cApp;
00221       cApp.calculate(posState, negState);
00222       if( !cApp.status() ) continue;
00223       float dca = fabs( cApp.distance() );
00224       GlobalPoint cxPt = cApp.crossingPoint();
00225 
00226       if (dca < 0. || dca > tkDCACut) continue;
00227       if (sqrt( cxPt.x()*cxPt.x() + cxPt.y()*cxPt.y() ) > 120. 
00228           || std::abs(cxPt.z()) > 300.) continue;
00229 
00230       // Get trajectory states for the tracks at POCA for later cuts
00231       TrajectoryStateClosestToPoint posTSCP = 
00232         posTransTkPtr->trajectoryStateClosestToPoint( cxPt );
00233       TrajectoryStateClosestToPoint negTSCP =
00234         negTransTkPtr->trajectoryStateClosestToPoint( cxPt );
00235 
00236       if( !posTSCP.isValid() || !negTSCP.isValid() ) continue;
00237 
00238 
00239       /*double posESq = posTSCP.momentum().mag2() + piMassSquared;
00240       double negESq = negTSCP.momentum().mag2() + piMassSquared;
00241       double posE = sqrt(posESq);
00242       double negE = sqrt(negESq);
00243       double totalE = posE + negE;*/
00244       double totalE = sqrt( posTSCP.momentum().mag2() + piMassSquared ) +
00245                       sqrt( negTSCP.momentum().mag2() + piMassSquared );
00246       double totalESq = totalE*totalE;
00247       double totalPSq =
00248         ( posTSCP.momentum() + negTSCP.momentum() ).mag2();
00249       double mass = sqrt( totalESq - totalPSq);
00250 
00251       //mPiPiMassOut << mass << std::endl;
00252 
00253       if( mass > mPiPiCut ) continue;
00254 
00255       // Create the vertex fitter object and vertex the tracks
00256       TransientVertex theRecoVertex;
00257       if(vtxFitter == std::string("KalmanVertexFitter")) {
00258         KalmanVertexFitter theKalmanFitter(useRefTrax == 0 ? false : true);
00259         theRecoVertex = theKalmanFitter.vertex(transTracks);
00260       }
00261       else if (vtxFitter == std::string("AdaptiveVertexFitter")) {
00262         useRefTrax = false;
00263         AdaptiveVertexFitter theAdaptiveFitter;
00264         theRecoVertex = theAdaptiveFitter.vertex(transTracks);
00265       }
00266     
00267       // If the vertex is valid, make a VertexCompositeCandidate with it
00268 
00269       if( !theRecoVertex.isValid() || theRecoVertex.totalChiSquared() < 0. ) {
00270         continue;
00271       }
00272 
00273       // Create reco::Vertex object for use in creating the Candidate
00274       reco::Vertex theVtx = theRecoVertex;
00275       // Create and fill vector of refitted TransientTracks
00276       //  (iff they've been created by the KVF)
00277       std::vector<TransientTrack> refittedTrax;
00278       if( theRecoVertex.hasRefittedTracks() ) {
00279         refittedTrax = theRecoVertex.refittedTracks();
00280       }
00281 
00282       // Do post-fit cuts if specified in config file.
00283 
00284       // Find the vertex d0 and its error
00285 
00286       typedef ROOT::Math::SMatrix<double, 3, 3, ROOT::Math::MatRepSym<double, 3> > SMatrixSym3D;
00287       typedef ROOT::Math::SVector<double, 3> SVector3;
00288 
00289       GlobalPoint vtxPos(theVtx.x(), theVtx.y(), theVtx.z());
00290 
00291       GlobalPoint beamSpotPos(theBeamSpotHandle->position().x(),
00292                               theBeamSpotHandle->position().y(),
00293                               theBeamSpotHandle->position().z());
00294 
00295       SMatrixSym3D totalCov = theBeamSpotHandle->rotatedCovariance3D() + theVtx.covariance();
00296       SVector3 distanceVector(vtxPos.x() - beamSpotPos.x(),
00297                               vtxPos.y() - beamSpotPos.y(),
00298                               0.);//so that we get radial values only, 
00299                                   //since z beamSpot uncertainty is huge
00300 
00301       double rVtxMag = ROOT::Math::Mag(distanceVector);
00302       double sigmaRvtxMag = sqrt(ROOT::Math::Similarity(totalCov, distanceVector)) / rVtxMag;
00303       
00304       // The methods innerOk() and innerPosition() require TrackExtra, which
00305       // is only available in the RECO data tier, not AOD. Setting innerHitPosCut
00306       // to -1 avoids this problem and allows to run on AOD.
00307       if( innerHitPosCut > 0. && positiveTrackRef->innerOk() ) {
00308         reco::Vertex::Point posTkHitPos = positiveTrackRef->innerPosition();
00309         double posTkHitPosD2 = 
00310           (posTkHitPos.x()-beamSpotPos.x())*(posTkHitPos.x()-beamSpotPos.x()) +
00311           (posTkHitPos.y()-beamSpotPos.y())*(posTkHitPos.y()-beamSpotPos.y());
00312         if( sqrt( posTkHitPosD2 ) < ( rVtxMag - sigmaRvtxMag*innerHitPosCut )
00313             ) {
00314           continue;
00315         }
00316       }
00317       if( innerHitPosCut > 0. && negativeTrackRef->innerOk() ) {
00318         reco::Vertex::Point negTkHitPos = negativeTrackRef->innerPosition();
00319         double negTkHitPosD2 = 
00320           (negTkHitPos.x()-beamSpotPos.x())*(negTkHitPos.x()-beamSpotPos.x()) +
00321           (negTkHitPos.y()-beamSpotPos.y())*(negTkHitPos.y()-beamSpotPos.y());
00322         if( sqrt( negTkHitPosD2 ) < ( rVtxMag - sigmaRvtxMag*innerHitPosCut )
00323             ) {
00324           continue;
00325         }
00326       }
00327       
00328       if( theVtx.normalizedChi2() > chi2Cut ||
00329           rVtxMag < rVtxCut ||
00330           rVtxMag / sigmaRvtxMag < vtxSigCut ) {
00331         continue;
00332       }
00333 
00334       // Cuts finished, now we create the candidates and push them back into the collections.
00335       
00336       std::auto_ptr<TrajectoryStateClosestToPoint> trajPlus;
00337       std::auto_ptr<TrajectoryStateClosestToPoint> trajMins;
00338 
00339       if( useRefTrax && refittedTrax.size() > 1 ) {
00340         // Need an iterator over the refitted tracks for below
00341         std::vector<TransientTrack>::iterator traxIter = refittedTrax.begin(),
00342           traxEnd = refittedTrax.end();
00343 
00344         // TransientTrack objects to hold the positive and negative
00345         //  refitted tracks
00346         TransientTrack* thePositiveRefTrack = 0;
00347         TransientTrack* theNegativeRefTrack = 0;
00348         
00349         for( ; traxIter != traxEnd; ++traxIter) {
00350           if( traxIter->track().charge() > 0. ) {
00351             thePositiveRefTrack = &*traxIter;
00352           }
00353           else if (traxIter->track().charge() < 0.) {
00354             theNegativeRefTrack = &*traxIter;
00355           }
00356         }
00357         if (thePositiveRefTrack == 0 || theNegativeRefTrack == 0) continue;
00358         trajPlus.reset(new TrajectoryStateClosestToPoint(thePositiveRefTrack->trajectoryStateClosestToPoint(vtxPos)));
00359         trajMins.reset(new TrajectoryStateClosestToPoint(theNegativeRefTrack->trajectoryStateClosestToPoint(vtxPos)));
00360       }
00361       else {
00362         trajPlus.reset(new TrajectoryStateClosestToPoint(posTransTkPtr->trajectoryStateClosestToPoint(vtxPos)));
00363         trajMins.reset(new TrajectoryStateClosestToPoint(negTransTkPtr->trajectoryStateClosestToPoint(vtxPos)));
00364 
00365       }
00366 
00367       if( trajPlus.get() == 0 || trajMins.get() == 0 || !trajPlus->isValid() || !trajMins->isValid() ) continue;
00368 
00369       posTransTkPtr = negTransTkPtr = 0;
00370 
00371       GlobalVector positiveP(trajPlus->momentum());
00372       GlobalVector negativeP(trajMins->momentum());
00373       GlobalVector totalP(positiveP + negativeP);
00374 
00375       //cleanup stuff we don't need anymore
00376       trajPlus.reset();
00377       trajMins.reset();
00378 
00379       // calculate total energy of V0 3 ways:
00380       //  Assume it's a kShort, a Lambda, or a LambdaBar.
00381       double piPlusE = sqrt( positiveP.mag2() + piMassSquared );
00382       double piMinusE = sqrt( negativeP.mag2() + piMassSquared );
00383       double protonE = sqrt( positiveP.mag2() + protonMassSquared );
00384       double antiProtonE = sqrt( negativeP.mag2() + protonMassSquared );
00385       double kShortETot = piPlusE + piMinusE;
00386       double lambdaEtot = protonE + piMinusE;
00387       double lambdaBarEtot = antiProtonE + piPlusE;
00388 
00389       using namespace reco;
00390 
00391       // Create momentum 4-vectors for the 3 candidate types
00392       const Particle::LorentzVector kShortP4(totalP.x(), 
00393                                              totalP.y(), totalP.z(), 
00394                                              kShortETot);
00395       const Particle::LorentzVector lambdaP4(totalP.x(), 
00396                                              totalP.y(), totalP.z(), 
00397                                              lambdaEtot);
00398       const Particle::LorentzVector lambdaBarP4(totalP.x(), 
00399                                                 totalP.y(), totalP.z(), 
00400                                                 lambdaBarEtot);
00401 
00402       Particle::Point vtx(theVtx.x(), theVtx.y(), theVtx.z());
00403       const Vertex::CovarianceMatrix vtxCov(theVtx.covariance());
00404       double vtxChi2(theVtx.chi2());
00405       double vtxNdof(theVtx.ndof());
00406 
00407       // Create the VertexCompositeCandidate object that will be stored in the Event
00408       VertexCompositeCandidate* theKshort = 0;
00409       VertexCompositeCandidate* theLambda = 0;
00410       VertexCompositeCandidate* theLambdaBar = 0;
00411 
00412       if( doKshorts ) {
00413         theKshort = new VertexCompositeCandidate(0, kShortP4, vtx, vtxCov, vtxChi2, vtxNdof);
00414       }
00415       if( doLambdas ) {
00416         if( positiveP.mag() > negativeP.mag() ) {
00417           theLambda = 
00418             new VertexCompositeCandidate(0, lambdaP4, vtx, vtxCov, vtxChi2, vtxNdof);
00419         }
00420         else {
00421           theLambdaBar = 
00422             new VertexCompositeCandidate(0, lambdaBarP4, vtx, vtxCov, vtxChi2, vtxNdof);
00423         }
00424       }
00425 
00426       // Create daughter candidates for the VertexCompositeCandidates
00427       RecoChargedCandidate 
00428         thePiPlusCand(1, Particle::LorentzVector(positiveP.x(), 
00429                                                  positiveP.y(), positiveP.z(),
00430                                                  piPlusE), vtx);
00431       thePiPlusCand.setTrack(positiveTrackRef);
00432       
00433       RecoChargedCandidate
00434         thePiMinusCand(-1, Particle::LorentzVector(negativeP.x(), 
00435                                                    negativeP.y(), negativeP.z(),
00436                                                    piMinusE), vtx);
00437       thePiMinusCand.setTrack(negativeTrackRef);
00438       
00439  
00440       RecoChargedCandidate
00441         theProtonCand(1, Particle::LorentzVector(positiveP.x(),
00442                                                  positiveP.y(), positiveP.z(),
00443                                                  protonE), vtx);
00444       theProtonCand.setTrack(positiveTrackRef);
00445 
00446       RecoChargedCandidate
00447         theAntiProtonCand(-1, Particle::LorentzVector(negativeP.x(),
00448                                                       negativeP.y(), negativeP.z(),
00449                                                       antiProtonE), vtx);
00450       theAntiProtonCand.setTrack(negativeTrackRef);
00451 
00452 
00453       AddFourMomenta addp4;
00454       // Store the daughter Candidates in the VertexCompositeCandidates 
00455       //    if they pass mass cuts
00456       if( doKshorts ) {
00457         theKshort->addDaughter(thePiPlusCand);
00458         theKshort->addDaughter(thePiMinusCand);
00459         theKshort->setPdgId(310);
00460         addp4.set( *theKshort );
00461         if( theKshort->mass() < kShortMass + kShortMassCut &&
00462             theKshort->mass() > kShortMass - kShortMassCut ) {
00463           theKshorts.push_back( *theKshort );
00464         }
00465       }
00466       
00467       if( doLambdas && theLambda ) {
00468         theLambda->addDaughter(theProtonCand);
00469         theLambda->addDaughter(thePiMinusCand);
00470         theLambda->setPdgId(3122);
00471         addp4.set( *theLambda );
00472         if( theLambda->mass() < lambdaMass + lambdaMassCut &&
00473             theLambda->mass() > lambdaMass - lambdaMassCut ) {
00474           theLambdas.push_back( *theLambda );
00475         }
00476       }
00477       else if ( doLambdas && theLambdaBar ) {
00478         theLambdaBar->addDaughter(theAntiProtonCand);
00479         theLambdaBar->addDaughter(thePiPlusCand);
00480         theLambdaBar->setPdgId(-3122);
00481         addp4.set( *theLambdaBar );
00482         if( theLambdaBar->mass() < lambdaMass + lambdaMassCut &&
00483             theLambdaBar->mass() > lambdaMass - lambdaMassCut ) {
00484           theLambdas.push_back( *theLambdaBar );
00485         }
00486       }
00487 
00488       if(theKshort) delete theKshort;
00489       if(theLambda) delete theLambda;
00490       if(theLambdaBar) delete theLambdaBar;
00491       theKshort = theLambda = theLambdaBar = 0;
00492 
00493     }
00494   }
00495 }
00496 
00497 // Get methods
00498 const reco::VertexCompositeCandidateCollection& V0Fitter::getKshorts() const {
00499   return theKshorts;
00500 }
00501 
00502 const reco::VertexCompositeCandidateCollection& V0Fitter::getLambdas() const {
00503   return theLambdas;
00504 }
00505 
00506 
00507 // Experimental
00508 double V0Fitter::findV0MassError(const GlobalPoint &vtxPos, const std::vector<reco::TransientTrack> &dauTracks) { 
00509   return -1.;
00510 }
00511 
00512 /*
00513 double V0Fitter::findV0MassError(const GlobalPoint &vtxPos, std::vector<reco::TransientTrack> dauTracks) {
00514   // Returns -99999. if trajectory states fail at vertex position
00515 
00516   // Load positive track trajectory at vertex into vector, then negative track
00517   std::vector<TrajectoryStateClosestToPoint> sortedTrajStatesAtVtx;
00518   for( unsigned int ndx = 0; ndx < dauTracks.size(); ndx++ ) {
00519     if( dauTracks[ndx].trajectoryStateClosestToPoint(vtxPos).isValid() ) {
00520       std::cout << "From TSCP: " 
00521                 << dauTracks[ndx].trajectoryStateClosestToPoint(vtxPos).perigeeParameters().transverseCurvature()
00522                 << "; From Track: " << dauTracks[ndx].track().qoverp() << std::endl;
00523     }
00524     if( sortedTrajStatesAtVtx.size() == 0 ) {
00525       if( dauTracks[ndx].charge() > 0 ) {
00526         sortedTrajStatesAtVtx.push_back( dauTracks[ndx].trajectoryStateClosestToPoint(vtxPos) );
00527       }
00528       else {
00529         sortedTrajStatesAtVtx.push_back( dauTracks[ndx].trajectoryStateClosestToPoint(vtxPos) );
00530       }
00531     }
00532   }
00533   std::vector<PerigeeTrajectoryParameters> param;
00534   std::vector<PerigeeTrajectoryError> paramError;
00535   std::vector<GlobalVector> momenta;
00536 
00537   for( unsigned int ndx2 = 0; ndx2 < sortedTrajStatesAtVtx.size(); ndx2++ ) {
00538     if( sortedTrajStatesAtVtx[ndx2].isValid() ) {
00539       param.push_back( sortedTrajStatesAtVtx[ndx2].perigeeParameters() );
00540       paramError.push_back( sortedTrajStatesAtVtx[ndx2].perigeeError() );
00541       momenta.push_back( sortedTrajStatesAtVtx[ndx2].momentum() );
00542     }
00543     else return -99999.;
00544   }
00545   return 0;
00546 }
00547 */
00548 
00549 
00550