CMS 3D CMS Logo

All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
List of all members | Public Member Functions | Protected Member Functions | Private Member Functions | Private Attributes
SiPixelRecHitsValid Class Reference

#include <SiPixelRecHitsValid.h>

Inheritance diagram for SiPixelRecHitsValid:
DQMEDAnalyzer edm::stream::EDProducer< edm::GlobalCache< DQMEDAnalyzerGlobalCache >, edm::EndRunProducer, edm::EndLuminosityBlockProducer, edm::Accumulator >

Public Member Functions

 SiPixelRecHitsValid (const edm::ParameterSet &conf)
 
 ~SiPixelRecHitsValid () override
 
- Public Member Functions inherited from DQMEDAnalyzer
void accumulate (edm::Event const &event, edm::EventSetup const &setup) final
 
void beginLuminosityBlock (edm::LuminosityBlock const &lumi, edm::EventSetup const &setup) final
 
void beginRun (edm::Run const &run, edm::EventSetup const &setup) final
 
void beginStream (edm::StreamID id) final
 
virtual void dqmBeginRun (edm::Run const &, edm::EventSetup const &)
 
 DQMEDAnalyzer ()
 
void endLuminosityBlock (edm::LuminosityBlock const &lumi, edm::EventSetup const &setup) final
 
void endRun (edm::Run const &run, edm::EventSetup const &setup) final
 
virtual bool getCanSaveByLumi ()
 
- Public Member Functions inherited from edm::stream::EDProducer< edm::GlobalCache< DQMEDAnalyzerGlobalCache >, edm::EndRunProducer, edm::EndLuminosityBlockProducer, edm::Accumulator >
 EDProducer ()=default
 
 EDProducer (const EDProducer &)=delete
 
bool hasAbilityToProduceInBeginLumis () const final
 
bool hasAbilityToProduceInBeginProcessBlocks () const final
 
bool hasAbilityToProduceInBeginRuns () const final
 
bool hasAbilityToProduceInEndLumis () const final
 
bool hasAbilityToProduceInEndProcessBlocks () const final
 
bool hasAbilityToProduceInEndRuns () const final
 
const EDProduceroperator= (const EDProducer &)=delete
 

Protected Member Functions

void analyze (const edm::Event &e, const edm::EventSetup &c) override
 
void bookHistograms (DQMStore::IBooker &ibooker, const edm::Run &run, const edm::EventSetup &es) override
 
- Protected Member Functions inherited from DQMEDAnalyzer
uint64_t meId () const
 

Private Member Functions

void fillBarrel (const SiPixelRecHit &, const PSimHit &, DetId, const PixelGeomDetUnit *, const TrackerTopology *tTopo)
 
void fillForward (const SiPixelRecHit &, const PSimHit &, DetId, const PixelGeomDetUnit *, const TrackerTopology *tTopo)
 

Private Attributes

MonitorElementclustChargeDisk1Plaquettes [7]
 
MonitorElementclustChargeDisk2Plaquettes [7]
 
MonitorElementclustChargeLayer1Modules [8]
 
MonitorElementclustChargeLayer2Modules [8]
 
MonitorElementclustChargeLayer3Modules [8]
 
MonitorElementclustXSizeDisk1Plaquettes [7]
 
MonitorElementclustXSizeDisk2Plaquettes [7]
 
MonitorElementclustXSizeLayer [3]
 
MonitorElementclustYSizeDisk1Plaquettes [7]
 
MonitorElementclustYSizeDisk2Plaquettes [7]
 
MonitorElementclustYSizeModule [8]
 
MonitorElementrecHitBunchB
 
MonitorElementrecHitBunchF
 
MonitorElementrecHitEventB
 
MonitorElementrecHitEventF
 
MonitorElementrecHitNsimHitDisk1
 
MonitorElementrecHitNsimHitDisk2
 
MonitorElementrecHitNsimHitLayer [3]
 
MonitorElementrecHitXFullModules
 
MonitorElementrecHitXHalfModules
 
MonitorElementrecHitXPlaquetteSize1
 
MonitorElementrecHitXPlaquetteSize2
 
MonitorElementrecHitXPullAllB
 
MonitorElementrecHitXPullAllF
 
MonitorElementrecHitXPullDisk1Plaquettes [7]
 
MonitorElementrecHitXPullDisk2Plaquettes [7]
 
MonitorElementrecHitXPullFlippedLadderLayers [3]
 
MonitorElementrecHitXPullNonFlippedLadderLayers [3]
 
MonitorElementrecHitXResAllB
 
MonitorElementrecHitXResAllF
 
MonitorElementrecHitXResDisk1Plaquettes [7]
 
MonitorElementrecHitXResDisk2Plaquettes [7]
 
MonitorElementrecHitXResFlippedLadderLayers [3]
 
MonitorElementrecHitXResNonFlippedLadderLayers [3]
 
MonitorElementrecHitYAllModules
 
MonitorElementrecHitYPlaquetteSize2
 
MonitorElementrecHitYPlaquetteSize3
 
MonitorElementrecHitYPlaquetteSize4
 
MonitorElementrecHitYPlaquetteSize5
 
MonitorElementrecHitYPullAllB
 
MonitorElementrecHitYPullAllF
 
MonitorElementrecHitYPullDisk1Plaquettes [7]
 
MonitorElementrecHitYPullDisk2Plaquettes [7]
 
MonitorElementrecHitYPullLayer1Modules [8]
 
MonitorElementrecHitYPullLayer2Modules [8]
 
MonitorElementrecHitYPullLayer3Modules [8]
 
MonitorElementrecHitYResAllB
 
MonitorElementrecHitYResAllF
 
MonitorElementrecHitYResDisk1Plaquettes [7]
 
MonitorElementrecHitYResDisk2Plaquettes [7]
 
MonitorElementrecHitYResLayer1Modules [8]
 
MonitorElementrecHitYResLayer2Modules [8]
 
MonitorElementrecHitYResLayer3Modules [8]
 
edm::EDGetTokenT< SiPixelRecHitCollectionsiPixelRecHitCollectionToken_
 
const edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecordtGeomEsToken
 
TrackerHitAssociator::Config trackerHitAssociatorConfig_
 
const edm::ESGetToken< TrackerTopology, TrackerTopologyRcdtTopoEsToken
 

Additional Inherited Members

- Public Types inherited from DQMEDAnalyzer
typedef dqm::reco::DQMStore DQMStore
 
typedef dqm::reco::MonitorElement MonitorElement
 
- Public Types inherited from edm::stream::EDProducer< edm::GlobalCache< DQMEDAnalyzerGlobalCache >, edm::EndRunProducer, edm::EndLuminosityBlockProducer, edm::Accumulator >
using CacheTypes = CacheContexts< T... >
 
using GlobalCache = typename CacheTypes::GlobalCache
 
using HasAbility = AbilityChecker< T... >
 
using InputProcessBlockCache = typename CacheTypes::InputProcessBlockCache
 
using LuminosityBlockCache = typename CacheTypes::LuminosityBlockCache
 
using LuminosityBlockContext = LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCache >
 
using LuminosityBlockSummaryCache = typename CacheTypes::LuminosityBlockSummaryCache
 
using RunCache = typename CacheTypes::RunCache
 
using RunContext = RunContextT< RunCache, GlobalCache >
 
using RunSummaryCache = typename CacheTypes::RunSummaryCache
 
- Static Public Member Functions inherited from DQMEDAnalyzer
static void globalEndJob (DQMEDAnalyzerGlobalCache const *)
 
static void globalEndLuminosityBlockProduce (edm::LuminosityBlock &lumi, edm::EventSetup const &setup, LuminosityBlockContext const *context)
 
static void globalEndRunProduce (edm::Run &run, edm::EventSetup const &setup, RunContext const *context)
 
static std::unique_ptr< DQMEDAnalyzerGlobalCacheinitializeGlobalCache (edm::ParameterSet const &)
 
- Protected Attributes inherited from DQMEDAnalyzer
edm::EDPutTokenT< DQMTokenlumiToken_
 
edm::EDPutTokenT< DQMTokenrunToken_
 
unsigned int streamId_
 

Detailed Description

File: SiPixelRecHitsValid.h

Author
Jason Shaev, JHU Created: 6/7/06

