File indexing completed on 2024-04-06 12:04:05
0001 #ifndef EgammaCandidates_Conversion_h
0002 #define EgammaCandidates_Conversion_h
0003
0004
0005
0006
0007
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
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
0082 Conversion* clone() const;
0083
0084 reco::CaloClusterPtrVector caloCluster() const { return caloCluster_; }
0085
0086 std::vector<edm::RefToBase<reco::Track> > const& tracks() const;
0087
0088 const reco::Vertex& conversionVertex() const { return theConversionVertex_; }
0089
0090 bool isConverted() const;
0091
0092 unsigned int nTracks() const { return tracks().size(); }
0093
0094 double MVAout() const { return theMVAout_; }
0095
0096 std::vector<float> const oneLegMVA() { return theOneLegMVA_; }
0097
0098 double pairInvariantMass() const;
0099
0100 double pairCotThetaSeparation() const;
0101
0102 math::XYZVectorF pairMomentum() const;
0103
0104 math::XYZTLorentzVectorF refittedPair4Momentum() const;
0105
0106 math::XYZVectorF refittedPairMomentum() const;
0107
0108
0109
0110 double EoverP() const;
0111
0112
0113
0114 double EoverPrefittedTracks() const;
0115
0116 double distOfMinimumApproach() const { return theMinDistOfApproach_; }
0117
0118 double dPhiTracksAtVtx() const;
0119
0120 double dPhiTracksAtEcal() const;
0121
0122 double dEtaTracksAtEcal() const;
0123
0124
0125
0126
0127
0128 double dxy(const math::XYZPoint& myBeamSpot = math::XYZPoint()) const;
0129
0130 double dz(const math::XYZPoint& myBeamSpot = math::XYZPoint()) const;
0131
0132 double lxy(const math::XYZPoint& myBeamSpot = math::XYZPoint()) const;
0133
0134 double lz(const math::XYZPoint& myBeamSpot = math::XYZPoint()) const;
0135
0136 double zOfPrimaryVertexFromTracks(const math::XYZPoint& myBeamSpot = math::XYZPoint()) const {
0137 return dz(myBeamSpot) + myBeamSpot.z();
0138 }
0139
0140
0141
0142 const std::vector<math::XYZPointF>& ecalImpactPosition() const { return thePositionAtEcal_; }
0143
0144 const std::vector<reco::CaloClusterPtr>& bcMatchingWithTracks() const { return theMatchingBCs_; }
0145
0146 std::vector<double> tracksSigned_d0() const;
0147
0148 const std::vector<math::XYZPointF>& tracksInnerPosition() const { return theTrackInnerPosition_; }
0149
0150 const std::vector<math::XYZVectorF>& tracksPout() const { return theTrackPout_; }
0151
0152 const std::vector<math::XYZVectorF>& tracksPin() const { return theTrackPin_; }
0153
0154 const std::vector<uint8_t>& nHitsBeforeVtx() const { return nHitsBeforeVtx_; }
0155
0156 const std::vector<Measurement1DFloat>& dlClosestHitToVtx() const { return dlClosestHitToVtx_; }
0157
0158 uint8_t nSharedHits() const { return nSharedHits_; }
0159
0160
0161 void setMVAout(const float& mva) { theMVAout_ = mva; }
0162
0163 void setOneLegMVA(const std::vector<float>& mva) { theOneLegMVA_ = mva; }
0164
0165 void setMatchingSuperCluster(const reco::CaloClusterPtrVector& sc) { caloCluster_ = sc; }
0166
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
0183 reco::CaloClusterPtrVector caloCluster_;
0184
0185 std::vector<edm::RefToBase<reco::Track> > trackToBaseRefs_;
0186
0187 std::vector<math::XYZPointF> thePositionAtEcal_;
0188
0189 reco::Vertex theConversionVertex_;
0190
0191 std::vector<reco::CaloClusterPtr> theMatchingBCs_;
0192
0193 std::vector<math::XYZPointF> theTrackInnerPosition_;
0194
0195 std::vector<math::XYZVectorF> theTrackPin_;
0196
0197 std::vector<math::XYZVectorF> theTrackPout_;
0198
0199 std::vector<uint8_t> nHitsBeforeVtx_;
0200
0201 std::vector<Measurement1DFloat> dlClosestHitToVtx_;
0202
0203 std::vector<float> theOneLegMVA_;
0204
0205 float theMinDistOfApproach_;
0206
0207 float theMVAout_;
0208 uint16_t qualityMask_;
0209
0210 uint8_t nSharedHits_;
0211
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)
0241 qualityMask_ |= (1 << q);
0242 else
0243 qualityMask_ &= (~(1 << q));
0244 }
0245
0246 }
0247
0248 #endif