CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_2_9/src/DQM/SiPixelMonitorRecHit/src/SiPixelRecHitModule.cc

Go to the documentation of this file.
00001 #include "DQM/SiPixelMonitorRecHit/interface/SiPixelRecHitModule.h"
00002 #include "DQMServices/Core/interface/DQMStore.h"
00003 #include "DQM/SiPixelCommon/interface/SiPixelHistogramId.h"
00005 #include "FWCore/ServiceRegistry/interface/Service.h"
00006 // STL
00007 #include <vector>
00008 #include <memory>
00009 #include <string>
00010 #include <iostream>
00011 #include <stdlib.h>
00012 
00013 // Data Formats
00014 #include "DataFormats/SiPixelDetId/interface/PixelBarrelName.h"
00015 #include "DataFormats/SiPixelDetId/interface/PixelEndcapName.h"
00016 #include "DataFormats/DetId/interface/DetId.h"
00017 #include "DataFormats/SiPixelDetId/interface/PixelSubdetector.h"
00018 //
00019 // Constructors
00020 //
00021 SiPixelRecHitModule::SiPixelRecHitModule() : id_(0) { }
00023 SiPixelRecHitModule::SiPixelRecHitModule(const uint32_t& id) : 
00024   id_(id)
00025 { 
00026 }
00027 
00028 //
00029 // Destructor
00030 //
00031 SiPixelRecHitModule::~SiPixelRecHitModule() {}
00032 //
00033 // Book histograms
00034 //
00035 void SiPixelRecHitModule::book(const edm::ParameterSet& iConfig, int type, 
00036                                bool twoD, bool reducedSet) {
00037 
00038   bool barrel = DetId(id_).subdetId() == static_cast<int>(PixelSubdetector::PixelBarrel);
00039   bool endcap = DetId(id_).subdetId() == static_cast<int>(PixelSubdetector::PixelEndcap);
00040   bool isHalfModule = false;
00041   if(barrel){
00042     isHalfModule = PixelBarrelName(DetId(id_)).isHalfModule(); 
00043   }
00044 
00045   std::string hid;
00046   // Get collection name and instantiate Histo Id builder
00047   edm::InputTag src = iConfig.getParameter<edm::InputTag>( "src" );
00048   // Get DQM interface
00049   DQMStore* theDMBE = edm::Service<DQMStore>().operator->();
00050 
00051 
00052   if(type==0){
00053     SiPixelHistogramId* theHistogramId = new SiPixelHistogramId( src.label() );
00054         if(!reducedSet)
00055         {
00056     if(twoD){
00057       // XYPosition
00058       hid = theHistogramId->setHistoId("xypos",id_);
00059       meXYPos_ = theDMBE->book2D(hid,"XY Position",100,-1.,1,100,-4,4);
00060       meXYPos_->setAxisTitle("X Position",1);
00061       meXYPos_->setAxisTitle("Y Position",2);
00062     }
00063     else{
00064       // projections of XYPosition
00065       hid = theHistogramId->setHistoId("xypos",id_);
00066       meXYPos_px_ = theDMBE->book1D(hid+"_px","X Position",100,-1.,1);
00067       meXYPos_px_->setAxisTitle("X Position",1);
00068       meXYPos_py_ = theDMBE->book1D(hid+"_py","Y Position",100,-4,4);
00069       meXYPos_py_->setAxisTitle("Y Position",1);
00070     }
00071         }
00072     hid = theHistogramId->setHistoId("ClustX",id_);
00073     meClustX_ = theDMBE->book1D(hid, "RecHit X size", 10, 0., 10.);
00074     meClustX_->setAxisTitle("RecHit size X dimension", 1);
00075     hid = theHistogramId->setHistoId("ClustY",id_);
00076     meClustY_ = theDMBE->book1D(hid, "RecHit Y size", 15, 0., 15.);
00077     meClustY_->setAxisTitle("RecHit size Y dimension", 1); 
00078 
00079     hid = theHistogramId->setHistoId("ErrorX",id_);
00080     meErrorX_ = theDMBE->book1D(hid, "RecHit error X", 100,0.,0.02);
00081     meErrorX_->setAxisTitle("RecHit error X", 1);
00082     hid = theHistogramId->setHistoId("ErrorY",id_);
00083     meErrorY_ = theDMBE->book1D(hid, "RecHit error Y", 100,0.,0.02);
00084     meErrorY_->setAxisTitle("RecHit error Y", 1);
00085 
00086     //Removed to save offline memory
00087     //hid = theHistogramId->setHistoId("nRecHits",id_);
00088     //menRecHits_ = theDMBE->book1D(hid, "# of rechits in this module", 8, 0, 8);
00089     //menRecHits_->setAxisTitle("number of rechits",1);  
00090     delete theHistogramId;
00091   }
00092 
00093   if(type==1 && barrel){
00094     uint32_t DBladder = PixelBarrelName(DetId(id_)).ladderName();
00095     char sladder[80]; sprintf(sladder,"Ladder_%02i",DBladder);
00096     hid = src.label() + "_" + sladder;
00097     if(isHalfModule) hid += "H";
00098     else hid += "F";
00099         if(!reducedSet)
00100         {
00101     if(twoD){
00102       meXYPosLad_ = theDMBE->book2D("xypos_" + hid,"XY Position",100,-1.,1,100,-4,4);
00103       meXYPosLad_->setAxisTitle("X Position",1);
00104       meXYPosLad_->setAxisTitle("Y Position",2);
00105     }
00106     else{
00107       // projections of XYPosition
00108       meXYPosLad_px_ = theDMBE->book1D("xypos_"+hid+"_px","X Position",100,-1.,1);
00109       meXYPosLad_px_->setAxisTitle("X Position",1);
00110       meXYPosLad_py_ = theDMBE->book1D("xypos_"+hid+"_py","Y Position",100,-4,4);
00111       meXYPosLad_py_->setAxisTitle("Y Position",1);
00112     }
00113         }
00114     meClustXLad_ = theDMBE->book1D("ClustX_" +hid, "RecHit X size", 10, 0., 10.);
00115     meClustXLad_->setAxisTitle("RecHit size X dimension", 1);
00116     meClustYLad_ = theDMBE->book1D("ClustY_" +hid,"RecHit Y size", 15, 0.,15.);
00117     meClustYLad_->setAxisTitle("RecHit size Y dimension", 1);
00118     meErrorXLad_ = theDMBE->book1D("ErrorX_"+hid, "RecHit error X", 100,0.,0.02);
00119     meErrorXLad_->setAxisTitle("RecHit error X", 1);
00120     meErrorYLad_ = theDMBE->book1D("ErrorY_"+hid, "RecHit error Y", 100,0.,0.02);
00121     meErrorYLad_->setAxisTitle("RecHit error Y", 1);
00122     menRecHitsLad_ = theDMBE->book1D("nRecHits_"+hid, "# of rechits in this module", 8, 0, 8);
00123     menRecHitsLad_->setAxisTitle("number of rechits",1);
00124 
00125   }
00126 
00127   if(type==2 && barrel){
00128     
00129     uint32_t DBlayer = PixelBarrelName(DetId(id_)).layerName();
00130     char slayer[80]; sprintf(slayer,"Layer_%i",DBlayer);
00131     hid = src.label() + "_" + slayer;
00132     
00133         if(!reducedSet)
00134         {
00135     if(twoD){
00136       meXYPosLay_ = theDMBE->book2D("xypos_" + hid,"XY Position",100,-1.,1,100,-4,4);
00137       meXYPosLay_->setAxisTitle("X Position",1);
00138       meXYPosLay_->setAxisTitle("Y Position",2);
00139     }
00140     else{
00141       // projections of XYPosition
00142       meXYPosLay_px_ = theDMBE->book1D("xypos_"+hid+"_px","X Position",100,-1.,1);
00143       meXYPosLay_px_->setAxisTitle("X Position",1);
00144       meXYPosLay_py_ = theDMBE->book1D("xypos_"+hid+"_py","Y Position",100,-4,4);
00145       meXYPosLay_py_->setAxisTitle("Y Position",1);
00146     }
00147         }
00148 
00149     meClustXLay_ = theDMBE->book1D("ClustX_" +hid, "RecHit X size", 10, 0., 10.);
00150     meClustXLay_->setAxisTitle("RecHit size X dimension", 1);
00151     meClustYLay_ = theDMBE->book1D("ClustY_" +hid,"RecHit Y size", 15, 0.,15.);
00152     meClustYLay_->setAxisTitle("RecHit size Y dimension", 1);
00153     meErrorXLay_ = theDMBE->book1D("ErrorX_"+hid, "RecHit error X", 100,0.,0.02);
00154     meErrorXLay_->setAxisTitle("RecHit error X", 1);
00155     meErrorYLay_ = theDMBE->book1D("ErrorY_"+hid, "RecHit error Y", 100,0.,0.02);
00156     meErrorYLay_->setAxisTitle("RecHit error Y", 1);
00157     menRecHitsLay_ = theDMBE->book1D("nRecHits_"+hid, "# of rechits in this module", 8, 0, 8);
00158     menRecHitsLay_->setAxisTitle("number of rechits",1);
00159 
00160   }
00161 
00162   if(type==3 && barrel){
00163     uint32_t DBmodule = PixelBarrelName(DetId(id_)).moduleName();
00164     char smodule[80]; sprintf(smodule,"Ring_%i",DBmodule);
00165     hid = src.label() + "_" + smodule;
00166     
00167         if(!reducedSet)
00168         {
00169     if(twoD){
00170       meXYPosPhi_ = theDMBE->book2D("xypos_" + hid,"XY Position",100,-1.,1,100,-4,4);
00171       meXYPosPhi_->setAxisTitle("X Position",1);
00172       meXYPosPhi_->setAxisTitle("Y Position",2);
00173     }
00174     else{
00175       // projections of XYPosition
00176       meXYPosPhi_px_ = theDMBE->book1D("xypos_"+hid+"_px","X Position",100,-1.,1);
00177       meXYPosPhi_px_->setAxisTitle("X Position",1);
00178       meXYPosPhi_py_ = theDMBE->book1D("xypos_"+hid+"_py","Y Position",100,-4,4);
00179       meXYPosPhi_py_->setAxisTitle("Y Position",1);
00180     }
00181         }
00182     meClustXPhi_ = theDMBE->book1D("ClustX_" +hid, "RecHit X size", 10, 0., 10.);
00183     meClustXPhi_->setAxisTitle("RecHit size X dimension", 1);
00184     meClustYPhi_ = theDMBE->book1D("ClustY_" +hid,"RecHit Y size", 15, 0.,15.);
00185     meClustYPhi_->setAxisTitle("RecHit size Y dimension", 1);
00186     meErrorXPhi_ = theDMBE->book1D("ErrorX_"+hid, "RecHit error X", 100,0.,0.02);
00187     meErrorXPhi_->setAxisTitle("RecHit error X", 1);
00188     meErrorYPhi_ = theDMBE->book1D("ErrorY_"+hid, "RecHit error Y", 100,0.,0.02);
00189     meErrorYPhi_->setAxisTitle("RecHit error Y", 1);
00190     menRecHitsPhi_ = theDMBE->book1D("nRecHits_"+hid, "# of rechits in this module", 8, 0, 8);
00191     menRecHitsPhi_->setAxisTitle("number of rechits",1);
00192 
00193   }
00194 
00195   if(type==4 && endcap){
00196     uint32_t blade= PixelEndcapName(DetId(id_)).bladeName();
00197     
00198     char sblade[80]; sprintf(sblade, "Blade_%02i",blade);
00199     hid = src.label() + "_" + sblade;
00200 //     meXYPosBlade_ = theDMBE->book2D("xypos_" + hid,"XY Position",100,-1.,1,100,-4,4);
00201 //     meXYPosBlade_->setAxisTitle("X Position",1);
00202 //     meXYPosBlade_->setAxisTitle("Y Position",2);
00203 
00204     meClustXBlade_ = theDMBE->book1D("ClustX_" +hid, "RecHit X size", 10, 0., 10.);
00205     meClustXBlade_->setAxisTitle("RecHit size X dimension", 1);
00206     meClustYBlade_ = theDMBE->book1D("ClustY_" +hid,"RecHit Y size", 15, 0.,15.);
00207     meClustYBlade_->setAxisTitle("RecHit size Y dimension", 1);
00208     meErrorXBlade_ = theDMBE->book1D("ErrorX_"+hid, "RecHit error X", 100,0.,0.02);
00209     meErrorXBlade_->setAxisTitle("RecHit error X", 1);
00210     meErrorYBlade_ = theDMBE->book1D("ErrorY_"+hid, "RecHit error Y", 100,0.,0.02);
00211     meErrorYBlade_->setAxisTitle("RecHit error Y", 1);
00212     menRecHitsBlade_ = theDMBE->book1D("nRecHits_"+hid, "# of rechits in this module", 8, 0, 8);
00213     menRecHitsBlade_->setAxisTitle("number of rechits",1);
00214 
00215   }
00216   if(type==5 && endcap){
00217     uint32_t disk = PixelEndcapName(DetId(id_)).diskName();
00218     
00219     char sdisk[80]; sprintf(sdisk, "Disk_%i",disk);
00220     hid = src.label() + "_" + sdisk;
00221 //     meXYPosDisk_ = theDMBE->book2D("xypos_" + hid,"XY Position",100,-1.,1,100,-4,4);
00222 //     meXYPosDisk_->setAxisTitle("X Position",1);
00223 //     meXYPosDisk_->setAxisTitle("Y Position",2);
00224 
00225     meClustXDisk_ = theDMBE->book1D("ClustX_" +hid, "RecHit X size", 10, 0., 10.);
00226     meClustXDisk_->setAxisTitle("RecHit size X dimension", 1);
00227     meClustYDisk_ = theDMBE->book1D("ClustY_" +hid,"RecHit Y size", 15, 0.,15.);
00228     meClustYDisk_->setAxisTitle("RecHit size Y dimension", 1);
00229     meErrorXDisk_ = theDMBE->book1D("ErrorX_"+hid, "RecHit error X", 100,0.,0.02);
00230     meErrorXDisk_->setAxisTitle("RecHit error X", 1);
00231     meErrorYDisk_ = theDMBE->book1D("ErrorY_"+hid, "RecHit error Y", 100,0.,0.02);
00232     meErrorYDisk_->setAxisTitle("RecHit error Y", 1);
00233     menRecHitsDisk_ = theDMBE->book1D("nRecHits_"+hid, "# of rechits in this module", 8, 0, 8);
00234     menRecHitsDisk_->setAxisTitle("number of rechits",1);
00235 
00236   }
00237 
00238   if(type==6 && endcap){
00239     uint32_t panel= PixelEndcapName(DetId(id_)).pannelName();
00240     uint32_t module= PixelEndcapName(DetId(id_)).plaquetteName();
00241     char slab[80]; sprintf(slab, "Panel_%i_Ring_%i",panel, module);
00242     hid = src.label() + "_" + slab;
00243     
00244         if(!reducedSet)
00245         {
00246     if(twoD){
00247       meXYPosRing_ = theDMBE->book2D("xypos_" + hid,"XY Position",100,-1.,1,100,-4,4);
00248       meXYPosRing_->setAxisTitle("X Position",1);
00249       meXYPosRing_->setAxisTitle("Y Position",2);
00250     }
00251     else{
00252       // projections of XYPosition
00253       meXYPosRing_px_ = theDMBE->book1D("xypos_"+hid+"_px","X Position",100,-1.,1);
00254       meXYPosRing_px_->setAxisTitle("X Position",1);
00255       meXYPosRing_py_ = theDMBE->book1D("xypos_"+hid+"_py","Y Position",100,-4,4);
00256       meXYPosRing_py_->setAxisTitle("Y Position",1);
00257     }
00258         }
00259     meClustXRing_ = theDMBE->book1D("ClustX_" +hid, "RecHit X size", 10, 0., 10.);
00260     meClustXRing_->setAxisTitle("RecHit size X dimension", 1);
00261     meClustYRing_ = theDMBE->book1D("ClustY_" +hid,"RecHit Y size", 15, 0.,15.);
00262     meClustYRing_->setAxisTitle("RecHit size Y dimension", 1);
00263     meErrorXRing_ = theDMBE->book1D("ErrorX_"+hid, "RecHit error X", 100,0.,0.02);
00264     meErrorXRing_->setAxisTitle("RecHit error X", 1);
00265     meErrorYRing_ = theDMBE->book1D("ErrorY_"+hid, "RecHit error Y", 100,0.,0.02);
00266     meErrorYRing_->setAxisTitle("RecHit error Y", 1);
00267     menRecHitsRing_ = theDMBE->book1D("nRecHits_"+hid, "# of rechits in this module", 8, 0, 8);
00268     menRecHitsRing_->setAxisTitle("number of rechits",1);
00269 
00270   }
00271 
00272 }
00273 //
00274 // Fill histograms
00275 //
00276 void SiPixelRecHitModule::fill(const float& rechit_x, const float& rechit_y, 
00277                                const int& sizeX, const int& sizeY, 
00278                                const float& lerr_x, const float& lerr_y, 
00279                                bool modon, bool ladon, bool layon, bool phion, 
00280                                bool bladeon, bool diskon, bool ringon, 
00281                                bool twoD, bool reducedSet) {
00282 
00283   bool barrel = DetId(id_).subdetId() == static_cast<int>(PixelSubdetector::PixelBarrel);
00284   bool endcap = DetId(id_).subdetId() == static_cast<int>(PixelSubdetector::PixelEndcap);
00285 
00286   //std::cout << rechit_x << " " << rechit_y << " " << sizeX << " " << sizeY << std::endl;
00287   if(modon){
00288     meClustX_->Fill(sizeX);
00289     meClustY_->Fill(sizeY);
00290     meErrorX_->Fill(lerr_x);
00291     meErrorY_->Fill(lerr_y);  
00292         if(!reducedSet)
00293         {
00294     if(twoD) meXYPos_->Fill(rechit_x, rechit_y);
00295     else {
00296       meXYPos_px_->Fill(rechit_x); 
00297       meXYPos_py_->Fill(rechit_y);
00298     }
00299         }
00300   }
00301   //std::cout<<"number of detector units="<<numberOfDetUnits<<std::endl;
00302 
00303   if(ladon && barrel){
00304     meClustXLad_->Fill(sizeX);
00305     meClustYLad_->Fill(sizeY);
00306     meErrorXLad_->Fill(lerr_x);
00307     meErrorYLad_->Fill(lerr_y);  
00308         if(!reducedSet)
00309         {
00310     if(twoD) meXYPosLad_->Fill(rechit_x, rechit_y);
00311     else{
00312       meXYPosLad_px_->Fill(rechit_x); 
00313       meXYPosLad_py_->Fill(rechit_y);
00314     }
00315         }
00316   }
00317 
00318   if(layon && barrel){
00319     meClustXLay_->Fill(sizeX);
00320     meClustYLay_->Fill(sizeY);
00321     meErrorXLay_->Fill(lerr_x);
00322     meErrorYLay_->Fill(lerr_y); 
00323         if(!reducedSet)
00324         {
00325     if(twoD) meXYPosLay_->Fill(rechit_x, rechit_y);
00326     else{
00327       meXYPosLay_px_->Fill(rechit_x); 
00328       meXYPosLay_py_->Fill(rechit_y);
00329     }
00330         }
00331   }
00332 
00333   if(phion && barrel){
00334     meClustXPhi_->Fill(sizeX);
00335     meClustYPhi_->Fill(sizeY);
00336     meErrorXPhi_->Fill(lerr_x);
00337     meErrorYPhi_->Fill(lerr_y); 
00338     if(!reducedSet)
00339         {
00340     if(twoD) meXYPosPhi_->Fill(rechit_x, rechit_y);
00341     else{
00342       meXYPosPhi_px_->Fill(rechit_x); 
00343       meXYPosPhi_py_->Fill(rechit_y);
00344     }
00345         }       
00346   }
00347 
00348   if(bladeon && endcap){
00349     //meXYPosBlade_->Fill(rechit_x, rechit_y);
00350     meClustXBlade_->Fill(sizeX);
00351     meClustYBlade_->Fill(sizeY);
00352     meErrorXBlade_->Fill(lerr_x);
00353     meErrorYBlade_->Fill(lerr_y); 
00354   }
00355 
00356   if(diskon && endcap){
00357     //meXYPosDisk_->Fill(rechit_x, rechit_y);
00358     meClustXDisk_->Fill(sizeX);
00359     meClustYDisk_->Fill(sizeY);
00360     meErrorXDisk_->Fill(lerr_x);
00361     meErrorYDisk_->Fill(lerr_y); 
00362   }
00363 
00364   if(ringon && endcap){
00365     meClustXRing_->Fill(sizeX);
00366     meClustYRing_->Fill(sizeY);
00367     meErrorXRing_->Fill(lerr_x);
00368     meErrorYRing_->Fill(lerr_y); 
00369         if(!reducedSet)
00370         {
00371     if(twoD) meXYPosRing_->Fill(rechit_x, rechit_y);
00372     else{
00373       meXYPosRing_px_->Fill(rechit_x); 
00374       meXYPosRing_py_->Fill(rechit_y);
00375     }
00376         }       
00377   }
00378 }
00379 
00380 void SiPixelRecHitModule::nfill(const int& nrec, bool modon, bool ladon, bool layon, bool phion, bool bladeon, bool diskon, bool ringon) {
00381   
00382   bool barrel = DetId(id_).subdetId() == static_cast<int>(PixelSubdetector::PixelBarrel);
00383   bool endcap = DetId(id_).subdetId() == static_cast<int>(PixelSubdetector::PixelEndcap);
00384 
00385   //if(modon) menRecHits_->Fill(nrec);
00386   //barrel
00387   if(ladon && barrel) menRecHitsLad_->Fill(nrec);
00388   if(layon && barrel) menRecHitsLay_->Fill(nrec);
00389   if(phion && barrel) menRecHitsPhi_->Fill(nrec);
00390   //endcap
00391   if(bladeon && endcap) menRecHitsBlade_->Fill(nrec);
00392   if(diskon && endcap) menRecHitsDisk_->Fill(nrec);
00393   if(ringon && endcap) menRecHitsRing_->Fill(nrec);
00394 }