CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Classes | Public Member Functions | Private Attributes
TrackerValidationVariables Class Reference

#include <TrackerValidationVariables.h>

Classes

struct  AVHitStruct
 
struct  AVTrackStruct
 

Public Member Functions

void fillHitQuantities (const Trajectory *trajectory, std::vector< AVHitStruct > &v_avhitout)
 
void fillHitQuantities (const edm::Event &, std::vector< AVHitStruct > &v_avhitout)
 
void fillTrackQuantities (const edm::Event &, std::vector< AVTrackStruct > &v_avtrackout)
 
 TrackerValidationVariables ()
 
 TrackerValidationVariables (const edm::EventSetup &, const edm::ParameterSet &)
 
 ~TrackerValidationVariables ()
 

Private Attributes

const edm::ParameterSet conf_
 
edm::ESHandle< MagneticFieldmagneticField_
 
edm::ESHandle< TrackerGeometrytkGeom_
 

Detailed Description

Definition at line 16 of file TrackerValidationVariables.h.

Constructor & Destructor Documentation

TrackerValidationVariables::TrackerValidationVariables ( )

Definition at line 45 of file TrackerValidationVariables.cc.

46 {
47 
48 }
TrackerValidationVariables::TrackerValidationVariables ( const edm::EventSetup es,
const edm::ParameterSet iSetup 
)

Definition at line 50 of file TrackerValidationVariables.cc.

References edm::EventSetup::get(), magneticField_, and tkGeom_.

TrackerValidationVariables::~TrackerValidationVariables ( )

Definition at line 57 of file TrackerValidationVariables.cc.

58 {
59 
60 }

Member Function Documentation

void TrackerValidationVariables::fillHitQuantities ( const Trajectory trajectory,
std::vector< AVHitStruct > &  v_avhitout 
)

Definition at line 63 of file TrackerValidationVariables.cc.

References RadialStripTopology::angularWidth(), funct::cos(), SiPixelRawToDigiRegional_cfi::deltaPhi, RadialStripTopology::detHeight(), TrackerValidationVariables::AVHitStruct::eta, PV3DBase< T, PVType, FrameType >::eta(), edm::hlt::Exception, F(), TrajectoryStateOnSurface::globalDirection(), TrackerValidationVariables::AVHitStruct::inside, TrackerValidationVariables::AVHitStruct::localAlpha, TrackerValidationVariables::AVHitStruct::localBeta, TrajectoryStateOnSurface::localDirection(), TrajectoryStateOnSurface::localError(), RadialStripTopology::localPosition(), TrajectoryStateOnSurface::localPosition(), RadialStripTopology::localStripLength(), TrackerValidationVariables::AVHitStruct::localX, TrackerValidationVariables::AVHitStruct::localXnorm, TrackerValidationVariables::AVHitStruct::localY, TrackerValidationVariables::AVHitStruct::localYnorm, RadialStripTopology::measurementError(), RadialStripTopology::measurementPosition(), Trajectory::measurements(), NULL, RadialStripTopology::originToIntersection(), PV3DBase< T, PVType, FrameType >::perp(), TrackerValidationVariables::AVHitStruct::phi, PV3DBase< T, PVType, FrameType >::phi(), PixelSubdetector::PixelBarrel, PixelSubdetector::PixelEndcap, LocalTrajectoryError::positionError(), funct::pow(), TrackerValidationVariables::AVHitStruct::rawDetId, DetId::rawId(), TrackerValidationVariables::AVHitStruct::resErrX, TrackerValidationVariables::AVHitStruct::resErrY, TrackerValidationVariables::AVHitStruct::resX, TrackerValidationVariables::AVHitStruct::resXatTrkY, TrackerValidationVariables::AVHitStruct::resXprime, TrackerValidationVariables::AVHitStruct::resXprimeErr, TrackerValidationVariables::AVHitStruct::resY, TrackerValidationVariables::AVHitStruct::resYprime, TrackerValidationVariables::AVHitStruct::resYprimeErr, funct::sin(), mathSSE::sqrt(), RadialStripTopology::stripAngle(), DetId::subdetId(), StripSubdetector::TEC, StripSubdetector::TIB, StripSubdetector::TID, StripSubdetector::TOB, Surface::toGlobal(), GeomDetType::topology(), DetId::Tracker, GeomDet::type(), MeasurementError::uu(), MeasurementError::vv(), create_public_lumi_plots::width, PV2DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::x(), LocalError::xx(), PV2DBase< T, PVType, FrameType >::y(), PV3DBase< T, PVType, FrameType >::y(), LocalError::yy(), and PV3DBase< T, PVType, FrameType >::z().

