Back to home page

sPhenix code displayed by LXR

 
 

    


File indexing completed on 2025-08-05 08:11:10

0001 // -- c++ includes --
0002 #include <filesystem>
0003 #include <fstream>
0004 #include <iomanip>
0005 #include <iostream>
0006 #include <string>
0007 
0008 // -- root includes --
0009 #include <TFile.h>
0010 #include <TProfile2D.h>
0011 #include <TLatex.h>
0012 
0013 #include <calobase/TowerInfoDefs.h>
0014 #include <cdbobjects/CDBTTree.h>
0015 
0016 using std::cerr;
0017 using std::cout;
0018 using std::endl;
0019 using std::max;
0020 using std::min;
0021 using std::ofstream;
0022 using std::pair;
0023 using std::string;
0024 using std::stringstream;
0025 using std::to_string;
0026 using std::vector;
0027 
0028 namespace myAnalysis
0029 {
0030 
0031   void histToCaloCDBTree(string outputfile, string fieldName, Int_t icalo, TProfile2D *hist);
0032   void analyze(const string &output);
0033 
0034   // utils
0035   pair<string, string> getRunDataset(const string &input);
0036   void makeDir(const string &output);
0037   Int_t readHists(const string &input);
0038 
0039   pair<string, string> run_dataset;
0040 
0041   TProfile2D *h_CaloFittingQA_cemc_etaphi_ZScrosscalib = nullptr;
0042   TProfile2D *h_CaloFittingQA_ihcal_etaphi_ZScrosscalib = nullptr;
0043   TProfile2D *h_CaloFittingQA_ohcal_etaphi_ZScrosscalib = nullptr;
0044 
0045   UInt_t cemc_bins_eta = 96;
0046   UInt_t cemc_bins_phi = 256;
0047   UInt_t hcal_bins_eta = 24;
0048   UInt_t hcal_bins_phi = 64;
0049 }  // namespace myAnalysis
0050 
0051 void myAnalysis::makeDir(const string &output)
0052 {
0053   if (std::filesystem::exists(output))
0054   {
0055     std::cout << "Directory '" << output << "' already exists." << std::endl;
0056   }
0057   else
0058   {
0059     try
0060     {
0061       std::filesystem::create_directory(output);
0062       std::cout << "Directory '" << output << "' created successfully." << std::endl;
0063     }
0064     catch (const std::filesystem::filesystem_error &e)
0065     {
0066       // Handle other potential errors
0067       std::cerr << "Error creating directory: " << e.what() << std::endl;
0068     }
0069   }
0070 }
0071 
0072 pair<string, string> myAnalysis::getRunDataset(const string &input)
0073 {
0074   string basename = std::filesystem::path(input).filename();
0075   string run = basename.substr(0, basename.find("_"));
0076   string dataset = basename.substr(basename.find("_") + 1, basename.size() - basename.find("_") - 6);
0077   return make_pair(run, dataset);
0078 }
0079 
0080 Int_t myAnalysis::readHists(const string &input)
0081 {
0082   // Create an input stream
0083   std::ifstream file(input);
0084 
0085   // Check if the file was successfully opened
0086   if (!file.is_open())
0087   {
0088     cerr << "Failed to open file list: " << input << endl;
0089     return 1;
0090   }
0091 
0092   run_dataset = getRunDataset(input);
0093 
0094   cout << "Reading Hists" << endl;
0095   cout << "======================================" << endl;
0096 
0097   delete h_CaloFittingQA_cemc_etaphi_ZScrosscalib;
0098   delete h_CaloFittingQA_ihcal_etaphi_ZScrosscalib;
0099   delete h_CaloFittingQA_ohcal_etaphi_ZScrosscalib;
0100 
0101   h_CaloFittingQA_cemc_etaphi_ZScrosscalib = new TProfile2D("cemc_etaphi_ZScrosscalib","", cemc_bins_eta, 0, cemc_bins_eta, cemc_bins_phi, 0, cemc_bins_phi);
0102   h_CaloFittingQA_ihcal_etaphi_ZScrosscalib = new TProfile2D("ihcal_etaphi_ZScrosscalib","", hcal_bins_eta, 0, hcal_bins_eta, hcal_bins_phi, 0, hcal_bins_phi);
0103   h_CaloFittingQA_ohcal_etaphi_ZScrosscalib = new TProfile2D("ohcal_etaphi_ZScrosscalib","", hcal_bins_eta, 0, hcal_bins_eta, hcal_bins_phi, 0, hcal_bins_phi);
0104 
0105   string line;
0106   while (std::getline(file, line))
0107   {
0108     cout << "Reading File: " << line << endl;
0109     auto tf = TFile::Open(line.c_str());
0110 
0111     auto h = static_cast<TProfile2D*>(tf->Get("h_CaloFittingQA_cemc_etaphi_ZScrosscalib"));
0112 
0113     if (h) h_CaloFittingQA_cemc_etaphi_ZScrosscalib->Add(h);
0114 
0115     h = static_cast<TProfile2D*>(tf->Get("h_CaloFittingQA_ihcal_etaphi_ZScrosscalib"));
0116 
0117     if (h) h_CaloFittingQA_ihcal_etaphi_ZScrosscalib->Add(h);
0118 
0119     h = static_cast<TProfile2D*>(tf->Get("h_CaloFittingQA_ohcal_etaphi_ZScrosscalib"));
0120 
0121     if (h) h_CaloFittingQA_ohcal_etaphi_ZScrosscalib->Add(h);
0122 
0123     tf->Close();
0124   }
0125 
0126   // Close the file
0127   file.close();
0128 
0129   return 0;
0130 }
0131 
0132 void myAnalysis::histToCaloCDBTree(string outputfile, string fieldName, Int_t icalo, TProfile2D *hist)
0133 {
0134   Int_t neta, nphi;
0135 
0136   if (icalo != 0 && icalo != 1) return;
0137 
0138   if (icalo == 0)
0139   {
0140     neta = 96;
0141     nphi = 256;
0142   }
0143   if (icalo == 1)
0144   {
0145     neta = 24;
0146     nphi = 64;
0147   }
0148 
0149   CDBTTree *cdbttree = new CDBTTree(outputfile);
0150 
0151   Float_t mean = 0;
0152   Int_t count = 0;
0153 
0154   for (Int_t ie = 0; ie < neta; ie++)
0155   {
0156     for (Int_t ip = 0; ip < nphi; ip++)
0157     {
0158       UInt_t key;
0159       if (icalo == 0) key = TowerInfoDefs::encode_emcal(ie, ip);
0160       if (icalo == 1) key = TowerInfoDefs::encode_hcal(ie, ip);
0161       Float_t binval = hist->GetBinContent(ie + 1, ip + 1);
0162       Float_t val = 1.0;
0163       if (binval < 1.2 && binval > 0.0) { // clean up of cross calib factors for towers with bit flip behavior 
0164         val = 1.0/binval; // flip cross calib factor for use in CaloTowerCalib
0165       }
0166       cdbttree->SetFloatValue(key, fieldName, val);
0167       mean += val;
0168       count++;
0169     }
0170   }
0171 
0172   cout << "Writing " << outputfile.c_str() << "   with mean=" << mean / count << endl;
0173   cdbttree->Commit();
0174   cdbttree->WriteCDBTTree();
0175   // cdbttree->Print();
0176   delete cdbttree;
0177 }
0178 
0179 void myAnalysis::analyze(const string &output)
0180 {
0181   stringstream t;
0182 
0183   t.str("");
0184 
0185   string run = run_dataset.first;
0186   string dataset = run_dataset.second;
0187 
0188   cout << "Processing: Run: " << run << ", Dataset: " << dataset << endl;
0189 
0190   t << output << "/" << run << "_" << dataset;
0191 
0192   makeDir(output);
0193   makeDir(t.str());
0194 
0195   string detector = "CEMC";
0196   // ZScrosscalib
0197   string payloadName = t.str() + "/" + detector + "_ZSCrossCalib" + "_" + dataset + "_" + run + ".root";
0198   if (h_CaloFittingQA_cemc_etaphi_ZScrosscalib) histToCaloCDBTree(payloadName, "ratio", 0, h_CaloFittingQA_cemc_etaphi_ZScrosscalib);
0199 
0200   detector = "HCALIN";
0201   // ZScrosscalib
0202   payloadName = t.str() + "/" + detector + "_ZSCrossCalib" + "_" + dataset + "_" + run + ".root";
0203   if (h_CaloFittingQA_ihcal_etaphi_ZScrosscalib) histToCaloCDBTree(payloadName, "ratio", 1, h_CaloFittingQA_ihcal_etaphi_ZScrosscalib);
0204 
0205   detector = "HCALOUT";
0206   // ZScrosscalib
0207   payloadName = t.str() + "/" + detector + "_ZSCrossCalib" + "_" + dataset + "_" + run + ".root";
0208   if (h_CaloFittingQA_ohcal_etaphi_ZScrosscalib) histToCaloCDBTree(payloadName, "ratio", 1, h_CaloFittingQA_ohcal_etaphi_ZScrosscalib);
0209 }
0210 
0211 void genFittingStatus(const string &input, const string &output = "output")
0212 {
0213   cout << "#############################" << endl;
0214   cout << "Run Parameters" << endl;
0215   cout << "input: " << input << endl;
0216   cout << "output: " << output << endl;
0217   cout << "#############################" << endl;
0218 
0219   // merges individal qa into one per run
0220   myAnalysis::readHists(input);
0221   myAnalysis::analyze(output);
0222 }
0223 
0224 #ifndef __CINT__
0225 Int_t main(Int_t argc, const char* const argv[])
0226 {
0227   if (argc < 2 || argc > 3)
0228   {
0229     cout << "usage: ./genFittingStatus input [output]" << endl;
0230     cout << "input: input list" << endl;
0231     cout << "output: output directory" << endl;
0232     return 1;
0233   }
0234 
0235   string output = "output";
0236 
0237   if (argc >= 3)
0238   {
0239     output = argv[2];
0240   }
0241 
0242   genFittingStatus(argv[1], output);
0243 
0244   cout << "======================================" << endl;
0245   cout << "done" << endl;
0246   return 0;
0247 }
0248 #endif