Definition at line 29 of file SiPixelRecHitsValid.h.

Constructor & Destructor Documentation

◆ SiPixelRecHitsValid()

SiPixelRecHitsValid::SiPixelRecHitsValid ( const edm::ParameterSet conf)

Definition at line 46 of file SiPixelRecHitsValid.cc.

49  trackerHitAssociatorConfig_(ps, consumesCollector()),
50  siPixelRecHitCollectionToken_(consumes<SiPixelRecHitCollection>(ps.getParameter<edm::InputTag>("src"))) {}
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
const edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > tGeomEsToken
TrackerHitAssociator::Config trackerHitAssociatorConfig_
edm::EDGetTokenT< SiPixelRecHitCollection > siPixelRecHitCollectionToken_
const edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > tTopoEsToken

◆ ~SiPixelRecHitsValid()

SiPixelRecHitsValid::~SiPixelRecHitsValid ( )
overridedefault

Member Function Documentation

◆ analyze()

void SiPixelRecHitsValid::analyze ( const edm::Event e,
const edm::EventSetup c 
)
overrideprotectedvirtual

Reimplemented from DQMEDAnalyzer.

Definition at line 275 of file SiPixelRecHitsValid.cc.

References TrackerHitAssociator::associateHit(), edmNew::DetSet< T >::begin(), gather_cfg::cout, TrackerGeometry::dets(), MillePedeFileConverter_cfg::e, edmNew::DetSet< T >::end(), edmNew::DetSetVector< T >::end(), dqm::impl::MonitorElement::Fill(), HcalObjRepresent::Fill(), fillBarrel(), fillForward(), edmNew::DetSetVector< T >::find(), edm::EventSetup::getData(), mps_fire::i, TrackerGeometry::idToDet(), visualization-live-secondInstance_cfg::m, muonTagProbeFilters_cff::matched, TrackerTopology::pxbLayer(), TrackerTopology::pxfDisk(), recHitNsimHitDisk1, recHitNsimHitDisk2, recHitNsimHitLayer, siPixelRecHitCollectionToken_, mathSSE::sqrt(), DetId::subdetId(), tGeomEsToken, trackerHitAssociatorConfig_, tTopoEsToken, PV3DBase< T, PVType, FrameType >::x(), and PV3DBase< T, PVType, FrameType >::y().

275  {
276  //Retrieve tracker topology from geometry
277  const TrackerTopology* tTopo = &es.getData(tTopoEsToken);
278 
279  edm::LogInfo("EventInfo") << " Run = " << e.id().run() << " Event = " << e.id().event();
280  if (e.id().event() % 1000 == 0)
281  std::cout << " Run = " << e.id().run() << " Event = " << e.id().event() << std::endl;
282 
283  //Get RecHits
285  e.getByToken(siPixelRecHitCollectionToken_, recHitColl);
286 
287  //Get event setup
288  const TrackerGeometry* theTracker = &es.getData(tGeomEsToken);
289 
291 
292  //iterate over detunits
293  //for (TrackerGeometry::DetContainer::const_iterator it = geom->dets().begin(); it != geom->dets().end(); it++) {
294  for (const auto& it : theTracker->dets()) {
295  DetId detId = it->geographicalId();
296  unsigned int subid = detId.subdetId();
297 
298  if (!((subid == 1) || (subid == 2)))
299  continue;
300 
301  const PixelGeomDetUnit* theGeomDet = dynamic_cast<const PixelGeomDetUnit*>(theTracker->idToDet(detId));
302 
303  SiPixelRecHitCollection::const_iterator pixeldet = recHitColl->find(detId);
304  if (pixeldet == recHitColl->end())
305  continue;
306  SiPixelRecHitCollection::DetSet pixelrechitRange = *pixeldet;
307  SiPixelRecHitCollection::DetSet::const_iterator pixelrechitRangeIteratorBegin = pixelrechitRange.begin();
308  SiPixelRecHitCollection::DetSet::const_iterator pixelrechitRangeIteratorEnd = pixelrechitRange.end();
309  SiPixelRecHitCollection::DetSet::const_iterator pixeliter = pixelrechitRangeIteratorBegin;
310  std::vector<PSimHit> matched;
311 
312  //----Loop over rechits for this detId
313  for (; pixeliter != pixelrechitRangeIteratorEnd; pixeliter++) {
314  matched.clear();
315  matched = associate.associateHit(*pixeliter);
316 
317  if (!matched.empty()) {
318  float closest = 9999.9;
319  std::vector<PSimHit>::const_iterator closestit = matched.begin();
320  LocalPoint lp = pixeliter->localPosition();
321  float rechit_x = lp.x();
322  float rechit_y = lp.y();
323 
324  //loop over sim hits and fill closet
325  for (std::vector<PSimHit>::const_iterator m = matched.begin(); m < matched.end(); m++) {
326  float sim_x1 = (*m).entryPoint().x();
327  float sim_x2 = (*m).exitPoint().x();
328  float sim_xpos = 0.5 * (sim_x1 + sim_x2);
329 
330  float sim_y1 = (*m).entryPoint().y();
331  float sim_y2 = (*m).exitPoint().y();
332  float sim_ypos = 0.5 * (sim_y1 + sim_y2);
333 
334  float x_res = fabs(sim_xpos - rechit_x);
335  float y_res = fabs(sim_ypos - rechit_y);
336 
337  float dist = sqrt(x_res * x_res + y_res * y_res);
338 
339  if (dist < closest) {
340  closest = dist;
341  closestit = m;
342  }
343  } // end sim hit loop
344 
345  if (subid == 1) { //<----------barrel
346  fillBarrel(*pixeliter, *closestit, detId, theGeomDet, tTopo);
347  } // end barrel
348  if (subid == 2) { // <-------forward
349  fillForward(*pixeliter, *closestit, detId, theGeomDet, tTopo);
350  }
351 
352  } // end matched emtpy
353 
354  int NsimHit = matched.size();
355  if (subid == 1) { //<----------barrel
356  for (unsigned int i = 0; i < 3; i++)
357  if (tTopo->pxbLayer(detId) == i + 1)
358  recHitNsimHitLayer[i]->Fill(NsimHit);
359  } // end barrel
360  if (subid == 2) { // <-------forward
361  if (tTopo->pxfDisk(detId) == 1)
362  recHitNsimHitDisk1->Fill(NsimHit);
363  else
364  recHitNsimHitDisk2->Fill(NsimHit);
365  }
366  } // <-----end rechit loop
367  } // <------ end detunit loop
368 }
MonitorElement * recHitNsimHitDisk1
unsigned int pxbLayer(const DetId &id) const
const edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > tGeomEsToken
MonitorElement * recHitNsimHitLayer[3]
const DetContainer & dets() const override
Returm a vector of all GeomDet (including all GeomDetUnits)
const_iterator end(bool update=false) const
void Fill(long long x)
T x() const
Definition: PV3DBase.h:59
T y() const
Definition: PV3DBase.h:60
TrackerHitAssociator::Config trackerHitAssociatorConfig_
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
T sqrt(T t)
Definition: SSEVec.h:19
edm::EDGetTokenT< SiPixelRecHitCollection > siPixelRecHitCollectionToken_
unsigned int pxfDisk(const DetId &id) const
void fillForward(const SiPixelRecHit &, const PSimHit &, DetId, const PixelGeomDetUnit *, const TrackerTopology *tTopo)
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:48
const TrackerGeomDet * idToDet(DetId) const override
Log< level::Info, false > LogInfo
const edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > tTopoEsToken
boost::transform_iterator< IterHelp, const_IdIter > const_iterator
Definition: DetId.h:17
iterator end()
Definition: DetSetNew.h:56
const_iterator find(id_type i, bool update=false) const
MonitorElement * recHitNsimHitDisk2
void fillBarrel(const SiPixelRecHit &, const PSimHit &, DetId, const PixelGeomDetUnit *, const TrackerTopology *tTopo)
iterator begin()
Definition: DetSetNew.h:54

◆ bookHistograms()

void SiPixelRecHitsValid::bookHistograms ( DQMStore::IBooker ibooker,
const edm::Run run,
const edm::EventSetup es 
)
overrideprotectedvirtual

