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 
26  applyVertexCut_=iConfig.getUntrackedParameter<bool>("VertexCut",true);
27 
28  offlinePrimaryVerticesToken_= consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("vertices"));
29 
30 }
31 
33  if( !checktrigger(iEvent,iSetup,DCS) ) return;
34 
36  iSetup.get<TrackerDigiGeometryRecord>().get(tracker);
37  assert(tracker.isValid());
38 
40  if(applyVertexCut_) {
41  iEvent.getByToken(offlinePrimaryVerticesToken_, vertices);
42  if (!vertices.isValid() || vertices->empty()) return;
43  }
44 
45  std::vector<TrackerValidationVariables::AVTrackStruct> vtracks;
46 
47  validator.fillTrackQuantities(iEvent, iSetup,
48  // tell the validator to only look at good tracks
49  [&](const reco::Track& track) -> bool {
50  return (!applyVertexCut_ || (track.pt() > 0.75
51  && std::abs( track.dxy(vertices->at(0).position()) ) < 5*track.dxyError())) ;
52  }, vtracks);
53 
54  for (auto& track : vtracks) {
55  for (auto& it : track.hits) {
56  auto id = DetId(it.rawDetId);
57  auto isPixel = id.subdetId() == 1 || id.subdetId() == 2;
58  if (!isPixel) continue;
59 
60  histo[RESIDUAL_X].fill(it.resXprime, id, &iEvent);
61  histo[RESIDUAL_Y].fill(it.resYprime, id, &iEvent);
62  }
63  }
64 
65 }
66 
68 
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:519
double dxyError() const
error on dxy
Definition: TrackBase.h:796
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
bool checktrigger(const edm::Event &iEvent, const edm::EventSetup &iSetup, const unsigned trgidx) const
int iEvent
Definition: GenABIO.cc:230
double pt() const
track transverse momentum
Definition: TrackBase.h:621
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)
void analyze(const edm::Event &, const edm::EventSetup &) override
edm::EDGetTokenT< reco::VertexCollection > offlinePrimaryVerticesToken_
const T & get() const
Definition: EventSetup.h:58
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:591