Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:33:46

0001 // SiPixelRecHitsValid.cc
0002 // Description: see SiPixelRecHitsValid.h
0003 // Author: Jason Shaev, JHU
0004 // Created 6/7/06
0005 //
0006 // G. Giurgiu, JHU (ggiurgiu@pha.jhu.edu)
0007 //             added pull distributions (12/27/06)
0008 //--------------------------------
0009 
0010 #include "Validation/TrackerRecHits/interface/SiPixelRecHitsValid.h"
0011 
0012 #include "DataFormats/Common/interface/DetSetVector.h"
0013 #include "DataFormats/Common/interface/Handle.h"
0014 #include "DataFormats/Common/interface/OwnVector.h"
0015 #include "DataFormats/Common/interface/Ref.h"
0016 #include "DataFormats/DetId/interface/DetId.h"
0017 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
0018 #include "DataFormats/GeometryVector/interface/LocalPoint.h"
0019 #include "DataFormats/SiPixelCluster/interface/SiPixelCluster.h"
0020 #include "DataFormats/SiPixelDetId/interface/PixelSubdetector.h"
0021 #include "DataFormats/TrackerCommon/interface/TrackerTopology.h"
0022 #include "DataFormats/TrackerRecHit2D/interface/SiPixelRecHit.h"
0023 
0024 #include "DQMServices/Core/interface/DQMStore.h"
0025 
0026 #include "FWCore/Framework/interface/Event.h"
0027 #include "FWCore/Framework/interface/EventSetup.h"
0028 #include "FWCore/Framework/interface/ESHandle.h"
0029 #include "FWCore/Framework/interface/MakerMacros.h"
0030 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0031 #include "FWCore/ServiceRegistry/interface/Service.h"
0032 
0033 #include "Geometry/CommonDetUnit/interface/GeomDetType.h"
0034 #include "Geometry/CommonDetUnit/interface/GeomDet.h"
0035 #include "Geometry/CommonTopologies/interface/PixelTopology.h"
0036 #include "Geometry/CommonDetUnit/interface/PixelGeomDetType.h"
0037 #include "Geometry/CommonDetUnit/interface/PixelGeomDetUnit.h"
0038 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
0039 #include "Geometry/TrackerNumberingBuilder/interface/GeometricDet.h"
0040 
0041 #include "SimDataFormats/TrackingHit/interface/PSimHit.h"
0042 #include "SimDataFormats/TrackingHit/interface/PSimHitContainer.h"
0043 
0044 #include <cmath>
0045 
0046 SiPixelRecHitsValid::SiPixelRecHitsValid(const edm::ParameterSet& ps)
0047     : tTopoEsToken(esConsumes()),
0048       tGeomEsToken(esConsumes()),
0049       trackerHitAssociatorConfig_(ps, consumesCollector()),
0050       siPixelRecHitCollectionToken_(consumes<SiPixelRecHitCollection>(ps.getParameter<edm::InputTag>("src"))) {}
0051 
0052 SiPixelRecHitsValid::~SiPixelRecHitsValid() = default;
0053 
0054 void SiPixelRecHitsValid::bookHistograms(DQMStore::IBooker& ibooker, const edm::Run& run, const edm::EventSetup& es) {
0055   ibooker.setCurrentFolder("TrackerRecHitsV/TrackerRecHits/Pixel/clustBPIX");
0056 
0057   Char_t histo[200];
0058 
0059   // ---------------------------------------------------------------
0060   // All histograms that depend on plaquette number have 7 indexes.
0061   // The first 4 (0-3) correspond to Panel 1 plaquettes 1-4.
0062   // The last 3 (4-6) correspond to Panel 2 plaquettes 1-3.
0063   // ---------------------------------------------------------------
0064 
0065   //Cluster y-size by module number for barrel
0066   for (int i = 0; i < 8; i++) {
0067     sprintf(histo, "Clust_y_size_Module%d", i + 1);
0068     clustYSizeModule[i] = ibooker.book1D(histo, "Cluster y-size by Module", 20, 0.5, 20.5);
0069   }  // end for
0070 
0071   //Cluster x-size by layer for barrel
0072   for (int i = 0; i < 3; i++) {
0073     sprintf(histo, "Clust_x_size_Layer%d", i + 1);
0074     clustXSizeLayer[i] = ibooker.book1D(histo, "Cluster x-size by Layer", 20, 0.5, 20.5);
0075   }  // end for
0076 
0077   //Cluster charge by module for 3 layers of barrel
0078   for (int i = 0; i < 8; i++) {
0079     //Cluster charge by module for Layer1
0080     sprintf(histo, "Clust_charge_Layer1_Module%d", i + 1);
0081     clustChargeLayer1Modules[i] = ibooker.book1D(histo, "Cluster charge Layer 1 by Module", 50, 0., 200000.);
0082 
0083     //Cluster charge by module for Layer2
0084     sprintf(histo, "Clust_charge_Layer2_Module%d", i + 1);
0085     clustChargeLayer2Modules[i] = ibooker.book1D(histo, "Cluster charge Layer 2 by Module", 50, 0., 200000.);
0086 
0087     //Cluster charge by module for Layer3
0088     sprintf(histo, "Clust_charge_Layer3_Module%d", i + 1);
0089     clustChargeLayer3Modules[i] = ibooker.book1D(histo, "Cluster charge Layer 3 by Module", 50, 0., 200000.);
0090   }  // end for
0091 
0092   ibooker.setCurrentFolder("TrackerRecHitsV/TrackerRecHits/Pixel/clustFPIX");
0093   //Cluster x-size, y-size, and charge by plaquette for Disks in Forward
0094   for (int i = 0; i < 7; i++) {
0095     //Cluster x-size for Disk1 by Plaquette
0096     sprintf(histo, "Clust_x_size_Disk1_Plaquette%d", i + 1);
0097     clustXSizeDisk1Plaquettes[i] = ibooker.book1D(histo, "Cluster X-size for Disk1 by Plaquette", 20, 0.5, 20.5);
0098 
0099     //Cluster x-size for Disk2 by Plaquette
0100     sprintf(histo, "Clust_x_size_Disk2_Plaquette%d", i + 1);
0101     clustXSizeDisk2Plaquettes[i] = ibooker.book1D(histo, "Cluster X-size for Disk2 by Plaquette", 20, 0.5, 20.5);
0102 
0103     //Cluster y-size for Disk1 by Plaquette
0104     sprintf(histo, "Clust_y_size_Disk1_Plaquette%d", i + 1);
0105     clustYSizeDisk1Plaquettes[i] = ibooker.book1D(histo, "Cluster Y-size for Disk1 by Plaquette", 20, 0.5, 20.5);
0106 
0107     //Cluster y-size for Disk2 by Plaquette
0108     sprintf(histo, "Clust_y_size_Disk2_Plaquette%d", i + 1);
0109     clustYSizeDisk2Plaquettes[i] = ibooker.book1D(histo, "Cluster Y-size for Disk2 by Plaquette", 20, 0.5, 20.5);
0110 
0111     //Cluster charge for Disk1 by Plaquette
0112     sprintf(histo, "Clust_charge_Disk1_Plaquette%d", i + 1);
0113     clustChargeDisk1Plaquettes[i] = ibooker.book1D(histo, "Cluster charge for Disk1 by Plaquette", 50, 0., 200000.);
0114 
0115     //Cluster charge for Disk2 by Plaquette
0116     sprintf(histo, "Clust_charge_Disk2_Plaquette%d", i + 1);
0117     clustChargeDisk2Plaquettes[i] = ibooker.book1D(histo, "Cluster charge for Disk2 by Plaquette", 50, 0., 200000.);
0118   }  // end for
0119 
0120   ibooker.setCurrentFolder("TrackerRecHitsV/TrackerRecHits/Pixel/recHitBPIX");
0121   //RecHit Bunch crossing all barrel hits
0122   recHitBunchB = ibooker.book1D("RecHit_Bunch_Barrel", "RecHit Bunch Crossing, Barrel", 20, -10., 10.);
0123 
0124   //RecHit Event, in-time bunch, all barrel hits
0125   recHitEventB = ibooker.book1D("RecHit_Event_Barrel", "RecHit Event (in-time bunch), Barrel", 100, 0., 100.);
0126 
0127   //RecHit X Resolution all barrel hits
0128   recHitXResAllB = ibooker.book1D("RecHit_xres_b_All", "RecHit X Res All Modules in Barrel", 100, -200., 200.);
0129 
0130   //RecHit Y Resolution all barrel hits
0131   recHitYResAllB = ibooker.book1D("RecHit_yres_b_All", "RecHit Y Res All Modules in Barrel", 100, -200., 200.);
0132 
0133   //RecHit X distribution for full modules for barrel
0134   recHitXFullModules = ibooker.book1D("RecHit_x_FullModules", "RecHit X distribution for full modules", 100, -2., 2.);
0135 
0136   //RecHit X distribution for half modules for barrel
0137   recHitXHalfModules = ibooker.book1D("RecHit_x_HalfModules", "RecHit X distribution for half modules", 100, -1., 1.);
0138 
0139   //RecHit Y distribution all modules for barrel
0140   recHitYAllModules = ibooker.book1D("RecHit_y_AllModules", "RecHit Y distribution for all modules", 100, -4., 4.);
0141 
0142   //RecHit X resolution for flipped and unflipped ladders by layer for barrel
0143   for (int i = 0; i < 3; i++) {
0144     //RecHit no. of matched simHits all ladders by layer
0145     sprintf(histo, "RecHit_NsimHit_Layer%d", i + 1);
0146     recHitNsimHitLayer[i] = ibooker.book1D(histo, "RecHit Number of simHits by Layer", 30, 0., 30.);
0147 
0148     //RecHit X resolution for flipped ladders by layer
0149     sprintf(histo, "RecHit_XRes_FlippedLadder_Layer%d", i + 1);
0150     recHitXResFlippedLadderLayers[i] = ibooker.book1D(histo, "RecHit XRes Flipped Ladders by Layer", 100, -200., 200.);
0151 
0152     //RecHit X resolution for unflipped ladders by layer
0153     sprintf(histo, "RecHit_XRes_UnFlippedLadder_Layer%d", i + 1);
0154     recHitXResNonFlippedLadderLayers[i] =
0155         ibooker.book1D(histo, "RecHit XRes NonFlipped Ladders by Layer", 100, -200., 200.);
0156   }  // end for
0157 
0158   //RecHit Y resolutions for layers by module for barrel
0159   for (int i = 0; i < 8; i++) {
0160     //Rec Hit Y resolution by module for Layer1
0161     sprintf(histo, "RecHit_YRes_Layer1_Module%d", i + 1);
0162     recHitYResLayer1Modules[i] = ibooker.book1D(histo, "RecHit YRes Layer1 by module", 100, -200., 200.);
0163 
0164     //RecHit Y resolution by module for Layer2
0165     sprintf(histo, "RecHit_YRes_Layer2_Module%d", i + 1);
0166     recHitYResLayer2Modules[i] = ibooker.book1D(histo, "RecHit YRes Layer2 by module", 100, -200., 200.);
0167 
0168     //RecHit Y resolution by module for Layer3
0169     sprintf(histo, "RecHit_YRes_Layer3_Module%d", i + 1);
0170     recHitYResLayer3Modules[i] = ibooker.book1D(histo, "RecHit YRes Layer3 by module", 100, -200., 200.);
0171   }  // end for
0172 
0173   ibooker.setCurrentFolder("TrackerRecHitsV/TrackerRecHits/Pixel/recHitFPIX");
0174   //RecHit Bunch crossing all plaquettes
0175   recHitBunchF = ibooker.book1D("RecHit_Bunch_Forward", "RecHit Bunch Crossing, Forward", 20, -10., 10.);
0176 
0177   //RecHit Event, in-time bunch, all plaquettes
0178   recHitEventF = ibooker.book1D("RecHit_Event_Forward", "RecHit Event (in-time bunch), Forward", 100, 0., 100.);
0179 
0180   //RecHit No. of simHits, by disk
0181   recHitNsimHitDisk1 = ibooker.book1D("RecHit_NsimHit_Disk1", "RecHit Number of simHits, Disk1", 30, 0., 30.);
0182   recHitNsimHitDisk2 = ibooker.book1D("RecHit_NsimHit_Disk2", "RecHit Number of simHits, Disk2", 30, 0., 30.);
0183 
0184   //RecHit X resolution all plaquettes
0185   recHitXResAllF = ibooker.book1D("RecHit_xres_f_All", "RecHit X Res All in Forward", 100, -200., 200.);
0186 
0187   //RecHit Y resolution all plaquettes
0188   recHitYResAllF = ibooker.book1D("RecHit_yres_f_All", "RecHit Y Res All in Forward", 100, -200., 200.);
0189 
0190   //RecHit X distribution for plaquette with x-size 1 in forward
0191   recHitXPlaquetteSize1 =
0192       ibooker.book1D("RecHit_x_Plaquette_xsize1", "RecHit X Distribution for plaquette x-size1", 100, -2., 2.);
0193 
0194   //RecHit X distribution for plaquette with x-size 2 in forward
0195   recHitXPlaquetteSize2 =
0196       ibooker.book1D("RecHit_x_Plaquette_xsize2", "RecHit X Distribution for plaquette x-size2", 100, -2., 2.);
0197 
0198   //RecHit Y distribution for plaquette with y-size 2 in forward
0199   recHitYPlaquetteSize2 =
0200       ibooker.book1D("RecHit_y_Plaquette_ysize2", "RecHit Y Distribution for plaquette y-size2", 100, -4., 4.);
0201 
0202   //RecHit Y distribution for plaquette with y-size 3 in forward
0203   recHitYPlaquetteSize3 =
0204       ibooker.book1D("RecHit_y_Plaquette_ysize3", "RecHit Y Distribution for plaquette y-size3", 100, -4., 4.);
0205 
0206   //RecHit Y distribution for plaquette with y-size 4 in forward
0207   recHitYPlaquetteSize4 =
0208       ibooker.book1D("RecHit_y_Plaquette_ysize4", "RecHit Y Distribution for plaquette y-size4", 100, -4., 4.);
0209 
0210   //RecHit Y distribution for plaquette with y-size 5 in forward
0211   recHitYPlaquetteSize5 =
0212       ibooker.book1D("RecHit_y_Plaquette_ysize5", "RecHit Y Distribution for plaquette y-size5", 100, -4., 4.);
0213 
0214   //X and Y resolutions for both disks by plaquette in forward
0215   for (int i = 0; i < 7; i++) {
0216     //X resolution for Disk1 by plaquette
0217     sprintf(histo, "RecHit_XRes_Disk1_Plaquette%d", i + 1);
0218     recHitXResDisk1Plaquettes[i] = ibooker.book1D(histo, "RecHit XRes Disk1 by plaquette", 100, -200., 200.);
0219     //X resolution for Disk2 by plaquette
0220     sprintf(histo, "RecHit_XRes_Disk2_Plaquette%d", i + 1);
0221     recHitXResDisk2Plaquettes[i] = ibooker.book1D(histo, "RecHit XRes Disk2 by plaquette", 100, -200., 200.);
0222 
0223     //Y resolution for Disk1 by plaquette
0224     sprintf(histo, "RecHit_YRes_Disk1_Plaquette%d", i + 1);
0225     recHitYResDisk1Plaquettes[i] = ibooker.book1D(histo, "RecHit YRes Disk1 by plaquette", 100, -200., 200.);
0226     //Y resolution for Disk2 by plaquette
0227     sprintf(histo, "RecHit_YRes_Disk2_Plaquette%d", i + 1);
0228     recHitYResDisk2Plaquettes[i] = ibooker.book1D(histo, "RecHit YRes Disk2 by plaquette", 100, -200., 200.);
0229   }
0230 
0231   ibooker.setCurrentFolder("TrackerRecHitsV/TrackerRecHits/Pixel/recHitPullsBPIX");
0232   recHitXPullAllB = ibooker.book1D("RecHit_xres_b_All", "RecHit X Pull All Modules in Barrel", 100, -10.0, 10.0);
0233   recHitYPullAllB = ibooker.book1D("RecHit_yres_b_All", "RecHit Y Pull All Modules in Barrel", 100, -10.0, 10.0);
0234 
0235   for (int i = 0; i < 3; i++) {
0236     sprintf(histo, "RecHit_XPull_FlippedLadder_Layer%d", i + 1);
0237     recHitXPullFlippedLadderLayers[i] =
0238         ibooker.book1D(histo, "RecHit XPull Flipped Ladders by Layer", 100, -10.0, 10.0);
0239 
0240     sprintf(histo, "RecHit_XPull_UnFlippedLadder_Layer%d", i + 1);
0241     recHitXPullNonFlippedLadderLayers[i] =
0242         ibooker.book1D(histo, "RecHit XPull NonFlipped Ladders by Layer", 100, -10.0, 10.0);
0243   }
0244 
0245   for (int i = 0; i < 8; i++) {
0246     sprintf(histo, "RecHit_YPull_Layer1_Module%d", i + 1);
0247     recHitYPullLayer1Modules[i] = ibooker.book1D(histo, "RecHit YPull Layer1 by module", 100, -10.0, 10.0);
0248 
0249     sprintf(histo, "RecHit_YPull_Layer2_Module%d", i + 1);
0250     recHitYPullLayer2Modules[i] = ibooker.book1D(histo, "RecHit YPull Layer2 by module", 100, -10.0, 10.0);
0251 
0252     sprintf(histo, "RecHit_YPull_Layer3_Module%d", i + 1);
0253     recHitYPullLayer3Modules[i] = ibooker.book1D(histo, "RecHit YPull Layer3 by module", 100, -10.0, 10.0);
0254   }
0255 
0256   ibooker.setCurrentFolder("TrackerRecHitsV/TrackerRecHits/Pixel/recHitPullsFPIX");
0257   recHitXPullAllF = ibooker.book1D("RecHit_XPull_f_All", "RecHit X Pull All in Forward", 100, -10.0, 10.0);
0258 
0259   recHitYPullAllF = ibooker.book1D("RecHit_YPull_f_All", "RecHit Y Pull All in Forward", 100, -10.0, 10.0);
0260 
0261   for (int i = 0; i < 7; i++) {
0262     sprintf(histo, "RecHit_XPull_Disk1_Plaquette%d", i + 1);
0263     recHitXPullDisk1Plaquettes[i] = ibooker.book1D(histo, "RecHit XPull Disk1 by plaquette", 100, -10.0, 10.0);
0264     sprintf(histo, "RecHit_XPull_Disk2_Plaquette%d", i + 1);
0265     recHitXPullDisk2Plaquettes[i] = ibooker.book1D(histo, "RecHit XPull Disk2 by plaquette", 100, -10.0, 10.0);
0266 
0267     sprintf(histo, "RecHit_YPull_Disk1_Plaquette%d", i + 1);
0268     recHitYPullDisk1Plaquettes[i] = ibooker.book1D(histo, "RecHit YPull Disk1 by plaquette", 100, -10.0, 10.0);
0269 
0270     sprintf(histo, "RecHit_YPull_Disk2_Plaquette%d", i + 1);
0271     recHitYPullDisk2Plaquettes[i] = ibooker.book1D(histo, "RecHit YPull Disk2 by plaquette", 100, -10.0, 10.0);
0272   }
0273 }
0274 
0275 void SiPixelRecHitsValid::analyze(const edm::Event& e, const edm::EventSetup& es) {
0276   //Retrieve tracker topology from geometry
0277   const TrackerTopology* tTopo = &es.getData(tTopoEsToken);
0278 
0279   edm::LogInfo("EventInfo") << " Run = " << e.id().run() << " Event = " << e.id().event();
0280   if (e.id().event() % 1000 == 0)
0281     std::cout << " Run = " << e.id().run() << " Event = " << e.id().event() << std::endl;
0282 
0283   //Get RecHits
0284   edm::Handle<SiPixelRecHitCollection> recHitColl;
0285   e.getByToken(siPixelRecHitCollectionToken_, recHitColl);
0286 
0287   //Get event setup
0288   const TrackerGeometry* theTracker = &es.getData(tGeomEsToken);
0289 
0290   TrackerHitAssociator associate(e, trackerHitAssociatorConfig_);
0291 
0292   //iterate over detunits
0293   //for (TrackerGeometry::DetContainer::const_iterator it = geom->dets().begin(); it != geom->dets().end(); it++) {
0294   for (const auto& it : theTracker->dets()) {
0295     DetId detId = it->geographicalId();
0296     unsigned int subid = detId.subdetId();
0297 
0298     if (!((subid == 1) || (subid == 2)))
0299       continue;
0300 
0301     const PixelGeomDetUnit* theGeomDet = dynamic_cast<const PixelGeomDetUnit*>(theTracker->idToDet(detId));
0302 
0303     SiPixelRecHitCollection::const_iterator pixeldet = recHitColl->find(detId);
0304     if (pixeldet == recHitColl->end())
0305       continue;
0306     SiPixelRecHitCollection::DetSet pixelrechitRange = *pixeldet;
0307     SiPixelRecHitCollection::DetSet::const_iterator pixelrechitRangeIteratorBegin = pixelrechitRange.begin();
0308     SiPixelRecHitCollection::DetSet::const_iterator pixelrechitRangeIteratorEnd = pixelrechitRange.end();
0309     SiPixelRecHitCollection::DetSet::const_iterator pixeliter = pixelrechitRangeIteratorBegin;
0310     std::vector<PSimHit> matched;
0311 
0312     //----Loop over rechits for this detId
0313     for (; pixeliter != pixelrechitRangeIteratorEnd; pixeliter++) {
0314       matched.clear();
0315       matched = associate.associateHit(*pixeliter);
0316 
0317       if (!matched.empty()) {
0318         float closest = 9999.9;
0319         std::vector<PSimHit>::const_iterator closestit = matched.begin();
0320         LocalPoint lp = pixeliter->localPosition();
0321         float rechit_x = lp.x();
0322         float rechit_y = lp.y();
0323 
0324         //loop over sim hits and fill closet
0325         for (std::vector<PSimHit>::const_iterator m = matched.begin(); m < matched.end(); m++) {
0326           float sim_x1 = (*m).entryPoint().x();
0327           float sim_x2 = (*m).exitPoint().x();
0328           float sim_xpos = 0.5 * (sim_x1 + sim_x2);
0329 
0330           float sim_y1 = (*m).entryPoint().y();
0331           float sim_y2 = (*m).exitPoint().y();
0332           float sim_ypos = 0.5 * (sim_y1 + sim_y2);
0333 
0334           float x_res = fabs(sim_xpos - rechit_x);
0335           float y_res = fabs(sim_ypos - rechit_y);
0336 
0337           float dist = sqrt(x_res * x_res + y_res * y_res);
0338 
0339           if (dist < closest) {
0340             closest = dist;
0341             closestit = m;
0342           }
0343         }  // end sim hit loop
0344 
0345         if (subid == 1) {  //<----------barrel
0346           fillBarrel(*pixeliter, *closestit, detId, theGeomDet, tTopo);
0347         }                  // end barrel
0348         if (subid == 2) {  // <-------forward
0349           fillForward(*pixeliter, *closestit, detId, theGeomDet, tTopo);
0350         }
0351 
0352       }  // end matched emtpy
0353 
0354       int NsimHit = matched.size();
0355       if (subid == 1) {  //<----------barrel
0356         for (unsigned int i = 0; i < 3; i++)
0357           if (tTopo->pxbLayer(detId) == i + 1)
0358             recHitNsimHitLayer[i]->Fill(NsimHit);
0359       }                  // end barrel
0360       if (subid == 2) {  // <-------forward
0361         if (tTopo->pxfDisk(detId) == 1)
0362           recHitNsimHitDisk1->Fill(NsimHit);
0363         else
0364           recHitNsimHitDisk2->Fill(NsimHit);
0365       }
0366     }  // <-----end rechit loop
0367   }    // <------ end detunit loop
0368 }
0369 
0370 void SiPixelRecHitsValid::fillBarrel(const SiPixelRecHit& recHit,
0371                                      const PSimHit& simHit,
0372                                      DetId detId,
0373                                      const PixelGeomDetUnit* theGeomDet,
0374                                      const TrackerTopology* tTopo) {
0375   const float cmtomicron = 10000.0;
0376 
0377   int bunch = simHit.eventId().bunchCrossing();
0378   int event = simHit.eventId().event();
0379 
0380   recHitBunchB->Fill(bunch);
0381   if (bunch == 0)
0382     recHitEventB->Fill(event);
0383 
0384   LocalPoint lp = recHit.localPosition();
0385   float lp_y = lp.y();
0386   float lp_x = lp.x();
0387 
0388   LocalError lerr = recHit.localPositionError();
0389   float lerr_x = sqrt(lerr.xx());
0390   float lerr_y = sqrt(lerr.yy());
0391 
0392   recHitYAllModules->Fill(lp_y);
0393 
0394   float sim_x1 = simHit.entryPoint().x();
0395   float sim_x2 = simHit.exitPoint().x();
0396   float sim_xpos = 0.5 * (sim_x1 + sim_x2);
0397   float res_x = (lp.x() - sim_xpos) * cmtomicron;
0398 
0399   recHitXResAllB->Fill(res_x);
0400 
0401   float sim_y1 = simHit.entryPoint().y();
0402   float sim_y2 = simHit.exitPoint().y();
0403   float sim_ypos = 0.5 * (sim_y1 + sim_y2);
0404   float res_y = (lp.y() - sim_ypos) * cmtomicron;
0405 
0406   recHitYResAllB->Fill(res_y);
0407 
0408   float pull_x = (lp_x - sim_xpos) / lerr_x;
0409   float pull_y = (lp_y - sim_ypos) / lerr_y;
0410 
0411   recHitXPullAllB->Fill(pull_x);
0412   recHitYPullAllB->Fill(pull_y);
0413 
0414   int rows = theGeomDet->specificTopology().nrows();
0415 
0416   if (rows == 160) {
0417     recHitXFullModules->Fill(lp_x);
0418   } else if (rows == 80) {
0419     recHitXHalfModules->Fill(lp_x);
0420   }
0421 
0422   float tmp1 = theGeomDet->surface().toGlobal(Local3DPoint(0., 0., 0.)).perp();
0423   float tmp2 = theGeomDet->surface().toGlobal(Local3DPoint(0., 0., 1.)).perp();
0424 
0425   if (tmp2 < tmp1) {  // flipped
0426     for (unsigned int i = 0; i < 3; i++) {
0427       if (tTopo->pxbLayer(detId) == i + 1) {
0428         recHitXResFlippedLadderLayers[i]->Fill(res_x);
0429         recHitXPullFlippedLadderLayers[i]->Fill(pull_x);
0430       }
0431     }
0432   } else {
0433     for (unsigned int i = 0; i < 3; i++) {
0434       if (tTopo->pxbLayer(detId) == i + 1) {
0435         recHitXResNonFlippedLadderLayers[i]->Fill(res_x);
0436         recHitXPullNonFlippedLadderLayers[i]->Fill(pull_x);
0437       }
0438     }
0439   }
0440 
0441   //get cluster
0442   SiPixelRecHit::ClusterRef const& clust = recHit.cluster();
0443 
0444   // fill module dependent info
0445   for (unsigned int i = 0; i < 8; i++) {
0446     if (tTopo->pxbModule(detId) == i + 1) {
0447       int sizeY = (*clust).sizeY();
0448       clustYSizeModule[i]->Fill(sizeY);
0449 
0450       if (tTopo->pxbLayer(detId) == 1) {
0451         float charge = (*clust).charge();
0452         clustChargeLayer1Modules[i]->Fill(charge);
0453         recHitYResLayer1Modules[i]->Fill(res_y);
0454         recHitYPullLayer1Modules[i]->Fill(pull_y);
0455       } else if (tTopo->pxbLayer(detId) == 2) {
0456         float charge = (*clust).charge();
0457         clustChargeLayer2Modules[i]->Fill(charge);
0458         recHitYResLayer2Modules[i]->Fill(res_y);
0459         recHitYPullLayer2Modules[i]->Fill(pull_y);
0460       } else if (tTopo->pxbLayer(detId) == 3) {
0461         float charge = (*clust).charge();
0462         clustChargeLayer3Modules[i]->Fill(charge);
0463         recHitYResLayer3Modules[i]->Fill(res_y);
0464         recHitYPullLayer3Modules[i]->Fill(pull_y);
0465       }
0466     }
0467   }
0468   int sizeX = (*clust).sizeX();
0469   if (tTopo->pxbLayer(detId) == 1)
0470     clustXSizeLayer[0]->Fill(sizeX);
0471   if (tTopo->pxbLayer(detId) == 2)
0472     clustXSizeLayer[1]->Fill(sizeX);
0473   if (tTopo->pxbLayer(detId) == 3)
0474     clustXSizeLayer[2]->Fill(sizeX);
0475 }
0476 
0477 void SiPixelRecHitsValid::fillForward(const SiPixelRecHit& recHit,
0478                                       const PSimHit& simHit,
0479                                       DetId detId,
0480                                       const PixelGeomDetUnit* theGeomDet,
0481                                       const TrackerTopology* tTopo) {
0482   int rows = theGeomDet->specificTopology().nrows();
0483   int cols = theGeomDet->specificTopology().ncolumns();
0484 
0485   const float cmtomicron = 10000.0;
0486 
0487   int bunch = simHit.eventId().bunchCrossing();
0488   int event = simHit.eventId().event();
0489 
0490   recHitBunchF->Fill(bunch);
0491   if (bunch == 0)
0492     recHitEventF->Fill(event);
0493 
0494   LocalPoint lp = recHit.localPosition();
0495   float lp_x = lp.x();
0496   float lp_y = lp.y();
0497 
0498   LocalError lerr = recHit.localPositionError();
0499   float lerr_x = sqrt(lerr.xx());
0500   float lerr_y = sqrt(lerr.yy());
0501 
0502   float sim_x1 = simHit.entryPoint().x();
0503   float sim_x2 = simHit.exitPoint().x();
0504   float sim_xpos = 0.5 * (sim_x1 + sim_x2);
0505 
0506   float sim_y1 = simHit.entryPoint().y();
0507   float sim_y2 = simHit.exitPoint().y();
0508   float sim_ypos = 0.5 * (sim_y1 + sim_y2);
0509 
0510   float pull_x = (lp_x - sim_xpos) / lerr_x;
0511   float pull_y = (lp_y - sim_ypos) / lerr_y;
0512 
0513   if (rows == 80) {
0514     recHitXPlaquetteSize1->Fill(lp_x);
0515   } else if (rows == 160) {
0516     recHitXPlaquetteSize2->Fill(lp_x);
0517   }
0518 
0519   if (cols == 104) {
0520     recHitYPlaquetteSize2->Fill(lp_y);
0521   } else if (cols == 156) {
0522     recHitYPlaquetteSize3->Fill(lp_y);
0523   } else if (cols == 208) {
0524     recHitYPlaquetteSize4->Fill(lp_y);
0525   } else if (cols == 260) {
0526     recHitYPlaquetteSize5->Fill(lp_y);
0527   }
0528 
0529   float res_x = (lp.x() - sim_xpos) * cmtomicron;
0530 
0531   recHitXResAllF->Fill(res_x);
0532   recHitXPullAllF->Fill(pull_x);
0533 
0534   float res_y = (lp.y() - sim_ypos) * cmtomicron;
0535 
0536   recHitYPullAllF->Fill(pull_y);
0537 
0538   // get cluster
0539   SiPixelRecHit::ClusterRef const& clust = recHit.cluster();
0540 
0541   // fill plaquette dependent info
0542   for (unsigned int i = 0; i < 7; i++) {
0543     if (tTopo->pxfModule(detId) == i + 1) {
0544       if (tTopo->pxfDisk(detId) == 1) {
0545         int sizeX = (*clust).sizeX();
0546         clustXSizeDisk1Plaquettes[i]->Fill(sizeX);
0547 
0548         int sizeY = (*clust).sizeY();
0549         clustYSizeDisk1Plaquettes[i]->Fill(sizeY);
0550 
0551         float charge = (*clust).charge();
0552         clustChargeDisk1Plaquettes[i]->Fill(charge);
0553 
0554         recHitXResDisk1Plaquettes[i]->Fill(res_x);
0555         recHitYResDisk1Plaquettes[i]->Fill(res_y);
0556 
0557         recHitXPullDisk1Plaquettes[i]->Fill(pull_x);
0558         recHitYPullDisk1Plaquettes[i]->Fill(pull_y);
0559       } else {
0560         int sizeX = (*clust).sizeX();
0561         clustXSizeDisk2Plaquettes[i]->Fill(sizeX);
0562 
0563         int sizeY = (*clust).sizeY();
0564         clustYSizeDisk2Plaquettes[i]->Fill(sizeY);
0565 
0566         float charge = (*clust).charge();
0567         clustChargeDisk2Plaquettes[i]->Fill(charge);
0568 
0569         recHitXResDisk2Plaquettes[i]->Fill(res_x);
0570         recHitYResDisk2Plaquettes[i]->Fill(res_y);
0571 
0572         recHitXPullDisk2Plaquettes[i]->Fill(pull_x);
0573         recHitYPullDisk2Plaquettes[i]->Fill(pull_y);
0574 
0575       }  // end else
0576     }    // end if module
0577     else if (tTopo->pxfPanel(detId) == 2 && (tTopo->pxfModule(detId) + 4) == i + 1) {
0578       if (tTopo->pxfDisk(detId) == 1) {
0579         int sizeX = (*clust).sizeX();
0580         clustXSizeDisk1Plaquettes[i]->Fill(sizeX);
0581 
0582         int sizeY = (*clust).sizeY();
0583         clustYSizeDisk1Plaquettes[i]->Fill(sizeY);
0584 
0585         float charge = (*clust).charge();
0586         clustChargeDisk1Plaquettes[i]->Fill(charge);
0587 
0588         recHitXResDisk1Plaquettes[i]->Fill(res_x);
0589         recHitYResDisk1Plaquettes[i]->Fill(res_y);
0590 
0591         recHitXPullDisk1Plaquettes[i]->Fill(pull_x);
0592         recHitYPullDisk1Plaquettes[i]->Fill(pull_y);
0593       } else {
0594         int sizeX = (*clust).sizeX();
0595         clustXSizeDisk2Plaquettes[i]->Fill(sizeX);
0596 
0597         int sizeY = (*clust).sizeY();
0598         clustYSizeDisk2Plaquettes[i]->Fill(sizeY);
0599 
0600         float charge = (*clust).charge();
0601         clustChargeDisk2Plaquettes[i]->Fill(charge);
0602 
0603         recHitXResDisk2Plaquettes[i]->Fill(res_x);
0604         recHitYResDisk2Plaquettes[i]->Fill(res_y);
0605 
0606         recHitXPullDisk2Plaquettes[i]->Fill(pull_x);
0607         recHitYPullDisk2Plaquettes[i]->Fill(pull_y);
0608 
0609       }  // end else
0610     }    // end else
0611   }      // end for
0612 }
0613 DEFINE_FWK_MODULE(SiPixelRecHitsValid);