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 
13 
22 
23 namespace {
24 
25 class SiPixelPhase1TrackResiduals final : public SiPixelPhase1Base {
26  enum {
27  RESIDUAL_X,
28  RESIDUAL_Y
29  };
30 
31  public:
32  explicit SiPixelPhase1TrackResiduals(const edm::ParameterSet& conf);
33  void analyze(const edm::Event&, const edm::EventSetup&) override;
34 
35  private:
37  edm::EDGetTokenT<reco::VertexCollection> offlinePrimaryVerticesToken_;
38 
39  bool applyVertexCut_;
40 };
41 
42 SiPixelPhase1TrackResiduals::SiPixelPhase1TrackResiduals(const edm::ParameterSet& iConfig) :
43  SiPixelPhase1Base(iConfig),
44  validator(iConfig, consumesCollector())
45 {
46 
47  applyVertexCut_=iConfig.getUntrackedParameter<bool>("VertexCut",true);
48 
49  offlinePrimaryVerticesToken_= consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("vertices"));
50 
51 }
52 
54  if( !checktrigger(iEvent,iSetup,DCS) ) return;
55 
57  iSetup.get<TrackerDigiGeometryRecord>().get(tracker);
58  assert(tracker.isValid());
59 
61  if(applyVertexCut_) {
62  iEvent.getByToken(offlinePrimaryVerticesToken_, vertices);
63  if (!vertices.isValid() || vertices->empty()) return;
64  }
65 
66  std::vector<TrackerValidationVariables::AVTrackStruct> vtracks;
67 
68  validator.fillTrackQuantities(iEvent, iSetup,
69  // tell the validator to only look at good tracks
70  [&](const reco::Track& track) -> bool {
71  return (!applyVertexCut_ || (track.pt() > 0.75
72  && std::abs( track.dxy(vertices->at(0).position()) ) < 5*track.dxyError())) ;
73  }, vtracks);
74 
75  for (auto& track : vtracks) {
76  for (auto& it : track.hits) {
77  auto id = DetId(it.rawDetId);
78  auto isPixel = id.subdetId() == 1 || id.subdetId() == 2;
79  if (!isPixel) continue;
80 
81  histo[RESIDUAL_X].fill(it.resXprime, id, &iEvent);
82  histo[RESIDUAL_Y].fill(it.resYprime, id, &iEvent);
83  }
84  }
85 
86 }
87 
88 } // namespace
89 
90 DEFINE_FWK_MODULE(SiPixelPhase1TrackResiduals);
91 
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
virtual example_stream void analyze(const edm::Event &, const edm::EventSetup &) override
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
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
void analyze(edm::Event const &e, edm::EventSetup const &eSetup) override=0
Definition: DetId.h:18
const T & get() const
Definition: EventSetup.h:59
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