CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TrackerValidationVariables.cc
Go to the documentation of this file.
1 
11 
17 
20 
28 
30 
40 
42 
44 
45 #include "TMath.h"
46 
47 #include <string>
48 
50 {
51 
52 }
53 
56 {
57  trajCollectionToken_ = iC.consumes<std::vector<Trajectory> >(edm::InputTag(config.getParameter<std::string>("trajectoryInput")));
59 }
60 
62 {
63 
64 }
65 
66 void
67 TrackerValidationVariables::fillHitQuantities(const Trajectory* trajectory, std::vector<AVHitStruct> & v_avhitout)
68 {
69  TrajectoryStateCombiner tsoscomb;
70 
71  const std::vector<TrajectoryMeasurement> &tmColl = trajectory->measurements();
72  for(std::vector<TrajectoryMeasurement>::const_iterator itTraj = tmColl.begin();
73  itTraj != tmColl.end();
74  ++itTraj) {
75 
76  if (!itTraj->updatedState().isValid()) continue;
77 
78  TrajectoryStateOnSurface tsos = tsoscomb( itTraj->forwardPredictedState(), itTraj->backwardPredictedState() );
79  if(!tsos.isValid()) continue;
81 
82  if(!hit->isValid() || hit->geographicalId().det() != DetId::Tracker) continue;
83 
84  AVHitStruct hitStruct;
85  const DetId& hit_detId = hit->geographicalId();
86  unsigned int IntRawDetID = (hit_detId.rawId());
87  unsigned int IntSubDetID = (hit_detId.subdetId());
88 
89  if(IntSubDetID == 0) continue;
90 
91  //first calculate residuals in cartesian coordinates in the local module coordinate system
92 
93  LocalPoint lPHit = hit->localPosition();
94  LocalPoint lPTrk = tsos.localPosition();
95  LocalVector lVTrk = tsos.localDirection();
96 
97  hitStruct.localAlpha = atan2(lVTrk.x(), lVTrk.z()); // wrt. normal tg(alpha)=x/z
98  hitStruct.localBeta = atan2(lVTrk.y(), lVTrk.z()); // wrt. normal tg(beta)= y/z
99 
100  LocalError errHit = hit->localPositionError();
101  // no need to add APE to hitError anymore
102  // AlgebraicROOTObject<2>::SymMatrix mat = asSMatrix<2>(hit->parametersError());
103  // LocalError errHit = LocalError( mat(0,0),mat(0,1),mat(1,1) );
104  LocalError errTrk = tsos.localError().positionError();
105 
106  //check for negative error values: track error can have negative value, if matrix inversion fails (very rare case)
107  //hit error should always give positive values
108  if(errHit.xx()<0. || errHit.yy()<0. || errTrk.xx()<0. || errTrk.yy()<0.){
109  edm::LogError("TrackerValidationVariables") << "@SUB=TrackerValidationVariables::fillHitQuantities"
110  << "One of the squared error methods gives negative result"
111  << "\n\terrHit.xx()\terrHit.yy()\terrTrk.xx()\terrTrk.yy()"
112  << "\n\t" << errHit.xx()
113  << "\t" << errHit.yy()
114  << "\t" << errTrk.xx()
115  << "\t" << errTrk.yy();
116  continue;
117  }
118 
119  align::LocalVector res = lPTrk - lPHit;
120 
121  float resXErr = std::sqrt( errHit.xx() + errTrk.xx() );
122  float resYErr = std::sqrt( errHit.yy() + errTrk.yy() );
123 
124  hitStruct.resX = res.x();
125  hitStruct.resY = res.y();
126  hitStruct.resErrX = resXErr;
127  hitStruct.resErrY = resYErr;
128 
129  // hitStruct.localX = lPhit.x();
130  // hitStruct.localY = lPhit.y();
131  // EM: use predictions for local coordinates
132  hitStruct.localX = lPTrk.x();
133  hitStruct.localY = lPTrk.y();
134 
135  // now calculate residuals taking global orientation of modules and radial topology in TID/TEC into account
136  float resXprime(999.F), resYprime(999.F);
137  float resXatTrkY(999.F);
138  float resXprimeErr(999.F), resYprimeErr(999.F);
139 
140  if(hit->detUnit()){ // is it a single physical module?
141  const GeomDetUnit& detUnit = *(hit->detUnit());
142  float uOrientation(-999.F), vOrientation(-999.F);
143  float resXTopol(999.F), resYTopol(999.F);
144  float resXatTrkYTopol(999.F);
145 
146  const Surface& surface = hit->detUnit()->surface();
147  const BoundPlane& boundplane = hit->detUnit()->surface();
148  const Bounds& bound = boundplane.bounds();
149 
150  float length = 0;
151  float width = 0;
152 
153  LocalPoint lPModule(0.,0.,0.), lUDirection(1.,0.,0.), lVDirection(0.,1.,0.);
154  GlobalPoint gPModule = surface.toGlobal(lPModule),
155  gUDirection = surface.toGlobal(lUDirection),
156  gVDirection = surface.toGlobal(lVDirection);
157 
158  if (IntSubDetID == PixelSubdetector::PixelBarrel ||
159  IntSubDetID == StripSubdetector::TIB ||
160  IntSubDetID == StripSubdetector::TOB) {
161 
162  uOrientation = deltaPhi(gUDirection.barePhi(),gPModule.barePhi()) >= 0. ? +1.F : -1.F;
163  vOrientation = gVDirection.z() - gPModule.z() >= 0 ? +1.F : -1.F;
164  resXTopol = res.x();
165  resXatTrkYTopol = res.x();
166  resYTopol = res.y();
167  resXprimeErr = resXErr;
168  resYprimeErr = resYErr;
169 
170  const RectangularPlaneBounds *rectangularBound = dynamic_cast<const RectangularPlaneBounds*>(&bound);
171  if (rectangularBound!=NULL) {
172  hitStruct.inside = rectangularBound->inside(lPTrk);
173  length = rectangularBound->length();
174  width = rectangularBound->width();
175  hitStruct.localXnorm = 2*hitStruct.localX/width;
176  hitStruct.localYnorm = 2*hitStruct.localY/length;
177  } else {
178  throw cms::Exception("Geometry Error")
179  << "[TrackerValidationVariables] Cannot cast bounds to RectangularPlaneBounds as expected for TPB, TIB and TOB";
180  }
181 
182  } else if (IntSubDetID == PixelSubdetector::PixelEndcap) {
183 
184  uOrientation = gUDirection.perp() - gPModule.perp() >= 0 ? +1.F : -1.F;
185  vOrientation = deltaPhi(gVDirection.barePhi(),gPModule.barePhi()) >= 0. ? +1.F : -1.F;
186  resXTopol = res.x();
187  resXatTrkYTopol = res.x();
188  resYTopol = res.y();
189  resXprimeErr = resXErr;
190  resYprimeErr = resYErr;
191 
192  const RectangularPlaneBounds *rectangularBound = dynamic_cast<const RectangularPlaneBounds*>(&bound);
193  if (rectangularBound!=NULL) {
194  hitStruct.inside = rectangularBound->inside(lPTrk);
195  length = rectangularBound->length();
196  width = rectangularBound->width();
197  hitStruct.localXnorm = 2*hitStruct.localX/width;
198  hitStruct.localYnorm = 2*hitStruct.localY/length;
199  } else {
200  throw cms::Exception("Geometry Error")
201  << "[TrackerValidationVariables] Cannot cast bounds to RectangularPlaneBounds as expected for TPE";
202  }
203 
204  } else if (IntSubDetID == StripSubdetector::TID ||
205  IntSubDetID == StripSubdetector::TEC) {
206 
207  uOrientation = deltaPhi(gUDirection.barePhi(),gPModule.barePhi()) >= 0. ? +1.F : -1.F;
208  vOrientation = gVDirection.perp() - gPModule.perp() >= 0. ? +1.F : -1.F;
209 
210  if (!dynamic_cast<const RadialStripTopology*>(&detUnit.type().topology()))continue;
211  const RadialStripTopology& topol = dynamic_cast<const RadialStripTopology&>(detUnit.type().topology());
212 
213  MeasurementPoint measHitPos = topol.measurementPosition(lPHit);
214  MeasurementPoint measTrkPos = topol.measurementPosition(lPTrk);
215 
216  MeasurementError measHitErr = topol.measurementError(lPHit,errHit);
217  MeasurementError measTrkErr = topol.measurementError(lPTrk,errTrk);
218 
219  if (measHitErr.uu()<0. ||
220  measHitErr.vv()<0. ||
221  measTrkErr.uu()<0. ||
222  measTrkErr.vv()<0.){
223  edm::LogError("TrackerValidationVariables") << "@SUB=TrackerValidationVariables::fillHitQuantities"
224  << "One of the squared error methods gives negative result"
225  << "\n\tmeasHitErr.uu()\tmeasHitErr.vv()\tmeasTrkErr.uu()\tmeasTrkErr.vv()"
226  << "\n\t" << measHitErr.uu()
227  << "\t" << measHitErr.vv()
228  << "\t" << measTrkErr.uu()
229  << "\t" << measTrkErr.vv();
230  continue;
231  }
232 
233  float localStripLengthHit = topol.localStripLength(lPHit);
234  float localStripLengthTrk = topol.localStripLength(lPTrk);
235  float phiHit = topol.stripAngle(measHitPos.x());
236  float phiTrk = topol.stripAngle(measTrkPos.x());
237  float r_0 = topol.originToIntersection();
238 
239  resXTopol = (phiTrk-phiHit)*r_0;
240 // resXTopol = (tan(phiTrk)-tan(phiHit))*r_0;
241 
242  LocalPoint LocalHitPosCor= topol.localPosition(MeasurementPoint(measHitPos.x(), measTrkPos.y()));
243  resXatTrkYTopol = lPTrk.x() - LocalHitPosCor.x();
244 
245 
246  //resYTopol = measTrkPos.y()*localStripLengthTrk - measHitPos.y()*localStripLengthHit;
247  float cosPhiHit(cos(phiHit)), cosPhiTrk(cos(phiTrk)),
248  sinPhiHit(sin(phiHit)), sinPhiTrk(sin(phiTrk));
249  float l_0 = r_0 - topol.detHeight()/2;
250  resYTopol = measTrkPos.y()*localStripLengthTrk - measHitPos.y()*localStripLengthHit + l_0*(1/cosPhiTrk - 1/cosPhiHit);
251 
252  resXprimeErr = std::sqrt(measHitErr.uu()+measTrkErr.uu())*topol.angularWidth()*r_0;
253  //resYprimeErr = std::sqrt(measHitErr.vv()*localStripLengthHit*localStripLengthHit + measTrkErr.vv()*localStripLengthTrk*localStripLengthTrk);
254  float helpSummand = l_0*l_0*topol.angularWidth()*topol.angularWidth()*(sinPhiHit*sinPhiHit/pow(cosPhiHit,4)*measHitErr.uu()
255  + sinPhiTrk*sinPhiTrk/pow(cosPhiTrk,4)*measTrkErr.uu() );
256  resYprimeErr = std::sqrt(measHitErr.vv()*localStripLengthHit*localStripLengthHit
257  + measTrkErr.vv()*localStripLengthTrk*localStripLengthTrk + helpSummand );
258 
259 
260  const TrapezoidalPlaneBounds *trapezoidalBound = dynamic_cast < const TrapezoidalPlaneBounds * >(& bound);
261  if (trapezoidalBound!=NULL) {
262  hitStruct.inside = trapezoidalBound->inside(lPTrk);
263  length = trapezoidalBound->length();
264  width = trapezoidalBound->width();
265  //float widthAtHalfLength = trapezoidalBound->widthAtHalfLength();
266 
267  // int yAxisOrientation=trapezoidalBound->yAxisOrientation();
268  // for trapezoidal shape modules, scale with as function of local y coordinate
269  // float widthAtlocalY=width-(1-yAxisOrientation*2*lPTrk.y()/length)*(width-widthAtHalfLength);
270  // hitStruct.localXnorm = 2*hitStruct.localX/widthAtlocalY;
271  hitStruct.localXnorm = 2*hitStruct.localX/width;
272  hitStruct.localYnorm = 2*hitStruct.localY/length;
273  } else {
274  throw cms::Exception("Geometry Error")
275  << "[TrackerValidationVariables] Cannot cast bounds to TrapezoidalPlaneBounds as expected for TID and TEC";
276  }
277 
278  } else {
279  edm::LogWarning("TrackerValidationVariables") << "@SUB=TrackerValidationVariables::fillHitQuantities"
280  << "No valid tracker subdetector " << IntSubDetID;
281  continue;
282  }
283 
284  resXprime = resXTopol*uOrientation;
285  resXatTrkY = resXatTrkYTopol;
286  resYprime = resYTopol*vOrientation;
287 
288  } else { // not a detUnit, so must be a virtual 2D-Module
289  // FIXME: at present only for det units residuals are calculated and filled in the hitStruct
290  // But in principle this method should also be useable for the gluedDets (2D modules in TIB, TID, TOB, TEC)
291  // In this case, only orientation should be taken into account for primeResiduals, but not the radial topology
292  // At present, default values (999.F) are given out
293  }
294 
295  hitStruct.resXprime = resXprime;
296  hitStruct.resXatTrkY = resXatTrkY;
297  hitStruct.resYprime = resYprime;
298  hitStruct.resXprimeErr = resXprimeErr;
299  hitStruct.resYprimeErr = resYprimeErr;
300 
301  hitStruct.rawDetId = IntRawDetID;
302  hitStruct.phi = tsos.globalDirection().phi();
303  hitStruct.eta = tsos.globalDirection().eta();
304 
305  v_avhitout.push_back(hitStruct);
306  }
307 }
308 
309 void
311  const edm::EventSetup& eventSetup,
312  std::vector<AVTrackStruct> & v_avtrackout)
313 {
314  fillTrackQuantities(event,
315  eventSetup,
316  [](const reco::Track&) -> bool { return true; },
317  v_avtrackout);
318 }
319 
320 void
322  const edm::EventSetup& eventSetup,
323  std::function<bool(const reco::Track&)> trackFilter,
324  std::vector<AVTrackStruct> & v_avtrackout)
325 {
327  eventSetup.get<IdealMagneticFieldRecord>().get(magneticField);
328 
330  event.getByToken(trajTracksToken_, TrajTracksMap);
331  LogDebug("TrackerValidationVariables") << "TrajTrack collection size " << TrajTracksMap->size();
332 
333  const Trajectory* trajectory;
334  const reco::Track* track;
335 
336  for ( TrajTrackAssociationCollection::const_iterator iPair = TrajTracksMap->begin();
337  iPair != TrajTracksMap->end();
338  ++iPair) {
339 
340  trajectory = &(*(*iPair).key);
341  track = &(*(*iPair).val);
342 
343  if (!trackFilter(*track)) continue;
344 
345  AVTrackStruct trackStruct;
346 
347  trackStruct.p = track->p();
348  trackStruct.pt = track->pt();
349  trackStruct.ptError = track->ptError();
350  trackStruct.px = track->px();
351  trackStruct.py = track->py();
352  trackStruct.pz = track->pz();
353  trackStruct.eta = track->eta();
354  trackStruct.phi = track->phi();
355  trackStruct.chi2 = track->chi2();
356  trackStruct.chi2Prob= TMath::Prob(track->chi2(),track->ndof());
357  trackStruct.normchi2 = track->normalizedChi2();
358  GlobalPoint gPoint(track->vx(), track->vy(), track->vz());
359  double theLocalMagFieldInInverseGeV = magneticField->inInverseGeV(gPoint).z();
360  trackStruct.kappa = -track->charge()*theLocalMagFieldInInverseGeV/track->pt();
361  trackStruct.charge = track->charge();
362  trackStruct.d0 = track->d0();
363  trackStruct.dz = track->dz();
364  trackStruct.numberOfValidHits = track->numberOfValidHits();
365  trackStruct.numberOfLostHits = track->numberOfLostHits();
366 
367  fillHitQuantities(trajectory, trackStruct.hits);
368 
369  v_avtrackout.push_back(trackStruct);
370  }
371 }
372 
#define LogDebug(id)
double p() const
momentum vector magnitude
Definition: TrackBase.h:610
GlobalPoint toGlobal(const Point2DBase< Scalar, LocalTag > lp) const
Definition: Surface.h:106
T getParameter(std::string const &) const
float xx() const
Definition: LocalError.h:24
float vv() const
T y() const
Definition: PV2DBase.h:46
T perp() const
Definition: PV3DBase.h:72
virtual float localStripLength(const LocalPoint &) const =0
double d0() const
dxy parameter in perigee convention (d0 = -dxy)
Definition: TrackBase.h:592
virtual const GeomDetType & type() const
Definition: GeomDet.cc:85
double normalizedChi2() const
chi-squared divided by n.d.o.f. (or chi-squared * 1e6 if n.d.o.f. is zero)
Definition: TrackBase.h:556
LocalVector localDirection() const
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
T y() const
Definition: PV3DBase.h:63
double phi() const
azimuthal angle of momentum vector
Definition: TrackBase.h:640
unsigned short numberOfLostHits() const
number of cases where track crossed a layer without getting a hit.
Definition: TrackBase.h:821
#define NULL
Definition: scimark2.h:8
void fillHitQuantities(const Trajectory *trajectory, std::vector< AVHitStruct > &v_avhitout)
virtual LocalPoint localPosition(float strip) const =0
double px() const
x coordinate of momentum vector
Definition: TrackBase.h:622
virtual float detHeight() const =0
LocalError positionError() const
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
T barePhi() const
Definition: PV3DBase.h:68
virtual bool inside(const Local2DPoint &p) const
DataContainer const & measurements() const
Definition: Trajectory.h:250
Measurement2DPoint MeasurementPoint
Measurement points are two-dimensional by default.
virtual MeasurementError measurementError(const LocalPoint &, const LocalError &) const =0
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:646
float yy() const
Definition: LocalError.h:26
double chi2() const
chi-squared of the fit
Definition: TrackBase.h:544
edm::EDGetTokenT< std::vector< Trajectory > > trajCollectionToken_
std::shared_ptr< TrackingRecHit const > ConstRecHitPointer
double ndof() const
number of degrees of freedom of the fit
Definition: TrackBase.h:550
T sqrt(T t)
Definition: SSEVec.h:18
helper::RootFunctionHelper< F, args >::root_function function(F &f)
Definition: rootFunction.h:14
double pt() const
track transverse momentum
Definition: TrackBase.h:616
T z() const
Definition: PV3DBase.h:64
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
double ptError() const
error on Pt (set to 1000 TeV if charge==0 for safety)
Definition: TrackBase.h:758
float uu() const
virtual MeasurementPoint measurementPosition(const LocalPoint &) const =0
virtual float angularWidth() const =0
unsigned short numberOfValidHits() const
number of valid hits found
Definition: TrackBase.h:815
const LocalTrajectoryError & localError() const
virtual float width() const
Width along local X.
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:37
double pz() const
z coordinate of momentum vector
Definition: TrackBase.h:634
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:604
tuple trackFilter
Definition: valSkim_cff.py:14
double vz() const
z coordinate of the reference point on track
Definition: TrackBase.h:664
virtual float stripAngle(float strip) const =0
Definition: DetId.h:18
void fillTrackQuantities(const edm::Event &, const edm::EventSetup &, std::vector< AVTrackStruct > &v_avtrackout)
edm::EDGetTokenT< TrajTrackAssociationCollection > trajTracksToken_
const T & get() const
Definition: EventSetup.h:56
virtual float originToIntersection() const =0
virtual const Topology & topology() const =0
T eta() const
Definition: PV3DBase.h:76
double vy() const
y coordinate of the reference point on track
Definition: TrackBase.h:658
virtual float length() const
Lenght along local Y.
int charge() const
track electric charge
Definition: TrackBase.h:562
Definition: Bounds.h:22
static uInt32 F(BLOWFISH_CTX *ctx, uInt32 x)
Definition: blowfish.cc:281
T x() const
Definition: PV2DBase.h:45
T x() const
Definition: PV3DBase.h:62
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
double py() const
y coordinate of momentum vector
Definition: TrackBase.h:628
double vx() const
x coordinate of the reference point on track
Definition: TrackBase.h:652
GlobalVector globalDirection() const