Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:08:24

0001 // Package:    SiPixelMonitorTrack
0002 // Class:      SiPixelTrackResidualModule
0003 //
0004 // class SiPixelTrackResidualModule SiPixelTrackResidualModule.cc
0005 //       DQM/SiPixelMonitorTrack/src/SiPixelTrackResidualModule.cc
0006 //
0007 // Description: SiPixel hit-to-track residual data quality monitoring modules
0008 // Implementation: prototype -> improved -> never final - end of the 1st step
0009 //
0010 // Original Author: Shan-Huei Chuang
0011 //         Created: Fri Mar 23 18:41:42 CET 2007
0012 
0013 #include <iostream>
0014 #include <string>
0015 
0016 #include "DQM/SiPixelCommon/interface/SiPixelHistogramId.h"
0017 #include "DQM/SiPixelMonitorTrack/interface/SiPixelTrackResidualModule.h"
0018 
0019 // Data Formats
0020 #include "DataFormats/DetId/interface/DetId.h"
0021 #include "DataFormats/TrackerCommon/interface/PixelBarrelName.h"
0022 #include "DataFormats/SiPixelDetId/interface/PixelBarrelNameUpgrade.h"
0023 #include "DataFormats/TrackerCommon/interface/PixelEndcapName.h"
0024 #include "DataFormats/SiPixelDetId/interface/PixelEndcapNameUpgrade.h"
0025 #include "DataFormats/SiPixelDetId/interface/PixelSubdetector.h"
0026 
0027 using namespace std;
0028 using namespace edm;
0029 
0030 SiPixelTrackResidualModule::SiPixelTrackResidualModule() : id_(0) { bBookTracks = true; }
0031 
0032 SiPixelTrackResidualModule::SiPixelTrackResidualModule(uint32_t id) : id_(id) { bBookTracks = true; }
0033 
0034 SiPixelTrackResidualModule::~SiPixelTrackResidualModule() {}
0035 
0036 void SiPixelTrackResidualModule::book(const edm::ParameterSet &iConfig,
0037                                       const TrackerTopology *pTT,
0038                                       DQMStore::IBooker &iBooker,
0039                                       bool reducedSet,
0040                                       int type,
0041                                       bool isUpgrade) {
0042   bool barrel = DetId(id_).subdetId() == static_cast<int>(PixelSubdetector::PixelBarrel);
0043   bool endcap = DetId(id_).subdetId() == static_cast<int>(PixelSubdetector::PixelEndcap);
0044   bool isHalfModule = false;
0045   if (barrel) {
0046     isHalfModule = PixelBarrelName(DetId(id_), pTT, isUpgrade).isHalfModule();
0047   }
0048 
0049   edm::InputTag src = iConfig.getParameter<edm::InputTag>("src");
0050   std::string hisID;
0051 
0052   if (type == 0) {
0053     SiPixelHistogramId *theHistogramId = new SiPixelHistogramId(src.label());
0054     hisID = theHistogramId->setHistoId("residualX", id_);
0055     meResidualX_ = iBooker.book1D(hisID, "Hit-to-Track Residual in r-phi", 100, -150, 150);
0056     meResidualX_->setAxisTitle("hit-to-track residual in r-phi (um)", 1);
0057     hisID = theHistogramId->setHistoId("residualY", id_);
0058     meResidualY_ = iBooker.book1D(hisID, "Hit-to-Track Residual in Z", 100, -300, 300);
0059     meResidualY_->setAxisTitle("hit-to-track residual in z (um)", 1);
0060     // Number of clusters
0061     hisID = theHistogramId->setHistoId("nclusters_OnTrack", id_);
0062     meNClusters_onTrack_ = iBooker.book1D(hisID, "Number of Clusters (on Track)", 10, 0., 10.);
0063     meNClusters_onTrack_->setAxisTitle("Number of Clusters on Track", 1);
0064     // Total cluster charge in ke
0065     hisID = theHistogramId->setHistoId("charge_OnTrack", id_);
0066     meCharge_onTrack_ = iBooker.book1D(hisID, "Normalized Cluster charge (on Track)", 100, 0., 200.);
0067     meCharge_onTrack_->setAxisTitle("Charge [kilo electrons]", 1);
0068     // Total cluster size (in pixels)
0069     hisID = theHistogramId->setHistoId("size_OnTrack", id_);
0070     meSize_onTrack_ = iBooker.book1D(hisID, "Total cluster size (on Track)", 30, 0., 30.);
0071     meSize_onTrack_->setAxisTitle("Cluster size [number of pixels]", 1);
0072     // Number of clusters
0073     hisID = theHistogramId->setHistoId("nclusters_OffTrack", id_);
0074     meNClusters_offTrack_ = iBooker.book1D(hisID, "Number of Clusters (off Track)", 35, 0., 35.);
0075     meNClusters_offTrack_->setAxisTitle("Number of Clusters off Track", 1);
0076     // Total cluster charge in ke
0077     hisID = theHistogramId->setHistoId("charge_OffTrack", id_);
0078     meCharge_offTrack_ = iBooker.book1D(hisID, "Cluster charge (off Track)", 100, 0., 200.);
0079     meCharge_offTrack_->setAxisTitle("Charge [kilo electrons]", 1);
0080     // Total cluster size (in pixels)
0081     hisID = theHistogramId->setHistoId("size_OffTrack", id_);
0082     meSize_offTrack_ = iBooker.book1D(hisID, "Total cluster size (off Track)", 30, 0., 30.);
0083     meSize_offTrack_->setAxisTitle("Cluster size [number of pixels]", 1);
0084     if (!reducedSet) {
0085       // Cluster width on the x-axis
0086       hisID = theHistogramId->setHistoId("sizeX_OnTrack", id_);
0087       meSizeX_onTrack_ = iBooker.book1D(hisID, "Cluster x-width (rows) (on Track)", 10, 0., 10.);
0088       meSizeX_onTrack_->setAxisTitle("Cluster x-size [rows]", 1);
0089       // Cluster width on the y-axis
0090       hisID = theHistogramId->setHistoId("sizeY_OnTrack", id_);
0091       meSizeY_onTrack_ = iBooker.book1D(hisID, "Cluster y-width (columns) (on Track)", 15, 0., 15.);
0092       meSizeY_onTrack_->setAxisTitle("Cluster y-size [columns]", 1);
0093       // Cluster width on the x-axis
0094       hisID = theHistogramId->setHistoId("sizeX_OffTrack", id_);
0095       meSizeX_offTrack_ = iBooker.book1D(hisID, "Cluster x-width (rows) (off Track)", 10, 0., 10.);
0096       meSizeX_offTrack_->setAxisTitle("Cluster x-size [rows]", 1);
0097       // Cluster width on the y-axis
0098       hisID = theHistogramId->setHistoId("sizeY_OffTrack", id_);
0099       meSizeY_offTrack_ = iBooker.book1D(hisID, "Cluster y-width (columns) (off Track)", 15, 0., 15.);
0100       meSizeY_offTrack_->setAxisTitle("Cluster y-size [columns]", 1);
0101     }
0102     delete theHistogramId;
0103   }
0104 
0105   if (type == 1 && barrel) {
0106     uint32_t DBladder;
0107     DBladder = PixelBarrelName(DetId(id_), pTT, isUpgrade).ladderName();
0108     char sladder[80];
0109     sprintf(sladder, "Ladder_%02i", DBladder);
0110     hisID = src.label() + "_" + sladder;
0111     if (isHalfModule)
0112       hisID += "H";
0113     else
0114       hisID += "F";
0115     meResidualXLad_ = iBooker.book1D("residualX_" + hisID, "Hit-to-Track Residual in r-phi", 100, -150, 150);
0116     meResidualXLad_->setAxisTitle("hit-to-track residual in r-phi (um)", 1);
0117     meResidualYLad_ = iBooker.book1D("residualY_" + hisID, "Hit-to-Track Residual in Z", 100, -300, 300);
0118     meResidualYLad_->setAxisTitle("hit-to-track residual in z (um)", 1);
0119     // Number of clusters
0120     meNClusters_onTrackLad_ =
0121         iBooker.book1D("nclusters_OnTrack_" + hisID, "Number of Clusters (on Track)", 10, 0., 10.);
0122     meNClusters_onTrackLad_->setAxisTitle("Number of Clusters on Track", 1);
0123     // Total cluster charge in MeV
0124     meCharge_onTrackLad_ =
0125         iBooker.book1D("charge_OnTrack_" + hisID, "Normalized Cluster charge (on Track)", 100, 0., 200.);
0126     meCharge_onTrackLad_->setAxisTitle("Charge [kilo electrons]", 1);
0127     // Total cluster size (in pixels)
0128     meSize_onTrackLad_ = iBooker.book1D("size_OnTrack_" + hisID, "Total cluster size (on Track)", 30, 0., 30.);
0129     meSize_onTrackLad_->setAxisTitle("Cluster size [number of pixels]", 1);
0130     // Number of clusters
0131     meNClusters_offTrackLad_ =
0132         iBooker.book1D("nclusters_OffTrack_" + hisID, "Number of Clusters (off Track)", 35, 0., 35.);
0133     meNClusters_offTrackLad_->setAxisTitle("Number of Clusters off Track", 1);
0134     // Total cluster charge in MeV
0135     meCharge_offTrackLad_ = iBooker.book1D("charge_OffTrack_" + hisID, "Cluster charge (off Track)", 100, 0., 200.);
0136     meCharge_offTrackLad_->setAxisTitle("Charge [kilo electrons]", 1);
0137     // Total cluster size (in pixels)
0138     meSize_offTrackLad_ = iBooker.book1D("size_OffTrack_" + hisID, "Total cluster size (off Track)", 30, 0., 30.);
0139     meSize_offTrackLad_->setAxisTitle("Cluster size [number of pixels]", 1);
0140     if (!reducedSet) {
0141       // Cluster width on the x-axis
0142       meSizeX_offTrackLad_ =
0143           iBooker.book1D("sizeX_OffTrack_" + hisID, "Cluster x-width (rows) (off Track)", 10, 0., 10.);
0144       meSizeX_offTrackLad_->setAxisTitle("Cluster x-size [rows]", 1);
0145       // Cluster width on the y-axis
0146       meSizeY_offTrackLad_ =
0147           iBooker.book1D("sizeY_OffTrack_" + hisID, "Cluster y-width (columns) (off Track)", 15, 0., 15.);
0148       meSizeY_offTrackLad_->setAxisTitle("Cluster y-size [columns]", 1);
0149       // Cluster width on the x-axis
0150       meSizeX_onTrackLad_ = iBooker.book1D("sizeX_OnTrack_" + hisID, "Cluster x-width (rows) (on Track)", 10, 0., 10.);
0151       meSizeX_onTrackLad_->setAxisTitle("Cluster x-size [rows]", 1);
0152       // Cluster width on the y-axis
0153       meSizeY_onTrackLad_ =
0154           iBooker.book1D("sizeY_OnTrack_" + hisID, "Cluster y-width (columns) (on Track)", 15, 0., 15.);
0155       meSizeY_onTrackLad_->setAxisTitle("Cluster y-size [columns]", 1);
0156     }
0157   }
0158 
0159   if (type == 2 && barrel) {
0160     uint32_t DBlayer;
0161     DBlayer = PixelBarrelName(DetId(id_), pTT, isUpgrade).layerName();
0162     char slayer[80];
0163     sprintf(slayer, "Layer_%i", DBlayer);
0164     hisID = src.label() + "_" + slayer;
0165     meResidualXLay_ = iBooker.book1D("residualX_" + hisID, "Hit-to-Track Residual in r-phi", 100, -150, 150);
0166     meResidualXLay_->setAxisTitle("hit-to-track residual in r-phi (um)", 1);
0167     meResidualYLay_ = iBooker.book1D("residualY_" + hisID, "Hit-to-Track Residual in Z", 100, -300, 300);
0168     meResidualYLay_->setAxisTitle("hit-to-track residual in z (um)", 1);
0169     // Number of clusters
0170     meNClusters_onTrackLay_ =
0171         iBooker.book1D("nclusters_OnTrack_" + hisID, "Number of Clusters (on Track)", 10, 0., 10.);
0172     meNClusters_onTrackLay_->setAxisTitle("Number of Clusters on Track", 1);
0173     // Total cluster charge in MeV
0174     meCharge_onTrackLay_ =
0175         iBooker.book1D("charge_OnTrack_" + hisID, "Normalized Cluster charge (on Track)", 100, 0., 200.);
0176     meCharge_onTrackLay_->setAxisTitle("Charge [kilo electrons]", 1);
0177     // Total cluster size (in pixels)
0178     meSize_onTrackLay_ = iBooker.book1D("size_OnTrack_" + hisID, "Total cluster size (on Track)", 30, 0., 30.);
0179     meSize_onTrackLay_->setAxisTitle("Cluster size [number of pixels]", 1);
0180     // Number of clusters
0181     meNClusters_offTrackLay_ =
0182         iBooker.book1D("nclusters_OffTrack_" + hisID, "Number of Clusters (off Track)", 35, 0., 35.);
0183     meNClusters_offTrackLay_->setAxisTitle("Number of Clusters off Track", 1);
0184     // Total cluster charge in MeV
0185     meCharge_offTrackLay_ = iBooker.book1D("charge_OffTrack_" + hisID, "Cluster charge (off Track)", 100, 0., 200.);
0186     meCharge_offTrackLay_->setAxisTitle("Charge [kilo electrons]", 1);
0187     // Total cluster size (in pixels)
0188     meSize_offTrackLay_ = iBooker.book1D("size_OffTrack_" + hisID, "Total cluster size (off Track)", 30, 0., 30.);
0189     meSize_offTrackLay_->setAxisTitle("Cluster size [number of pixels]", 1);
0190     if (!reducedSet) {
0191       // Cluster width on the x-axis
0192       meSizeX_onTrackLay_ = iBooker.book1D("sizeX_OnTrack_" + hisID, "Cluster x-width (rows) (on Track)", 10, 0., 10.);
0193       meSizeX_onTrackLay_->setAxisTitle("Cluster x-size [rows]", 1);
0194       // Cluster width on the y-axis
0195       meSizeY_onTrackLay_ =
0196           iBooker.book1D("sizeY_OnTrack_" + hisID, "Cluster y-width (columns) (on Track)", 15, 0., 15.);
0197       meSizeY_onTrackLay_->setAxisTitle("Cluster y-size [columns]", 1);
0198       // Cluster width on the x-axis
0199       meSizeX_offTrackLay_ =
0200           iBooker.book1D("sizeX_OffTrack_" + hisID, "Cluster x-width (rows) (off Track)", 10, 0., 10.);
0201       meSizeX_offTrackLay_->setAxisTitle("Cluster x-size [rows]", 1);
0202       // Cluster width on the y-axis
0203       meSizeY_offTrackLay_ =
0204           iBooker.book1D("sizeY_OffTrack_" + hisID, "Cluster y-width (columns) (off Track)", 15, 0., 15.);
0205       meSizeY_offTrackLay_->setAxisTitle("Cluster y-size [columns]", 1);
0206     }
0207   }
0208 
0209   if (type == 3 && barrel) {
0210     uint32_t DBmodule;
0211     DBmodule = PixelBarrelName(DetId(id_), pTT, isUpgrade).moduleName();
0212     char smodule[80];
0213     sprintf(smodule, "Ring_%i", DBmodule);
0214     hisID = src.label() + "_" + smodule;
0215     meResidualXPhi_ = iBooker.book1D("residualX_" + hisID, "Hit-to-Track Residual in r-phi", 100, -150, 150);
0216     meResidualXPhi_->setAxisTitle("hit-to-track residual in r-phi (um)", 1);
0217     meResidualYPhi_ = iBooker.book1D("residualY_" + hisID, "Hit-to-Track Residual in Z", 100, -300, 300);
0218     meResidualYPhi_->setAxisTitle("hit-to-track residual in z (um)", 1);
0219     // Number of clusters
0220     meNClusters_onTrackPhi_ =
0221         iBooker.book1D("nclusters_OnTrack_" + hisID, "Number of Clusters (on Track)", 10, 0., 10.);
0222     meNClusters_onTrackPhi_->setAxisTitle("Number of Clusters on Track", 1);
0223     // Total cluster charge in MeV
0224     meCharge_onTrackPhi_ =
0225         iBooker.book1D("charge_OnTrack_" + hisID, "Normalized Cluster charge (on Track)", 100, 0., 200.);
0226     meCharge_onTrackPhi_->setAxisTitle("Charge [kilo electrons]", 1);
0227     // Total cluster size (in pixels)
0228     meSize_onTrackPhi_ = iBooker.book1D("size_OnTrack_" + hisID, "Total cluster size (on Track)", 30, 0., 30.);
0229     meSize_onTrackPhi_->setAxisTitle("Cluster size [number of pixels]", 1);
0230     // Number of clusters
0231     meNClusters_offTrackPhi_ =
0232         iBooker.book1D("nclusters_OffTrack_" + hisID, "Number of Clusters (off Track)", 35, 0., 35.);
0233     meNClusters_offTrackPhi_->setAxisTitle("Number of Clusters off Track", 1);
0234     // Total cluster charge in MeV
0235     meCharge_offTrackPhi_ = iBooker.book1D("charge_OffTrack_" + hisID, "Cluster charge (off Track)", 100, 0., 200.);
0236     meCharge_offTrackPhi_->setAxisTitle("Charge [kilo electrons]", 1);
0237     // Total cluster size (in pixels)
0238     meSize_offTrackPhi_ = iBooker.book1D("size_OffTrack_" + hisID, "Total cluster size (off Track)", 30, 0., 30.);
0239     meSize_offTrackPhi_->setAxisTitle("Cluster size [number of pixels]", 1);
0240     if (!reducedSet) {
0241       // Cluster width on the x-axis
0242       meSizeX_onTrackPhi_ = iBooker.book1D("sizeX_OnTrack_" + hisID, "Cluster x-width (rows) (on Track)", 10, 0., 10.);
0243       meSizeX_onTrackPhi_->setAxisTitle("Cluster x-size [rows]", 1);
0244       // Cluster width on the y-axis
0245       meSizeY_onTrackPhi_ =
0246           iBooker.book1D("sizeY_OnTrack_" + hisID, "Cluster y-width (columns) (on Track)", 15, 0., 15.);
0247       meSizeY_onTrackPhi_->setAxisTitle("Cluster y-size [columns]", 1);
0248       // Cluster width on the x-axis
0249       meSizeX_offTrackPhi_ =
0250           iBooker.book1D("sizeX_OffTrack_" + hisID, "Cluster x-width (rows) (off Track)", 10, 0., 10.);
0251       meSizeX_offTrackPhi_->setAxisTitle("Cluster x-size [rows]", 1);
0252       // Cluster width on the y-axis
0253       meSizeY_offTrackPhi_ =
0254           iBooker.book1D("sizeY_OffTrack_" + hisID, "Cluster y-width (columns) (off Track)", 15, 0., 15.);
0255       meSizeY_offTrackPhi_->setAxisTitle("Cluster y-size [columns]", 1);
0256     }
0257   }
0258 
0259   if (type == 4 && endcap) {
0260     uint32_t blade;
0261     blade = PixelEndcapName(DetId(id_), pTT, isUpgrade).bladeName();
0262     char sblade[80];
0263     sprintf(sblade, "Blade_%02i", blade);
0264     hisID = src.label() + "_" + sblade;
0265     meResidualXBlade_ = iBooker.book1D("residualX_" + hisID, "Hit-to-Track Residual in r-phi", 100, -150, 150);
0266     meResidualXBlade_->setAxisTitle("hit-to-track residual in r-phi (um)", 1);
0267     meResidualYBlade_ = iBooker.book1D("residualY_" + hisID, "Hit-to-Track Residual in Z", 100, -300, 300);
0268     meResidualYBlade_->setAxisTitle("hit-to-track residual in z (um)", 1);
0269     // Number of clusters
0270     meNClusters_onTrackBlade_ =
0271         iBooker.book1D("nclusters_OnTrack_" + hisID, "Number of Clusters (on Track)", 10, 0., 10.);
0272     meNClusters_onTrackBlade_->setAxisTitle("Number of Clusters on Track", 1);
0273     // Total cluster charge in MeV
0274     meCharge_onTrackBlade_ =
0275         iBooker.book1D("charge_OnTrack_" + hisID, "Normalized Cluster charge (on Track)", 100, 0., 200.);
0276     meCharge_onTrackBlade_->setAxisTitle("Charge [kilo electrons]", 1);
0277     // Total cluster size (in pixels)
0278     meSize_onTrackBlade_ = iBooker.book1D("size_OnTrack_" + hisID, "Total cluster size (on Track)", 30, 0., 30.);
0279     meSize_onTrackBlade_->setAxisTitle("Cluster size [number of pixels]", 1);
0280     // Number of clusters
0281     meNClusters_offTrackBlade_ =
0282         iBooker.book1D("nclusters_OffTrack_" + hisID, "Number of Clusters (off Track)", 35, 0., 35.);
0283     meNClusters_offTrackBlade_->setAxisTitle("Number of Clusters off Track", 1);
0284     // Total cluster charge in MeV
0285     meCharge_offTrackBlade_ = iBooker.book1D("charge_OffTrack_" + hisID, "Cluster charge (off Track)", 100, 0., 200.);
0286     meCharge_offTrackBlade_->setAxisTitle("Charge [kilo electrons]", 1);
0287     // Total cluster size (in pixels)
0288     meSize_offTrackBlade_ = iBooker.book1D("size_OffTrack_" + hisID, "Total cluster size (off Track)", 30, 0., 30.);
0289     meSize_offTrackBlade_->setAxisTitle("Cluster size [number of pixels]", 1);
0290     if (!reducedSet) {
0291       // Cluster width on the x-axis
0292       meSizeX_onTrackBlade_ =
0293           iBooker.book1D("sizeX_OnTrack_" + hisID, "Cluster x-width (rows) (on Track)", 10, 0., 10.);
0294       meSizeX_onTrackBlade_->setAxisTitle("Cluster x-size [rows]", 1);
0295       // Cluster width on the y-axis
0296       meSizeY_onTrackBlade_ =
0297           iBooker.book1D("sizeY_OnTrack_" + hisID, "Cluster y-width (columns) (on Track)", 15, 0., 15.);
0298       meSizeY_onTrackBlade_->setAxisTitle("Cluster y-size [columns]", 1);
0299       // Cluster width on the x-axis
0300       meSizeX_offTrackBlade_ =
0301           iBooker.book1D("sizeX_OffTrack_" + hisID, "Cluster x-width (rows) (off Track)", 10, 0., 10.);
0302       meSizeX_offTrackBlade_->setAxisTitle("Cluster x-size [rows]", 1);
0303       // Cluster width on the y-axis
0304       meSizeY_offTrackBlade_ =
0305           iBooker.book1D("sizeY_OffTrack_" + hisID, "Cluster y-width (columns) (off Track)", 15, 0., 15.);
0306       meSizeY_offTrackBlade_->setAxisTitle("Cluster y-size [columns]", 1);
0307     }
0308   }
0309 
0310   if (type == 5 && endcap) {
0311     uint32_t disk;
0312     disk = PixelEndcapName(DetId(id_), pTT, isUpgrade).diskName();
0313 
0314     char sdisk[80];
0315     sprintf(sdisk, "Disk_%i", disk);
0316     hisID = src.label() + "_" + sdisk;
0317     meResidualXDisk_ = iBooker.book1D("residualX_" + hisID, "Hit-to-Track Residual in r-phi", 100, -150, 150);
0318     meResidualXDisk_->setAxisTitle("hit-to-track residual in r-phi (um)", 1);
0319     meResidualYDisk_ = iBooker.book1D("residualY_" + hisID, "Hit-to-Track Residual in Z", 100, -300, 300);
0320     meResidualYDisk_->setAxisTitle("hit-to-track residual in z (um)", 1);
0321     // Number of clusters
0322     meNClusters_onTrackDisk_ =
0323         iBooker.book1D("nclusters_OnTrack_" + hisID, "Number of Clusters (on Track)", 10, 0., 10.);
0324     meNClusters_onTrackDisk_->setAxisTitle("Number of Clusters on Track", 1);
0325     // Total cluster charge in MeV
0326     meCharge_onTrackDisk_ =
0327         iBooker.book1D("charge_OnTrack_" + hisID, "Normalized Cluster charge (on Track)", 100, 0., 200.);
0328     meCharge_onTrackDisk_->setAxisTitle("Charge [kilo electrons]", 1);
0329     // Total cluster size (in pixels)
0330     meSize_onTrackDisk_ = iBooker.book1D("size_OnTrack_" + hisID, "Total cluster size (on Track)", 30, 0., 30.);
0331     meSize_onTrackDisk_->setAxisTitle("Cluster size [number of pixels]", 1);
0332     // Number of clusters
0333     meNClusters_offTrackDisk_ =
0334         iBooker.book1D("nclusters_OffTrack_" + hisID, "Number of Clusters (off Track)", 35, 0., 35.);
0335     meNClusters_offTrackDisk_->setAxisTitle("Number of Clusters off Track", 1);
0336     // Total cluster charge in MeV
0337     meCharge_offTrackDisk_ = iBooker.book1D("charge_OffTrack_" + hisID, "Cluster charge (off Track)", 100, 0., 200.);
0338     meCharge_offTrackDisk_->setAxisTitle("Charge [kilo electrons]", 1);
0339     // Total cluster size (in pixels)
0340     meSize_offTrackDisk_ = iBooker.book1D("size_OffTrack_" + hisID, "Total cluster size (off Track)", 30, 0., 30.);
0341     meSize_offTrackDisk_->setAxisTitle("Cluster size [number of pixels]", 1);
0342     if (!reducedSet) {
0343       // Cluster width on the x-axis
0344       meSizeX_onTrackDisk_ = iBooker.book1D("sizeX_OnTrack_" + hisID, "Cluster x-width (rows) (on Track)", 10, 0., 10.);
0345       meSizeX_onTrackDisk_->setAxisTitle("Cluster x-size [rows]", 1);
0346       // Cluster width on the y-axis
0347       meSizeY_onTrackDisk_ =
0348           iBooker.book1D("sizeY_OnTrack_" + hisID, "Cluster y-width (columns) (on Track)", 15, 0., 15.);
0349       meSizeY_onTrackDisk_->setAxisTitle("Cluster y-size [columns]", 1);
0350       // Cluster width on the x-axis
0351       meSizeX_offTrackDisk_ =
0352           iBooker.book1D("sizeX_OffTrack_" + hisID, "Cluster x-width (rows) (off Track)", 10, 0., 10.);
0353       meSizeX_offTrackDisk_->setAxisTitle("Cluster x-size [rows]", 1);
0354       // Cluster width on the y-axis
0355       meSizeY_offTrackDisk_ =
0356           iBooker.book1D("sizeY_OffTrack_" + hisID, "Cluster y-width (columns) (off Track)", 15, 0., 15.);
0357       meSizeY_offTrackDisk_->setAxisTitle("Cluster y-size [columns]", 1);
0358     }
0359   }
0360 
0361   if (type == 6 && endcap) {
0362     uint32_t panel;
0363     uint32_t module;
0364     panel = PixelEndcapName(DetId(id_), pTT, isUpgrade).pannelName();
0365     module = PixelEndcapName(DetId(id_), pTT, isUpgrade).plaquetteName();
0366 
0367     char slab[80];
0368     sprintf(slab, "Panel_%i_Ring_%i", panel, module);
0369     hisID = src.label() + "_" + slab;
0370     meResidualXRing_ = iBooker.book1D("residualX_" + hisID, "Hit-to-Track Residual in r-phi", 100, -150, 150);
0371     meResidualXRing_->setAxisTitle("hit-to-track residual in r-phi (um)", 1);
0372     meResidualYRing_ = iBooker.book1D("residualY_" + hisID, "Hit-to-Track Residual in Z", 100, -300, 300);
0373     meResidualYRing_->setAxisTitle("hit-to-track residual in z (um)", 1);
0374     // Number of clusters
0375     meNClusters_onTrackRing_ =
0376         iBooker.book1D("nclusters_OnTrack_" + hisID, "Number of Clusters (on Track)", 10, 0., 10.);
0377     meNClusters_onTrackRing_->setAxisTitle("Number of Clusters on Track", 1);
0378     // Total cluster charge in MeV
0379     meCharge_onTrackRing_ =
0380         iBooker.book1D("charge_OnTrack_" + hisID, "Normalized Cluster charge (on Track)", 100, 0., 200.);
0381     meCharge_onTrackRing_->setAxisTitle("Charge [kilo electrons]", 1);
0382     // Total cluster size (in pixels)
0383     meSize_onTrackRing_ = iBooker.book1D("size_OnTrack_" + hisID, "Total cluster size (on Track)", 30, 0., 30.);
0384     meSize_onTrackRing_->setAxisTitle("Cluster size [number of pixels]", 1);
0385     // Number of clusters
0386     meNClusters_offTrackRing_ =
0387         iBooker.book1D("nclusters_OffTrack_" + hisID, "Number of Clusters (off Track)", 35, 0., 35.);
0388     meNClusters_offTrackRing_->setAxisTitle("Number of Clusters off Track", 1);
0389     // Total cluster charge in MeV
0390     meCharge_offTrackRing_ = iBooker.book1D("charge_OffTrack_" + hisID, "Cluster charge (off Track)", 100, 0., 200.);
0391     meCharge_offTrackRing_->setAxisTitle("Charge [kilo electrons]", 1);
0392     // Total cluster size (in pixels)
0393     meSize_offTrackRing_ = iBooker.book1D("size_OffTrack_" + hisID, "Total cluster size (off Track)", 30, 0., 30.);
0394     meSize_offTrackRing_->setAxisTitle("Cluster size [number of pixels]", 1);
0395     if (!reducedSet) {
0396       // Cluster width on the x-axis
0397       meSizeX_onTrackRing_ = iBooker.book1D("sizeX_OnTrack_" + hisID, "Cluster x-width (rows) (on Track)", 10, 0., 10.);
0398       meSizeX_onTrackRing_->setAxisTitle("Cluster x-size [rows]", 1);
0399       // Cluster width on the y-axis
0400       meSizeY_onTrackRing_ =
0401           iBooker.book1D("sizeY_OnTrack_" + hisID, "Cluster y-width (columns) (on Track)", 15, 0., 15.);
0402       meSizeY_onTrackRing_->setAxisTitle("Cluster y-size [columns]", 1);
0403       // Cluster width on the x-axis
0404       meSizeX_offTrackRing_ =
0405           iBooker.book1D("sizeX_OffTrack_" + hisID, "Cluster x-width (rows) (off Track)", 10, 0., 10.);
0406       meSizeX_offTrackRing_->setAxisTitle("Cluster x-size [rows]", 1);
0407       // Cluster width on the y-axis
0408       meSizeY_offTrackRing_ =
0409           iBooker.book1D("sizeY_OffTrack_" + hisID, "Cluster y-width (columns) (off Track)", 15, 0., 15.);
0410       meSizeY_offTrackRing_->setAxisTitle("Cluster y-size [columns]", 1);
0411     }
0412   }
0413 }
0414 
0415 void SiPixelTrackResidualModule::fill(const Measurement2DVector &residual,
0416                                       bool reducedSet,
0417                                       bool modon,
0418                                       bool ladon,
0419                                       bool layon,
0420                                       bool phion,
0421                                       bool bladeon,
0422                                       bool diskon,
0423                                       bool ringon) {
0424   bool barrel = DetId(id_).subdetId() == static_cast<int>(PixelSubdetector::PixelBarrel);
0425   bool endcap = DetId(id_).subdetId() == static_cast<int>(PixelSubdetector::PixelEndcap);
0426 
0427   if (modon) {
0428     (meResidualX_)->Fill(residual.x());
0429     (meResidualY_)->Fill(residual.y());
0430   }
0431 
0432   if (ladon && barrel) {
0433     (meResidualXLad_)->Fill(residual.x());
0434     (meResidualYLad_)->Fill(residual.y());
0435   }
0436 
0437   if (layon && barrel) {
0438     (meResidualXLay_)->Fill(residual.x());
0439     (meResidualYLay_)->Fill(residual.y());
0440   }
0441   if (phion && barrel) {
0442     (meResidualXPhi_)->Fill(residual.x());
0443     (meResidualYPhi_)->Fill(residual.y());
0444   }
0445 
0446   if (bladeon && endcap) {
0447     (meResidualXBlade_)->Fill(residual.x());
0448     (meResidualYBlade_)->Fill(residual.y());
0449   }
0450 
0451   if (diskon && endcap) {
0452     (meResidualXDisk_)->Fill(residual.x());
0453     (meResidualYDisk_)->Fill(residual.y());
0454   }
0455 
0456   if (ringon && endcap) {
0457     (meResidualXRing_)->Fill(residual.x());
0458     (meResidualYRing_)->Fill(residual.y());
0459   }
0460 }
0461 
0462 void SiPixelTrackResidualModule::fill(const SiPixelCluster &clust,
0463                                       bool onTrack,
0464                                       double corrCharge,
0465                                       bool reducedSet,
0466                                       bool modon,
0467                                       bool ladon,
0468                                       bool layon,
0469                                       bool phion,
0470                                       bool bladeon,
0471                                       bool diskon,
0472                                       bool ringon) {
0473   bool barrel = DetId(id_).subdetId() == static_cast<int>(PixelSubdetector::PixelBarrel);
0474   bool endcap = DetId(id_).subdetId() == static_cast<int>(PixelSubdetector::PixelEndcap);
0475 
0476   float charge = 0.001 * (clust.charge());  // total charge of cluster
0477   if (onTrack)
0478     charge = corrCharge;      // corrected cluster charge
0479   int size = clust.size();    // total size of cluster (in pixels)
0480   int sizeX = clust.sizeX();  // size of cluster in x-clustrection
0481   int sizeY = clust.sizeY();  // size of cluster in y-direction
0482 
0483   if (onTrack) {
0484     if (modon) {
0485       (meCharge_onTrack_)->Fill((float)charge);
0486       (meSize_onTrack_)->Fill((int)size);
0487       if (!reducedSet) {
0488         (meSizeX_onTrack_)->Fill((int)sizeX);
0489         (meSizeY_onTrack_)->Fill((int)sizeY);
0490       }
0491     }
0492     if (barrel && ladon) {
0493       (meCharge_onTrackLad_)->Fill((float)charge);
0494       (meSize_onTrackLad_)->Fill((int)size);
0495       if (!reducedSet) {
0496         (meSizeX_onTrackLad_)->Fill((int)sizeX);
0497         (meSizeY_onTrackLad_)->Fill((int)sizeY);
0498       }
0499     }
0500     if (barrel && layon) {
0501       (meCharge_onTrackLay_)->Fill((float)charge);
0502       (meSize_onTrackLay_)->Fill((int)size);
0503       if (!reducedSet) {
0504         (meSizeX_onTrackLay_)->Fill((int)sizeX);
0505         (meSizeY_onTrackLay_)->Fill((int)sizeY);
0506       }
0507     }
0508     if (barrel && phion) {
0509       (meCharge_onTrackPhi_)->Fill((float)charge);
0510       (meSize_onTrackPhi_)->Fill((int)size);
0511       if (!reducedSet) {
0512         (meSizeX_onTrackPhi_)->Fill((int)sizeX);
0513         (meSizeY_onTrackPhi_)->Fill((int)sizeY);
0514       }
0515     }
0516     if (endcap && bladeon) {
0517       (meCharge_onTrackBlade_)->Fill((float)charge);
0518       (meSize_onTrackBlade_)->Fill((int)size);
0519       if (!reducedSet) {
0520         (meSizeX_onTrackBlade_)->Fill((int)sizeX);
0521         (meSizeY_onTrackBlade_)->Fill((int)sizeY);
0522       }
0523     }
0524     if (endcap && diskon) {
0525       (meCharge_onTrackDisk_)->Fill((float)charge);
0526       (meSize_onTrackDisk_)->Fill((int)size);
0527       if (!reducedSet) {
0528         (meSizeX_onTrackDisk_)->Fill((int)sizeX);
0529         (meSizeY_onTrackDisk_)->Fill((int)sizeY);
0530       }
0531     }
0532     if (endcap && ringon) {
0533       (meCharge_onTrackRing_)->Fill((float)charge);
0534       (meSize_onTrackRing_)->Fill((int)size);
0535       if (!reducedSet) {
0536         (meSizeX_onTrackRing_)->Fill((int)sizeX);
0537         (meSizeY_onTrackRing_)->Fill((int)sizeY);
0538       }
0539     }
0540   }
0541 
0542   if (!onTrack) {
0543     if (modon) {
0544       (meCharge_offTrack_)->Fill((float)charge);
0545       (meSize_offTrack_)->Fill((int)size);
0546       if (!reducedSet) {
0547         (meSizeX_offTrack_)->Fill((int)sizeX);
0548         (meSizeY_offTrack_)->Fill((int)sizeY);
0549       }
0550     }
0551     if (barrel && ladon) {
0552       (meCharge_offTrackLad_)->Fill((float)charge);
0553       (meSize_offTrackLad_)->Fill((int)size);
0554       if (!reducedSet) {
0555         (meSizeX_offTrackLad_)->Fill((int)sizeX);
0556         (meSizeY_offTrackLad_)->Fill((int)sizeY);
0557       }
0558     }
0559     if (barrel && layon) {
0560       (meCharge_offTrackLay_)->Fill((float)charge);
0561       (meSize_offTrackLay_)->Fill((int)size);
0562       if (!reducedSet) {
0563         (meSizeX_offTrackLay_)->Fill((int)sizeX);
0564         (meSizeY_offTrackLay_)->Fill((int)sizeY);
0565       }
0566     }
0567     if (barrel && phion) {
0568       (meCharge_offTrackPhi_)->Fill((float)charge);
0569       (meSize_offTrackPhi_)->Fill((int)size);
0570       if (!reducedSet) {
0571         (meSizeX_offTrackPhi_)->Fill((int)sizeX);
0572         (meSizeY_offTrackPhi_)->Fill((int)sizeY);
0573       }
0574     }
0575     if (endcap && bladeon) {
0576       (meCharge_offTrackBlade_)->Fill((float)charge);
0577       (meSize_offTrackBlade_)->Fill((int)size);
0578       if (!reducedSet) {
0579         (meSizeX_offTrackBlade_)->Fill((int)sizeX);
0580         (meSizeY_offTrackBlade_)->Fill((int)sizeY);
0581       }
0582     }
0583     if (endcap && diskon) {
0584       (meCharge_offTrackDisk_)->Fill((float)charge);
0585       (meSize_offTrackDisk_)->Fill((int)size);
0586       if (!reducedSet) {
0587         (meSizeX_offTrackDisk_)->Fill((int)sizeX);
0588         (meSizeY_offTrackDisk_)->Fill((int)sizeY);
0589       }
0590     }
0591     if (endcap && ringon) {
0592       (meCharge_offTrackRing_)->Fill((float)charge);
0593       (meSize_offTrackRing_)->Fill((int)size);
0594       if (!reducedSet) {
0595         (meSizeX_offTrackRing_)->Fill((int)sizeX);
0596         (meSizeY_offTrackRing_)->Fill((int)sizeY);
0597       }
0598     }
0599   }
0600 }
0601 
0602 void SiPixelTrackResidualModule::nfill(int onTrack,
0603                                        int offTrack,
0604                                        bool reducedSet,
0605                                        bool modon,
0606                                        bool ladon,
0607                                        bool layon,
0608                                        bool phion,
0609                                        bool bladeon,
0610                                        bool diskon,
0611                                        bool ringon) {
0612   bool barrel = DetId(id_).subdetId() == static_cast<int>(PixelSubdetector::PixelBarrel);
0613   bool endcap = DetId(id_).subdetId() == static_cast<int>(PixelSubdetector::PixelEndcap);
0614   bool fillOn = false;
0615   if (onTrack > 0)
0616     fillOn = true;
0617   bool fillOff = false;
0618   if (offTrack > 0)
0619     fillOff = true;
0620 
0621   if (modon) {
0622     if (fillOn)
0623       meNClusters_onTrack_->Fill(onTrack);
0624     if (fillOff)
0625       meNClusters_offTrack_->Fill(offTrack);
0626   }
0627   if (ladon && barrel) {
0628     if (fillOn)
0629       meNClusters_onTrackLad_->Fill(onTrack);
0630     if (fillOff)
0631       meNClusters_offTrackLad_->Fill(offTrack);
0632   }
0633 
0634   if (layon && barrel) {
0635     if (fillOn)
0636       meNClusters_onTrackLay_->Fill(onTrack);
0637     if (fillOff)
0638       meNClusters_offTrackLay_->Fill(offTrack);
0639   }
0640   if (phion && barrel) {
0641     if (fillOn)
0642       meNClusters_onTrackPhi_->Fill(onTrack);
0643     if (fillOff)
0644       meNClusters_offTrackPhi_->Fill(offTrack);
0645   }
0646   if (bladeon && endcap) {
0647     if (fillOn)
0648       meNClusters_onTrackBlade_->Fill(onTrack);
0649     if (fillOff)
0650       meNClusters_offTrackBlade_->Fill(offTrack);
0651   }
0652   if (diskon && endcap) {
0653     if (fillOn)
0654       meNClusters_onTrackDisk_->Fill(onTrack);
0655     if (fillOff)
0656       meNClusters_offTrackDisk_->Fill(offTrack);
0657   }
0658   if (ringon && endcap) {
0659     if (fillOn)
0660       meNClusters_onTrackRing_->Fill(onTrack);
0661     if (fillOff)
0662       meNClusters_offTrackRing_->Fill(offTrack);
0663   }
0664 }