Back to home page

sPhenix code displayed by LXR

 
 

    


File indexing completed on 2025-08-05 08:15:16

0001 #include "../CommonTools.h"
0002 
0003 void initializeModules(std::map<std::string,std::vector<double>>& modules);
0004 void collectData(std::map<std::string,std::vector<double>>& modules, 
0005          std::string& filename);
0006 void plotTimes(std::map<std::string, std::vector<double>>& modules, std::string& filename);
0007 
0008 
0009 
0010 /*
0011  * This root macro analyzes the output from parseTimers.csh
0012  * The input is a filename to the text file created from parseTimers.csh, e.g.
0013  * root -l -b -q AnalyzeTime.C'("time.txt")'
0014  * The output is a root file which contains the time of each tracking module
0015  * The output also produces a .png with each module's individual contribution
0016  * to the total time. The output rootfile has a total time graph which
0017  * sums the contributions
0018  */
0019 void AnalyzeTime(std::string infile, std::string outfile)
0020 {
0021 
0022   std::map<std::string,std::vector<double>> modules;
0023   initializeModules(modules);
0024  
0025   collectData(modules,infile);
0026     
0027   plotTimes(modules, outfile);
0028 }
0029 
0030 void plotTimes(std::map<std::string, std::vector<double>>& modules,
0031            std::string& filename)
0032 {
0033   ostringstream canname;
0034   int nbins = 200;
0035   double bins[nbins+1];
0036   double stepsize = 1.07;
0037   bins[0] = 1.;
0038 
0039   for(int i=1;i <nbins+1; i++)
0040     {
0041       bins[i]=bins[i-1]*stepsize;
0042     }
0043 
0044 
0045   std::vector<double> totalTime;
0046   const int totentries = modules.at("MvtxClusterizer").size();
0047  
0048   for(int i=0; i<totentries; i++)
0049     { totalTime.push_back(0); }
0050 
0051   TFile *outfile = new TFile(filename.c_str(),"recreate");
0052   for(const auto& [moduleName, values] : modules)
0053     {
0054       canname.str("");
0055       canname<<"can_"<<moduleName.c_str();
0056       TCanvas *can = new TCanvas(canname.str().c_str(),canname.str().c_str(),
0057                  200,200,800,600);
0058       //gStyle->SetOptStat(0);
0059       TH1F *histo = new TH1F(moduleName.c_str(),moduleName.c_str(),
0060                  nbins,bins);
0061       
0062       for(int i=0; i<values.size(); i++)
0063     {
0064       histo->Fill(values.at(i));
0065       if(i<totentries)
0066         totalTime.at(i) += values.at(i);
0067     }
0068 
0069       histo->GetXaxis()->SetTitle("time [ms]");
0070       can->cd();
0071  
0072       gPad->SetLogx();
0073       
0074       histo->Draw();
0075       outfile->cd();
0076       histo->Write();
0077       myText(0.2,0.96,kBlack,moduleName.c_str());
0078       canname<<".png";
0079       //can->Print(canname.str().c_str());
0080     }
0081 
0082   TH1F *totalTimes = new TH1F("totalTime",";time [ms]",nbins,bins);
0083   for(const auto& val : totalTime)
0084     {  totalTimes->Fill(val); }
0085 
0086   totalTimes->Write();
0087   outfile->Write();
0088   outfile->Close();
0089 
0090 }
0091 
0092 void collectData(std::map<std::string,std::vector<double>>& modules, 
0093          std::string& filename)
0094 {
0095 
0096   ifstream file;
0097   file.open(filename.c_str());
0098   std::string line;
0099 
0100   if(file.is_open()) {
0101     while(getline(file,line)) {
0102       for(auto& [module, values] : modules) {
0103     if(line.find(module) != std::string::npos) {
0104       if(line.find("TOP") != std::string::npos) {
0105         const auto pos = line.find("(ms):  ");
0106         const auto val = line.substr(pos+7,11);
0107         values.push_back(atof(val.c_str()));
0108       }
0109       else if(line.find("Acts seed time") != std::string::npos) {
0110         const auto apos = line.find("time ");
0111         const auto aval = line.substr(apos+5,6);
0112         values.push_back(atof(aval.c_str()));
0113       }
0114     }
0115       }
0116     }
0117   }  
0118 }
0119 
0120 
0121 void initializeModules(std::map<std::string,std::vector<double>>& modules)
0122 {
0123   std::vector<double> mvtx, intt, tpc, tpclean, silseed, caseed, preprop, kfprop, circlefit, clusterMover, tpcresiduals, secondtrkfit, trackmatch, trackfit, vtxfinder, vertexprop, actsseed, mm, mmmatch, dzcor;
0124   modules.insert(std::make_pair("MvtxClusterizer",mvtx));
0125   modules.insert(std::make_pair("InttClusterizer",intt));
0126   modules.insert(std::make_pair("TpcClusterizer",tpc));
0127   modules.insert(std::make_pair("MicromegasClusterizer",mm));
0128   modules.insert(std::make_pair("PHMicromegasTpcTrackMatching",mmmatch));
0129   modules.insert(std::make_pair("TpcClusterCleaner",tpclean));
0130   modules.insert(std::make_pair("PHActsSiliconSeeding",silseed));
0131   modules.insert(std::make_pair("Acts seed time",actsseed));
0132   modules.insert(std::make_pair("PHCASeeding",caseed));
0133   modules.insert(std::make_pair("PrePropagatorPHTpcTrackSeedCircleFit",preprop));
0134   modules.insert(std::make_pair("PHSimpleKFProp",kfprop));
0135   modules.insert(std::make_pair("PHTpcTrackSeedCircleFit",circlefit));
0136   modules.insert(std::make_pair("PHSiliconTpcTrackMatching",trackmatch));
0137   modules.insert(std::make_pair("PHActsFirstTrkFitter",trackfit));
0138   modules.insert(std::make_pair("PHSimpleVertexFinder",vtxfinder));
0139   modules.insert(std::make_pair("PHActsVertexPropagator",vertexprop));
0140   modules.insert(std::make_pair("PHTpcClusterMover",clusterMover));
0141   modules.insert(std::make_pair("PHTpcResiduals",tpcresiduals));
0142   modules.insert(std::make_pair("PHActsSecondTrkFitter",secondtrkfit));
0143   modules.insert(std::make_pair("PHTpcDeltaZCorrection",dzcor));
0144 }
0145