Implements DQMEDAnalyzer.

Definition at line 54 of file SiPixelRecHitsValid.cc.

References dqm::implementation::IBooker::book1D(), clustChargeDisk1Plaquettes, clustChargeDisk2Plaquettes, clustChargeLayer1Modules, clustChargeLayer2Modules, clustChargeLayer3Modules, clustXSizeDisk1Plaquettes, clustXSizeDisk2Plaquettes, clustXSizeLayer, clustYSizeDisk1Plaquettes, clustYSizeDisk2Plaquettes, clustYSizeModule, timingPdfMaker::histo, mps_fire::i, recHitBunchB, recHitBunchF, recHitEventB, recHitEventF, recHitNsimHitDisk1, recHitNsimHitDisk2, recHitNsimHitLayer, recHitXFullModules, recHitXHalfModules, recHitXPlaquetteSize1, recHitXPlaquetteSize2, recHitXPullAllB, recHitXPullAllF, recHitXPullDisk1Plaquettes, recHitXPullDisk2Plaquettes, recHitXPullFlippedLadderLayers, recHitXPullNonFlippedLadderLayers, recHitXResAllB, recHitXResAllF, recHitXResDisk1Plaquettes, recHitXResDisk2Plaquettes, recHitXResFlippedLadderLayers, recHitXResNonFlippedLadderLayers, recHitYAllModules, recHitYPlaquetteSize2, recHitYPlaquetteSize3, recHitYPlaquetteSize4, recHitYPlaquetteSize5, recHitYPullAllB, recHitYPullAllF, recHitYPullDisk1Plaquettes, recHitYPullDisk2Plaquettes, recHitYPullLayer1Modules, recHitYPullLayer2Modules, recHitYPullLayer3Modules, recHitYResAllB, recHitYResAllF, recHitYResDisk1Plaquettes, recHitYResDisk2Plaquettes, recHitYResLayer1Modules, recHitYResLayer2Modules, recHitYResLayer3Modules, and dqm::implementation::NavigatorBase::setCurrentFolder().

