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  offlinePrimaryVerticesToken_ = consumes<reco::VertexCollection>(std::string("offlinePrimaryVertices"));
26 
27  applyVertexCut_=iConfig.getUntrackedParameter<bool>("VertexCut",true);
28 }
29 
31 
33  iSetup.get<TrackerDigiGeometryRecord>().get(tracker);
34  assert(tracker.isValid());
35 
37  iEvent.getByToken(offlinePrimaryVerticesToken_, vertices);
38 
39  if (applyVertexCut_ && (!vertices.isValid() || vertices->size() == 0)) return;
40 
41  std::vector<TrackerValidationVariables::AVTrackStruct> vtracks;
42 
43  validator.fillTrackQuantities(iEvent, iSetup,
44  // tell the validator to only look at good tracks
45  [&](const reco::Track& track) -> bool {
46  return (!applyVertexCut_ || (track.pt() > 0.75
47  && std::abs( track.dxy(vertices->at(0).position()) ) < 5*track.dxyError())) ;
48  }, vtracks);
49 
50  for (auto& track : vtracks) {
51  for (auto& it : track.hits) {
52  auto id = DetId(it.rawDetId);
53  auto isPixel = id.subdetId() == 1 || id.subdetId() == 2;
54  if (!isPixel) continue;
55 
56  //TO BE UPDATED WITH VINCENZO STUFF
57  /*
58  const PixelGeomDetUnit* geomdetunit = dynamic_cast<const PixelGeomDetUnit*> ( tracker->idToDet(id) );
59  const PixelTopology& topol = geomdetunit->specificTopology();
60 
61  float lpx=it.localX;
62  float lpy=it.localY;
63  LocalPoint lp(lpx,lpy);
64 
65  MeasurementPoint mp = topol.measurementPosition(lp);
66  int row = (int) mp.x();
67  int col = (int) mp.y();
68  */
69 
70  histo[RESIDUAL_X].fill(it.resXprime, id, &iEvent);
71  histo[RESIDUAL_Y].fill(it.resYprime, id, &iEvent);
72  }
73  }
74 
75 }
76 
78 
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:796
#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: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)
edm::EDGetTokenT< reco::VertexCollection > offlinePrimaryVerticesToken_
const T & get() const
Definition: EventSetup.h:55
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