Back to home page

sPhenix code displayed by LXR

 
 

    


File indexing completed on 2025-12-18 09:17:48

0001 #include "InttCalib.h"
0002 
0003 #include <intt/InttMapping.h>
0004 
0005 #include <cdbobjects/CDBTTree.h>
0006 
0007 #include <ffarawobjects/Gl1Packet.h>
0008 #include <ffarawobjects/InttRawHit.h>
0009 #include <ffarawobjects/InttRawHitContainer.h>
0010 
0011 #include <qautils/QAHistManagerDef.h>
0012 
0013 #include <fun4all/Fun4AllHistoManager.h>
0014 #include <fun4all/Fun4AllReturnCodes.h>
0015 
0016 #include <phool/getClass.h>
0017 #include <phool/phool.h>
0018 
0019 #include <TCanvas.h>
0020 #include <TF1.h>
0021 #include <TFile.h>
0022 #include <TH1.h>
0023 #include <TLine.h>
0024 #include <TPolyLine.h>
0025 #include <TStyle.h>
0026 #include <TText.h>
0027 #include <TTree.h>
0028 
0029 #include <algorithm>
0030 #include <cmath>
0031 #include <cstddef>  // for size_t
0032 #include <cstdint>  // for uint64_t
0033 #include <format>
0034 #include <iostream>
0035 #include <limits>
0036 #include <memory>
0037 #include <sstream>  // for basic_ostringstream
0038 
0039 InttCalib::InttCalib(const std::string &name)
0040   : SubsysReco(name)
0041 {
0042 }
0043 
0044 InttCalib::~InttCalib()
0045 {
0046 // only delete histograms which were not registered with the HistoManager
0047 // otherwise the code will segfault in the histomanager cleanup
0048   for (auto &hist : m_hist_fee)
0049   {
0050     delete hist;
0051   }
0052   for (auto &fit : m_fit_fee)
0053   {
0054     delete fit;
0055   }
0056   for (auto &hist : m_hist_half)
0057   {
0058     delete hist;
0059   }
0060 }
0061 
0062 int InttCalib::InitRun(PHCompositeNode * /*unused*/)
0063 {
0064   m_evts = 0;
0065   for (auto const &raw : InttNameSpace::AllRawDataChannels())
0066   {
0067     for (int bco = 0; bco < 129; ++bco)
0068     {
0069       m_hitmap[raw.felix_server][raw.felix_channel][raw.chip][raw.channel][bco] = 0;
0070     }
0071   }
0072 
0073   m_do_nothing = false;
0074   if (Verbosity() > 0)
0075   {
0076     std::cout << "INITRUNEND" << std::endl;
0077   }
0078   return Fun4AllReturnCodes::EVENT_OK;
0079 }
0080 
0081 int InttCalib::process_event(PHCompositeNode *top_node)
0082 {
0083   if (m_do_nothing)
0084   {
0085     return Fun4AllReturnCodes::EVENT_OK;
0086   }
0087   Gl1Packet *gl1 = findNode::getClass<Gl1Packet>(top_node, "GL1RAWHIT");
0088   uint64_t gl1_bco = (gl1 != nullptr) ? (gl1->getBCO() & 0xFFFFFFFFFFU)
0089                                       : std::numeric_limits<uint64_t>::max();
0090 
0091   if (!gl1)
0092   {
0093     std::cout << PHWHERE << "\n"
0094               << "\tCould not get 'GL1RAWHIT' from node tree" << std::endl;
0095     if (m_streaming)
0096     {
0097       std::cout << "\tRunmode is Streaming \n"
0098                 << "\tModule will do nothing" << std::endl;
0099       m_do_nothing = true;
0100       return Fun4AllReturnCodes::EVENT_OK;
0101     }
0102 
0103     std::cout << "\tRunmode is Triggered \n"
0104               << "\tModule will process without GL1" << std::endl;
0105   }
0106   InttRawHitContainer *intt_raw_hit_container =
0107       findNode::getClass<InttRawHitContainer>(top_node, m_rawhit_container_name);
0108   if (!intt_raw_hit_container)
0109   {
0110     std::cout << PHWHERE << "\n"
0111               << "\tCould not get 'INTTRAWHIT' from node tree\n"
0112               << "\tModule will do nothing" << std::endl;
0113     m_do_nothing = true;
0114     // gSystem->Exit(1);
0115     // exit(1);
0116     return Fun4AllReturnCodes::EVENT_OK;
0117   }
0118 
0119   for (size_t n = 0, N = intt_raw_hit_container->get_nhits(); n < N; ++n)
0120   {
0121     InttRawHit *intt_raw_hit = intt_raw_hit_container->get_hit(n);
0122     if (!intt_raw_hit)
0123     {
0124       continue;
0125     }
0126 
0127     InttNameSpace::RawData_s raw{
0128         .felix_server = intt_raw_hit->get_packetid() - 3001,  //
0129         .felix_channel = intt_raw_hit->get_fee(),             //
0130         .chip = (intt_raw_hit->get_chip_id() + 25) % 26,      //
0131         .channel = intt_raw_hit->get_channel_id(),            //
0132     };
0133 
0134     // Trigger mode BCO Offset
0135     // old convention
0136     // int bco_diff = ((intt_raw_hit->get_bco() & 0x7fU) -
0137     // intt_raw_hit->get_FPHX_BCO() + 128) % 128; int bco_diff =
0138     // (intt_raw_hit->get_FPHX_BCO() - (intt_raw_hit->get_bco() & 0x7fU) + 128)
0139     // % 128;
0140     int bco_diff =
0141         (m_streaming) ?
0142                       // Streaming mode BCO Offset : Hit BCO(FPHX BCO + INTT
0143                       // BCO) - GL1 BCO
0144             intt_raw_hit->get_FPHX_BCO() + intt_raw_hit->get_bco() - gl1_bco
0145                       :
0146                       // Trigger mode BCO Offset
0147             (intt_raw_hit->get_FPHX_BCO() - (intt_raw_hit->get_bco() & 0x7fU) +
0148              128) %
0149                 128;
0150     // Not use the hit from abort gap region
0151     if (m_streaming && ((intt_raw_hit->get_FPHX_BCO() > 116) ||
0152                         (intt_raw_hit->get_FPHX_BCO() < 5)))
0153     {
0154       bco_diff = -999;
0155     }
0156 
0157     if (bco_diff > -1)
0158     {
0159       ++m_hitmap[raw.felix_server][raw.felix_channel][raw.chip][raw.channel][bco_diff];
0160     }
0161     ++m_hitmap[raw.felix_server][raw.felix_channel][raw.chip][raw.channel][128];
0162   }
0163 
0164   ++m_evts;
0165   if ((Verbosity() > 1) && ((m_evts % 1000) == 0))
0166   {
0167     std::cout << "event: " << m_evts << std::endl;
0168   }
0169   if (m_evts == m_evts_bco && m_evts_bco != 0)
0170   {
0171     ConfigureBcoMap();
0172     MakeBcoMapCdb();
0173     MakeBcoMapPng();
0174     m_do_make_bco = false;
0175   }
0176   return Fun4AllReturnCodes::EVENT_OK;
0177 }
0178 
0179 int InttCalib::EndRun(int const run_number)
0180 {
0181   if (m_do_nothing)
0182   {
0183     std::cout << PHWHERE << "\n"
0184               << "\tMember 'm_do_nothing' set\n"
0185               << "\tDoing nothing" << std::endl;
0186     return Fun4AllReturnCodes::EVENT_OK;
0187   }
0188 
0189   m_run_num = run_number;
0190   if (m_do_fee)
0191   {
0192     ConfigureHotMap_fee();
0193     MakeHotMapCdb_fee();
0194     MakeHotMapROOT_fee();
0195   }
0196   else
0197   {
0198     ConfigureHotMap_v3();
0199     MakeHotMapCdb_v3();
0200     MakeHotMapPng_v3();
0201   }
0202   if (m_do_make_bco)
0203   {
0204     ConfigureBcoMap();
0205     MakeBcoMapCdb();
0206     MakeBcoMapPng();
0207   }
0208 
0209   return Fun4AllReturnCodes::EVENT_OK;
0210 }
0211 
0212 int InttCalib::ConfigureHotMap_fee()
0213 {
0214   std::map<double, int> hitrate_pdf[m_MAX_LADDER]{};
0215   std::string name[m_MAX_LADDER];
0216   std::string title[m_MAX_LADDER];
0217   for (int i = 0; i < m_MAX_LADDER; ++i)
0218   {
0219     // name[i] = std::format("intt{:01d}", (i / 4));
0220     name[i] = std::format("h_InttCalib_intt{:01d}_fee{}", (i / 14), (i % 14));
0221     title[i] = name[i];
0222   }
0223 
0224   for (auto const &raw : InttNameSpace::AllRawDataChannels())
0225   {
0226     double hitrate =
0227         m_hitmap[raw.felix_server][raw.felix_channel][raw.chip][raw.channel][128] / m_evts;
0228     InttNameSpace::Offline_s ofl = InttNameSpace::ToOffline(raw);
0229 
0230     int index = GetFeeIndex(raw, ofl);
0231     adjust_hitrate(ofl, hitrate);
0232 
0233     ++hitrate_pdf[index][hitrate];
0234   }
0235   std::vector<double> middle_keys(m_MAX_LADDER);
0236   for (int i = 0; i < m_MAX_LADDER; ++i)
0237   {
0238     size_t map_size = hitrate_pdf[i].size();
0239 
0240     size_t mid_index = map_size / 2;
0241     double middle_key = 0.;
0242 
0243     auto it = hitrate_pdf[i].begin();
0244     std::advance(it, mid_index);
0245     middle_key = it->first;
0246     middle_keys[i] = middle_key;
0247   }
0248   std::sort(middle_keys.begin(), middle_keys.end());
0249   double global_maxbin = 5 * (middle_keys[56] + middle_keys[57]) / 2.0;
0250   for (int i = 0; i < m_MAX_LADDER; ++i)
0251   {
0252     ConfigureHist_v3(m_hist_fee[i], m_fit_fee[i], global_maxbin, hitrate_pdf[i], name[i], title[i]);
0253 
0254     double mean = m_fit_fee[i]->GetParameter(1);
0255     double sigma = m_fit_fee[i]->GetParameter(2);
0256     m_mean_fee[i] = mean;
0257     m_sigma_fee[i] = sigma;
0258 
0259     m_min_fee[i] = mean - m_NUM_SIGMA_COLD * sigma;
0260     if (m_min_fee[i] <= 0)
0261     {
0262       m_min_fee[i] = -999;
0263     }
0264     m_max_fee[i] = mean + m_NUM_SIGMA_HOT * sigma;
0265   }
0266 
0267   return Fun4AllReturnCodes::EVENT_OK;
0268 }
0269 
0270 int InttCalib::ConfigureHotMap_v3()
0271 {
0272   std::map<double, int> hitrate_pdf[m_MAX_INDEX]{};
0273   std::string name[m_MAX_INDEX];
0274   std::string title[m_MAX_INDEX];
0275   for (int i = 0; i < m_MAX_INDEX; ++i)
0276   {
0277     // name[i] = std::format("intt{:01d}", % (i / 4));
0278     name[i] = std::format("h_InttCalib_intt{:01d}", i);
0279     title[i] = name[i];
0280   }
0281 
0282   for (auto const &raw : InttNameSpace::AllRawDataChannels())
0283   {
0284     if (m_FELIX_TARGET != -1 && m_FELIX_TARGET != raw.felix_server)
0285     {
0286       continue;
0287     }
0288     double hitrate =
0289         m_hitmap[raw.felix_server][raw.felix_channel][raw.chip][raw.channel][128] / m_evts;
0290     InttNameSpace::Offline_s ofl = InttNameSpace::ToOffline(raw);
0291 
0292     int index = GetIndex(raw, ofl);
0293     adjust_hitrate(ofl, hitrate);
0294 
0295     ++hitrate_pdf[index][hitrate];
0296   }
0297   std::vector<double> middle_keys(8);
0298   double global_maxbin = 0;
0299   for (int i = 0; i < m_MAX_INDEX; ++i)
0300   {
0301     if (m_FELIX_TARGET != -1 && m_FELIX_TARGET != i)
0302     {
0303       continue;
0304     }
0305     size_t map_size = hitrate_pdf[i].size();
0306 
0307     size_t mid_index = map_size / 2;
0308     double middle_key = 0.;
0309 
0310     auto it = hitrate_pdf[i].begin();
0311     std::advance(it, mid_index);
0312     middle_key = it->first;
0313     middle_keys[i] = middle_key;
0314     global_maxbin = 5 * middle_key;
0315   }
0316   std::sort(middle_keys.begin(), middle_keys.end());
0317   if (m_FELIX_TARGET == -1)
0318   {
0319     global_maxbin = 5 * (middle_keys[3] + middle_keys[4]) / 2.0;
0320   }
0321   for (int i = 0; i < m_MAX_INDEX; ++i)
0322   {
0323     if (m_FELIX_TARGET != -1 && m_FELIX_TARGET != i)
0324     {
0325       continue;
0326     }
0327     ConfigureHist_v3(m_hist[i], m_fit[i], global_maxbin, hitrate_pdf[i], name[i], title[i]);
0328     QAHistManagerDef::getHistoManager()->registerHisto(m_hist[i]);
0329     QAHistManagerDef::getHistoManager()->registerHisto(m_fit[i]);
0330     int nBins = m_hist[i]->GetNbinsX();
0331     double xMin = m_hist[i]->GetXaxis()->GetXmin();
0332     double xMax = m_hist[i]->GetXaxis()->GetXmax();
0333 
0334     m_hist_half[i] =
0335         new TH1D(std::format("h_InttCalib_half_hist_{}", i).c_str(),
0336                  "New Histogram with Same Binning", nBins, xMin, xMax);
0337 
0338     double mean = m_fit[i]->GetParameter(1);
0339     double sigma = m_fit[i]->GetParameter(2);
0340     m_mean[i] = mean;
0341     m_sigma[i] = sigma;
0342     // m_min[i] = mean - m_NUM_SIGMA * sigma;
0343 
0344     m_min[i] = mean - m_NUM_SIGMA_COLD * sigma;
0345     if (m_min[i] <= 0)
0346     {
0347       m_min[i] = -999;
0348     }
0349     m_max[i] = mean + m_NUM_SIGMA_HOT * sigma;
0350   }
0351 
0352   return Fun4AllReturnCodes::EVENT_OK;
0353 }
0354 
0355 int InttCalib::MakeHotMapCdb_fee()
0356 {
0357   if (m_hotmap_cdb_file.empty())
0358   {
0359     return Fun4AllReturnCodes::EVENT_OK;
0360   }
0361 
0362   CDBTTree *cdbttree = new CDBTTree(m_hotmap_cdb_file);
0363   int size = 0;
0364   for (auto const &raw : InttNameSpace::AllRawDataChannels())
0365   {
0366     double hitrate =
0367         m_hitmap[raw.felix_server][raw.felix_channel][raw.chip][raw.channel][128] / m_evts;
0368     InttNameSpace::Offline_s ofl = InttNameSpace::ToOffline(raw);
0369 
0370     int index = GetFeeIndex(raw, ofl);
0371     adjust_hitrate(ofl, hitrate);
0372     // Example Part How to add channel masking mannually.
0373     //  if( (raw.pid-3001 == 2) && (raw.felix_channel == 9) && (raw.chip == 15))
0374     //  {
0375     //    cdbttree->SetIntValue(size, "felix_server", raw.felix_server);
0376     //    cdbttree->SetIntValue(size, "felix_channel", raw.felix_channel);
0377     //    cdbttree->SetIntValue(size, "chip", raw.chip);
0378     //    cdbttree->SetIntValue(size, "channel", raw.channel);
0379     //    cdbttree->SetIntValue(size, "flag", 4);
0380     //    ++size;
0381     //    continue;
0382     //  }
0383     // End of Example Part
0384     if (hitrate == 0)
0385     {  // dead channel
0386       cdbttree->SetIntValue(size, "felix_server", raw.felix_server);
0387       cdbttree->SetIntValue(size, "felix_channel", raw.felix_channel);
0388       cdbttree->SetIntValue(size, "chip", raw.chip);
0389       cdbttree->SetIntValue(size, "channel", raw.channel);
0390       cdbttree->SetIntValue(size, "flag", 1);
0391       ++size;
0392       continue;
0393     }
0394     if (m_min_fee[index] < hitrate && hitrate < m_max_fee[index])
0395     {  // good channel
0396       continue;
0397     }
0398     if (hitrate > m_max_fee[index])
0399     {  // hot channel
0400       cdbttree->SetIntValue(size, "felix_server", raw.felix_server);
0401       cdbttree->SetIntValue(size, "felix_channel", raw.felix_channel);
0402       cdbttree->SetIntValue(size, "chip", raw.chip);
0403       cdbttree->SetIntValue(size, "channel", raw.channel);
0404       cdbttree->SetIntValue(size, "flag", 8);
0405       ++size;
0406       continue;
0407     }
0408     if (hitrate < m_min_fee[index])
0409     {  // cold channel
0410       cdbttree->SetIntValue(size, "felix_server", raw.felix_server);
0411       cdbttree->SetIntValue(size, "felix_channel", raw.felix_channel);
0412       cdbttree->SetIntValue(size, "chip", raw.chip);
0413       cdbttree->SetIntValue(size, "channel", raw.channel);
0414       cdbttree->SetIntValue(size, "flag", 4);
0415       ++size;
0416       continue;
0417     }
0418   }
0419   cdbttree->SetSingleIntValue("size", size);
0420   cdbttree->SetSingleIntValue("event", m_evts);
0421   for (int i = 0; i < m_MAX_LADDER; i++)
0422   {
0423     std::string meanname = "mean" + std::to_string(i);
0424     std::string sigmaname = "sigma" + std::to_string(i);
0425     cdbttree->SetSingleDoubleValue(meanname, m_mean_fee[i]);
0426     cdbttree->SetSingleDoubleValue(sigmaname, m_sigma_fee[i]);
0427   }
0428 
0429   cdbttree->Commit();
0430   cdbttree->CommitSingle();
0431   cdbttree->WriteCDBTTree();
0432 
0433   return Fun4AllReturnCodes::EVENT_OK;
0434 }
0435 
0436 int InttCalib::MakeHotMapCdb_v3()
0437 {
0438   if (m_hotmap_cdb_file.empty())
0439   {
0440     return Fun4AllReturnCodes::EVENT_OK;
0441   }
0442   auto cdbttree = std::make_unique<CDBTTree>(m_hotmap_cdb_file);
0443 
0444 //  CDBTTree *cdbttree = new CDBTTree(m_hotmap_cdb_file);
0445   int size = 0;
0446   for (auto const &raw : InttNameSpace::AllRawDataChannels())
0447   {
0448     if (m_FELIX_TARGET != -1 && m_FELIX_TARGET != raw.felix_server)
0449     {
0450       continue;
0451     }
0452     double hitrate =
0453         m_hitmap[raw.felix_server][raw.felix_channel][raw.chip][raw.channel][128] / m_evts;
0454     InttNameSpace::Offline_s ofl = InttNameSpace::ToOffline(raw);
0455 
0456     int index = GetIndex(raw, ofl);
0457     adjust_hitrate(ofl, hitrate);
0458     if (m_half_min[index] < hitrate && hitrate < m_half_max[index])
0459     {
0460       m_hitmap_half[raw.felix_server][raw.felix_channel][raw.chip]++;
0461       // if (m_hitmap_half[raw.felix_server][raw.felix_channel][raw.chip] > 100)
0462       // {
0463       //   cdbttree->SetIntValue(size, "felix_server", raw.felix_server);
0464       //   cdbttree->SetIntValue(size, "felix_channel", raw.felix_channel);
0465       //   cdbttree->SetIntValue(size, "chip", raw.chip);
0466       //   cdbttree->SetIntValue(size, "channel", raw.channel);
0467       //   cdbttree->SetIntValue(size, "flag", 2);
0468       //   ++size;
0469       //   continue;
0470       // }
0471     }
0472 
0473     if (hitrate == 0)
0474     {  // dead channel
0475       cdbttree->SetIntValue(size, "felix_server", raw.felix_server);
0476       cdbttree->SetIntValue(size, "felix_channel", raw.felix_channel);
0477       cdbttree->SetIntValue(size, "chip", raw.chip);
0478       cdbttree->SetIntValue(size, "channel", raw.channel);
0479       cdbttree->SetIntValue(size, "flag", 1);
0480       ++size;
0481       continue;
0482     }
0483     if (m_min[index] < hitrate && hitrate < m_max[index])
0484     {  // good channel
0485       continue;
0486     }
0487     if (hitrate > m_max[index])
0488     {  // hot channel
0489       cdbttree->SetIntValue(size, "felix_server", raw.felix_server);
0490       cdbttree->SetIntValue(size, "felix_channel", raw.felix_channel);
0491       cdbttree->SetIntValue(size, "chip", raw.chip);
0492       cdbttree->SetIntValue(size, "channel", raw.channel);
0493       cdbttree->SetIntValue(size, "flag", 8);
0494       ++size;
0495       continue;
0496     }
0497     if (hitrate < m_min[index])
0498     {  // cold channel
0499       cdbttree->SetIntValue(size, "felix_server", raw.felix_server);
0500       cdbttree->SetIntValue(size, "felix_channel", raw.felix_channel);
0501       cdbttree->SetIntValue(size, "chip", raw.chip);
0502       cdbttree->SetIntValue(size, "channel", raw.channel);
0503       cdbttree->SetIntValue(size, "flag", 4);
0504       ++size;
0505       continue;
0506     }
0507   }
0508   cdbttree->SetSingleIntValue("size", size);
0509   cdbttree->SetSingleIntValue("event", m_evts);
0510   for (int i = 0; i < 8; i++)
0511   {
0512     std::string meanname = "mean" + std::to_string(i);
0513     std::string sigmaname = "sigma" + std::to_string(i);
0514     cdbttree->SetSingleDoubleValue(meanname, m_mean[i]);
0515     cdbttree->SetSingleDoubleValue(sigmaname, m_sigma[i]);
0516   }
0517 
0518   cdbttree->Commit();
0519   cdbttree->CommitSingle();
0520   cdbttree->WriteCDBTTree();
0521   return Fun4AllReturnCodes::EVENT_OK;
0522 }
0523 
0524 int InttCalib::MakeHotMapROOT_fee()
0525 {
0526   const int rows = 8;
0527   const int cols = 14;
0528   const int numHists = rows * cols;
0529   TCanvas *c1 = new TCanvas("c1", "Histograms with Fits", 2000, 1400);
0530   c1->Divide(cols, rows);
0531 
0532   gStyle->SetOptStat(0);
0533   gStyle->SetPadTopMargin(0.01);
0534   gStyle->SetPadBottomMargin(0.01);
0535   gStyle->SetPadLeftMargin(0.01);
0536   gStyle->SetPadRightMargin(0.01);
0537 
0538   for (int i = 0; i < numHists; ++i)
0539   {
0540     c1->cd(i + 1);
0541     m_hist_fee[i]->Draw();
0542     m_fit_fee[i]->Draw("same");
0543   }
0544   if (!m_hotmap_png_file.empty())
0545   {
0546     c1->SaveAs(m_hotmap_png_file.c_str());
0547   }
0548   return Fun4AllReturnCodes::EVENT_OK;
0549 }
0550 
0551 int InttCalib::MakeHotMapPng_v3()
0552 {
0553   // canvas
0554   gStyle->SetOptStat(0);
0555   TCanvas *cnvs = new TCanvas("hitrate_cnvs",  //
0556                               "hitrate_cnvs",  //
0557                               1280, 720        //
0558   );
0559   double n_hot = 0;
0560   double n_cold = 0;
0561   double n_dead = 0;
0562   double n_total = 0;
0563   for (auto const &raw : InttNameSpace::AllRawDataChannels())
0564   {
0565     double hitrate =
0566         m_hitmap[raw.felix_server][raw.felix_channel][raw.chip][raw.channel][128] / m_evts;
0567     InttNameSpace::Offline_s ofl = InttNameSpace::ToOffline(raw);
0568 
0569     int index = GetIndex(raw, ofl);
0570     adjust_hitrate(ofl, hitrate);
0571 
0572     if (!(m_min[index] < hitrate))
0573     {
0574       ++n_cold;
0575     }
0576     if (!(hitrate < m_max[index]))
0577     {
0578       ++n_hot;
0579     }
0580     if (hitrate == 0)
0581     {
0582       ++n_dead;
0583     }
0584     ++n_total;
0585     if (m_hitmap_half[raw.felix_server][raw.felix_channel][raw.chip] > 100)
0586     {
0587       if(m_hist_half[raw.felix_server])
0588       {
0589       m_hist_half[raw.felix_server]->Fill(hitrate);
0590       }
0591     }
0592   }
0593   for (int j = 0; j < 8; ++j)
0594   {
0595     std::string name = std::format("hist_pad_{:01d}", j);
0596 
0597     cnvs->cd();
0598     TPad *hist_pad = new TPad(  //
0599         name.c_str(),
0600         name.c_str(),  //
0601                        // (j % 4 + 0.0) / 4.0 * 0.9 + 0.0, (1.0 - j / 4) / 2.0 *
0602                        // 0.9 + 0.1, // (j % 4 + 1.0) / 4.0 * 0.9 + 0.0, (2.0 - j
0603                        // / 4) / 2.0 * 0.9 + 0.1  //
0604                        // NOLINTNEXTLINE(bugprone-integer-division)
0605         ((j % 4 + 0.0) / 4.0 * 1.0) + 0.0,
0606         // NOLINTNEXTLINE(bugprone-integer-division)
0607         ((1.0 - j / 4) / 2.0 * 0.9) + 0.1,  //
0608                                             // NOLINTNEXTLINE(bugprone-integer-division)
0609         ((j % 4 + 1.0) / 4.0 * 1.0) + 0.0,
0610         // NOLINTNEXTLINE(bugprone-integer-division)
0611         ((2.0 - j / 4) / 2.0 * 0.9) + 0.1  //
0612     );
0613 
0614     hist_pad->SetFillStyle(4000);
0615     hist_pad->Range(0.0, 0.0, 1.0, 1.0);
0616     hist_pad->SetLogy();
0617     hist_pad->Draw();
0618 
0619     hist_pad->cd();
0620     double x_max = 0;
0621     double y_max = 0;
0622     // for(int i = j * 4; i < (j + 1) * 4; ++i)
0623     for (int i = j; i < j + 1; ++i)
0624     {
0625       // m_hist[i]->SetLineColor(GetFeeColor(i - j * 4));
0626       if (!m_hist[i])
0627       {
0628         continue;
0629       }
0630       m_hist[i]->SetLineColor(kBlack);
0631       m_hist[i]->SetLineWidth(2);
0632       m_hist_half[i]->SetLineColor(kRed);
0633       m_hist_half[i]->SetLineWidth(3);
0634 
0635       // m_fit[i]->SetLineColor(GetFeeColor(i - j * 4));
0636       m_fit[i]->SetLineColor(kBlue);
0637       m_fit[i]->SetLineWidth(2);
0638 
0639       double temp_max;
0640 
0641       temp_max = m_hist[i]->GetBinContent(m_hist[i]->GetMaximumBin());
0642       y_max = std::max(y_max, temp_max);
0643 
0644       temp_max = m_hist[i]->GetXaxis()->GetBinLowEdge(
0645           m_hist[i]->GetXaxis()->GetNbins() - 1);
0646       temp_max += m_hist[i]->GetXaxis()->GetBinWidth(
0647           m_hist[i]->GetXaxis()->GetNbins() - 1);
0648       x_max = std::max(x_max, temp_max);
0649     }
0650     y_max *= 10;
0651 
0652     for (int i = j; i < j + 1; ++i)
0653     {
0654       if (!m_hist[i])
0655       {
0656         continue;
0657       }
0658       m_hist[i]->GetXaxis()->SetRangeUser(0, x_max);
0659       m_hist[i]->GetYaxis()->SetRangeUser(1, y_max);
0660       m_hist[i]->Draw("same");
0661       // histogram for half entry chips. Not used for now
0662       // m_hist_half[i]->Draw("same");
0663       m_fit[i]->Draw("same");
0664 
0665       TLine line;
0666       line.SetLineColor(kRed);
0667       line.SetLineWidth(2);
0668       line.DrawLine(m_max[i], 0, m_max[i], y_max);
0669 
0670       TLine line2;
0671       line2.SetLineColor(kBlue);
0672       line2.DrawLine(m_min[i], 0, m_min[i], y_max);
0673     }
0674   }
0675   cnvs->cd();
0676   TPad *legend_pad = new TPad("legend_pad", "legend_pad", 0.9, 0.1, 1.0, 1.0);
0677   legend_pad->SetFillStyle(4000);
0678   legend_pad->Range(0.0, 0.0, 1.0, 1.0);
0679   legend_pad->Draw(); // this adds the tpad to the canvas, so it gets deleted when the tcanvas is deleted
0680   legend_pad->cd();
0681   for (int i = 0; i < 4; ++i)
0682   {
0683     TText text;
0684     text.SetTextColor(GetFeeColor(i));
0685     text.SetTextAlign(22);
0686     text.SetTextSize(0.15);
0687     if (!m_hist[i])
0688     {
0689       continue;
0690     }
0691     std::string title = m_hist[i]->GetName();
0692     text.DrawText(0.5, (2.0 * i + 1.0) / (2.0 * 4), title.substr(6, 7).c_str());
0693   }
0694 
0695   cnvs->cd();
0696   TPad *caption_pad =
0697       new TPad("caption_pad", "caption_pad", 0.0, 0.0, 1.0, 0.1);
0698   caption_pad->SetFillStyle(4000);
0699   caption_pad->Range(0.0, 0.0, 1.0, 1.0);
0700   caption_pad->Draw();
0701 
0702   caption_pad->cd();
0703   TText caption;
0704   caption.SetTextColor(kBlack);
0705   caption.SetTextAlign(22);
0706   caption.SetTextSize(0.25);
0707   caption.DrawText(0.5, 0.75, std::format("Run: {:08d} Events: {}", m_run_num, m_evts).c_str());
0708   caption.DrawText(0.5, 0.50,
0709                    std::format("Fraction Cold: {:.3f}% Fraction Dead: {:.3f}%", (n_cold * 100 / n_total), (n_dead * 100 / n_total)).c_str());
0710   caption.DrawText(0.5, 0.25, std::format("Fraction Hot: {:.3f}%", (n_hot * 100 / n_total)).c_str());
0711 
0712   cnvs->Update();
0713   cnvs->Show();
0714   if (!m_hotmap_png_file.empty())
0715   {
0716     cnvs->SaveAs(m_hotmap_png_file.c_str());
0717   }
0718   delete cnvs;
0719   return Fun4AllReturnCodes::EVENT_OK;
0720 }
0721 
0722 int InttCalib::ConfigureBcoMap()
0723 {
0724   for (int felix = 0; felix < 8; felix++)
0725   {
0726     for (int fee = 0; fee < 14; fee++)
0727     {
0728       if (m_FELIX_TARGET != -1 && m_FELIX_TARGET != felix)
0729       {
0730         continue;
0731       }
0732       // Find chip with highest total hits for this pid/fee across all channels
0733       int chp_most = 0;
0734       int max_hits = 0;
0735       for (int chp = 0; chp < 26; chp++)
0736       {
0737         int chip_total = 0;
0738         // Sum hits across all channels for this chip
0739         for (int chan = 0; chan < 128; chan++)
0740         {
0741           chip_total += m_hitmap[felix][fee][chp][chan][128];
0742         }
0743         if (chip_total > max_hits)
0744         {
0745           max_hits = chip_total;
0746           chp_most = chp;
0747         }
0748       }
0749       for (int chp = 0; chp < 26; chp++)
0750       {
0751         // Sum hits across all channels for this chip
0752         for (int chan = 0; chan < 128; chan++)
0753         {
0754           if (chp != chp_most)
0755           {
0756             for (int bco = 0; bco < 128; ++bco)
0757             {
0758               m_bcorates_fee[felix][fee][bco] +=
0759                   m_hitmap[felix][fee][chp][chan][bco];
0760             }
0761           }
0762         }
0763       }
0764     }
0765   }
0766 
0767   for (int felix = 0; felix < 8; felix++)
0768   {
0769     if (m_FELIX_TARGET != -1 && m_FELIX_TARGET != felix)
0770     {
0771       continue;
0772     }
0773     for (int fee = 0; fee < 14; fee++)
0774     {
0775       int max_counts = 0;
0776       int bco_peak = 0;
0777       for (int bco = 0; bco < 128; bco++)
0778       {
0779         if (max_counts < m_bcorates_fee[felix][fee][bco])
0780         {
0781           bco_peak = bco;
0782           max_counts = m_bcorates_fee[felix][fee][bco];
0783         }
0784       }
0785 
0786       if (max_counts < 50)  // if max_count is less than 50(masked ladder but
0787                             // somethimes it has few hits), set bco_peak as -1
0788       {
0789         bco_peak = -1;
0790       }
0791       m_bcopeaks_fee[felix][fee] = bco_peak;
0792     }
0793   }
0794 
0795   return Fun4AllReturnCodes::EVENT_OK;
0796 }
0797 
0798 int InttCalib::MakeBcoMapCdb()
0799 {
0800   if (m_bcomap_cdb_file.empty())
0801   {
0802     return Fun4AllReturnCodes::EVENT_OK;
0803   }
0804 
0805   auto cdbttree = std::make_unique<CDBTTree>(m_bcomap_cdb_file);
0806 
0807   int size = 0;
0808   std::vector<int> bco_temp_container;
0809   bco_temp_container.clear();
0810   for (int felix = 0; felix < 8; felix++)
0811   {
0812     if (m_FELIX_TARGET != -1 && m_FELIX_TARGET != felix)
0813     {
0814       continue;
0815     }
0816     std::string name_bco_peak = "h_InttCalib_BCOOffSet_INTT" + std::to_string(felix);
0817     m_bco_peak[felix] = new TH1I(name_bco_peak.c_str(), name_bco_peak.c_str(), 14, 0, 14);
0818     for (int fee = 0; fee < 14; fee++)
0819     {
0820       int bco = m_bcopeaks_fee[felix][fee];
0821       cdbttree->SetIntValue(size, "felix_server", felix);
0822       cdbttree->SetIntValue(size, "felix_channel", fee);
0823       cdbttree->SetIntValue(size, "bco_diff", bco);
0824       m_bco_peak[felix]->SetBinContent(fee + 1, bco);
0825       bco_temp_container.push_back(bco);
0826       ++size;
0827     }
0828     QAHistManagerDef::getHistoManager()->registerHisto(m_bco_peak[felix]);
0829   }
0830 
0831   cdbttree->SetSingleIntValue("size", size);
0832   std::pair<double, double> stats =
0833       CalculateStandardDeviation(bco_temp_container);
0834   m_bco_mean = stats.first;
0835   m_bco_stdDev = stats.second;
0836   int runmode = m_streaming ? 1 : 0;  //
0837   cdbttree->SetSingleDoubleValue("StdDev", m_bco_stdDev);
0838   cdbttree->SetSingleIntValue("runmode", runmode);
0839   cdbttree->SetSingleIntValue("events", m_evts);
0840   cdbttree->Commit();
0841   cdbttree->CommitSingle();
0842   cdbttree->WriteCDBTTree();
0843   bco_temp_container.clear();
0844   return Fun4AllReturnCodes::EVENT_OK;
0845 }
0846 
0847 int InttCalib::MakeBcoMapPng()
0848 {
0849   // Canvas
0850   gStyle->SetOptStat(0);
0851   TCanvas *bco_cnvs = new TCanvas(  //
0852       "bco_cnvs", "bco_cnvs",       //
0853       1280, 720                     //
0854   );
0855   bco_cnvs->Draw();
0856 
0857   TH1 *bco_hist[112] = {};
0858   for (int i = 0; i < 16; ++i)
0859   {
0860     std::string tpadtitle = std::format("bco_pad_{:02d}", i);
0861     bco_cnvs->cd();
0862     TPad *bco_pdf_pad = new TPad(  //
0863         tpadtitle.c_str(), tpadtitle.c_str(),
0864         // NOLINTNEXTLINE(bugprone-integer-division)
0865         ((i % 4 + 0.0) / 4.0 * 0.9) + 0.0, ((3.0 - i / 4) / 4.0 * 0.9) + 0.1,
0866         // NOLINTNEXTLINE(bugprone-integer-division)
0867         ((i % 4 + 1.0) / 4.0 * 0.9) + 0.0, ((4.0 - i / 4) / 4.0 * 0.9) + 0.1);
0868     bco_pdf_pad->SetFillStyle(4000);  // transparent
0869     bco_pdf_pad->SetLogy();
0870     bco_pdf_pad->SetLeftMargin(0.15);
0871     bco_pdf_pad->SetRightMargin(0.05);
0872     if ((i / 4) % 2)
0873     {
0874       bco_pdf_pad->SetTopMargin(0.0);
0875       bco_pdf_pad->SetBottomMargin(0.15);
0876     }
0877     else
0878     {
0879       bco_pdf_pad->SetTopMargin(0.15);
0880       bco_pdf_pad->SetBottomMargin(0.0);
0881     }
0882     bco_pdf_pad->Range(0.0, 0.0, 1.0, 1.0);
0883     bco_pdf_pad->Draw();
0884 
0885     int felix = (i % 4) + (4 * (i / 8));
0886     int fee_index_start = (i / 4) % 2 ? 0 : 7;
0887     int fee_index_end = (i / 4) % 2 ? 7 : 14;
0888     double max = 0;
0889     for (int fee = fee_index_start; fee < fee_index_end; fee++)
0890     {
0891       int h = (felix * 14) + fee;
0892       bco_pdf_pad->cd();
0893       std::string htitle = std::format("h_InttCalib_bco_hist_{:03d}", h);
0894       // float min_bin = (m_streaming) ? -127.5 : -0.5;
0895       float min_bin = -0.5;
0896       float max_bin = 127.5;
0897       // int nbins = (m_streaming) ? 256 : 128;
0898       int nbins = 128;
0899       bco_hist[h] = new TH1D(              //
0900           htitle.c_str(), htitle.c_str(),  //
0901           nbins, min_bin, max_bin          //
0902                                            // 128, -0.5, 127.5                //
0903       );
0904       bco_hist[h]->SetTitle(std::format(";BCO Difference;intt{:01d} ({:01d} - {:02d})", felix, fee, fee).c_str());
0905       bco_hist[h]->GetYaxis()->CenterTitle();
0906       bco_hist[h]->GetYaxis()->SetTitleSize(0.12);
0907       bco_hist[h]->GetYaxis()->SetTitleOffset(0.6);
0908       bco_hist[h]->GetYaxis()->SetLabelSize(0.07);
0909       bco_hist[h]->GetXaxis()->SetTitleSize(0.10);
0910       bco_hist[h]->GetXaxis()->SetTitleOffset(0.6);
0911       bco_hist[h]->GetXaxis()->SetLabelSize(0.07);
0912       bco_hist[h]->SetLineColor(GetFeeColor(fee));
0913       bco_hist[h]->Draw("same");
0914 
0915       for (int bco = 0; bco < 128; ++bco)
0916       {
0917         bco_hist[h]->SetBinContent(bco + 1, m_bcorates_fee[felix][fee][bco]);
0918         max = std::max<double>(max, m_bcorates_fee[felix][fee][bco]);
0919       }
0920     }
0921     //  for (InttNameSpace::RawData_s raw = raw_begin; raw <= raw_end; ++raw)
0922     //  {
0923     //    int h = (raw.felix_server) * 14 + raw.felix_channel;
0924     //    bco_hist[h]->GetYaxis()->SetRangeUser(0.1, 10 * max);
0925     //  }
0926   }
0927 
0928   // Legend
0929   bco_cnvs->cd();
0930   TPad *legend_pad = new TPad("legend_pad", "legend_pad",  //
0931                               0.9, 0.1, 1.0, 1.0           //
0932   );
0933   legend_pad->SetFillStyle(4000);  // transparent
0934   legend_pad->Range(0.0, 0.0, 1.0, 1.0);
0935   legend_pad->Draw();
0936   legend_pad->cd();
0937 
0938   for (int fee = 0; fee < 7; ++fee)
0939   {
0940     TText legend_text;
0941     legend_text.SetTextAlign(22);
0942     legend_text.SetTextColor(kBlack);
0943     legend_text.SetTextSize(0.1);
0944     legend_text.DrawText(
0945         0.6, (2.0 * fee + 1.0) / 14.0,
0946         std::format("FCh {:01d}, {:02d}", fee, (fee + 7)).c_str());
0947 
0948     double x[4] = {-1.0, +1.0, +1.0, -1.0};
0949     double y[4] = {-1.0, -1.0, +1.0, +1.0};
0950     for (int i = 0; i < 4; ++i)
0951     {
0952       x[i] *= 0.1;
0953       x[i] += 0.2;
0954 
0955       y[i] *= 0.008;
0956       y[i] += (2.0 * fee + 1.0) / 14.0;
0957     }
0958 
0959     TPolyLine box;
0960     box.SetFillColor(GetFeeColor(fee));
0961     box.SetLineColor(kBlack);
0962     box.SetLineWidth(1);
0963     box.DrawPolyLine(4, x, y, "f");
0964   }
0965 
0966   // Caption
0967   bco_cnvs->cd();
0968   TPad *caption_pad = new TPad("caption_pad", "caption_pad",  //
0969                                0.0, 0.0, 1.0, 0.1             //
0970   );
0971   caption_pad->SetFillStyle(4000);  // transparent
0972   caption_pad->Range(0.0, 0.0, 1.0, 1.0);
0973   caption_pad->Draw();
0974 
0975   caption_pad->cd();
0976   TText caption_text;
0977   caption_text.SetTextAlign(22);
0978   caption_text.SetTextSize(0.40);
0979   caption_text.SetTextColor(kBlack);
0980   caption_text.DrawText(
0981       0.5, 0.75,
0982       std::format("Run: {:08d} Events: {} BCO StdDev: {:.2f} BCO Offset: {:.2f}",
0983                   m_run_num, m_evts, m_bco_stdDev, m_bco_mean)
0984           .c_str());
0985 
0986   TText caption_tag;
0987   caption_tag.SetTextAlign(22);
0988   caption_tag.SetTextSize(0.40);
0989   if (m_evts < 1000)
0990   {
0991     caption_tag.SetTextColor(kRed);
0992     caption_tag.DrawText(0.5, 0.3, "LOW STATSTICS events < 1000!!");
0993   }
0994   else if (m_bco_stdDev != 0)
0995   {
0996     caption_tag.SetTextColor(kRed);
0997     caption_tag.DrawText(0.5, 0.3, "FEE MISALIGNED!!");
0998   }
0999   else
1000   {
1001     caption_tag.SetTextColor(kBlue);
1002     caption_tag.DrawText(0.5, 0.3, "GOOD");
1003   }
1004 
1005   std::map<int, std::vector<int>> maskedLadders;  // To store fee -> fc mappings
1006 
1007   // Iterate through the bco_hist array
1008   for (int i = 0; i < 112; ++i)
1009   {
1010     if (bco_hist[i]->GetBinContent(bco_hist[i]->GetMaximumBin()) < 5)
1011     {
1012       int sev = i / 14;
1013       int fee = i % 14;
1014       // Add fc to the corresponding fee in the map
1015       maskedLadders[sev].push_back(fee);
1016     }
1017   }
1018 
1019   // Now, prepare the text to add to the canvas
1020   std::ostringstream maskedText;
1021   maskedText << "Masked Ladders: ";
1022 
1023   // Iterate through the map and create the output text
1024   for (const auto &entry : maskedLadders)
1025   {
1026     maskedText << "INTT " << entry.first << "(FC ";  // Fee numbering starts at 1
1027     for (auto j = 0U; j < entry.second.size(); ++j)
1028     {
1029       maskedText << entry.second[j];
1030       if (j != entry.second.size() - 1)
1031       {
1032         maskedText << ",";  // Add commas between FC values
1033       }
1034     }
1035     maskedText << ") ";
1036   }
1037 
1038   // Print or draw the text on the canvas
1039   TText caption_masked;
1040   caption_masked.SetTextAlign(32);
1041   caption_masked.SetTextSize(0.20);
1042   caption_masked.DrawText(0.9, 0.3, maskedText.str().c_str());
1043   bco_cnvs->Update();
1044   bco_cnvs->Show();
1045   if (!m_bcomap_png_file.empty())
1046   {
1047     bco_cnvs->SaveAs(m_bcomap_png_file.c_str());
1048   }
1049 
1050   delete bco_cnvs;
1051 
1052   return Fun4AllReturnCodes::EVENT_OK;
1053 }
1054 
1055 int InttCalib::SaveHitrates()
1056 {
1057   TFile *file = TFile::Open(m_hotmap_png_file.c_str(), "RECREATE");
1058   if (!file)
1059   {
1060     std::cerr << "\n"
1061               << PHWHERE << "\n\tfile\n"
1062               << std::endl;
1063     return 1;
1064   }
1065   file->cd();
1066   TTree *tree = new TTree("hitrate_tree", "hitrate_tree");
1067   tree->SetDirectory(file);
1068 
1069   int pid{0};
1070   int fee{0};
1071   int chp{0};
1072   int chn{0};
1073   tree->Branch("pid", &pid);
1074   tree->Branch("fee", &fee);
1075   tree->Branch("chp", &chp);
1076   tree->Branch("chn", &chn);
1077 
1078   double hitrate;
1079   tree->Branch("hitrate", &hitrate);
1080 
1081   for (auto const &raw : InttNameSpace::AllRawDataChannels())
1082   {
1083     pid = raw.felix_server + 3001;
1084     fee = raw.felix_channel;
1085     chp = raw.chip;
1086     chn = raw.channel;
1087 
1088     hitrate = m_hitmap[raw.felix_server][raw.felix_channel][raw.chip][raw.channel][128] / m_evts;
1089     tree->Fill();
1090   }
1091 
1092   tree->Write();
1093   file->Write();
1094   file->Close();
1095 
1096   return 0;
1097 }
1098 
1099 int InttCalib::LoadHitrates()
1100 {
1101   if (m_hotmap_png_file.empty())
1102   {
1103     std::cout << PHWHERE << "No TFile with hitrates specified with SetHotMapPngFile(filename)" << std::endl;
1104   }
1105   TFile *file = TFile::Open(m_hotmap_png_file.c_str(), "READ");
1106   if (!file)
1107   {
1108     std::cerr << "\n"
1109               << PHWHERE << "\n\tfile\n"
1110               << std::endl;
1111     return 1;
1112   }
1113 
1114   TTree *tree = dynamic_cast<TTree *>(file->Get("hitrate_tree"));
1115   if (!tree)
1116   {
1117     std::cerr << "\n"
1118               << PHWHERE << "\n\ttree\n"
1119               << std::endl;
1120     return 1;
1121   }
1122 
1123   int pid{0};
1124   int fee{0};
1125   int chp{0};
1126   int chn{0};
1127   tree->SetBranchAddress("pid", &pid);
1128   tree->SetBranchAddress("fee", &fee);
1129   tree->SetBranchAddress("chp", &chp);
1130   tree->SetBranchAddress("chn", &chn);
1131 
1132   double hitrate;
1133   tree->SetBranchAddress("hitrate", &hitrate);
1134 
1135   for (Int_t n = 0, N = tree->GetEntriesFast(); n < N; ++n)
1136   {
1137     tree->GetEntry(n);
1138     InttNameSpace::RawData_s raw{
1139         .felix_server = pid + 3001,
1140         .felix_channel = fee,
1141         .chip = chp,
1142         .channel = chn,
1143     };
1144     m_hitmap[raw.felix_server][raw.felix_channel][raw.chip][raw.channel][128] = hitrate;
1145   }
1146 
1147   m_evts = 1.0;
1148 
1149   return 0;
1150 }
1151 int InttCalib::ConfigureHist_v3(TH1 *&hist, TF1 *&fit, double maxbin,
1152                                 std::map<double, int> const &hitrate_map,
1153                                 std::string const &name,
1154                                 std::string const &title)
1155 {
1156   size_t map_size = hitrate_map.size();
1157   double global_middle = maxbin / 5;  // Global middle value is maxbin / 5
1158   size_t mid_index = map_size / 2;
1159   double middle_key = 0.;
1160 
1161   auto it = hitrate_map.begin();
1162   std::advance(it, mid_index);
1163   middle_key = it->first;
1164 
1165   // Make hist
1166   delete hist;
1167   hist = new TH1D(                                        //
1168       (name + " hitrates").c_str(),                       //
1169       title.c_str(),                                      //
1170       100, std::next(hitrate_map.begin())->first, maxbin  //
1171   );
1172 
1173   TH1 *fit_hist = new TH1D(                               //
1174       (name + " fit_hitrates").c_str(),                   //
1175       (title + " (fit)").c_str(),                         //
1176       100, std::next(hitrate_map.begin())->first, maxbin  //
1177   );
1178 
1179   for (auto const &[hitrate, count] : hitrate_map)
1180   {
1181     for (int i = 0; i < count; ++i)
1182     {
1183       hist->Fill(hitrate);
1184       if (hitrate / global_middle > 0.7)
1185       {
1186         fit_hist->Fill(hitrate);
1187       }
1188     }
1189   }
1190 
1191   delete fit;
1192   fit = new TF1(                //
1193       (name + "_fit").c_str(),  //
1194       "gaus",                   //
1195       middle_key / 10, maxbin   //
1196   );
1197 
1198   if (Verbosity())
1199   {
1200     fit_hist->Fit(fit, "R");  // range, no-draw, log likelihood
1201   }
1202   else
1203   {
1204     fit_hist->Fit(fit, "QR");  // range, no-draw, log likelihood, quiet
1205   }
1206 
1207   std::cout << "Global middle: " << global_middle << std::endl;
1208   delete fit_hist;
1209   return 0;
1210 }
1211 
1212 int InttCalib::adjust_hitrate(InttNameSpace::Offline_s const &ofl,
1213                               double &hitrate)
1214 {
1215   hitrate /= (ofl.ladder_z % 2) ? 2.0 : 1.6;  // normalize by sensor length
1216   if (ofl.layer < 5)
1217   {
1218     hitrate /= (10.005 / 7.4994);  // Inner = 7.4994, Outer = 10.005
1219   }
1220 
1221   // InttSurveyMap::val_t transform;
1222   // if (m_survey.GetStripTransform(ofl, transform))
1223   // {
1224   //   std::cout << PHWHERE << "\n"
1225   //             << "InttSurveyMap::GetStripTransform failed for\n"
1226   //             << ofl
1227   //             << std::endl;
1228   //   return 1;
1229   // }
1230   // Eigen::Vector3d v = transform.translation() - m_vertex;
1231   // Eigen::Vector3d n{transform.linear()(0, 1), transform.linear()(1, 1),
1232   // transform.linear()(2, 1)}; hitrate *= v.squaredNorm(); hitrate /= (-1.0 *
1233   // n.dot(v) / v.norm());
1234   return 0;
1235 }
1236 
1237 int InttCalib::GetIndex(InttNameSpace::RawData_s const &raw,
1238                         InttNameSpace::Offline_s const & /*ofl*/)
1239 {
1240   // int index = 0;
1241   // index += (raw.felix_server) * 4;
1242   // index += (ofl.layer < 5) ? 0 : 2; // +2 for outer
1243   // index += (ofl.ladder_z % 2) ? 0 : 1; // +1 for type B
1244   // return index;
1245   return raw.felix_server;
1246 }
1247 
1248 int InttCalib::GetFeeIndex(InttNameSpace::RawData_s const &raw,
1249                            InttNameSpace::Offline_s const & /*ofl*/)
1250 {
1251   int index = 0;
1252   index = (raw.felix_server) * 14 + raw.felix_channel;
1253   return index;
1254 }
1255 std::pair<double, double>
1256 InttCalib::CalculateStandardDeviation(const std::vector<int> &data)
1257 {
1258   if (data.empty())
1259   {
1260     return std::make_pair(-1, -1);
1261   }
1262 
1263   double sum = 0.0;
1264   int n_masked_ladder = 0;
1265   for (int i : data)
1266   {
1267     if (i == -1)
1268     {
1269       n_masked_ladder++;
1270       continue;
1271     }
1272     sum += i;
1273   }
1274   double mean = sum / static_cast<double>(data.size() - n_masked_ladder);
1275   double sumSquaredDiffs = 0.0;
1276   //  int count = 0;
1277   for (int i : data)
1278   {
1279     if (i == -1)
1280     {
1281       continue;  // do not include maksed ladder for std calculation
1282     }
1283     sumSquaredDiffs += (i - mean) * (i - mean);
1284     //    count++;
1285   }
1286   double stddev = std::sqrt(sumSquaredDiffs /
1287                             static_cast<double>(data.size() - n_masked_ladder));
1288   // return stddev;
1289   return std::make_pair(mean, stddev);
1290 }
1291 
1292 Color_t InttCalib::GetFeeColor(int fee)
1293 {
1294   switch (fee % 7)
1295   {
1296   case 1:
1297     return kRed;
1298   case 2:
1299     return kGreen;
1300   case 3:
1301     return kBlue;
1302   case 4:
1303     return kOrange;
1304   case 5:
1305     return kMagenta;
1306   case 6:
1307     return kCyan;
1308   default:  // 0
1309     break;
1310   }
1311   return kBlack;
1312 }