54  {
55  ibooker.setCurrentFolder("TrackerRecHitsV/TrackerRecHits/Pixel/clustBPIX");
56 
57  Char_t histo[200];
58 
59  // ---------------------------------------------------------------
60  // All histograms that depend on plaquette number have 7 indexes.
61  // The first 4 (0-3) correspond to Panel 1 plaquettes 1-4.
62  // The last 3 (4-6) correspond to Panel 2 plaquettes 1-3.
63  // ---------------------------------------------------------------
64 
65  //Cluster y-size by module number for barrel
66  for (int i = 0; i < 8; i++) {
67  sprintf(histo, "Clust_y_size_Module%d", i + 1);
68  clustYSizeModule[i] = ibooker.book1D(histo, "Cluster y-size by Module", 20, 0.5, 20.5);
69  } // end for
70 
71  //Cluster x-size by layer for barrel
72  for (int i = 0; i < 3; i++) {
73  sprintf(histo, "Clust_x_size_Layer%d", i + 1);
74  clustXSizeLayer[i] = ibooker.book1D(histo, "Cluster x-size by Layer", 20, 0.5, 20.5);
75  } // end for
76 
77  //Cluster charge by module for 3 layers of barrel
78  for (int i = 0; i < 8; i++) {
79  //Cluster charge by module for Layer1
80  sprintf(histo, "Clust_charge_Layer1_Module%d", i + 1);
81  clustChargeLayer1Modules[i] = ibooker.book1D(histo, "Cluster charge Layer 1 by Module", 50, 0., 200000.);
82 
83  //Cluster charge by module for Layer2
84  sprintf(histo, "Clust_charge_Layer2_Module%d", i + 1);
85  clustChargeLayer2Modules[i] = ibooker.book1D(histo, "Cluster charge Layer 2 by Module", 50, 0., 200000.);
86 
87  //Cluster charge by module for Layer3
88  sprintf(histo, "Clust_charge_Layer3_Module%d", i + 1);
89  clustChargeLayer3Modules[i] = ibooker.book1D(histo, "Cluster charge Layer 3 by Module", 50, 0., 200000.);
90  } // end for
91 
92  ibooker.setCurrentFolder("TrackerRecHitsV/TrackerRecHits/Pixel/clustFPIX");
93  //Cluster x-size, y-size, and charge by plaquette for Disks in Forward
94  for (int i = 0; i < 7; i++) {
95  //Cluster x-size for Disk1 by Plaquette
96  sprintf(histo, "Clust_x_size_Disk1_Plaquette%d", i + 1);
97  clustXSizeDisk1Plaquettes[i] = ibooker.book1D(histo, "Cluster X-size for Disk1 by Plaquette", 20, 0.5, 20.5);
98 
99  //Cluster x-size for Disk2 by Plaquette
100  sprintf(histo, "Clust_x_size_Disk2_Plaquette%d", i + 1);
101  clustXSizeDisk2Plaquettes[i] = ibooker.book1D(histo, "Cluster X-size for Disk2 by Plaquette", 20, 0.5, 20.5);
102 
103  //Cluster y-size for Disk1 by Plaquette
104  sprintf(histo, "Clust_y_size_Disk1_Plaquette%d", i + 1);
105  clustYSizeDisk1Plaquettes[i] = ibooker.book1D(histo, "Cluster Y-size for Disk1 by Plaquette", 20, 0.5, 20.5);
106 
107  //Cluster y-size for Disk2 by Plaquette
108  sprintf(histo, "Clust_y_size_Disk2_Plaquette%d", i + 1);
109  clustYSizeDisk2Plaquettes[i] = ibooker.book1D(histo, "Cluster Y-size for Disk2 by Plaquette", 20, 0.5, 20.5);
110 
111  //Cluster charge for Disk1 by Plaquette
112  sprintf(histo, "Clust_charge_Disk1_Plaquette%d", i + 1);
113  clustChargeDisk1Plaquettes[i] = ibooker.book1D(histo, "Cluster charge for Disk1 by Plaquette", 50, 0., 200000.);
114 
115  //Cluster charge for Disk2 by Plaquette
116  sprintf(histo, "Clust_charge_Disk2_Plaquette%d", i + 1);
117  clustChargeDisk2Plaquettes[i] = ibooker.book1D(histo, "Cluster charge for Disk2 by Plaquette", 50, 0., 200000.);
118  } // end for
119 
120  ibooker.setCurrentFolder("TrackerRecHitsV/TrackerRecHits/Pixel/recHitBPIX");
121  //RecHit Bunch crossing all barrel hits
122  recHitBunchB = ibooker.book1D("RecHit_Bunch_Barrel", "RecHit Bunch Crossing, Barrel", 20, -10., 10.);
123 
124  //RecHit Event, in-time bunch, all barrel hits
125  recHitEventB = ibooker.book1D("RecHit_Event_Barrel", "RecHit Event (in-time bunch), Barrel", 100, 0., 100.);
126 
127  //RecHit X Resolution all barrel hits
128  recHitXResAllB = ibooker.book1D("RecHit_xres_b_All", "RecHit X Res All Modules in Barrel", 100, -200., 200.);
129 
130  //RecHit Y Resolution all barrel hits
131  recHitYResAllB = ibooker.book1D("RecHit_yres_b_All", "RecHit Y Res All Modules in Barrel", 100, -200., 200.);
132 
133  //RecHit X distribution for full modules for barrel
134  recHitXFullModules = ibooker.book1D("RecHit_x_FullModules", "RecHit X distribution for full modules", 100, -2., 2.);
135 
136  //RecHit X distribution for half modules for barrel
137  recHitXHalfModules = ibooker.book1D("RecHit_x_HalfModules", "RecHit X distribution for half modules", 100, -1., 1.);
138 
139  //RecHit Y distribution all modules for barrel
140  recHitYAllModules = ibooker.book1D("RecHit_y_AllModules", "RecHit Y distribution for all modules", 100, -4., 4.);
141 
142  //RecHit X resolution for flipped and unflipped ladders by layer for barrel
143  for (int i = 0; i < 3; i++) {
144  //RecHit no. of matched simHits all ladders by layer
145  sprintf(histo, "RecHit_NsimHit_Layer%d", i + 1);
146  recHitNsimHitLayer[i] = ibooker.book1D(histo, "RecHit Number of simHits by Layer", 30, 0., 30.);
147 
148  //RecHit X resolution for flipped ladders by layer
149  sprintf(histo, "RecHit_XRes_FlippedLadder_Layer%d", i + 1);
150  recHitXResFlippedLadderLayers[i] = ibooker.book1D(histo, "RecHit XRes Flipped Ladders by Layer", 100, -200., 200.);
151 
152  //RecHit X resolution for unflipped ladders by layer
153  sprintf(histo, "RecHit_XRes_UnFlippedLadder_Layer%d", i + 1);
155  ibooker.book1D(histo, "RecHit XRes NonFlipped Ladders by Layer", 100, -200., 200.);
156  } // end for
157 
158  //RecHit Y resolutions for layers by module for barrel
159  for (int i = 0; i < 8; i++) {
160  //Rec Hit Y resolution by module for Layer1
161  sprintf(histo, "RecHit_YRes_Layer1_Module%d", i + 1);
162  recHitYResLayer1Modules[i] = ibooker.book1D(histo, "RecHit YRes Layer1 by module", 100, -200., 200.);
163 
164  //RecHit Y resolution by module for Layer2
165  sprintf(histo, "RecHit_YRes_Layer2_Module%d", i + 1);
166  recHitYResLayer2Modules[i] = ibooker.book1D(histo, "RecHit YRes Layer2 by module", 100, -200., 200.);
167 
168  //RecHit Y resolution by module for Layer3
169  sprintf(histo, "RecHit_YRes_Layer3_Module%d", i + 1);
170  recHitYResLayer3Modules[i] = ibooker.book1D(histo, "RecHit YRes Layer3 by module", 100, -200., 200.);
171  } // end for
172 
173  ibooker.setCurrentFolder("TrackerRecHitsV/TrackerRecHits/Pixel/recHitFPIX");
174  //RecHit Bunch crossing all plaquettes
175  recHitBunchF = ibooker.book1D("RecHit_Bunch_Forward", "RecHit Bunch Crossing, Forward", 20, -10., 10.);
176 
177  //RecHit Event, in-time bunch, all plaquettes
178  recHitEventF = ibooker.book1D("RecHit_Event_Forward", "RecHit Event (in-time bunch), Forward", 100, 0., 100.);
179 
180  //RecHit No. of simHits, by disk
181  recHitNsimHitDisk1 = ibooker.book1D("RecHit_NsimHit_Disk1", "RecHit Number of simHits, Disk1", 30, 0., 30.);
182  recHitNsimHitDisk2 = ibooker.book1D("RecHit_NsimHit_Disk2", "RecHit Number of simHits, Disk2", 30, 0., 30.);
183 
184  //RecHit X resolution all plaquettes
185  recHitXResAllF = ibooker.book1D("RecHit_xres_f_All", "RecHit X Res All in Forward", 100, -200., 200.);
186 
187  //RecHit Y resolution all plaquettes
188  recHitYResAllF = ibooker.book1D("RecHit_yres_f_All", "RecHit Y Res All in Forward", 100, -200., 200.);
189 
190  //RecHit X distribution for plaquette with x-size 1 in forward
192  ibooker.book1D("RecHit_x_Plaquette_xsize1", "RecHit X Distribution for plaquette x-size1", 100, -2., 2.);
193 
194  //RecHit X distribution for plaquette with x-size 2 in forward
196  ibooker.book1D("RecHit_x_Plaquette_xsize2", "RecHit X Distribution for plaquette x-size2", 100, -2., 2.);
197 
198  //RecHit Y distribution for plaquette with y-size 2 in forward
200  ibooker.book1D("RecHit_y_Plaquette_ysize2", "RecHit Y Distribution for plaquette y-size2", 100, -4., 4.);
201 
202  //RecHit Y distribution for plaquette with y-size 3 in forward
204  ibooker.book1D("RecHit_y_Plaquette_ysize3", "RecHit Y Distribution for plaquette y-size3", 100, -4., 4.);
205 
206  //RecHit Y distribution for plaquette with y-size 4 in forward
208  ibooker.book1D("RecHit_y_Plaquette_ysize4", "RecHit Y Distribution for plaquette y-size4", 100, -4., 4.);
209 
210  //RecHit Y distribution for plaquette with y-size 5 in forward
212  ibooker.book1D("RecHit_y_Plaquette_ysize5", "RecHit Y Distribution for plaquette y-size5", 100, -4., 4.);
213 
214  //X and Y resolutions for both disks by plaquette in forward
215  for (int i = 0; i < 7; i++) {
216  //X resolution for Disk1 by plaquette
217  sprintf(histo, "RecHit_XRes_Disk1_Plaquette%d", i + 1);
218  recHitXResDisk1Plaquettes[i] = ibooker.book1D(histo, "RecHit XRes Disk1 by plaquette", 100, -200., 200.);
219  //X resolution for Disk2 by plaquette
220  sprintf(histo, "RecHit_XRes_Disk2_Plaquette%d", i + 1);
221  recHitXResDisk2Plaquettes[i] = ibooker.book1D(histo, "RecHit XRes Disk2 by plaquette", 100, -200., 200.);
222 
223  //Y resolution for Disk1 by plaquette
224  sprintf(histo, "RecHit_YRes_Disk1_Plaquette%d", i + 1);
225  recHitYResDisk1Plaquettes[i] = ibooker.book1D(histo, "RecHit YRes Disk1 by plaquette", 100, -200., 200.);
226  //Y resolution for Disk2 by plaquette
227  sprintf(histo, "RecHit_YRes_Disk2_Plaquette%d", i + 1);
228  recHitYResDisk2Plaquettes[i] = ibooker.book1D(histo, "RecHit YRes Disk2 by plaquette", 100, -200., 200.);
229  }
230 
231  ibooker.setCurrentFolder("TrackerRecHitsV/TrackerRecHits/Pixel/recHitPullsBPIX");
232  recHitXPullAllB = ibooker.book1D("RecHit_xres_b_All", "RecHit X Pull All Modules in Barrel", 100, -10.0, 10.0);
233  recHitYPullAllB = ibooker.book1D("RecHit_yres_b_All", "RecHit Y Pull All Modules in Barrel", 100, -10.0, 10.0);
234 
235  for (int i = 0; i < 3; i++) {
236  sprintf(histo, "RecHit_XPull_FlippedLadder_Layer%d", i + 1);
238  ibooker.book1D(histo, "RecHit XPull Flipped Ladders by Layer", 100, -10.0, 10.0);
239 
240  sprintf(histo, "RecHit_XPull_UnFlippedLadder_Layer%d", i + 1);
242  ibooker.book1D(histo, "RecHit XPull NonFlipped Ladders by Layer", 100, -10.0, 10.0);
243  }
244 
245  for (int i = 0; i < 8; i++) {
246  sprintf(histo, "RecHit_YPull_Layer1_Module%d", i + 1);
247  recHitYPullLayer1Modules[i] = ibooker.book1D(histo, "RecHit YPull Layer1 by module", 100, -10.0, 10.0);
248 
249  sprintf(histo, "RecHit_YPull_Layer2_Module%d", i + 1);
250  recHitYPullLayer2Modules[i] = ibooker.book1D(histo, "RecHit YPull Layer2 by module", 100, -10.0, 10.0);
251 
252  sprintf(histo, "RecHit_YPull_Layer3_Module%d", i + 1);
253  recHitYPullLayer3Modules[i] = ibooker.book1D(histo, "RecHit YPull Layer3 by module", 100, -10.0, 10.0);
254  }
255 
256  ibooker.setCurrentFolder("TrackerRecHitsV/TrackerRecHits/Pixel/recHitPullsFPIX");
257  recHitXPullAllF = ibooker.book1D("RecHit_XPull_f_All", "RecHit X Pull All in Forward", 100, -10.0, 10.0);
258 
259  recHitYPullAllF = ibooker.book1D("RecHit_YPull_f_All", "RecHit Y Pull All in Forward", 100, -10.0, 10.0);
260 
261  for (int i = 0; i < 7; i++) {
262  sprintf(histo, "RecHit_XPull_Disk1_Plaquette%d", i + 1);
263  recHitXPullDisk1Plaquettes[i] = ibooker.book1D(histo, "RecHit XPull Disk1 by plaquette", 100, -10.0, 10.0);
264  sprintf(histo, "RecHit_XPull_Disk2_Plaquette%d", i + 1);
265  recHitXPullDisk2Plaquettes[i] = ibooker.book1D(histo, "RecHit XPull Disk2 by plaquette", 100, -10.0, 10.0);
266 
267  sprintf(histo, "RecHit_YPull_Disk1_Plaquette%d", i + 1);
268  recHitYPullDisk1Plaquettes[i] = ibooker.book1D(histo, "RecHit YPull Disk1 by plaquette", 100, -10.0, 10.0);
269 
270  sprintf(histo, "RecHit_YPull_Disk2_Plaquette%d", i + 1);
271  recHitYPullDisk2Plaquettes[i] = ibooker.book1D(histo, "RecHit YPull Disk2 by plaquette", 100, -10.0, 10.0);
272  }
273 }
MonitorElement * clustYSizeModule[8]
MonitorElement * clustXSizeDisk1Plaquettes[7]
MonitorElement * recHitNsimHitDisk1
MonitorElement * recHitXPullDisk2Plaquettes[7]
MonitorElement * recHitXResFlippedLadderLayers[3]
MonitorElement * recHitXFullModules
MonitorElement * recHitXResAllF
MonitorElement * recHitXHalfModules
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:32
MonitorElement * recHitYPullAllF
MonitorElement * recHitXResDisk1Plaquettes[7]
MonitorElement * recHitEventF
MonitorElement * clustXSizeDisk2Plaquettes[7]
MonitorElement * recHitYResAllB
MonitorElement * recHitYPullDisk2Plaquettes[7]
MonitorElement * recHitNsimHitLayer[3]
MonitorElement * clustXSizeLayer[3]
MonitorElement * clustChargeLayer1Modules[8]
MonitorElement * recHitYPlaquetteSize5
MonitorElement * recHitEventB
MonitorElement * recHitXPlaquetteSize1
MonitorElement * recHitXResDisk2Plaquettes[7]
MonitorElement * recHitXPullNonFlippedLadderLayers[3]
MonitorElement * recHitYPullDisk1Plaquettes[7]
MonitorElement * recHitXPullAllB
MonitorElement * recHitXResAllB
MonitorElement * recHitBunchB
MonitorElement * clustChargeLayer2Modules[8]
MonitorElement * recHitYResDisk2Plaquettes[7]
MonitorElement * recHitXResNonFlippedLadderLayers[3]
MonitorElement * recHitYResDisk1Plaquettes[7]
MonitorElement * recHitYResLayer2Modules[8]
MonitorElement * recHitYPlaquetteSize2
MonitorElement * recHitYPullLayer3Modules[8]
MonitorElement * clustYSizeDisk1Plaquettes[7]
MonitorElement * recHitBunchF
MonitorElement * recHitXPlaquetteSize2
MonitorElement * clustChargeLayer3Modules[8]
MonitorElement * recHitYResAllF
MonitorElement * recHitYResLayer3Modules[8]
MonitorElement * recHitYPullLayer2Modules[8]
MonitorElement * clustYSizeDisk2Plaquettes[7]
MonitorElement * recHitYPlaquetteSize4
MonitorElement * recHitXPullFlippedLadderLayers[3]
MonitorElement * recHitNsimHitDisk2
MonitorElement * recHitYResLayer1Modules[8]
MonitorElement * recHitYPullLayer1Modules[8]
MonitorElement * recHitXPullDisk1Plaquettes[7]
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
MonitorElement * recHitYPullAllB
MonitorElement * recHitXPullAllF
MonitorElement * recHitYPlaquetteSize3
MonitorElement * clustChargeDisk1Plaquettes[7]
MonitorElement * clustChargeDisk2Plaquettes[7]
MonitorElement * recHitYAllModules

