CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_0/src/DataFormats/EgammaCandidates/interface/Conversion.h

Go to the documentation of this file.
00001 #ifndef EgammaCandidates_Conversion_h
00002 #define EgammaCandidates_Conversion_h
00003 
00013 #include <bitset>
00014 #include "DataFormats/Math/interface/Point3D.h"
00015 #include "DataFormats/TrackReco/interface/TrackFwd.h" 
00016 #include "DataFormats/GeometryVector/interface/GlobalVector.h"
00017 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
00018 #include "DataFormats/GeometryCommonDetAlgo/interface/Measurement1DFloat.h"
00019 #include "DataFormats/VertexReco/interface/Vertex.h"
00020 #include "DataFormats/CaloRecHit/interface/CaloCluster.h"
00021 #include "DataFormats/CaloRecHit/interface/CaloClusterFwd.h"
00022 
00023 
00024 namespace reco {
00025     class Conversion  {
00026   public:
00027 
00028       enum ConversionAlgorithm {undefined=0, 
00029                                 ecalSeeded=1, 
00030                                 trackerOnly=2, 
00031                                 mixed=3, 
00032                                 pflow=4,
00033                                 algoSize=5}; 
00034 
00035       enum ConversionQuality {generalTracksOnly=0, 
00036                               arbitratedEcalSeeded=1, 
00037                               arbitratedMerged=2,
00038                               arbitratedMergedEcalGeneral=3,
00039                               highPurity=8, 
00040                               highEfficiency=9,
00041                               ecalMatched1Track=10,
00042                               ecalMatched2Track=11};
00043 
00044       static const std::string algoNames[];      
00045 
00046       // Default constructor
00047       Conversion();
00048       
00049       Conversion( const reco::CaloClusterPtrVector clu, 
00050                   const std::vector<edm::RefToBase<reco::Track> > tr,
00051                   const std::vector<math::XYZPointF> trackPositionAtEcal , 
00052                   const reco::Vertex  &  convVtx,
00053                   const std::vector<reco::CaloClusterPtr> & matchingBC,
00054                   const float DCA,        
00055                   const std::vector<math::XYZPointF> & innPoint,
00056                   const std::vector<math::XYZVectorF> & trackPin ,
00057                   const std::vector<math::XYZVectorF> & trackPout,
00058                   const std::vector<uint8_t> nHitsBeforeVtx,
00059                   const std::vector<Measurement1DFloat> & dlClosestHitToVtx,
00060                   uint8_t nSharedHits,
00061                   const float mva,
00062                   ConversionAlgorithm=undefined);
00063 
00064 
00065       Conversion( const reco::CaloClusterPtrVector clu, 
00066                   const std::vector<reco::TrackRef> tr,
00067                   const std::vector<math::XYZPointF> trackPositionAtEcal , 
00068                   const reco::Vertex  &  convVtx,
00069                   const std::vector<reco::CaloClusterPtr> & matchingBC,
00070                   const float DCA,        
00071                   const std::vector<math::XYZPointF> & innPoint,
00072                   const std::vector<math::XYZVectorF> & trackPin ,
00073                   const std::vector<math::XYZVectorF> & trackPout,
00074                   const float mva,
00075                   ConversionAlgorithm=undefined);
00076       
00077 
00078 
00079       
00080       Conversion( const reco::CaloClusterPtrVector clu, 
00081                   const std::vector<reco::TrackRef> tr,
00082                   const reco::Vertex  &  convVtx,
00083                   ConversionAlgorithm=undefined);
00084       
00085       Conversion( const reco::CaloClusterPtrVector clu, 
00086                   const std::vector<edm::RefToBase<reco::Track> > tr,
00087                   const reco::Vertex  &  convVtx,
00088                   ConversionAlgorithm=undefined);
00089       
00090            
00091       
00093       virtual ~Conversion();
00095       Conversion * clone() const;
00097       reco::CaloClusterPtrVector caloCluster() const {return caloCluster_ ;}
00099       std::vector<edm::RefToBase<reco::Track> > tracks() const ; 
00101       const reco::Vertex & conversionVertex() const  { return theConversionVertex_ ; }
00103       bool isConverted() const;
00105       unsigned int nTracks() const {return  tracks().size(); }
00107       double MVAout() const { return theMVAout_;}
00109       std::vector<float> const oneLegMVA() {return theOneLegMVA_;}
00111       double pairInvariantMass() const;
00113       double pairCotThetaSeparation() const;
00115       math::XYZVectorF  pairMomentum() const;
00117       math::XYZTLorentzVectorF   refittedPair4Momentum() const;
00119       math::XYZVectorF  refittedPairMomentum() const;
00123       double EoverP() const;
00127       double EoverPrefittedTracks() const;
00128       // Dist of minimum approach between tracks
00129       double distOfMinimumApproach() const {return  theMinDistOfApproach_;}
00130       // deltaPhi tracks at innermost point
00131       double dPhiTracksAtVtx() const;
00132       // deltaPhi tracks at ECAl
00133       double dPhiTracksAtEcal() const;
00134       // deltaEta tracks at ECAl
00135       double dEtaTracksAtEcal() const;
00136 
00137       //impact parameter and decay length computed with respect to given beamspot or vertex
00138       //computed from refittedPairMomentum
00139 
00140       //transverse impact parameter
00141       double dxy(const math::XYZPoint& myBeamSpot = math::XYZPoint()) const;
00142       //longitudinal impact parameter
00143       double dz(const math::XYZPoint& myBeamSpot = math::XYZPoint()) const;
00144       //transverse decay length
00145       double lxy(const math::XYZPoint& myBeamSpot = math::XYZPoint()) const;
00146       //longitudinal decay length
00147       double lz(const math::XYZPoint& myBeamSpot = math::XYZPoint()) const;
00148       //z position of intersection with beamspot in rz plane (possible tilt of beamspot is neglected)
00149       double zOfPrimaryVertexFromTracks(const math::XYZPoint& myBeamSpot = math::XYZPoint()) const { return dz(myBeamSpot) + myBeamSpot.z(); }
00150 
00153       const std::vector<math::XYZPointF> & ecalImpactPosition() const  {return thePositionAtEcal_;}
00154       //  pair of BC matching a posteriori the tracks
00155       const std::vector<reco::CaloClusterPtr>&  bcMatchingWithTracks() const { return theMatchingBCs_;}
00157       std::vector<double> tracksSigned_d0() const ;
00159       const std::vector<math::XYZPointF>& tracksInnerPosition() const {return theTrackInnerPosition_;}
00161       const std::vector<math::XYZVectorF>& tracksPout() const {return theTrackPout_;}
00163       const std::vector<math::XYZVectorF>& tracksPin() const  {return theTrackPin_;}
00165       const std::vector<uint8_t> &nHitsBeforeVtx() const { return nHitsBeforeVtx_; }
00167       const std::vector<Measurement1DFloat> &dlClosestHitToVtx() const { return dlClosestHitToVtx_; }
00169       uint8_t nSharedHits() const { return nSharedHits_; }
00170 
00171 
00173       void setMVAout(const float& mva) { theMVAout_=mva;}
00175       void setOneLegMVA(const std::vector<float>& mva) { theOneLegMVA_=mva;}
00176       // Set the ptr to the Super cluster if not set in the constructor 
00177       void setMatchingSuperCluster ( const  reco::CaloClusterPtrVector& sc) { caloCluster_= sc;}
00179       void setConversionAlgorithm(const ConversionAlgorithm a, bool set=true) { if (set) algorithm_=a; else algorithm_=undefined;}
00180       ConversionAlgorithm algo() const ;
00181       std::string algoName() const;
00182       static std::string algoName(ConversionAlgorithm );
00183       static ConversionAlgorithm  algoByName(const std::string &name);      
00184 
00185       bool quality(ConversionQuality q) const { return  (qualityMask_ & (1<<q))>>q; }
00186       void setQuality(ConversionQuality q, bool b);
00187 
00188 
00189       
00190     private:
00191       
00193       reco::CaloClusterPtrVector caloCluster_;
00195       std::vector<reco::TrackRef>  tracks_;
00197       mutable std::vector<edm::RefToBase<reco::Track> >  trackToBaseRefs_;
00199       std::vector<math::XYZPointF>  thePositionAtEcal_;
00201       reco::Vertex theConversionVertex_;
00203       std::vector<reco::CaloClusterPtr> theMatchingBCs_;
00205       float theMinDistOfApproach_;
00207       std::vector<math::XYZPointF> theTrackInnerPosition_;    
00209       std::vector<math::XYZVectorF> theTrackPin_;    
00211       std::vector<math::XYZVectorF> theTrackPout_;    
00213       std::vector<uint8_t> nHitsBeforeVtx_;
00215       std::vector<Measurement1DFloat> dlClosestHitToVtx_;
00217       uint8_t nSharedHits_;
00219       float theMVAout_;
00221       std::vector<float>  theOneLegMVA_;
00223       uint8_t algorithm_;
00224       uint16_t qualityMask_;
00225 
00226 
00227   };
00228 
00229     inline Conversion::ConversionAlgorithm Conversion::algo() const {
00230       return (ConversionAlgorithm) algorithm_;
00231     }
00232     
00233     
00234     inline std::string Conversion::algoName() const{
00235             
00236       switch(algorithm_)
00237         {
00238         case undefined: return "undefined";
00239         case ecalSeeded: return "ecalSeeded";
00240         case trackerOnly: return "trackerOnly";
00241         case mixed: return "mixed";
00242         case pflow: return "pflow";
00243 
00244         }
00245       return "undefined";
00246     }
00247 
00248     inline std::string Conversion::algoName(ConversionAlgorithm a){
00249       if(int(a) < int(algoSize) && int(a)>0) return algoNames[int(a)];
00250       return "undefined";
00251     }
00252 
00253     inline void Conversion::setQuality(ConversionQuality q, bool b){
00254       if (b)//regular OR if setting value to true
00255         qualityMask_ |= (1<<q) ;
00256       else // doing "half-XOR" if unsetting value
00257         qualityMask_ &= (~(1<<q));
00258 
00259     }
00260   
00261 }
00262 
00263 #endif