CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_5/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/TrackerCommon/interface/TrackerTopology.h"
00023 #include "Geometry/Records/interface/IdealGeometryRecord.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   //Retrieve tracker topology from geometry
00259   edm::ESHandle<TrackerTopology> tTopoHand;
00260   es.get<IdealGeometryRecord>().get(tTopoHand);
00261   const TrackerTopology *tTopo=tTopoHand.product();
00262 
00263   LogInfo("EventInfo") << " Run = " << e.id().run() << " Event = " << e.id().event();
00264   if ( (int) e.id().event() % 1000 == 0 )
00265     cout << " Run = " << e.id().run() << " Event = " << e.id().event() << endl;
00266   
00267   //Get RecHits
00268   edm::Handle<SiPixelRecHitCollection> recHitColl;
00269   e.getByLabel( src_, recHitColl);
00270   
00271   //Get event setup
00272   edm::ESHandle<TrackerGeometry> geom;
00273   es.get<TrackerDigiGeometryRecord>().get(geom); 
00274   const TrackerGeometry& theTracker(*geom);
00275   
00276   TrackerHitAssociator associate( e, conf_ ); 
00277   
00278   //iterate over detunits
00279   for (TrackerGeometry::DetContainer::const_iterator it = geom->dets().begin(); it != geom->dets().end(); it++) 
00280     {
00281       DetId detId = ((*it)->geographicalId());
00282       unsigned int subid=detId.subdetId();
00283       
00284       if (! ((subid==1) || (subid==2))) continue;
00285       
00286       const PixelGeomDetUnit * theGeomDet = dynamic_cast<const PixelGeomDetUnit*>(theTracker.idToDet(detId) );
00287       
00288       SiPixelRecHitCollection::const_iterator pixeldet = recHitColl->find(detId);
00289       if (pixeldet == recHitColl->end()) continue;
00290       SiPixelRecHitCollection::DetSet pixelrechitRange = *pixeldet;
00291       SiPixelRecHitCollection::DetSet::const_iterator pixelrechitRangeIteratorBegin = pixelrechitRange.begin();
00292       SiPixelRecHitCollection::DetSet::const_iterator pixelrechitRangeIteratorEnd   = pixelrechitRange.end();
00293       SiPixelRecHitCollection::DetSet::const_iterator pixeliter = pixelrechitRangeIteratorBegin;
00294       std::vector<PSimHit> matched;
00295       
00296       //----Loop over rechits for this detId
00297       for ( ; pixeliter != pixelrechitRangeIteratorEnd; pixeliter++) 
00298         {
00299           matched.clear();
00300           matched = associate.associateHit(*pixeliter);
00301           
00302           if ( !matched.empty() ) 
00303             {
00304               float closest = 9999.9;
00305               std::vector<PSimHit>::const_iterator closestit = matched.begin();
00306               LocalPoint lp = pixeliter->localPosition();
00307               float rechit_x = lp.x();
00308               float rechit_y = lp.y();
00309 
00310               //loop over sim hits and fill closet
00311               for (std::vector<PSimHit>::const_iterator m = matched.begin(); m<matched.end(); m++) 
00312                 {
00313                   float sim_x1 = (*m).entryPoint().x();
00314                   float sim_x2 = (*m).exitPoint().x();
00315                   float sim_xpos = 0.5*(sim_x1+sim_x2);
00316 
00317                   float sim_y1 = (*m).entryPoint().y();
00318                   float sim_y2 = (*m).exitPoint().y();
00319                   float sim_ypos = 0.5*(sim_y1+sim_y2);
00320                   
00321                   float x_res = fabs(sim_xpos - rechit_x);
00322                   float y_res = fabs(sim_ypos - rechit_y);
00323                   
00324                   float dist = sqrt(x_res*x_res + y_res*y_res);
00325 
00326                   if ( dist < closest ) 
00327                     {
00328                       closest = x_res;
00329                       closestit = m;
00330                     }
00331                 } // end sim hit loop
00332               
00333               if (subid==1) 
00334                 { //<----------barrel
00335                   fillBarrel(*pixeliter, *closestit, detId, theGeomDet,tTopo);  
00336                 } // end barrel
00337               if (subid==2) 
00338                 { // <-------forward
00339                   fillForward(*pixeliter, *closestit, detId, theGeomDet,tTopo);
00340                 }
00341               
00342             } // end matched emtpy
00343         } // <-----end rechit loop 
00344     } // <------ end detunit loop
00345 }
00346 
00347 void SiPixelRecHitsValid::fillBarrel(const SiPixelRecHit& recHit, const PSimHit& simHit, 
00348                                      DetId detId, const PixelGeomDetUnit* theGeomDet,
00349                                      const TrackerTopology *tTopo) 
00350 {
00351   const float cmtomicron = 10000.0; 
00352   
00353   LocalPoint lp = recHit.localPosition();
00354   float lp_y = lp.y();  
00355   float lp_x = lp.x();
00356 
00357   LocalError lerr = recHit.localPositionError();
00358   float lerr_x = sqrt(lerr.xx());
00359   float lerr_y = sqrt(lerr.yy());
00360   
00361   recHitYAllModules->Fill(lp_y);
00362   
00363   float sim_x1 = simHit.entryPoint().x();
00364   float sim_x2 = simHit.exitPoint().x();
00365   float sim_xpos = 0.5*(sim_x1 + sim_x2);
00366   float res_x = (lp.x() - sim_xpos)*cmtomicron;
00367   
00368   recHitXResAllB->Fill(res_x);
00369   
00370   float sim_y1 = simHit.entryPoint().y();
00371   float sim_y2 = simHit.exitPoint().y();
00372   float sim_ypos = 0.5*(sim_y1 + sim_y2);
00373   float res_y = (lp.y() - sim_ypos)*cmtomicron;
00374   
00375   recHitYResAllB->Fill(res_y);
00376   
00377   float pull_x = ( lp_x - sim_xpos ) / lerr_x;
00378   float pull_y = ( lp_y - sim_ypos ) / lerr_y;
00379 
00380   recHitXPullAllB->Fill(pull_x);  
00381   recHitYPullAllB->Fill(pull_y);
00382 
00383   int rows = theGeomDet->specificTopology().nrows();
00384   
00385   if (rows == 160) 
00386     {
00387       recHitXFullModules->Fill(lp_x);
00388     } 
00389   else if (rows == 80) 
00390     {
00391       recHitXHalfModules->Fill(lp_x);
00392     }
00393   
00394   float tmp1 = theGeomDet->surface().toGlobal(Local3DPoint(0.,0.,0.)).perp();
00395   float tmp2 = theGeomDet->surface().toGlobal(Local3DPoint(0.,0.,1.)).perp();
00396   
00397   if (tmp2<tmp1) 
00398     { // flipped
00399       for (unsigned int i=0; i<3; i++) 
00400         {
00401           if (tTopo->pxbLayer(detId) == i+1) 
00402             {
00403               recHitXResFlippedLadderLayers[i]->Fill(res_x);
00404               recHitXPullFlippedLadderLayers[i]->Fill(pull_x);
00405             }
00406         }
00407     } 
00408   else 
00409     {
00410       for (unsigned int i=0; i<3; i++) 
00411         {
00412           if (tTopo->pxbLayer(detId) == i+1) 
00413             {
00414               recHitXResNonFlippedLadderLayers[i]->Fill(res_x);
00415               recHitXPullNonFlippedLadderLayers[i]->Fill(pull_x);
00416             }
00417         }
00418     }
00419   
00420   //get cluster
00421   SiPixelRecHit::ClusterRef const& clust = recHit.cluster();
00422   
00423   // fill module dependent info
00424   for (unsigned int i=0; i<8; i++) 
00425     {
00426       if (tTopo->pxbModule(detId) == i+1) 
00427         {
00428           int sizeY = (*clust).sizeY();
00429           clustYSizeModule[i]->Fill(sizeY);
00430           
00431           if (tTopo->pxbLayer(detId) == 1) 
00432             {
00433               float charge = (*clust).charge();
00434               clustChargeLayer1Modules[i]->Fill(charge);
00435               recHitYResLayer1Modules[i]->Fill(res_y);
00436               recHitYPullLayer1Modules[i]->Fill(pull_y);
00437             }
00438           else if (tTopo->pxbLayer(detId) == 2) 
00439             {
00440               float charge = (*clust).charge();
00441               clustChargeLayer2Modules[i]->Fill(charge);
00442               recHitYResLayer2Modules[i]->Fill(res_y);
00443               recHitYPullLayer2Modules[i]->Fill(pull_y);
00444             }
00445           else if (tTopo->pxbLayer(detId) == 3) 
00446             {
00447               float charge = (*clust).charge();
00448               clustChargeLayer3Modules[i]->Fill(charge);
00449               recHitYResLayer3Modules[i]->Fill(res_y);
00450               recHitYPullLayer3Modules[i]->Fill(pull_y);
00451             }
00452         }
00453     }
00454   int sizeX = (*clust).sizeX();
00455   if (tTopo->pxbLayer(detId) == 1) clustXSizeLayer[0]->Fill(sizeX);
00456   if (tTopo->pxbLayer(detId) == 2) clustXSizeLayer[1]->Fill(sizeX);
00457   if (tTopo->pxbLayer(detId) == 3) clustXSizeLayer[2]->Fill(sizeX);
00458 }
00459 
00460 void SiPixelRecHitsValid::fillForward(const SiPixelRecHit & recHit, const PSimHit & simHit, 
00461                                       DetId detId,const PixelGeomDetUnit * theGeomDet,
00462                                       const TrackerTopology *tTopo) 
00463 {
00464   int rows = theGeomDet->specificTopology().nrows();
00465   int cols = theGeomDet->specificTopology().ncolumns();
00466   
00467   const float cmtomicron = 10000.0;
00468   
00469   LocalPoint lp = recHit.localPosition();
00470   float lp_x = lp.x();
00471   float lp_y = lp.y();
00472   
00473   LocalError lerr = recHit.localPositionError();
00474   float lerr_x = sqrt(lerr.xx());
00475   float lerr_y = sqrt(lerr.yy());
00476 
00477   float sim_x1 = simHit.entryPoint().x();
00478   float sim_x2 = simHit.exitPoint().x();
00479   float sim_xpos = 0.5*(sim_x1 + sim_x2);
00480   
00481   float sim_y1 = simHit.entryPoint().y();
00482   float sim_y2 = simHit.exitPoint().y();
00483   float sim_ypos = 0.5*(sim_y1 + sim_y2);
00484   
00485   float pull_x = ( lp_x - sim_xpos ) / lerr_x;
00486   float pull_y = ( lp_y - sim_ypos ) / lerr_y;
00487 
00488 
00489   if (rows == 80) 
00490     {
00491       recHitXPlaquetteSize1->Fill(lp_x);
00492     } 
00493   else if (rows == 160) 
00494     {
00495       recHitXPlaquetteSize2->Fill(lp_x);
00496     }
00497   
00498   if (cols == 104) 
00499     {
00500       recHitYPlaquetteSize2->Fill(lp_y);
00501     } 
00502   else if (cols == 156) 
00503     {
00504       recHitYPlaquetteSize3->Fill(lp_y);
00505     } 
00506   else if (cols == 208) 
00507     {
00508       recHitYPlaquetteSize4->Fill(lp_y);
00509     } 
00510   else if (cols == 260) 
00511     {
00512       recHitYPlaquetteSize5->Fill(lp_y);
00513     }
00514   
00515   float res_x = (lp.x() - sim_xpos)*cmtomicron;
00516   
00517   recHitXResAllF->Fill(res_x);
00518   recHitXPullAllF->Fill(pull_x);
00519 
00520   float res_y = (lp.y() - sim_ypos)*cmtomicron;
00521   
00522   recHitYPullAllF->Fill(pull_y);
00523   
00524   // get cluster
00525   SiPixelRecHit::ClusterRef const& clust = recHit.cluster();
00526   
00527   // fill plaquette dependent info
00528   for (unsigned int i=0; i<7; i++) 
00529     {
00530       if (tTopo->pxfModule(detId) == i+1) 
00531         {
00532           if (tTopo->pxfDisk(detId) == 1) 
00533             {
00534               int sizeX = (*clust).sizeX();
00535               clustXSizeDisk1Plaquettes[i]->Fill(sizeX);
00536               
00537               int sizeY = (*clust).sizeY();
00538               clustYSizeDisk1Plaquettes[i]->Fill(sizeY);
00539               
00540               float charge = (*clust).charge();
00541               clustChargeDisk1Plaquettes[i]->Fill(charge);
00542               
00543               recHitXResDisk1Plaquettes[i]->Fill(res_x);
00544               recHitYResDisk1Plaquettes[i]->Fill(res_y);
00545 
00546               recHitXPullDisk1Plaquettes[i]->Fill(pull_x);
00547               recHitYPullDisk1Plaquettes[i]->Fill(pull_y);
00548             }
00549           else 
00550             {
00551               int sizeX = (*clust).sizeX();
00552               clustXSizeDisk2Plaquettes[i]->Fill(sizeX);
00553               
00554               int sizeY = (*clust).sizeY();
00555               clustYSizeDisk2Plaquettes[i]->Fill(sizeY);
00556               
00557               float charge = (*clust).charge();
00558               clustChargeDisk2Plaquettes[i]->Fill(charge);
00559               
00560               recHitXResDisk2Plaquettes[i]->Fill(res_x);
00561               recHitYResDisk2Plaquettes[i]->Fill(res_y);
00562 
00563               recHitXPullDisk2Plaquettes[i]->Fill(pull_x);
00564               recHitYPullDisk2Plaquettes[i]->Fill(pull_y);
00565               
00566             } // end else
00567         } // end if module
00568       else if (tTopo->pxfPanel(detId) == 2 && (tTopo->pxfModule(detId)+4) == i+1) 
00569         {
00570           if (tTopo->pxfDisk(detId) == 1) 
00571             {
00572               int sizeX = (*clust).sizeX();
00573               clustXSizeDisk1Plaquettes[i]->Fill(sizeX);
00574               
00575               int sizeY = (*clust).sizeY();
00576               clustYSizeDisk1Plaquettes[i]->Fill(sizeY);
00577               
00578               float charge = (*clust).charge();
00579               clustChargeDisk1Plaquettes[i]->Fill(charge);
00580               
00581               recHitXResDisk1Plaquettes[i]->Fill(res_x);
00582               recHitYResDisk1Plaquettes[i]->Fill(res_y);
00583 
00584               recHitXPullDisk1Plaquettes[i]->Fill(pull_x);
00585               recHitYPullDisk1Plaquettes[i]->Fill(pull_y);
00586             }
00587           else 
00588             {
00589               int sizeX = (*clust).sizeX();
00590               clustXSizeDisk2Plaquettes[i]->Fill(sizeX);
00591               
00592               int sizeY = (*clust).sizeY();
00593               clustYSizeDisk2Plaquettes[i]->Fill(sizeY);
00594               
00595               float charge = (*clust).charge();
00596               clustChargeDisk2Plaquettes[i]->Fill(charge);
00597               
00598               recHitXResDisk2Plaquettes[i]->Fill(res_x);
00599               recHitYResDisk2Plaquettes[i]->Fill(res_y);
00600 
00601               recHitXPullDisk2Plaquettes[i]->Fill(pull_x);
00602               recHitYPullDisk2Plaquettes[i]->Fill(pull_y);
00603 
00604             } // end else
00605         } // end else
00606     } // end for
00607 }