CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_1/src/DataFormats/ParticleFlowReco/interface/PFDisplacedVertex.h

Go to the documentation of this file.
00001 #ifndef DataFormat_ParticleFlowReco_PFDisplacedVertex_h
00002 #define DataFormat_ParticleFlowReco_PFDisplacedVertex_h 
00003 
00004 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00005 #include "DataFormats/VertexReco/interface/Vertex.h"
00006 #include "DataFormats/Math/interface/LorentzVector.h"
00007 
00008 #include <vector>
00009 #include <string>
00010 #include <iostream>
00011 
00012 namespace reco {
00013 
00014   
00016 
00024   class PFDisplacedVertex : public Vertex{
00025 
00026   public:
00027 
00029     typedef std::pair <unsigned int, unsigned int> PFTrackHitInfo;
00030     typedef std::pair <PFTrackHitInfo, PFTrackHitInfo> PFTrackHitFullInfo;
00031 
00033     enum M_Hypo {
00034       M_CUSTOM,
00035       M_MASSLESS,
00036       M_PION,
00037       M_KAON,
00038       M_LAMBDA
00039     };
00040 
00041 
00046     enum VertexTrackType {
00047       T_NOT_FROM_VERTEX,
00048       T_TO_VERTEX,
00049       T_FROM_VERTEX,
00050       T_MERGED
00051     };
00052 
00055     enum VertexType {
00056       ANY = 0,
00057       FAKE = 1,
00058       LOOPER = 2,
00059       NUCL = 10,
00060       NUCL_LOOSE = 11,
00061       NUCL_KINK = 12,
00062       CONVERSION = 20,
00063       CONVERSION_LOOSE = 21,
00064       CONVERTED_BREMM = 22,
00065       K0_DECAY = 30,
00066       LAMBDA_DECAY = 31,
00067       LAMBDABAR_DECAY = 32,
00068       KPLUS_DECAY = 40,
00069       KMINUS_DECAY = 41,
00070       KPLUS_DECAY_LOOSE = 42,
00071       KMINUS_DECAY_LOOSE = 43,
00072       BSM_VERTEX = 100
00073     };
00074 
00075 
00077     PFDisplacedVertex();
00078 
00080     PFDisplacedVertex(reco::Vertex&);
00081 
00083     void addElement( const TrackBaseRef & r, const Track & refTrack, 
00084                      const PFTrackHitFullInfo& hitInfo , 
00085                      VertexTrackType trackType = T_NOT_FROM_VERTEX, float w=1.0 );
00086 
00088     void cleanTracks();
00089     
00091     void setVertexType(VertexType vertexType) {vertexType_ = vertexType;}
00092 
00095     void setPrimaryDirection(const math::XYZPoint& pvtx);
00096 
00098     VertexType vertexType(){return vertexType_;}
00099 
00100     const std::vector < PFTrackHitFullInfo > trackHitFullInfos() const
00101       {return trackHitFullInfos_;}
00102 
00103     const std::vector <VertexTrackType> trackTypes() const
00104       {return trackTypes_;}
00105 
00106 
00108 
00110     const bool isTherePrimaryTracks() const 
00111     {return isThereKindTracks(T_TO_VERTEX);}
00112 
00114     const bool isThereMergedTracks() const
00115     {return isThereKindTracks(T_MERGED);}
00116 
00118     const bool isThereSecondaryTracks() const
00119     {return isThereKindTracks(T_FROM_VERTEX);}
00120 
00122     const bool isThereNotFromVertexTracks() const
00123     {return isThereKindTracks(T_NOT_FROM_VERTEX);}
00124 
00125 
00126 
00127 
00129     const bool isPrimaryTrack(const reco::TrackBaseRef& originalTrack) const 
00130     {
00131       size_t itrk = trackPosition(originalTrack);
00132       return isTrack(itrk, T_TO_VERTEX);
00133     }
00134 
00136     const bool isSecondaryTrack(const reco::TrackBaseRef& originalTrack) const 
00137     {
00138       size_t itrk = trackPosition(originalTrack);
00139       return isTrack(itrk, T_FROM_VERTEX);
00140     }
00141 
00143     const bool isMergedTrack(const reco::TrackBaseRef& originalTrack) const 
00144     {
00145       size_t itrk = trackPosition(originalTrack);
00146       return isTrack(itrk, T_MERGED);
00147     }
00148 
00149     const PFTrackHitFullInfo trackHitFullInfo(const reco::TrackBaseRef& originalTrack) const{
00150       size_t itrk = trackPosition(originalTrack);
00151       return trackHitFullInfos_[itrk];
00152     }
00153 
00154 
00155 
00157     const bool isIncomingTrack(const reco::TrackBaseRef& originalTrack) const 
00158     {
00159       size_t itrk = trackPosition(originalTrack);
00160       return isTrack(itrk, T_MERGED) || isTrack(itrk, T_TO_VERTEX);
00161     }
00162  
00164     const bool isOutgoingTrack(const reco::TrackBaseRef& originalTrack) const 
00165     {
00166       size_t itrk = trackPosition(originalTrack);
00167       return isTrack(itrk, T_FROM_VERTEX);    
00168     }
00169 
00170 
00172     const int nPrimaryTracks() const
00173     {return nKindTracks(T_TO_VERTEX);}
00174 
00176     const int nMergedTracks() const
00177     {return nKindTracks(T_MERGED);}
00178 
00180     const int nSecondaryTracks() const
00181     {return nKindTracks(T_FROM_VERTEX);}
00182 
00184     const int nNotFromVertexTracks() const
00185     {return nKindTracks(T_NOT_FROM_VERTEX);}
00186 
00188     const int nTracks() const {return trackTypes_.size();}
00189 
00190     //    const reco::VertexTrackType vertexTrackType(reco::TrackBaseRef tkRef) const;
00191 
00195     const math::XYZTLorentzVector 
00196       secondaryMomentum(std::string massHypo = "PI", 
00197                         bool useRefitted = true, double mass = 0.0) const 
00198       {return momentum(massHypo, T_FROM_VERTEX, useRefitted, mass);}
00199 
00201     const math::XYZTLorentzVector 
00202       primaryMomentum(std::string massHypo = "PI", 
00203                       bool useRefitted = true, double mass = 0.0) const
00204       {return momentum(massHypo, T_TO_VERTEX, useRefitted, mass);}
00205 
00206 
00207 
00211     const math::XYZTLorentzVector 
00212       secondaryMomentum(M_Hypo massHypo, 
00213                         bool useRefitted = true, double mass = 0.0) const 
00214       {return momentum(massHypo, T_FROM_VERTEX, useRefitted, mass);}
00215 
00217     const math::XYZTLorentzVector 
00218       primaryMomentum(M_Hypo massHypo, 
00219                       bool useRefitted = true, double mass = 0.0) const
00220       {return momentum(massHypo, T_TO_VERTEX, useRefitted, mass);}
00221 
00222 
00223     void calcKinematics(){
00224       defaultPrimaryMomentum_ = momentum("PI", T_TO_VERTEX, false, 0.0);
00225       defaultSecondaryMomentum_ = momentum("PI", T_FROM_VERTEX, true, 0.0);
00226     }
00227 
00228 
00229     const double secondaryPt() const 
00230       {return defaultPrimaryMomentum_.Pt();}
00231 
00233     const double primaryPt() const
00234       {return defaultSecondaryMomentum_.Pt();}
00235 
00236 
00237 
00239     const int totalCharge() const;
00240 
00243     const double angle_io() const;
00244 
00245 
00247     const math::XYZVector primaryDirection() const;
00248 
00249 
00250 
00251     bool isFake() const { return vertexType_ == FAKE;}
00252     bool isLooper() const { return vertexType_ ==  LOOPER;}
00253     bool isNucl() const { return vertexType_ ==  NUCL;}
00254     bool isNucl_Loose() const { return vertexType_ ==  NUCL_LOOSE;}
00255     bool isNucl_Kink() const { return vertexType_ ==   NUCL_KINK;}
00256     bool isConv() const { return vertexType_ ==   CONVERSION;}
00257     bool isConv_Loose() const { return vertexType_ ==   CONVERSION_LOOSE;}
00258     bool isConvertedBremm() const { return vertexType_ ==   CONVERTED_BREMM;}
00259     bool isK0() const { return vertexType_ ==   K0_DECAY;}
00260     bool isLambda() const { return vertexType_ ==   LAMBDA_DECAY;}
00261     bool isLambdaBar() const { return vertexType_ ==   LAMBDABAR_DECAY;}
00262     bool isKplus() const { return vertexType_ ==   KPLUS_DECAY;}
00263     bool isKminus() const { return vertexType_ ==   KMINUS_DECAY;}
00264     bool isKplus_Loose() const { return vertexType_ ==   KPLUS_DECAY_LOOSE;}
00265     bool isKminus_Loose() const { return vertexType_ ==   KMINUS_DECAY_LOOSE;}
00266     bool isBSM() const { return vertexType_ ==   BSM_VERTEX;}
00267 
00268 
00269     std::string nameVertexType() const;
00270 
00272     void Dump(std::ostream& out = std::cout) const;
00273 
00274   private:
00275 
00277 
00279     const bool isThereKindTracks(VertexTrackType) const;
00280 
00282     const int nKindTracks(VertexTrackType) const;
00283 
00285     const  math::XYZTLorentzVector momentum(std::string, 
00286                                             VertexTrackType,
00287                                             bool, double mass) const;
00288 
00290     const  math::XYZTLorentzVector momentum(M_Hypo massHypo, 
00291                                             VertexTrackType,
00292                                             bool, double mass) const;
00293 
00295     const double getMass2(M_Hypo, double) const;
00296 
00297     const size_t trackPosition(const reco::TrackBaseRef& originalTrack) const;
00298 
00299     const bool isTrack(size_t itrk, VertexTrackType T) const {
00300       return  trackTypes_[itrk] == T;
00301     }
00302 
00303 
00305 
00307     VertexType vertexType_;
00308 
00310     std::vector < VertexTrackType > trackTypes_;
00311 
00313     std::vector < PFTrackHitFullInfo > trackHitFullInfos_;
00314 
00315     math::XYZVector primaryDirection_;
00316 
00317     math::XYZTLorentzVector defaultPrimaryMomentum_;
00318     math::XYZTLorentzVector defaultSecondaryMomentum_;
00319 
00320     
00321   };
00322 }
00323 
00324 #endif
00325 
00326 
00327