◆ fillBarrel()

void SiPixelRecHitsValid::fillBarrel ( const SiPixelRecHit recHit,
const PSimHit simHit,
DetId  detId,
const PixelGeomDetUnit theGeomDet,
const TrackerTopology tTopo 
)
private

Definition at line 370 of file SiPixelRecHitsValid.cc.

References ALCARECOTkAlJpsiMuMu_cff::charge, clustChargeLayer1Modules, clustChargeLayer2Modules, clustChargeLayer3Modules, clustXSizeLayer, clustYSizeModule, cmtomicron, dqm::impl::MonitorElement::Fill(), HcalObjRepresent::Fill(), mps_fire::i, PixelTopology::nrows(), PV3DBase< T, PVType, FrameType >::perp(), TrackerTopology::pxbLayer(), TrackerTopology::pxbModule(), rpcPointValidation_cfi::recHit, recHitBunchB, recHitEventB, recHitXFullModules, recHitXHalfModules, recHitXPullAllB, recHitXPullFlippedLadderLayers, recHitXPullNonFlippedLadderLayers, recHitXResAllB, recHitXResFlippedLadderLayers, recHitXResNonFlippedLadderLayers, recHitYAllModules, recHitYPullAllB, recHitYPullLayer1Modules, recHitYPullLayer2Modules, recHitYPullLayer3Modules, recHitYResAllB, recHitYResLayer1Modules, recHitYResLayer2Modules, recHitYResLayer3Modules, postprocess-scan-build::rows, rpcPointValidation_cfi::simHit, RecoTauValidation_cfi::sizeX, RecoTauValidation_cfi::sizeY, PixelGeomDetUnit::specificTopology(), mathSSE::sqrt(), GeomDet::surface(), Surface::toGlobal(), PV3DBase< T, PVType, FrameType >::x(), LocalError::xx(), PV3DBase< T, PVType, FrameType >::y(), and LocalError::yy().

Referenced by analyze().

