Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:09:04

0001 #include "DQM/TrackerRemapper/interface/Phase1PixelROCMaps.h"
0002 
0003 using modBins = std::vector<std::pair<int, int>>;
0004 using rocBins = std::vector<std::tuple<int, int, int>>;
0005 
0006 // find detector coordinates for filling
0007 /*--------------------------------------------------------------------*/
0008 DetCoordinates Phase1PixelROCMaps::findDetCoordinates(const uint32_t& t_detid)
0009 /*--------------------------------------------------------------------*/
0010 {
0011   DetCoordinates coord;
0012 
0013   auto myDetId = DetId(t_detid);
0014   int subid = DetId(t_detid).subdetId();
0015 
0016   if (subid == PixelSubdetector::PixelBarrel) {
0017     coord.m_layer = m_trackerTopo.pxbLayer(myDetId);
0018     coord.m_s_ladder = this->signed_ladder(myDetId, true);
0019     coord.m_s_module = this->signed_module(myDetId, true);
0020 
0021     bool isFlipped = this->isBPixOuterLadder(myDetId, false);
0022     if ((coord.m_layer > 1 && coord.m_s_module < 0))
0023       isFlipped = !isFlipped;
0024 
0025     coord.m_isFlipped = isFlipped;
0026 
0027   }  // if it's barrel
0028   else if (subid == PixelSubdetector::PixelEndcap) {
0029     coord.m_ring = this->ring(myDetId, true);
0030     coord.m_s_blade = this->signed_blade(myDetId, true);
0031     coord.m_s_disk = this->signed_disk(myDetId, true);
0032     coord.m_panel = m_trackerTopo.pxfPanel(t_detid);
0033     coord.m_isFlipped = (coord.m_s_disk > 0) ? (coord.m_panel == 1) : (coord.m_panel == 2);
0034   }  // it it's endcap
0035   else {
0036     throw cms::Exception("LogicError") << "Unknown Pixel SubDet ID " << std::endl;
0037   }
0038 
0039   if (std::strcmp(m_option, kVerbose) == 0) {
0040     coord.printCoordinates();
0041   }
0042 
0043   return coord;
0044 }
0045 
0046 // overloaded method: mask entire module
0047 /*--------------------------------------------------------------------*/
0048 modBins Phase1PixelROCMaps::maskedBarrelRocsToBins(DetCoordinates coord)
0049 /*--------------------------------------------------------------------*/
0050 {
0051   modBins rocsToMask;
0052   int nlad = nlad_list[coord.m_layer - 1];
0053 
0054   int start_x = coord.m_s_module > 0 ? ((coord.m_s_module + 4) * 8) + 1 : ((4 - (std::abs(coord.m_s_module))) * 8) + 1;
0055   int start_y =
0056       coord.m_s_ladder > 0 ? ((coord.m_s_ladder + nlad) * 2) + 1 : ((nlad - (std::abs(coord.m_s_ladder))) * 2) + 1;
0057 
0058   int end_x = start_x + 7;
0059   int end_y = start_y + 1;
0060 
0061   for (int bin_x = 1; bin_x <= 72; bin_x++) {
0062     for (int bin_y = 1; bin_y <= (nlad * 4 + 2); bin_y++) {
0063       if (bin_x >= start_x && bin_x <= end_x && bin_y >= start_y && bin_y <= end_y) {
0064         rocsToMask.push_back(std::make_pair(bin_x, bin_y));
0065       }
0066     }
0067   }
0068   return rocsToMask;
0069 }
0070 
0071 // overloaded method: mask single ROCs
0072 /*--------------------------------------------------------------------*/
0073 rocBins Phase1PixelROCMaps::maskedBarrelRocsToBins(DetCoordinates coord, std::bitset<16> myRocs)
0074 /*--------------------------------------------------------------------*/
0075 {
0076   rocBins rocsToMask;
0077   int nlad = nlad_list[coord.m_layer - 1];
0078 
0079   int start_x = coord.m_s_module > 0 ? ((coord.m_s_module + 4) * 8) + 1 : ((4 - (std::abs(coord.m_s_module))) * 8) + 1;
0080   int start_y =
0081       coord.m_s_ladder > 0 ? ((coord.m_s_ladder + nlad) * 2) + 1 : ((nlad - (std::abs(coord.m_s_ladder))) * 2) + 1;
0082 
0083   int roc0_x = ((coord.m_layer == 1) || (coord.m_layer > 1 && coord.m_s_module > 0)) ? start_x + 7 : start_x;
0084   int roc0_y = start_y - 1;
0085 
0086   size_t idx = 0;
0087   while (idx < myRocs.size()) {
0088     if (myRocs.test(idx)) {
0089       //////////////////////////////////////////////////////////////////////////////////////
0090       //                                |                     //
0091       // In BPix Layer1 & module > 0 in L2,3,4  |   In BPix Layer 2,3,4 module < 0        //
0092       //                                        |                     //
0093       // ROCs are ordered in the following      |   ROCs are ordered in the following     //
0094       // fashion for unflipped modules          |   fashion for unflipped modules         //
0095       //                        |                         //
0096       // +---+---+---+---+---+---+---+---+      |   +---+---+---+---+---+---+---+---+     //
0097       // | 8 |9  |10 |11 |12 |13 |14 |15 |      |   |15 |14 |13 |12 |11 |10 | 9 | 8 |     //
0098       // +---+---+---+---+---+---+---+---+      |   +---+---+---+---+---+---+---+---+     //
0099       // | 7 | 6 | 5 | 4 | 3 | 2 | 1 | 0 |      |   | 0 | 1 | 2 | 3 | 4 | 5 | 6 | 7 |     //
0100       // +---+---+---+---+---+---+---+---+      |   +---+---+---+---+---+---+---+---+     //
0101       //                    |                     //
0102       // if the module is flipped the ordering  |   if the module is flipped the ordering //
0103       // is reveresed                           |   is reversed                           //
0104       //                    |                                         //
0105       // +---+---+---+---+---+---+---+---+      |   +---+---+---+---+---+---+---+---+     //
0106       // | 7 | 6 | 5 | 4 | 3 | 2 | 1 | 0 |      |   | 0 | 1 | 2 | 3 | 4 | 5 | 6 | 7 |     //
0107       // +---+---+---+---+---+---+---+---+      |   +---+---+---+---+---+---+---+---+     //
0108       // | 8 | 9 |10 |11 |12 |13 |14 |15 |      |   |15 |14 |13 |12 |11 |10 | 9 | 8 |     //
0109       // +---+---+---+---+---+---+---+---+      |   +---+---+---+---+---+---+---+---+     //
0110       //////////////////////////////////////////////////////////////////////////////////////
0111 
0112       int roc_x(0), roc_y(0);
0113 
0114       if ((coord.m_layer == 1) || (coord.m_layer > 1 && coord.m_s_module > 0)) {
0115         if (!coord.m_isFlipped) {
0116           roc_x = idx < 8 ? roc0_x - idx : (start_x - 8) + idx;
0117           roc_y = idx < 8 ? roc0_y + 1 : roc0_y + 2;
0118         } else {
0119           roc_x = idx < 8 ? roc0_x - idx : (start_x - 8) + idx;
0120           roc_y = idx < 8 ? roc0_y + 2 : roc0_y + 1;
0121         }
0122       } else {
0123         if (!coord.m_isFlipped) {
0124           roc_x = idx < 8 ? roc0_x + idx : (roc0_x + 7) - (idx - 8);
0125           roc_y = idx < 8 ? roc0_y + 1 : roc0_y + 2;
0126         } else {
0127           roc_x = idx < 8 ? roc0_x + idx : (roc0_x + 7) - (idx - 8);
0128           roc_y = idx < 8 ? roc0_y + 2 : roc0_y + 1;
0129         }
0130       }
0131       rocsToMask.push_back(std::make_tuple(roc_x, roc_y, idx));
0132     }
0133     ++idx;
0134   }
0135   return rocsToMask;
0136 }
0137 
0138 // overloaded method: mask entire module
0139 /*--------------------------------------------------------------------*/
0140 modBins Phase1PixelROCMaps::maskedForwardRocsToBins(DetCoordinates coord)
0141 /*--------------------------------------------------------------------*/
0142 {
0143   modBins rocsToMask;
0144   int nybins = nybins_list[coord.m_ring - 1];
0145 
0146   int start_x = coord.m_s_disk > 0 ? ((coord.m_s_disk + 3) * 8) + 1 : ((3 - (std::abs(coord.m_s_disk))) * 8) + 1;
0147   int start_y = coord.m_s_blade > 0 ? (nybins / 2) + (coord.m_s_blade * 4) - (coord.m_panel * 2) + 3
0148                                     : ((nybins / 2) - (std::abs(coord.m_s_blade) * 4) - coord.m_panel * 2) + 3;
0149 
0150   int end_x = start_x + 7;
0151   int end_y = start_y + 1;
0152 
0153   for (int bin_x = 1; bin_x <= 56; bin_x++) {
0154     for (int bin_y = 1; bin_y <= nybins; bin_y++) {
0155       if (bin_x >= start_x && bin_x <= end_x && bin_y >= start_y && bin_y <= end_y) {
0156         rocsToMask.push_back(std::make_pair(bin_x, bin_y));
0157       }
0158     }
0159   }
0160   return rocsToMask;
0161 }
0162 
0163 // overloaded method: mask single ROCs
0164 /*--------------------------------------------------------------------*/
0165 rocBins Phase1PixelROCMaps::maskedForwardRocsToBins(DetCoordinates coord, std::bitset<16> myRocs)
0166 /*--------------------------------------------------------------------*/
0167 {
0168   rocBins rocsToMask;
0169   int nybins = nybins_list[coord.m_ring - 1];
0170 
0171   int start_x = coord.m_s_disk > 0 ? ((coord.m_s_disk + 3) * 8) + 1 : ((3 - (std::abs(coord.m_s_disk))) * 8) + 1;
0172   int start_y = coord.m_s_blade > 0 ? (nybins / 2) + (coord.m_s_blade * 4) - (coord.m_panel * 2) + 3
0173                                     : ((nybins / 2) - (std::abs(coord.m_s_blade) * 4) - coord.m_panel * 2) + 3;
0174 
0175   int roc0_x = coord.m_s_disk > 0 ? start_x + 7 : start_x;
0176   int roc0_y = start_y - 1;
0177 
0178   size_t idx = 0;
0179   while (idx < myRocs.size()) {
0180     if (myRocs.test(idx)) {
0181       int roc_x(0), roc_y(0);
0182 
0183       //////////////////////////////////////////////////////////////////////////////////////
0184       //                                |                     //
0185       // In FPix + (Disk 1,2,3)                 |   In FPix - (Disk -1,-2,-3)             //
0186       //                                        |                     //
0187       // ROCs are ordered in the following      |   ROCs are ordered in the following     //
0188       // fashion for unflipped modules          |   fashion for unflipped modules         //
0189       //                    |                         //
0190       // +---+---+---+---+---+---+---+---+      |   +---+---+---+---+---+---+---+---+     //
0191       // | 8 |9  |10 |11 |12 |13 |14 |15 |      |   |15 |14 |13 |12 |11 |10 | 9 | 8 |     //
0192       // +---+---+---+---+---+---+---+---+      |   +---+---+---+---+---+---+---+---+     //
0193       // | 7 | 6 | 5 | 4 | 3 | 2 | 1 | 0 |      |   | 0 | 1 | 2 | 3 | 4 | 5 | 6 | 7 |     //
0194       // +---+---+---+---+---+---+---+---+      |   +---+---+---+---+---+---+---+---+     //
0195       //                    |                     //
0196       // if the module is flipped the ordering  |   if the module is flipped the ordering //
0197       // is reveresed                           |   is reversed                           //
0198       //                    |                                         //
0199       // +---+---+---+---+---+---+---+---+      |   +---+---+---+---+---+---+---+---+     //
0200       // | 7 | 6 | 5 | 4 | 3 | 2 | 1 | 0 |      |   | 0 | 1 | 2 | 3 | 4 | 5 | 6 | 7 |     //
0201       // +---+---+---+---+---+---+---+---+      |   +---+---+---+---+---+---+---+---+     //
0202       // | 8 | 9 |10 |11 |12 |13 |14 |15 |      |   |15 |14 |13 |12 |11 |10 | 9 | 8 |     //
0203       // +---+---+---+---+---+---+---+---+      |   +---+---+---+---+---+---+---+---+     //
0204       //////////////////////////////////////////////////////////////////////////////////////
0205 
0206       if (coord.m_s_disk > 0) {
0207         if (!coord.m_isFlipped) {
0208           roc_x = idx < 8 ? roc0_x - idx : (start_x - 8) + idx;
0209           roc_y = idx < 8 ? roc0_y + 1 : roc0_y + 2;
0210         } else {
0211           roc_x = idx < 8 ? roc0_x - idx : (start_x - 8) + idx;
0212           roc_y = idx < 8 ? roc0_y + 2 : roc0_y + 1;
0213         }
0214       } else {
0215         if (!coord.m_isFlipped) {
0216           roc_x = idx < 8 ? roc0_x + idx : (roc0_x + 7) - (idx - 8);
0217           roc_y = idx < 8 ? roc0_y + 1 : roc0_y + 2;
0218         } else {
0219           roc_x = idx < 8 ? roc0_x + idx : (roc0_x + 7) - (idx - 8);
0220           roc_y = idx < 8 ? roc0_y + 2 : roc0_y + 1;
0221         }
0222       }
0223 
0224       rocsToMask.push_back(std::make_tuple(roc_x, roc_y, idx));
0225     }
0226     ++idx;
0227   }
0228   return rocsToMask;
0229 }
0230 
0231 /*--------------------------------------------------------------------*/
0232 void Phase1PixelROCMaps::fillWholeModule(const uint32_t& detid, double value)
0233 /*--------------------------------------------------------------------*/
0234 {
0235   auto coord = findDetCoordinates(detid);
0236   auto rocsToMark = coord.isBarrel() ? this->maskedBarrelRocsToBins(coord) : this->maskedForwardRocsToBins(coord);
0237 
0238   if (coord.isBarrel()) {
0239     for (const auto& bin : rocsToMark) {
0240       double x = h_bpix_maps[coord.m_layer - 1]->GetXaxis()->GetBinCenter(bin.first);
0241       double y = h_bpix_maps[coord.m_layer - 1]->GetYaxis()->GetBinCenter(bin.second);
0242       h_bpix_maps[coord.m_layer - 1]->Fill(x, y, value);
0243     }
0244   } else {
0245     for (const auto& bin : rocsToMark) {
0246       double x = h_fpix_maps[coord.m_ring - 1]->GetXaxis()->GetBinCenter(bin.first);
0247       double y = h_fpix_maps[coord.m_ring - 1]->GetYaxis()->GetBinCenter(bin.second);
0248       h_fpix_maps[coord.m_ring - 1]->Fill(x, y, value);
0249     }
0250   }
0251   return;
0252 }
0253 
0254 /*--------------------------------------------------------------------*/
0255 void Phase1PixelROCMaps::fillSelectedRocs(const uint32_t& detid, const std::bitset<16>& theROCs, double value)
0256 /*--------------------------------------------------------------------*/
0257 {
0258   auto coord = findDetCoordinates(detid);
0259   auto rocsToMark =
0260       coord.isBarrel() ? this->maskedBarrelRocsToBins(coord, theROCs) : this->maskedForwardRocsToBins(coord, theROCs);
0261 
0262   if (coord.isBarrel()) {
0263     for (const auto& bin : rocsToMark) {
0264       double x = h_bpix_maps[coord.m_layer - 1]->GetXaxis()->GetBinCenter(std::get<0>(bin));
0265       double y = h_bpix_maps[coord.m_layer - 1]->GetYaxis()->GetBinCenter(std::get<1>(bin));
0266       h_bpix_maps[coord.m_layer - 1]->Fill(x, y, value);
0267     }
0268   } else {
0269     for (const auto& bin : rocsToMark) {
0270       double x = h_fpix_maps[coord.m_ring - 1]->GetXaxis()->GetBinCenter(std::get<0>(bin));
0271       double y = h_fpix_maps[coord.m_ring - 1]->GetYaxis()->GetBinCenter(std::get<1>(bin));
0272       h_fpix_maps[coord.m_ring - 1]->Fill(x, y, value);
0273     }
0274   }
0275 
0276   return;
0277 }
0278 
0279 /*--------------------------------------------------------------------*/
0280 void PixelROCMapHelper::draw_line(
0281     double x1, double x2, double y1, double y2, int width = 2, int style = 1, int color = 1)
0282 /*--------------------------------------------------------------------*/
0283 {
0284   TLine* l = new TLine(x1, y1, x2, y2);
0285   l->SetBit(kCanDelete);
0286   l->SetLineWidth(width);
0287   l->SetLineStyle(style);
0288   l->SetLineColor(color);
0289   l->Draw();
0290 }
0291 
0292 /*--------------------------------------------------------------------*/
0293 void PixelROCMapHelper::dress_plot(TPad*& canv,
0294                                    TH2* h,
0295                                    int lay,
0296                                    int ring = 0,
0297                                    int phase = 0,
0298                                    bool standard_palette = true,
0299                                    bool half_shift = true,
0300                                    bool mark_zero = true)
0301 /*--------------------------------------------------------------------*/
0302 {
0303   std::string s_title;
0304   const auto zAxisTitle = fmt::sprintf("%s", h->GetZaxis()->GetTitle());
0305 
0306   if (lay > 0) {
0307     canv->cd(lay);
0308     canv->cd(lay)->SetTopMargin(0.05);
0309     canv->cd(lay)->SetBottomMargin(0.07);
0310     canv->cd(lay)->SetLeftMargin(0.1);
0311     if (!zAxisTitle.empty()) {
0312       h->GetZaxis()->SetTitleOffset(1.3);
0313       h->GetZaxis()->CenterTitle(true);
0314       canv->cd(lay)->SetRightMargin(0.14);
0315     } else {
0316       canv->cd(lay)->SetRightMargin(0.11);
0317     }
0318     s_title = "Barrel Pixel Layer " + std::to_string(lay);
0319   } else {
0320     canv->cd(ring);
0321     canv->cd(ring)->SetTopMargin(0.05);
0322     canv->cd(ring)->SetBottomMargin(0.07);
0323     canv->cd(ring)->SetLeftMargin(0.1);
0324     if (!zAxisTitle.empty()) {
0325       h->GetZaxis()->SetTitleOffset(1.3);
0326       h->GetZaxis()->CenterTitle(true);
0327       canv->cd(ring)->SetRightMargin(0.14);
0328     } else {
0329       canv->cd(ring)->SetRightMargin(0.11);
0330     }
0331     if (ring > 4) {
0332       ring = ring - 4;
0333     }
0334     s_title = "Forward Pixel Ring " + std::to_string(ring);
0335   }
0336 
0337   if (standard_palette) {
0338     gStyle->SetPalette(1);
0339   } else {
0340     /*
0341     const Int_t NRGBs = 5;
0342     const Int_t NCont = 255;
0343 
0344     Double_t stops[NRGBs] = {0.00, 0.34, 0.61, 0.84, 1.00};
0345     Double_t red[NRGBs] = {0.00, 0.00, 0.87, 1.00, 0.51};
0346     Double_t green[NRGBs] = {0.00, 0.81, 1.00, 0.20, 0.00};
0347     Double_t blue[NRGBs] = {0.51, 1.00, 0.12, 0.00, 0.00};
0348     TColor::CreateGradientColorTable(NRGBs, stops, red, green, blue, NCont);
0349     gStyle->SetNumberContours(NCont);
0350     */
0351 
0352     // this is the fine gradient palette (blue to red)
0353     double max = h->GetMaximum();
0354     double min = h->GetMinimum();
0355     double val_white = 0.;
0356     double per_white = (max != min) ? ((val_white - min) / (max - min)) : 0.5;
0357 
0358     const int Number = 3;
0359     double Red[Number] = {0., 1., 1.};
0360     double Green[Number] = {0., 1., 0.};
0361     double Blue[Number] = {1., 1., 0.};
0362     double Stops[Number] = {0., per_white, 1.};
0363     int nb = 256;
0364     h->SetContour(nb);
0365     TColor::CreateGradientColorTable(Number, Stops, Red, Green, Blue, nb);
0366     // if max == min impose the range to be the same as it was a real diff
0367     if (max == min)
0368       h->GetZaxis()->SetRangeUser(-1., 1.);
0369   }
0370 
0371   h->SetMarkerSize(0.7);
0372   h->Draw("colz1");
0373 
0374   auto ltx = TLatex();
0375   ltx.SetTextFont(62);
0376   ltx.SetTextColor(1);
0377   ltx.SetTextSize(0.06);
0378   ltx.SetTextAlign(31);
0379   ltx.DrawLatexNDC(1 - gPad->GetRightMargin(), 1 - gPad->GetTopMargin() + 0.01, (s_title).c_str());
0380 
0381   // Draw Lines around modules
0382   if (lay > 0) {
0383     std::vector<std::vector<int>> nladder = {{10, 16, 22}, {6, 14, 22, 32}};
0384     int nlad = nladder[phase][lay - 1];
0385     for (int xsign = -1; xsign <= 1; xsign += 2)
0386       for (int ysign = -1; ysign <= 1; ysign += 2) {
0387         float xlow = xsign * (half_shift * 0.5);
0388         float xhigh = xsign * (half_shift * 0.5 + 4);
0389         float ylow = ysign * (half_shift * 0.5 + (phase == 0) * 0.5);
0390         float yhigh = ysign * (half_shift * 0.5 - (phase == 0) * 0.5 + nlad);
0391         // Outside box
0392         PixelROCMapHelper::draw_line(xlow, xhigh, ylow, ylow, 1);    // bottom
0393         PixelROCMapHelper::draw_line(xlow, xhigh, yhigh, yhigh, 1);  // top
0394         PixelROCMapHelper::draw_line(xlow, xlow, ylow, yhigh, 1);    // left
0395         PixelROCMapHelper::draw_line(xhigh, xhigh, ylow, yhigh, 1);  // right
0396         // Inner Horizontal lines
0397         for (int lad = 1; lad < nlad; ++lad) {
0398           float y = ysign * (lad + half_shift * 0.5);
0399           PixelROCMapHelper::draw_line(xlow, xhigh, y, y, 1);
0400         }
0401         for (int lad = 1; lad <= nlad; ++lad)
0402           if (!(phase == 0 && (lad == 1 || lad == nlad))) {
0403             float y = ysign * (lad + half_shift * 0.5 - 0.5);
0404             PixelROCMapHelper::draw_line(xlow, xhigh, y, y, 1, 3);
0405           }
0406         // Inner Vertical lines
0407         for (int mod = 1; mod < 4; ++mod) {
0408           float x = xsign * (mod + half_shift * 0.5);
0409           PixelROCMapHelper::draw_line(x, x, ylow, yhigh, 1);
0410         }
0411         // Make a BOX around ROC 0
0412         // Phase 0 - ladder +1 is always non-flipped
0413         // Phase 1 - ladder +1 is always     flipped
0414         if (mark_zero) {
0415           for (int mod = 1; mod <= 4; ++mod)
0416             for (int lad = 1; lad <= nlad; ++lad) {
0417               bool flipped = ysign == 1 ? lad % 2 == 0 : lad % 2 == 1;
0418               if (phase == 1)
0419                 flipped = !flipped;
0420               int roc0_orientation = flipped ? -1 : 1;
0421               if (xsign == -1)
0422                 roc0_orientation *= -1;
0423               if (ysign == -1)
0424                 roc0_orientation *= -1;
0425               float x1 = xsign * (mod + half_shift * 0.5);
0426               float x2 = xsign * (mod + half_shift * 0.5 - 1. / 8);
0427               float y1 = ysign * (lad + half_shift * 0.5 - 0.5);
0428               float y2 = ysign * (lad + half_shift * 0.5 - 0.5 + roc0_orientation * 1. / 2);
0429               if (!(phase == 0 && (lad == 1 || lad == nlad) && xsign == -1)) {
0430                 if (lay == 1 && xsign <= -1) {
0431                   float x1 = xsign * ((mod - 1) + half_shift * 0.5);
0432                   float x2 = xsign * ((mod - 1) + half_shift * 0.5 + 1. / 8);
0433                   float y1 = ysign * (lad + half_shift * 0.5 - 0.5 + roc0_orientation);
0434                   float y2 = ysign * (lad + half_shift * 0.5 - 0.5 + roc0_orientation * 3. / 2);
0435                   PixelROCMapHelper::draw_line(x1, x2, y1, y1, 1);
0436                   PixelROCMapHelper::draw_line(x2, x2, y1, y2, 1);
0437                 } else {
0438                   PixelROCMapHelper::draw_line(x1, x2, y1, y1, 1);
0439                   //PixelROCMapHelper::draw_line(x1, x2, y2, y2, 1);
0440                   //PixelROCMapHelper::draw_line(x1, x1, y1, y2, 1);
0441                   PixelROCMapHelper::draw_line(x2, x2, y1, y2, 1);
0442                 }
0443               }
0444             }
0445         }
0446       }
0447   } else {
0448     // FPIX
0449     for (int dsk = 1, ndsk = 2 + (phase == 1); dsk <= ndsk; ++dsk) {
0450       for (int xsign = -1; xsign <= 1; xsign += 2)
0451         for (int ysign = -1; ysign <= 1; ysign += 2) {
0452           if (phase == 0) {
0453             int first_roc = 3, nbin = 16;
0454             for (int bld = 1, nbld = 12; bld <= nbld; ++bld) {
0455               // Horizontal lines
0456               for (int plq = 1, nplq = 7; plq <= nplq; ++plq) {
0457                 float xlow =
0458                     xsign * (half_shift * 0.5 + dsk - 1 + (first_roc - 3 + 2 * plq + (plq == 1)) / (float)nbin);
0459                 float xhigh =
0460                     xsign * (half_shift * 0.5 + dsk - 1 + (first_roc - 3 + 2 * (plq + 1) - (plq == 7)) / (float)nbin);
0461                 float ylow = ysign * (half_shift * 0.5 + (bld - 0.5) - (2 + plq / 2) * 0.1);
0462                 float yhigh = ysign * (half_shift * 0.5 + (bld - 0.5) + (2 + plq / 2) * 0.1);
0463                 PixelROCMapHelper::draw_line(xlow, xhigh, ylow, ylow, 1);    // bottom
0464                 PixelROCMapHelper::draw_line(xlow, xhigh, yhigh, yhigh, 1);  // top
0465               }
0466               // Vertical lines
0467               for (int plq = 1, nplq = 7 + 1; plq <= nplq; ++plq) {
0468                 float x = xsign * (half_shift * 0.5 + dsk - 1 +
0469                                    (first_roc - 3 + 2 * plq + (plq == 1) - (plq == 8)) / (float)nbin);
0470                 float ylow = ysign * (half_shift * 0.5 + (bld - 0.5) - (2 + (plq - (plq == 8)) / 2) * 0.1);
0471                 float yhigh = ysign * (half_shift * 0.5 + (bld - 0.5) + (2 + (plq - (plq == 8)) / 2) * 0.1);
0472                 PixelROCMapHelper::draw_line(x, x, ylow, yhigh, 1);
0473               }
0474               // Panel 2 has dashed mid-plane
0475               for (int plq = 2, nplq = 6; plq <= nplq; ++plq)
0476                 if (plq % 2 == 0) {
0477                   float x = xsign * (half_shift * 0.5 + dsk - 1 +
0478                                      (first_roc - 3 + 2 * plq + (plq == 1) - (plq == 8) + 1) / (float)nbin);
0479                   float ylow = ysign * (half_shift * 0.5 + (bld - 0.5) - (2 + (plq - (plq == 8)) / 2) * 0.1);
0480                   float yhigh = ysign * (half_shift * 0.5 + (bld - 0.5) + (2 + (plq - (plq == 8)) / 2) * 0.1);
0481                   PixelROCMapHelper::draw_line(x, x, ylow, yhigh, 1, 2);
0482                 }
0483               // Make a BOX around ROC 0
0484               for (int plq = 1, nplq = 7; plq <= nplq; ++plq) {
0485                 float x1 = xsign * (half_shift * 0.5 + dsk - 1 + (first_roc - 3 + 2 * plq + (plq == 1)) / (float)nbin);
0486                 float x2 =
0487                     xsign * (half_shift * 0.5 + dsk - 1 + (first_roc - 3 + 2 * plq + (plq == 1) + 1) / (float)nbin);
0488                 int sign = xsign * ysign * ((plq % 2) ? 1 : -1);
0489                 float y1 = ysign * (half_shift * 0.5 + (bld - 0.5) + sign * (2 + plq / 2) * 0.1);
0490                 float y2 = ysign * (half_shift * 0.5 + (bld - 0.5) + sign * (plq / 2) * 0.1);
0491                 //PixelROCMapHelper::draw_line(x1, x2, y1, y1, 1);
0492                 PixelROCMapHelper::draw_line(x1, x2, y2, y2, 1);
0493                 //PixelROCMapHelper::draw_line(x1, x1, y1, y2, 1);
0494                 PixelROCMapHelper::draw_line(x2, x2, y1, y2, 1);
0495               }
0496             }
0497           } else if (phase == 1) {
0498             if (ring == 0) {  // both
0499               for (int ring = 1; ring <= 2; ++ring)
0500                 for (int bld = 1, nbld = 5 + ring * 6; bld <= nbld; ++bld) {
0501                   float scale = (ring == 1) ? 1.5 : 1;
0502                   Color_t p1_color = 1, p2_color = 1;
0503                   // Horizontal lines
0504                   // Panel 2 has dashed mid-plane
0505                   float x1 = xsign * (half_shift * 0.5 + dsk - 1 + (ring - 1) * 0.5);
0506                   float x2 = xsign * (half_shift * 0.5 + dsk - 1 + ring * 0.5);
0507                   int sign = ysign;
0508                   float y1 = ysign * (half_shift * 0.5 - 0.5 + scale * bld + sign * 0.5);
0509                   //float yp1_mid = ysign * (half_shift*0.5 - 0.5 + scale*bld + sign*0.25);
0510                   float y2 = ysign * (half_shift * 0.5 - 0.5 + scale * bld);
0511                   float yp2_mid = ysign * (half_shift * 0.5 - 0.5 + scale * bld - sign * 0.25);
0512                   float y3 = ysign * (half_shift * 0.5 - 0.5 + scale * bld - sign * 0.5);
0513                   PixelROCMapHelper::draw_line(x1, x2, y1, y1, 1, 1, p1_color);
0514                   //PixelROCMapHelper::draw_line(x1, x2, yp1_mid, yp1_mid, 1, 3);
0515                   PixelROCMapHelper::draw_line(x1, x2, y2, y2, 1, 1, p1_color);
0516                   PixelROCMapHelper::draw_line(x1, x2, yp2_mid, yp2_mid, 1, 2);
0517                   PixelROCMapHelper::draw_line(x1, x2, y3, y3, 1, 1, p2_color);
0518                   // Vertical lines
0519                   float x = xsign * (half_shift * 0.5 + dsk - 1 + (ring - 1) * 0.5);
0520                   PixelROCMapHelper::draw_line(x, x, y1, y2, 1, 1, p1_color);
0521                   PixelROCMapHelper::draw_line(x, x, y2, y3, 1, 1, p2_color);
0522                   if (ring == 2) {
0523                     //PixelROCMapHelper::draw_line(x,  x,  y2,  y3, 1, 1, p1_color);
0524                     x = xsign * (half_shift * 0.5 + dsk);
0525                     PixelROCMapHelper::draw_line(x, x, y1, y2, 1, 1, p1_color);
0526                     PixelROCMapHelper::draw_line(x, x, y2, y3, 1, 1, p2_color);
0527                   }
0528                   // Make a BOX around ROC 0
0529                   x1 = xsign * (half_shift * 0.5 + dsk - 1 + ring * 0.5 - 1 / 16.);
0530                   x2 = xsign * (half_shift * 0.5 + dsk - 1 + ring * 0.5);
0531                   float y1_p1 = ysign * (half_shift * 0.5 - 0.5 + scale * bld + sign * 0.25);
0532                   float y2_p1 = ysign * (half_shift * 0.5 - 0.5 + scale * bld + sign * 0.25 + xsign * ysign * 0.25);
0533                   PixelROCMapHelper::draw_line(x1, x2, y1_p1, y1_p1, 1);
0534                   //PixelROCMapHelper::draw_line(x1, x2, y2_p1, y2_p1, 1);
0535                   PixelROCMapHelper::draw_line(x1, x1, y1_p1, y2_p1, 1);
0536                   //PixelROCMapHelper::draw_line(x2, x2, y1_p1, y2_p1, 1);
0537                   float y1_p2 = ysign * (half_shift * 0.5 - 0.5 + scale * bld - sign * 0.25);
0538                   float y2_p2 = ysign * (half_shift * 0.5 - 0.5 + scale * bld - sign * 0.25 - xsign * ysign * 0.25);
0539                   PixelROCMapHelper::draw_line(x1, x2, y1_p2, y1_p2, 1);
0540                   //PixelROCMapHelper::draw_line(x1, x2, y2_p2, y2_p2, 1);
0541                   PixelROCMapHelper::draw_line(x1, x1, y1_p2, y2_p2, 1);
0542                   //PixelROCMapHelper::draw_line(x2, x2, y1_p2, y2_p2, 1);
0543                 }
0544             } else {  // only one ring, 1 or 2
0545               for (int bld = 1, nbld = 5 + ring * 6; bld <= nbld; ++bld) {
0546                 Color_t p1_color = 1, p2_color = 1;
0547                 // Horizontal lines
0548                 // Panel 2 has dashed mid-plane
0549                 float x1 = xsign * (half_shift * 0.5 + dsk - 1);
0550                 float x2 = xsign * (half_shift * 0.5 + dsk);
0551                 int sign = ysign;
0552                 float y1 = ysign * (half_shift * 0.5 - 0.5 + bld + sign * 0.5);
0553                 //float yp1_mid = ysign * (half_shift*0.5 - 0.5 + bld + sign*0.25);
0554                 float y2 = ysign * (half_shift * 0.5 - 0.5 + bld);
0555                 float yp2_mid = ysign * (half_shift * 0.5 - 0.5 + bld - sign * 0.25);
0556                 float y3 = ysign * (half_shift * 0.5 - 0.5 + bld - sign * 0.5);
0557                 PixelROCMapHelper::draw_line(x1, x2, y1, y1, 1, 1, p1_color);
0558                 //PixelROCMapHelper::draw_line(x1, x2, yp1_mid, yp1_mid, 1, 3);
0559                 PixelROCMapHelper::draw_line(x1, x2, y2, y2, 1, 1, p1_color);
0560                 PixelROCMapHelper::draw_line(x1, x2, yp2_mid, yp2_mid, 1, 2);
0561                 PixelROCMapHelper::draw_line(x1, x2, y3, y3, 1, 1, p2_color);
0562                 // Vertical lines
0563                 float x = xsign * (half_shift * 0.5 + dsk - 1);
0564                 PixelROCMapHelper::draw_line(x, x, y1, y2, 1, 1, p1_color);
0565                 PixelROCMapHelper::draw_line(x, x, y2, y3, 1, 1, p2_color);
0566                 if (ring == 2) {
0567                   //PixelROCMapHelper::draw_line(x,  x,  y2,  y3, 1, 1, p1_color);
0568                   x = xsign * (half_shift * 0.5 + dsk);
0569                   PixelROCMapHelper::draw_line(x, x, y1, y2, 1, 1, p1_color);
0570                   PixelROCMapHelper::draw_line(x, x, y2, y3, 1, 1, p2_color);
0571                 }
0572                 // Make a BOX around ROC 0
0573                 x1 = xsign * (half_shift * 0.5 + dsk - 1 / 8.);
0574                 x2 = xsign * (half_shift * 0.5 + dsk);
0575                 float y1_p1 = ysign * (half_shift * 0.5 - 0.5 + bld + sign * 0.25);
0576                 float y2_p1 = ysign * (half_shift * 0.5 - 0.5 + bld + sign * 0.25 + xsign * ysign * 0.25);
0577                 PixelROCMapHelper::draw_line(x1, x2, y1_p1, y1_p1, 1);
0578                 //PixelROCMapHelper::draw_line(x1, x2, y2_p1, y2_p1, 1);
0579                 PixelROCMapHelper::draw_line(x1, x1, y1_p1, y2_p1, 1);
0580                 //PixelROCMapHelper::draw_line(x2, x2, y1_p1, y2_p1, 1);
0581                 float y1_p2 = ysign * (half_shift * 0.5 - 0.5 + bld - sign * 0.25);
0582                 float y2_p2 = ysign * (half_shift * 0.5 - 0.5 + bld - sign * 0.25 - xsign * ysign * 0.25);
0583                 PixelROCMapHelper::draw_line(x1, x2, y1_p2, y1_p2, 1);
0584                 //PixelROCMapHelper::draw_line(x1, x2, y2_p2, y2_p2, 1);
0585                 PixelROCMapHelper::draw_line(x1, x1, y1_p2, y2_p2, 1);
0586                 //PixelROCMapHelper::draw_line(x2, x2, y1_p2, y2_p2, 1);
0587               }
0588             }
0589           }
0590         }
0591     }
0592     // Special shifted "rebin" for Phase 0
0593     // Y axis should always have at least half-roc granularity because
0594     // there are half-ROC size shifts implemented in the coordinates
0595     // To remove this and show full ROC granularity
0596     // We merge bin contents in each pair of bins corresponding to one ROC
0597     // TODO: make sure this works for Profiles
0598     if (phase == 0 && h->GetNbinsY() == 250 && h->GetNbinsX() == 80) {
0599       int nentries = h->GetEntries();
0600       for (int binx = 1; binx <= 80; ++binx) {
0601         double sum = 0;
0602         for (int biny = 1; biny <= 250; ++biny) {
0603           bool odd_nrocy = (binx - 1 < 40) != (((binx - 1) / 4) % 2);
0604           if (biny % 2 == odd_nrocy)
0605             sum += h->GetBinContent(binx, biny);
0606           else {
0607             sum += h->GetBinContent(binx, biny);
0608             if (sum) {
0609               h->SetBinContent(binx, biny, sum);
0610               h->SetBinContent(binx, biny - 1, sum);
0611             }
0612             sum = 0;
0613           }
0614         }
0615       }
0616       h->SetEntries(nentries);
0617     }
0618   }
0619 }
0620 
0621 /*--------------------------------------------------------------------*/
0622 void Phase1PixelROCMaps::drawBarrelMaps(TCanvas& canvas, const std::string& text)
0623 /*--------------------------------------------------------------------*/
0624 {
0625   canvas.cd();
0626   canvas.Modified();
0627 
0628   auto topPad = new TPad("pad1", "upper pad", 0.005, 0.96, 0.995, 0.995);
0629   topPad->Draw();
0630   topPad->cd();
0631   auto ltx = TLatex();
0632   ltx.SetTextFont(62);
0633 
0634   std::size_t found = text.find("Delta");
0635   if (found != std::string::npos) {
0636     ltx.SetTextSize(0.7);
0637   } else {
0638     ltx.SetTextSize(1.);
0639   }
0640   ltx.DrawLatexNDC(0.02, 0.3, text.c_str());
0641 
0642   canvas.cd();
0643   auto bottomPad = new TPad("pad2", "lower pad", 0.005, 0.005, 0.995, 0.955);
0644   bottomPad->Draw();
0645   bottomPad->cd();
0646   bottomPad->Divide(2, 2);
0647   for (unsigned int lay = 1; lay <= n_layers; lay++) {
0648     h_bpix_maps[lay - 1]->SetStats(false);
0649     PixelROCMapHelper::dress_plot(bottomPad, h_bpix_maps[lay - 1].get(), lay, 0, 1, found == std::string::npos);
0650   }
0651 }
0652 
0653 /*--------------------------------------------------------------------*/
0654 void Phase1PixelROCMaps::drawForwardMaps(TCanvas& canvas, const std::string& text)
0655 /*--------------------------------------------------------------------*/
0656 {
0657   canvas.cd();
0658   canvas.Modified();
0659 
0660   auto topPad = new TPad("pad1", "upper pad", 0.005, 0.94, 0.995, 0.995);
0661   topPad->Draw();
0662   topPad->cd();
0663   auto ltx = TLatex();
0664   ltx.SetTextFont(62);
0665 
0666   std::size_t found = text.find("Delta");
0667   if (found != std::string::npos) {
0668     ltx.SetTextSize(0.7);
0669   } else {
0670     ltx.SetTextSize(1.);
0671   }
0672   ltx.DrawLatexNDC(0.02, 0.3, text.c_str());
0673 
0674   canvas.cd();
0675   auto bottomPad = new TPad("pad2", "lower pad", 0.005, 0.005, 0.995, 0.935);
0676   bottomPad->Draw();
0677   bottomPad->cd();
0678   bottomPad->Divide(2, 1);
0679   for (unsigned int ring = 1; ring <= n_rings; ring++) {
0680     h_fpix_maps[ring - 1]->SetStats(false);
0681     PixelROCMapHelper::dress_plot(bottomPad, h_fpix_maps[ring - 1].get(), 0, ring, 1, found == std::string::npos);
0682   }
0683 }
0684 
0685 /*--------------------------------------------------------------------*/
0686 void Phase1PixelROCMaps::drawMaps(TCanvas& canvas, const std::string& text)
0687 /*--------------------------------------------------------------------*/
0688 {
0689   canvas.cd();
0690   canvas.Modified();
0691 
0692   auto topPad = new TPad("pad1", "upper pad", 0.005, 0.97, 0.995, 0.995);
0693   topPad->Draw();
0694   topPad->cd();
0695   auto ltx = TLatex();
0696   ltx.SetTextFont(62);
0697 
0698   std::size_t found = text.find("Delta");
0699   if (found != std::string::npos) {
0700     ltx.SetTextSize(0.7);
0701   } else {
0702     ltx.SetTextSize(1.);
0703   }
0704   ltx.DrawLatexNDC(0.02, 0.2, text.c_str());
0705 
0706   canvas.cd();
0707   auto bottomPad = new TPad("pad2", "lower pad", 0.005, 0.005, 0.995, 0.97);
0708   bottomPad->Draw();
0709   bottomPad->cd();
0710   bottomPad->Divide(2, 3);
0711 
0712   // dress the plots
0713   for (unsigned int lay = 1; lay <= n_layers; lay++) {
0714     h_bpix_maps[lay - 1]->SetStats(false);
0715     PixelROCMapHelper::dress_plot(bottomPad, h_bpix_maps[lay - 1].get(), lay, 0, 1, found == std::string::npos);
0716   }
0717 
0718   bottomPad->Update();
0719   bottomPad->Modified();
0720   bottomPad->cd();
0721 
0722   for (unsigned int ring = 1; ring <= n_rings; ring++) {
0723     h_fpix_maps[ring - 1]->SetStats(false);
0724     PixelROCMapHelper::dress_plot(
0725         bottomPad, h_fpix_maps[ring - 1].get(), 0, n_layers + ring, 1, found == std::string::npos);
0726   }
0727 }