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 {
25  RESIDUAL_X,
26  RESIDUAL_Y,
27  RESONEDGE_X,
28  RESONEDGE_Y,
29  RESOTHERBAD_X,
30  RESOTHERBAD_Y,
31  NormRes_X,
32  NormRes_Y,
33  DRNR_X,
34  DRNR_Y
35  };
36 
37  public:
38  explicit SiPixelPhase1TrackResiduals(const edm::ParameterSet& conf);
39  void analyze(const edm::Event&, const edm::EventSetup&) override;
40 
41  private:
43  edm::EDGetTokenT<reco::VertexCollection> offlinePrimaryVerticesToken_;
44 
46 
47  bool applyVertexCut_;
48  };
49 
50  SiPixelPhase1TrackResiduals::SiPixelPhase1TrackResiduals(const edm::ParameterSet& iConfig)
51  : SiPixelPhase1Base(iConfig), validator(iConfig, consumesCollector()) {
52  applyVertexCut_ = iConfig.getUntrackedParameter<bool>("VertexCut", true);
53 
54  offlinePrimaryVerticesToken_ = consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("vertices"));
55 
56  trackerGeomToken_ = esConsumes<TrackerGeometry, TrackerDigiGeometryRecord>();
57  }
58 
60  if (!checktrigger(iEvent, iSetup, DCS))
61  return;
62 
63  edm::ESHandle<TrackerGeometry> tracker = iSetup.getHandle(trackerGeomToken_);
64  assert(tracker.isValid());
65 
67  if (applyVertexCut_) {
68  iEvent.getByToken(offlinePrimaryVerticesToken_, vertices);
69  if (!vertices.isValid() || vertices->empty())
70  return;
71  }
72 
73  std::vector<TrackerValidationVariables::AVTrackStruct> vtracks;
74 
75  validator.fillTrackQuantities(
76  iEvent,
77  iSetup,
78  // tell the validator to only look at good tracks
79  [&](const reco::Track& track) -> bool {
80  return (!applyVertexCut_ ||
81  (track.pt() > 0.75 && std::abs(track.dxy(vertices->at(0).position())) < 5 * track.dxyError()));
82  },
83  vtracks);
84 
85  for (auto& track : vtracks) {
86  for (auto& it : track.hits) {
87  auto id = DetId(it.rawDetId);
88  auto isPixel = id.subdetId() == 1 || id.subdetId() == 2;
89  if (!isPixel)
90  continue;
91 
92  histo[RESIDUAL_X].fill(it.resXprime, id, &iEvent);
93  histo[RESIDUAL_Y].fill(it.resYprime, id, &iEvent);
94 
95  if (it.resXprimeErr != 0 && it.resYprimeErr != 0) {
96  histo[NormRes_X].fill(it.resXprime / it.resXprimeErr, id, &iEvent);
97  histo[NormRes_Y].fill(it.resYprime / it.resYprimeErr, id, &iEvent);
98  histo[DRNR_X].fill(it.resXprime / it.resXprimeErr, id, &iEvent);
99  histo[DRNR_Y].fill(it.resYprime / it.resYprimeErr, id, &iEvent);
100  }
101 
102  if (it.isOnEdgePixel) {
103  histo[RESONEDGE_X].fill(it.resXprime, id, &iEvent);
104  histo[RESONEDGE_Y].fill(it.resYprime, id, &iEvent);
105  }
106 
107  if (it.isOtherBadPixel) {
108  histo[RESOTHERBAD_X].fill(it.resXprime, id, &iEvent);
109  histo[RESOTHERBAD_Y].fill(it.resYprime, id, &iEvent);
110  }
111  }
112  }
113  }
114 
115 } // namespace
116 
117 DEFINE_FWK_MODULE(SiPixelPhase1TrackResiduals);
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
assert(be >=bs)
example_stream void analyze(const edm::Event &, const edm::EventSetup &) override
T getUntrackedParameter(std::string const &, T const &) const
int iEvent
Definition: GenABIO.cc:224
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:130
Definition: DetId.h:17
void analyze(edm::Event const &e, edm::EventSetup const &) override=0
bool isPixel(HitType hitType)