Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:04:05

0001 #ifndef EgammaCandidates_Conversion_h
0002 #define EgammaCandidates_Conversion_h
0003 /** \class reco::Conversion Conversion.h DataFormats/EgammaCandidates/interface/Conversion.h
0004  *
0005  * 
0006  *
0007  * \author N.Marinelli  University of Notre Dame, US
0008  *
0009  *
0010  */
0011 
0012 #include <bitset>
0013 #include "DataFormats/Math/interface/Point3D.h"
0014 #include "DataFormats/TrackReco/interface/TrackFwd.h"
0015 #include "DataFormats/GeometryVector/interface/GlobalVector.h"
0016 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
0017 #include "DataFormats/GeometryCommonDetAlgo/interface/Measurement1DFloat.h"
0018 #include "DataFormats/VertexReco/interface/Vertex.h"
0019 #include "DataFormats/CaloRecHit/interface/CaloCluster.h"
0020 #include "DataFormats/CaloRecHit/interface/CaloClusterFwd.h"
0021 
0022 namespace reco {
0023   class Conversion {
0024   public:
0025     enum ConversionAlgorithm { undefined = 0, ecalSeeded = 1, trackerOnly = 2, mixed = 3, pflow = 4, algoSize = 5 };
0026 
0027     enum ConversionQuality {
0028       generalTracksOnly = 0,
0029       arbitratedEcalSeeded = 1,
0030       arbitratedMerged = 2,
0031       arbitratedMergedEcalGeneral = 3,
0032       gsfTracksOpenOnly = 4,
0033       highPurity = 8,
0034       highEfficiency = 9,
0035       ecalMatched1Track = 10,
0036       ecalMatched2Track = 11
0037     };
0038 
0039     static const std::string algoNames[];
0040 
0041     // Default constructor
0042     Conversion();
0043 
0044     Conversion(const reco::CaloClusterPtrVector& clu,
0045                const std::vector<edm::RefToBase<reco::Track> >& tr,
0046                const std::vector<math::XYZPointF>& trackPositionAtEcal,
0047                const reco::Vertex& convVtx,
0048                const std::vector<reco::CaloClusterPtr>& matchingBC,
0049                const float DCA,
0050                const std::vector<math::XYZPointF>& innPoint,
0051                const std::vector<math::XYZVectorF>& trackPin,
0052                const std::vector<math::XYZVectorF>& trackPout,
0053                const std::vector<uint8_t>& nHitsBeforeVtx,
0054                const std::vector<Measurement1DFloat>& dlClosestHitToVtx,
0055                uint8_t nSharedHits,
0056                const float mva,
0057                ConversionAlgorithm = undefined);
0058 
0059     Conversion(const reco::CaloClusterPtrVector& clu,
0060                const std::vector<reco::TrackRef>& tr,
0061                const std::vector<math::XYZPointF>& trackPositionAtEcal,
0062                const reco::Vertex& convVtx,
0063                const std::vector<reco::CaloClusterPtr>& matchingBC,
0064                const float DCA,
0065                const std::vector<math::XYZPointF>& innPoint,
0066                const std::vector<math::XYZVectorF>& trackPin,
0067                const std::vector<math::XYZVectorF>& trackPout,
0068                const float mva,
0069                ConversionAlgorithm = undefined);
0070 
0071     Conversion(const reco::CaloClusterPtrVector& clu,
0072                const std::vector<reco::TrackRef>& tr,
0073                const reco::Vertex& convVtx,
0074                ConversionAlgorithm = undefined);
0075 
0076     Conversion(const reco::CaloClusterPtrVector& clu,
0077                const std::vector<edm::RefToBase<reco::Track> >& tr,
0078                const reco::Vertex& convVtx,
0079                ConversionAlgorithm = undefined);
0080 
0081     /// returns a clone of the candidate
0082     Conversion* clone() const;
0083     /// Pointer to CaloCluster (foe Egamma Conversions it points to  a SuperCluster)
0084     reco::CaloClusterPtrVector caloCluster() const { return caloCluster_; }
0085     /// vector of track to base references
0086     std::vector<edm::RefToBase<reco::Track> > const& tracks() const;
0087     /// returns  the reco conversion vertex
0088     const reco::Vertex& conversionVertex() const { return theConversionVertex_; }
0089     /// Bool flagging objects having track size >0
0090     bool isConverted() const;
0091     /// Number of tracks= 0,1,2
0092     unsigned int nTracks() const { return tracks().size(); }
0093     /// get the value  of the TMVA output
0094     double MVAout() const { return theMVAout_; }
0095     /// get the MVS output from PF for one leg conversions
0096     std::vector<float> const oneLegMVA() { return theOneLegMVA_; }
0097     /// if nTracks=2 returns the pair invariant mass. Original tracks are used here
0098     double pairInvariantMass() const;
0099     /// Delta cot(Theta) where Theta is the angle in the (y,z) plane between the two tracks. Original tracks are used
0100     double pairCotThetaSeparation() const;
0101     /// Conversion tracks momentum from the tracks inner momentum
0102     math::XYZVectorF pairMomentum() const;
0103     /// Conversion track pair 4-momentum from the tracks refitted with vertex constraint
0104     math::XYZTLorentzVectorF refittedPair4Momentum() const;
0105     /// Conversion tracks momentum from the tracks refitted with vertex constraint
0106     math::XYZVectorF refittedPairMomentum() const;
0107     /// Super Cluster energy divided by track pair momentum if Standard seeding method. If a pointer to two (or more clusters)
0108     /// is stored in the conversion, this method returns the energy sum of clusters divided by the track pair momentum
0109     /// Track innermost momentum is used here
0110     double EoverP() const;
0111     /// Super Cluster energy divided by track pair momentum if Standard seeing method. If a pointer to two (or more clusters)
0112     /// is stored in the conversion, this method returns the energy sum of clusters divided by the track pair momentum
0113     ///  Track momentum refitted with vertex constraint is used
0114     double EoverPrefittedTracks() const;
0115     // Dist of minimum approach between tracks
0116     double distOfMinimumApproach() const { return theMinDistOfApproach_; }
0117     // deltaPhi tracks at innermost point
0118     double dPhiTracksAtVtx() const;
0119     // deltaPhi tracks at ECAl
0120     double dPhiTracksAtEcal() const;
0121     // deltaEta tracks at ECAl
0122     double dEtaTracksAtEcal() const;
0123 
0124     //impact parameter and decay length computed with respect to given beamspot or vertex
0125     //computed from refittedPairMomentum
0126 
0127     //transverse impact parameter
0128     double dxy(const math::XYZPoint& myBeamSpot = math::XYZPoint()) const;
0129     //longitudinal impact parameter
0130     double dz(const math::XYZPoint& myBeamSpot = math::XYZPoint()) const;
0131     //transverse decay length
0132     double lxy(const math::XYZPoint& myBeamSpot = math::XYZPoint()) const;
0133     //longitudinal decay length
0134     double lz(const math::XYZPoint& myBeamSpot = math::XYZPoint()) const;
0135     //z position of intersection with beamspot in rz plane (possible tilt of beamspot is neglected)
0136     double zOfPrimaryVertexFromTracks(const math::XYZPoint& myBeamSpot = math::XYZPoint()) const {
0137       return dz(myBeamSpot) + myBeamSpot.z();
0138     }
0139 
0140     ///// The following are variables provided per each track
0141     /// positions of the track extrapolation at the ECAL front face
0142     const std::vector<math::XYZPointF>& ecalImpactPosition() const { return thePositionAtEcal_; }
0143     //  pair of BC matching a posteriori the tracks
0144     const std::vector<reco::CaloClusterPtr>& bcMatchingWithTracks() const { return theMatchingBCs_; }
0145     /// signed transverse impact parameter for each track
0146     std::vector<double> tracksSigned_d0() const;
0147     /// Vector containing the position of the innermost hit of each track
0148     const std::vector<math::XYZPointF>& tracksInnerPosition() const { return theTrackInnerPosition_; }
0149     /// Vector of track momentum measured at the outermost hit
0150     const std::vector<math::XYZVectorF>& tracksPout() const { return theTrackPout_; }
0151     /// Vector of track momentum measured at the innermost hit
0152     const std::vector<math::XYZVectorF>& tracksPin() const { return theTrackPin_; }
0153     ///Vector of the number of hits before the vertex along each track trajector
0154     const std::vector<uint8_t>& nHitsBeforeVtx() const { return nHitsBeforeVtx_; }
0155     ///Vector of signed decay length with uncertainty from nearest hit on track to the conversion vtx positions
0156     const std::vector<Measurement1DFloat>& dlClosestHitToVtx() const { return dlClosestHitToVtx_; }
0157     ///number of shared hits btw the two track
0158     uint8_t nSharedHits() const { return nSharedHits_; }
0159 
0160     /// set the value  of the TMVA output
0161     void setMVAout(const float& mva) { theMVAout_ = mva; }
0162     /// set the MVS output from PF for one leg conversions
0163     void setOneLegMVA(const std::vector<float>& mva) { theOneLegMVA_ = mva; }
0164     // Set the ptr to the Super cluster if not set in the constructor
0165     void setMatchingSuperCluster(const reco::CaloClusterPtrVector& sc) { caloCluster_ = sc; }
0166     /// Conversion Track algorithm/provenance
0167     void setConversionAlgorithm(const ConversionAlgorithm a, bool set = true) {
0168       if (set)
0169         algorithm_ = a;
0170       else
0171         algorithm_ = undefined;
0172     }
0173     ConversionAlgorithm algo() const;
0174     std::string algoName() const;
0175     static std::string algoName(ConversionAlgorithm);
0176     static ConversionAlgorithm algoByName(const std::string& name);
0177 
0178     bool quality(ConversionQuality q) const { return (qualityMask_ & (1 << q)) >> q; }
0179     void setQuality(ConversionQuality q, bool b);
0180 
0181   private:
0182     /// vector pointer to a/multiple seed CaloCluster(s)
0183     reco::CaloClusterPtrVector caloCluster_;
0184     /// vector Track RefToBase
0185     std::vector<edm::RefToBase<reco::Track> > trackToBaseRefs_;
0186     /// position at the ECAl surface of the track extrapolation
0187     std::vector<math::XYZPointF> thePositionAtEcal_;
0188     /// Fitted Kalman conversion vertex
0189     reco::Vertex theConversionVertex_;
0190     /// Clusters mathing the tracks (these are not the seeds)
0191     std::vector<reco::CaloClusterPtr> theMatchingBCs_;
0192     /// P_in of tracks
0193     std::vector<math::XYZPointF> theTrackInnerPosition_;
0194     /// P_in of tracks
0195     std::vector<math::XYZVectorF> theTrackPin_;
0196     /// P_out of tracks
0197     std::vector<math::XYZVectorF> theTrackPout_;
0198     ///number of hits before the vertex on each trackerOnly
0199     std::vector<uint8_t> nHitsBeforeVtx_;
0200     ///signed decay length and uncertainty from nearest hit on track to conversion vertex
0201     std::vector<Measurement1DFloat> dlClosestHitToVtx_;
0202     /// vectors of TMVA outputs from pflow for one leg conversions
0203     std::vector<float> theOneLegMVA_;
0204     /// Distance of min approach of the two tracks
0205     float theMinDistOfApproach_;
0206     /// TMVA output
0207     float theMVAout_;
0208     uint16_t qualityMask_;
0209     ///number of shared hits between tracks
0210     uint8_t nSharedHits_;
0211     /// conversion algorithm/provenance
0212     uint8_t algorithm_;
0213   };
0214 
0215   inline Conversion::ConversionAlgorithm Conversion::algo() const { return (ConversionAlgorithm)algorithm_; }
0216 
0217   inline std::string Conversion::algoName() const {
0218     switch (algorithm_) {
0219       case undefined:
0220         return "undefined";
0221       case ecalSeeded:
0222         return "ecalSeeded";
0223       case trackerOnly:
0224         return "trackerOnly";
0225       case mixed:
0226         return "mixed";
0227       case pflow:
0228         return "pflow";
0229     }
0230     return "undefined";
0231   }
0232 
0233   inline std::string Conversion::algoName(ConversionAlgorithm a) {
0234     if (int(a) < int(algoSize) && int(a) > 0)
0235       return algoNames[int(a)];
0236     return "undefined";
0237   }
0238 
0239   inline void Conversion::setQuality(ConversionQuality q, bool b) {
0240     if (b)  //regular OR if setting value to true
0241       qualityMask_ |= (1 << q);
0242     else  // doing "half-XOR" if unsetting value
0243       qualityMask_ &= (~(1 << q));
0244   }
0245 
0246 }  // namespace reco
0247 
0248 #endif