CMS 3D CMS Logo

SiPixelPhase1RecHitsV.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: SiPixelPhase1RecHitsV
4 // Class: SiPixelPhase1RecHitsV
5 //
6 
7 // Original Author: Marcel Schneider
8 // Additional Authors: Alexander Morton - modifying code for validation use
9 
12 
15 
18 
20  SiPixelPhase1Base(iConfig),
21  trackerHitAssociatorConfig_( iConfig, consumesCollector() ),
22  srcToken_ ( consumes<SiPixelRecHitCollection>(iConfig.getParameter<edm::InputTag>("src")) )
23 {}
24 
27  iEvent.getByToken(srcToken_, input);
28  if (!input.isValid()) return;
29 
31 
33  for (it = input->begin(); it != input->end(); ++it) {
34  auto id = DetId(it->detId());
35 
36  for(SiPixelRecHit const& rechit : *it) {
37  SiPixelRecHit::ClusterRef const& clust = rechit.cluster();
38 
39  std::vector<PSimHit> associateSimHit;
40  associateSimHit = associate.associateHit(rechit);
41  std::vector<PSimHit>::const_iterator closestIt = associateSimHit.begin();
42 
43  LocalPoint lp = rechit.localPosition();
44  float rechit_x = lp.x();
45  float rechit_y = lp.y();
46 
47  LocalError lerr = rechit.localPositionError();
48  float lerr_x = sqrt(lerr.xx());
49  float lerr_y = sqrt(lerr.yy());
50 
51  // loop over associated sim hits and find the closest
52  if ( !associateSimHit.empty() ) {
53  float closestSimHit = 9999.9;
54 
55  for (std::vector<PSimHit>::const_iterator m = associateSimHit.begin(); m < associateSimHit.end(); m++) {
56  float sim_x1 ( (*m).entryPoint().x() ), sim_x2 ( (*m).exitPoint().x() ), sim_xpos ( 0.5*(sim_x1+sim_x2) );
57  float sim_y1 ( (*m).entryPoint().y() ), sim_y2 ( (*m).exitPoint().y() ), sim_ypos ( 0.5*(sim_y1+sim_y2) );
58 
59  float xres ( std::abs(sim_xpos - rechit_x) ), yres ( std::abs(sim_ypos - rechit_y) );
60  float dist = std::sqrt( xres*xres + yres*yres );
61 
62  if ( dist < closestSimHit ) {
63  closestSimHit = dist;
64  closestIt = m;
65  }
66  }
67  }
68 
69  // Sim Hit stuff
70 
71  if ( !associateSimHit.empty() ) {
72 
73  const PSimHit& simHit = *closestIt;
74 
75  int bunch = simHit.eventId().bunchCrossing();
76  int event = simHit.eventId().event();
77 
78  float sim_x1 ( simHit.entryPoint().x() ), sim_x2 ( simHit.exitPoint().x() ), sim_xpos ( 0.5*(sim_x1 + sim_x2) );
79  float sim_y1 ( simHit.entryPoint().y() ), sim_y2 ( simHit.exitPoint().y() ), sim_ypos ( 0.5*(sim_y1 + sim_y2) );
80 
81  float res_x = (rechit_x - sim_xpos) * 10000.0;
82  float res_y = (rechit_y - sim_ypos) * 10000.0;
83 
84  float pull_x = ( rechit_x - sim_xpos ) / lerr_x;
85  float pull_y = ( rechit_y - sim_ypos ) / lerr_y;
86 
87  // Now Plotting stuff
88  if ( bunch == 0 ) histo[IN_TIME_BUNCH].fill(bunch, id, &iEvent);
89  if ( bunch != 0 ) histo[OUT_TIME_BUNCH].fill(bunch, id, &iEvent);
90 
91  histo[NSIMHITS].fill(event, id, &iEvent);
92 
93  histo[RECHIT_X].fill(rechit_x, id, &iEvent);
94  histo[RECHIT_Y].fill(rechit_y, id, &iEvent);
95 
96  histo[RES_X].fill(res_x, id, &iEvent);
97  histo[RES_Y].fill(res_y, id, &iEvent);
98 
99  histo[ERROR_X].fill(lerr_x, id, &iEvent);
100  histo[ERROR_Y].fill(lerr_y, id, &iEvent);
101 
102  histo[PULL_X].fill(pull_x, id, &iEvent);
103  histo[PULL_Y].fill(pull_y, id, &iEvent);
104  }
105  }
106  }
107 }
108 
110 
float xx() const
Definition: LocalError.h:24
boost::transform_iterator< IterHelp, const_IdIter > const_iterator
const_iterator end(bool update=false) const
int event() const
get the contents of the subdetector field (should be protected?)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:460
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
T y() const
Definition: PV3DBase.h:63
static std::string const input
Definition: EdmProvDump.cc:44
int iEvent
Definition: GenABIO.cc:230
Local3DPoint exitPoint() const
Exit point in the local Det frame.
Definition: PSimHit.h:38
float yy() const
Definition: LocalError.h:26
T sqrt(T t)
Definition: SSEVec.h:18
int bunchCrossing() const
get the detector field from this detid
SiPixelPhase1RecHitsV(const edm::ParameterSet &conf)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
edm::EDGetTokenT< SiPixelRecHitCollection > srcToken_
EncodedEventId eventId() const
Definition: PSimHit.h:105
bool isValid() const
Definition: HandleBase.h:74
TrackerHitAssociator::Config trackerHitAssociatorConfig_
Definition: DetId.h:18
void analyze(const edm::Event &, const edm::EventSetup &)
std::vector< HistogramManager > histo
HLT enums.
std::vector< PSimHit > associateHit(const TrackingRecHit &thit) const
T x() const
Definition: PV3DBase.h:62
Local3DPoint entryPoint() const
Entry point in the local Det frame.
Definition: PSimHit.h:35
const_iterator begin(bool update=false) const
Definition: event.py:1
Our base class.
Definition: SiPixelRecHit.h:23