]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWGPP/analysisQA/processCFv2vsPt.C
Merge branch 'feature-movesplit'
[u/mrichter/AliRoot.git] / PWGPP / analysisQA / processCFv2vsPt.C
CommitLineData
948b38b2 1/***************************************************************
2 processCFv2vsPt.C
3 Post Processing of Flow QA task in Analysis QA train
4 Clone copy of v2vsPt_0d.C from Ante
5
6 Modification Done // sjena
7 To save file and putting ani unique naming convention
8
9***************************************************************/
10
11// Name of the merged, large statistics file obtained with the merging macros:
12TString mergedFileName = "AnalysisResults.root";
13
14// Final output file name holding correct final results for large statistics sample:
15TString outputFileName = "AnalysisResults.root";
16
17void processCFv2vsPt(const char *filePath = "AnalysisResults.root",
18 TString suffix="eps",
19 const char *outfile="CFv2vsPt_output.root") {
20 // Load flow libraries:
21 gSystem->Load("libPWGflowBase");
22
23 // Let's go:
24 TString mergedFileFullPathName(filePath);
25 TFile *mergedFile = NULL;
26 if(gSystem->AccessPathName(mergedFileFullPathName.Data(),kFileExists))
27 {
28 cout<<endl;
29 cout<<" WARNING: Couldn't find a file: "<<mergedFileFullPathName.Data()<<endl;
30 cout<<endl;
31 exit(0);
32 } else
33 {
34 // Create temporarily copy of <mergedFileName> if neccessary:
35 if(!(mergedFileName == outputFileName))
36 {
37 TSystemFile *fileTemp = new TSystemFile(mergedFileFullPathName.Data(),".");
38 fileTemp->Copy("mergedAnalysisResultsTemp.root");
39 delete fileTemp;
40 }
41 // Access <mergedFileName>:
42 mergedFile = TFile::Open(mergedFileFullPathName.Data(),"UPDATE");
43 }
44
45 // Accessing the file and updating it:
46 TList* mergedFileKeys = mergedFile->GetListOfKeys();
47 for(Int_t i=0; i<mergedFileKeys->GetEntries(); i++)
48 {
49 TDirectory* directory = dynamic_cast<TDirectory*>(mergedFile->Get(mergedFileKeys->At(i)->GetName()));
50 if (!directory) continue;
51 TList* listTemp = directory->GetListOfKeys();
52 if(!listTemp) continue;
53 for (Int_t icent=0; icent<listTemp->GetEntries(); icent++)
54 {
55 TList* list = dynamic_cast<TList*>(directory->Get(listTemp->At(icent)->GetName()));
56 if (!list) continue;
57 // QC:
58 if(TString(list->GetName()).Contains("QC"))
59 {
60 AliFlowAnalysisWithQCumulants* qc = new AliFlowAnalysisWithQCumulants();
61 qc->GetOutputHistograms(list);
62 Bool_t bApplyCorrectionForNUA = kFALSE; // apply correction for non-uniform acceptance
63 qc->GetIntFlowFlags()->SetBinContent(3,(Int_t)bApplyCorrectionForNUA);
64 qc->Finish();
65 directory->Add(list,kTRUE);
66 directory->Write(directory->GetName(),TObject::kSingleKey+TObject::kWriteDelete);
67 } // if(TString(list->GetName()).Contains("QC"))
68 } // for (Int_t icent=0; icent<listTemp->GetEntries(); icent++)
69 } // for(Int_t i=0; i<mergedFileKeys->GetEntries(); i++)
70
71 // Close the final output file:
72 delete mergedFile;
73
74 // Changing the final name, if neccessary:
75 if(!(mergedFileName == outputFileName))
76 {
77 TSystemFile *outputFileFinal = new TSystemFile(mergedFileName.Data(),".");
78 outputFileFinal->Rename(outputFileName.Data());
79 delete outputFileFinal;
80 TSystemFile *mergedFileFinal = new TSystemFile("mergedAnalysisResultsTemp.root",".");
81 mergedFileFinal->Rename(mergedFileName.Data());
82 delete mergedFileFinal;
83 } // end of if(!(mergedFileName == outputFileName))
84
85 // Finally, make some nice plots...:
86 // Access "v2 vs. pT" histogram from the specified output file:
5d3299b0 87 TFile *f = TFile::Open(filePath,"READ");
948b38b2 88 TDirectoryFile *df = (TDirectoryFile*)f->FindObjectAny("outputQCanalysisTPCstandalone");
5d3299b0 89 if(!df){df = (TDirectoryFile*)f->FindObjectAny("outputQCanalysis");} // for toy MC
90 if(!df){cout<<"df is NULL pointer!!!!"<<endl; exit(0);}
948b38b2 91 TList *l = NULL;
92 TList *listTemp = df->GetListOfKeys();
93 if(listTemp && listTemp->GetEntries() == 1)
94 {
95 TString listName = listTemp->At(0)->GetName();
96 df->GetObject(listName.Data(),l);
97 }
98 TList *diffFlowList = dynamic_cast<TList*> l->FindObject("Differential Flow");
5d3299b0 99 TList *diffFlowResultsList = dynamic_cast<TList*> diffFlowList->FindObject("Results");
100 TList *diffFlowDifFlowList = dynamic_cast<TList*> diffFlowResultsList->FindObject("Differential flow (POI, p_{T})");
948b38b2 101
102 // v2{2} vs. pT:
5d3299b0 103 TH1D *v22 = dynamic_cast<TH1D*> diffFlowDifFlowList->FindObject("fDiffFlow, POI, p_{T}, v'{2}");
104 if(!v22){cout<<"v22 is NULL"<<endl; exit(0);}
948b38b2 105 v22->SetMarkerStyle(kFullCircle);
106 v22->SetMarkerColor(kBlue);
107 v22->SetLineColor(kBlue);
108
109 // v2{4} vs. pT:
5d3299b0 110 TH1D *v24 = dynamic_cast<TH1D*> diffFlowDifFlowList->FindObject("fDiffFlow, POI, p_{T}, v'{4}");
111 if(!v24){cout<<"v24 is NULL"<<endl; exit(0);}
948b38b2 112 v24->SetMarkerStyle(kFullSquare);
113 v24->SetMarkerColor(kRed);
114 v24->SetLineColor(kRed);
115
116 // Legend:
117 TLegend *legend = new TLegend(0.15,0.67,0.37,0.87);
118 legend->SetFillStyle(0); // white legend background
119 legend->SetTextSize(0.04);
120 legend->SetBorderSize(0.0);
121 //legend->SetTextFont(22);
122 legend->SetMargin(0.1);
123 //legend->AddEntry(inputValues_g,"theoretical/input values","p");
124 legend->AddEntry(v22,"v_{2}{2}","p");
125 legend->AddEntry(v24,"v_{2}{4}","p");
126
127 // Style histogram:
128 TH1D *hist = new TH1D("","",10,0.,10.);
129 hist->SetStats(kFALSE);
130 hist->SetTitle("");
131 hist->GetXaxis()->SetTitle("p_{T}");
5d3299b0 132 Double_t maxFlow_y = v22->GetMaximum(0.75); // Hardwiring that max_value is <= 0.75 for aesthetic reasons
133 if(v24->GetMaximum(0.75)>maxFlow_y){maxFlow_y = v24->GetMaximum(0.75);}
134 hist->GetYaxis()->SetRangeUser(0.0,maxFlow_y);
948b38b2 135
136 // Final drawing:
137 // v2 vs. pT:
138 TCanvas *c1 = new TCanvas("c1","v2 vs. pT");
139 hist->Draw();
140 legend->Draw("same");
141 v22->Draw("psame");
142 v24->Draw("psame");
143
5d3299b0 144 c1->SaveAs(Form("fig_cf_flow_fDiffFlow.%s",suffix.Data()));
145
146 //===========================================================================================
948b38b2 147
148 // Bias from non-uniform azimuthal acceptance (NUA):
149 TList *intFlowList = dynamic_cast<TList*> l->FindObject("Integrated Flow");
150 intFlowList = dynamic_cast<TList*> intFlowList->FindObject("Results");
151 TH1D *fIntFlowDetectorBias = dynamic_cast<TH1D*> intFlowList->FindObject("fIntFlowDetectorBias");
152 fIntFlowDetectorBias->SetStats(kFALSE);
153 fIntFlowDetectorBias->SetTitle("Quantifying bias from non-uniform acceptance");
154 fIntFlowDetectorBias->SetFillColor(kBlue-10);
155 fIntFlowDetectorBias->GetXaxis()->SetRangeUser(0.,2.);
156 fIntFlowDetectorBias->GetYaxis()->SetRangeUser(0.8044,1.0544);
157 fIntFlowDetectorBias->SetBinError(1,0.); // stat. errors need to be reimplemented
158 fIntFlowDetectorBias->SetBinError(2,0.); // stat. errors need to be reimplemented
159
160 // Final drawing:
161 TCanvas *c2 = new TCanvas("c2","NUA");
162 fIntFlowDetectorBias->Draw();
163
164 c2->SaveAs(Form("fig_cf_flow_IntFlowDetectorBias.%s",suffix.Data()));
165
5d3299b0 166 //===========================================================================================
167
168 // Differential (vs. pT) cumulants:
169 if(!diffFlowResultsList){cout<<"diffFlowResultsList is NULL"<<endl; exit(0);}
170 TList *diffFlowDifCumulantsList = dynamic_cast<TList*> diffFlowResultsList->FindObject("Differential Q-cumulants (POI, p_{T})");
171 if(!diffFlowDifCumulantsList){cout<<"diffFlowDifCumulantsList is NULL"<<endl; exit(0);}
172
173 // d2{2} vs. pT:
174 TH1D *d22 = dynamic_cast<TH1D*> diffFlowDifCumulantsList->FindObject("fDiffFlowCumulants, POI, p_{T}, QC{2'}");
175 if(!d22){cout<<"d22 is NULL"<<endl; exit(0);}
176 d22->SetMarkerStyle(kFullCircle);
177 d22->SetMarkerColor(kBlue);
178 d22->SetLineColor(kBlue);
179 d22->Scale(10); // Rescale this histogram by 10^2
180
181 // d2{4} vs. pT:
182 TH1D *d24 = dynamic_cast<TH1D*> diffFlowDifCumulantsList->FindObject("fDiffFlowCumulants, POI, p_{T}, QC{4'}");
183 if(!d24){cout<<"d24 is NULL"<<endl; exit(0);}
184 d24->SetMarkerStyle(kFullSquare);
185 d24->SetMarkerColor(kRed);
186 d24->SetLineColor(kRed);
187 d24->Scale(1000); // Rescale this histogram by 10^4
188
189 // Legend:
190 TLegend *legend_c = new TLegend(0.15,0.67,0.37,0.87);
191 legend_c->SetFillStyle(0); // white legend background
192 legend_c->SetTextSize(0.04);
193 legend_c->SetBorderSize(0.0);
194 //legend_c->SetTextFont(22);
195 legend_c->SetMargin(0.1);
196 legend_c->AddEntry(d22,"d_{2}{2} #times 10","p");
197 legend_c->AddEntry(d24,"d_{2}{4} #times 10^{3}","p");
198
199 // Style histogram:
200 TH1D *hist_c = new TH1D("","",10,0.,10.);
201 hist_c->SetStats(kFALSE);
202 hist_c->SetTitle("");
203 hist_c->GetXaxis()->SetTitle("p_{T}");
204 Double_t maxCumulant_y = d22->GetMaximum(0.75); // Hardwiring that max_value is <= 0.75 for aesthetic reasons
205 if(d24->GetMaximum(0.75)>maxCumulant_y){maxCumulant_y = d24->GetMaximum(0.75);}
206 Double_t minCumulant_y = d22->GetMinimum(-0.75); // Hardwiring that max_value is <= 0.75 for aesthetic reasons
207 if(d24->GetMinimum(-0.75)<minCumulant_y){minCumulant_y = d24->GetMinimum(-0.75);}
208 hist_c->GetYaxis()->SetRangeUser(minCumulant_y,maxCumulant_y);
209
210 // Final drawing:
211 // d2 vs. pT:
212 TCanvas *c3 = new TCanvas("c3","differential Q-cumulants vs. pT");
213 hist_c->Draw();
214 legend_c->Draw("same");
215 d22->Draw("psame");
216 d24->Draw("psame");
217
218 c3->SaveAs(Form("fig_cf_flow_fDiffFlowCumulants.%s",suffix.Data()));
219
220 //f->Close();
948b38b2 221
222 cout<<endl;
223
224 //Added by sjena
225 TFile *fout = TFile::Open(outfile,"UPDATE");
226 fout->ls();
227
228 TDirectoryFile *cdd = NULL;
229 cdd = (TDirectoryFile*)fout->Get("CF");
230 if(!cdd) {
231 Printf("Warning: CF <dir> doesn't exist, creating a new one");
232 cdd = (TDirectoryFile*)fout->mkdir("CF");
233 }
234 cdd->cd();
235 cdd->ls();
236
237
238 v22->SetName(Form("fig_cf_flow_22", v22->GetName()));
239 // v22->SetName(Form("fig_cf_flow_fDiffFlow22"));
240 v22->Write();
241
242 v24->SetName(Form("fig_cf_flow_24", v24->GetName()));
243 // v24->SetName(Form("fig_cf_flow_fDiffFlow24"));
244 v24->Write();
245
246 fIntFlowDetectorBias->SetName(Form("fig_cf_flow_IntFlowDetectorBias", fIntFlowDetectorBias->GetName()));
247 fIntFlowDetectorBias->Write();
248
249 fout->cd();
250 fout->Close();
251
252
253
254} // void void v2vsPt_0d()
255