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 
20 
21 namespace {
22 
23  class SiPixelPhase1TrackResiduals final : public SiPixelPhase1Base {
24  enum { RESIDUAL_X, RESIDUAL_Y, RESONEDGE_X, RESONEDGE_Y, RESOTHERBAD_X, RESOTHERBAD_Y };
25 
26  public:
27  explicit SiPixelPhase1TrackResiduals(const edm::ParameterSet& conf);
28  void analyze(const edm::Event&, const edm::EventSetup&) override;
29 
30  private:
32  edm::EDGetTokenT<reco::VertexCollection> offlinePrimaryVerticesToken_;
33 
35 
36  bool applyVertexCut_;
37  };
38 
39  SiPixelPhase1TrackResiduals::SiPixelPhase1TrackResiduals(const edm::ParameterSet& iConfig)
40  : SiPixelPhase1Base(iConfig), validator(iConfig, consumesCollector()) {
41  applyVertexCut_ = iConfig.getUntrackedParameter<bool>("VertexCut", true);
42 
43  offlinePrimaryVerticesToken_ = consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("vertices"));
44 
45  trackerGeomToken_ = esConsumes<TrackerGeometry, TrackerDigiGeometryRecord>();
46  }
47 
49  if (!checktrigger(iEvent, iSetup, DCS))
50  return;
51 
52  edm::ESHandle<TrackerGeometry> tracker = iSetup.getHandle(trackerGeomToken_);
53  assert(tracker.isValid());
54 
56  if (applyVertexCut_) {
57  iEvent.getByToken(offlinePrimaryVerticesToken_, vertices);
58  if (!vertices.isValid() || vertices->empty())
59  return;
60  }
61 
62  std::vector<TrackerValidationVariables::AVTrackStruct> vtracks;
63 
64  validator.fillTrackQuantities(
65  iEvent,
66  iSetup,
67  // tell the validator to only look at good tracks
68  [&](const reco::Track& track) -> bool {
69  return (!applyVertexCut_ ||
70  (track.pt() > 0.75 && std::abs(track.dxy(vertices->at(0).position())) < 5 * track.dxyError()));
71  },
72  vtracks);
73 
74  for (auto& track : vtracks) {
75  for (auto& it : track.hits) {
76  auto id = DetId(it.rawDetId);
77  auto isPixel = id.subdetId() == 1 || id.subdetId() == 2;
78  if (!isPixel)
79  continue;
80 
81  histo[RESIDUAL_X].fill(it.resXprime, id, &iEvent);
82  histo[RESIDUAL_Y].fill(it.resYprime, id, &iEvent);
83 
84  if (it.isOnEdgePixel) {
85  histo[RESONEDGE_X].fill(it.resXprime, id, &iEvent);
86  histo[RESONEDGE_Y].fill(it.resYprime, id, &iEvent);
87  }
88 
89  if (it.isOtherBadPixel) {
90  histo[RESOTHERBAD_X].fill(it.resXprime, id, &iEvent);
91  histo[RESOTHERBAD_Y].fill(it.resYprime, id, &iEvent);
92  }
93  }
94  }
95  }
96 
97 } // namespace
98 
99 DEFINE_FWK_MODULE(SiPixelPhase1TrackResiduals);
TrackerValidationVariables.h
HLT_FULL_cff.track
track
Definition: HLT_FULL_cff.py:11713
TrackerGeometry.h
PixelTopology.h
edm::EDGetTokenT< reco::VertexCollection >
timingPdfMaker.histo
histo
Definition: timingPdfMaker.py:279
cms::cuda::assert
assert(be >=bs)
TrackerValidationVariables
Definition: TrackerValidationVariables.h:22
edm::ParameterSet::getUntrackedParameter
T getUntrackedParameter(std::string const &, T const &) const
edm::Handle< reco::VertexCollection >
SiPixelPhase1Base
Definition: SiPixelPhase1Base.h:46
DetId
Definition: DetId.h:17
MakerMacros.h
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
reco::Track
Definition: Track.h:27
edm::ESHandle< TrackerGeometry >
SiPixelRecHit.h
Vertex.h
PbPb_ZMuSkimMuonDPG_cff.tracker
tracker
Definition: PbPb_ZMuSkimMuonDPG_cff.py:60
TrackerDigiGeometryRecord.h
edm::ParameterSet
Definition: ParameterSet.h:47
fastTrackerRecHitType::isPixel
bool isPixel(HitType hitType)
Definition: FastTrackerRecHit.h:37
CertificationClient_cfi.DCS
DCS
Definition: CertificationClient_cfi.py:9
iEvent
int iEvent
Definition: GenABIO.cc:224
edm::EventSetup::getHandle
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:155
analyze
example_stream void analyze(const edm::Event &, const edm::EventSetup &) override
edm::EventSetup
Definition: EventSetup.h:58
edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord >
SiPixelPhase1Base::analyze
void analyze(edm::Event const &e, edm::EventSetup const &) override=0
VertexFwd.h
SiPixelPhase1Base.h
PixelGeomDetUnit.h
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
edm::Event
Definition: Event.h:73
edm::InputTag
Definition: InputTag.h:15
pwdgSkimBPark_cfi.vertices
vertices
Definition: pwdgSkimBPark_cfi.py:7