374  {
375  const float cmtomicron = 10000.0;
376 
377  int bunch = simHit.eventId().bunchCrossing();
378  int event = simHit.eventId().event();
379 
380  recHitBunchB->Fill(bunch);
381  if (bunch == 0)
383 
384  LocalPoint lp = recHit.localPosition();
385  float lp_y = lp.y();
386  float lp_x = lp.x();
387 
388  LocalError lerr = recHit.localPositionError();
389  float lerr_x = sqrt(lerr.xx());
390  float lerr_y = sqrt(lerr.yy());
391 
392  recHitYAllModules->Fill(lp_y);
393 
394  float sim_x1 = simHit.entryPoint().x();
395  float sim_x2 = simHit.exitPoint().x();
396  float sim_xpos = 0.5 * (sim_x1 + sim_x2);
397  float res_x = (lp.x() - sim_xpos) * cmtomicron;
398 
399  recHitXResAllB->Fill(res_x);
400 
401  float sim_y1 = simHit.entryPoint().y();
402  float sim_y2 = simHit.exitPoint().y();
403  float sim_ypos = 0.5 * (sim_y1 + sim_y2);
404  float res_y = (lp.y() - sim_ypos) * cmtomicron;
405 
406  recHitYResAllB->Fill(res_y);
407 
408  float pull_x = (lp_x - sim_xpos) / lerr_x;
409  float pull_y = (lp_y - sim_ypos) / lerr_y;
410 
411  recHitXPullAllB->Fill(pull_x);
412  recHitYPullAllB->Fill(pull_y);
413 
414  int rows = theGeomDet->specificTopology().nrows();
415 
416  if (rows == 160) {
417  recHitXFullModules->Fill(lp_x);
418  } else if (rows == 80) {
419  recHitXHalfModules->Fill(lp_x);
420  }
421 
422  float tmp1 = theGeomDet->surface().toGlobal(Local3DPoint(0., 0., 0.)).perp();
423  float tmp2 = theGeomDet->surface().toGlobal(Local3DPoint(0., 0., 1.)).perp();
424 
425  if (tmp2 < tmp1) { // flipped
426  for (unsigned int i = 0; i < 3; i++) {
427  if (tTopo->pxbLayer(detId) == i + 1) {
430  }
431  }
432  } else {
433  for (unsigned int i = 0; i < 3; i++) {
434  if (tTopo->pxbLayer(detId) == i + 1) {
437  }
438  }
439  }
440 
441  //get cluster
442  SiPixelRecHit::ClusterRef const& clust = recHit.cluster();
443 
444  // fill module dependent info
445  for (unsigned int i = 0; i < 8; i++) {
446  if (tTopo->pxbModule(detId) == i + 1) {
447  int sizeY = (*clust).sizeY();
449 
450  if (tTopo->pxbLayer(detId) == 1) {
451  float charge = (*clust).charge();
453  recHitYResLayer1Modules[i]->Fill(res_y);
454  recHitYPullLayer1Modules[i]->Fill(pull_y);
455  } else if (tTopo->pxbLayer(detId) == 2) {
456  float charge = (*clust).charge();
458  recHitYResLayer2Modules[i]->Fill(res_y);
459  recHitYPullLayer2Modules[i]->Fill(pull_y);
460  } else if (tTopo->pxbLayer(detId) == 3) {
461  float charge = (*clust).charge();
463  recHitYResLayer3Modules[i]->Fill(res_y);
464  recHitYPullLayer3Modules[i]->Fill(pull_y);
465  }
466  }
467  }
468  int sizeX = (*clust).sizeX();
469  if (tTopo->pxbLayer(detId) == 1)
471  if (tTopo->pxbLayer(detId) == 2)
473  if (tTopo->pxbLayer(detId) == 3)
475 }
MonitorElement * clustYSizeModule[8]
unsigned int pxbLayer(const DetId &id) const
MonitorElement * recHitXResFlippedLadderLayers[3]
T perp() const
Definition: PV3DBase.h:69
MonitorElement * recHitXFullModules
MonitorElement * recHitXHalfModules
const float cmtomicron
MonitorElement * recHitYResAllB
MonitorElement * clustXSizeLayer[3]
MonitorElement * clustChargeLayer1Modules[8]
virtual int nrows() const =0
MonitorElement * recHitEventB
MonitorElement * recHitXPullNonFlippedLadderLayers[3]
MonitorElement * recHitXPullAllB
MonitorElement * recHitXResAllB
void Fill(long long x)
MonitorElement * recHitBunchB
float yy() const
Definition: LocalError.h:24
T x() const
Definition: PV3DBase.h:59
T y() const
Definition: PV3DBase.h:60
MonitorElement * clustChargeLayer2Modules[8]
MonitorElement * recHitXResNonFlippedLadderLayers[3]
MonitorElement * recHitYResLayer2Modules[8]
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
T sqrt(T t)
Definition: SSEVec.h:19
MonitorElement * recHitYPullLayer3Modules[8]
GlobalPoint toGlobal(const Point2DBase< Scalar, LocalTag > lp) const
Definition: Surface.h:79
MonitorElement * clustChargeLayer3Modules[8]
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:37
MonitorElement * recHitYResLayer3Modules[8]
MonitorElement * recHitYPullLayer2Modules[8]
MonitorElement * recHitXPullFlippedLadderLayers[3]
virtual const PixelTopology & specificTopology() const
Returns a reference to the pixel proxy topology.
MonitorElement * recHitYResLayer1Modules[8]
MonitorElement * recHitYPullLayer1Modules[8]
unsigned int pxbModule(const DetId &id) const
MonitorElement * recHitYPullAllB
float xx() const
Definition: LocalError.h:22
Definition: event.py:1
MonitorElement * recHitYAllModules
Point3DBase< float, LocalTag > Local3DPoint
Definition: LocalPoint.h:9

◆ fillForward()

void SiPixelRecHitsValid::fillForward ( const SiPixelRecHit recHit,
const PSimHit simHit,
DetId  detId,
const PixelGeomDetUnit theGeomDet,
const TrackerTopology tTopo 
)
private

Definition at line 477 of file SiPixelRecHitsValid.cc.

References ALCARECOTkAlJpsiMuMu_cff::charge, clustChargeDisk1Plaquettes, clustChargeDisk2Plaquettes, clustXSizeDisk1Plaquettes, clustXSizeDisk2Plaquettes, clustYSizeDisk1Plaquettes, clustYSizeDisk2Plaquettes, cmtomicron, dqm::impl::MonitorElement::Fill(), mps_fire::i, PixelTopology::ncolumns(), PixelTopology::nrows(), TrackerTopology::pxfDisk(), TrackerTopology::pxfModule(), TrackerTopology::pxfPanel(), rpcPointValidation_cfi::recHit, recHitBunchF, recHitEventF, recHitXPlaquetteSize1, recHitXPlaquetteSize2, recHitXPullAllF, recHitXPullDisk1Plaquettes, recHitXPullDisk2Plaquettes, recHitXResAllF, recHitXResDisk1Plaquettes, recHitXResDisk2Plaquettes, recHitYPlaquetteSize2, recHitYPlaquetteSize3, recHitYPlaquetteSize4, recHitYPlaquetteSize5, recHitYPullAllF, recHitYPullDisk1Plaquettes, recHitYPullDisk2Plaquettes, recHitYResDisk1Plaquettes, recHitYResDisk2Plaquettes, postprocess-scan-build::rows, rpcPointValidation_cfi::simHit, RecoTauValidation_cfi::sizeX, RecoTauValidation_cfi::sizeY, PixelGeomDetUnit::specificTopology(), mathSSE::sqrt(), PV3DBase< T, PVType, FrameType >::x(), LocalError::xx(), PV3DBase< T, PVType, FrameType >::y(), and LocalError::yy().

Referenced by analyze().

481  {
482  int rows = theGeomDet->specificTopology().nrows();
483  int cols = theGeomDet->specificTopology().ncolumns();
484 
485  const float cmtomicron = 10000.0;
486 
487  int bunch = simHit.eventId().bunchCrossing();
488  int event = simHit.eventId().event();
489 
490  recHitBunchF->Fill(bunch);
491  if (bunch == 0)
493 
494  LocalPoint lp = recHit.localPosition();
495  float lp_x = lp.x();
496  float lp_y = lp.y();
497 
498  LocalError lerr = recHit.localPositionError();
499  float lerr_x = sqrt(lerr.xx());
500  float lerr_y = sqrt(lerr.yy());
501 
502  float sim_x1 = simHit.entryPoint().x();
503  float sim_x2 = simHit.exitPoint().x();
504  float sim_xpos = 0.5 * (sim_x1 + sim_x2);
505 
506  float sim_y1 = simHit.entryPoint().y();
507  float sim_y2 = simHit.exitPoint().y();
508  float sim_ypos = 0.5 * (sim_y1 + sim_y2);
509 
510  float pull_x = (lp_x - sim_xpos) / lerr_x;
511  float pull_y = (lp_y - sim_ypos) / lerr_y;
512 
513  if (rows == 80) {
515  } else if (rows == 160) {
517  }
518 
519  if (cols == 104) {
521  } else if (cols == 156) {
523  } else if (cols == 208) {
525  } else if (cols == 260) {
527  }
528 
529  float res_x = (lp.x() - sim_xpos) * cmtomicron;
530 
531  recHitXResAllF->Fill(res_x);
532  recHitXPullAllF->Fill(pull_x);
533 
534  float res_y = (lp.y() - sim_ypos) * cmtomicron;
535 
536  recHitYPullAllF->Fill(pull_y);
537 
538  // get cluster
539  SiPixelRecHit::ClusterRef const& clust = recHit.cluster();
540 
541  // fill plaquette dependent info
542  for (unsigned int i = 0; i < 7; i++) {
543  if (tTopo->pxfModule(detId) == i + 1) {
544  if (tTopo->pxfDisk(detId) == 1) {
545  int sizeX = (*clust).sizeX();
547 
548  int sizeY = (*clust).sizeY();
550 
551  float charge = (*clust).charge();
553 
556 
559  } else {
560  int sizeX = (*clust).sizeX();
562 
563  int sizeY = (*clust).sizeY();
565 
566  float charge = (*clust).charge();
568 
571 
574 
575  } // end else
576  } // end if module
577  else if (tTopo->pxfPanel(detId) == 2 && (tTopo->pxfModule(detId) + 4) == i + 1) {
578  if (tTopo->pxfDisk(detId) == 1) {
579  int sizeX = (*clust).sizeX();
581 
582  int sizeY = (*clust).sizeY();
584 
585  float charge = (*clust).charge();
587 
590 
593  } else {
594  int sizeX = (*clust).sizeX();
596 
597  int sizeY = (*clust).sizeY();
599 
600  float charge = (*clust).charge();
602 
605 
608 
609  } // end else
610  } // end else
611  } // end for
612 }
MonitorElement * clustXSizeDisk1Plaquettes[7]
MonitorElement * recHitXPullDisk2Plaquettes[7]
MonitorElement * recHitXResAllF
virtual int ncolumns() const =0
MonitorElement * recHitYPullAllF
const float cmtomicron
MonitorElement * recHitXResDisk1Plaquettes[7]
MonitorElement * recHitEventF
MonitorElement * clustXSizeDisk2Plaquettes[7]
unsigned int pxfModule(const DetId &id) const
MonitorElement * recHitYPullDisk2Plaquettes[7]
virtual int nrows() const =0
MonitorElement * recHitYPlaquetteSize5
MonitorElement * recHitXPlaquetteSize1
MonitorElement * recHitXResDisk2Plaquettes[7]
MonitorElement * recHitYPullDisk1Plaquettes[7]
void Fill(long long x)
float yy() const
Definition: LocalError.h:24
T x() const
Definition: PV3DBase.h:59
T y() const
Definition: PV3DBase.h:60
MonitorElement * recHitYResDisk2Plaquettes[7]
MonitorElement * recHitYResDisk1Plaquettes[7]
MonitorElement * recHitYPlaquetteSize2
T sqrt(T t)
Definition: SSEVec.h:19
MonitorElement * clustYSizeDisk1Plaquettes[7]
unsigned int pxfDisk(const DetId &id) const
MonitorElement * recHitBunchF
MonitorElement * recHitXPlaquetteSize2
unsigned int pxfPanel(const DetId &id) const
MonitorElement * clustYSizeDisk2Plaquettes[7]
MonitorElement * recHitYPlaquetteSize4
virtual const PixelTopology & specificTopology() const
Returns a reference to the pixel proxy topology.
MonitorElement * recHitXPullDisk1Plaquettes[7]
MonitorElement * recHitXPullAllF
MonitorElement * recHitYPlaquetteSize3
MonitorElement * clustChargeDisk1Plaquettes[7]
float xx() const
Definition: LocalError.h:22
MonitorElement * clustChargeDisk2Plaquettes[7]
Definition: event.py:1

