CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_1/src/RecoParticleFlow/PFTracking/src/PFDisplacedVertexHelper.cc

Go to the documentation of this file.
00001 #include "RecoParticleFlow/PFTracking/interface/PFDisplacedVertexHelper.h"
00002 
00003 #include "TVector3.h"
00004 #include "DataFormats/VertexReco/interface/Vertex.h"
00005 
00006 #include "TMath.h"
00007 
00008 using namespace std;
00009 using namespace reco;
00010 
00011 
00012 const double PFDisplacedVertexHelper::pion_mass2 = 0.0194;
00013 const double PFDisplacedVertexHelper::muon_mass2 = 0.106*0.106;
00014 const double PFDisplacedVertexHelper::proton_mass2 = 0.938*0.938;
00015 
00016 //for debug only 
00017 //#define PFLOW_DEBUG
00018 
00019 PFDisplacedVertexHelper::PFDisplacedVertexHelper() : 
00020   tracksSelector_(),
00021   vertexIdentifier_(),      
00022   pvtx_(math::XYZPoint(0,0,0)) {}
00023 
00024 PFDisplacedVertexHelper::~PFDisplacedVertexHelper() {}
00025 
00026 void PFDisplacedVertexHelper::setPrimaryVertex(
00027                                                edm::Handle< reco::VertexCollection > mainVertexHandle, 
00028                                                edm::Handle< reco::BeamSpot > beamSpotHandle){
00029 
00030   const math::XYZPoint beamSpot = beamSpotHandle.isValid() ? 
00031     math::XYZPoint(beamSpotHandle->x0(), beamSpotHandle->y0(), beamSpotHandle->z0()) : 
00032     math::XYZPoint(0, 0, 0);
00033 
00034   // The primary vertex is taken from the refitted list, 
00035   // if does not exist from the average offline beam spot position  
00036   // if does not exist (0,0,0) is used
00037   pvtx_ = mainVertexHandle.isValid() ? 
00038     math::XYZPoint(mainVertexHandle->begin()->x(), 
00039                    mainVertexHandle->begin()->y(), 
00040                    mainVertexHandle->begin()->z()) :
00041     beamSpot;
00042 }
00043 
00044 bool 
00045 PFDisplacedVertexHelper::isTrackSelected(const reco::Track& trk, 
00046                                          const reco::PFDisplacedVertex::VertexTrackType vertexTrackType) const {
00047     
00048 
00049   if (!tracksSelector_.selectTracks()) return true;
00050 
00051   bool isGoodTrack = false;
00052     
00053   bool isHighPurity = trk.quality( trk.qualityByName(tracksSelector_.quality().data()) );
00054         
00055   double nChi2  = trk.normalizedChi2(); 
00056   double pt = trk.pt();
00057   int nHits = trk.numberOfValidHits();
00058 
00059   bool bIsPrimary = 
00060     (
00061      (vertexTrackType == reco::PFDisplacedVertex::T_TO_VERTEX)
00062      ||
00063      (vertexTrackType == reco::PFDisplacedVertex::T_MERGED)
00064      );  
00065     
00066   
00067   if (bIsPrimary) {
00068     // Primary or merged track selection
00069     isGoodTrack = 
00070       ( ( nChi2 > tracksSelector_.nChi2_min()
00071           && nChi2 <  tracksSelector_.nChi2_max())
00072         || isHighPurity )
00073       && pt >  tracksSelector_.pt_min();
00074   } else {
00075     // Secondary tracks selection
00076     int nOuterHits = trk.trackerExpectedHitsOuter().numberOfHits();
00077 
00078     double dxy = trk.dxy(pvtx_);
00079       
00080     isGoodTrack =   
00081       nChi2 <  tracksSelector_.nChi2_max() 
00082       && pt >  tracksSelector_.pt_min() 
00083       && fabs(dxy) >  tracksSelector_.dxy_min() 
00084       && nHits >=  tracksSelector_.nHits_min() 
00085       && nOuterHits <=  tracksSelector_.nOuterHits_max();
00086       
00087   }
00088     
00089   return isGoodTrack;
00090       
00091 }
00092 
00093 reco::PFDisplacedVertex::VertexType 
00094 PFDisplacedVertexHelper::identifyVertex(const reco::PFDisplacedVertex& v) const{
00095 
00096   if (!vertexIdentifier_.identifyVertices()) return PFDisplacedVertex::ANY;
00097 
00098   PFDisplacedVertex::M_Hypo massElec = PFDisplacedVertex::M_MASSLESS;
00099   PFDisplacedVertex::M_Hypo massPion = PFDisplacedVertex::M_PION;
00100 
00101   math::XYZTLorentzVector mom_ee = v.secondaryMomentum(massElec, true);
00102   math::XYZTLorentzVector mom_pipi = v.secondaryMomentum(massPion, true);
00103 
00104   // ===== (1) Identify fake and looper vertices ===== //
00105 
00106   double ang = v.angle_io();
00107   double pt_ee = mom_ee.Pt();
00108   double eta_vtx = v.position().eta();
00109 
00110   //cout << "Angle = " << ang << endl;
00111 
00112   bool bDirectionFake = ang > vertexIdentifier_.angle_max();
00113   bool bLowPt = pt_ee < vertexIdentifier_.pt_min();
00114   bool bLooperEta = fabs(eta_vtx) < vertexIdentifier_.looper_eta_max();
00115 
00116   bool isFake = bDirectionFake ||(bLowPt && !bLooperEta); 
00117   bool isLooper = !bDirectionFake && bLowPt && bLooperEta;
00118 
00119   if (isFake) return  PFDisplacedVertex::FAKE;
00120   if (isLooper) return PFDisplacedVertex::LOOPER;
00121 
00122   // ===== (2) Identify Decays and Conversions ===== //
00123 
00124   int c1 = v.originalTrack(v.refittedTracks()[0])->charge();
00125   int c2 = v.originalTrack(v.refittedTracks()[1])->charge();
00126   double mass_ee = mom_ee.M();
00127 
00128   int nTracks = v.nTracks();
00129   int nSecondaryTracks = v.nSecondaryTracks();
00130   bool bPrimaryTracks = v.isTherePrimaryTracks();
00131   bool bMergedTracks = v.isThereMergedTracks();
00132 
00133   bool bPair = (nTracks == nSecondaryTracks) && (nTracks == 2);
00134   bool bOpposite = (c1*c2 < -0.1);
00135   bool bDirection = ang < vertexIdentifier_.angle_V0Conv_max();
00136   bool bConvMass = mass_ee < vertexIdentifier_.mConv_max();
00137 
00138   bool bV0Conv = bPair && bOpposite && bDirection;
00139 
00140   // If the basic configuration of conversions and V0 decays is respected
00141   // pair of secondary track with opposite charge and going in the right direction
00142   // the selection is then based on mass limits
00143   if (bV0Conv){
00144 
00145     // == (2.1) Identify Conversions == //
00146 
00147     bool isConv = bConvMass;
00148 
00149     if (isConv) return PFDisplacedVertex::CONVERSION_LOOSE;
00150  
00151     // == (2.2) Identify K0 == //
00152 
00153     double mass_pipi = mom_pipi.M();
00154 
00155     bool bK0Mass = 
00156       mass_pipi < vertexIdentifier_.mK0_max() 
00157       &&  mass_pipi > vertexIdentifier_.mK0_min();
00158 
00159     bool isK0 =  bK0Mass;
00160 
00161     if (isK0) return PFDisplacedVertex::K0_DECAY;
00162 
00163     // == (2.3) Identify Lambda == //
00164 
00165     int lambdaKind = lambdaCP(v); 
00166     
00167 
00168     bool isLambda = (lambdaKind == 1);
00169     bool isLambdaBar = (lambdaKind == -1);
00170 
00171     if (isLambda) return PFDisplacedVertex::LAMBDA_DECAY;
00172     if (isLambdaBar) return PFDisplacedVertex::LAMBDABAR_DECAY;
00173 
00174   }
00175 
00176   // == (2.4) Identify K- and K+ ==
00177   bool bK = 
00178     (nSecondaryTracks == 1) && bPrimaryTracks && !bMergedTracks
00179     && !bOpposite;
00180 
00181   if(bK){
00182 
00183     bool bKMass = isKaonMass(v);
00184  
00185     bool isKPlus = bKMass && c1 > 0;
00186     bool isKMinus = bKMass && c1 < 0;
00187 
00188     if(isKMinus) return PFDisplacedVertex::KMINUS_DECAY_LOOSE;
00189     if(isKPlus) return PFDisplacedVertex::KPLUS_DECAY_LOOSE;
00190 
00191   }
00192 
00193   // ===== (3) Identify Nuclears, Kinks and Remaining Fakes ===== //
00194 
00195   math::XYZTLorentzVector mom_prim = v.primaryMomentum(massPion, true);
00196   
00197   double p_prim = mom_prim.P();
00198   double p_sec = mom_pipi.P();
00199   double pt_prim = mom_prim.Pt();
00200   
00201   bool bLog = false;
00202   if (p_prim > 0.05 && p_sec > 0.05)  bLog = log10(p_prim/p_sec) > vertexIdentifier_.logPrimSec_min();
00203   bool bPtMin = pt_prim >  vertexIdentifier_.pt_kink_min();
00204 
00205   // A vertex with at least 3 tracks is considered as high purity nuclear interaction
00206   // the only exception is K- decay into 3 prongs. To be studied.
00207   bool isNuclearHighPurity = nTracks > 2 && mass_ee > vertexIdentifier_.mNucl_min(); 
00208   bool isFakeHighPurity = nTracks > 2 && mass_ee < vertexIdentifier_.mNucl_min();
00209   // Two secondary tracks with some minimal tracks angular opening are still accepted
00210   // as nuclear interactions
00211   bool isNuclearLowPurity = 
00212     (nTracks == nSecondaryTracks)
00213     && (nTracks == 2)
00214     && mass_ee > vertexIdentifier_.mNucl_min();
00215 
00216   bool isFakeNucl = 
00217     (nTracks == nSecondaryTracks)
00218     && (nTracks == 2)
00219     && mass_ee < vertexIdentifier_.mNucl_min();
00220 
00221   // Kinks: 1 primary + 1 secondary is accepted only if the primary tracks 
00222   // has more energy than the secondary and primary have some minimal pT 
00223   // to produce a nuclear interaction
00224   bool isNuclearKink = 
00225     (nSecondaryTracks == 1) && bPrimaryTracks && !bMergedTracks 
00226     && bLog && bPtMin;
00227 
00228   // Here some loopers may hide appearing in Particle Isolation plots. To be studied...
00229   bool isFakeKink = 
00230     ( (nSecondaryTracks == 1) && bMergedTracks && !bPrimaryTracks )
00231     ||
00232     ( (nSecondaryTracks == 1) && bPrimaryTracks && !bMergedTracks 
00233     && (!bLog || !bPtMin) );
00234 
00235   if (isNuclearHighPurity) return PFDisplacedVertex::NUCL;
00236   if (isNuclearLowPurity) return PFDisplacedVertex::NUCL_LOOSE;
00237   if (isFakeKink || isFakeNucl || isFakeHighPurity) return  PFDisplacedVertex::FAKE;
00238   if (isNuclearKink)  return  PFDisplacedVertex::NUCL_KINK;
00239   
00240 
00241   return PFDisplacedVertex::ANY;
00242 
00243 }
00244 
00245 
00246 int PFDisplacedVertexHelper::lambdaCP(const PFDisplacedVertex& v) const {
00247 
00248   int lambdaCP = 0;
00249 
00250   vector <Track> refittedTracks = v.refittedTracks();
00251 
00252   math::XYZTLorentzVector totalMomentumDcaRefit_lambda;
00253   math::XYZTLorentzVector totalMomentumDcaRefit_lambdabar;
00254 
00255 
00256   reco::Track tMomentumDcaRefit_0 = refittedTracks[0];
00257   reco::Track tMomentumDcaRefit_1 = refittedTracks[1];
00258 
00259   double mass2_0 = 0, mass2_1 = 0;
00260 
00261   int c1 = v.originalTrack(v.refittedTracks()[1])->charge();
00262 
00263   // --------------------------- lambda --------------------
00264 
00265 
00266   if (c1 > 0.1) mass2_0 = pion_mass2, mass2_1 = proton_mass2;
00267   else mass2_0 = proton_mass2, mass2_1 = pion_mass2;
00268       
00269   double E = sqrt(tMomentumDcaRefit_0.p()*tMomentumDcaRefit_0.p() + mass2_0);
00270 
00271   math::XYZTLorentzVector momentumDcaRefit_0(tMomentumDcaRefit_0.px(), tMomentumDcaRefit_0.py(), 
00272                                     tMomentumDcaRefit_0.pz(), E); 
00273 
00274 
00275 
00276   E = sqrt(tMomentumDcaRefit_1.p()*tMomentumDcaRefit_1.p() + mass2_1);
00277 
00278   math::XYZTLorentzVector momentumDcaRefit_1(tMomentumDcaRefit_1.px(), tMomentumDcaRefit_1.py(), 
00279                                     tMomentumDcaRefit_1.pz(), E); 
00280 
00281 
00282   totalMomentumDcaRefit_lambda = momentumDcaRefit_0 + momentumDcaRefit_1;
00283 
00284 
00285 
00286   // --------------------------- anti - lambda --------------------
00287 
00288   if (c1 > 0.1) mass2_1 = pion_mass2, mass2_0 = proton_mass2;
00289   else mass2_1 = proton_mass2, mass2_0 = pion_mass2;
00290   
00291       
00292   E = sqrt(tMomentumDcaRefit_0.p()*tMomentumDcaRefit_0.p() + mass2_0);
00293 
00294   math::XYZTLorentzVector momentumDcaRefit_01(tMomentumDcaRefit_0.px(), tMomentumDcaRefit_0.py(), 
00295                                     tMomentumDcaRefit_0.pz(), E); 
00296 
00297 
00298 
00299   E = sqrt(tMomentumDcaRefit_1.p()*tMomentumDcaRefit_1.p() + mass2_1);
00300 
00301   math::XYZTLorentzVector momentumDcaRefit_11(tMomentumDcaRefit_1.px(), tMomentumDcaRefit_1.py(), 
00302                                     tMomentumDcaRefit_1.pz(), E); 
00303 
00304 
00305   totalMomentumDcaRefit_lambdabar = momentumDcaRefit_01 + momentumDcaRefit_11;
00306 
00307   double mass_l    = totalMomentumDcaRefit_lambda.M();
00308   double mass_lbar = totalMomentumDcaRefit_lambdabar.M();
00309 
00310   //  cout << "mass_l = " << mass_l << " mass_lbar = " <<  mass_lbar << endl;
00311 
00312 
00313   if (
00314       mass_l < mass_lbar 
00315       && mass_l > vertexIdentifier_.mLambda_min()
00316       && mass_l < vertexIdentifier_.mLambda_max()) 
00317     lambdaCP = 1;
00318   else if (mass_lbar < mass_l 
00319            && mass_lbar > vertexIdentifier_.mLambda_min()
00320            && mass_lbar < vertexIdentifier_.mLambda_max()) 
00321     lambdaCP = -1;
00322   else 
00323     lambdaCP = 0;
00324 
00325   //cout << "lambdaCP = " << lambdaCP << endl;
00326 
00327   return lambdaCP;
00328 
00329 }
00330 
00331 
00332 
00333 bool PFDisplacedVertexHelper::isKaonMass(const PFDisplacedVertex& v) const {
00334 
00335   math::XYZVector  trkInit = v.refittedTracks()[1].momentum(), 
00336     trkFinal = v.refittedTracks()[0].momentum();
00337 
00338   if (v.trackTypes()[0] == PFDisplacedVertex::T_TO_VERTEX)
00339     trkInit = v.refittedTracks()[0].momentum(),
00340       trkFinal =  v.refittedTracks()[1].momentum();
00341 
00342 
00343     math::XYZVector trkNeutre(trkInit.x()-trkFinal.x(),  trkInit.y()-trkFinal.y(),
00344                           trkInit.z()-trkFinal.z());
00345 
00346     double Ec = sqrt(muon_mass2 + trkFinal.Mag2());
00347     double En = sqrt(0*0        + trkNeutre.Mag2());
00348 
00349 
00350 
00351     math::XYZTLorentzVectorD trkMuNu(trkInit.x(), trkInit.y(), trkInit.z(), Ec+En);
00352     double massMuNu = trkMuNu.M();
00353 
00354     bool bKMass = 
00355       massMuNu > vertexIdentifier_.mK_min() 
00356       && massMuNu < vertexIdentifier_.mK_max(); 
00357 
00358     return bKMass; 
00359 
00360 }
00361 
00362 
00363 
00364 void PFDisplacedVertexHelper::Dump(std::ostream& out) const {
00365     tracksSelector_.Dump();
00366     vertexIdentifier_.Dump();
00367     out << " pvtx_ = " << pvtx_ << std::endl;
00368   }