8 #include "boost/foreach.hpp"
11 : tracks_token_( consumes<edm::
View<
reco::Track> >( iConfig.getParameter<edm::
InputTag>(
"Tracks") )),
13 Suffix ( iConfig.getParameter<std::
string>(
"Suffix") )
16 produces <std::vector<double> > (
Prefix +
"chi2" +
Suffix );
17 produces <std::vector<double> > (
Prefix +
"ndof" +
Suffix );
18 produces <std::vector<double> > (
Prefix +
"chi2ndof" +
Suffix );
19 produces <std::vector<float> > (
Prefix +
"charge" +
Suffix );
20 produces <std::vector<float> > (
Prefix +
"momentum" +
Suffix );
21 produces <std::vector<float> > (
Prefix +
"pt" +
Suffix );
22 produces <std::vector<float> > (
Prefix +
"pterr" +
Suffix );
23 produces <std::vector<unsigned int> > (
Prefix +
"hitsvalid" +
Suffix );
24 produces <std::vector<unsigned int> > (
Prefix +
"hitslost" +
Suffix );
25 produces <std::vector<double> > (
Prefix +
"theta" +
Suffix );
26 produces <std::vector<double> > (
Prefix +
"thetaerr" +
Suffix );
27 produces <std::vector<double> > (
Prefix +
"phi" +
Suffix );
28 produces <std::vector<double> > (
Prefix +
"phierr" +
Suffix );
29 produces <std::vector<double> > (
Prefix +
"eta" +
Suffix );
30 produces <std::vector<double> > (
Prefix +
"etaerr" +
Suffix );
31 produces <std::vector<double> > (
Prefix +
"dxy" +
Suffix );
32 produces <std::vector<double> > (
Prefix +
"dxyerr" +
Suffix );
33 produces <std::vector<double> > (
Prefix +
"dsz" +
Suffix );
34 produces <std::vector<double> > (
Prefix +
"dszerr" +
Suffix );
35 produces <std::vector<double> > (
Prefix +
"qoverp" +
Suffix );
36 produces <std::vector<double> > (
Prefix +
"qoverperr" +
Suffix );
37 produces <std::vector<double> > (
Prefix +
"vx" +
Suffix );
38 produces <std::vector<double> > (
Prefix +
"vy" +
Suffix );
39 produces <std::vector<double> > (
Prefix +
"vz" +
Suffix );
40 produces <std::vector<int> > (
Prefix +
"algo" +
Suffix );
45 std::auto_ptr<unsigned int>
number (
new unsigned int(0) );
46 std::auto_ptr<std::vector<double> >
chi2 (
new std::vector<double>() );
47 std::auto_ptr<std::vector<double> >
ndof (
new std::vector<double>() );
48 std::auto_ptr<std::vector<double> > chi2ndof (
new std::vector<double>() );
49 std::auto_ptr<std::vector<float> >
charge (
new std::vector<float>() );
50 std::auto_ptr<std::vector<float> > momentum (
new std::vector<float>() );
51 std::auto_ptr<std::vector<float> >
pt (
new std::vector<float>() );
52 std::auto_ptr<std::vector<float> > pterr (
new std::vector<float>() );
53 std::auto_ptr<std::vector<unsigned int> > hitsvalid (
new std::vector<unsigned int>() );
54 std::auto_ptr<std::vector<unsigned int> > hitslost (
new std::vector<unsigned int>() );
55 std::auto_ptr<std::vector<double> >
theta (
new std::vector<double>() );
56 std::auto_ptr<std::vector<double> > thetaerr (
new std::vector<double>() );
57 std::auto_ptr<std::vector<double> >
phi (
new std::vector<double>() );
58 std::auto_ptr<std::vector<double> > phierr (
new std::vector<double>() );
59 std::auto_ptr<std::vector<double> >
eta (
new std::vector<double>() );
60 std::auto_ptr<std::vector<double> > etaerr (
new std::vector<double>() );
61 std::auto_ptr<std::vector<double> > dxy (
new std::vector<double>() );
62 std::auto_ptr<std::vector<double> > dxyerr (
new std::vector<double>() );
63 std::auto_ptr<std::vector<double> > dsz (
new std::vector<double>() );
64 std::auto_ptr<std::vector<double> > dszerr (
new std::vector<double>() );
65 std::auto_ptr<std::vector<double> > qoverp (
new std::vector<double>() );
66 std::auto_ptr<std::vector<double> > qoverperr (
new std::vector<double>() );
67 std::auto_ptr<std::vector<double> > vx (
new std::vector<double>() );
68 std::auto_ptr<std::vector<double> > vy (
new std::vector<double>() );
69 std::auto_ptr<std::vector<double> > vz (
new std::vector<double>() );
70 std::auto_ptr<std::vector<int> >
algo (
new std::vector<int>() );
74 *number = tracks->size();
76 chi2->push_back( track.
chi2() );
77 ndof->push_back( track.
ndof() );
78 chi2ndof->push_back( track.
chi2()/track.
ndof() );
79 charge->push_back( track.
charge() );
80 momentum->push_back( track.
p() );
81 pt->push_back( track.
pt() );
82 pterr->push_back( track.
ptError() );
85 theta->push_back( track.
theta() );
87 phi->push_back( track.
phi() );
88 phierr->push_back( track.
phiError() );
89 eta->push_back( track.
eta() );
90 etaerr->push_back( track.
etaError() );
91 dxy->push_back( track.
dxy() );
92 dxyerr->push_back( track.
dxyError() );
93 dsz->push_back( track.
dsz() );
94 dszerr->push_back( track.
dszError() );
95 qoverp->push_back( track.
qoverp() );
97 vx->push_back( track.
vx() );
98 vy->push_back( track.
vy() );
99 vz->push_back( track.
vz() );
100 algo->push_back( (
int) track.
algo() );
double qoverp() const
q / p
double p() const
momentum vector magnitude
bool getByToken(EDGetToken token, Handle< PROD > &result) const
double theta() const
polar angle
double dxyError() const
error on dxy
Geom::Theta< T > theta() const
double etaError() const
error on eta
double phi() const
azimuthal angle of momentum vector
unsigned short numberOfLostHits() const
number of cases where track crossed a layer without getting a hit.
void produce(edm::Event &, const edm::EventSetup &)
double dsz() const
dsz parameter (THIS IS NOT the SZ impact parameter to (0,0,0) if refPoint is far from (0...
const edm::EDGetTokenT< edm::View< reco::Track > > tracks_token_
TrackAlgorithm algo() const
double dszError() const
error on dsz
double eta() const
pseudorapidity of momentum vector
double chi2() const
chi-squared of the fit
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
double ndof() const
number of degrees of freedom of the fit
double pt() const
track transverse momentum
double ptError() const
error on Pt (set to 1000 TeV if charge==0 for safety)
double phiError() const
error on phi
unsigned short numberOfValidHits() const
number of valid hits found
ShallowTracksProducer(const edm::ParameterSet &)
double qoverpError() const
error on signed transverse curvature
double vz() const
z coordinate of the reference point on track
double vy() const
y coordinate of the reference point on track
int charge() const
track electric charge
double dxy() const
dxy parameter. (This is the transverse impact parameter w.r.t. to (0,0,0) ONLY if refPoint is close t...
double vx() const
x coordinate of the reference point on track
double thetaError() const
error on theta