Referenced by MonitorTrackResiduals::analyze(), fillHitQuantities(), and fillTrackQuantities().

64 {
65  TrajectoryStateCombiner tsoscomb;
66 
67  const std::vector<TrajectoryMeasurement> &tmColl = trajectory->measurements();
68  for(std::vector<TrajectoryMeasurement>::const_iterator itTraj = tmColl.begin();
69  itTraj != tmColl.end();
70  ++itTraj) {
71 
72  if (!itTraj->updatedState().isValid()) continue;
73 
74  TrajectoryStateOnSurface tsos = tsoscomb( itTraj->forwardPredictedState(), itTraj->backwardPredictedState() );
76 
77  if(!hit->isValid() || hit->geographicalId().det() != DetId::Tracker) continue;
78 
79  AVHitStruct hitStruct;
80  const DetId& hit_detId = hit->geographicalId();
81  unsigned int IntRawDetID = (hit_detId.rawId());
82  unsigned int IntSubDetID = (hit_detId.subdetId());
83 
84  if(IntSubDetID == 0) continue;
85 
86  //first calculate residuals in cartesian coordinates in the local module coordinate system
87 
88  LocalPoint lPHit = hit->localPosition();
89  LocalPoint lPTrk = tsos.localPosition();
90  LocalVector lVTrk = tsos.localDirection();
91 
92  hitStruct.localAlpha = atan2(lVTrk.x(), lVTrk.z()); // wrt. normal tg(alpha)=x/z
93  hitStruct.localBeta = atan2(lVTrk.y(), lVTrk.z()); // wrt. normal tg(beta)= y/z
94 
95  LocalError errHit = hit->localPositionError();
96  // no need to add APE to hitError anymore
97  // AlgebraicROOTObject<2>::SymMatrix mat = asSMatrix<2>(hit->parametersError());
98  // LocalError errHit = LocalError( mat(0,0),mat(0,1),mat(1,1) );
99  LocalError errTrk = tsos.localError().positionError();
100 
101  //check for negative error values: track error can have negative value, if matrix inversion fails (very rare case)
102  //hit error should always give positive values
103  if(errHit.xx()<0. || errHit.yy()<0. || errTrk.xx()<0. || errTrk.yy()<0.){
104  edm::LogError("TrackerValidationVariables") << "@SUB=TrackerValidationVariables::fillHitQuantities"
105  << "One of the squared error methods gives negative result"
106  << "\n\terrHit.xx()\terrHit.yy()\terrTrk.xx()\terrTrk.yy()"
107  << "\n\t" << errHit.xx()
108  << "\t" << errHit.yy()
109  << "\t" << errTrk.xx()
110  << "\t" << errTrk.yy();
111  continue;
112  }
113 
114  align::LocalVector res = lPTrk - lPHit;
115 
116  float resXErr = std::sqrt( errHit.xx() + errTrk.xx() );
117  float resYErr = std::sqrt( errHit.yy() + errTrk.yy() );
118 
119  hitStruct.resX = res.x();
120  hitStruct.resY = res.y();
121  hitStruct.resErrX = resXErr;
122  hitStruct.resErrY = resYErr;
123 
124  // hitStruct.localX = lPhit.x();
125  // hitStruct.localY = lPhit.y();
126  // EM: use predictions for local coordinates
127  hitStruct.localX = lPTrk.x();
128  hitStruct.localY = lPTrk.y();
129 
130  // now calculate residuals taking global orientation of modules and radial topology in TID/TEC into account
131  float resXprime(999.F), resYprime(999.F);
132  float resXatTrkY(999.F);
133  float resXprimeErr(999.F), resYprimeErr(999.F);
134 
135  if(hit->detUnit()){ // is it a single physical module?
136  const GeomDetUnit& detUnit = *(hit->detUnit());
137  float uOrientation(-999.F), vOrientation(-999.F);
138  float resXTopol(999.F), resYTopol(999.F);
139  float resXatTrkYTopol(999.F);
140 
141  const Surface& surface = hit->detUnit()->surface();
142  const BoundPlane& boundplane = hit->detUnit()->surface();
143  const Bounds& bound = boundplane.bounds();
144 
145  float length = 0;
146  float width = 0;
147 
148  LocalPoint lPModule(0.,0.,0.), lUDirection(1.,0.,0.), lVDirection(0.,1.,0.);
149  GlobalPoint gPModule = surface.toGlobal(lPModule),
150  gUDirection = surface.toGlobal(lUDirection),
151  gVDirection = surface.toGlobal(lVDirection);
152 
153  if (IntSubDetID == PixelSubdetector::PixelBarrel ||
154  IntSubDetID == StripSubdetector::TIB ||
155  IntSubDetID == StripSubdetector::TOB) {
156 
157  uOrientation = deltaPhi(gUDirection.phi(),gPModule.phi()) >= 0. ? +1.F : -1.F;
158  vOrientation = gVDirection.z() - gPModule.z() >= 0 ? +1.F : -1.F;
159  resXTopol = res.x();
160  resXatTrkYTopol = res.x();
161  resYTopol = res.y();
162  resXprimeErr = resXErr;
163  resYprimeErr = resYErr;
164 
165  const RectangularPlaneBounds *rectangularBound = dynamic_cast<const RectangularPlaneBounds*>(&bound);
166  if (rectangularBound!=NULL) {
167  hitStruct.inside = rectangularBound->inside(lPTrk);
168  length = rectangularBound->length();
169  width = rectangularBound->width();
170  hitStruct.localXnorm = 2*hitStruct.localX/width;
171  hitStruct.localYnorm = 2*hitStruct.localY/length;
172  } else {
173  throw cms::Exception("Geometry Error")
174  << "[TrackerValidationVariables] Cannot cast bounds to RectangularPlaneBounds as expected for TPB, TIB and TOB";
175  }
176 
177  } else if (IntSubDetID == PixelSubdetector::PixelEndcap) {
178 
179  uOrientation = gUDirection.perp() - gPModule.perp() >= 0 ? +1.F : -1.F;
180  vOrientation = deltaPhi(gVDirection.phi(),gPModule.phi()) >= 0. ? +1.F : -1.F;
181  resXTopol = res.x();
182  resXatTrkYTopol = res.x();
183  resYTopol = res.y();
184  resXprimeErr = resXErr;
185  resYprimeErr = resYErr;
186 
187  const RectangularPlaneBounds *rectangularBound = dynamic_cast<const RectangularPlaneBounds*>(&bound);
188  if (rectangularBound!=NULL) {
189  hitStruct.inside = rectangularBound->inside(lPTrk);
190  length = rectangularBound->length();
191  width = rectangularBound->width();
192  hitStruct.localXnorm = 2*hitStruct.localX/width;
193  hitStruct.localYnorm = 2*hitStruct.localY/length;
194  } else {
195  throw cms::Exception("Geometry Error")
196  << "[TrackerValidationVariables] Cannot cast bounds to RectangularPlaneBounds as expected for TPE";
197  }
198 
199  } else if (IntSubDetID == StripSubdetector::TID ||
200  IntSubDetID == StripSubdetector::TEC) {
201 
202  uOrientation = deltaPhi(gUDirection.phi(),gPModule.phi()) >= 0. ? +1.F : -1.F;
203  vOrientation = gVDirection.perp() - gPModule.perp() >= 0. ? +1.F : -1.F;
204 
205  if (!dynamic_cast<const RadialStripTopology*>(&detUnit.type().topology()))continue;
206  const RadialStripTopology& topol = dynamic_cast<const RadialStripTopology&>(detUnit.type().topology());
207 
208  MeasurementPoint measHitPos = topol.measurementPosition(lPHit);
209  MeasurementPoint measTrkPos = topol.measurementPosition(lPTrk);
210 
211  MeasurementError measHitErr = topol.measurementError(lPHit,errHit);
212  MeasurementError measTrkErr = topol.measurementError(lPTrk,errTrk);
213 
214  if (measHitErr.uu()<0. ||
215  measHitErr.vv()<0. ||
216  measTrkErr.uu()<0. ||
217  measTrkErr.vv()<0.){
218  edm::LogError("TrackerValidationVariables") << "@SUB=TrackerValidationVariables::fillHitQuantities"
219  << "One of the squared error methods gives negative result"
220  << "\n\tmeasHitErr.uu()\tmeasHitErr.vv()\tmeasTrkErr.uu()\tmeasTrkErr.vv()"
221  << "\n\t" << measHitErr.uu()
222  << "\t" << measHitErr.vv()
223  << "\t" << measTrkErr.uu()
224  << "\t" << measTrkErr.vv();
225  continue;
226  }
227 
228  float localStripLengthHit = topol.localStripLength(lPHit);
229  float localStripLengthTrk = topol.localStripLength(lPTrk);
230  float phiHit = topol.stripAngle(measHitPos.x());
231  float phiTrk = topol.stripAngle(measTrkPos.x());
232  float r_0 = topol.originToIntersection();
233 
234  resXTopol = (phiTrk-phiHit)*r_0;
235 // resXTopol = (tan(phiTrk)-tan(phiHit))*r_0;
236 
237  LocalPoint LocalHitPosCor= topol.localPosition(MeasurementPoint(measHitPos.x(), measTrkPos.y()));
238  resXatTrkYTopol = lPTrk.x() - LocalHitPosCor.x();
239 
240 
241  //resYTopol = measTrkPos.y()*localStripLengthTrk - measHitPos.y()*localStripLengthHit;
242  float cosPhiHit(cos(phiHit)), cosPhiTrk(cos(phiTrk)),
243  sinPhiHit(sin(phiHit)), sinPhiTrk(sin(phiTrk));
244  float l_0 = r_0 - topol.detHeight()/2;
245  resYTopol = measTrkPos.y()*localStripLengthTrk - measHitPos.y()*localStripLengthHit + l_0*(1/cosPhiTrk - 1/cosPhiHit);
246 
247  resXprimeErr = std::sqrt(measHitErr.uu()+measTrkErr.uu())*topol.angularWidth()*r_0;
248  //resYprimeErr = std::sqrt(measHitErr.vv()*localStripLengthHit*localStripLengthHit + measTrkErr.vv()*localStripLengthTrk*localStripLengthTrk);
249  float helpSummand = l_0*l_0*topol.angularWidth()*topol.angularWidth()*(sinPhiHit*sinPhiHit/pow(cosPhiHit,4)*measHitErr.uu()
250  + sinPhiTrk*sinPhiTrk/pow(cosPhiTrk,4)*measTrkErr.uu() );
251  resYprimeErr = std::sqrt(measHitErr.vv()*localStripLengthHit*localStripLengthHit
252  + measTrkErr.vv()*localStripLengthTrk*localStripLengthTrk + helpSummand );
253 
254 
255  const TrapezoidalPlaneBounds *trapezoidalBound = dynamic_cast < const TrapezoidalPlaneBounds * >(& bound);
256  if (trapezoidalBound!=NULL) {
257  hitStruct.inside = trapezoidalBound->inside(lPTrk);
258  length = trapezoidalBound->length();
259  width = trapezoidalBound->width();
260  //float widthAtHalfLength = trapezoidalBound->widthAtHalfLength();
261 
262  // int yAxisOrientation=trapezoidalBound->yAxisOrientation();
263  // for trapezoidal shape modules, scale with as function of local y coordinate
264  // float widthAtlocalY=width-(1-yAxisOrientation*2*lPTrk.y()/length)*(width-widthAtHalfLength);
265  // hitStruct.localXnorm = 2*hitStruct.localX/widthAtlocalY;
266  hitStruct.localXnorm = 2*hitStruct.localX/width;
267  hitStruct.localYnorm = 2*hitStruct.localY/length;
268  } else {
269  throw cms::Exception("Geometry Error")
270  << "[TrackerValidationVariables] Cannot cast bounds to TrapezoidalPlaneBounds as expected for TID and TEC";
271  }
272 
273  } else {
274  edm::LogWarning("TrackerValidationVariables") << "@SUB=TrackerValidationVariables::fillHitQuantities"
275  << "No valid tracker subdetector " << IntSubDetID;
276  continue;
277  }
278 
279  resXprime = resXTopol*uOrientation;
280  resXatTrkY = resXatTrkYTopol;
281  resYprime = resYTopol*vOrientation;
282 
283  } else { // not a detUnit, so must be a virtual 2D-Module
284  // FIXME: at present only for det units residuals are calculated and filled in the hitStruct
285  // But in principle this method should also be useable for the gluedDets (2D modules in TIB, TID, TOB, TEC)
286  // In this case, only orientation should be taken into account for primeResiduals, but not the radial topology
287  // At present, default values (999.F) are given out
288  }
289 
290  hitStruct.resXprime = resXprime;
291  hitStruct.resXatTrkY = resXatTrkY;
292  hitStruct.resYprime = resYprime;
293  hitStruct.resXprimeErr = resXprimeErr;
294  hitStruct.resYprimeErr = resYprimeErr;
295 
296  hitStruct.rawDetId = IntRawDetID;
297  hitStruct.phi = tsos.globalDirection().phi();
298  hitStruct.eta = tsos.globalDirection().eta();
299 
300  v_avhitout.push_back(hitStruct);
301  }
302 }
GlobalPoint toGlobal(const Point2DBase< Scalar, LocalTag > lp) const
Definition: Surface.h:114
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
virtual const GeomDetType & type() const
Definition: GeomDet.cc:90
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
#define NULL
Definition: scimark2.h:8
virtual LocalPoint localPosition(float strip) const =0
virtual float detHeight() const =0
LocalError positionError() const
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
DataContainer const & measurements() const
Definition: Trajectory.h:203
virtual bool inside(const Local2DPoint &p) const
Measurement2DPoint MeasurementPoint
Measurement points are two-dimensional by default.
virtual MeasurementError measurementError(const LocalPoint &, const LocalError &) const =0
float yy() const
Definition: LocalError.h:26
std::shared_ptr< TrackingRecHit const > ConstRecHitPointer
T sqrt(T t)
Definition: SSEVec.h:48
T z() const
Definition: PV3DBase.h:64
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
float uu() const
virtual MeasurementPoint measurementPosition(const LocalPoint &) const =0
virtual float angularWidth() const =0
const LocalTrajectoryError & localError() const
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:37
virtual float stripAngle(float strip) const =0
Definition: DetId.h:18
virtual float originToIntersection() const =0
virtual const Topology & topology() const =0
T eta() const
Definition: PV3DBase.h:76
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
GlobalVector globalDirection() const
void TrackerValidationVariables::fillHitQuantities ( const edm::Event event,
std::vector< AVHitStruct > &  v_avhitout 
)

