CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
SiPixelPhase1TrackResiduals.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: SiPixelPhase1TrackResiduals
4 // Class: SiPixelPhase1TrackResiduals
5 //
6 
7 // Original Author: Marcel Schneider
8 
11 
14 
15 
17  SiPixelPhase1Base(iConfig),
18  validator(iConfig, consumesCollector())
19 {
20  offlinePrimaryVerticesToken_ = consumes<reco::VertexCollection>(std::string("offlinePrimaryVertices"));
21 }
22 
24 
26  iEvent.getByToken(offlinePrimaryVerticesToken_, vertices);
27  if (!vertices.isValid() || vertices->size() == 0) return;
28  const auto primaryVertex = vertices->at(0);
29 
30  std::vector<TrackerValidationVariables::AVTrackStruct> vtracks;
31  validator.fillTrackQuantities(iEvent, iSetup,
32  // tell the validator to only look at good tracks
33  [&](const reco::Track& track) -> bool {
34  return track.pt() > 0.75
35  && std::abs( track.dxy(primaryVertex.position()) ) < 5*track.dxyError();
36  }, vtracks);
37 
38  for (auto& track : vtracks) {
39  for (auto& it : track.hits) {
40  auto id = DetId(it.rawDetId);
41  auto isPixel = id.subdetId() == 1 || id.subdetId() == 2;
42  if (!isPixel) continue;
43 
44  histo[RESIDUAL_X].fill(it.resXprime, id, &iEvent);
45  histo[RESIDUAL_Y].fill(it.resYprime, id, &iEvent);
46  }
47  }
48 
49 }
50 
52 
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
double dxyError() const
error on dxy
Definition: TrackBase.h:791
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
void analyze(const edm::Event &, const edm::EventSetup &)
int iEvent
Definition: GenABIO.cc:230
double pt() const
track transverse momentum
Definition: TrackBase.h:616
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
SiPixelPhase1TrackResiduals(const edm::ParameterSet &conf)
bool isValid() const
Definition: HandleBase.h:75
TrackerValidationVariables validator
Definition: DetId.h:18
void fillTrackQuantities(const edm::Event &, const edm::EventSetup &, std::vector< AVTrackStruct > &v_avtrackout)
edm::EDGetTokenT< reco::VertexCollection > offlinePrimaryVerticesToken_
std::vector< HistogramManager > histo
bool isPixel(HitType hitType)
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:586