CMS 3D CMS Logo

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