Member Data Documentation

◆ clustChargeDisk1Plaquettes

MonitorElement* SiPixelRecHitsValid::clustChargeDisk1Plaquettes[7]
private

Definition at line 61 of file SiPixelRecHitsValid.h.

Referenced by bookHistograms(), and fillForward().

◆ clustChargeDisk2Plaquettes

MonitorElement* SiPixelRecHitsValid::clustChargeDisk2Plaquettes[7]
private

Definition at line 62 of file SiPixelRecHitsValid.h.

Referenced by bookHistograms(), and fillForward().

◆ clustChargeLayer1Modules

MonitorElement* SiPixelRecHitsValid::clustChargeLayer1Modules[8]
private

Definition at line 52 of file SiPixelRecHitsValid.h.

Referenced by bookHistograms(), and fillBarrel().

◆ clustChargeLayer2Modules

MonitorElement* SiPixelRecHitsValid::clustChargeLayer2Modules[8]
private

Definition at line 53 of file SiPixelRecHitsValid.h.

Referenced by bookHistograms(), and fillBarrel().

◆ clustChargeLayer3Modules

MonitorElement* SiPixelRecHitsValid::clustChargeLayer3Modules[8]
private

Definition at line 54 of file SiPixelRecHitsValid.h.

Referenced by bookHistograms(), and fillBarrel().

◆ clustXSizeDisk1Plaquettes

MonitorElement* SiPixelRecHitsValid::clustXSizeDisk1Plaquettes[7]
private

Definition at line 57 of file SiPixelRecHitsValid.h.

Referenced by bookHistograms(), and fillForward().

◆ clustXSizeDisk2Plaquettes

MonitorElement* SiPixelRecHitsValid::clustXSizeDisk2Plaquettes[7]
private

Definition at line 58 of file SiPixelRecHitsValid.h.

Referenced by bookHistograms(), and fillForward().

◆ clustXSizeLayer

MonitorElement* SiPixelRecHitsValid::clustXSizeLayer[3]
private

Definition at line 51 of file SiPixelRecHitsValid.h.

Referenced by bookHistograms(), and fillBarrel().

◆ clustYSizeDisk1Plaquettes

MonitorElement* SiPixelRecHitsValid::clustYSizeDisk1Plaquettes[7]
private

Definition at line 59 of file SiPixelRecHitsValid.h.

Referenced by bookHistograms(), and fillForward().

◆ clustYSizeDisk2Plaquettes

MonitorElement* SiPixelRecHitsValid::clustYSizeDisk2Plaquettes[7]
private

Definition at line 60 of file SiPixelRecHitsValid.h.

Referenced by bookHistograms(), and fillForward().

◆ clustYSizeModule

MonitorElement* SiPixelRecHitsValid::clustYSizeModule[8]
private

Definition at line 50 of file SiPixelRecHitsValid.h.

Referenced by bookHistograms(), and fillBarrel().

◆ recHitBunchB

MonitorElement* SiPixelRecHitsValid::recHitBunchB
private

Definition at line 75 of file SiPixelRecHitsValid.h.

Referenced by bookHistograms(), and fillBarrel().

◆ recHitBunchF

MonitorElement* SiPixelRecHitsValid::recHitBunchF
private

Definition at line 92 of file SiPixelRecHitsValid.h.

Referenced by bookHistograms(), and fillForward().

◆ recHitEventB

MonitorElement* SiPixelRecHitsValid::recHitEventB
private

Definition at line 76 of file SiPixelRecHitsValid.h.

Referenced by bookHistograms(), and fillBarrel().

◆ recHitEventF

MonitorElement* SiPixelRecHitsValid::recHitEventF
private

Definition at line 93 of file SiPixelRecHitsValid.h.

Referenced by bookHistograms(), and fillForward().

◆ recHitNsimHitDisk1

MonitorElement* SiPixelRecHitsValid::recHitNsimHitDisk1
private

Definition at line 94 of file SiPixelRecHitsValid.h.

Referenced by analyze(), and bookHistograms().

◆ recHitNsimHitDisk2

MonitorElement* SiPixelRecHitsValid::recHitNsimHitDisk2
private

Definition at line 95 of file SiPixelRecHitsValid.h.

Referenced by analyze(), and bookHistograms().

◆ recHitNsimHitLayer

MonitorElement* SiPixelRecHitsValid::recHitNsimHitLayer[3]
private

Definition at line 77 of file SiPixelRecHitsValid.h.

Referenced by analyze(), and bookHistograms().

◆ recHitXFullModules

MonitorElement* SiPixelRecHitsValid::recHitXFullModules
private

Definition at line 67 of file SiPixelRecHitsValid.h.

Referenced by bookHistograms(), and fillBarrel().

◆ recHitXHalfModules

MonitorElement* SiPixelRecHitsValid::recHitXHalfModules
private

Definition at line 68 of file SiPixelRecHitsValid.h.

Referenced by bookHistograms(), and fillBarrel().

◆ recHitXPlaquetteSize1

MonitorElement* SiPixelRecHitsValid::recHitXPlaquetteSize1
private

Definition at line 82 of file SiPixelRecHitsValid.h.

Referenced by bookHistograms(), and fillForward().

◆ recHitXPlaquetteSize2

MonitorElement* SiPixelRecHitsValid::recHitXPlaquetteSize2
private

Definition at line 83 of file SiPixelRecHitsValid.h.

Referenced by bookHistograms(), and fillForward().

◆ recHitXPullAllB

MonitorElement* SiPixelRecHitsValid::recHitXPullAllB
private

Definition at line 99 of file SiPixelRecHitsValid.h.

Referenced by bookHistograms(), and fillBarrel().

◆ recHitXPullAllF

MonitorElement* SiPixelRecHitsValid::recHitXPullAllF
private

Definition at line 109 of file SiPixelRecHitsValid.h.

Referenced by bookHistograms(), and fillForward().

◆ recHitXPullDisk1Plaquettes

