CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PixelTrackBuilder.cc
Go to the documentation of this file.
2 
3 
7 
12 
17 
18 #include <sstream>
19 using namespace std;
20 using namespace reco;
21 
22 
23 template <class T> T inline sqr( T t) {return t*t;}
24 
25 namespace {
27  const Measurement1D & pt,
28  const Measurement1D & phi,
29  const Measurement1D & cotTheta,
30  const Measurement1D & tip,
31  const Measurement1D & zip,
32  float chi2,
33  int charge)
34  {
35  ostringstream str;
36  str <<"\t pt: " << pt.value() <<"+/-"<<pt.error()
37  <<"\t phi: " << phi.value() <<"+/-"<<phi.error()
38  <<"\t cot: " << cotTheta.value() <<"+/-"<<cotTheta.error()
39  <<"\t tip: " << tip.value() <<"+/-"<<tip.error()
40  <<"\t zip: " << zip.value() <<"+/-"<<zip.error()
41  <<"\t charge: " << charge;
42  return str.str();
43  }
44 
45  std::string print(const reco::Track & track, const GlobalPoint & origin)
46  {
47 
48  math::XYZPoint bs(origin.x(), origin.y(), origin.z());
49 
50  Measurement1D phi( track.phi(), track.phiError());
51 
52  float theta = track.theta();
53  float cosTheta = cos(theta);
54  float sinTheta = sin(theta);
55  float errLambda2 = sqr( track.lambdaError() );
56  Measurement1D cotTheta(cosTheta/sinTheta, sqrt(errLambda2)/sqr(sinTheta));
57 
58  float pt_v = track.pt();
59  float errInvP2 = sqr(track.qoverpError());
60  float covIPtTheta = track.covariance(TrackBase::i_qoverp, TrackBase::i_lambda);
61  float errInvPt2 = ( errInvP2
62  + sqr(cosTheta/pt_v)*errLambda2
63  + 2*(cosTheta/pt_v)*covIPtTheta
64  ) / sqr(sinTheta);
65  Measurement1D pt(pt_v, sqr(pt_v)*sqrt(errInvPt2));
66 
67  Measurement1D tip(track.dxy(bs), track.d0Error());
68 
69  Measurement1D zip(track.dz(bs), track.dzError());
70 
71  return print(pt, phi, cotTheta, tip, zip, track.chi2(), track.charge());
72  }
73 
75  {
76  // TrajectoryStateOnSurface state(bstate);
77  float pt_v = state.globalMomentum().perp();
78  float phi_v = state.globalMomentum().phi();
79  float theta_v = state.globalMomentum().theta();
80 
81  CurvilinearTrajectoryError curv = state.curvilinearError();
82  float errPhi2 = curv.matrix()(3,3);
83  float errLambda2 = curv.matrix()(2,2);
84  float errInvP2 = curv.matrix()(1,1);
85  float covIPtTheta = curv.matrix()(1,2);
86  float cosTheta = cos(theta_v);
87  float sinTheta = sin(theta_v);
88  float errInvPt2 = ( errInvP2
89  + sqr(cosTheta/pt_v)*errLambda2
90  + 2*(cosTheta/pt_v)*covIPtTheta) / sqr(sinTheta);
91  float errCotTheta = sqrt(errLambda2)/sqr(sinTheta) ;
92  Measurement1D pt(pt_v, sqr(pt_v) * sqrt(errInvPt2));
93  Measurement1D phi(phi_v, sqrt(errPhi2) );
94  Measurement1D cotTheta(cosTheta/sinTheta, errCotTheta);
95 
96  float zip_v = state.globalPosition().z();
97  float zip_e = sqrt( state.localError().matrix()(4,4));
98  Measurement1D zip(zip_v, zip_e);
99 
100  float tip_v = state.localPosition().x();
101  int tip_sign = (state.localMomentum().y()*cotTheta.value() > 0) ? -1 : 1;
102  float tip_e = sqrt( state.localError().matrix()(3,3) );
103  Measurement1D tip( tip_sign*tip_v, tip_e);
104 
105  return print(pt, phi, cotTheta, tip, zip, 0., state.charge());
106  }
107 
108 
109  inline void checkState(const BasicTrajectoryStateOnSurface & bstate, const MagneticField* mf, const GlobalPoint & origin)
110  {
111  TrajectoryStateOnSurface state(bstate.clone());
112 
113  LogTrace("")<<" *** PixelTrackBuilder::checkState: ";
114  LogTrace("")<<"INPUT, ROTATION" << endl<<state.surface().rotation();
115  LogTrace("")<<"INPUT, TSOS:"<<endl<<state;
116 
118  TrajectoryStateOnSurface test= tipe.extrapolate(state, origin);
119  LogTrace("")<<"CHECK-1 ROTATION" << endl<<"\n"<<test.surface().rotation();
120  LogTrace("")<<"CHECK-1 TSOS" << endl<<test;
121 
122  TSCPBuilderNoMaterial tscpBuilder;
124  tscpBuilder(*(state.freeState()), origin);
125  FreeTrajectoryState fs = tscp.theState();
126  LogTrace("")<<"CHECK-2 FTS: " << fs;
127  }
128 
129 }
130 
131 
133  const Measurement1D & pt,
134  const Measurement1D & phi,
135  const Measurement1D & cotTheta,
136  const Measurement1D & tip,
137  const Measurement1D & zip,
138  float chi2,
139  int charge,
140  const std::vector<const TrackingRecHit* >& hits,
141  const MagneticField * mf,
142  const GlobalPoint & origin) const
143 {
144 
145  LogDebug("PixelTrackBuilder::build");
146  LogTrace("")<<"reconstructed TRIPLET kinematics:\n"<<print(pt,phi,cotTheta,tip,zip,chi2,charge);
147 
148  double sinTheta = 1/std::sqrt(1+sqr(cotTheta.value()));
149  double cosTheta = cotTheta.value()*sinTheta;
150  int tipSign = tip.value() > 0 ? 1 : -1;
151 
153  double invPtErr = 1./sqr(pt.value()) * pt.error();
154  m(0,0) = sqr(sinTheta) * (
155  sqr(invPtErr)
156  + sqr(cotTheta.error()/pt.value()*cosTheta * sinTheta)
157  );
158  m(0,2) = sqr( cotTheta.error()) * cosTheta * sqr(sinTheta) / pt.value();
159  m(1,1) = sqr( phi.error() );
160  m(2,2) = sqr( cotTheta.error());
161  m(3,3) = sqr( tip.error() );
162  m(4,4) = sqr( zip.error() );
164 
166  LocalPoint(tipSign*tip.value(), -tipSign*zip.value(), 0),
167  LocalVector(0., -tipSign*pt.value()*cotTheta.value(), pt.value()),
168  charge);
169 
170 
171  float sp = std::sin(phi.value());
172  float cp = std::cos(phi.value());
174  sp*tipSign, -cp*tipSign, 0,
175  0 , 0, -tipSign,
176  cp , sp , 0);
177 
178  // BTSOS hold Surface in a shared pointer and will be autodeleted when BTSOS goes out of scope...
179  // to avoid memory churn we allocate it locally and just avoid it be deleted by refcount...
180  Plane impPointPlane(origin, rot);
181  // (twice just to be sure!)
182  impPointPlane.addReference(); impPointPlane.addReference();
183  // use Base (to avoid a useless new)
184  BasicTrajectoryStateOnSurface impactPointState( lpar , error, impPointPlane, mf);
185 
186  //checkState(impactPointState,mf);
187  LogTrace("")<<"constructed TSOS :\n"<<print(impactPointState);
188 
189  int ndof = 2*hits.size()-5;
190  GlobalPoint vv = impactPointState.globalPosition();
191  math::XYZPoint pos( vv.x(), vv.y(), vv.z() );
192  GlobalVector pp = impactPointState.globalMomentum();
193  math::XYZVector mom( pp.x(), pp.y(), pp.z() );
194 
195  reco::Track * track = new reco::Track( chi2, ndof, pos, mom,
196  impactPointState.charge(), impactPointState.curvilinearError());
197 
198  LogTrace("") <<"RECONSTRUCTED TRACK (0,0,0):\n"<< print(*track,GlobalPoint(0,0,0))<<std::endl;
199  LogTrace("") <<"RECONSTRUCTED TRACK "<<origin<<"\n"<< print(*track,origin)<<std::endl;
200 
201  return track;
202 }
203 
#define LogDebug(id)
double d0Error() const
error on d0
Definition: TrackBase.h:209
Local3DVector LocalVector
Definition: LocalVector.h:12
tuple pp
Definition: createTree.py:15
const FreeTrajectoryState & theState() const
double theta() const
polar angle
Definition: TrackBase.h:115
std::string print(const Track &, edm::Verbosity=edm::Concise)
Track print utility.
Definition: print.cc:8
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
Geom::Theta< T > theta() const
T y() const
Definition: PV3DBase.h:63
double error() const
Definition: Measurement1D.h:30
double phi() const
azimuthal angle of momentum vector
Definition: TrackBase.h:137
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepSym< double, 5 > > AlgebraicSymMatrix55
double charge(const std::vector< uint8_t > &Ampls)
Definition: Plane.h:17
CovarianceMatrix covariance() const
return track covariance matrix
Definition: TrackBase.h:180
tuple zip
Definition: archive.py:476
const SurfaceType & surface() const
double chi2() const
chi-squared of the fit
Definition: TrackBase.h:105
T sqrt(T t)
Definition: SSEVec.h:48
double pt() const
track transverse momentum
Definition: TrackBase.h:129
T z() const
Definition: PV3DBase.h:64
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
reco::Track * build(const Measurement1D &pt, const Measurement1D &phi, const Measurement1D &cotTheta, const Measurement1D &tip, const Measurement1D &zip, float chi2, int charge, const std::vector< const TrackingRecHit * > &hits, const MagneticField *mf, const GlobalPoint &reference=GlobalPoint(0, 0, 0)) const
double phiError() const
error on phi
Definition: TrackBase.h:205
#define LogTrace(id)
double qoverpError() const
error on signed transverse curvature
Definition: TrackBase.h:190
double dz() const
dz parameter (= dsz/cos(lambda)). This is the track z0 w.r.t (0,0,0) only if the refPoint is close to...
Definition: TrackBase.h:125
double dzError() const
error on dz
Definition: TrackBase.h:213
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
BasicSingleTrajectoryState BasicTrajectoryStateOnSurface
double value() const
Definition: Measurement1D.h:28
void addReference() const
double lambdaError() const
error on lambda
Definition: TrackBase.h:201
const AlgebraicSymMatrix55 & matrix() const
Square< F >::type sqr(const F &f)
Definition: Square.h:13
Local3DPoint LocalPoint
Definition: LocalPoint.h:11
const RotationType & rotation() const
int charge() const
track electric charge
Definition: TrackBase.h:111
double dxy() const
dxy parameter. (This is the transverse impact parameter w.r.t. to (0,0,0) ONLY if refPoint is close t...
Definition: TrackBase.h:119
long double T
T x() const
Definition: PV3DBase.h:62
Definition: DDAxes.h:10