CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TemplatedSecondaryVertexTagInfo.h
Go to the documentation of this file.
1 #ifndef DataFormats_BTauReco_TemplatedSecondaryVertexTagInfo_h
2 #define DataFormats_BTauReco_TemplatedSecondaryVertexTagInfo_h
3 
4 #include <utility>
5 #include <vector>
6 
18 #include <functional>
19 #include <ext/functional>
20 #include <algorithm>
21 
25 
26 
27 
28 
29 namespace reco {
30 namespace btag{
31 inline float weight(const reco::TrackRef & t, const reco::Vertex &v) {return v.trackWeight(t);}
32 inline float weight(const reco::CandidatePtr & c, const reco::VertexCompositePtrCandidate &v) {return std::find(v.daughterPtrVector().begin(),v.daughterPtrVector().end(),c)!=v.daughterPtrVector().end();}
33 
34  struct TrackData {
35  enum Status {
39  };
40 
41  inline bool usedForVertexFit() const
42  { return (int)svStatus >= (int)trackUsedForVertexFit; }
43  inline bool associatedToVertex() const
44  { return (int)svStatus >= (int)trackAssociatedToVertex; }
45  inline bool associatedToVertex(unsigned int index) const
46  { return (int)svStatus == (int)trackAssociatedToVertex + (int)index; }
47 
49  };
50  typedef std::pair<unsigned int, TrackData> IndexedTrackData;
51 
52 }
53 template <class IPTI, class VTX>
55  public:
58 
59  struct VertexData {
60  VTX vertex;
63 
64  // Used by ROOT storage
66  };
67 
68  struct TrackFinder {
69  TrackFinder(const typename IPTI::input_container &tracks,
70  const typename IPTI::input_container::value_type &track) :
71  tracks(tracks), track(track) {}
72 
73  bool operator () (const IndexedTrackData &idt)
74  { return tracks[idt.first] == track; }
75 
76  const typename IPTI::input_container &tracks;
78  };
79 
81  bool operator () (const IndexedTrackData &idt)
82  { return idt.second.associatedToVertex(); }
83  };
84 
87  index(index) {}
88 
89  bool operator () (const IndexedTrackData &idt)
90  { return idt.second.associatedToVertex(index); }
91 
92  unsigned int index;
93  };
94 
95 
96 
97  typedef typename IPTI::input_container input_container;
98 
101 
103  const std::vector<IndexedTrackData> &trackData,
104  const std::vector<VertexData> &svData,
105  unsigned int vertexCandidates,
106  const edm::Ref<std::vector<IPTI> >&);
107 
109  virtual TemplatedSecondaryVertexTagInfo * clone(void) const {
110  return new TemplatedSecondaryVertexTagInfo(*this);
111  }
112 
114  { return m_trackIPTagInfoRef; }
115 
116  virtual edm::RefToBase<Jet> jet(void) const
117  { return m_trackIPTagInfoRef->jet(); }
118 
119 // virtual input_container ipTracks(void) const
120 // { return m_trackIPTagInfoRef->tracks(); }
121 
122 //AR TODO is it needed?
123 // const JetTracksAssociationRef &jtaRef(void) const
124 // { return m_trackIPTagInfoRef->jtaRef(); }
125 
126  const VTX &secondaryVertex(unsigned int index) const
127  { return m_svData[index].vertex; }
128 
129  unsigned int nSelectedTracks() const { return m_trackData.size(); }
130  unsigned int nVertexTracks() const;
131  unsigned int nVertexTracks(unsigned int index) const;
132  unsigned int nVertices() const { return m_svData.size(); }
133  unsigned int nVertexCandidates() const { return m_vertexCandidates; }
134 
137  input_container vertexTracks(unsigned int index) const;
138 
139  typename input_container::value_type track(unsigned int index) const;
140  unsigned int findTrack(const typename input_container::value_type &track) const;
141 
142  const TrackData &trackData(unsigned int index) const;
143  const TrackData &trackData(const typename input_container::value_type &track) const;
144 
145  const reco::btag::TrackIPData &trackIPData(unsigned int index) const;
147 
148  float trackWeight(unsigned int svIndex, unsigned int trackindex) const;
149  float trackWeight(unsigned int svIndex, const typename input_container::value_type &track) const;
150 
152  flightDistance(unsigned int index, bool in2d = false) const
153  { return in2d ? m_svData[index].dist2d : m_svData[index].dist3d; }
154  const GlobalVector &flightDirection(unsigned int index) const
155  { return m_svData[index].direction; }
156 
157  virtual TaggingVariableList taggingVariables() const;
158 
159  // Used by ROOT storage
161 
162  private:
164  std::vector<VertexData> m_svData;
165  unsigned int m_vertexCandidates;
166 
167  edm::Ref<std::vector<IPTI> > m_trackIPTagInfoRef;
168 };
169 
170 
172  const std::vector<IndexedTrackData> &trackData,
173  const std::vector<VertexData> &svData,
174  unsigned int vertexCandidates,
175  const edm::Ref<std::vector<IPTI> > & trackIPTagInfoRef) :
176  m_trackData(trackData),
177  m_svData(svData),
178  m_vertexCandidates(vertexCandidates),
179  m_trackIPTagInfoRef(trackIPTagInfoRef)
180 {
181 }
182 
183 template<class IPTI,class VTX> unsigned int TemplatedSecondaryVertexTagInfo<IPTI,VTX>::nVertexTracks() const
184 {
185  return std::count_if(m_trackData.begin(), m_trackData.end(),
187 }
188 
189 template<class IPTI,class VTX> unsigned int TemplatedSecondaryVertexTagInfo<IPTI,VTX>::nVertexTracks(unsigned int index) const
190 {
191  return std::count_if(m_trackData.begin(), m_trackData.end(),
193 }
194 
195 template<class IPTI,class VTX>
197 {
198  input_container trackRefs;
199  const input_container &trackIPTrackRefs =
200  m_trackIPTagInfoRef->selectedTracks();
201 
202  for(typename std::vector<typename reco::TemplatedSecondaryVertexTagInfo<IPTI,VTX>::IndexedTrackData>::const_iterator iter =
203  m_trackData.begin(); iter != m_trackData.end(); iter++)
204 
205  trackRefs.push_back(trackIPTrackRefs[iter->first]);
206 
207  return trackRefs;
208 }
209 
210 template<class IPTI,class VTX>
212 {
213  input_container trackRefs;
214  const input_container &trackIPTrackRefs =
215  m_trackIPTagInfoRef->selectedTracks();
216 
217  for(typename std::vector<typename reco::TemplatedSecondaryVertexTagInfo<IPTI,VTX>::IndexedTrackData>::const_iterator iter =
218  m_trackData.begin(); iter != m_trackData.end(); iter++)
219 
220  if (iter->second.associatedToVertex())
221  trackRefs.push_back(trackIPTrackRefs[iter->first]);
222 
223  return trackRefs;
224 }
225 
226 template<class IPTI,class VTX>
228 {
229  input_container trackRefs;
230  const input_container &trackIPTrackRefs =
231  m_trackIPTagInfoRef->selectedTracks();
232 
233  for(typename std::vector<typename reco::TemplatedSecondaryVertexTagInfo<IPTI,VTX>::IndexedTrackData>::const_iterator iter =
234  m_trackData.begin(); iter != m_trackData.end(); iter++)
235 
236  if (iter->second.associatedToVertex(index))
237  trackRefs.push_back(trackIPTrackRefs[iter->first]);
238 
239  return trackRefs;
240 }
241 
243 {
244  return m_trackIPTagInfoRef->selectedTracks()[m_trackData[index].first];
245 }
246 
247 template<class IPTI,class VTX> unsigned int TemplatedSecondaryVertexTagInfo<IPTI,VTX>::findTrack(const typename input_container::value_type &track) const
248 {
250  std::find_if(m_trackData.begin(), m_trackData.end(),
251  TrackFinder(m_trackIPTagInfoRef->selectedTracks(),
252  track));
253 
254  if (pos == m_trackData.end())
256  << "Track not found in "
257  " TemplatedSecondaryVertexTagInfo<IPTI,VTX>::findTrack." << std::endl;
258 
259  return pos - m_trackData.begin();
260 }
261 
262 template<class IPTI,class VTX>
264 {
265  return m_trackData[index].second;
266 }
267 
268 template<class IPTI,class VTX>
270 {
271  return m_trackData[findTrack(track)].second;
272 }
273 
274 template<class IPTI,class VTX>
276 {
277  return m_trackIPTagInfoRef->impactParameterData()[
278  m_trackData[index].first];
279 }
280 
281 template<class IPTI,class VTX>
283 {
284  return trackIPData(findTrack(track));
285 }
286 
287 template<class IPTI,class VTX> float TemplatedSecondaryVertexTagInfo<IPTI,VTX>::trackWeight(unsigned int svIndex,
288  const typename input_container::value_type &track) const
289 {
290  return reco::btag::weight(track,m_svData[svIndex].vertex);
291 }
292 
293 template<class IPTI,class VTX> float TemplatedSecondaryVertexTagInfo<IPTI,VTX>::trackWeight(unsigned int svIndex,
294  unsigned int trackIndex) const
295 {
296  return trackWeight(svIndex, track(trackIndex));
297 }
298 
300 {
301  TaggingVariableList vars;
302 
303  for(typename std::vector<typename TemplatedSecondaryVertexTagInfo<IPTI,VTX>::VertexData>::const_iterator iter = m_svData.begin();
304  iter != m_svData.end(); iter++) {
306  iter->dist2d.value(), true);
308  iter->dist2d.significance(), true);
310  iter->dist3d.value(), true);
312  iter->dist3d.significance(), true);
313 
315  Geom::deltaR(iter->direction, jet()->momentum()), true);
316  }
317 
320 
321  vars.finalize();
322  return vars;
323 }
324 
325 }
326 #endif // DataFormats_BTauReco_TemplatedSecondaryVertexTagInfo_h
Measurement1D flightDistance(unsigned int index, bool in2d=false) const
const VTX & secondaryVertex(unsigned int index) const
float weight(const reco::TrackRef &t, const reco::Vertex &v)
#define CMS_CLASS_VERSION(_version_)
Definition: classes.h:31
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
const edm::Ref< std::vector< IPTI > > & trackIPTagInfoRef() const
unsigned int findTrack(const typename input_container::value_type &track) const
float trackWeight(unsigned int svIndex, unsigned int trackindex) const
float trackWeight(const TrackBaseRef &r) const
returns the weight with which a Track has contributed to the vertex-fit.
const reco::btag::TrackIPData & trackIPData(unsigned int index) const
const TrackData & trackData(unsigned int index) const
#define private
Definition: FWEveView.cc:22
Container::value_type value_type
input_container::value_type track(unsigned int index) const
const GlobalVector & flightDirection(unsigned int index) const
tuple btag
Definition: BTagSF.py:15
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< float >, ROOT::Math::GlobalCoordinateSystemTag > GlobalVector
vector in glovbal coordinate system
Definition: Vector3D.h:27
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
TrackFinder(const typename IPTI::input_container &tracks, const typename IPTI::input_container::value_type &track)
string const
Definition: compareJSON.py:14
bool associatedToVertex(unsigned int index) const
virtual TemplatedSecondaryVertexTagInfo * clone(void) const
clone
const daughters & daughterPtrVector() const
references to daughtes
virtual TaggingVariableList taggingVariables() const
returns a description of the extended informations in a TaggingVariableList
std::pair< unsigned int, TrackData > IndexedTrackData
def template
Definition: svgfig.py:520
void insert(const TaggingVariable &variable, bool delayed=false)
virtual edm::RefToBase< Jet > jet(void) const
returns a polymorphic reference to the tagged jet