]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGPP/analysisQA/processFemtoQA.C
New macro for Femto QA: sjena
[u/mrichter/AliRoot.git] / PWGPP / analysisQA / processFemtoQA.C
1 /***************************************************************
2  processFemtoQA.C
3  Post Processing of Femto QA task in Analysis QA train
4  Author: Maciej Szymanski, maszyman@cern.ch
5 ***************************************************************/
6
7 enum collidingSystem {
8   PbPb,
9   pPb,
10   pp
11 };
12
13 Double_t calculateNormalizationFactor(TH1D *num,TH1D *den, Double_t qlo=0.3,Double_t qhi=0.4)
14 {
15   Double_t binlo = num->GetXaxis()->FindFixBin(qlo);
16   Double_t binhi = num->GetXaxis()->FindFixBin(qhi);
17   Double_t integralNum = num->Integral(binlo, binhi);
18   Double_t integralDen = den->Integral(binlo, binhi);
19   return integralDen / integralNum;
20 }
21
22 void processFemtoQA(const char *filePath = "AnalysisResults.root",
23                     const char *listname = "femtolist",
24                     const char *suffix = "png",
25                     enum collidingSystem system = PbPb) {
26
27   TFile *_file = new TFile(filePath,"read");
28   TList* _femtolist = (TList*)_file->Get(Form("PWG2FEMTO/%s",listname));
29
30   // 1D pion correlation function low kT
31   TH1D* numCFlowkT = (TH1D*)_femtolist->FindObject("NumcqinvpimtpcM0kT0");
32   TH1D* denCFlowkT = (TH1D*)_femtolist->FindObject("DencqinvpimtpcM0kT0");
33   Double_t norm = calculateNormalizationFactor(numCFlowkT,denCFlowkT,0.4,0.45 );
34   numCFlowkT->Divide(denCFlowkT);
35   numCFlowkT->Scale(norm);
36   numCFlowkT->SetXTitle("q_{inv} (GeV/c)");
37   numCFlowkT->SetYTitle("C(q_{inv})");
38   numCFlowkT->GetXaxis()->SetRangeUser(0,0.5);
39   numCFlowkT->GetYaxis()->SetRangeUser(0.5,2.5);
40
41   // 1D pion correlation function high kT
42   TH1D* numCFhighkT = (TH1D*)_femtolist->FindObject("NumcqinvpimtpcM0kT3");
43   TH1D* denCFhighkt = (TH1D*)_femtolist->FindObject("DencqinvpimtpcM0kT3");
44   Double_t norm = calculateNormalizationFactor(numCFhighkT,denCFhighkt,0.4,0.45 );
45   numCFhighkT->Divide(denCFhighkt);
46   numCFhighkT->Scale(norm);
47   numCFhighkT->SetXTitle("q_{inv} (GeV/c)");
48   numCFhighkT->SetYTitle("C(q_{inv})");
49   numCFhighkT->SetTitle(Form());
50   numCFhighkT->GetXaxis()->SetRangeUser(0,0.5);
51   numCFhighkT->GetYaxis()->SetRangeUser(0.5,2.5);
52
53   // delta eta - delta phi* low kT
54   TH2D* numPhiEtalowkT = (TH2D*)_femtolist->FindObject("NumRadDPhistarEtapimtpcM0kT0");
55   TH2D* denPhiEtalowkT = (TH2D*)_femtolist->FindObject("DenRadDPhistarEtapimtpcM0kT0");
56   numPhiEtalowkT->Divide(denPhiEtalowkT);
57   numPhiEtalowkT->SetXTitle("#Delta #phi*");
58   numPhiEtalowkT->SetYTitle("#Delta #eta");
59   numPhiEtalowkT->SetTitle(Form());
60
61   // delta eta - delta phi* high kT
62   TH2D* numPhiEtahighkT = (TH2D*)_femtolist->FindObject("NumRadDPhistarEtapimtpcM0kT3");
63   TH2D* denPhiEtahighkt = (TH2D*)_femtolist->FindObject("DenRadDPhistarEtapimtpcM0kT3");
64   numPhiEtahighkT->Divide(denPhiEtahighkt);
65   numPhiEtahighkT->SetXTitle("#Delta #phi*");
66   numPhiEtahighkT->SetYTitle("#Delta #eta");
67   numPhiEtahighkT->SetTitle(Form());
68
69   // qinv vs. separation at TPC entrance low kT
70   TH2D* numQinvEtpclowkT = (TH2D*)_femtolist->FindObject("NumDTPCPhistarEtapimtpcM0kT0");
71   TH2D* denQinvEtpclowkT = (TH2D*)_femtolist->FindObject("DenDTPCPhistarEtapimtpcM0kT0");
72   numQinvEtpclowkT->Divide(denQinvEtpclowkT);
73   numQinvEtpclowkT->SetXTitle("q_{inv} (GeV/c)");
74   numQinvEtpclowkT->SetYTitle("separation at TPC entrance");
75   numQinvEtpclowkT->SetTitle(Form());
76
77   // qinv vs. separation at TPC entrance high kT
78   TH2D* numQinvEtpchighkT = (TH2D*)_femtolist->FindObject("NumDTPCPhistarEtapimtpcM0kT3");
79   TH2D* denQinvEtpchighkt = (TH2D*)_femtolist->FindObject("DenDTPCPhistarEtapimtpcM0kT3");
80   numQinvEtpchighkT->Divide(denQinvEtpchighkt);
81   numQinvEtpchighkT->SetXTitle("q_{inv} (GeV/c)");
82   numQinvEtpchighkT->SetYTitle("separation at TPC entrance");
83   numQinvEtpchighkT->SetTitle(Form());
84
85   if ( system == PbPb ) {
86     numCFlowkT->SetTitle("#pi^{-}#pi^{-} 0-10%, 0.2 < k_{T} < 0.3 GeV/c");
87     numCFhighkT->SetTitle("#pi^{-}#pi^{-} 0-10%, 0.6 < k_{T} < 0.7 GeV/c");
88   }
89   else if ( system == pPb ) {
90     numCFlowkT->SetTitle("#pi^{-}#pi^{-} 0-20%, 0.2 < k_{T} < 0.3 GeV/c");
91     numCFhighkT->SetTitle("#pi^{-}#pi^{-} 0-20%, 0.6 < k_{T} < 0.7 GeV/c");
92   }
93   else if ( system == pp ) {
94     numCFlowkT->SetTitle("#pi^{-}#pi^{-} N_{ch} 50-150, 0.2 < k_{T} < 0.3 GeV/c");
95     numCFhighkT->SetTitle("#pi^{-}#pi^{-} N_{ch} 50-150, 0.6 < k_{T} < 0.7 GeV/c");
96   }
97
98   gStyle->SetOptStat(0);
99   TCanvas* _can = new TCanvas("Femto QA","Femto QA");
100   _can->Divide(2,3);
101   _can->cd(1);
102   numCFlowkT->Draw();
103   _can->cd(2);
104   numCFhighkT->Draw();
105   _can->cd(3);
106   numPhiEtalowkT->Draw("colz");
107   _can->cd(4);
108   numPhiEtahighkT->Draw("colz");
109   _can->cd(5);
110   numQinvEtpclowkT->Draw("colz");
111   _can->cd(6);
112   numQinvEtpchighkT->Draw("colz");
113
114   _can->SaveAs(Form("fig_cf_FemtoQA.%s",suffix));
115
116   _file->Close();
117
118   delete _file;
119   delete _can;
120 }