CMS 3D CMS Logo

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 
20 
22  SiPixelPhase1Base(iConfig),
23  validator(iConfig, consumesCollector())
24 {
25  applyVertexCut_=iConfig.getUntrackedParameter<bool>("VertexCut",true);
26  offlinePrimaryVerticesToken_= consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("vertices"));
27 }
28 
30 
32  iSetup.get<TrackerDigiGeometryRecord>().get(tracker);
33  assert(tracker.isValid());
34 
36 
37  if(applyVertexCut_) {
38  iEvent.getByToken(offlinePrimaryVerticesToken_, vertices);
39  if (!vertices.isValid() || vertices->empty()) return;
40  }
41 
42  std::vector<TrackerValidationVariables::AVTrackStruct> vtracks;
43 
44  validator.fillTrackQuantities(iEvent, iSetup,
45  // tell the validator to only look at good tracks
46  [&](const reco::Track& track) -> bool {
47  return (!applyVertexCut_ || (track.pt() > 0.75
48  && std::abs( track.dxy(vertices->at(0).position()) ) < 5*track.dxyError())) ;
49  }, vtracks);
50 
51  for (auto& track : vtracks) {
52  for (auto& it : track.hits) {
53  auto id = DetId(it.rawDetId);
54  auto isPixel = id.subdetId() == 1 || id.subdetId() == 2;
55  if (!isPixel) continue;
56 
57  histo[RESIDUAL_X].fill(it.resXprime, id, &iEvent);
58  histo[RESIDUAL_Y].fill(it.resYprime, id, &iEvent);
59  }
60  }
61 
62 }
63 
65 
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:460
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)
TrackerValidationVariables validator
Definition: DetId.h:18
void fillTrackQuantities(const edm::Event &, const edm::EventSetup &, std::vector< AVTrackStruct > &v_avtrackout)
edm::EDGetTokenT< reco::VertexCollection > offlinePrimaryVerticesToken_
const T & get() const
Definition: EventSetup.h:56
std::vector< HistogramManager > histo
bool isPixel(HitType hitType)
bool isValid() const
Definition: ESHandle.h:47
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