]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWG2/FLOW/macros/runFlowAnalysisOnTheFly.C
added fHarmonic setter and getter + bug fix for RP's and POI's last in in pt and...
[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
63911185 111int runFlowAnalysisOnTheFly(Int_t mode=mLocal, Int_t nEvts=440)
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
226 // qc->SetCalculate2DFlow(kFALSE); // default
ce4a88f5 227 qc->SetEvaluateNestedLoopsForIntFlow(kFALSE);
228 qc->SetEvaluateNestedLoopsForDiffFlow(kFALSE);
229 qc->Init();
93ff27bd 230 }
231
232 // GFC = Generating Function Cumulants
233 if(GFC) {
234 AliFlowAnalysisWithCumulants* gfc = new AliFlowAnalysisWithCumulants();
7fb984f6 235 if(listWithWeights) gfc->SetWeightsList(listWithWeights);
236 if(usePhiWeights) gfc->SetUsePhiWeights(usePhiWeights);
237 if(usePtWeights) gfc->SetUsePtWeights(usePtWeights);
238 if(useEtaWeights) gfc->SetUseEtaWeights(useEtaWeights);
ce4a88f5 239 gfc->Init();
93ff27bd 240 }
241
242 // FQD = Fitting q-distribution
243 if(FQD) {
ce4a88f5 244 AliFlowAnalysisWithFittingQDistribution* fqd = new AliFlowAnalysisWithFittingQDistribution();
7fb984f6 245 if(listWithWeights) fqd->SetWeightsList(listWithWeights);
246 if(usePhiWeights) fqd->SetUsePhiWeights(usePhiWeights);
ce4a88f5 247 fqd->Init();
93ff27bd 248 }
249
250 // SP = Scalar Product
251 if(SP) {
252 AliFlowAnalysisWithScalarProduct* sp = new AliFlowAnalysisWithScalarProduct();
253 sp->Init();
254 }
255
256 // LYZ1 = Lee-Yang Zeroes first run
c741f5d0 257 if(LYZ1SUM) {
258 AliFlowAnalysisWithLeeYangZeros* lyz1sum = new AliFlowAnalysisWithLeeYangZeros();
259 lyz1sum->SetFirstRun(kTRUE);
260 lyz1sum->SetUseSum(kTRUE);
261 lyz1sum->Init();
262 }
263 if(LYZ1PROD) {
264 AliFlowAnalysisWithLeeYangZeros* lyz1prod = new AliFlowAnalysisWithLeeYangZeros();
265 lyz1prod->SetFirstRun(kTRUE);
266 lyz1prod->SetUseSum(kFALSE);
267 lyz1prod->Init();
93ff27bd 268 }
93ff27bd 269 // LYZ2 = Lee-Yang Zeroes second run
c741f5d0 270 if(LYZ2SUM) {
271 AliFlowAnalysisWithLeeYangZeros* lyz2sum = new AliFlowAnalysisWithLeeYangZeros();
93ff27bd 272 // read the input file from the first run
c741f5d0 273 TString inputFileNameLYZ2SUM = "outputLYZ1SUManalysis.root" ;
274 TFile* inputFileLYZ2SUM = new TFile(inputFileNameLYZ2SUM.Data(),"READ");
275 if(!inputFileLYZ2SUM || inputFileLYZ2SUM->IsZombie()) {
5b40431d 276 cerr << " ERROR: To run LYZ2SUM you need the output file from LYZ1SUM. This file is not there! Please run LYZ1SUM first." << endl ;
93ff27bd 277 break;
278 }
279 else {
c741f5d0 280 TList* inputListLYZ2SUM = (TList*)inputFileLYZ2SUM->Get("cobjLYZ1SUM");
281 if (!inputListLYZ2SUM) {cout<<"Input list LYZ2SUM is NULL pointer!"<<endl; break;}
93ff27bd 282 else {
c741f5d0 283 cout<<"LYZ2SUM input file/list read..."<<endl;
284 lyz2sum->SetFirstRunList(inputListLYZ2SUM);
285 lyz2sum->SetFirstRun(kFALSE);
286 lyz2sum->SetUseSum(kTRUE);
287 lyz2sum->Init();
93ff27bd 288 }
289 }
290 }
c741f5d0 291 if(LYZ2PROD) {
292 AliFlowAnalysisWithLeeYangZeros* lyz2prod = new AliFlowAnalysisWithLeeYangZeros();
293 // read the input file from the first run
294 TString inputFileNameLYZ2PROD = "outputLYZ1PRODanalysis.root" ;
295 TFile* inputFileLYZ2PROD = new TFile(inputFileNameLYZ2PROD.Data(),"READ");
296 if(!inputFileLYZ2PROD || inputFileLYZ2PROD->IsZombie()) {
5b40431d 297 cerr << " ERROR: To run LYZ2PROD you need the output file from LYZ1PROD. This file is not there! Please run LYZ1PROD first." << endl ;
c741f5d0 298 break;
299 }
300 else {
301 TList* inputListLYZ2PROD = (TList*)inputFileLYZ2PROD->Get("cobjLYZ1PROD");
302 if (!inputListLYZ2PROD) {cout<<"Input list LYZ2PROD is NULL pointer!"<<endl; break;}
303 else {
304 cout<<"LYZ2PROD input file/list read..."<<endl;
305 lyz2prod->SetFirstRunList(inputListLYZ2PROD);
306 lyz2prod->SetFirstRun(kFALSE);
307 lyz2prod->SetUseSum(kFALSE);
308 lyz2prod->Init();
309 }
310 }
311 }
312
93ff27bd 313 // LYZEP = Lee-Yang Zeroes event plane
314 if(LYZEP) {
315 AliFlowLYZEventPlane* ep = new AliFlowLYZEventPlane() ;
316 AliFlowAnalysisWithLYZEventPlane* lyzep = new AliFlowAnalysisWithLYZEventPlane();
317 // read the input file from the second lyz run
c741f5d0 318 TString inputFileNameLYZEP = "outputLYZ2SUManalysis.root" ;
93ff27bd 319 TFile* inputFileLYZEP = new TFile(inputFileNameLYZEP.Data(),"READ");
320 if(!inputFileLYZEP || inputFileLYZEP->IsZombie()) {
5b40431d 321 cerr << " ERROR: To run LYZEP you need the output file from LYZ2SUM. This file is not there! Please run LYZ2SUM first." << endl ;
93ff27bd 322 break;
323 }
324 else {
c741f5d0 325 TList* inputListLYZEP = (TList*)inputFileLYZEP->Get("cobjLYZ2SUM");
326 if (!inputListLYZEP) {cout<<"Input list LYZEP is NULL pointer!"<<endl; break;}
93ff27bd 327 else {
328 cout<<"LYZEP input file/list read..."<<endl;
329 ep ->SetSecondRunList(inputListLYZEP);
330 lyzep->SetSecondRunList(inputListLYZEP);
331 ep ->Init();
332 lyzep->Init();
333 }
334 }
335 }
336 //---------------------------------------------------------------------------------------
ce4a88f5 337
338 //---------------------------------------------------------------------------------------
339 // Initialize all the flow methods for additional analysis with different settings/aims
340 // Label each setting/aim with different label !!!!
341
342 // GFC:
343 TString gfcDefaultName = "outputGFCanalysis";
344 // 1.) GFC analysis for elliptic flow with r0 = 1.5:
345 AliFlowAnalysisWithCumulants *gfc_1;
346 TString gfcAnalysisLabels_1 = "_r0_1.5"; // all histograms and output file name will have this label
347 TString gfcOutputFileName_1;
348 gfcOutputFileName_1 = gfcDefaultName.Data();
349 gfcOutputFileName_1 += gfcAnalysisLabels_1.Data();
350 gfcOutputFileName_1 += ".root";
351 if(GFC_Additional_Analysis)
352 {
353 gfc_1 = new AliFlowAnalysisWithCumulants();
354 //gfc_1->SetAnalysisLabel(gfcAnalysisLabels_1.Data());
355 if(listWithWeights) gfc_1->SetWeightsList(listWithWeights);
356 if(usePhiWeights) gfc_1->SetUsePhiWeights(usePhiWeights);
357 if(usePtWeights) gfc_1->SetUsePtWeights(usePtWeights);
358 if(useEtaWeights) gfc_1->SetUseEtaWeights(useEtaWeights);
359 gfc_1->Init();
360 }
361
362 // QC:
363 TString qcDefaultName = "outputQCanalysis";
364 // 1.) QC analysis for directed flow:
365 AliFlowAnalysisWithQCumulants *qc_1;
366 TString qcAnalysisLabels_1 = "_v1"; // all histograms and output file name will have this label
367 TString qcOutputFileName_1;
368 qcOutputFileName_1 = qcDefaultName.Data();
369 qcOutputFileName_1 += qcAnalysisLabels_1.Data();
370 qcOutputFileName_1 += ".root";
371 if(QC_Additional_Analysis)
372 {
373 qc_1 = new AliFlowAnalysisWithQCumulants();
374 //qc_1->SetAnalysisLabel(qcAnalysisLabels_1->Data());
375 if(listWithWeights) qc_1->SetWeightsList(listWithWeights);
376 if(usePhiWeights) qc_1->SetUsePhiWeights(usePhiWeights);
377 qc_1->Init();
378 }
379
380 // FQD:
381 TString fqdDefaultName = "outputFQDanalysis";
382 // 1.) FQD fitting with fixed sigma:
383 AliFlowAnalysisWithFittingQDistribution *fqd_1;
384 TString fqdAnalysisLabels_1 = "_fixedSigma"; // all histograms and output file name will have this label
385 TString fqdOutputFileName_1;
386 fqdOutputFileName_1 = fqdDefaultName.Data();
387 fqdOutputFileName_1 += fqdAnalysisLabels_1.Data();
388 fqdOutputFileName_1 += ".root";
389 if(FQD_Additional_Analysis)
390 {
391 fqd_1 = new AliFlowAnalysisWithFittingQDistribution();
392 //fqd_1->SetAnalysisLabel(fqdAnalysisLabels_1->Data());
393 if(listWithWeights) fqd_1->SetWeightsList(listWithWeights);
394 if(usePhiWeights) fqd_1->SetUsePhiWeights(usePhiWeights);
395 fqd_1->Init();
396 }
397 //---------------------------------------------------------------------------------------
93ff27bd 398
e8c7e66c 399 // set the global event parameters:
aaa62cf0 400 eventMakerOnTheFly->SetNoOfLoops(iLoops);
678f8252 401
402 if(bMultDistrOfRPsIsGauss)
403 {
404 eventMakerOnTheFly->SetMultDistrOfRPsIsGauss(bMultDistrOfRPsIsGauss);
405 eventMakerOnTheFly->SetMultiplicityOfRP(iMultiplicityOfRP);
406 eventMakerOnTheFly->SetMultiplicitySpreadOfRP(dMultiplicitySpreadOfRP);
407 } else
408 {
409 eventMakerOnTheFly->SetMultDistrOfRPsIsGauss(bMultDistrOfRPsIsGauss);
410 eventMakerOnTheFly->SetMinMultOfRP(iMinMultOfRP);
411 eventMakerOnTheFly->SetMaxMultOfRP(iMaxMultOfRP);
412 }
413
cf90787f 414 eventMakerOnTheFly->SetTemperatureOfRP(dTemperatureOfRP);
4d603348 415
e8c7e66c 416 eventMakerOnTheFly->SetV1RP(dV1RP);
417 eventMakerOnTheFly->SetV1SpreadRP(dV1SpreadRP);
4d603348 418 eventMakerOnTheFly->SetV4RP(dV4RP);
419 eventMakerOnTheFly->SetV4SpreadRP(dV4SpreadRP);
420
421 // constant harmonic V2:
422 if(bConstantHarmonics)
423 {
424 eventMakerOnTheFly->SetUseConstantHarmonics(bConstantHarmonics);
55381065 425 if(bV2DistrOfRPsIsGauss)
426 {
427 eventMakerOnTheFly->SetV2DistrOfRPsIsGauss(bV2DistrOfRPsIsGauss);
428 eventMakerOnTheFly->SetV2RP(dV2RP);
429 eventMakerOnTheFly->SetV2SpreadRP(dV2SpreadRP);
430 } else
431 {
432 eventMakerOnTheFly->SetV2DistrOfRPsIsGauss(bV2DistrOfRPsIsGauss);
433 eventMakerOnTheFly->SetMinV2RP(dMinV2RP);
434 eventMakerOnTheFly->SetMaxV2RP(dMaxV2RP);
435 }
4d603348 436 }
e23cbd39 437
4d603348 438 // (pt,eta) dependent harmonic V2:
439 if(!bConstantHarmonics)
440 {
441 eventMakerOnTheFly->SetUseConstantHarmonics(bConstantHarmonics);
442 eventMakerOnTheFly->SetV2RPMax(dV2RPMax);
443 eventMakerOnTheFly->SetPtCutOff(dPtCutOff);
444 }
e23cbd39 445
446 // non-uniform acceptance:
447 if(!uniformAcceptance)
448 {
449 eventMakerOnTheFly->SetFirstSectorPhiMin(phimin1);
450 eventMakerOnTheFly->SetFirstSectorPhiMax(phimax1);
451 eventMakerOnTheFly->SetFirstSectorProbability(p1);
452 eventMakerOnTheFly->SetSecondSectorPhiMin(phimin2);
453 eventMakerOnTheFly->SetSecondSectorPhiMax(phimax2);
454 eventMakerOnTheFly->SetSecondSectorProbability(p2);
455 }
4d603348 456
93ff27bd 457 //---------------------------------------------------------------------------------------
458 // create and analyze events 'on the fly':
459
93ff27bd 460 for(Int_t i=0;i<nEvts;i++) {
461 // creating the event with above settings:
93ff27bd 462 AliFlowEventSimple *event = eventMakerOnTheFly->CreateEventOnTheFly();
93ff27bd 463
464 // analyzing the created event 'on the fly':
93ff27bd 465 // do flow analysis for various methods:
c741f5d0 466 if(MCEP) mcep->Make(event);
467 if(QC) qc->Make(event);
468 if(GFC) gfc->Make(event);
469 if(FQD) fqd->Make(event);
470 if(LYZ1SUM) lyz1sum->Make(event);
471 if(LYZ1PROD)lyz1prod->Make(event);
472 if(LYZ2SUM) lyz2sum->Make(event);
473 if(LYZ2PROD)lyz2prod->Make(event);
474 if(LYZEP) lyzep->Make(event,ep);
475 if(SP) sp->Make(event);
93ff27bd 476
ce4a88f5 477 if(GFC_Additional_Analysis)
478 {
479 // r0 = 1.5:
480 gfc_1->Make(event);
481 }
482 if(QC_Additional_Analysis)
483 {
484 // v1:
485 qc_1->Make(event);
486 }
487 if(FQD_Additional_Analysis)
488 {
489 // fixed sigma:
490 fqd_1->Make(event);
491 }
492
93ff27bd 493 delete event;
494 } // end of for(Int_t i=0;i<nEvts;i++)
495 //---------------------------------------------------------------------------------------
496
497
498
499 //---------------------------------------------------------------------------------------
ce4a88f5 500 // calculating and storing the final results of default flow analysis:
c741f5d0 501 if(MCEP) {mcep->Finish(); mcep->WriteHistograms("outputMCEPanalysis.root");}
502 if(SP) {sp->Finish(); sp->WriteHistograms("outputSPanalysis.root");}
503 if(QC) {qc->Finish(); qc->WriteHistograms("outputQCanalysis.root");}
504 if(GFC) {gfc->Finish(); gfc->WriteHistograms("outputGFCanalysis.root");}
505 if(FQD) {fqd->Finish(); fqd->WriteHistograms("outputFQDanalysis.root");}
506 if(LYZ1SUM) {lyz1sum->Finish(); lyz1sum->WriteHistograms("outputLYZ1SUManalysis.root");}
507 if(LYZ1PROD){lyz1prod->Finish();lyz1prod->WriteHistograms("outputLYZ1PRODanalysis.root");}
508 if(LYZ2SUM) {lyz2sum->Finish(); lyz2sum->WriteHistograms("outputLYZ2SUManalysis.root");}
509 if(LYZ2PROD){lyz2prod->Finish();lyz2prod->WriteHistograms("outputLYZ2PRODanalysis.root");}
510 if(LYZEP) {lyzep->Finish(); lyzep->WriteHistograms("outputLYZEPanalysis.root");}
93ff27bd 511 //---------------------------------------------------------------------------------------
512
ce4a88f5 513 //---------------------------------------------------------------------------------------
514 // calculating and storing the final results of flow analysis with different settings/aims:
515 if(GFC_Additional_Analysis)
516 {
517 // r0 = 1.5:
518 gfc_1->Finish();
519 gfc_1->WriteHistograms(gfcOutputFileName_1.Data());
520 }
521 if(QC_Additional_Analysis)
522 {
523 // v1:
524 qc_1->Finish();
525 qc_1->WriteHistograms(qcOutputFileName_1.Data());
526 }
527 if(FQD_Additional_Analysis)
528 {
529 // fixed sigma:
530 fqd_1->Finish();
531 fqd_1->WriteHistograms(fqdOutputFileName_1.Data());
532 }
533 //---------------------------------------------------------------------------------------
93ff27bd 534
535 cout<<endl;
536 cout<<endl;
537 cout<<" ---- LANDED SUCCESSFULLY ---- "<<endl;
538 cout<<endl;
539
540 timer.Stop();
541 cout << endl;
542 timer.Print();
543}
544
545void SetupPar(char* pararchivename)
546{
547 //Load par files, create analysis libraries
548 //For testing, if par file already decompressed and modified
549 //classes then do not decompress.
550
551 TString cdir(Form("%s", gSystem->WorkingDirectory() )) ;
552 TString parpar(Form("%s.par", pararchivename)) ;
553 if ( gSystem->AccessPathName(parpar.Data()) ) {
554 gSystem->ChangeDirectory(gSystem->Getenv("ALICE_ROOT")) ;
555 TString processline(Form(".! make %s", parpar.Data())) ;
556 gROOT->ProcessLine(processline.Data()) ;
557 gSystem->ChangeDirectory(cdir) ;
558 processline = Form(".! mv /tmp/%s .", parpar.Data()) ;
559 gROOT->ProcessLine(processline.Data()) ;
560 }
561 if ( gSystem->AccessPathName(pararchivename) ) {
562 TString processline = Form(".! tar xvzf %s",parpar.Data()) ;
563 gROOT->ProcessLine(processline.Data());
564 }
565
566 TString ocwd = gSystem->WorkingDirectory();
567 gSystem->ChangeDirectory(pararchivename);
568
569 // check for BUILD.sh and execute
570 if (!gSystem->AccessPathName("PROOF-INF/BUILD.sh")) {
571 printf("*******************************\n");
572 printf("*** Building PAR archive ***\n");
573 cout<<pararchivename<<endl;
574 printf("*******************************\n");
575
576 if (gSystem->Exec("PROOF-INF/BUILD.sh")) {
577 Error("runProcess","Cannot Build the PAR Archive! - Abort!");
578 return -1;
579 }
580 }
581 // check for SETUP.C and execute
582 if (!gSystem->AccessPathName("PROOF-INF/SETUP.C")) {
583 printf("*******************************\n");
584 printf("*** Setup PAR archive ***\n");
585 cout<<pararchivename<<endl;
586 printf("*******************************\n");
587 gROOT->Macro("PROOF-INF/SETUP.C");
588 }
589
590 gSystem->ChangeDirectory(ocwd.Data());
591 printf("Current dir: %s\n", ocwd.Data());
592}
593
594void LoadLibraries(const anaModes mode) {
595
596 //--------------------------------------
597 // Load the needed libraries most of them already loaded by aliroot
598 //--------------------------------------
599 gSystem->Load("libTree.so");
600 gSystem->Load("libGeom.so");
601 gSystem->Load("libVMC.so");
602 gSystem->Load("libXMLIO.so");
603 gSystem->Load("libPhysics.so");
604
605 //----------------------------------------------------------
606 // >>>>>>>>>>> Local mode <<<<<<<<<<<<<<
607 //----------------------------------------------------------
608 if (mode==mLocal) {
609 //--------------------------------------------------------
610 // If you want to use already compiled libraries
611 // in the aliroot distribution
612 //--------------------------------------------------------
613 gSystem->Load("libSTEERBase");
614 gSystem->Load("libESD");
615 gSystem->Load("libAOD");
616 gSystem->Load("libANALYSIS");
617 gSystem->Load("libANALYSISalice");
618 gSystem->Load("libCORRFW.so");
619 cerr<<"libCORRFW.so loaded..."<<endl;
620 gSystem->Load("libPWG2flowCommon.so");
621 cerr<<"libPWG2flowCommon.so loaded..."<<endl;
622 gSystem->Load("libPWG2flowTasks.so");
623 cerr<<"libPWG2flowTasks.so loaded..."<<endl;
624 }
625
626 else if (mode == mLocalPAR) {
627 //--------------------------------------------------------
628 //If you want to use root and par files from aliroot
629 //--------------------------------------------------------
630 //If you want to use root and par files from aliroot
631 //--------------------------------------------------------
632 SetupPar("STEERBase");
633 SetupPar("ESD");
634 SetupPar("AOD");
635 SetupPar("ANALYSIS");
636 SetupPar("ANALYSISalice");
637 SetupPar("PWG2AOD");
638 SetupPar("CORRFW");
639 SetupPar("PWG2flowCommon");
640 cerr<<"PWG2flowCommon.par loaded..."<<endl;
641 SetupPar("PWG2flowTasks");
642 cerr<<"PWG2flowTasks.par loaded..."<<endl;
643 }
644
645 //---------------------------------------------------------
646 // <<<<<<<<<< Source mode >>>>>>>>>>>>
647 //---------------------------------------------------------
648 else if (mode==mLocalSource) {
649
650 // In root inline compile
651
652
653 // Constants
654 gROOT->LoadMacro("AliFlowCommon/AliFlowCommonConstants.cxx+");
655 gROOT->LoadMacro("AliFlowCommon/AliFlowLYZConstants.cxx+");
656 gROOT->LoadMacro("AliFlowCommon/AliFlowCumuConstants.cxx+");
657
658 // Flow event
659 gROOT->LoadMacro("AliFlowCommon/AliFlowVector.cxx+");
660 gROOT->LoadMacro("AliFlowCommon/AliFlowTrackSimple.cxx+");
661 gROOT->LoadMacro("AliFlowCommon/AliFlowEventSimple.cxx+");
662
663 // Cuts
664 gROOT->LoadMacro("AliFlowCommon/AliFlowTrackSimpleCuts.cxx+");
665
666 // Output histosgrams
667 gROOT->LoadMacro("AliFlowCommon/AliFlowCommonHist.cxx+");
668 gROOT->LoadMacro("AliFlowCommon/AliFlowCommonHistResults.cxx+");
669 gROOT->LoadMacro("AliFlowCommon/AliFlowLYZHist1.cxx+");
670 gROOT->LoadMacro("AliFlowCommon/AliFlowLYZHist2.cxx+");
671
672 // Functions needed for various methods
673 gROOT->LoadMacro("AliFlowCommon/AliCumulantsFunctions.cxx+");
93ff27bd 674 gROOT->LoadMacro("AliFlowCommon/AliFlowLYZEventPlane.cxx+");
675
676 // Flow Analysis code for various methods
677 gROOT->LoadMacro("AliFlowCommon/AliFlowAnalysisWithMCEventPlane.cxx+");
678 gROOT->LoadMacro("AliFlowCommon/AliFlowAnalysisWithScalarProduct.cxx+");
679 gROOT->LoadMacro("AliFlowCommon/AliFlowAnalysisWithLYZEventPlane.cxx+");
680 gROOT->LoadMacro("AliFlowCommon/AliFlowAnalysisWithLeeYangZeros.cxx+");
681 gROOT->LoadMacro("AliFlowCommon/AliFlowAnalysisWithCumulants.cxx+");
682 gROOT->LoadMacro("AliFlowCommon/AliFlowAnalysisWithQCumulants.cxx+");
ce4a88f5 683 gROOT->LoadMacro("AliFlowCommon/AliFlowAnalysisWithFittingQDistribution.cxx+");
93ff27bd 684
a47dc7e4 685 // Class to fill the FlowEvent on the fly (generate Monte Carlo events)
686 gROOT->LoadMacro("AliFlowCommon/AliFlowEventSimpleMakerOnTheFly.cxx+");
93ff27bd 687
688 cout << "finished loading macros!" << endl;
689
690 }
691
692}
693
694