Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-03-17 11:22:39

0001 #include "StackValidation.hh"
0002 
0003 StackValidation::StackValidation(const TString& label,
0004                                  const TString& extra,
0005                                  const Bool_t cmsswComp,
0006                                  const TString& suite)
0007     : label(label), extra(extra), cmsswComp(cmsswComp), suite(suite) {
0008   // setup style for plotting
0009   setupStyle();
0010 
0011   // setup suite enum
0012   setupSUITEEnum(suite);
0013 
0014   // setup build options : true for isBenchmark-type plots, and include cmssw in plots if simval
0015   setupBuilds(false, !cmsswComp);
0016 
0017   // set legend y height
0018   y1 = 0.80;
0019   y2 = y1 + nbuilds * 0.04;  // full + CMSSW = 0.04*5 + 0.8 = 1.0
0020 
0021   // open the files
0022   files.resize(nbuilds);
0023   for (auto b = 0U; b < nbuilds; b++) {
0024     const auto& build = builds[b];
0025     auto& file = files[b];
0026 
0027     file = TFile::Open("validation_" + label + "_" + build.name + extra + "/plots.root");
0028   }
0029 
0030   // setup ref
0031   setupRef(cmsswComp);
0032 
0033   // setup rates
0034   setupRates(cmsswComp);
0035 
0036   // setup ptcuts
0037   setupPtCuts();
0038 }
0039 
0040 StackValidation::~StackValidation() {
0041   for (auto& file : files)
0042     delete file;
0043 }
0044 
0045 void StackValidation::MakeValidationStacks() {
0046   StackValidation::MakeRatioStacks("build");
0047   StackValidation::MakeKinematicDiffStacks("build");
0048   StackValidation::MakeQualityStacks("build");
0049 
0050   if (cmsswComp) {
0051     StackValidation::MakeRatioStacks("fit");
0052     StackValidation::MakeKinematicDiffStacks("fit");
0053     StackValidation::MakeQualityStacks("fit");
0054   }
0055 }
0056 
0057 void StackValidation::MakeRatioStacks(const TString& trk) {
0058   // kinematic variables to plot
0059   std::vector<TString> vars = {"pt", "eta", "phi", "nLayers"};
0060   const UInt_t nvars = vars.size();
0061 
0062   // indices for loops match PlotValidation.cpp
0063   for (auto l = 0U; l < nrates; l++) {
0064     const auto& rate = rates[l];
0065 
0066     for (auto k = 0U; k < nptcuts; k++) {
0067       const auto& ptcut = ptcuts[k];
0068 
0069       for (auto i = 0U; i < nvars; i++) {
0070         const auto& var = vars[i];
0071 
0072         auto canv = new TCanvas();
0073         canv->cd();
0074 
0075         auto leg = new TLegend(0.85, y1, 1.0, y2);
0076 
0077         // tmp axis titles, not sure why ROOT is deleting them
0078         TString xtitle = "";
0079         TString ytitle = "";
0080 
0081         std::vector<TGraphAsymmErrors*> graphs(nbuilds);
0082 
0083         for (auto b = 0U; b < nbuilds; b++) {
0084           const auto& build = builds[b];
0085           auto& file = files[b];
0086           auto& graph = graphs[b];
0087 
0088           graph = ((TEfficiency*)file->Get(rate.dir + refdir + "/" + rate.rate + "_" + rate.sORr + "_" + var + "_" +
0089                                            trk + "_pt" + ptcut))
0090                       ->CreateGraph();
0091           graph->SetLineColor(build.color);
0092           graph->SetMarkerColor(build.color);
0093 
0094           // store tmp titles
0095           if (b == 0) {
0096             xtitle = graph->GetXaxis()->GetTitle();
0097             ytitle = graph->GetYaxis()->GetTitle();
0098           }
0099 
0100           graph->Draw(b > 0 ? "PZ SAME" : "APZ");
0101 
0102           if (!rate.rate.Contains("ineff", TString::kExact) && !rate.rate.Contains("dr", TString::kExact))
0103             graph->GetYaxis()->SetRangeUser(0.0, 1.05);
0104           else
0105             graph->GetYaxis()->SetRangeUser(0.0, 0.25);
0106 
0107           leg->AddEntry(graph, build.label.Data(), "LEP");
0108         }
0109 
0110         // print standard plot for every rate/variable
0111         leg->Draw("SAME");
0112         canv->SaveAs(label + "_" + rate.rate + "_" + var + "_" + trk + "_pt" + ptcut + extra + ".png");
0113 
0114         // zoom in on pt range
0115         if (i == 0) {
0116           std::vector<TGraphAsymmErrors*> zoomgraphs(nbuilds);
0117           for (auto b = 0U; b < nbuilds; b++) {
0118             auto& graph = graphs[b];
0119             auto& zoomgraph = zoomgraphs[b];
0120 
0121             zoomgraph = (TGraphAsymmErrors*)graph->Clone(Form("%s_zoom", graph->GetName()));
0122             zoomgraph->GetXaxis()->SetRangeUser(0, 10);
0123             zoomgraph->Draw(b > 0 ? "PZ SAME" : "APZ");
0124           }
0125 
0126           leg->Draw("SAME");
0127           canv->SaveAs(label + "_" + rate.rate + "_" + var + "_zoom_" + trk + "_pt" + ptcut + extra + ".png");
0128 
0129           for (auto& zoomgraph : zoomgraphs)
0130             delete zoomgraph;
0131         }
0132 
0133         // make logx plots for pt: causes LOTS of weird effects... workarounds for now
0134         if (i == 0) {
0135           canv->SetLogx(1);
0136 
0137           // apparently logx removes titles and ranges???
0138           for (auto b = 0U; b < nbuilds; b++) {
0139             auto& graph = graphs[b];
0140             graph->GetXaxis()->SetRangeUser(0.01, graph->GetXaxis()->GetBinUpEdge(graph->GetXaxis()->GetNbins()));
0141 
0142             if (!rate.rate.Contains("ineff", TString::kExact) && !rate.rate.Contains("dr", TString::kExact))
0143               graph->GetYaxis()->SetRangeUser(0.0, 1.05);
0144             else
0145               graph->GetYaxis()->SetRangeUser(0.0, 0.25);
0146 
0147             graph->GetXaxis()->SetTitle(xtitle);
0148             graph->GetYaxis()->SetTitle(ytitle);
0149 
0150             graph->Draw(b > 0 ? "PZ SAME" : "APZ");
0151           }
0152 
0153           leg->Draw("SAME");
0154           canv->SaveAs(label + "_" + rate.rate + "_" + var + "_logx_" + trk + "_pt" + ptcut + extra + ".png");
0155         }
0156 
0157         delete leg;
0158         for (auto& graph : graphs)
0159           delete graph;
0160         delete canv;
0161       }
0162     }
0163   }
0164 }
0165 
0166 void StackValidation::MakeKinematicDiffStacks(const TString& trk) {
0167   // variables to plot
0168   std::vector<TString> diffs = {"nHits", "invpt", "eta", "phi"};
0169   const UInt_t ndiffs = diffs.size();
0170 
0171   // diffferent reco collections
0172   std::vector<TString> colls = {"allmatch", "bestmatch"};
0173   const UInt_t ncolls = colls.size();
0174 
0175   // indices for loops match PlotValidation.cpp
0176   for (auto o = 0U; o < ncolls; o++) {
0177     const auto& coll = colls[o];
0178 
0179     for (auto p = 0U; p < ndiffs; p++) {
0180       const auto& diff = diffs[p];
0181 
0182       for (auto k = 0U; k < nptcuts; k++) {
0183         const auto& ptcut = ptcuts[k];
0184 
0185         const Bool_t isLogy = true;
0186         auto canv = new TCanvas();
0187         canv->cd();
0188         canv->SetLogy(isLogy);
0189 
0190         auto leg = new TLegend(0.85, y1, 1.0, y2);
0191 
0192         // tmp min/max
0193         Double_t min = 1e9;
0194         Double_t max = -1e9;
0195 
0196         std::vector<TH1F*> hists(nbuilds);
0197         for (auto b = 0U; b < nbuilds; b++) {
0198           const auto& build = builds[b];
0199           auto& file = files[b];
0200           auto& hist = hists[b];
0201 
0202           hist = (TH1F*)file->Get("kindiffs" + refdir + "/h_d" + diff + "_" + coll + "_" + trk + "_pt" + ptcut);
0203           hist->SetLineColor(build.color);
0204           hist->SetMarkerColor(build.color);
0205 
0206           hist->Scale(1.f / hist->Integral());
0207           hist->GetYaxis()->SetTitle("Fraction of Tracks");
0208 
0209           GetMinMaxHist(hist, min, max);
0210         }
0211 
0212         for (auto b = 0U; b < nbuilds; b++) {
0213           const auto& build = builds[b];
0214           auto& hist = hists[b];
0215 
0216           SetMinMaxHist(hist, min, max, isLogy);
0217           hist->Draw(b > 0 ? "EP SAME" : "EP");
0218 
0219           const TString mean = Form("%4.1f", hist->GetMean());
0220           leg->AddEntry(hist, build.label + " " + " [#mu = " + mean + "]", "LEP");
0221         }
0222 
0223         leg->Draw("SAME");
0224         canv->SaveAs(label + "_" + coll + "_d" + diff + "_" + trk + "_pt" + ptcut + extra + ".png");
0225 
0226         delete leg;
0227         for (auto& hist : hists)
0228           delete hist;
0229         delete canv;
0230       }  // end pt cut loop
0231     }    // end var loop
0232   }      // end coll loop
0233 }
0234 
0235 void StackValidation::MakeQualityStacks(const TString& trk) {
0236   // diffferent reco collections
0237   std::vector<TString> colls = {"allreco", "fake", "allmatch", "bestmatch"};
0238   const UInt_t ncolls = colls.size();
0239 
0240   // quality plots to use: nHits/track and track score
0241   std::vector<TString> quals = {"nHits", "score"};
0242   const UInt_t nquals = quals.size();
0243 
0244   // indices for loops match PlotValidation.cpp
0245   for (auto o = 0U; o < ncolls; o++) {
0246     const auto& coll = colls[o];
0247 
0248     for (auto k = 0U; k < nptcuts; k++) {
0249       const auto& ptcut = ptcuts[k];
0250 
0251       for (auto n = 0U; n < nquals; n++) {
0252         const auto& qual = quals[n];
0253 
0254         const Bool_t isLogy = true;
0255         auto canv = new TCanvas();
0256         canv->cd();
0257         canv->SetLogy(isLogy);
0258 
0259         auto leg = new TLegend(0.85, y1, 1.0, y2);
0260 
0261         // tmp min/max
0262         Double_t min = 1e9;
0263         Double_t max = -1e9;
0264 
0265         std::vector<TH1F*> hists(nbuilds);
0266         for (auto b = 0U; b < nbuilds; b++) {
0267           const auto& build = builds[b];
0268           auto& file = files[b];
0269           auto& hist = hists[b];
0270 
0271           hist = (TH1F*)file->Get("quality" + refdir + "/h_" + qual + "_" + coll + "_" + trk + "_pt" + ptcut);
0272           hist->SetLineColor(build.color);
0273           hist->SetMarkerColor(build.color);
0274 
0275           hist->Scale(1.f / hist->Integral());
0276           hist->GetYaxis()->SetTitle("Fraction of Tracks");
0277 
0278           GetMinMaxHist(hist, min, max);
0279         }
0280 
0281         for (auto b = 0U; b < nbuilds; b++) {
0282           const auto& build = builds[b];
0283           auto& hist = hists[b];
0284 
0285           SetMinMaxHist(hist, min, max, isLogy);
0286           hist->Draw(b > 0 ? "EP SAME" : "EP");
0287 
0288           const TString mean = Form("%4.1f", hist->GetMean());
0289           leg->AddEntry(hist, build.label + " " + " [#mu = " + mean + "]", "LEP");
0290         }
0291 
0292         leg->Draw("SAME");
0293         canv->SaveAs(label + "_" + coll + "_" + qual + "_" + trk + "_pt" + ptcut + extra + ".png");
0294 
0295         delete leg;
0296         for (auto& hist : hists)
0297           delete hist;
0298         delete canv;
0299 
0300       }  // end loop over quality variable
0301     }    // end pt cut loop
0302   }      // end coll loop
0303 }