Back to home page

sPhenix code displayed by LXR

 
 

    


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

0001 #include <cstdlib>
0002 
0003 float TickSize = 0.03;
0004 float AxisTitleSize = 0.05;
0005 float AxisLabelSize = 0.045;
0006 float LeftMargin = 0.15;
0007 float RightMargin = 0.05;
0008 float TopMargin = 0.08;
0009 float BottomMargin = 0.13;
0010 float textscale = 2.6;
0011 
0012 std::vector<std::string> v_color{"#99cc99", "#6699cc", "#f2777a"};
0013 
0014 void recotkl_matchEff()
0015 {
0016     // directory to save the output plots; if the directory does not exist, create it using
0017     std::string plotdir = "recotkl_matchEff";
0018     system(Form("mkdir -p %s", plotdir.c_str()));
0019 
0020     ROOT::EnableImplicitMT();
0021     ROOT::RDataFrame df("minitree", "/sphenix/tg/tg01/hf/hjheng/ppg02/minitree/TrackletMinitree_Sim_SIMPLE_ana466_20250214/dRcut999p0_NominalVtxZ_RandomClusSet0_clusAdcCutSet0_clusPhiSizeCutSet0/minitree_00*.root");
0022 
0023     auto h_recotklraw_dR = df.Histo1D({"h_recotklraw_dR", "h_recotklraw_dR:Reco-tracklet #DeltaR:Entries", 100, 0, 1}, "recotklraw_dR");
0024     auto h_recotklraw_dPhi = df.Histo1D({"h_recotklraw_dPhi", "h_recotklraw_dPhi:Reco-tracklet #Delta#phi:Entries", 200, -1, 1}, "recotklraw_dphi");
0025     auto h_recotklraw_dEta = df.Histo1D({"h_recotklraw_dEta", "h_recotklraw_dEta:Reco-tracklet #Delta#eta:Entries", 200, -1, 1}, "recotklraw_deta");
0026     auto h_recotklraw_dR_zoomin = df.Histo1D({"h_recotklraw_dR_zoomin", "h_recotklraw_dR_zoomin:Reco-tracklet #DeltaR:Entries", 100, 0, 0.2}, "recotklraw_dR");
0027     auto h_recotklraw_dPhi_zoomin = df.Histo1D({"h_recotklraw_dPhi_zoomin", "h_recotklraw_dPhi_zoomin:Reco-tracklet #Delta#phi:Entries", 200, -0.1, 0.1}, "recotklraw_dphi");
0028     auto h_recotklraw_dEta_zoomin = df.Histo1D({"h_recotklraw_dEta_zoomin", "h_recotklraw_dEta_zoomin:Reco-tracklet #Delta#eta:Entries", 200, -0.2, 0.2}, "recotklraw_deta");
0029     auto h_recotklraw_DCA3d = df.Histo1D({"h_recotklraw_DCA3d", "h_recotklraw_DCA3d:Reco-tracklet DCA (w.r.t event vertex) [cm]:Entries", 100, 0, 5}, "recotklraw_dca3dvtx");
0030     auto h_genhadron_pt = df.Histo1D({"h_genhadron_pt", "h_genhadron_pt:Gen-hadron p_{T}:Entries", 100, 0, 10}, "GenHadron_Pt");
0031     auto h_genhadron_eta = df.Histo1D({"h_genhadron_eta", "h_genhadron_eta:Gen-hadron #eta:Entries", 105, -1.05, 1.05}, "GenHadron_eta");
0032     auto h_genhadron_phi = df.Histo1D({"h_genhadron_phi", "h_genhadron_phi:Gen-hadron #phi:Entries", 128, -3.2, 3.2}, "GenHadron_phi");
0033     // define a new column for the tracklets with the corresponding value in recotkl_isMatchedGChadron is true
0034     auto df_with_matched = df.Define("recotkl_dR_matched_posTrackID",
0035                                      [](const std::vector<float> &recotklraw_dR, const std::vector<int> &recotkl_clus1_matchedtrackID, const std::vector<int> &recotkl_clus2_matchedtrackID)
0036                                      {
0037                                          std::vector<float> recotkl_dR_matched;
0038                                          for (size_t i = 0; i < recotklraw_dR.size(); i++)
0039                                          {
0040                                              if (recotkl_clus1_matchedtrackID[i] == recotkl_clus2_matchedtrackID[i] && recotkl_clus1_matchedtrackID[i] > 0)
0041                                              {
0042                                                  recotkl_dR_matched.push_back(recotklraw_dR[i]);
0043                                              }
0044                                          }
0045                                          return recotkl_dR_matched;
0046                                      },
0047                                      {"recotklraw_dR", "recotkl_clus1_matchedtrackID", "recotkl_clus2_matchedtrackID"})
0048                                .Define("recotkl_dR_matched_negTrackID",
0049                                        [](const std::vector<float> &recotklraw_dR, const std::vector<int> &recotkl_clus1_matchedtrackID, const std::vector<int> &recotkl_clus2_matchedtrackID)
0050                                        {
0051                                            std::vector<float> recotkl_dR_matched;
0052                                            for (size_t i = 0; i < recotklraw_dR.size(); i++)
0053                                            {
0054                                                if (recotkl_clus1_matchedtrackID[i] == recotkl_clus2_matchedtrackID[i] && recotkl_clus1_matchedtrackID[i] < 0)
0055                                                {
0056                                                    recotkl_dR_matched.push_back(recotklraw_dR[i]);
0057                                                }
0058                                            }
0059                                            return recotkl_dR_matched;
0060                                        },
0061                                        {"recotklraw_dR", "recotkl_clus1_matchedtrackID", "recotkl_clus2_matchedtrackID"})
0062                                .Define("recotkl_dR_notmatched",
0063                                        [](const std::vector<float> &recotklraw_dR, const std::vector<int> &recotkl_clus1_matchedtrackID, const std::vector<int> &recotkl_clus2_matchedtrackID)
0064                                        {
0065                                            std::vector<float> recotkl_dR_notmatched;
0066                                            for (size_t i = 0; i < recotklraw_dR.size(); i++)
0067                                            {
0068                                                if (recotkl_clus1_matchedtrackID[i] != recotkl_clus2_matchedtrackID[i])
0069                                                {
0070                                                    recotkl_dR_notmatched.push_back(recotklraw_dR[i]);
0071                                                }
0072                                            }
0073                                            return recotkl_dR_notmatched;
0074                                        },
0075                                        {"recotklraw_dR", "recotkl_clus1_matchedtrackID", "recotkl_clus2_matchedtrackID"})
0076                                .Define("recotkl_dPhi_matched_posTrackID",
0077                                        [](const std::vector<float> &recotklraw_dphi, const std::vector<int> &recotkl_clus1_matchedtrackID, const std::vector<int> &recotkl_clus2_matchedtrackID)
0078                                        {
0079                                            std::vector<float> recotkl_dPhi_matched;
0080                                            for (size_t i = 0; i < recotklraw_dphi.size(); i++)
0081                                            {
0082                                                if (recotkl_clus1_matchedtrackID[i] == recotkl_clus2_matchedtrackID[i] && recotkl_clus1_matchedtrackID[i] > 0)
0083                                                {
0084                                                    recotkl_dPhi_matched.push_back(recotklraw_dphi[i]);
0085                                                }
0086                                            }
0087                                            return recotkl_dPhi_matched;
0088                                        },
0089                                        {"recotklraw_dphi", "recotkl_clus1_matchedtrackID", "recotkl_clus2_matchedtrackID"})
0090                                .Define("recotkl_dPhi_matched_negTrackID",
0091                                        [](const std::vector<float> &recotklraw_dphi, const std::vector<int> &recotkl_clus1_matchedtrackID, const std::vector<int> &recotkl_clus2_matchedtrackID)
0092                                        {
0093                                            std::vector<float> recotkl_dPhi_matched;
0094                                            for (size_t i = 0; i < recotklraw_dphi.size(); i++)
0095                                            {
0096                                                if (recotkl_clus1_matchedtrackID[i] == recotkl_clus2_matchedtrackID[i] && recotkl_clus1_matchedtrackID[i] < 0)
0097                                                {
0098                                                    recotkl_dPhi_matched.push_back(recotklraw_dphi[i]);
0099                                                }
0100                                            }
0101                                            return recotkl_dPhi_matched;
0102                                        },
0103                                        {"recotklraw_dphi", "recotkl_clus1_matchedtrackID", "recotkl_clus2_matchedtrackID"})
0104                                .Define("recotkl_dPhi_notmatched",
0105                                        [](const std::vector<float> &recotklraw_dphi, const std::vector<int> &recotkl_clus1_matchedtrackID, const std::vector<int> &recotkl_clus2_matchedtrackID)
0106                                        {
0107                                            std::vector<float> recotkl_dPhi_notmatched;
0108                                            for (size_t i = 0; i < recotklraw_dphi.size(); i++)
0109                                            {
0110                                                if (recotkl_clus1_matchedtrackID[i] != recotkl_clus2_matchedtrackID[i])
0111                                                {
0112                                                    recotkl_dPhi_notmatched.push_back(recotklraw_dphi[i]);
0113                                                }
0114                                            }
0115                                            return recotkl_dPhi_notmatched;
0116                                        },
0117                                        {"recotklraw_dphi", "recotkl_clus1_matchedtrackID", "recotkl_clus2_matchedtrackID"})
0118                                .Define("recotkl_dEta_matched_posTrackID",
0119                                        [](const std::vector<float> &recotklraw_deta, const std::vector<int> &recotkl_clus1_matchedtrackID, const std::vector<int> &recotkl_clus2_matchedtrackID)
0120                                        {
0121                                            std::vector<float> recotkl_dEta_matched;
0122                                            for (size_t i = 0; i < recotklraw_deta.size(); i++)
0123                                            {
0124                                                if (recotkl_clus1_matchedtrackID[i] == recotkl_clus2_matchedtrackID[i] && recotkl_clus1_matchedtrackID[i] > 0)
0125                                                {
0126                                                    recotkl_dEta_matched.push_back(recotklraw_deta[i]);
0127                                                }
0128                                            }
0129                                            return recotkl_dEta_matched;
0130                                        },
0131                                        {"recotklraw_deta", "recotkl_clus1_matchedtrackID", "recotkl_clus2_matchedtrackID"})
0132                                .Define("recotkl_dEta_matched_negTrackID",
0133                                        [](const std::vector<float> &recotklraw_deta, const std::vector<int> &recotkl_clus1_matchedtrackID, const std::vector<int> &recotkl_clus2_matchedtrackID)
0134                                        {
0135                                            std::vector<float> recotkl_dEta_matched;
0136                                            for (size_t i = 0; i < recotklraw_deta.size(); i++)
0137                                            {
0138                                                if (recotkl_clus1_matchedtrackID[i] == recotkl_clus2_matchedtrackID[i] && recotkl_clus1_matchedtrackID[i] < 0)
0139                                                {
0140                                                    recotkl_dEta_matched.push_back(recotklraw_deta[i]);
0141                                                }
0142                                            }
0143                                            return recotkl_dEta_matched;
0144                                        },
0145                                        {"recotklraw_deta", "recotkl_clus1_matchedtrackID", "recotkl_clus2_matchedtrackID"})
0146                                .Define("recotkl_dEta_notmatched",
0147                                        [](const std::vector<float> &recotklraw_deta, const std::vector<int> &recotkl_clus1_matchedtrackID, const std::vector<int> &recotkl_clus2_matchedtrackID)
0148                                        {
0149                                            std::vector<float> recotkl_dEta_notmatched;
0150                                            for (size_t i = 0; i < recotklraw_deta.size(); i++)
0151                                            {
0152                                                if (recotkl_clus1_matchedtrackID[i] != recotkl_clus2_matchedtrackID[i])
0153                                                {
0154                                                    recotkl_dEta_notmatched.push_back(recotklraw_deta[i]);
0155                                                }
0156                                            }
0157                                            return recotkl_dEta_notmatched;
0158                                        },
0159                                        {"recotklraw_deta", "recotkl_clus1_matchedtrackID", "recotkl_clus2_matchedtrackID"})
0160                                .Define("recotklraw_dca3dvtx_posTrackID",
0161                                        [](const std::vector<float> &recotklraw_dca3dvtx, const std::vector<int> &recotkl_clus1_matchedtrackID, const std::vector<int> &recotkl_clus2_matchedtrackID)
0162                                        {
0163                                            std::vector<float> recotkl_dca3dvtx;
0164                                            for (size_t i = 0; i < recotklraw_dca3dvtx.size(); i++)
0165                                            {
0166                                                if (recotkl_clus1_matchedtrackID[i] == recotkl_clus2_matchedtrackID[i] && recotkl_clus1_matchedtrackID[i] > 0)
0167                                                {
0168                                                    recotkl_dca3dvtx.push_back(recotklraw_dca3dvtx[i]);
0169                                                }
0170                                            }
0171                                            return recotkl_dca3dvtx;
0172                                        },
0173                                        {"recotklraw_dca3dvtx", "recotkl_clus1_matchedtrackID", "recotkl_clus2_matchedtrackID"})
0174                                .Define("recotklraw_dca3dvtx_negTrackID",
0175                                        [](const std::vector<float> &recotklraw_dca3dvtx, const std::vector<int> &recotkl_clus1_matchedtrackID, const std::vector<int> &recotkl_clus2_matchedtrackID)
0176                                        {
0177                                            std::vector<float> recotkl_dca3dvtx;
0178                                            for (size_t i = 0; i < recotklraw_dca3dvtx.size(); i++)
0179                                            {
0180                                                if (recotkl_clus1_matchedtrackID[i] == recotkl_clus2_matchedtrackID[i] && recotkl_clus1_matchedtrackID[i] < 0)
0181                                                {
0182                                                    recotkl_dca3dvtx.push_back(recotklraw_dca3dvtx[i]);
0183                                                }
0184                                            }
0185                                            return recotkl_dca3dvtx;
0186                                        },
0187                                        {"recotklraw_dca3dvtx", "recotkl_clus1_matchedtrackID", "recotkl_clus2_matchedtrackID"})
0188                                .Define("recotklraw_dca3dvtx_notmatched",
0189                                        [](const std::vector<float> &recotklraw_dca3dvtx, const std::vector<int> &recotkl_clus1_matchedtrackID, const std::vector<int> &recotkl_clus2_matchedtrackID)
0190                                        {
0191                                            std::vector<float> recotkl_dca3dvtx;
0192                                            for (size_t i = 0; i < recotklraw_dca3dvtx.size(); i++)
0193                                            {
0194                                                if (recotkl_clus1_matchedtrackID[i] != recotkl_clus2_matchedtrackID[i])
0195                                                {
0196                                                    recotkl_dca3dvtx.push_back(recotklraw_dca3dvtx[i]);
0197                                                }
0198                                            }
0199                                            return recotkl_dca3dvtx;
0200                                        },
0201                                        {"recotklraw_dca3dvtx", "recotkl_clus1_matchedtrackID", "recotkl_clus2_matchedtrackID"})
0202                                .Define("genhadron_pt_matched",
0203                                        [](const std::vector<float> &GenHadron_Pt, const std::vector<bool> &GenHadron_IsRecotkl)
0204                                        {
0205                                            std::vector<float> genhadron_pt_matched;
0206                                            for (size_t i = 0; i < GenHadron_Pt.size(); i++)
0207                                            {
0208                                                if (GenHadron_IsRecotkl[i])
0209                                                {
0210                                                    genhadron_pt_matched.push_back(GenHadron_Pt[i]);
0211                                                }
0212                                            }
0213                                            return genhadron_pt_matched;
0214                                        },
0215                                        {"GenHadron_Pt", "GenHadron_IsRecotkl"})
0216                                .Define("genhadron_eta_matched",
0217                                        [](const std::vector<float> &GenHadron_eta, const std::vector<bool> &GenHadron_IsRecotkl)
0218                                        {
0219                                            std::vector<float> genhadron_eta_matched;
0220                                            for (size_t i = 0; i < GenHadron_eta.size(); i++)
0221                                            {
0222                                                if (GenHadron_IsRecotkl[i])
0223                                                {
0224                                                    genhadron_eta_matched.push_back(GenHadron_eta[i]);
0225                                                }
0226                                            }
0227                                            return genhadron_eta_matched;
0228                                        },
0229                                        {"GenHadron_eta", "GenHadron_IsRecotkl"})
0230                                .Define("genhadron_phi_matched",
0231                                        [](const std::vector<float> &GenHadron_phi, const std::vector<bool> &GenHadron_IsRecotkl)
0232                                        {
0233                                            std::vector<float> genhadron_phi_matched;
0234                                            for (size_t i = 0; i < GenHadron_phi.size(); i++)
0235                                            {
0236                                                if (GenHadron_IsRecotkl[i])
0237                                                {
0238                                                    genhadron_phi_matched.push_back(GenHadron_phi[i]);
0239                                                }
0240                                            }
0241                                            return genhadron_phi_matched;
0242                                        },
0243                                        {"GenHadron_phi", "GenHadron_IsRecotkl"});
0244 
0245     auto h_recotklraw_dR_matched_posTrackID = df_with_matched.Histo1D({"h_recotklraw_dR_matched_posTrackID", "h_recotklraw_dR_matched_posTrackID:Tracklet #DeltaR:Entries", 100, 0, 1}, "recotkl_dR_matched_posTrackID");
0246     auto h_recotklraw_dPhi_matched_posTrackID = df_with_matched.Histo1D({"h_recotklraw_dPhi_matched_posTrackID", "h_recotklraw_dPhi_matched_posTrackID:Tracklet #Delta#phi:Entries", 200, -1, 1}, "recotkl_dPhi_matched_posTrackID");
0247     auto h_recotklraw_dEta_matched_posTrackID = df_with_matched.Histo1D({"h_recotklraw_dEta_matched_posTrackID", "h_recotklraw_dEta_matched_posTrackID:Tracklet #Delta#eta:Entries", 200, -1, 1}, "recotkl_dEta_matched_posTrackID");
0248     auto h_recotklraw_dR_matched_negTrackID = df_with_matched.Histo1D({"h_recotklraw_dR_matched_negTrackID", "h_recotklraw_dR_matched_negTrackID:Tracklet #DeltaR:Entries", 100, 0, 1}, "recotkl_dR_matched_negTrackID");
0249     auto h_recotklraw_dPhi_matched_negTrackID = df_with_matched.Histo1D({"h_recotklraw_dPhi_matched_negTrackID", "h_recotklraw_dPhi_matched_negTrackID:Tracklet #Delta#phi:Entries", 200, -1, 1}, "recotkl_dPhi_matched_negTrackID");
0250     auto h_recotklraw_dEta_matched_negTrackID = df_with_matched.Histo1D({"h_recotklraw_dEta_matched_negTrackID", "h_recotklraw_dEta_matched_negTrackID:Tracklet #Delta#eta:Entries", 200, -1, 1}, "recotkl_dEta_matched_negTrackID");
0251     auto h_recotklraw_dR_notmatched = df_with_matched.Histo1D({"h_recotklraw_dR_notmatched", "h_recotklraw_dR_notmatched:Tracklet #DeltaR:Entries", 100, 0, 1}, "recotkl_dR_notmatched");
0252     auto h_recotklraw_dPhi_notmatched = df_with_matched.Histo1D({"h_recotklraw_dPhi_notmatched", "h_recotklraw_dPhi_notmatched:Tracklet #Delta#phi:Entries", 200, -1, 1}, "recotkl_dPhi_notmatched");
0253     auto h_recotklraw_dEta_notmatched = df_with_matched.Histo1D({"h_recotklraw_dEta_notmatched", "h_recotklraw_dEta_notmatched:Tracklet #Delta#eta:Entries", 200, -1, 1}, "recotkl_dEta_notmatched");
0254     auto h_recotklraw_dR_zoomin_matched_posTrackID = df_with_matched.Histo1D({"h_recotklraw_dR_zoomin_matched_posTrackID", "h_recotklraw_dR_zoomin_matched_posTrackID:Tracklet #DeltaR:Entries", 100, 0, 0.2}, "recotkl_dR_matched_posTrackID");
0255     auto h_recotklraw_dPhi_zoomin_matched_posTrackID = df_with_matched.Histo1D({"h_recotklraw_dPhi_zoomin_matched_posTrackID", "h_recotklraw_dPhi_zoomin_matched_posTrackID:Tracklet #Delta#phi:Entries", 200, -0.1, 0.1}, "recotkl_dPhi_matched_posTrackID");
0256     auto h_recotklraw_dEta_zoomin_matched_posTrackID = df_with_matched.Histo1D({"h_recotklraw_dEta_zoomin_matched_posTrackID", "h_recotklraw_dEta_zoomin_matched_posTrackID:Tracklet #Delta#eta:Entries", 200, -0.2, 0.2}, "recotkl_dEta_matched_posTrackID");
0257     auto h_recotklraw_dR_zoomin_matched_negTrackID = df_with_matched.Histo1D({"h_recotklraw_dR_zoomin_matched_negTrackID", "h_recotklraw_dR_zoomin_matched_negTrackID:Tracklet #DeltaR:Entries", 100, 0, 0.2}, "recotkl_dR_matched_negTrackID");
0258     auto h_recotklraw_dPhi_zoomin_matched_negTrackID = df_with_matched.Histo1D({"h_recotklraw_dPhi_zoomin_matched_negTrackID", "h_recotklraw_dPhi_zoomin_matched_negTrackID:Tracklet #Delta#phi:Entries", 200, -0.1, 0.1}, "recotkl_dPhi_matched_negTrackID");
0259     auto h_recotklraw_dEta_zoomin_matched_negTrackID = df_with_matched.Histo1D({"h_recotklraw_dEta_zoomin_matched_negTrackID", "h_recotklraw_dEta_zoomin_matched_negTrackID:Tracklet #Delta#eta:Entries", 200, -0.2, 0.2}, "recotkl_dEta_matched_negTrackID");
0260     auto h_recotklraw_dR_zoomin_notmatched = df_with_matched.Histo1D({"h_recotklraw_dR_zoomin_notmatched", "h_recotklraw_dR_zoomin_notmatched:Tracklet #DeltaR:Entries", 100, 0, 0.2}, "recotkl_dR_notmatched");
0261     auto h_recotklraw_dPhi_zoomin_notmatched = df_with_matched.Histo1D({"h_recotklraw_dPhi_zoomin_notmatched", "h_recotklraw_dPhi_zoomin_notmatched:Tracklet #Delta#phi:Entries", 200, -0.1, 0.1}, "recotkl_dPhi_notmatched");
0262     auto h_recotklraw_dEta_zoomin_notmatched = df_with_matched.Histo1D({"h_recotklraw_dEta_zoomin_notmatched", "h_recotklraw_dEta_zoomin_notmatched:Tracklet #Delta#eta:Entries", 200, -0.2, 0.2}, "recotkl_dEta_notmatched");
0263     auto h_recotklraw_DCA3d_matched_posTrackID = df_with_matched.Histo1D({"h_recotklraw_DCA3d_matched_posTrackID", "h_recotklraw_DCA3d_matched_posTrackID:Tracklet DCA (w.r.t event vertex) [cm]:Entries", 100, 0, 5}, "recotklraw_dca3dvtx_posTrackID");
0264     auto h_recotklraw_DCA3d_matched_negTrackID = df_with_matched.Histo1D({"h_recotklraw_DCA3d_matched_negTrackID", "h_recotklraw_DCA3d_matched_negTrackID:Tracklet DCA (w.r.t event vertex) [cm]:Entries", 100, 0, 5}, "recotklraw_dca3dvtx_negTrackID");
0265     auto h_recotklraw_DCA3d_notmatched = df_with_matched.Histo1D({"h_recotklraw_DCA3d_notmatched", "h_recotklraw_DCA3d_notmatched:Tracklet DCA (w.r.t event vertex) [cm]:Entries", 100, 0, 5}, "recotklraw_dca3dvtx_notmatched");
0266 
0267     auto h_genhadron_pt_matched = df_with_matched.Histo1D({"h_genhadron_pt_matched", "h_genhadron_pt_matched:Gen-hadron p_{T}:Entries", 100, 0, 10}, "genhadron_pt_matched");
0268     auto h_genhadron_eta_matched = df_with_matched.Histo1D({"h_genhadron_eta_matched", "h_genhadron_eta_matched:Gen-hadron #eta:Entries", 105, -1.05, 1.05}, "genhadron_eta_matched");
0269     auto h_genhadron_phi_matched = df_with_matched.Histo1D({"h_genhadron_phi_matched", "h_genhadron_phi_matched:Gen-hadron #phi:Entries", 128, -3.2, 3.2}, "genhadron_phi_matched");
0270 
0271     // draw the histograms on a canvas; split the canvas into two pads, the top pad to draw the two histograms and the bottom pad to draw the ratio plot
0272     auto drawComparison = [&](TH1 *h_all, const std::vector<TH1 *> &vec_h_stack, bool logy, float ymaxscale, const std::string &xTitle, const std::vector<std::string> &v_leg, const std::string &plotName)
0273     {
0274         // if xTitle contains "tracklet" set leg_all to "All tracklets", if "genhadron" set leg_all to "All gen-hadrons"
0275         std::string leg_all = "";
0276         if (xTitle.find("tracklet") != std::string::npos)
0277         {
0278             leg_all = "All tracklets";
0279         }
0280         else if (xTitle.find("hadron") != std::string::npos)
0281         {
0282             leg_all = "All gen-hadrons";
0283         }
0284 
0285         TCanvas *c = new TCanvas("c", "c", 800, 700);
0286         TPad *pad1 = new TPad("pad1", "pad1", 0, 0.3, 1, 1);
0287         pad1->SetBottomMargin(0);
0288         pad1->SetTopMargin(0.28);
0289         pad1->SetRightMargin(RightMargin);
0290         pad1->Draw();
0291         TPad *pad2 = new TPad("pad2", "pad2", 0, 0, 1, 0.3);
0292         pad2->SetRightMargin(RightMargin);
0293         pad2->SetTopMargin(0);
0294         pad2->SetBottomMargin(0.31);
0295         pad2->Draw();
0296         pad1->cd();
0297         pad1->SetLogy(logy);
0298         float ymin = 1E10;
0299         THStack *hs = new THStack("hs", "");
0300         for (size_t i = 0; i < vec_h_stack.size(); ++i)
0301         {
0302             vec_h_stack[i]->SetLineColor(TColor::GetColor(v_color[i].c_str()));
0303             vec_h_stack[i]->SetFillColorAlpha(TColor::GetColor(v_color[i].c_str()), 0.5);
0304             hs->Add(vec_h_stack[i]);
0305 
0306             if (vec_h_stack[i]->GetMinimum(0) < ymin)
0307             {
0308                 ymin = vec_h_stack[i]->GetMinimum(0);
0309             }
0310         }
0311         hs->Draw();
0312         gPad->Modified();
0313         gPad->Update();
0314         hs->GetXaxis()->SetTitleSize(0);
0315         hs->GetXaxis()->SetLabelSize(0);
0316         // hs->GetXaxis()->SetTitle(xTitle.c_str());
0317         hs->SetMaximum(h_all->GetMaximum() * ymaxscale);
0318         hs->SetMinimum((logy) ? ymin * 0.9 : 0);
0319         hs->GetYaxis()->SetTitle("Counts");
0320         hs->GetYaxis()->SetTitleOffset(1.4);
0321         hs->GetYaxis()->SetTitleSize(AxisTitleSize);
0322         hs->GetYaxis()->SetLabelSize(AxisLabelSize);
0323         gPad->Modified();
0324         gPad->Update();
0325         h_all->Draw("hist same");
0326 
0327         TLegend *leg = new TLegend(pad1->GetLeftMargin(),           //
0328                                    1 - pad1->GetTopMargin() + 0.05, //
0329                                    pad1->GetLeftMargin() + 0.25,     //
0330                                    0.97                             //
0331         );
0332         leg->AddEntry(h_all, leg_all.c_str(), "l");
0333         for (size_t i = 0; i < vec_h_stack.size(); ++i)
0334         {
0335             leg->AddEntry(vec_h_stack[i], v_leg[i].c_str(), "lf");
0336         }
0337         leg->SetBorderSize(0);
0338         leg->SetFillStyle(0);
0339         leg->SetTextSize(0.045);
0340         leg->Draw();
0341 
0342         pad2->cd();
0343         pad2->SetLogy(logy);
0344         float ratio_min = 999.;
0345         THStack *hs_ratio = new THStack("hs_ratio", "");
0346         for (size_t i = 0; i < vec_h_stack.size(); ++i)
0347         {
0348             TH1 *h_ratio = (TH1 *)vec_h_stack[i]->Clone(Form("h_ratio_%zu", i));
0349             h_ratio->Divide(h_all);
0350             h_ratio->SetFillColorAlpha(TColor::GetColor(v_color[i].c_str()), 0.5);
0351             h_ratio->SetLineColor(TColor::GetColor(v_color[i].c_str()));
0352             h_ratio->SetLineWidth(2);
0353             hs_ratio->Add(h_ratio);
0354 
0355             if (h_ratio->GetMinimum(0) < ratio_min)
0356             {
0357                 ratio_min = h_ratio->GetMinimum(0);
0358             }
0359         }
0360         hs_ratio->Draw();
0361         // reset the style of the ratio plot
0362         gPad->Modified();
0363         gPad->Update();
0364         hs_ratio->GetYaxis()->SetTitle("Ratio");
0365         hs_ratio->SetMaximum(1.1);
0366         hs_ratio->SetMinimum((logy) ? ratio_min * 0.9 : 0);
0367         hs_ratio->GetYaxis()->SetTitleSize(AxisTitleSize * textscale);
0368         hs_ratio->GetYaxis()->SetTitleOffset(0.5);
0369         hs_ratio->GetYaxis()->SetLabelSize(AxisLabelSize * textscale);
0370         hs_ratio->GetXaxis()->SetTitle(xTitle.c_str());
0371         hs_ratio->GetXaxis()->SetTitleSize(AxisTitleSize * textscale);
0372         hs_ratio->GetXaxis()->SetTitleOffset(1.1);
0373         hs_ratio->GetXaxis()->SetLabelSize(AxisLabelSize * textscale);
0374         hs_ratio->Draw();
0375         TLine *line = new TLine(hs_ratio->GetXaxis()->GetXmin(), 1, hs_ratio->GetXaxis()->GetXmax(), 1);
0376         line->SetLineStyle(2);
0377         line->Draw();
0378 
0379         c->SaveAs(Form("%s/%s.pdf", plotdir.c_str(), plotName.c_str()));
0380         c->SaveAs(Form("%s/%s.png", plotdir.c_str(), plotName.c_str()));
0381     };
0382 
0383     std::vector<std::string> hleg{"Matched to G4P with positive track ID", "Matched to G4P with negative track ID", "Not matched"};
0384     drawComparison(h_recotklraw_dR.GetPtr(), {h_recotklraw_dR_matched_posTrackID.GetPtr(), h_recotklraw_dR_matched_negTrackID.GetPtr(), h_recotklraw_dR_notmatched.GetPtr()}, true, 1, "Reco-tracklet #DeltaR", hleg, "recotkl_matchEff_dR");
0385     drawComparison(h_recotklraw_dPhi.GetPtr(), {h_recotklraw_dPhi_matched_posTrackID.GetPtr(), h_recotklraw_dPhi_matched_negTrackID.GetPtr(), h_recotklraw_dPhi_notmatched.GetPtr()}, true, 1, "Reco-tracklet #Delta#phi", hleg, "recotkl_matchEff_dPhi");
0386     drawComparison(h_recotklraw_dEta.GetPtr(), {h_recotklraw_dEta_matched_posTrackID.GetPtr(), h_recotklraw_dEta_matched_negTrackID.GetPtr(), h_recotklraw_dEta_notmatched.GetPtr()}, true, 1, "Reco-tracklet #Delta#eta", hleg, "recotkl_matchEff_dEta");
0387     drawComparison(h_recotklraw_dR_zoomin.GetPtr(), {h_recotklraw_dR_zoomin_matched_posTrackID.GetPtr(), h_recotklraw_dR_zoomin_matched_negTrackID.GetPtr(), h_recotklraw_dR_zoomin_notmatched.GetPtr()}, true, 1, "Reco-tracklet #DeltaR", hleg, "recotkl_matchEff_dR_zoomin");
0388     drawComparison(h_recotklraw_dPhi_zoomin.GetPtr(), {h_recotklraw_dPhi_zoomin_matched_posTrackID.GetPtr(), h_recotklraw_dPhi_zoomin_matched_negTrackID.GetPtr(), h_recotklraw_dPhi_zoomin_notmatched.GetPtr()}, true, 1, "Reco-tracklet #Delta#phi", hleg, "recotkl_matchEff_dPhi_zoomin");
0389     drawComparison(h_recotklraw_dEta_zoomin.GetPtr(), {h_recotklraw_dEta_zoomin_matched_posTrackID.GetPtr(), h_recotklraw_dEta_zoomin_matched_negTrackID.GetPtr(), h_recotklraw_dEta_zoomin_notmatched.GetPtr()}, true, 1, "Reco-tracklet #Delta#eta", hleg, "recotkl_matchEff_dEta_zoomin");
0390     drawComparison(h_recotklraw_DCA3d.GetPtr(), {h_recotklraw_DCA3d_matched_posTrackID.GetPtr(), h_recotklraw_DCA3d_matched_negTrackID.GetPtr(), h_recotklraw_DCA3d_notmatched.GetPtr()}, false, 1.1, "Reco-tracklet DCA_{3D} (w.r.t event vertex) [cm]", hleg, "recotkl_matchEff_DCA3d");
0391 
0392     drawComparison(h_genhadron_pt.GetPtr(), {h_genhadron_pt_matched.GetPtr()}, false, 1.03, "Gen-hadron p_{T} [GeV]", {"Matched to reco-tracklets"}, "recotkl_matchEff_genhadron_pt");
0393     drawComparison(h_genhadron_eta.GetPtr(), {h_genhadron_eta_matched.GetPtr()}, false, 1.03, "Gen-hadron #eta", {"Matched to reco-tracklets"}, "recotkl_matchEff_genhadron_eta");
0394     drawComparison(h_genhadron_phi.GetPtr(), {h_genhadron_phi_matched.GetPtr()}, false, 1.03, "Gen-hadron #phi [Radian]", {"Matched to reco-tracklets"}, "recotkl_matchEff_genhadron_phi");
0395 }