Definition at line 351 of file TrackerValidationVariables.cc.

References conf_, fillHitQuantities(), edm::ParameterSet::getParameter(), LogDebug, and AlCaHLTBitMon_QueryRunRegistry::string.

352 {
353  edm::Handle<std::vector<Trajectory> > trajCollectionHandle;
354  event.getByLabel(conf_.getParameter<std::string>("trajectoryInput"), trajCollectionHandle);
355 
356  LogDebug("TrackerValidationVariables") << "trajColl->size(): " << trajCollectionHandle->size() ;
357 
358  for (std::vector<Trajectory>::const_iterator it = trajCollectionHandle->begin(), itEnd = trajCollectionHandle->end();
359  it!=itEnd;
360  ++it) {
361 
362  fillHitQuantities(&(*it), v_avhitout);
363  }
364 }
#define LogDebug(id)
T getParameter(std::string const &) const
void fillHitQuantities(const Trajectory *trajectory, std::vector< AVHitStruct > &v_avhitout)
void TrackerValidationVariables::fillTrackQuantities ( const edm::Event event,
std::vector< AVTrackStruct > &  v_avtrackout 
)

Definition at line 305 of file TrackerValidationVariables.cc.

References TrackerValidationVariables::AVTrackStruct::charge, reco::TrackBase::charge(), TrackerValidationVariables::AVTrackStruct::chi2, reco::TrackBase::chi2(), TrackerValidationVariables::AVTrackStruct::chi2Prob, conf_, edm::AssociationMap< edm::OneToOne< std::vector< Trajectory >, reco::TrackCollection, unsigned short > >::const_iterator, TrackerValidationVariables::AVTrackStruct::d0, reco::TrackBase::d0(), TrackerValidationVariables::AVTrackStruct::dz, reco::TrackBase::dz(), TrackerValidationVariables::AVTrackStruct::eta, reco::TrackBase::eta(), fillHitQuantities(), edm::ParameterSet::getParameter(), TrackerValidationVariables::AVTrackStruct::hits, TrackerValidationVariables::AVTrackStruct::kappa, LogDebug, magneticField_, reco::TrackBase::ndof(), reco::TrackBase::normalizedChi2(), TrackerValidationVariables::AVTrackStruct::normchi2, TrackerValidationVariables::AVTrackStruct::numberOfLostHits, reco::TrackBase::numberOfLostHits(), TrackerValidationVariables::AVTrackStruct::numberOfValidHits, reco::TrackBase::numberOfValidHits(), TrackerValidationVariables::AVTrackStruct::p, reco::TrackBase::p(), TrackerValidationVariables::AVTrackStruct::phi, reco::TrackBase::phi(), TrackerValidationVariables::AVTrackStruct::pt, reco::TrackBase::pt(), TrackerValidationVariables::AVTrackStruct::ptError, reco::TrackBase::ptError(), TrackerValidationVariables::AVTrackStruct::px, reco::TrackBase::px(), TrackerValidationVariables::AVTrackStruct::py, reco::TrackBase::py(), TrackerValidationVariables::AVTrackStruct::pz, reco::TrackBase::pz(), reco::TrackBase::vx(), reco::TrackBase::vy(), and reco::TrackBase::vz().

