CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch9/src/Validation/TrackerRecHits/plugins/SiPixelRecHitsValid.cc

Go to the documentation of this file.
00001 // SiPixelRecHitsValid.cc
00002 // Description: see SiPixelRecHitsValid.h
00003 // Author: Jason Shaev, JHU
00004 // Created 6/7/06
00005 //
00006 // G. Giurgiu, JHU (ggiurgiu@pha.jhu.edu)
00007 //             added pull distributions (12/27/06)
00008 //--------------------------------
00009 
00010 #include "Validation/TrackerRecHits/interface/SiPixelRecHitsValid.h"
00011 
00012 #include "SimDataFormats/TrackingHit/interface/PSimHitContainer.h"
00013 
00014 #include "DataFormats/GeometryVector/interface/LocalPoint.h"
00015 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
00016 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
00017 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
00018 #include "Geometry/CommonDetUnit/interface/GeomDetUnit.h"
00019 
00020 #include "DataFormats/TrackerRecHit2D/interface/SiPixelRecHit.h"
00021 #include "DataFormats/Common/interface/OwnVector.h"
00022 #include "DataFormats/SiPixelDetId/interface/PXBDetId.h"
00023 #include "DataFormats/SiPixelDetId/interface/PXFDetId.h"
00024 #include "DataFormats/SiPixelDetId/interface/PixelSubdetector.h"
00025 #include "DataFormats/Common/interface/Ref.h"
00026 #include "DataFormats/Common/interface/DetSetVector.h"
00027 
00028 #include "SimTracker/TrackerHitAssociation/interface/TrackerHitAssociator.h"
00029 
00030 #include <math.h>
00031 #include "DQMServices/Core/interface/DQMStore.h"
00032 
00033 using namespace std;
00034 using namespace edm;
00035 
00036 
00037 SiPixelRecHitsValid::SiPixelRecHitsValid(const ParameterSet& ps): 
00038   dbe_(0), 
00039   conf_(ps),
00040   src_( ps.getParameter<edm::InputTag>( "src" ) ) 
00041 {
00042   outputFile_ = ps.getUntrackedParameter<string>("outputFile", "pixelrechitshisto.root");
00043   dbe_ = Service<DQMStore>().operator->();
00044   //dbe_->showDirStructure();
00045   dbe_->setCurrentFolder("TrackerRecHitsV/TrackerRecHits/Pixel/clustBPIX");
00046   
00047   Char_t histo[200];
00048   
00049   // ---------------------------------------------------------------
00050   // All histograms that depend on plaquette number have 7 indexes.
00051   // The first 4 (0-3) correspond to Panel 1 plaquettes 1-4.
00052   // The last 3 (4-6) correspond to Panel 2 plaquettes 1-3.
00053   // ---------------------------------------------------------------
00054   
00055   //Cluster y-size by module number for barrel
00056   for (int i=0; i<8; i++) {
00057     sprintf(histo, "Clust_y_size_Module%d", i+1);
00058     clustYSizeModule[i] = dbe_->book1D(histo,"Cluster y-size by Module", 20, 0.5, 20.5); 
00059   } // end for
00060   
00061   //Cluster x-size by layer for barrel
00062   for (int i=0; i<3; i++) {
00063     sprintf(histo, "Clust_x_size_Layer%d", i+1);
00064     clustXSizeLayer[i] = dbe_->book1D(histo,"Cluster x-size by Layer", 20, 0.5, 20.5);
00065   } // end for
00066   
00067   //Cluster charge by module for 3 layers of barrel
00068   for (int i=0; i<8; i++) {
00069     //Cluster charge by module for Layer1
00070     sprintf(histo, "Clust_charge_Layer1_Module%d", i+1);
00071     clustChargeLayer1Modules[i] = dbe_->book1D(histo, "Cluster charge Layer 1 by Module", 50, 0., 200000.);
00072     
00073     //Cluster charge by module for Layer2
00074     sprintf(histo, "Clust_charge_Layer2_Module%d", i+1);
00075     clustChargeLayer2Modules[i] = dbe_->book1D(histo, "Cluster charge Layer 2 by Module", 50, 0., 200000.);
00076     
00077     //Cluster charge by module for Layer3
00078     sprintf(histo, "Clust_charge_Layer3_Module%d", i+1);
00079     clustChargeLayer3Modules[i] = dbe_->book1D(histo, "Cluster charge Layer 3 by Module",50, 0., 200000.);      
00080   } // end for
00081   
00082   dbe_->setCurrentFolder("TrackerRecHitsV/TrackerRecHits/Pixel/clustFPIX");
00083   //Cluster x-size, y-size, and charge by plaquette for Disks in Forward
00084   for (int i=0; i<7; i++) {
00085     //Cluster x-size for Disk1 by Plaquette
00086     sprintf(histo, "Clust_x_size_Disk1_Plaquette%d", i+1);
00087     clustXSizeDisk1Plaquettes[i] = dbe_->book1D(histo, "Cluster X-size for Disk1 by Plaquette", 20, 0.5, 20.5);
00088     
00089     //Cluster x-size for Disk2 by Plaquette
00090     sprintf(histo, "Clust_x_size_Disk2_Plaquette%d", i+1);
00091     clustXSizeDisk2Plaquettes[i] = dbe_->book1D(histo, "Cluster X-size for Disk2 by Plaquette", 20, 0.5, 20.5);
00092     
00093     //Cluster y-size for Disk1 by Plaquette
00094     sprintf(histo, "Clust_y_size_Disk1_Plaquette%d", i+1);
00095     clustYSizeDisk1Plaquettes[i] = dbe_->book1D(histo, "Cluster Y-size for Disk1 by Plaquette", 20, 0.5, 20.5);
00096     
00097     //Cluster y-size for Disk2 by Plaquette
00098     sprintf(histo, "Clust_y_size_Disk2_Plaquette%d", i+1);
00099     clustYSizeDisk2Plaquettes[i] = dbe_->book1D(histo, "Cluster Y-size for Disk2 by Plaquette", 20, 0.5, 20.5);
00100     
00101     //Cluster charge for Disk1 by Plaquette
00102     sprintf(histo, "Clust_charge_Disk1_Plaquette%d", i+1);
00103     clustChargeDisk1Plaquettes[i] = dbe_->book1D(histo, "Cluster charge for Disk1 by Plaquette", 50, 0., 200000.);
00104     
00105     //Cluster charge for Disk2 by Plaquette
00106     sprintf(histo, "Clust_charge_Disk2_Plaquette%d", i+1);
00107     clustChargeDisk2Plaquettes[i] = dbe_->book1D(histo, "Cluster charge for Disk2 by Plaquette", 50, 0., 200000.);
00108   } // end for
00109   
00110 
00111 
00112   dbe_->setCurrentFolder("TrackerRecHitsV/TrackerRecHits/Pixel/recHitBPIX");
00113   //RecHit X Resolution all barrel hits
00114   recHitXResAllB = dbe_->book1D("RecHit_xres_b_All","RecHit X Res All Modules in Barrel", 100, -200., 200.);
00115   
00116   //RecHit Y Resolution all barrel hits
00117   recHitYResAllB = dbe_->book1D("RecHit_yres_b_All","RecHit Y Res All Modules in Barrel", 100, -200., 200.);
00118   
00119   //RecHit X distribution for full modules for barrel
00120   recHitXFullModules = dbe_->book1D("RecHit_x_FullModules", "RecHit X distribution for full modules", 100,-2., 2.);
00121   
00122   //RecHit X distribution for half modules for barrel
00123   recHitXHalfModules = dbe_->book1D("RecHit_x_HalfModules", "RecHit X distribution for half modules", 100, -1., 1.);
00124   
00125   //RecHit Y distribution all modules for barrel
00126   recHitYAllModules = dbe_->book1D("RecHit_y_AllModules", "RecHit Y distribution for all modules", 100, -4., 4.);
00127   
00128   //RecHit X resolution for flipped and unflipped ladders by layer for barrel
00129   for (int i=0; i<3; i++) {
00130     //RecHit X resolution for flipped ladders by layer
00131     sprintf(histo, "RecHit_XRes_FlippedLadder_Layer%d", i+1);
00132     recHitXResFlippedLadderLayers[i] = dbe_->book1D(histo, "RecHit XRes Flipped Ladders by Layer", 100, -200., 200.);
00133     
00134     //RecHit X resolution for unflipped ladders by layer
00135     sprintf(histo, "RecHit_XRes_UnFlippedLadder_Layer%d", i+1);
00136     recHitXResNonFlippedLadderLayers[i] = dbe_->book1D(histo, "RecHit XRes NonFlipped Ladders by Layer", 100, -200., 200.);
00137   } // end for
00138   
00139   //RecHit Y resolutions for layers by module for barrel
00140   for (int i=0; i<8; i++) {
00141     //Rec Hit Y resolution by module for Layer1
00142     sprintf(histo, "RecHit_YRes_Layer1_Module%d", i+1);
00143     recHitYResLayer1Modules[i] = dbe_->book1D(histo, "RecHit YRes Layer1 by module", 100, -200., 200.);
00144     
00145     //RecHit Y resolution by module for Layer2
00146     sprintf(histo, "RecHit_YRes_Layer2_Module%d", i+1);
00147     recHitYResLayer2Modules[i] = dbe_->book1D(histo, "RecHit YRes Layer2 by module", 100, -200., 200.);
00148     
00149     //RecHit Y resolution by module for Layer3
00150     sprintf(histo, "RecHit_YRes_Layer3_Module%d", i+1);
00151     recHitYResLayer3Modules[i] = dbe_->book1D(histo, "RecHit YRes Layer3 by module", 100, -200., 200.); 
00152   } // end for
00153   
00154   dbe_->setCurrentFolder("TrackerRecHitsV/TrackerRecHits/Pixel/recHitFPIX");
00155   //RecHit X resolution all plaquettes
00156   recHitXResAllF = dbe_->book1D("RecHit_xres_f_All", "RecHit X Res All in Forward", 100, -200., 200.);
00157   
00158   //RecHit Y resolution all plaquettes
00159   recHitYResAllF = dbe_->book1D("RecHit_yres_f_All", "RecHit Y Res All in Forward", 100, -200., 200.);
00160   
00161   //RecHit X distribution for plaquette with x-size 1 in forward
00162   recHitXPlaquetteSize1 = dbe_->book1D("RecHit_x_Plaquette_xsize1", "RecHit X Distribution for plaquette x-size1", 100, -2., 2.);
00163   
00164   //RecHit X distribution for plaquette with x-size 2 in forward
00165   recHitXPlaquetteSize2 = dbe_->book1D("RecHit_x_Plaquette_xsize2", "RecHit X Distribution for plaquette x-size2", 100, -2., 2.);
00166   
00167   //RecHit Y distribution for plaquette with y-size 2 in forward
00168   recHitYPlaquetteSize2 = dbe_->book1D("RecHit_y_Plaquette_ysize2", "RecHit Y Distribution for plaquette y-size2", 100, -4., 4.);
00169   
00170   //RecHit Y distribution for plaquette with y-size 3 in forward
00171   recHitYPlaquetteSize3 = dbe_->book1D("RecHit_y_Plaquette_ysize3", "RecHit Y Distribution for plaquette y-size3", 100, -4., 4.);
00172   
00173   //RecHit Y distribution for plaquette with y-size 4 in forward
00174   recHitYPlaquetteSize4 = dbe_->book1D("RecHit_y_Plaquette_ysize4", "RecHit Y Distribution for plaquette y-size4", 100, -4., 4.);
00175   
00176   //RecHit Y distribution for plaquette with y-size 5 in forward
00177   recHitYPlaquetteSize5 = dbe_->book1D("RecHit_y_Plaquette_ysize5", "RecHit Y Distribution for plaquette y-size5", 100, -4., 4.);
00178   
00179   //X and Y resolutions for both disks by plaquette in forward
00180   for (int i=0; i<7; i++) {
00181     //X resolution for Disk1 by plaquette
00182     sprintf(histo, "RecHit_XRes_Disk1_Plaquette%d", i+1);
00183     recHitXResDisk1Plaquettes[i] = dbe_->book1D(histo, "RecHit XRes Disk1 by plaquette", 100, -200., 200.); 
00184     //X resolution for Disk2 by plaquette
00185     sprintf(histo, "RecHit_XRes_Disk2_Plaquette%d", i+1);
00186     recHitXResDisk2Plaquettes[i] = dbe_->book1D(histo, "RecHit XRes Disk2 by plaquette", 100, -200., 200.);  
00187     
00188     //Y resolution for Disk1 by plaquette
00189     sprintf(histo, "RecHit_YRes_Disk1_Plaquette%d", i+1);
00190     recHitYResDisk1Plaquettes[i] = dbe_->book1D(histo, "RecHit YRes Disk1 by plaquette", 100, -200., 200.);
00191     //Y resolution for Disk2 by plaquette
00192     sprintf(histo, "RecHit_YRes_Disk2_Plaquette%d", i+1);
00193     recHitYResDisk2Plaquettes[i] = dbe_->book1D(histo, "RecHit YRes Disk2 by plaquette", 100, -200., 200.);
00194     
00195   }
00196 
00197 
00198   dbe_->setCurrentFolder("TrackerRecHitsV/TrackerRecHits/Pixel/recHitPullsBPIX");
00199   recHitXPullAllB        = dbe_->book1D("RecHit_xres_b_All"       , "RecHit X Pull All Modules in Barrel"        , 100, -10.0, 10.0);
00200   recHitYPullAllB        = dbe_->book1D("RecHit_yres_b_All"       , "RecHit Y Pull All Modules in Barrel"        , 100, -10.0, 10.0);
00201 
00202   for (int i=0; i<3; i++) 
00203     {
00204       sprintf(histo, "RecHit_XPull_FlippedLadder_Layer%d", i+1);
00205       recHitXPullFlippedLadderLayers[i] = dbe_->book1D(histo, "RecHit XPull Flipped Ladders by Layer", 100, -10.0, 10.0);
00206       
00207       sprintf(histo, "RecHit_XPull_UnFlippedLadder_Layer%d", i+1);
00208       recHitXPullNonFlippedLadderLayers[i] = dbe_->book1D(histo, "RecHit XPull NonFlipped Ladders by Layer", 100, -10.0, 10.0);
00209     }
00210   
00211   for (int i=0; i<8; i++) 
00212     {
00213       sprintf(histo, "RecHit_YPull_Layer1_Module%d", i+1);
00214       recHitYPullLayer1Modules[i] = dbe_->book1D(histo, "RecHit YPull Layer1 by module", 100, -10.0, 10.0);
00215       
00216       sprintf(histo, "RecHit_YPull_Layer2_Module%d", i+1);
00217       recHitYPullLayer2Modules[i] = dbe_->book1D(histo, "RecHit YPull Layer2 by module", 100, -10.0, 10.0);
00218       
00219       sprintf(histo, "RecHit_YPull_Layer3_Module%d", i+1);
00220       recHitYPullLayer3Modules[i] = dbe_->book1D(histo, "RecHit YPull Layer3 by module", 100, -10.0, 10.0); 
00221     }
00222   
00223   dbe_->setCurrentFolder("TrackerRecHitsV/TrackerRecHits/Pixel/recHitPullsFPIX");
00224   recHitXPullAllF = dbe_->book1D("RecHit_XPull_f_All", "RecHit X Pull All in Forward", 100, -10.0, 10.0);
00225   
00226   recHitYPullAllF = dbe_->book1D("RecHit_YPull_f_All", "RecHit Y Pull All in Forward", 100, -10.0, 10.0);
00227   
00228   for (int i=0; i<7; i++) 
00229     {
00230       sprintf(histo, "RecHit_XPull_Disk1_Plaquette%d", i+1);
00231       recHitXPullDisk1Plaquettes[i] = dbe_->book1D(histo, "RecHit XPull Disk1 by plaquette", 100, -10.0, 10.0); 
00232       sprintf(histo, "RecHit_XPull_Disk2_Plaquette%d", i+1);
00233       recHitXPullDisk2Plaquettes[i] = dbe_->book1D(histo, "RecHit XPull Disk2 by plaquette", 100, -10.0, 10.0);  
00234       
00235       sprintf(histo, "RecHit_YPull_Disk1_Plaquette%d", i+1);
00236       recHitYPullDisk1Plaquettes[i] = dbe_->book1D(histo, "RecHit YPull Disk1 by plaquette", 100, -10.0, 10.0);
00237       
00238       sprintf(histo, "RecHit_YPull_Disk2_Plaquette%d", i+1);
00239       recHitYPullDisk2Plaquettes[i] = dbe_->book1D(histo, "RecHit YPull Disk2 by plaquette", 100, -10.0, 10.0);
00240     }
00241 
00242 }
00243 
00244 SiPixelRecHitsValid::~SiPixelRecHitsValid() {
00245 }
00246 
00247 void SiPixelRecHitsValid::beginJob() {
00248   
00249 }
00250 
00251 void SiPixelRecHitsValid::endJob() {
00252   if ( outputFile_.size() != 0 && dbe_ ) dbe_->save(outputFile_);
00253 }
00254 
00255 void SiPixelRecHitsValid::analyze(const edm::Event& e, const edm::EventSetup& es) 
00256 {
00257   
00258   LogInfo("EventInfo") << " Run = " << e.id().run() << " Event = " << e.id().event();
00259   if ( (int) e.id().event() % 1000 == 0 )
00260     cout << " Run = " << e.id().run() << " Event = " << e.id().event() << endl;
00261   
00262   //Get RecHits
00263   edm::Handle<SiPixelRecHitCollection> recHitColl;
00264   e.getByLabel( src_, recHitColl);
00265   
00266   //Get event setup
00267   edm::ESHandle<TrackerGeometry> geom;
00268   es.get<TrackerDigiGeometryRecord>().get(geom); 
00269   const TrackerGeometry& theTracker(*geom);
00270   
00271   TrackerHitAssociator associate( e, conf_ ); 
00272   
00273   //iterate over detunits
00274   for (TrackerGeometry::DetContainer::const_iterator it = geom->dets().begin(); it != geom->dets().end(); it++) 
00275     {
00276       DetId detId = ((*it)->geographicalId());
00277       unsigned int subid=detId.subdetId();
00278       
00279       if (! ((subid==1) || (subid==2))) continue;
00280       
00281       const PixelGeomDetUnit * theGeomDet = dynamic_cast<const PixelGeomDetUnit*>(theTracker.idToDet(detId) );
00282       
00283       SiPixelRecHitCollection::const_iterator pixeldet = recHitColl->find(detId);
00284       if (pixeldet == recHitColl->end()) continue;
00285       SiPixelRecHitCollection::DetSet pixelrechitRange = *pixeldet;
00286       SiPixelRecHitCollection::DetSet::const_iterator pixelrechitRangeIteratorBegin = pixelrechitRange.begin();
00287       SiPixelRecHitCollection::DetSet::const_iterator pixelrechitRangeIteratorEnd   = pixelrechitRange.end();
00288       SiPixelRecHitCollection::DetSet::const_iterator pixeliter = pixelrechitRangeIteratorBegin;
00289       std::vector<PSimHit> matched;
00290       
00291       //----Loop over rechits for this detId
00292       for ( ; pixeliter != pixelrechitRangeIteratorEnd; pixeliter++) 
00293         {
00294           matched.clear();
00295           matched = associate.associateHit(*pixeliter);
00296           
00297           if ( !matched.empty() ) 
00298             {
00299               float closest = 9999.9;
00300               std::vector<PSimHit>::const_iterator closestit = matched.begin();
00301               LocalPoint lp = pixeliter->localPosition();
00302               float rechit_x = lp.x();
00303               float rechit_y = lp.y();
00304 
00305               //loop over sim hits and fill closet
00306               for (std::vector<PSimHit>::const_iterator m = matched.begin(); m<matched.end(); m++) 
00307                 {
00308                   float sim_x1 = (*m).entryPoint().x();
00309                   float sim_x2 = (*m).exitPoint().x();
00310                   float sim_xpos = 0.5*(sim_x1+sim_x2);
00311 
00312                   float sim_y1 = (*m).entryPoint().y();
00313                   float sim_y2 = (*m).exitPoint().y();
00314                   float sim_ypos = 0.5*(sim_y1+sim_y2);
00315                   
00316                   float x_res = fabs(sim_xpos - rechit_x);
00317                   float y_res = fabs(sim_ypos - rechit_y);
00318                   
00319                   float dist = sqrt(x_res*x_res + y_res*y_res);
00320 
00321                   if ( dist < closest ) 
00322                     {
00323                       closest = x_res;
00324                       closestit = m;
00325                     }
00326                 } // end sim hit loop
00327               
00328               if (subid==1) 
00329                 { //<----------barrel
00330                   fillBarrel(*pixeliter, *closestit, detId, theGeomDet);        
00331                 } // end barrel
00332               if (subid==2) 
00333                 { // <-------forward
00334                   fillForward(*pixeliter, *closestit, detId, theGeomDet);
00335                 }
00336               
00337             } // end matched emtpy
00338         } // <-----end rechit loop 
00339     } // <------ end detunit loop
00340 }
00341 
00342 void SiPixelRecHitsValid::fillBarrel(const SiPixelRecHit& recHit, const PSimHit& simHit, 
00343                                      DetId detId, const PixelGeomDetUnit* theGeomDet) 
00344 {
00345   const float cmtomicron = 10000.0; 
00346   
00347   LocalPoint lp = recHit.localPosition();
00348   float lp_y = lp.y();  
00349   float lp_x = lp.x();
00350 
00351   LocalError lerr = recHit.localPositionError();
00352   float lerr_x = sqrt(lerr.xx());
00353   float lerr_y = sqrt(lerr.yy());
00354   
00355   recHitYAllModules->Fill(lp_y);
00356   
00357   float sim_x1 = simHit.entryPoint().x();
00358   float sim_x2 = simHit.exitPoint().x();
00359   float sim_xpos = 0.5*(sim_x1 + sim_x2);
00360   float res_x = (lp.x() - sim_xpos)*cmtomicron;
00361   
00362   recHitXResAllB->Fill(res_x);
00363   
00364   float sim_y1 = simHit.entryPoint().y();
00365   float sim_y2 = simHit.exitPoint().y();
00366   float sim_ypos = 0.5*(sim_y1 + sim_y2);
00367   float res_y = (lp.y() - sim_ypos)*cmtomicron;
00368   
00369   recHitYResAllB->Fill(res_y);
00370   
00371   float pull_x = ( lp_x - sim_xpos ) / lerr_x;
00372   float pull_y = ( lp_y - sim_ypos ) / lerr_y;
00373 
00374   recHitXPullAllB->Fill(pull_x);  
00375   recHitYPullAllB->Fill(pull_y);
00376 
00377   int rows = theGeomDet->specificTopology().nrows();
00378   
00379   if (rows == 160) 
00380     {
00381       recHitXFullModules->Fill(lp_x);
00382     } 
00383   else if (rows == 80) 
00384     {
00385       recHitXHalfModules->Fill(lp_x);
00386     }
00387   
00388   float tmp1 = theGeomDet->surface().toGlobal(Local3DPoint(0.,0.,0.)).perp();
00389   float tmp2 = theGeomDet->surface().toGlobal(Local3DPoint(0.,0.,1.)).perp();
00390   
00391   if (tmp2<tmp1) 
00392     { // flipped
00393       for (unsigned int i=0; i<3; i++) 
00394         {
00395           if (PXBDetId(detId).layer() == i+1) 
00396             {
00397               recHitXResFlippedLadderLayers[i]->Fill(res_x);
00398               recHitXPullFlippedLadderLayers[i]->Fill(pull_x);
00399             }
00400         }
00401     } 
00402   else 
00403     {
00404       for (unsigned int i=0; i<3; i++) 
00405         {
00406           if (PXBDetId(detId).layer() == i+1) 
00407             {
00408               recHitXResNonFlippedLadderLayers[i]->Fill(res_x);
00409               recHitXPullNonFlippedLadderLayers[i]->Fill(pull_x);
00410             }
00411         }
00412     }
00413   
00414   //get cluster
00415   SiPixelRecHit::ClusterRef const& clust = recHit.cluster();
00416   
00417   // fill module dependent info
00418   for (unsigned int i=0; i<8; i++) 
00419     {
00420       if (PXBDetId(detId).module() == i+1) 
00421         {
00422           int sizeY = (*clust).sizeY();
00423           clustYSizeModule[i]->Fill(sizeY);
00424           
00425           if (PXBDetId(detId).layer() == 1) 
00426             {
00427               float charge = (*clust).charge();
00428               clustChargeLayer1Modules[i]->Fill(charge);
00429               recHitYResLayer1Modules[i]->Fill(res_y);
00430               recHitYPullLayer1Modules[i]->Fill(pull_y);
00431             }
00432           else if (PXBDetId(detId).layer() == 2) 
00433             {
00434               float charge = (*clust).charge();
00435               clustChargeLayer2Modules[i]->Fill(charge);
00436               recHitYResLayer2Modules[i]->Fill(res_y);
00437               recHitYPullLayer2Modules[i]->Fill(pull_y);
00438             }
00439           else if (PXBDetId(detId).layer() == 3) 
00440             {
00441               float charge = (*clust).charge();
00442               clustChargeLayer3Modules[i]->Fill(charge);
00443               recHitYResLayer3Modules[i]->Fill(res_y);
00444               recHitYPullLayer3Modules[i]->Fill(pull_y);
00445             }
00446         }
00447     }
00448   int sizeX = (*clust).sizeX();
00449   if (PXBDetId(detId).layer() == 1) clustXSizeLayer[0]->Fill(sizeX);
00450   if (PXBDetId(detId).layer() == 2) clustXSizeLayer[1]->Fill(sizeX);
00451   if (PXBDetId(detId).layer() == 3) clustXSizeLayer[2]->Fill(sizeX);
00452 }
00453 
00454 void SiPixelRecHitsValid::fillForward(const SiPixelRecHit & recHit, const PSimHit & simHit, 
00455                                       DetId detId,const PixelGeomDetUnit * theGeomDet ) 
00456 {
00457   int rows = theGeomDet->specificTopology().nrows();
00458   int cols = theGeomDet->specificTopology().ncolumns();
00459   
00460   const float cmtomicron = 10000.0;
00461   
00462   LocalPoint lp = recHit.localPosition();
00463   float lp_x = lp.x();
00464   float lp_y = lp.y();
00465   
00466   LocalError lerr = recHit.localPositionError();
00467   float lerr_x = sqrt(lerr.xx());
00468   float lerr_y = sqrt(lerr.yy());
00469 
00470   float sim_x1 = simHit.entryPoint().x();
00471   float sim_x2 = simHit.exitPoint().x();
00472   float sim_xpos = 0.5*(sim_x1 + sim_x2);
00473   
00474   float sim_y1 = simHit.entryPoint().y();
00475   float sim_y2 = simHit.exitPoint().y();
00476   float sim_ypos = 0.5*(sim_y1 + sim_y2);
00477   
00478   float pull_x = ( lp_x - sim_xpos ) / lerr_x;
00479   float pull_y = ( lp_y - sim_ypos ) / lerr_y;
00480 
00481 
00482   if (rows == 80) 
00483     {
00484       recHitXPlaquetteSize1->Fill(lp_x);
00485     } 
00486   else if (rows == 160) 
00487     {
00488       recHitXPlaquetteSize2->Fill(lp_x);
00489     }
00490   
00491   if (cols == 104) 
00492     {
00493       recHitYPlaquetteSize2->Fill(lp_y);
00494     } 
00495   else if (cols == 156) 
00496     {
00497       recHitYPlaquetteSize3->Fill(lp_y);
00498     } 
00499   else if (cols == 208) 
00500     {
00501       recHitYPlaquetteSize4->Fill(lp_y);
00502     } 
00503   else if (cols == 260) 
00504     {
00505       recHitYPlaquetteSize5->Fill(lp_y);
00506     }
00507   
00508   float res_x = (lp.x() - sim_xpos)*cmtomicron;
00509   
00510   recHitXResAllF->Fill(res_x);
00511   recHitXPullAllF->Fill(pull_x);
00512 
00513   float res_y = (lp.y() - sim_ypos)*cmtomicron;
00514   
00515   recHitYPullAllF->Fill(pull_y);
00516   
00517   // get cluster
00518   SiPixelRecHit::ClusterRef const& clust = recHit.cluster();
00519   
00520   // fill plaquette dependent info
00521   for (unsigned int i=0; i<7; i++) 
00522     {
00523       if (PXFDetId(detId).module() == i+1) 
00524         {
00525           if (PXFDetId(detId).disk() == 1) 
00526             {
00527               int sizeX = (*clust).sizeX();
00528               clustXSizeDisk1Plaquettes[i]->Fill(sizeX);
00529               
00530               int sizeY = (*clust).sizeY();
00531               clustYSizeDisk1Plaquettes[i]->Fill(sizeY);
00532               
00533               float charge = (*clust).charge();
00534               clustChargeDisk1Plaquettes[i]->Fill(charge);
00535               
00536               recHitXResDisk1Plaquettes[i]->Fill(res_x);
00537               recHitYResDisk1Plaquettes[i]->Fill(res_y);
00538 
00539               recHitXPullDisk1Plaquettes[i]->Fill(pull_x);
00540               recHitYPullDisk1Plaquettes[i]->Fill(pull_y);
00541             }
00542           else 
00543             {
00544               int sizeX = (*clust).sizeX();
00545               clustXSizeDisk2Plaquettes[i]->Fill(sizeX);
00546               
00547               int sizeY = (*clust).sizeY();
00548               clustYSizeDisk2Plaquettes[i]->Fill(sizeY);
00549               
00550               float charge = (*clust).charge();
00551               clustChargeDisk2Plaquettes[i]->Fill(charge);
00552               
00553               recHitXResDisk2Plaquettes[i]->Fill(res_x);
00554               recHitYResDisk2Plaquettes[i]->Fill(res_y);
00555 
00556               recHitXPullDisk2Plaquettes[i]->Fill(pull_x);
00557               recHitYPullDisk2Plaquettes[i]->Fill(pull_y);
00558               
00559             } // end else
00560         } // end if module
00561       else if (PXFDetId(detId).panel() == 2 && (PXFDetId(detId).module()+4) == i+1) 
00562         {
00563           if (PXFDetId(detId).disk() == 1) 
00564             {
00565               int sizeX = (*clust).sizeX();
00566               clustXSizeDisk1Plaquettes[i]->Fill(sizeX);
00567               
00568               int sizeY = (*clust).sizeY();
00569               clustYSizeDisk1Plaquettes[i]->Fill(sizeY);
00570               
00571               float charge = (*clust).charge();
00572               clustChargeDisk1Plaquettes[i]->Fill(charge);
00573               
00574               recHitXResDisk1Plaquettes[i]->Fill(res_x);
00575               recHitYResDisk1Plaquettes[i]->Fill(res_y);
00576 
00577               recHitXPullDisk1Plaquettes[i]->Fill(pull_x);
00578               recHitYPullDisk1Plaquettes[i]->Fill(pull_y);
00579             }
00580           else 
00581             {
00582               int sizeX = (*clust).sizeX();
00583               clustXSizeDisk2Plaquettes[i]->Fill(sizeX);
00584               
00585               int sizeY = (*clust).sizeY();
00586               clustYSizeDisk2Plaquettes[i]->Fill(sizeY);
00587               
00588               float charge = (*clust).charge();
00589               clustChargeDisk2Plaquettes[i]->Fill(charge);
00590               
00591               recHitXResDisk2Plaquettes[i]->Fill(res_x);
00592               recHitYResDisk2Plaquettes[i]->Fill(res_y);
00593 
00594               recHitXPullDisk2Plaquettes[i]->Fill(pull_x);
00595               recHitYPullDisk2Plaquettes[i]->Fill(pull_y);
00596 
00597             } // end else
00598         } // end else
00599     } // end for
00600 }