Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 // -*- C++ -*-
0002 //
0003 // Package:    SiPixelPhase1ResidualsExtra
0004 // Class:      SiPixelPhase1ResidualsExtra
0005 //
0006 /**\class 
0007 
0008  Description: Create the Phsae 1 pixel DRnR plots
0009 
0010  Implementation:
0011      <Notes on implementation>
0012 */
0013 //
0014 // Original Author:  Alessandro Rossi
0015 //         Created:  25th May 2021
0016 //
0017 //
0018 #include "DQM/SiPixelPhase1Track/interface/SiPixelPhase1ResidualsExtra.h"
0019 // Framework
0020 #include "FWCore/ServiceRegistry/interface/Service.h"
0021 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0022 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0023 // DQM Framework
0024 #include "DQM/SiPixelCommon/interface/SiPixelFolderOrganizer.h"
0025 #include "DQMServices/Core/interface/DQMStore.h"
0026 // Geometry
0027 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
0028 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
0029 #include "Geometry/CommonDetUnit/interface/PixelGeomDetUnit.h"
0030 #include "Geometry/CommonTopologies/interface/PixelTopology.h"
0031 #include "Geometry/Records/interface/TrackerTopologyRcd.h"
0032 // DataFormats
0033 #include "DataFormats/DetId/interface/DetId.h"
0034 #include "DataFormats/SiPixelDetId/interface/PixelSubdetector.h"
0035 #include "DataFormats/TrackerCommon/interface/PixelBarrelName.h"
0036 #include "DataFormats/SiPixelDetId/interface/PixelBarrelNameUpgrade.h"
0037 #include "DataFormats/TrackerCommon/interface/PixelEndcapName.h"
0038 #include "DataFormats/SiPixelDetId/interface/PixelEndcapNameUpgrade.h"
0039 //
0040 #include <string>
0041 #include <cstdlib>
0042 #include <iostream>
0043 #include <fstream>
0044 #include <sstream>
0045 
0046 using namespace std;
0047 using namespace edm;
0048 
0049 SiPixelPhase1ResidualsExtra::SiPixelPhase1ResidualsExtra(const edm::ParameterSet& iConfig)
0050     : DQMEDHarvester(iConfig), conf_(iConfig) {
0051   LogInfo("PixelDQM") << "SiPixelPhase1ResidualsExtra::SiPixelPhase1ResidualsExtra: Got DQM BackEnd interface" << endl;
0052   topFolderName_ = conf_.getParameter<std::string>("TopFolderName");
0053   minHits_ = conf_.getParameter<int>("MinHits");
0054 }
0055 
0056 SiPixelPhase1ResidualsExtra::~SiPixelPhase1ResidualsExtra() {
0057   // do anything here that needs to be done at desctruction time
0058   // (e.g. close files, deallocate resources etc.)
0059   LogInfo("PixelDQM") << "SiPixelPhase1ResidualsExtra::~SiPixelPhase1ResidualsExtra: Destructor" << endl;
0060 }
0061 
0062 void SiPixelPhase1ResidualsExtra::beginRun(edm::Run const& run, edm::EventSetup const& eSetup) {}
0063 
0064 void SiPixelPhase1ResidualsExtra::dqmEndJob(DQMStore::IBooker& iBooker, DQMStore::IGetter& iGetter) {
0065   bookMEs(iBooker);
0066   fillMEs(iBooker, iGetter);
0067 }
0068 
0069 //------------------------------------------------------------------
0070 // Used to book the MEs
0071 //------------------------------------------------------------------
0072 void SiPixelPhase1ResidualsExtra::bookMEs(DQMStore::IBooker& iBooker) {
0073   iBooker.cd();
0074 
0075   //New residual plots for the PXBarrel separated by inner and outer modules per layer
0076   iBooker.setCurrentFolder(topFolderName_ + "/PXBarrel");
0077 
0078   for (std::string layer : {"1", "2", "3", "4"}) {
0079     float mean_range = 20.;
0080     float rms_range = 200.;
0081     if (layer == "1") {
0082       mean_range = 50.;
0083       rms_range = 1000.;
0084     }
0085     residuals_["residual_mean_x_Inner_PXLayer_" + layer] =
0086         iBooker.book1D("residual_mean_x_Inner_PXLayer_" + layer,
0087                        "Mean of Track Residuals X Inner Modules for Layer " + layer + ";mean(x_rec-x_pred)[#mum]",
0088                        100,
0089                        -1 * mean_range,
0090                        mean_range);
0091     residuals_["residual_mean_x_Outer_PXLayer_" + layer] =
0092         iBooker.book1D("residual_mean_x_Outer_PXLayer_" + layer,
0093                        "Mean of Track Residuals X Outer Modules for Layer " + layer + ";mean(x_rec-x_pred)[#mum]",
0094                        100,
0095                        -1 * mean_range,
0096                        mean_range);
0097     residuals_["residual_mean_y_Inner_PXLayer_" + layer] =
0098         iBooker.book1D("residual_mean_y_Inner_PXLayer_" + layer,
0099                        "Mean of Track Residuals Y Inner Modules for Layer " + layer + ";mean(y_rec-y_pred)[#mum]",
0100                        100,
0101                        -1 * mean_range,
0102                        mean_range);
0103     residuals_["residual_mean_y_Outer_PXLayer_" + layer] =
0104         iBooker.book1D("residual_mean_y_Outer_PXLayer_" + layer,
0105                        "Mean of Track Residuals Y Outer Modules for Layer " + layer + ";mean(y_rec-y_pred)[#mum]",
0106                        100,
0107                        -1 * mean_range,
0108                        mean_range);
0109 
0110     residuals_["residual_rms_x_Inner_PXLayer_" + layer] =
0111         iBooker.book1D("residual_rms_x_Inner_PXLayer_" + layer,
0112                        "RMS of Track Residuals X Inner Modules for Layer " + layer + ";rms(x_rec-x_pred)[#mum]",
0113                        100,
0114                        0.,
0115                        rms_range);
0116     residuals_["residual_rms_x_Outer_PXLayer_" + layer] =
0117         iBooker.book1D("residual_rms_x_Outer_PXLayer_" + layer,
0118                        "RMS of Track Residuals X Outer Modules for Layer " + layer + ";rms(x_rec-x_pred)[#mum]",
0119                        100,
0120                        0.,
0121                        rms_range);
0122     residuals_["residual_rms_y_Inner_PXLayer_" + layer] =
0123         iBooker.book1D("residual_rms_y_Inner_PXLayer_" + layer,
0124                        "RMS of Track Residuals Y Inner Modules for Layer " + layer + ";rms(y_rec-y_pred)[#mum]",
0125                        100,
0126                        0.,
0127                        rms_range);
0128     residuals_["residual_rms_y_Outer_PXLayer_" + layer] =
0129         iBooker.book1D("residual_rms_y_Outer_PXLayer_" + layer,
0130                        "RMS of Track Residuals Y Outer Modules for Layer " + layer + ";rms(y_rec-y_pred)[#mum]",
0131                        100,
0132                        0.,
0133                        rms_range);
0134     ///Normalized Resiuduals Plots
0135     DRnR_["NormRes_mean_x_Inner_PXLayer_" + layer] = iBooker.book1D(
0136         "NormRes_mean_x_Inner_PXLayer_" + layer,
0137         "Mean of Normalized Track Residuals X Inner Modules for Layer " + layer + ";mean((x_rec-x_pred)/x_err)",
0138         100,
0139         -1 * 1,
0140         1);
0141     DRnR_["NormRes_mean_x_Outer_PXLayer_" + layer] = iBooker.book1D(
0142         "NormRes_mean_x_Outer_PXLayer_" + layer,
0143         "Mean of Normalized Track Residuals X Outer Modules for Layer " + layer + ";mean((x_rec-x_pred)/x_err)",
0144         100,
0145         -1 * 1,
0146         1);
0147     DRnR_["NormRes_mean_y_Inner_PXLayer_" + layer] = iBooker.book1D(
0148         "NormRes_mean_y_Inner_PXLayer_" + layer,
0149         "Mean of Normalized Track Residuals Y Inner Modules for Layer " + layer + ";mean((y_rec-y_pred)/y_err)",
0150         100,
0151         -1 * 1,
0152         1);
0153     DRnR_["NormRes_mean_y_Outer_PXLayer_" + layer] = iBooker.book1D(
0154         "NormRes_mean_y_Outer_PXLayer_" + layer,
0155         "Mean of Normalized Track Residuals Y Outer Modules for Layer " + layer + ";mean((y_rec-y_pred)/y_err)",
0156         100,
0157         -1 * 1,
0158         1);
0159 
0160     DRnR_["DRnR_x_Inner_PXLayer_" + layer] = iBooker.book1D(
0161         "DRnR_x_Inner_PXLayer_" + layer,
0162         "RMS of Normalized Track Residuals X Inner Modules for Layer " + layer + ";rms((x_rec-x_pred)/x_err)",
0163         100,
0164         0.,
0165         2);
0166     DRnR_["DRnR_x_Outer_PXLayer_" + layer] = iBooker.book1D(
0167         "DRnR_x_Outer_PXLayer_" + layer,
0168         "RMS of Normalized Track Residuals X Outer Modules for Layer " + layer + ";rms((x_rec-x_pred)/x_err)",
0169         100,
0170         0.,
0171         2);
0172     DRnR_["DRnR_y_Inner_PXLayer_" + layer] = iBooker.book1D(
0173         "DRnR_y_Inner_PXLayer_" + layer,
0174         "RMS of Normalized Track Residuals Y Inner Modules for Layer " + layer + ";rms((y_rec-y_pred)/y_err)",
0175         100,
0176         0.,
0177         2);
0178     DRnR_["DRnR_y_Outer_PXLayer_" + layer] = iBooker.book1D(
0179         "DRnR_y_Outer_PXLayer_" + layer,
0180         "RMS of Normalized Track Residuals Y Outer Modules for Layer " + layer + ";rms((y_rec-y_pred)/y_err)",
0181         100,
0182         0.,
0183         2);
0184   }
0185 
0186   //New residual plots for the PXForward separated by inner and outer modules
0187   iBooker.setCurrentFolder(topFolderName_ + "/PXForward");
0188 
0189   residuals_["residual_mean_x_Inner"] = iBooker.book1D(
0190       "residual_mean_x_Inner", "Mean of Track Residuals X Inner Modules;mean(x_rec-x_pred)[#mum]", 100, -20., 20.);
0191   residuals_["residual_mean_x_Outer"] = iBooker.book1D(
0192       "residual_mean_x_Outer", "Mean of Track Residuals X Outer Modules;mean(x_rec-x_pred)[#mum]", 100, -20., 20.);
0193   residuals_["residual_mean_y_Inner"] = iBooker.book1D(
0194       "residual_mean_y_Inner", "Mean of Track Residuals Y Inner Modules;mean(y_rec-y_pred)[#mum]", 100, -20., 20.);
0195   residuals_["residual_mean_y_Outer"] = iBooker.book1D(
0196       "residual_mean_y_Outer", "Mean of Track Residuals Y Outer Modules;mean(y_rec-y_pred)[#mum]", 100, -20., 20.);
0197 
0198   residuals_["residual_rms_x_Inner"] = iBooker.book1D(
0199       "residual_rms_x_Inner", "RMS of Track Residuals X Inner Modules;rms(x_rec-x_pred)[#mum]", 100, 0., 200.);
0200   residuals_["residual_rms_x_Outer"] = iBooker.book1D(
0201       "residual_rms_x_Outer", "RMS of Track Residuals X Outer Modules;rms(x_rec-x_pred)[#mum]", 100, 0., 200.);
0202   residuals_["residual_rms_y_Inner"] = iBooker.book1D(
0203       "residual_rms_y_Inner", "RMS of Track Residuals Y Inner Modules;rms(y_rec-y_pred)[#mum]", 100, 0., 200.);
0204   residuals_["residual_rms_y_Outer"] = iBooker.book1D(
0205       "residual_rms_y_Outer", "RMS of Track Residuals Y Outer Modules;rms(y_rec-y_pred)[#mum]", 100, 0., 200.);
0206   //Normalize Residuals inner/outer
0207   DRnR_["NormRes_mean_x_Inner"] =
0208       iBooker.book1D("NormRes_mean_x_Inner",
0209                      "Mean of Normalized Track Residuals X Inner Modules;mean((x_rec-x_pred)/x_err)",
0210                      100,
0211                      -1.,
0212                      1.);
0213   DRnR_["NormRes_mean_x_Outer"] =
0214       iBooker.book1D("NormRes_mean_x_Outer",
0215                      "Mean of Normalized Track Residuals X Outer Modules;mean((x_rec-x_pred)/x_err)",
0216                      100,
0217                      -1.,
0218                      1.);
0219   DRnR_["NormRes_mean_y_Inner"] =
0220       iBooker.book1D("NormRes_mean_y_Inner",
0221                      "Mean of Normalized Track Residuals Y Inner Modules;mean((y_rec-y_pred)/y_err)",
0222                      100,
0223                      -1.,
0224                      1.);
0225   DRnR_["NormRes_mean_y_Outer"] =
0226       iBooker.book1D("NormRes_mean_y_Outer",
0227                      "Mean of Normalized Track Residuals Y Outer Modules;mean((y_rec-y_pred)/y_err)",
0228                      100,
0229                      -1.,
0230                      1.);
0231 
0232   DRnR_["DRnR_x_Inner"] = iBooker.book1D(
0233       "DRnR_x_Inner", "RMS of Normalized Track Residuals X Inner Modules;rms((x_rec-x_pred)/x_err)", 100, 0., 2.);
0234   DRnR_["DRnR_x_Outer"] = iBooker.book1D(
0235       "DRnR_x_Outer", "RMS of Normalized Track Residuals X Outer Modules;rms((x_rec-x_pred)/x_err)", 100, 0., 2.);
0236   DRnR_["DRnR_y_Inner"] = iBooker.book1D(
0237       "DRnR_y_Inner", "RMS of Normalized Track Residuals Y Inner Modules;rms((y_rec-y_pred)/y_err)", 100, 0., 2.);
0238   DRnR_["DRnR_y_Outer"] = iBooker.book1D(
0239       "DRnR_y_Outer", "RMS of Normalized Track Residuals Y Outer Modules;rms((y_rec-y_pred)/y_err)", 100, 0., 2.);
0240 
0241   //New residual plots for the PXForward separated by positive and negative side
0242   iBooker.setCurrentFolder(topFolderName_ + "/PXForward");
0243 
0244   residuals_["residual_mean_x_pos"] = iBooker.book1D(
0245       "residual_mean_x_pos", "Mean of Track Residuals X pos. Side;mean(x_rec-x_pred)[#mum]", 100, -20., 20.);
0246   residuals_["residual_mean_x_neg"] = iBooker.book1D(
0247       "residual_mean_x_neg", "Mean of Track Residuals X neg. Side;mean(x_rec-x_pred)[#mum]", 100, -20., 20.);
0248   residuals_["residual_mean_y_pos"] = iBooker.book1D(
0249       "residual_mean_y_pos", "Mean of Track Residuals Y pos. Side;mean(y_rec-y_pred)[#mum]", 100, -20., 20.);
0250   residuals_["residual_mean_y_neg"] = iBooker.book1D(
0251       "residual_mean_y_neg", "Mean of Track Residuals Y neg. Side;mean(y_rec-y_pred)[#mum]", 100, -20., 20.);
0252 
0253   residuals_["residual_rms_x_pos"] =
0254       iBooker.book1D("residual_rms_x_pos", "RMS of Track Residuals X pos. Side;rms(x_rec-x_pred)[#mum]", 100, 0., 200.);
0255   residuals_["residual_rms_x_neg"] =
0256       iBooker.book1D("residual_rms_x_neg", "RMS of Track Residuals X neg. Side;rms(x_rec-x_pred)[#mum]", 100, 0., 200.);
0257   residuals_["residual_rms_y_pos"] =
0258       iBooker.book1D("residual_rms_y_pos", "RMS of Track Residuals Y pos. Side;rms(y_rec-y_pred)[#mum]", 100, 0., 200.);
0259   residuals_["residual_rms_y_neg"] =
0260       iBooker.book1D("residual_rms_y_neg", "RMS of Track Residuals Y neg. Side;rms(y_rec-y_pred)[#mum]", 100, 0., 200.);
0261   //Normalized Residuals pos/neg
0262   DRnR_["NormRes_mean_x_pos"] = iBooker.book1D(
0263       "NormRes_mean_x_pos", "Mean of Normalized Track Residuals X pos. Side;mean((x_rec-x_pred)/x_err)", 100, -1., 1.);
0264   DRnR_["NormRes_mean_x_neg"] = iBooker.book1D(
0265       "NormRes_mean_x_neg", "Mean of Normalized Track Residuals X neg. Side;mean((x_rec-x_pred)/x_err)", 100, -1., 1.);
0266   DRnR_["NormRes_mean_y_pos"] = iBooker.book1D(
0267       "NormRes_mean_y_pos", "Mean of Normalized Track Residuals Y pos. Side;mean((y_rec-y_pred)/y_err)", 100, -1., 1.);
0268   DRnR_["NormRes_mean_y_neg"] = iBooker.book1D(
0269       "NormRes_mean_y_neg", "Mean of Normalized Track Residuals Y neg. Side;mean((y_rec-y_pred)/y_err)", 100, -1., 1.);
0270 
0271   DRnR_["DRnR_x_pos"] = iBooker.book1D(
0272       "DRnR_x_pos", "RMS of Normalized Track Residuals X pos. Side;rms((x_rec-x_pred)/x_err)", 100, 0., 2.);
0273   DRnR_["DRnR_x_neg"] = iBooker.book1D(
0274       "DRnR_x_neg", "RMS of Normalized Track Residuals X neg. Side;rms((x_rec-x_pred)/x_err)", 100, 0., 2.);
0275   DRnR_["DRnR_y_pos"] = iBooker.book1D(
0276       "DRnR_y_pos", "RMS of Normalized Track Residuals Y pos. Side;rms((y_rec-y_pred)/y_err)", 100, 0., 2.);
0277   DRnR_["DRnR_y_neg"] = iBooker.book1D(
0278       "DRnR_y_neg", "RMS of Normalized Track Residuals Y neg. Side;rms((y_rec-y_pred)/y_err)", 100, 0., 2.);
0279 
0280   //Reset the iBooker
0281   iBooker.setCurrentFolder("PixelPhase1/");
0282 }
0283 
0284 //------------------------------------------------------------------
0285 // Fill the MEs
0286 //------------------------------------------------------------------
0287 void SiPixelPhase1ResidualsExtra::fillMEs(DQMStore::IBooker& iBooker, DQMStore::IGetter& iGetter) {
0288   //Fill additional residuals plots
0289   //PXBarrel
0290 
0291   //constexpr int minHits_ = 30;  //Miniminal number of hits needed for module to be filled in histograms
0292 
0293   for (std::string layer : {"1", "2", "3", "4"}) {
0294     MonitorElement* me_x =
0295         iGetter.get("PixelPhase1/Tracks/PXBarrel/residual_x_per_SignedModule_per_SignedLadder_PXLayer_" + layer);
0296     MonitorElement* me_y =
0297         iGetter.get("PixelPhase1/Tracks/PXBarrel/residual_y_per_SignedModule_per_SignedLadder_PXLayer_" + layer);
0298     MonitorElement* me2_x = iGetter.get(
0299         "PixelPhase1/Tracks/ResidualsExtra/PXBarrel/DRnR_x_per_SignedModule_per_SignedLadder_PXLayer_" + layer);
0300     MonitorElement* me2_y = iGetter.get(
0301         "PixelPhase1/Tracks/ResidualsExtra/PXBarrel/DRnR_y_per_SignedModule_per_SignedLadder_PXLayer_" + layer);
0302 
0303     if (me_x == nullptr || me_y == nullptr || me2_x == nullptr || me2_y == nullptr) {
0304       edm::LogWarning("SiPixelPhase1ResidualsExtra")
0305           << "Residuals plots for Pixel BPIX Layer" << layer << " not found. Skipping ResidualsExtra plots generation.";
0306       continue;
0307     }
0308 
0309     for (int i = 1; i <= me_x->getNbinsY(); i++) {
0310       if (i == (me_x->getNbinsY() / 2 + 1))
0311         continue;  //Middle bin of y axis is empty
0312 
0313       if (i <= me_x->getNbinsY() / 2) {
0314         bool iAmInner = (i % 2 == 0);  //Check whether current ladder is inner or outer ladder
0315         if (iAmInner) {
0316           for (int j : {1, 2, 3, 4, 6, 7, 8, 9}) {
0317             if (me_x->getBinEntries(j, i) < minHits_)  //Fill only if number of hits is above threshold
0318               continue;
0319             residuals_["residual_mean_x_Inner_PXLayer_" + layer]->Fill(me_x->getBinContent(j, i) * 1e4);
0320             residuals_["residual_mean_y_Inner_PXLayer_" + layer]->Fill(me_y->getBinContent(j, i) * 1e4);
0321             residuals_["residual_rms_x_Inner_PXLayer_" + layer]->Fill(me_x->getBinError(j, i) * 1e4);
0322             residuals_["residual_rms_y_Inner_PXLayer_" + layer]->Fill(me_y->getBinError(j, i) * 1e4);
0323             DRnR_["NormRes_mean_x_Inner_PXLayer_" + layer]->Fill(me2_x->getBinContent(j, i));
0324             DRnR_["NormRes_mean_y_Inner_PXLayer_" + layer]->Fill(me2_y->getBinContent(j, i));
0325             DRnR_["DRnR_x_Inner_PXLayer_" + layer]->Fill(me2_x->getBinError(j, i));
0326             DRnR_["DRnR_y_Inner_PXLayer_" + layer]->Fill(me2_y->getBinError(j, i));
0327           }
0328         } else {
0329           for (int j : {1, 2, 3, 4, 6, 7, 8, 9}) {
0330             if (me_x->getBinEntries(j, i) < minHits_)  //Fill only if number of hits is above threshold
0331               continue;
0332             residuals_["residual_mean_x_Outer_PXLayer_" + layer]->Fill(me_x->getBinContent(j, i) * 1e4);
0333             residuals_["residual_mean_y_Outer_PXLayer_" + layer]->Fill(me_y->getBinContent(j, i) * 1e4);
0334             residuals_["residual_rms_x_Outer_PXLayer_" + layer]->Fill(me_x->getBinError(j, i) * 1e4);
0335             residuals_["residual_rms_y_Outer_PXLayer_" + layer]->Fill(me_y->getBinError(j, i) * 1e4);
0336             DRnR_["NormRes_mean_x_Outer_PXLayer_" + layer]->Fill(me2_x->getBinContent(j, i));
0337             DRnR_["NormRes_mean_y_Outer_PXLayer_" + layer]->Fill(me2_y->getBinContent(j, i));
0338             DRnR_["DRnR_x_Outer_PXLayer_" + layer]->Fill(me2_x->getBinError(j, i));
0339             DRnR_["DRnR_y_Outer_PXLayer_" + layer]->Fill(me2_y->getBinError(j, i));
0340           }
0341         }
0342       } else {
0343         bool iAmInner = (i % 2 == 1);  //Check whether current ladder is inner or outer ladder
0344         if (iAmInner) {
0345           for (int j : {1, 2, 3, 4, 6, 7, 8, 9}) {
0346             if (me_x->getBinEntries(j, i) < minHits_)  //Fill only if number of hits is above threshold
0347               continue;
0348             residuals_["residual_mean_x_Inner_PXLayer_" + layer]->Fill(me_x->getBinContent(j, i) * 1e4);
0349             residuals_["residual_mean_y_Inner_PXLayer_" + layer]->Fill(me_y->getBinContent(j, i) * 1e4);
0350             residuals_["residual_rms_x_Inner_PXLayer_" + layer]->Fill(me_x->getBinError(j, i) * 1e4);
0351             residuals_["residual_rms_y_Inner_PXLayer_" + layer]->Fill(me_y->getBinError(j, i) * 1e4);
0352             DRnR_["NormRes_mean_x_Inner_PXLayer_" + layer]->Fill(me2_x->getBinContent(j, i));
0353             DRnR_["NormRes_mean_y_Inner_PXLayer_" + layer]->Fill(me2_y->getBinContent(j, i));
0354             DRnR_["DRnR_x_Inner_PXLayer_" + layer]->Fill(me2_x->getBinError(j, i));
0355             DRnR_["DRnR_y_Inner_PXLayer_" + layer]->Fill(me2_y->getBinError(j, i));
0356           }
0357         } else {
0358           for (int j : {1, 2, 3, 4, 6, 7, 8, 9}) {
0359             if (me_x->getBinEntries(j, i) < minHits_)  //Fill only if number of hits is above threshold
0360               continue;
0361             residuals_["residual_mean_x_Outer_PXLayer_" + layer]->Fill(me_x->getBinContent(j, i) * 1e4);
0362             residuals_["residual_mean_y_Outer_PXLayer_" + layer]->Fill(me_y->getBinContent(j, i) * 1e4);
0363             residuals_["residual_rms_x_Outer_PXLayer_" + layer]->Fill(me_x->getBinError(j, i) * 1e4);
0364             residuals_["residual_rms_y_Outer_PXLayer_" + layer]->Fill(me_y->getBinError(j, i) * 1e4);
0365             DRnR_["NormRes_mean_x_Outer_PXLayer_" + layer]->Fill(me2_x->getBinContent(j, i));
0366             DRnR_["NormRes_mean_y_Outer_PXLayer_" + layer]->Fill(me2_y->getBinContent(j, i));
0367             DRnR_["DRnR_x_Outer_PXLayer_" + layer]->Fill(me2_x->getBinError(j, i));
0368             DRnR_["DRnR_y_Outer_PXLayer_" + layer]->Fill(me2_y->getBinError(j, i));
0369           }
0370         }
0371       }
0372       for (int j : {1, 2, 3, 4, 6, 7, 8, 9}) {
0373         if (me_x->getBinEntries(j, i) < minHits_) {
0374           me2_x->setBinContent(j, i, 0);
0375           me2_y->setBinContent(j, i, 0);
0376           me2_x->setBinEntries(me2_x->getBin(j, i), 0);
0377           me2_y->setBinEntries(me2_y->getBin(j, i), 0);
0378         } else {
0379           me2_x->setBinContent(j, i, me2_x->getBinError(j, i));
0380           me2_y->setBinContent(j, i, me2_y->getBinError(j, i));
0381           me2_x->setBinEntries(me2_x->getBin(j, i), 1);
0382           me2_y->setBinEntries(me2_y->getBin(j, i), 1);
0383         }
0384       }
0385     }
0386   }
0387 
0388   //PXForward separating outer and inner modules as well as positive and negative side
0389   for (std::string ring : {"1", "2"}) {
0390     MonitorElement* me_x =
0391         iGetter.get("PixelPhase1/Tracks/PXForward/residual_x_per_PXDisk_per_SignedBladePanel_PXRing_" + ring);
0392     MonitorElement* me_y =
0393         iGetter.get("PixelPhase1/Tracks/PXForward/residual_y_per_PXDisk_per_SignedBladePanel_PXRing_" + ring);
0394     MonitorElement* me2_x = iGetter.get(
0395         "PixelPhase1/Tracks/ResidualsExtra/PXForward/DRnR_x_per_PXDisk_per_SignedBladePanel_PXRing_" + ring);
0396     MonitorElement* me2_y = iGetter.get(
0397         "PixelPhase1/Tracks/ResidualsExtra/PXForward/DRnR_y_per_PXDisk_per_SignedBladePanel_PXRing_" + ring);
0398 
0399     if (me_x == nullptr || me_y == nullptr || me2_x == nullptr || me2_y == nullptr) {
0400       edm::LogWarning("SiPixelPhase1ResidualsExtra")
0401           << "Residuals plots for Pixel FPIX Ring" << ring << " not found. Skipping ResidualsExtra plots generation.";
0402       continue;
0403     }
0404 
0405     bool posSide = false;
0406     for (int j = 1; j <= me_x->getNbinsX(); j++) {
0407       if (j == 4)
0408         continue;  //fourth x-bin in profile plots is empty
0409 
0410       if (j == 5)
0411         posSide = true;  //change to postive side
0412 
0413       for (int i = 1; i <= me_x->getNbinsY(); i++) {
0414         if (i == me_x->getNbinsY() / 2)
0415           continue;  //Middle bins of y axis is empty
0416         if (i == (me_x->getNbinsY() / 2) + 1)
0417           continue;
0418         if (me_x->getBinEntries(j, i) >= minHits_) {  //Fill only if number of hits is above threshold
0419 
0420           bool iAmInner = (i % 2 == 0);  //separate inner and outer modules
0421           if (iAmInner) {
0422             residuals_["residual_mean_x_Inner"]->Fill(me_x->getBinContent(j, i) * 1e4);
0423             residuals_["residual_mean_y_Inner"]->Fill(me_y->getBinContent(j, i) * 1e4);
0424             residuals_["residual_rms_x_Inner"]->Fill(me_x->getBinError(j, i) * 1e4);
0425             residuals_["residual_rms_y_Inner"]->Fill(me_y->getBinError(j, i) * 1e4);
0426             DRnR_["NormRes_mean_x_Inner"]->Fill(me2_x->getBinContent(j, i));
0427             DRnR_["NormRes_mean_y_Inner"]->Fill(me2_y->getBinContent(j, i));
0428             DRnR_["DRnR_x_Inner"]->Fill(me2_x->getBinError(j, i));
0429             DRnR_["DRnR_y_Inner"]->Fill(me2_y->getBinError(j, i));
0430           } else {
0431             residuals_["residual_mean_x_Outer"]->Fill(me_x->getBinContent(j, i) * 1e4);
0432             residuals_["residual_mean_y_Outer"]->Fill(me_y->getBinContent(j, i) * 1e4);
0433             residuals_["residual_rms_x_Outer"]->Fill(me_x->getBinError(j, i) * 1e4);
0434             residuals_["residual_rms_y_Outer"]->Fill(me_y->getBinError(j, i) * 1e4);
0435             DRnR_["NormRes_mean_x_Outer"]->Fill(me2_x->getBinContent(j, i));
0436             DRnR_["NormRes_mean_y_Outer"]->Fill(me2_y->getBinContent(j, i));
0437             DRnR_["DRnR_x_Outer"]->Fill(me2_x->getBinError(j, i));
0438             DRnR_["DRnR_y_Outer"]->Fill(me2_y->getBinError(j, i));
0439           }
0440 
0441           if (!posSide) {  //separate postive and negative side
0442             residuals_["residual_mean_x_neg"]->Fill(me_x->getBinContent(j, i) * 1e4);
0443             residuals_["residual_mean_y_neg"]->Fill(me_y->getBinContent(j, i) * 1e4);
0444             residuals_["residual_rms_x_neg"]->Fill(me_x->getBinError(j, i) * 1e4);
0445             residuals_["residual_rms_y_neg"]->Fill(me_y->getBinError(j, i) * 1e4);
0446             DRnR_["NormRes_mean_x_neg"]->Fill(me2_x->getBinContent(j, i));
0447             DRnR_["NormRes_mean_y_neg"]->Fill(me2_y->getBinContent(j, i));
0448             DRnR_["DRnR_x_neg"]->Fill(me2_x->getBinError(j, i));
0449             DRnR_["DRnR_y_neg"]->Fill(me2_y->getBinError(j, i));
0450           } else {
0451             residuals_["residual_mean_x_pos"]->Fill(me_x->getBinContent(j, i) * 1e4);
0452             residuals_["residual_mean_y_pos"]->Fill(me_y->getBinContent(j, i) * 1e4);
0453             residuals_["residual_rms_x_pos"]->Fill(me_x->getBinError(j, i) * 1e4);
0454             residuals_["residual_rms_y_pos"]->Fill(me_y->getBinError(j, i) * 1e4);
0455             DRnR_["NormRes_mean_x_pos"]->Fill(me2_x->getBinContent(j, i));
0456             DRnR_["NormRes_mean_y_pos"]->Fill(me2_y->getBinContent(j, i));
0457             DRnR_["DRnR_x_pos"]->Fill(me2_x->getBinError(j, i));
0458             DRnR_["DRnR_y_pos"]->Fill(me2_y->getBinError(j, i));
0459           }
0460           me2_x->setBinContent(j, i, me2_x->getBinError(j, i));
0461           me2_y->setBinContent(j, i, me2_y->getBinError(j, i));
0462           me2_x->setBinEntries(me2_x->getBin(j, i), 1);
0463           me2_y->setBinEntries(me2_y->getBin(j, i), 1);
0464         } else {
0465           me2_x->setBinContent(j, i, 0);
0466           me2_y->setBinContent(j, i, 0);
0467           me2_x->setBinEntries(me2_x->getBin(j, i), 0);
0468           me2_y->setBinEntries(me2_y->getBin(j, i), 0);
0469         }
0470       }
0471     }
0472   }
0473 }
0474 
0475 //define this as a plug-in
0476 DEFINE_FWK_MODULE(SiPixelPhase1ResidualsExtra);