Referenced by TrackerOfflineValidation::analyze().

306 {
307  edm::InputTag TrjTrackTag = conf_.getParameter<edm::InputTag>("Tracks");
309  event.getByLabel(TrjTrackTag, TrajTracksMap);
310  LogDebug("TrackerValidationVariables") << "TrajTrack collection size " << TrajTracksMap->size();
311 
312  const Trajectory* trajectory;
313  const reco::Track* track;
314 
315  for ( TrajTrackAssociationCollection::const_iterator iPair = TrajTracksMap->begin();
316  iPair != TrajTracksMap->end();
317  ++iPair) {
318 
319  trajectory = &(*(*iPair).key);
320  track = &(*(*iPair).val);
321 
322  AVTrackStruct trackStruct;
323 
324  trackStruct.p = track->p();
325  trackStruct.pt = track->pt();
326  trackStruct.ptError = track->ptError();
327  trackStruct.px = track->px();
328  trackStruct.py = track->py();
329  trackStruct.pz = track->pz();
330  trackStruct.eta = track->eta();
331  trackStruct.phi = track->phi();
332  trackStruct.chi2 = track->chi2();
333  trackStruct.chi2Prob= TMath::Prob(track->chi2(),track->ndof());
334  trackStruct.normchi2 = track->normalizedChi2();
335  GlobalPoint gPoint(track->vx(), track->vy(), track->vz());
336  double theLocalMagFieldInInverseGeV = magneticField_->inInverseGeV(gPoint).z();
337  trackStruct.kappa = -track->charge()*theLocalMagFieldInInverseGeV/track->pt();
338  trackStruct.charge = track->charge();
339  trackStruct.d0 = track->d0();
340  trackStruct.dz = track->dz();
341  trackStruct.numberOfValidHits = track->numberOfValidHits();
342  trackStruct.numberOfLostHits = track->numberOfLostHits();
343 
344  fillHitQuantities(trajectory, trackStruct.hits);
345 
346  v_avtrackout.push_back(trackStruct);
347  }
348 }
#define LogDebug(id)
double p() const
momentum vector magnitude
Definition: TrackBase.h:591
T getParameter(std::string const &) const
double d0() const
dxy parameter in perigee convention (d0 = -dxy)
Definition: TrackBase.h:573
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:537
double phi() const
azimuthal angle of momentum vector
Definition: TrackBase.h:621
unsigned short numberOfLostHits() const
number of cases where track crossed a layer without getting a hit.
Definition: TrackBase.h:802
void fillHitQuantities(const Trajectory *trajectory, std::vector< AVHitStruct > &v_avhitout)
double px() const
x coordinate of momentum vector
Definition: TrackBase.h:603
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:627
edm::ESHandle< MagneticField > magneticField_
double chi2() const
chi-squared of the fit
Definition: TrackBase.h:525
double ndof() const
number of degrees of freedom of the fit
Definition: TrackBase.h:531
double pt() const
track transverse momentum
Definition: TrackBase.h:597
double ptError() const
error on Pt (set to 1000 TeV if charge==0 for safety)
Definition: TrackBase.h:739
unsigned short numberOfValidHits() const
number of valid hits found
Definition: TrackBase.h:796
double pz() const
z coordinate of momentum vector
Definition: TrackBase.h:615
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:585
double vz() const
z coordinate of the reference point on track
Definition: TrackBase.h:645
double vy() const
y coordinate of the reference point on track
Definition: TrackBase.h:639
int charge() const
track electric charge
Definition: TrackBase.h:543
double py() const
y coordinate of momentum vector
Definition: TrackBase.h:609
double vx() const
x coordinate of the reference point on track
Definition: TrackBase.h:633

Member Data Documentation

const edm::ParameterSet TrackerValidationVariables::conf_
private

Definition at line 83 of file TrackerValidationVariables.h.

Referenced by fillHitQuantities(), and fillTrackQuantities().

edm::ESHandle<MagneticField> TrackerValidationVariables::magneticField_
private

Definition at line 85 of file TrackerValidationVariables.h.

Referenced by fillTrackQuantities(), and TrackerValidationVariables().

edm::ESHandle<TrackerGeometry> TrackerValidationVariables::tkGeom_
private

Definition at line 84 of file TrackerValidationVariables.h.

Referenced by TrackerValidationVariables().