MonitorElement* SiPixelRecHitsValid::recHitXPullDisk1Plaquettes[7]
private

Definition at line 112 of file SiPixelRecHitsValid.h.

Referenced by bookHistograms(), and fillForward().

◆ recHitXPullDisk2Plaquettes

MonitorElement* SiPixelRecHitsValid::recHitXPullDisk2Plaquettes[7]
private

Definition at line 113 of file SiPixelRecHitsValid.h.

Referenced by bookHistograms(), and fillForward().

◆ recHitXPullFlippedLadderLayers

MonitorElement* SiPixelRecHitsValid::recHitXPullFlippedLadderLayers[3]
private

Definition at line 102 of file SiPixelRecHitsValid.h.

Referenced by bookHistograms(), and fillBarrel().

◆ recHitXPullNonFlippedLadderLayers

MonitorElement* SiPixelRecHitsValid::recHitXPullNonFlippedLadderLayers[3]
private

Definition at line 103 of file SiPixelRecHitsValid.h.

Referenced by bookHistograms(), and fillBarrel().

◆ recHitXResAllB

MonitorElement* SiPixelRecHitsValid::recHitXResAllB
private

Definition at line 65 of file SiPixelRecHitsValid.h.

Referenced by bookHistograms(), and fillBarrel().

◆ recHitXResAllF

MonitorElement* SiPixelRecHitsValid::recHitXResAllF
private

Definition at line 80 of file SiPixelRecHitsValid.h.

Referenced by bookHistograms(), and fillForward().

◆ recHitXResDisk1Plaquettes

MonitorElement* SiPixelRecHitsValid::recHitXResDisk1Plaquettes[7]
private

Definition at line 88 of file SiPixelRecHitsValid.h.

Referenced by bookHistograms(), and fillForward().

◆ recHitXResDisk2Plaquettes

MonitorElement* SiPixelRecHitsValid::recHitXResDisk2Plaquettes[7]
private

Definition at line 89 of file SiPixelRecHitsValid.h.

Referenced by bookHistograms(), and fillForward().

◆ recHitXResFlippedLadderLayers

MonitorElement* SiPixelRecHitsValid::recHitXResFlippedLadderLayers[3]
private

Definition at line 70 of file SiPixelRecHitsValid.h.

Referenced by bookHistograms(), and fillBarrel().

◆ recHitXResNonFlippedLadderLayers

MonitorElement* SiPixelRecHitsValid::recHitXResNonFlippedLadderLayers[3]
private

Definition at line 71 of file SiPixelRecHitsValid.h.

Referenced by bookHistograms(), and fillBarrel().

◆ recHitYAllModules

MonitorElement* SiPixelRecHitsValid::recHitYAllModules
private

Definition at line 69 of file SiPixelRecHitsValid.h.

Referenced by bookHistograms(), and fillBarrel().

◆ recHitYPlaquetteSize2

MonitorElement* SiPixelRecHitsValid::recHitYPlaquetteSize2
private

Definition at line 84 of file SiPixelRecHitsValid.h.

Referenced by bookHistograms(), and fillForward().

◆ recHitYPlaquetteSize3

MonitorElement* SiPixelRecHitsValid::recHitYPlaquetteSize3
private

Definition at line 85 of file SiPixelRecHitsValid.h.

Referenced by bookHistograms(), and fillForward().

◆ recHitYPlaquetteSize4

MonitorElement* SiPixelRecHitsValid::recHitYPlaquetteSize4
private

Definition at line 86 of file SiPixelRecHitsValid.h.

Referenced by bookHistograms(), and fillForward().

◆ recHitYPlaquetteSize5

MonitorElement* SiPixelRecHitsValid::recHitYPlaquetteSize5
private

Definition at line 87 of file SiPixelRecHitsValid.h.

Referenced by bookHistograms(), and fillForward().

◆ recHitYPullAllB

MonitorElement* SiPixelRecHitsValid::recHitYPullAllB
private

Definition at line 100 of file SiPixelRecHitsValid.h.

Referenced by bookHistograms(), and fillBarrel().

◆ recHitYPullAllF

MonitorElement* SiPixelRecHitsValid::recHitYPullAllF
private

Definition at line 110 of file SiPixelRecHitsValid.h.

Referenced by bookHistograms(), and fillForward().

◆ recHitYPullDisk1Plaquettes

MonitorElement* SiPixelRecHitsValid::recHitYPullDisk1Plaquettes[7]
private

Definition at line 114 of file SiPixelRecHitsValid.h.

Referenced by bookHistograms(), and fillForward().

◆ recHitYPullDisk2Plaquettes

MonitorElement* SiPixelRecHitsValid::recHitYPullDisk2Plaquettes[7]
private

Definition at line 115 of file SiPixelRecHitsValid.h.

Referenced by bookHistograms(), and fillForward().

◆ recHitYPullLayer1Modules

MonitorElement* SiPixelRecHitsValid::recHitYPullLayer1Modules[8]
private

Definition at line 104 of file SiPixelRecHitsValid.h.

Referenced by bookHistograms(), and fillBarrel().

◆ recHitYPullLayer2Modules

MonitorElement* SiPixelRecHitsValid::recHitYPullLayer2Modules[8]
private

Definition at line 105 of file SiPixelRecHitsValid.h.

Referenced by bookHistograms(), and fillBarrel().

◆ recHitYPullLayer3Modules

MonitorElement* SiPixelRecHitsValid::recHitYPullLayer3Modules[8]
private

Definition at line 106 of file SiPixelRecHitsValid.h.

Referenced by bookHistograms(), and fillBarrel().

◆ recHitYResAllB

MonitorElement* SiPixelRecHitsValid::recHitYResAllB
private

Definition at line 66 of file SiPixelRecHitsValid.h.

Referenced by bookHistograms(), and fillBarrel().

◆ recHitYResAllF

MonitorElement* SiPixelRecHitsValid::recHitYResAllF
private

Definition at line 81 of file SiPixelRecHitsValid.h.

Referenced by bookHistograms().

◆ recHitYResDisk1Plaquettes

MonitorElement* SiPixelRecHitsValid::recHitYResDisk1Plaquettes[7]
private

Definition at line 90 of file SiPixelRecHitsValid.h.

Referenced by bookHistograms(), and fillForward().

◆ recHitYResDisk2Plaquettes

MonitorElement* SiPixelRecHitsValid::recHitYResDisk2Plaquettes[7]
private

Definition at line 91 of file SiPixelRecHitsValid.h.

Referenced by bookHistograms(), and fillForward().

◆ recHitYResLayer1Modules

MonitorElement* SiPixelRecHitsValid::recHitYResLayer1Modules[8]
private

Definition at line 72 of file SiPixelRecHitsValid.h.

Referenced by bookHistograms(), and fillBarrel().

◆ recHitYResLayer2Modules

MonitorElement* SiPixelRecHitsValid::recHitYResLayer2Modules[8]
private

Definition at line 73 of file SiPixelRecHitsValid.h.

Referenced by bookHistograms(), and fillBarrel().

◆ recHitYResLayer3Modules

MonitorElement* SiPixelRecHitsValid::recHitYResLayer3Modules[8]
private

Definition at line 74 of file SiPixelRecHitsValid.h.

Referenced by bookHistograms(), and fillBarrel().

◆ siPixelRecHitCollectionToken_

edm::EDGetTokenT<SiPixelRecHitCollection> SiPixelRecHitsValid::siPixelRecHitCollectionToken_
private

Definition at line 118 of file SiPixelRecHitsValid.h.

Referenced by analyze().

◆ tGeomEsToken

const edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> SiPixelRecHitsValid::tGeomEsToken
private

Definition at line 47 of file SiPixelRecHitsValid.h.

Referenced by analyze().

◆ trackerHitAssociatorConfig_

TrackerHitAssociator::Config SiPixelRecHitsValid::trackerHitAssociatorConfig_
private

Definition at line 117 of file SiPixelRecHitsValid.h.

Referenced by analyze().

◆ tTopoEsToken

const edm::ESGetToken<TrackerTopology, TrackerTopologyRcd> SiPixelRecHitsValid::tTopoEsToken
private

Definition at line 46 of file SiPixelRecHitsValid.h.

Referenced by analyze().