Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-06-30 07:55:32

0001 // Code to compare the tracking performances: Truth seeding vs real seeding
0002 // Shyam Kumar; shyam.kumar@ba.infn.it; shyam055119@gmail.com
0003 
0004 #include "TGraphErrors.h"
0005 #include "TF1.h"
0006 #include "TRandom.h"
0007 #include "TCanvas.h"
0008 #include "TLegend.h"
0009 #include "TMath.h"
0010 
0011 void draw_req_DCA(double etamin, double etamax, double xmin=0., double xmax=0.);
0012 void doCompare_truth_real_widebins_dcaT(TString particle = "pi-",double etamin=-1.0, double etamax=1.0, Bool_t drawreq=1, TString extra_legend = "") // name = p, pt for getting p or pt dependence fitted results
0013 {
0014 
0015 //=== style of the plot=========
0016    gStyle->SetPalette(1);
0017    gStyle->SetOptTitle(1);
0018    gStyle->SetTitleOffset(1.0,"XY");
0019    gStyle->SetTitleSize(.04,"XY");
0020    gStyle->SetLabelSize(.04,"XY");
0021    gStyle->SetHistLineWidth(2);
0022    gStyle->SetOptFit(1);
0023    gStyle->SetOptStat(1);
0024   
0025    const Int_t nptbins = 10;
0026    double pt[nptbins] ={0.2, 0.3, 0.5,1.0, 1.5, 2.0, 5.0, 8.0, 10., 15.0};
0027    Double_t variation = 0.1; // 10 % variation 
0028    std::vector<double> momV_truth, momV_real;
0029    std::vector<double> dcaTresolV_truth, err_dcaTresolV_truth, dcaTresolV_real, err_dcaTresolV_real;
0030    momV_truth.clear(); momV_real.clear(); dcaTresolV_truth.clear(); err_dcaTresolV_truth.clear(); dcaTresolV_real.clear(); err_dcaTresolV_real.clear();
0031    
0032    
0033    TString symbolname = "";
0034    if (particle == "pi-") symbolname = "#pi^{-}"; 
0035    else symbolname = particle;
0036 
0037    
0038    TF1 *f1=new TF1("f1","FitPointingAngle",0.,30.0,2);
0039    f1->SetParLimits(0,0.,50000);    
0040    f1->SetParLimits(1,0.,50000);    
0041   
0042    
0043    TCanvas *c_dcaxy = new TCanvas("c_dcaxy","c_dcaxy",1400,1000);
0044    c_dcaxy->SetMargin(0.10, 0.05 ,0.1,0.05);
0045    c_dcaxy->SetGridy();
0046  
0047    //Reading the root file
0048     TFile *fDCA_truth, *fDCA_real;
0049     TGraphErrors *gr_dcaT_truth, *gr_dcaT_real, *gr_dcaZ_truth, *gr_dcaZ_real;
0050      
0051     TMultiGraph *mgDCAT; 
0052     TLegend *lDCAT; 
0053      mgDCAT = new TMultiGraph("mgDCAT",";p_{T} (GeV/c); #sigma_{DCA_{T}} (#mum)");
0054      lDCAT = new TLegend(0.65,0.80,0.90,0.93);
0055      lDCAT->SetTextSize(0.03);
0056      lDCAT->SetBorderSize(0);
0057      lDCAT->SetHeader(extra_legend.Data(), "C");
0058      lDCAT->AddEntry((TObject*)0, Form("%s, %1.1f < #eta < %1.1f", symbolname.Data(), etamin, etamax), "C");
0059       
0060       fDCA_truth = TFile::Open(Form("./truthseed/%s/dca/final_hist_dca_truthseed.root",particle.Data()));
0061       fDCA_real = TFile::Open(Form("./realseed/%s/dca/final_hist_dca_realseed.root",particle.Data()));
0062      
0063      // Truth seeding histograms
0064       TH3D *hist_d0xy_truth = (TH3D*) fDCA_truth->Get("h_d0xy_3d");
0065       TH3D *hist_d0xy_real = (TH3D*) fDCA_real->Get("h_d0xy_3d");
0066       
0067       // d0xy calculation for truth/real (binning are same in both cases)
0068      Int_t etamin_bin = hist_d0xy_truth->GetYaxis()->FindBin(etamin+0.0001);
0069      Int_t etamax_bin = hist_d0xy_truth->GetYaxis()->FindBin(etamax-0.0001);
0070     
0071      TF1 *func_truth = new TF1("func_truth","gaus",-0.5,0.5);
0072      TF1 *func_real = new TF1("func_real","gaus",-0.5,0.5);
0073     
0074     for(int iptbin=0; iptbin<nptbins; ++iptbin){
0075         
0076    TCanvas *cp = new TCanvas("cp","cp",1400,1000);
0077    cp->SetMargin(0.10, 0.05 ,0.1,0.07);
0078      
0079     double ptmin = (1.0-variation)*pt[iptbin]; // 10% range  
0080     double ptmax = (1.0+variation)*pt[iptbin]; // 10% range 
0081   
0082     Int_t ptmin_bin = hist_d0xy_truth->GetZaxis()->FindBin(ptmin+0.0001);
0083     Int_t ptmax_bin = hist_d0xy_truth->GetZaxis()->FindBin(ptmax-0.0001);
0084         
0085     TH1D *histd0xy_truth_1d = (TH1D*)hist_d0xy_truth->ProjectionX(Form("histd0xy_truth_eta%1.1f_%1.1f_pt%1.1f_%1.1f",etamin,etamax,ptmin,ptmax),etamin_bin,etamax_bin,ptmin_bin,ptmax_bin,"o");
0086     histd0xy_truth_1d->SetTitle(Form("d0_{xy} (truth): %1.1f <#eta< %1.1f && %1.2f <p_{T}< %1.2f",etamin,etamax,ptmin,ptmax));   
0087     histd0xy_truth_1d->SetName(Form("eta_%1.1f_%1.1f_d0xy_truth_pt_%1.1f",etamin,etamax,pt[iptbin]));  
0088    
0089   // if (histd0xy_truth_1d->GetEntries()<100) continue;   
0090    double mu_truth = histd0xy_truth_1d->GetMean(); 
0091    double sigma_truth = histd0xy_truth_1d->GetStdDev();
0092    func_truth->SetRange(mu_truth-2.0*sigma_truth,mu_truth+2.0*sigma_truth); // fit with in 2 sigma range
0093    histd0xy_truth_1d->Fit(func_truth,"NR+");
0094    mu_truth = func_truth->GetParameter(2);
0095    sigma_truth = func_truth->GetParError(2);
0096    func_truth->SetRange(mu_truth-2.0*sigma_truth,mu_truth+2.0*sigma_truth); // fit with in 2 sigma range
0097    histd0xy_truth_1d->Fit(func_truth,"R+");
0098    float truth_par2 = func_truth->GetParameter(2)*10000; // cm to mum 10000 factor
0099    float truth_par2_err = func_truth->GetParError(2)*10000;
0100    momV_truth.push_back(pt[iptbin]);
0101    dcaTresolV_truth.push_back(truth_par2);
0102    err_dcaTresolV_truth.push_back(truth_par2_err);
0103  
0104     TH1D *histd0xy_real_1d = (TH1D*)hist_d0xy_real->ProjectionX(Form("histd0xy_real_eta%1.1f_%1.1f_pt%1.1f_%1.1f",etamin,etamax,ptmin,ptmax),etamin_bin,etamax_bin,ptmin_bin,ptmax_bin,"o");
0105     histd0xy_real_1d->SetTitle(Form("d0_{xy} (real): %1.1f <#eta< %1.1f && %1.2f <p_{T}< %1.2f",etamin,etamax,ptmin,ptmax));   
0106     histd0xy_real_1d->SetName(Form("eta_%1.1f_%1.1f_d0xy_real_pt_%1.1f",etamin,etamax,pt[iptbin])); 
0107    
0108   // if (histd0xy_real_1d->GetEntries()<100) continue;    
0109    double mu_real = histd0xy_real_1d->GetMean(); 
0110    double sigma_real = histd0xy_real_1d->GetStdDev();
0111    func_real->SetRange(mu_real-2.0*sigma_real,mu_real+2.0*sigma_real); // fit with in 2 sigma range
0112    histd0xy_real_1d->Fit(func_real,"NR+");
0113    mu_real = func_real->GetParameter(1); 
0114    sigma_real = func_real->GetParameter(2);
0115    func_real->SetRange(mu_real-2.0*sigma_real,mu_real+2.0*sigma_real); // fit with in 2 sigma range
0116    histd0xy_real_1d->Fit(func_real,"R+");
0117    float real_par2 = func_real->GetParameter(2)*10000; // cm to mum 10000 factor
0118    float real_par2_err = func_real->GetParError(2)*10000;
0119    momV_real.push_back(pt[iptbin]);
0120    dcaTresolV_real.push_back(real_par2);
0121    err_dcaTresolV_real.push_back(real_par2_err);
0122    
0123    cp->cd();
0124    histd0xy_truth_1d->Draw();
0125    cp->SaveAs(Form("Debug_Plots/truth/%s/dca/truth_dcaxy_resol_mom%1.1f_%1.1f_eta_%1.1f.png",particle.Data(),pt[iptbin],etamin,etamax));
0126    cp->Clear();
0127    cp->cd();
0128    histd0xy_real_1d->Draw();
0129    cp->SaveAs(Form("Debug_Plots/real/%s/dca/real_dcaxy_resol_mom%1.1f_%1.1f_eta_%1.1f.png",particle.Data(),pt[iptbin],etamin,etamax));
0130    }   // ptbin
0131       
0132     const int size_truth = momV_truth.size();
0133     double pt_truth[size_truth], err_pt_truth[size_truth], sigma_dcaxy_truth[size_truth], err_sigma_dcaxy_truth[size_truth]; 
0134     
0135     for (int i=0; i<size_truth; i++){
0136     pt_truth[i] = momV_truth.at(i);
0137     sigma_dcaxy_truth[i] = dcaTresolV_truth.at(i);
0138     err_sigma_dcaxy_truth[i] = err_dcaTresolV_truth.at(i);
0139     err_pt_truth[i] = 0.;
0140     }
0141     
0142      const int size_real = momV_real.size();
0143     double pt_real[size_real], err_pt_real[size_real], sigma_dcaxy_real[size_real], err_sigma_dcaxy_real[size_real]; 
0144     
0145     for (int i=0; i<size_real; i++){
0146     pt_real[i] = momV_real.at(i);
0147      sigma_dcaxy_real[i] = dcaTresolV_real.at(i);
0148     err_sigma_dcaxy_real[i] = err_dcaTresolV_real.at(i);
0149     err_pt_real[i] = 0.;
0150     }
0151 
0152      TFile *fout = new TFile(Form("Final_Results/%s/dca/dcaxy_resol_%1.1f_eta_%1.1f.root",particle.Data(),etamin,etamax),"recreate");
0153     TGraphErrors *gr1 = new TGraphErrors(size_truth,pt_truth,sigma_dcaxy_truth,err_pt_truth,err_sigma_dcaxy_truth);
0154      gr1->SetName("gr_truthseed");
0155     gr1->SetMarkerStyle(25);
0156     gr1->SetMarkerColor(kMagenta);
0157     gr1->SetMarkerSize(2.0);
0158     gr1->SetTitle(";p_{T} (GeV/c); #sigma_{DCA_{T}} (#mum)");
0159     gr1->GetXaxis()->CenterTitle();
0160     gr1->GetYaxis()->CenterTitle();
0161     
0162      TGraphErrors *gr2 = new TGraphErrors(size_real,pt_real,sigma_dcaxy_real,err_pt_real,err_sigma_dcaxy_real);
0163      gr2->SetName("gr_realseed");
0164     gr2->SetMarkerStyle(34);
0165     gr2->SetMarkerColor(kRed);
0166     gr2->SetMarkerSize(2.0);
0167     gr2->SetTitle(";p_{T} (GeV/c); #sigma_{DCA_{T}} (#mum)");
0168     gr2->GetXaxis()->CenterTitle();
0169     gr2->GetYaxis()->CenterTitle();
0170     
0171 //  mgDCAT->Add(gr1);
0172     mgDCAT->Add(gr2);
0173     c_dcaxy->cd();
0174 //  c_dcaxy->SetLogy();
0175      mgDCAT->GetXaxis()->SetRangeUser(0.18,mgDCAT->GetXaxis()->GetXmax()); 
0176        mgDCAT->GetYaxis()->SetRangeUser(0.0,1.50*TMath::MaxElement(gr2->GetN(),gr2->GetY())); 
0177      mgDCAT->Draw("AP");
0178     // lDCAT->AddEntry(gr1,"Truth Seeding");
0179      lDCAT->AddEntry(gr2,"Realistic Seeding");
0180      lDCAT->Draw("same");
0181      draw_req_DCA(etamin,etamax,0.,mgDCAT->GetXaxis()->GetXmax());
0182      c_dcaxy->SaveAs(Form("Final_Results/%s/dca/dcaxy_resol_%1.1f_eta_%1.1f.png",particle.Data(),etamin,etamax));
0183  
0184      fout->cd();
0185      mgDCAT->SetName(Form("dcaxy_resol_%1.1f_eta_%1.1f",etamin,etamax));
0186      mgDCAT->Write();
0187      fout->Close();
0188 }
0189 
0190 //From Yellow report from section 11.2.2
0191 //===Fit Pointing Resolution
0192 float FitPointingAngle(Double_t *x, Double_t *par)
0193 {
0194   float func = sqrt((par[0]*par[0])/(x[0]*x[0])+par[1]*par[1]);
0195   return func;
0196 }
0197 
0198 void draw_req_DCA(double etamin, double etamax, double xmin=0., double xmax=0.)
0199 {
0200   TF1 *PWGReq_DCA2D;
0201   if (etamin >= -3.5 && etamax <= -2.5) PWGReq_DCA2D = new TF1("PWGReq_DCA2D", "TMath::Sqrt((30./x)^2+40.0^2)",xmin,xmax);
0202   else if (etamin >= -2.5 && etamax <= -1.0) PWGReq_DCA2D = new TF1("PWGReq_DCA2D", "TMath::Sqrt((30./x)^2+20.0^2)",xmin,xmax);
0203   else if (etamin >= -1.0 && etamax <= 1.0) PWGReq_DCA2D = new TF1("PWGReq_DCA2D", "TMath::Sqrt((20./x)^2+5.0^2)",xmin,xmax);
0204   else if (etamin >= 1.0 && etamax <= 2.5) PWGReq_DCA2D = new TF1("PWGReq_DCA2D", "TMath::Sqrt((30./x)^2+20.0^2)",xmin,xmax);
0205   else if (etamin >= 2.5 && etamax <= 3.5) PWGReq_DCA2D = new TF1("PWGReq_DCA2D", "TMath::Sqrt((30./x)^2+40.0^2)",xmin,xmax);
0206   else return;
0207   PWGReq_DCA2D->SetLineStyle(7);
0208   PWGReq_DCA2D->SetLineWidth(3.0);
0209   PWGReq_DCA2D->SetLineColor(kBlue);
0210   PWGReq_DCA2D->Draw("same");
0211         
0212   TLegend *l= new TLegend(0.70,0.75,0.90,0.80);
0213   l->SetTextSize(0.03);
0214   l->SetBorderSize(0);
0215   l->AddEntry(PWGReq_DCA2D,"PWGReq","l");
0216   l->Draw("same");
0217 }