1 const Int_t numberOfCentralityBins = 9;
2 TString centralityArray[numberOfCentralityBins] = {"0-5","5-10","10-20","20-30","30-40","40-50","50-60","60-70","70-80"};
3 Double_t gMinCentrality[numberOfCentralityBins] = {0.,5.,10.,20.,30.,40.,50.,60.,70.};
4 Double_t gMaxCentrality[numberOfCentralityBins] = {5.,10.,20.,30.,40.,50.,60.,70.,80.};
5 TString gAnalysisType[7] = {"y","eta","qlong","qout","qside","qinv","phi"};
7 const Int_t gRebin = 1;
8 void drawCorrelationFunctionPsi(const char* filename = "AnalysisResults.root",
9 Double_t psiMin = 0., Double_t psiMax = 7.5) {
10 //Macro that draws the correlation functions from the balance function
11 //analysis vs the reaction plane
12 //Author: Panos.Christakoglou@nikhef.nl
13 //Load the PWG2ebye library
14 gSystem->Load("libANALYSIS.so");
15 gSystem->Load("libANALYSISalice.so");
16 gSystem->Load("libEventMixing.so");
17 gSystem->Load("libCORRFW.so");
18 gSystem->Load("libPWGTools.so");
19 gSystem->Load("libPWGCFebye.so");
21 //Prepare the objects and return them
22 TList *list = GetListOfObjects(filename);
24 Printf("The TList object was not created");
28 draw(list,psiMin,psiMax);
31 //______________________________________________________//
32 TList *GetListOfObjects(const char* filename) {
33 //Get the TList objects (QA, bf, bf shuffled)
36 TList *listBFShuffling = 0x0;
39 TFile *f = TFile::Open(filename);
40 if((!f)||(!f->IsOpen())) {
41 Printf("The file %s is not found. Aborting...",filename);
46 TDirectoryFile *dir = dynamic_cast<TDirectoryFile *>(f->Get("PWGCFEbyE.outputBalanceFunctionAnalysis"));
48 Printf("The TDirectoryFile is not found. Aborting...",filename);
53 TString listBFName = "listBF_0-100_V0M_vZ10.0_DCAxy-1.0_DCAz-1.0_Pt0.3-5.0_Eta-0.8-0.8_Chi-1.0_nClus-1_Bit1_withCentralTrigger";
54 listBF = dynamic_cast<TList *>(dir->Get(listBFName.Data()));
58 TString histoName = "fHistPV0M";
59 AliTHn *fHistP = dynamic_cast<AliTHn *>(listBF->FindObject(histoName.Data()));
61 Printf("fHistP %s not found!!!",histoName.Data());
64 fHistP->FillParent(); fHistP->DeleteContainers();
66 histoName = "fHistNV0M";
67 AliTHn *fHistN = dynamic_cast<AliTHn *>(listBF->FindObject(histoName.Data()));
69 Printf("fHistN %s not found!!!",histoName.Data());
72 fHistN->FillParent(); fHistN->DeleteContainers();
74 histoName = "fHistPNV0M";
75 AliTHn *fHistPN = dynamic_cast<AliTHn *>(listBF->FindObject(histoName.Data()));
77 Printf("fHistPN %s not found!!!",histoName.Data());
80 fHistPN->FillParent(); fHistPN->DeleteContainers();
82 histoName = "fHistNPV0M";
83 AliTHn *fHistNP = dynamic_cast<AliTHn *>(listBF->FindObject(histoName.Data()));
85 Printf("fHistNP %s not found!!!",histoName.Data());
88 fHistNP->FillParent(); fHistNP->DeleteContainers();
90 histoName = "fHistPPV0M";
91 AliTHn *fHistPP = dynamic_cast<AliTHn *>(listBF->FindObject(histoName.Data()));
93 Printf("fHistPP %s not found!!!",histoName.Data());
96 fHistPP->FillParent(); fHistPP->DeleteContainers();
98 histoName = "fHistNNV0M";
99 AliTHn *fHistNN = dynamic_cast<AliTHn *>(listBF->FindObject(histoName.Data()));
101 Printf("fHistNN %s not found!!!",histoName.Data());
104 fHistNN->FillParent(); fHistNN->DeleteContainers();
109 //______________________________________________________//
110 void draw(TList *list, Double_t psiMin, Double_t psiMax) {
111 //Draws the correlation functions for every centrality bin
112 //(+-), (-+), (++), (--)
113 gROOT->LoadMacro("~/SetPlotStyle.C");
115 gStyle->SetPalette(1,0);
124 hP = (AliTHn*) list->FindObject("fHistPV0M");
125 hN = (AliTHn*) list->FindObject("fHistNV0M");
126 hPN = (AliTHn*) list->FindObject("fHistPNV0M");
127 hNP = (AliTHn*) list->FindObject("fHistNPV0M");
128 hPP = (AliTHn*) list->FindObject("fHistPPV0M");
129 hNN = (AliTHn*) list->FindObject("fHistNNV0M");
131 //Create the AliBalancePsi object and fill it with the AliTHn objects
132 AliBalancePsi *b = new AliBalancePsi();
139 TH2D *gHistPN[numberOfCentralityBins];
140 TH2D *gHistNP[numberOfCentralityBins];
141 TH2D *gHistPP[numberOfCentralityBins];
142 TH2D *gHistNN[numberOfCentralityBins];
144 TCanvas *cPN[numberOfCentralityBins];
145 TCanvas *cNP[numberOfCentralityBins];
146 TCanvas *cPP[numberOfCentralityBins];
147 TCanvas *cNN[numberOfCentralityBins];
148 TString histoTitle, pngName;
149 //loop over the centrality bins
150 for(Int_t iCentralityBin = 0; iCentralityBin < numberOfCentralityBins; iCentralityBin++) {
151 histoTitle = "(+-) | Centrality: ";
152 histoTitle += centralityArray[iCentralityBin];
154 histoTitle += " | "; histoTitle += psiMin;
155 histoTitle += " < #phi - #Psi_{2} < "; histoTitle += psiMax;
156 gHistPN[iCentralityBin] = b->GetCorrelationFunctionPN(gMinCentrality[iCentralityBin],gMaxCentrality[iCentralityBin],psiMin,psiMax);
157 gHistPN[iCentralityBin]->GetYaxis()->SetTitleOffset(1.5);
158 gHistPN[iCentralityBin]->SetTitle(histoTitle.Data());
159 cPN[iCentralityBin] = new TCanvas(histoTitle.Data(),"",0,0,400,400);
160 cPN[iCentralityBin]->SetFillColor(10);
161 cPN[iCentralityBin]->SetHighLightColor(10);
162 gHistPN[iCentralityBin]->DrawCopy("lego2");
163 pngName = "DeltaPhiDeltaEta.Centrality";
164 pngName += centralityArray[iCentralityBin];
165 pngName += ".Psi"; pngName += psiMin; pngName += "To"; pngName += psiMax;
166 pngName += ".PositiveNegative.png";
167 cPN[iCentralityBin]->SaveAs(pngName.Data());
169 histoTitle = "(-+) | Centrality: ";
170 histoTitle += centralityArray[iCentralityBin];
172 histoTitle += " | "; histoTitle += psiMin;
173 histoTitle += " < #phi - #Psi_{2} < "; histoTitle += psiMax;
174 gHistNP[iCentralityBin] = b->GetCorrelationFunctionNP(gMinCentrality[iCentralityBin],gMaxCentrality[iCentralityBin],psiMin,psiMax);
175 gHistNP[iCentralityBin]->GetYaxis()->SetTitleOffset(1.5);
176 gHistNP[iCentralityBin]->SetTitle(histoTitle.Data());
177 cNP[iCentralityBin] = new TCanvas(histoTitle.Data(),"",400,0,400,400);
178 cNP[iCentralityBin]->SetFillColor(10);
179 cNP[iCentralityBin]->SetHighLightColor(10);
180 gHistNP[iCentralityBin]->DrawCopy("lego2");
181 pngName = "DeltaPhiDeltaEta.Centrality";
182 pngName += centralityArray[iCentralityBin];
183 pngName += ".Psi"; pngName += psiMin; pngName += "To"; pngName += psiMax;
184 pngName += ".NegativePositive.png";
185 cNP[iCentralityBin]->SaveAs(pngName.Data());
187 histoTitle = "(++) | Centrality: ";
188 histoTitle += centralityArray[iCentralityBin];
190 histoTitle += " | "; histoTitle += psiMin;
191 histoTitle += " < #phi - #Psi_{2} < "; histoTitle += psiMax;
192 gHistPP[iCentralityBin] = b->GetCorrelationFunctionPP(gMinCentrality[iCentralityBin],gMaxCentrality[iCentralityBin],psiMin,psiMax);
193 gHistPP[iCentralityBin]->GetYaxis()->SetTitleOffset(1.5);
194 gHistPP[iCentralityBin]->SetTitle(histoTitle.Data());
195 cPP[iCentralityBin] = new TCanvas(histoTitle.Data(),"",0,400,400,400);
196 cPP[iCentralityBin]->SetFillColor(10);
197 cPP[iCentralityBin]->SetHighLightColor(10);
198 gHistPP[iCentralityBin]->DrawCopy("lego2");
199 pngName = "DeltaPhiDeltaEta.Centrality";
200 pngName += centralityArray[iCentralityBin];
201 pngName += ".Psi"; pngName += psiMin; pngName += "To"; pngName += psiMax;
202 pngName += ".PositivePositive.png";
203 cPP[iCentralityBin]->SaveAs(pngName.Data());
205 histoTitle = "(--) | Centrality: ";
206 histoTitle += centralityArray[iCentralityBin];
208 histoTitle += " | "; histoTitle += psiMin;
209 histoTitle += " < #phi - #Psi_{2} < "; histoTitle += psiMax;
210 gHistNN[iCentralityBin] = b->GetCorrelationFunctionNN(gMinCentrality[iCentralityBin],gMaxCentrality[iCentralityBin],psiMin,psiMax);
211 gHistNN[iCentralityBin]->GetYaxis()->SetTitleOffset(1.5);
212 gHistNN[iCentralityBin]->SetTitle(histoTitle.Data());
213 cNN[iCentralityBin] = new TCanvas(histoTitle.Data(),"",400,400,400,400);
214 cNN[iCentralityBin]->SetFillColor(10);
215 cNN[iCentralityBin]->SetHighLightColor(10);
216 gHistNN[iCentralityBin]->DrawCopy("lego2");
217 pngName = "DeltaPhiDeltaEta.Centrality";
218 pngName += centralityArray[iCentralityBin];
219 pngName += ".Psi"; pngName += psiMin; pngName += "To"; pngName += psiMax;
220 pngName += ".NegativeNegative.png";
221 cNN[iCentralityBin]->SaveAs(pngName.Data());
222 }//end of loop over centralities