]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWG2/FLOW/macros/runFlowAnalysisOnTheFly.C
option to calculate cumulants with various multiplicity weights
[u/mrichter/AliRoot.git] / PWG2 / FLOW / macros / runFlowAnalysisOnTheFly.C
CommitLineData
93ff27bd 1#include "TStopwatch.h"
2#include "TObjArray"
3#include "Riostream.h"
4#include "TFile.h"
5
6//--------------------------------------------------------------------------------------
7// RUN SETTINGS
8// flow analysis method can be: (set to kTRUE or kFALSE)
c741f5d0 9Bool_t SP = kTRUE;
5b40431d 10Bool_t LYZ1SUM = kTRUE;
11Bool_t LYZ1PROD = kTRUE;
c741f5d0 12Bool_t LYZ2SUM = kFALSE;
13Bool_t LYZ2PROD = kFALSE;
5b40431d 14Bool_t LYZEP = kFALSE;
c741f5d0 15Bool_t GFC = kTRUE;
16Bool_t QC = kTRUE;
e23cbd39 17Bool_t FQD = kTRUE;
18Bool_t MCEP = kTRUE;
93ff27bd 19//--------------------------------------------------------------------------------------
20
7fb984f6 21// Weights
22// use weights for Q-vector:
23Bool_t usePhiWeights = kFALSE; // phi weights (correction for non-uniform azimuthal acceptance)
24Bool_t usePtWeights = kFALSE; // pt weights
25Bool_t useEtaWeights = kFALSE; // eta weights
26
ce4a88f5 27// Run same flow analysis method but with different settings/aims
28// You will have to label each setting/aim with your own label (see examples bellow):
29Bool_t GFC_Additional_Analysis = kFALSE;
30Bool_t QC_Additional_Analysis = kFALSE;
31Bool_t FQD_Additional_Analysis = kFALSE;
32
7fb984f6 33// Parameters for the simulation of events 'on the fly':
e23cbd39 34Bool_t bSameSeed = kFALSE; // use always the same seed for random generators.
c741f5d0 35 // usage of same seed (kTRUE) is relevant in two cases:
7fb984f6 36 // 1.) If you want to use LYZ method to calcualte differential flow;
37 // 2.) If you want to use phi weights for GFC, QC and FQD
678f8252 38
7fb984f6 39Bool_t bConstantHarmonics = kTRUE; // harmonics V1, V2, V4... are constant (kTRUE) or functions of pt and eta (kFALSE)
4d603348 40
aaa62cf0 41Int_t iLoops = 1; // number of times to use each track (to simulate nonflow)
4d603348 42
e23cbd39 43Bool_t bMultDistrOfRPsIsGauss = kTRUE; // 1.) if kTRUE = multiplicitiy of RPs is sampled e-b-e from Gaussian distribution with
678f8252 44 // mean = iMultiplicityOfRP and spread = dMultiplicitySpreadOfRP
45 // 2.) if kFALSE = multiplicitiy of RPs is sampled e-b-e uniformly from
46 // interval [iMinMultOfRP,iMaxMultOfRP]
47 // 3.) for a fixed multiplicity use Gaussian with zero spread or use uniform with iMinMult=iMaxMult
55381065 48
e23cbd39 49Bool_t bV2DistrOfRPsIsGauss = kTRUE; // 1.) if kTRUE = elliptic flow of RPs is sampled e-b-e from Gaussian distribution with
55381065 50 // mean = dV2RP and spread = dV2SpreadRP
51 // 2.) if kFALSE = elliptic flow of RPs is sampled e-b-e uniformly from
52 // interval [dMinV2RP,dMaxV2RP]
53 // 3.) for a fixed elliptic flow use Gaussian with zero spread or use uniform with dMinV2RP=dMaxV2RP
e23cbd39 54
55Bool_t uniformAcceptance = kTRUE; // 1.) if kTRUE = detectors has uniform azimuthal acceptance
56 // 2.) if kFALSE = you will simulate detector with non-uniform acceptance in one or two sectors.
57 // For each of two sectors you specify phi_min, phi_max and probability p. Then all particles
58 // going in direction phi_min < phi < phi_max will be taken with probability p. If p = 0, that
59 // sector is blocked. Set bellow phimin1, phimax1, p1 for the first sector and phimin2, phimax2, p2
60 // for the second sector. If you set phimin2 = phimax2 = p2 = 0, only first non-uniform sector is
61 // simulated.
62
678f8252 63Int_t iMultiplicityOfRP = 500; // mean multiplicity of RPs (if sampled from Gaussian)
64Double_t dMultiplicitySpreadOfRP = 0; // multiplicity spread of RPs (if sampled from Gaussian)
65Int_t iMinMultOfRP = 400; // minimal multiplicity of RPs (if sampled uniformly)
66Int_t iMaxMultOfRP = 600; // maximal multiplicity of RPs (if sampled uniformly)
67
cf90787f 68Double_t dTemperatureOfRP = 0.44; // 'temperature' of RPs in GeV/c (increase this parameter to get more high pt RPs)
93ff27bd 69
4d603348 70//......................................................................................
71// if you use (pt,eta) dependent harmonics (bConstantHarmonics = kFALSE):
72Double_t dPtCutOff = 2.0; // V2(pt) is linear up to pt = 2 GeV and for pt > 2 GeV it is constant: V2(pt) = dVRPMax
cf90787f 73Double_t dV2RPMax = 0.20; // maximum value of V2(pt) for pt >= 2GeV
4d603348 74//......................................................................................
75
76//......................................................................................
55381065 77// if you use constant harmonics (bConstantHarmonics = kTRUE) (i.e. no pt dependence):
78Double_t dV2RP = 0.05; // elliptic flow of RPs (if sampled from Gaussian)
79Double_t dV2SpreadRP = 0.0; // elliptic flow spread of RPs (if sampled from Gaussian)
80Double_t dMinV2RP = 0.04; // minimal elliptic flow of RPs (if sampled uniformly)
81Double_t dMaxV2RP = 0.06; // maximal elliptic flow of RPs (if sampled uniformly)
93ff27bd 82
e8c7e66c 83Double_t dV1RP = 0.0; // directed flow of RPs
84Double_t dV1SpreadRP = 0.0; // directed flow spread of RPs
93ff27bd 85
4d603348 86Double_t dV4RP = 0.0; // harmonic V4 of RPs (to be improved: name needed)
87Double_t dV4SpreadRP = 0.0; // harmonic V4's spread of RPs (to be improved: name needed)
88//......................................................................................
89
e23cbd39 90//......................................................................................
91// settings for non-uniform acceptance:
92// Remark: set the angles in degrees from interval [0,360] and probability from interval [0,1]
93
94// 1st non-uniform sector:
95Double_t phimin1 = 60; // first non-uniform sector starts at this azimuth
96Double_t phimax1 = 120; // first non-uniform sector ends at this azimuth
97Double_t p1 = 0.33; // e.g. if p1 = 0 all particles emitted in phimin1 < phi < phimax1 are blocked
98 // e.g. if p1 = 0.5 half of the particles emitted in phimin1 < phi < phimax1 are blocked
99
100// 2nd non-uniform sector (Remark: if you do NOT want to simulate this sector, set phimin2 = phimax2 = p2 = 0):
101Double_t phimin2 = 0.0; // second non-uniform sector starts at this azimuth (make sure phimin2 > phimax1 !!!!)
102Double_t phimax2 = 0.0; // second non-uniform sector ends at this azimuth
103Double_t p2 = 0.0;
104//......................................................................................
105
4d603348 106enum anaModes {mLocal,mLocalSource,mLocalPAR};
93ff27bd 107// mLocal: Analyze data on your computer using aliroot
108// mLocalPAR: Analyze data on your computer using root + PAR files
109// mLocalSource: Analyze data on your computer using root + source files
110
129292f8 111int runFlowAnalysisOnTheFly( Int_t nEvts=440, Int_t mode=mLocal)
93ff27bd 112{
113 TStopwatch timer;
114 timer.Start();
115
c741f5d0 116 if (LYZ1SUM && LYZ2SUM) {cout<<"WARNING: you cannot run LYZ1 and LYZ2 at the same time! LYZ2 needs the output from LYZ1. "<<endl; exit(); }
117 if (LYZ1PROD && LYZ2PROD) {cout<<"WARNING: you cannot run LYZ1 and LYZ2 at the same time! LYZ2 needs the output from LYZ1. "<<endl; exit(); }
118 if (LYZ2SUM && LYZEP) {cout<<"WARNING: you cannot run LYZ2 and LYZEP at the same time! LYZEP needs the output from LYZ2."<<endl; exit(); }
119 if (LYZ1SUM && LYZEP) {cout<<"WARNING: you cannot run LYZ1 and LYZEP at the same time! LYZEP needs the output from LYZ2."<<endl; exit(); }
e23cbd39 120
121 if(!uniformAcceptance && phimin1 > phimax1)
122 {
123 cout<<"WARNING: you must have phimin1 < phimax1 !!!!"<<endl;
124 break;
125 }
126
127 if (!uniformAcceptance && !((phimin2 == 0.) && (phimax2 == 0.) && (p2 == 0.)) && (phimin2 < phimax1 || phimin2 > phimax2))
128 {
129 cout<<"WARNING: you must have phimin2 > phimax1 and phimin2 < phimax2 !!!!"<<endl;
130 break;
131 }
132
133 if((phimin1 < 0 || phimin1 > 360) || (phimax1 < 0 || phimax1 > 360) ||
134 (phimin2 < 0 || phimin2 > 360) || (phimax2 < 0 || phimax2 > 360) )
135 {
136 cout<<"WARNING: you must take azimuthal angles from interval [0,360] !!!!"<<endl;
137 break;
138 }
139
140 if((p1 < 0 || p1 > 1) || (p2 < 0 || p2 > 1))
141 {
142 cout<<"WARNING: you must take p1 and p2 from interval [0,1] !!!!"<<endl;
143 break;
144 }
93ff27bd 145
146 cout<<endl;
147 cout<<endl;
148 cout<<" ---- ARE YOU READY TO FLY ? ---- "<<endl;
149 cout<<endl;
150
151 cout<<endl;
152 cout<<" ---- BEGIN FLOW ANALYSIS 'ON THE FLY' ---- "<<endl;
153 cout<<endl;
154 cout<<endl;
155
93ff27bd 156 LoadLibraries(mode);
157
4d603348 158 // Initialize the seed for random generator
159 UInt_t sseed = 0;
160
161 if(bSameSeed)
162 {
163 sseed = 44; // the default constant value for seed for random generators
164 }
5b8dccde 165
4d603348 166 if(!bSameSeed)
167 {
168 TTimeStamp dt;
169 sseed = dt.GetNanoSec()/1000;
170 }
171
e23cbd39 172 cout<<endl;
173 cout<<"Seed for the random generators is "<<sseed<<endl;
174 cout<<endl;
175
7fb984f6 176 //---------------------------------------------------------------------------------------
177 // If the weights are used:
178 TFile *fileWithWeights = NULL;
179 TList *listWithWeights = NULL;
180
181 if(usePhiWeights||usePtWeights||useEtaWeights) {
182 fileWithWeights = TFile::Open("weights.root","READ");
183 if(fileWithWeights) {
184 listWithWeights = (TList*)fileWithWeights->Get("weights");
185 }
186 else
187 {cout << " WARNING: the file <weights.root> with weights from the previous run was not found."<<endl;
188 break;
189 }
190 }
191
93ff27bd 192 //---------------------------------------------------------------------------------------
193 // Initialize the flowevent maker
5b8dccde 194 AliFlowEventSimpleMakerOnTheFly* eventMakerOnTheFly = new AliFlowEventSimpleMakerOnTheFly(sseed);
4d603348 195 eventMakerOnTheFly->Init();
93ff27bd 196
93ff27bd 197 //---------------------------------------------------------------------------------------
ce4a88f5 198 // Initialize all the flow methods for default analysis:
199 AliFlowAnalysisWithQCumulants *qc = NULL;
200 AliFlowAnalysisWithCumulants *gfc = NULL;
201 AliFlowAnalysisWithFittingQDistribution *fqd = NULL;
202 AliFlowAnalysisWithLeeYangZeros *lyz1sum = NULL;
203 AliFlowAnalysisWithLeeYangZeros *lyz1prod = NULL;
204 AliFlowAnalysisWithLeeYangZeros *lyz2sum = NULL;
205 AliFlowAnalysisWithLeeYangZeros *lyz2prod = NULL;
206 AliFlowAnalysisWithLYZEventPlane *lyzep = NULL;
207 AliFlowAnalysisWithScalarProduct *sp = NULL;
208 AliFlowAnalysisWithMCEventPlane *mcep = NULL;
93ff27bd 209
210 // MCEP = monte carlo event plane
211 if (MCEP) {
212 AliFlowAnalysisWithMCEventPlane *mcep = new AliFlowAnalysisWithMCEventPlane();
afa4af05 213 // mcep->SetHarmonic(2); // default is v2
93ff27bd 214 mcep->Init();
215 }
216
ce4a88f5 217 // QC = Q-cumulants
93ff27bd 218 if(QC) {
219 AliFlowAnalysisWithQCumulants* qc = new AliFlowAnalysisWithQCumulants();
7fb984f6 220 if(listWithWeights) qc->SetWeightsList(listWithWeights);
221 if(usePhiWeights) qc->SetUsePhiWeights(usePhiWeights);
222 if(usePtWeights) qc->SetUsePtWeights(usePtWeights);
223 if(useEtaWeights) qc->SetUseEtaWeights(useEtaWeights);
afa4af05 224 // qc->SetHarmonic(2); // default is v2
225 // qc->SetApplyCorrectionForNUA(kTRUE); // default
a5b7efd0 226 // qc->SetCalculate2DFlow(kFALSE); // default
227 // qc->SetMultiplicityWeight("combinations"); // default
228 // qc->SetMultiplicityWeight("unit");
229 // qc->SetMultiplicityWeight("multiplicity");
ce4a88f5 230 qc->Init();
93ff27bd 231 }
232
233 // GFC = Generating Function Cumulants
234 if(GFC) {
235 AliFlowAnalysisWithCumulants* gfc = new AliFlowAnalysisWithCumulants();
7fb984f6 236 if(listWithWeights) gfc->SetWeightsList(listWithWeights);
237 if(usePhiWeights) gfc->SetUsePhiWeights(usePhiWeights);
238 if(usePtWeights) gfc->SetUsePtWeights(usePtWeights);
239 if(useEtaWeights) gfc->SetUseEtaWeights(useEtaWeights);
ce4a88f5 240 gfc->Init();
93ff27bd 241 }
242
243 // FQD = Fitting q-distribution
244 if(FQD) {
ce4a88f5 245 AliFlowAnalysisWithFittingQDistribution* fqd = new AliFlowAnalysisWithFittingQDistribution();
7fb984f6 246 if(listWithWeights) fqd->SetWeightsList(listWithWeights);
247 if(usePhiWeights) fqd->SetUsePhiWeights(usePhiWeights);
ce4a88f5 248 fqd->Init();
93ff27bd 249 }
250
251 // SP = Scalar Product
252 if(SP) {
253 AliFlowAnalysisWithScalarProduct* sp = new AliFlowAnalysisWithScalarProduct();
254 sp->Init();
255 }
256
257 // LYZ1 = Lee-Yang Zeroes first run
c741f5d0 258 if(LYZ1SUM) {
259 AliFlowAnalysisWithLeeYangZeros* lyz1sum = new AliFlowAnalysisWithLeeYangZeros();
260 lyz1sum->SetFirstRun(kTRUE);
261 lyz1sum->SetUseSum(kTRUE);
262 lyz1sum->Init();
263 }
264 if(LYZ1PROD) {
265 AliFlowAnalysisWithLeeYangZeros* lyz1prod = new AliFlowAnalysisWithLeeYangZeros();
266 lyz1prod->SetFirstRun(kTRUE);
267 lyz1prod->SetUseSum(kFALSE);
268 lyz1prod->Init();
93ff27bd 269 }
93ff27bd 270 // LYZ2 = Lee-Yang Zeroes second run
c741f5d0 271 if(LYZ2SUM) {
272 AliFlowAnalysisWithLeeYangZeros* lyz2sum = new AliFlowAnalysisWithLeeYangZeros();
93ff27bd 273 // read the input file from the first run
c741f5d0 274 TString inputFileNameLYZ2SUM = "outputLYZ1SUManalysis.root" ;
275 TFile* inputFileLYZ2SUM = new TFile(inputFileNameLYZ2SUM.Data(),"READ");
276 if(!inputFileLYZ2SUM || inputFileLYZ2SUM->IsZombie()) {
5b40431d 277 cerr << " ERROR: To run LYZ2SUM you need the output file from LYZ1SUM. This file is not there! Please run LYZ1SUM first." << endl ;
93ff27bd 278 break;
279 }
280 else {
c741f5d0 281 TList* inputListLYZ2SUM = (TList*)inputFileLYZ2SUM->Get("cobjLYZ1SUM");
282 if (!inputListLYZ2SUM) {cout<<"Input list LYZ2SUM is NULL pointer!"<<endl; break;}
93ff27bd 283 else {
c741f5d0 284 cout<<"LYZ2SUM input file/list read..."<<endl;
285 lyz2sum->SetFirstRunList(inputListLYZ2SUM);
286 lyz2sum->SetFirstRun(kFALSE);
287 lyz2sum->SetUseSum(kTRUE);
288 lyz2sum->Init();
93ff27bd 289 }
290 }
291 }
c741f5d0 292 if(LYZ2PROD) {
293 AliFlowAnalysisWithLeeYangZeros* lyz2prod = new AliFlowAnalysisWithLeeYangZeros();
294 // read the input file from the first run
295 TString inputFileNameLYZ2PROD = "outputLYZ1PRODanalysis.root" ;
296 TFile* inputFileLYZ2PROD = new TFile(inputFileNameLYZ2PROD.Data(),"READ");
297 if(!inputFileLYZ2PROD || inputFileLYZ2PROD->IsZombie()) {
5b40431d 298 cerr << " ERROR: To run LYZ2PROD you need the output file from LYZ1PROD. This file is not there! Please run LYZ1PROD first." << endl ;
c741f5d0 299 break;
300 }
301 else {
302 TList* inputListLYZ2PROD = (TList*)inputFileLYZ2PROD->Get("cobjLYZ1PROD");
303 if (!inputListLYZ2PROD) {cout<<"Input list LYZ2PROD is NULL pointer!"<<endl; break;}
304 else {
305 cout<<"LYZ2PROD input file/list read..."<<endl;
306 lyz2prod->SetFirstRunList(inputListLYZ2PROD);
307 lyz2prod->SetFirstRun(kFALSE);
308 lyz2prod->SetUseSum(kFALSE);
309 lyz2prod->Init();
310 }
311 }
312 }
313
93ff27bd 314 // LYZEP = Lee-Yang Zeroes event plane
315 if(LYZEP) {
316 AliFlowLYZEventPlane* ep = new AliFlowLYZEventPlane() ;
317 AliFlowAnalysisWithLYZEventPlane* lyzep = new AliFlowAnalysisWithLYZEventPlane();
318 // read the input file from the second lyz run
c741f5d0 319 TString inputFileNameLYZEP = "outputLYZ2SUManalysis.root" ;
93ff27bd 320 TFile* inputFileLYZEP = new TFile(inputFileNameLYZEP.Data(),"READ");
321 if(!inputFileLYZEP || inputFileLYZEP->IsZombie()) {
5b40431d 322 cerr << " ERROR: To run LYZEP you need the output file from LYZ2SUM. This file is not there! Please run LYZ2SUM first." << endl ;
93ff27bd 323 break;
324 }
325 else {
c741f5d0 326 TList* inputListLYZEP = (TList*)inputFileLYZEP->Get("cobjLYZ2SUM");
327 if (!inputListLYZEP) {cout<<"Input list LYZEP is NULL pointer!"<<endl; break;}
93ff27bd 328 else {
329 cout<<"LYZEP input file/list read..."<<endl;
330 ep ->SetSecondRunList(inputListLYZEP);
331 lyzep->SetSecondRunList(inputListLYZEP);
332 ep ->Init();
333 lyzep->Init();
334 }
335 }
336 }
337 //---------------------------------------------------------------------------------------
ce4a88f5 338
339 //---------------------------------------------------------------------------------------
340 // Initialize all the flow methods for additional analysis with different settings/aims
341 // Label each setting/aim with different label !!!!
342
343 // GFC:
344 TString gfcDefaultName = "outputGFCanalysis";
345 // 1.) GFC analysis for elliptic flow with r0 = 1.5:
346 AliFlowAnalysisWithCumulants *gfc_1;
347 TString gfcAnalysisLabels_1 = "_r0_1.5"; // all histograms and output file name will have this label
348 TString gfcOutputFileName_1;
349 gfcOutputFileName_1 = gfcDefaultName.Data();
350 gfcOutputFileName_1 += gfcAnalysisLabels_1.Data();
351 gfcOutputFileName_1 += ".root";
352 if(GFC_Additional_Analysis)
353 {
354 gfc_1 = new AliFlowAnalysisWithCumulants();
355 //gfc_1->SetAnalysisLabel(gfcAnalysisLabels_1.Data());
356 if(listWithWeights) gfc_1->SetWeightsList(listWithWeights);
357 if(usePhiWeights) gfc_1->SetUsePhiWeights(usePhiWeights);
358 if(usePtWeights) gfc_1->SetUsePtWeights(usePtWeights);
359 if(useEtaWeights) gfc_1->SetUseEtaWeights(useEtaWeights);
360 gfc_1->Init();
361 }
362
363 // QC:
364 TString qcDefaultName = "outputQCanalysis";
365 // 1.) QC analysis for directed flow:
366 AliFlowAnalysisWithQCumulants *qc_1;
367 TString qcAnalysisLabels_1 = "_v1"; // all histograms and output file name will have this label
368 TString qcOutputFileName_1;
369 qcOutputFileName_1 = qcDefaultName.Data();
370 qcOutputFileName_1 += qcAnalysisLabels_1.Data();
371 qcOutputFileName_1 += ".root";
372 if(QC_Additional_Analysis)
373 {
374 qc_1 = new AliFlowAnalysisWithQCumulants();
375 //qc_1->SetAnalysisLabel(qcAnalysisLabels_1->Data());
376 if(listWithWeights) qc_1->SetWeightsList(listWithWeights);
377 if(usePhiWeights) qc_1->SetUsePhiWeights(usePhiWeights);
378 qc_1->Init();
379 }
380
381 // FQD:
382 TString fqdDefaultName = "outputFQDanalysis";
383 // 1.) FQD fitting with fixed sigma:
384 AliFlowAnalysisWithFittingQDistribution *fqd_1;
385 TString fqdAnalysisLabels_1 = "_fixedSigma"; // all histograms and output file name will have this label
386 TString fqdOutputFileName_1;
387 fqdOutputFileName_1 = fqdDefaultName.Data();
388 fqdOutputFileName_1 += fqdAnalysisLabels_1.Data();
389 fqdOutputFileName_1 += ".root";
390 if(FQD_Additional_Analysis)
391 {
392 fqd_1 = new AliFlowAnalysisWithFittingQDistribution();
393 //fqd_1->SetAnalysisLabel(fqdAnalysisLabels_1->Data());
394 if(listWithWeights) fqd_1->SetWeightsList(listWithWeights);
395 if(usePhiWeights) fqd_1->SetUsePhiWeights(usePhiWeights);
396 fqd_1->Init();
397 }
398 //---------------------------------------------------------------------------------------
93ff27bd 399
e8c7e66c 400 // set the global event parameters:
aaa62cf0 401 eventMakerOnTheFly->SetNoOfLoops(iLoops);
678f8252 402
403 if(bMultDistrOfRPsIsGauss)
404 {
405 eventMakerOnTheFly->SetMultDistrOfRPsIsGauss(bMultDistrOfRPsIsGauss);
406 eventMakerOnTheFly->SetMultiplicityOfRP(iMultiplicityOfRP);
407 eventMakerOnTheFly->SetMultiplicitySpreadOfRP(dMultiplicitySpreadOfRP);
408 } else
409 {
410 eventMakerOnTheFly->SetMultDistrOfRPsIsGauss(bMultDistrOfRPsIsGauss);
411 eventMakerOnTheFly->SetMinMultOfRP(iMinMultOfRP);
412 eventMakerOnTheFly->SetMaxMultOfRP(iMaxMultOfRP);
413 }
414
cf90787f 415 eventMakerOnTheFly->SetTemperatureOfRP(dTemperatureOfRP);
4d603348 416
e8c7e66c 417 eventMakerOnTheFly->SetV1RP(dV1RP);
418 eventMakerOnTheFly->SetV1SpreadRP(dV1SpreadRP);
4d603348 419 eventMakerOnTheFly->SetV4RP(dV4RP);
420 eventMakerOnTheFly->SetV4SpreadRP(dV4SpreadRP);
421
422 // constant harmonic V2:
423 if(bConstantHarmonics)
424 {
425 eventMakerOnTheFly->SetUseConstantHarmonics(bConstantHarmonics);
55381065 426 if(bV2DistrOfRPsIsGauss)
427 {
428 eventMakerOnTheFly->SetV2DistrOfRPsIsGauss(bV2DistrOfRPsIsGauss);
429 eventMakerOnTheFly->SetV2RP(dV2RP);
430 eventMakerOnTheFly->SetV2SpreadRP(dV2SpreadRP);
431 } else
432 {
433 eventMakerOnTheFly->SetV2DistrOfRPsIsGauss(bV2DistrOfRPsIsGauss);
434 eventMakerOnTheFly->SetMinV2RP(dMinV2RP);
435 eventMakerOnTheFly->SetMaxV2RP(dMaxV2RP);
436 }
4d603348 437 }
e23cbd39 438
4d603348 439 // (pt,eta) dependent harmonic V2:
440 if(!bConstantHarmonics)
441 {
442 eventMakerOnTheFly->SetUseConstantHarmonics(bConstantHarmonics);
443 eventMakerOnTheFly->SetV2RPMax(dV2RPMax);
444 eventMakerOnTheFly->SetPtCutOff(dPtCutOff);
445 }
e23cbd39 446
447 // non-uniform acceptance:
448 if(!uniformAcceptance)
449 {
450 eventMakerOnTheFly->SetFirstSectorPhiMin(phimin1);
451 eventMakerOnTheFly->SetFirstSectorPhiMax(phimax1);
452 eventMakerOnTheFly->SetFirstSectorProbability(p1);
453 eventMakerOnTheFly->SetSecondSectorPhiMin(phimin2);
454 eventMakerOnTheFly->SetSecondSectorPhiMax(phimax2);
455 eventMakerOnTheFly->SetSecondSectorProbability(p2);
456 }
4d603348 457
93ff27bd 458 //---------------------------------------------------------------------------------------
459 // create and analyze events 'on the fly':
460
93ff27bd 461 for(Int_t i=0;i<nEvts;i++) {
462 // creating the event with above settings:
93ff27bd 463 AliFlowEventSimple *event = eventMakerOnTheFly->CreateEventOnTheFly();
93ff27bd 464
465 // analyzing the created event 'on the fly':
93ff27bd 466 // do flow analysis for various methods:
c741f5d0 467 if(MCEP) mcep->Make(event);
468 if(QC) qc->Make(event);
469 if(GFC) gfc->Make(event);
470 if(FQD) fqd->Make(event);
471 if(LYZ1SUM) lyz1sum->Make(event);
472 if(LYZ1PROD)lyz1prod->Make(event);
473 if(LYZ2SUM) lyz2sum->Make(event);
474 if(LYZ2PROD)lyz2prod->Make(event);
475 if(LYZEP) lyzep->Make(event,ep);
476 if(SP) sp->Make(event);
93ff27bd 477
ce4a88f5 478 if(GFC_Additional_Analysis)
479 {
480 // r0 = 1.5:
481 gfc_1->Make(event);
482 }
483 if(QC_Additional_Analysis)
484 {
485 // v1:
486 qc_1->Make(event);
487 }
488 if(FQD_Additional_Analysis)
489 {
490 // fixed sigma:
491 fqd_1->Make(event);
492 }
493
93ff27bd 494 delete event;
495 } // end of for(Int_t i=0;i<nEvts;i++)
496 //---------------------------------------------------------------------------------------
497
498
499
500 //---------------------------------------------------------------------------------------
ce4a88f5 501 // calculating and storing the final results of default flow analysis:
c741f5d0 502 if(MCEP) {mcep->Finish(); mcep->WriteHistograms("outputMCEPanalysis.root");}
503 if(SP) {sp->Finish(); sp->WriteHistograms("outputSPanalysis.root");}
504 if(QC) {qc->Finish(); qc->WriteHistograms("outputQCanalysis.root");}
505 if(GFC) {gfc->Finish(); gfc->WriteHistograms("outputGFCanalysis.root");}
506 if(FQD) {fqd->Finish(); fqd->WriteHistograms("outputFQDanalysis.root");}
507 if(LYZ1SUM) {lyz1sum->Finish(); lyz1sum->WriteHistograms("outputLYZ1SUManalysis.root");}
508 if(LYZ1PROD){lyz1prod->Finish();lyz1prod->WriteHistograms("outputLYZ1PRODanalysis.root");}
509 if(LYZ2SUM) {lyz2sum->Finish(); lyz2sum->WriteHistograms("outputLYZ2SUManalysis.root");}
510 if(LYZ2PROD){lyz2prod->Finish();lyz2prod->WriteHistograms("outputLYZ2PRODanalysis.root");}
511 if(LYZEP) {lyzep->Finish(); lyzep->WriteHistograms("outputLYZEPanalysis.root");}
93ff27bd 512 //---------------------------------------------------------------------------------------
513
ce4a88f5 514 //---------------------------------------------------------------------------------------
515 // calculating and storing the final results of flow analysis with different settings/aims:
516 if(GFC_Additional_Analysis)
517 {
518 // r0 = 1.5:
519 gfc_1->Finish();
520 gfc_1->WriteHistograms(gfcOutputFileName_1.Data());
521 }
522 if(QC_Additional_Analysis)
523 {
524 // v1:
525 qc_1->Finish();
526 qc_1->WriteHistograms(qcOutputFileName_1.Data());
527 }
528 if(FQD_Additional_Analysis)
529 {
530 // fixed sigma:
531 fqd_1->Finish();
532 fqd_1->WriteHistograms(fqdOutputFileName_1.Data());
533 }
534 //---------------------------------------------------------------------------------------
93ff27bd 535
536 cout<<endl;
537 cout<<endl;
538 cout<<" ---- LANDED SUCCESSFULLY ---- "<<endl;
539 cout<<endl;
540
541 timer.Stop();
542 cout << endl;
543 timer.Print();
544}
545
546void SetupPar(char* pararchivename)
547{
548 //Load par files, create analysis libraries
549 //For testing, if par file already decompressed and modified
550 //classes then do not decompress.
551
552 TString cdir(Form("%s", gSystem->WorkingDirectory() )) ;
553 TString parpar(Form("%s.par", pararchivename)) ;
554 if ( gSystem->AccessPathName(parpar.Data()) ) {
555 gSystem->ChangeDirectory(gSystem->Getenv("ALICE_ROOT")) ;
556 TString processline(Form(".! make %s", parpar.Data())) ;
557 gROOT->ProcessLine(processline.Data()) ;
558 gSystem->ChangeDirectory(cdir) ;
559 processline = Form(".! mv /tmp/%s .", parpar.Data()) ;
560 gROOT->ProcessLine(processline.Data()) ;
561 }
562 if ( gSystem->AccessPathName(pararchivename) ) {
563 TString processline = Form(".! tar xvzf %s",parpar.Data()) ;
564 gROOT->ProcessLine(processline.Data());
565 }
566
567 TString ocwd = gSystem->WorkingDirectory();
568 gSystem->ChangeDirectory(pararchivename);
569
570 // check for BUILD.sh and execute
571 if (!gSystem->AccessPathName("PROOF-INF/BUILD.sh")) {
572 printf("*******************************\n");
573 printf("*** Building PAR archive ***\n");
574 cout<<pararchivename<<endl;
575 printf("*******************************\n");
576
577 if (gSystem->Exec("PROOF-INF/BUILD.sh")) {
578 Error("runProcess","Cannot Build the PAR Archive! - Abort!");
579 return -1;
580 }
581 }
582 // check for SETUP.C and execute
583 if (!gSystem->AccessPathName("PROOF-INF/SETUP.C")) {
584 printf("*******************************\n");
585 printf("*** Setup PAR archive ***\n");
586 cout<<pararchivename<<endl;
587 printf("*******************************\n");
588 gROOT->Macro("PROOF-INF/SETUP.C");
589 }
590
591 gSystem->ChangeDirectory(ocwd.Data());
592 printf("Current dir: %s\n", ocwd.Data());
593}
594
595void LoadLibraries(const anaModes mode) {
596
597 //--------------------------------------
598 // Load the needed libraries most of them already loaded by aliroot
599 //--------------------------------------
5d040cf3 600 gSystem->Load("libTree");
601 gSystem->Load("libGeom");
602 gSystem->Load("libVMC");
603 gSystem->Load("libXMLIO");
604 gSystem->Load("libPhysics");
93ff27bd 605
606 //----------------------------------------------------------
607 // >>>>>>>>>>> Local mode <<<<<<<<<<<<<<
608 //----------------------------------------------------------
609 if (mode==mLocal) {
610 //--------------------------------------------------------
611 // If you want to use already compiled libraries
612 // in the aliroot distribution
613 //--------------------------------------------------------
614 gSystem->Load("libSTEERBase");
615 gSystem->Load("libESD");
616 gSystem->Load("libAOD");
617 gSystem->Load("libANALYSIS");
618 gSystem->Load("libANALYSISalice");
5d040cf3 619 gSystem->Load("libCORRFW");
620 cerr<<"libCORRFW loaded..."<<endl;
621 gSystem->Load("libPWG2flowCommon");
622 cerr<<"libPWG2flowCommon loaded..."<<endl;
623 gSystem->Load("libPWG2flowTasks");
624 cerr<<"libPWG2flowTasks loaded..."<<endl;
93ff27bd 625 }
626
627 else if (mode == mLocalPAR) {
628 //--------------------------------------------------------
629 //If you want to use root and par files from aliroot
630 //--------------------------------------------------------
631 //If you want to use root and par files from aliroot
632 //--------------------------------------------------------
633 SetupPar("STEERBase");
634 SetupPar("ESD");
635 SetupPar("AOD");
636 SetupPar("ANALYSIS");
637 SetupPar("ANALYSISalice");
638 SetupPar("PWG2AOD");
639 SetupPar("CORRFW");
640 SetupPar("PWG2flowCommon");
641 cerr<<"PWG2flowCommon.par loaded..."<<endl;
642 SetupPar("PWG2flowTasks");
643 cerr<<"PWG2flowTasks.par loaded..."<<endl;
644 }
645
646 //---------------------------------------------------------
647 // <<<<<<<<<< Source mode >>>>>>>>>>>>
648 //---------------------------------------------------------
649 else if (mode==mLocalSource) {
650
651 // In root inline compile
652
653
654 // Constants
655 gROOT->LoadMacro("AliFlowCommon/AliFlowCommonConstants.cxx+");
656 gROOT->LoadMacro("AliFlowCommon/AliFlowLYZConstants.cxx+");
657 gROOT->LoadMacro("AliFlowCommon/AliFlowCumuConstants.cxx+");
658
659 // Flow event
660 gROOT->LoadMacro("AliFlowCommon/AliFlowVector.cxx+");
661 gROOT->LoadMacro("AliFlowCommon/AliFlowTrackSimple.cxx+");
662 gROOT->LoadMacro("AliFlowCommon/AliFlowEventSimple.cxx+");
663
664 // Cuts
665 gROOT->LoadMacro("AliFlowCommon/AliFlowTrackSimpleCuts.cxx+");
666
667 // Output histosgrams
668 gROOT->LoadMacro("AliFlowCommon/AliFlowCommonHist.cxx+");
669 gROOT->LoadMacro("AliFlowCommon/AliFlowCommonHistResults.cxx+");
670 gROOT->LoadMacro("AliFlowCommon/AliFlowLYZHist1.cxx+");
671 gROOT->LoadMacro("AliFlowCommon/AliFlowLYZHist2.cxx+");
672
673 // Functions needed for various methods
674 gROOT->LoadMacro("AliFlowCommon/AliCumulantsFunctions.cxx+");
93ff27bd 675 gROOT->LoadMacro("AliFlowCommon/AliFlowLYZEventPlane.cxx+");
676
677 // Flow Analysis code for various methods
678 gROOT->LoadMacro("AliFlowCommon/AliFlowAnalysisWithMCEventPlane.cxx+");
679 gROOT->LoadMacro("AliFlowCommon/AliFlowAnalysisWithScalarProduct.cxx+");
680 gROOT->LoadMacro("AliFlowCommon/AliFlowAnalysisWithLYZEventPlane.cxx+");
681 gROOT->LoadMacro("AliFlowCommon/AliFlowAnalysisWithLeeYangZeros.cxx+");
682 gROOT->LoadMacro("AliFlowCommon/AliFlowAnalysisWithCumulants.cxx+");
683 gROOT->LoadMacro("AliFlowCommon/AliFlowAnalysisWithQCumulants.cxx+");
ce4a88f5 684 gROOT->LoadMacro("AliFlowCommon/AliFlowAnalysisWithFittingQDistribution.cxx+");
93ff27bd 685
a47dc7e4 686 // Class to fill the FlowEvent on the fly (generate Monte Carlo events)
687 gROOT->LoadMacro("AliFlowCommon/AliFlowEventSimpleMakerOnTheFly.cxx+");
93ff27bd 688
689 cout << "finished loading macros!" << endl;
690
691 }
692
693}
694
695