File indexing completed on 2024-04-06 12:03:47
0001 #ifndef DataFormats_BTauReco_TemplatedSecondaryVertexTagInfo_h
0002 #define DataFormats_BTauReco_TemplatedSecondaryVertexTagInfo_h
0003
0004 #include <utility>
0005 #include <vector>
0006
0007 #include "DataFormats/Common/interface/CMS_CLASS_VERSION.h"
0008 #include "DataFormats/GeometryCommonDetAlgo/interface/Measurement1D.h"
0009 #include "DataFormats/TrackReco/interface/TrackFwd.h"
0010 #include "DataFormats/VertexReco/interface/Vertex.h"
0011 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
0012 #include "DataFormats/GeometryVector/interface/GlobalVector.h"
0013 #include "DataFormats/JetReco/interface/JetTracksAssociation.h"
0014 #include "DataFormats/BTauReco/interface/JTATagInfo.h"
0015 #include "DataFormats/BTauReco/interface/TrackIPTagInfo.h"
0016 #include "DataFormats/BTauReco/interface/CandIPTagInfo.h"
0017 #include "DataFormats/Candidate/interface/VertexCompositePtrCandidate.h"
0018 #include <functional>
0019 #include <algorithm>
0020
0021 #include "FWCore/Utilities/interface/EDMException.h"
0022 #include "DataFormats/GeometryVector/interface/VectorUtil.h"
0023 #include "DataFormats/BTauReco/interface/TaggingVariable.h"
0024
0025 namespace reco {
0026 namespace btag {
0027 inline float weight(const reco::TrackRef &t, const reco::Vertex &v) { return v.trackWeight(t); }
0028 inline float weight(const reco::CandidatePtr &c, const reco::VertexCompositePtrCandidate &v) {
0029 return std::find(v.daughterPtrVector().begin(), v.daughterPtrVector().end(), c) != v.daughterPtrVector().end();
0030 }
0031
0032 struct TrackData {
0033 static constexpr int trackSelected = 0;
0034 static constexpr int trackUsedForVertexFit = 1;
0035 static constexpr int trackAssociatedToVertex = 2;
0036 inline bool usedForVertexFit() const { return svStatus >= trackUsedForVertexFit; }
0037 inline bool associatedToVertex() const { return svStatus >= trackAssociatedToVertex; }
0038 inline bool associatedToVertex(unsigned int index) const {
0039 return svStatus == trackAssociatedToVertex + (int)index;
0040 }
0041 int svStatus;
0042 };
0043 typedef std::pair<unsigned int, TrackData> IndexedTrackData;
0044
0045 }
0046 template <class IPTI, class VTX>
0047 class TemplatedSecondaryVertexTagInfo : public BaseTagInfo {
0048 public:
0049 typedef reco::btag::TrackData TrackData;
0050 typedef reco::btag::IndexedTrackData IndexedTrackData;
0051
0052 struct VertexData {
0053 VTX vertex;
0054 Measurement1D dist1d, dist2d, dist3d;
0055 GlobalVector direction;
0056
0057
0058 CMS_CLASS_VERSION(12)
0059 };
0060
0061 struct TrackFinder {
0062 TrackFinder(const typename IPTI::input_container &tracks, const typename IPTI::input_container::value_type &track)
0063 : tracks(tracks), track(track) {}
0064
0065 bool operator()(const IndexedTrackData &idt) { return tracks[idt.first] == track; }
0066
0067 const typename IPTI::input_container &tracks;
0068 const typename IPTI::input_container::value_type &track;
0069 };
0070
0071 struct VertexTrackSelector {
0072 bool operator()(const IndexedTrackData &idt) { return idt.second.associatedToVertex(); }
0073 };
0074
0075 struct IndexedVertexTrackSelector {
0076 IndexedVertexTrackSelector(unsigned int index) : index(index) {}
0077
0078 bool operator()(const IndexedTrackData &idt) { return idt.second.associatedToVertex(index); }
0079
0080 unsigned int index;
0081 };
0082
0083 typedef typename IPTI::input_container input_container;
0084
0085 TemplatedSecondaryVertexTagInfo() {}
0086 ~TemplatedSecondaryVertexTagInfo() override {}
0087
0088 TemplatedSecondaryVertexTagInfo(const std::vector<IndexedTrackData> &trackData,
0089 const std::vector<VertexData> &svData,
0090 unsigned int vertexCandidates,
0091 const edm::Ref<std::vector<IPTI> > &);
0092
0093
0094 TemplatedSecondaryVertexTagInfo *clone(void) const override { return new TemplatedSecondaryVertexTagInfo(*this); }
0095
0096 const edm::Ref<std::vector<IPTI> > &trackIPTagInfoRef() const { return m_trackIPTagInfoRef; }
0097
0098 edm::RefToBase<Jet> jet(void) const override { return m_trackIPTagInfoRef->jet(); }
0099
0100
0101
0102
0103
0104
0105
0106
0107 const VTX &secondaryVertex(unsigned int index) const { return m_svData[index].vertex; }
0108
0109 unsigned int nSelectedTracks() const { return m_trackData.size(); }
0110 unsigned int nVertexTracks() const;
0111 unsigned int nVertexTracks(unsigned int index) const;
0112 unsigned int nVertices() const { return m_svData.size(); }
0113 unsigned int nVertexCandidates() const { return m_vertexCandidates; }
0114
0115 input_container selectedTracks() const;
0116 input_container vertexTracks() const;
0117 input_container vertexTracks(unsigned int index) const;
0118
0119 typename input_container::value_type track(unsigned int index) const;
0120 unsigned int findTrack(const typename input_container::value_type &track) const;
0121
0122 const TrackData &trackData(unsigned int index) const;
0123 const TrackData &trackData(const typename input_container::value_type &track) const;
0124
0125 const reco::btag::TrackIPData &trackIPData(unsigned int index) const;
0126 const reco::btag::TrackIPData &trackIPData(const typename input_container::value_type &track) const;
0127
0128 float trackWeight(unsigned int svIndex, unsigned int trackindex) const;
0129 float trackWeight(unsigned int svIndex, const typename input_container::value_type &track) const;
0130
0131 Measurement1D flightDistance(unsigned int index, int dim = 0) const {
0132 if (dim == 1)
0133 return m_svData[index].dist1d;
0134 else if (dim == 2)
0135 return m_svData[index].dist2d;
0136 else
0137 return m_svData[index].dist3d;
0138 }
0139 const GlobalVector &flightDirection(unsigned int index) const { return m_svData[index].direction; }
0140 TaggingVariableList taggingVariables() const override;
0141
0142
0143 CMS_CLASS_VERSION(11)
0144
0145 private:
0146 std::vector<IndexedTrackData> m_trackData;
0147 std::vector<VertexData> m_svData;
0148 unsigned int m_vertexCandidates;
0149
0150 edm::Ref<std::vector<IPTI> > m_trackIPTagInfoRef;
0151 };
0152
0153 template <class IPTI, class VTX>
0154 TemplatedSecondaryVertexTagInfo<IPTI, VTX>::TemplatedSecondaryVertexTagInfo(
0155 const std::vector<IndexedTrackData> &trackData,
0156 const std::vector<VertexData> &svData,
0157 unsigned int vertexCandidates,
0158 const edm::Ref<std::vector<IPTI> > &trackIPTagInfoRef)
0159 : m_trackData(trackData),
0160 m_svData(svData),
0161 m_vertexCandidates(vertexCandidates),
0162 m_trackIPTagInfoRef(trackIPTagInfoRef) {}
0163
0164 template <class IPTI, class VTX>
0165 unsigned int TemplatedSecondaryVertexTagInfo<IPTI, VTX>::nVertexTracks() const {
0166 return std::count_if(m_trackData.begin(), m_trackData.end(), VertexTrackSelector());
0167 }
0168
0169 template <class IPTI, class VTX>
0170 unsigned int TemplatedSecondaryVertexTagInfo<IPTI, VTX>::nVertexTracks(unsigned int index) const {
0171 return std::count_if(m_trackData.begin(), m_trackData.end(), IndexedVertexTrackSelector(index));
0172 }
0173
0174 template <class IPTI, class VTX>
0175 typename reco::TemplatedSecondaryVertexTagInfo<IPTI, VTX>::input_container
0176 TemplatedSecondaryVertexTagInfo<IPTI, VTX>::selectedTracks() const {
0177 input_container trackRefs;
0178 const input_container &trackIPTrackRefs = m_trackIPTagInfoRef->selectedTracks();
0179
0180 for (typename std::vector<
0181 typename reco::TemplatedSecondaryVertexTagInfo<IPTI, VTX>::IndexedTrackData>::const_iterator iter =
0182 m_trackData.begin();
0183 iter != m_trackData.end();
0184 iter++)
0185
0186 trackRefs.push_back(trackIPTrackRefs[iter->first]);
0187
0188 return trackRefs;
0189 }
0190
0191 template <class IPTI, class VTX>
0192 typename TemplatedSecondaryVertexTagInfo<IPTI, VTX>::input_container
0193 TemplatedSecondaryVertexTagInfo<IPTI, VTX>::vertexTracks() const {
0194 input_container trackRefs;
0195 const input_container &trackIPTrackRefs = m_trackIPTagInfoRef->selectedTracks();
0196
0197 for (typename std::vector<
0198 typename reco::TemplatedSecondaryVertexTagInfo<IPTI, VTX>::IndexedTrackData>::const_iterator iter =
0199 m_trackData.begin();
0200 iter != m_trackData.end();
0201 iter++)
0202
0203 if (iter->second.associatedToVertex())
0204 trackRefs.push_back(trackIPTrackRefs[iter->first]);
0205
0206 return trackRefs;
0207 }
0208
0209 template <class IPTI, class VTX>
0210 typename TemplatedSecondaryVertexTagInfo<IPTI, VTX>::input_container
0211 TemplatedSecondaryVertexTagInfo<IPTI, VTX>::vertexTracks(unsigned int index) const {
0212 input_container trackRefs;
0213 const input_container &trackIPTrackRefs = m_trackIPTagInfoRef->selectedTracks();
0214
0215 for (typename std::vector<
0216 typename reco::TemplatedSecondaryVertexTagInfo<IPTI, VTX>::IndexedTrackData>::const_iterator iter =
0217 m_trackData.begin();
0218 iter != m_trackData.end();
0219 iter++)
0220
0221 if (iter->second.associatedToVertex(index))
0222 trackRefs.push_back(trackIPTrackRefs[iter->first]);
0223
0224 return trackRefs;
0225 }
0226
0227 template <class IPTI, class VTX>
0228 typename TemplatedSecondaryVertexTagInfo<IPTI, VTX>::input_container::value_type
0229 TemplatedSecondaryVertexTagInfo<IPTI, VTX>::track(unsigned int index) const {
0230 return m_trackIPTagInfoRef->selectedTracks()[m_trackData[index].first];
0231 }
0232
0233 template <class IPTI, class VTX>
0234 unsigned int TemplatedSecondaryVertexTagInfo<IPTI, VTX>::findTrack(
0235 const typename input_container::value_type &track) const {
0236 typename std::vector<typename reco::TemplatedSecondaryVertexTagInfo<IPTI, VTX>::IndexedTrackData>::const_iterator
0237 pos = std::find_if(
0238 m_trackData.begin(), m_trackData.end(), TrackFinder(m_trackIPTagInfoRef->selectedTracks(), track));
0239
0240 if (pos == m_trackData.end())
0241 throw edm::Exception(edm::errors::InvalidReference) << "Track not found in "
0242 " TemplatedSecondaryVertexTagInfo<IPTI,VTX>::findTrack."
0243 << std::endl;
0244
0245 return pos - m_trackData.begin();
0246 }
0247
0248 template <class IPTI, class VTX>
0249 const typename TemplatedSecondaryVertexTagInfo<IPTI, VTX>::TrackData &
0250 TemplatedSecondaryVertexTagInfo<IPTI, VTX>::trackData(unsigned int index) const {
0251 return m_trackData[index].second;
0252 }
0253
0254 template <class IPTI, class VTX>
0255 const typename TemplatedSecondaryVertexTagInfo<IPTI, VTX>::TrackData &
0256 TemplatedSecondaryVertexTagInfo<IPTI, VTX>::trackData(const typename input_container::value_type &track) const {
0257 return m_trackData[findTrack(track)].second;
0258 }
0259
0260 template <class IPTI, class VTX>
0261 const reco::btag::TrackIPData &TemplatedSecondaryVertexTagInfo<IPTI, VTX>::trackIPData(unsigned int index) const {
0262 return m_trackIPTagInfoRef->impactParameterData()[m_trackData[index].first];
0263 }
0264
0265 template <class IPTI, class VTX>
0266 const reco::btag::TrackIPData &TemplatedSecondaryVertexTagInfo<IPTI, VTX>::trackIPData(
0267 const typename input_container::value_type &track) const {
0268 return trackIPData(findTrack(track));
0269 }
0270
0271 template <class IPTI, class VTX>
0272 float TemplatedSecondaryVertexTagInfo<IPTI, VTX>::trackWeight(
0273 unsigned int svIndex, const typename input_container::value_type &track) const {
0274 return reco::btag::weight(track, m_svData[svIndex].vertex);
0275 }
0276
0277 template <class IPTI, class VTX>
0278 float TemplatedSecondaryVertexTagInfo<IPTI, VTX>::trackWeight(unsigned int svIndex, unsigned int trackIndex) const {
0279 return trackWeight(svIndex, track(trackIndex));
0280 }
0281
0282 template <class IPTI, class VTX>
0283 TaggingVariableList TemplatedSecondaryVertexTagInfo<IPTI, VTX>::taggingVariables() const {
0284 TaggingVariableList vars;
0285
0286 for (typename std::vector<typename TemplatedSecondaryVertexTagInfo<IPTI, VTX>::VertexData>::const_iterator iter =
0287 m_svData.begin();
0288 iter != m_svData.end();
0289 iter++) {
0290 vars.insert(btau::flightDistance1dVal, iter->dist1d.value(), true);
0291 vars.insert(btau::flightDistance1dSig, iter->dist1d.significance(), true);
0292 vars.insert(btau::flightDistance2dVal, iter->dist2d.value(), true);
0293 vars.insert(btau::flightDistance2dSig, iter->dist2d.significance(), true);
0294 vars.insert(btau::flightDistance3dVal, iter->dist3d.value(), true);
0295 vars.insert(btau::flightDistance3dSig, iter->dist3d.significance(), true);
0296
0297 vars.insert(btau::vertexJetDeltaR, Geom::deltaR(iter->direction, jet()->momentum()), true);
0298 }
0299
0300 vars.insert(btau::jetNSecondaryVertices, m_vertexCandidates, true);
0301 vars.insert(btau::vertexNTracks, nVertexTracks(), true);
0302
0303 vars.finalize();
0304 return vars;
0305 }
0306
0307 }
0308 #endif