Main Page
Namespaces
Classes
Package Documentation
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Properties
Friends
Macros
Pages
src
RecoTracker
DebugTools
interface
GetTrackTrajInfo.h
Go to the documentation of this file.
1
#ifndef RecoTracker_DebugTools_GetTrackTrajInfo_H
2
#define RecoTracker_DebugTools_GetTrackTrajInfo_H
3
4
/*
5
* Determine the track trajectory and detLayer at each layer that the track produces a hit in.
6
* This info can then be used to get the coordinates and momentum vector of the track at each of these
7
* layers etc.
8
*
9
* Call function analyze() for each track you are interested in. See comments below for that function.
10
* From the "result" that it returns, you can do things such as result.detTSOS.globalPosition(),
11
* to get the estimated position at which the track intercepts the layer.
12
*
13
* N.B. This information is obtained by extrapolating the track trajectory from its point of closest
14
* approach to the beam-line. It is therefore approximate, and should not be used for hit resolution
15
* studies.
16
* If you are using RECO, you can get more precise results by refitting the track instead of using this
17
* class. However, this class will work even on AOD.
18
*
19
* N.B. Your _cfg.py must load RecoTracker.Configuration.RecoTracker_cff to use this.
20
*
21
* Author: Ian Tomalin
22
* Date: Oct. 2011
23
*/
24
25
#include <memory>
26
#include "
FWCore/Framework/interface/Frameworkfwd.h
"
27
#include "
FWCore/Framework/interface/Event.h
"
28
#include "
FWCore/ParameterSet/interface/ParameterSet.h
"
29
#include "
FWCore/Framework/interface/ESHandle.h
"
30
31
#include "
FWCore/Framework/interface/EventSetup.h
"
32
#include "
DataFormats/TrackReco/interface/Track.h
"
33
#include "
TrackingTools/TransientTrack/interface/TransientTrackBuilder.h
"
34
#include "
TrackingTools/TrajectoryState/interface/TrajectoryStateOnSurface.h
"
35
36
#include <vector>
37
38
class
DetLayer
;
39
40
class
GetTrackTrajInfo
{
41
42
public
:
43
44
// Used to return results.
45
struct
Result
{
46
// If this is false, the information in the struct for this hit is invalid, as the track trajectory
47
// did not cross this layer. (Should only happen in very rare cases).
48
bool
valid
;
49
// If this is false, then although the track trajectory intercepted the layer, it did not intercept
50
// a sensor inside the layer. (Can happen rarely, if the track scattered, for example).
51
// If it is false, the detTSOS is evaluated at the intercept with the layer, not with the sensor,
52
// so will be slightly less accurate.
53
bool
accurate
;
54
// This is the DetLayer returned by GeometricSearchTracker.
55
// You can cast it into a specific type, such as BarrelDetLayer, before using it.
56
const
DetLayer
*
detLayer
;
57
// This is the track trajectory evaluated at the sensor. You can use it to get the coordinates
58
// where the track crosses the sensor and its momentum vector at that point.
59
TrajectoryStateOnSurface
detTSOS
;
60
};
61
62
GetTrackTrajInfo
() {}
63
64
~GetTrackTrajInfo
() {}
65
66
// For each hit on the track, return the information listed in the struct Result above.
67
// (See comments for struct Results).
68
// There is a one-to-one correspondence between this vector and the hits returned by track::hitPattern().
69
// i.e. They are ordered by the order in which the track crossed them.
70
std::vector<Result>
analyze
(
const
edm::EventSetup
& iSetup,
const
reco::Track
& track);
71
72
private
:
73
// Create map indicating r/z values of all layers/disks.
74
void
init
(
const
edm::EventSetup
& iSetup);
75
76
private
:
77
// Makes TransientTracks needed for vertex fitting.
78
edm::ESHandle<TransientTrackBuilder>
trkTool_
;
79
};
80
81
#endif
TransientTrackBuilder.h
GetTrackTrajInfo::GetTrackTrajInfo
GetTrackTrajInfo()
Definition:
GetTrackTrajInfo.h:62
GetTrackTrajInfo::init
void init(const edm::EventSetup &iSetup)
Event.h
EventSetup.h
TrajectoryStateOnSurface
Definition:
TrajectoryStateOnSurface.h:15
GetTrackTrajInfo
Definition:
GetTrackTrajInfo.h:40
Frameworkfwd.h
ParameterSet.h
GetTrackTrajInfo::~GetTrackTrajInfo
~GetTrackTrajInfo()
Definition:
GetTrackTrajInfo.h:64
GetTrackTrajInfo::Result::detLayer
const DetLayer * detLayer
Definition:
GetTrackTrajInfo.h:56
edm::ESHandle< TransientTrackBuilder >
GetTrackTrajInfo::Result::accurate
bool accurate
Definition:
GetTrackTrajInfo.h:53
ESHandle.h
GetTrackTrajInfo::trkTool_
edm::ESHandle< TransientTrackBuilder > trkTool_
Definition:
GetTrackTrajInfo.h:78
GetTrackTrajInfo::Result::detTSOS
TrajectoryStateOnSurface detTSOS
Definition:
GetTrackTrajInfo.h:59
edm::EventSetup
Definition:
EventSetup.h:44
DetLayer
Definition:
DetLayer.h:24
GetTrackTrajInfo::analyze
std::vector< Result > analyze(const edm::EventSetup &iSetup, const reco::Track &track)
Definition:
GetTrackTrajInfo.cc:28
reco::Track
Definition:
Track.h:26
TrajectoryStateOnSurface.h
GetTrackTrajInfo::Result
Definition:
GetTrackTrajInfo.h:45
GetTrackTrajInfo::Result::valid
bool valid
Definition:
GetTrackTrajInfo.h:48
Track.h
Generated for CMSSW Reference Manual by
1.8.5