]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWGCF/FLOW/macros/runFlowAnalysisOnTheFly.C
from Redmer Bertens:
[u/mrichter/AliRoot.git] / PWGCF / FLOW / macros / runFlowAnalysisOnTheFly.C
CommitLineData
93417736 1// Settings for the simulation of events 'on the fly':
2// a) Determine how many events you want to create;
3// b) Set random or same seed for random generator;
4// c) Determine multiplicites of events;
5// d) Parametrize the phi distribution;
6// d1) Enable/disable uniform event-wise fluctuations of v2;
7// d2) Enable/diable pt dependence of v2;
8// e) Parametrize the pt distribution;
9// f) Determine how many times each sampled particle will be taken (simulating nonflow);
10// g) Configure detector's acceptance;
11// h) Decide which flow analysis methods you will use;
12// i) Define simple cuts for Reference Particle (RP) selection;
13// j) Define simple cuts for Particle of Interest (POI) selection;
14// k) Define the ranges for two subevents separated with eta gap (needed only for SP method);
15// l) Enable/disable usage of particle weights.
7fb984f6 16
93417736 17// a) Determine how many events you want to create:
18Int_t iNevts = 1000; // total statistics
7a01f4a7 19
93417736 20// b) Set random or same seed for random generator:
21Bool_t bSameSeed = kFALSE; // if kTRUE, the created events are the same when re-doing flow analysis 'on the fly'
e1911c19 22
93417736 23// c) Determine multiplicites of events:
24// Remark 1: Multiplicity M for each event is sampled uniformly from interval iMinMult <= M < iMaxMult;
25// Remark 2: For constant M of e.g. 500 for each event, set iMinMult = 500 and iMaxMult = 501.
26Int_t iMinMult = 500; // uniformly sampled multiplicity is >= iMinMult
27Int_t iMaxMult = 501; // uniformly sampled multiplicity is < iMaxMult
4d603348 28
93417736 29// d) Parametrize the phi distribution:
30// Remark 1: Hardwired is Fourier-like distribution f(phi) = (1/2pi)(1+sum_{n=1}^{4} 2v_n cos[n(phi-rp)]),
31// where reaction plane (rp) is sampled uniformly for each event from interval [0,2pi]
32Double_t dV1 = 0.0; // constant harmonic v1
33Double_t dV2 = 0.05; // constant harmonic v2
34Double_t dV3 = 0.0; // constant harmonic v3
35Double_t dV4 = 0.0; // constant harmonic v4
dd1d2f85 36Double_t dV5 = 0.0; // constant harmonic v5
37Double_t dV6 = 0.0; // constant harmonic v6
93417736 38// Remark 2: By default all harmonics are constant for each event and for each particle. However, for v2
39// the uniform event-wise fluctuations or pt dependence can be enabled:
40// d1) Enable/disable uniform event-wise fluctuations of v2:
41Bool_t bUniformFluctuationsV2 = kFALSE; // enable uniform event-wise flow fluctuations (set than also dMinV2 and dMaxV2 bellow)
42Double_t dMinV2 = 0.04; // lower boundary on v2, when bUniformFluctuationsV2 = kTRUE
43Double_t dMaxV2 = 0.06; // upper boundary on v2, when bUniformFluctuationsV2 = kTRUE
44// d2) Enable/disable pt dependence of v2:
45Bool_t bPtDependentV2 = kFALSE; // enable pt dependence of v2 (set then also dV2vsPtMax and dV2vsPtCutOff bellow)
46Double_t dV2vsPtCutOff = 2.0; // up to pt = dV2vsPtCutOff v2 is growing linearly as a function of pt
47Double_t dV2vsPtMax = 0.20; // for pt >= dV2vsPtCutOff, v2(pt) = dV2vsPtMax
93ff27bd 48
93417736 49// e) Parametrize the pt distribution:
50// Remark: Hardwired is Boltzmann distribution f(pt) = pt*exp[-sqrt(dMass^2+pt^2)/dT]
51Double_t dMass = 0.13957; // mass in GeV/c^2 (e.g. m_{pions} = 0.13957)
52Double_t dTemperature = 0.44; // "temperature" in GeV/c (increase this parameter to get more high pt particles)
598d292e 53
93417736 54// f) Determine how many times each sampled particle will be taken in the analysis (simulating nonflow):
55Int_t nTimes = 1; // e.g. for nTimes = 2, strong 2-particle nonflow correlations are introduced
7a01f4a7 56
93417736 57// g) Configure detector's acceptance:
58Bool_t uniformAcceptance = kTRUE; // if kTRUE: detectors has uniform azimuthal acceptance.
59 // if kFALSE: you will simulate detector with non-uniform acceptance in one or
60 // two sectors. For each sector you specify phiMin, phiMax and probability p.
61 // Then all particles emitted in direction phiMin < phi < phiMax will be taken
62 // with probability p. If p = 0, that sector is completely blocked. Set bellow
63 // phiMin1, phiMax1, p1 for the first sector and phiMin2, phiMax2, p2 for the second
64 // sector. If you set phiMin2 = phiMax2 = p2 = 0, only first non-uniform sector is
65 // simulated.
66// 1st non-uniform sector:
67Double_t phiMin1 = 60; // first non-uniform sector starts at this azimuth (in degrees)
68Double_t phiMax1 = 120; // first non-uniform sector ends at this azimuth (in degrees)
69Double_t p1 = 0.5; // probablitity that particles emitted in [phiMin1,phiMax1] are taken
70// 2nd non-uniform sector:
71Double_t phiMin2 = 0.; // first non-uniform sector starts at this azimuth (in degrees)
72Double_t phiMax2 = 0.; // first non-uniform sector ends at this azimuth (in degrees)
73Double_t p2 = 0.; // probablitity that particles emitted in [phiMin2,phiMax2] are taken
7a01f4a7 74
93417736 75// h) Decide which flow analysis methods you will use:
76Bool_t MCEP = kTRUE; // Monte Carlo Event Plane
77Bool_t SP = kTRUE; // Scalar Product (a.k.a 'flow analysis with eta gaps')
78Bool_t GFC = kTRUE; // Generating Function Cumulants
79Bool_t QC = kTRUE; // Q-cumulants
80Bool_t FQD = kTRUE; // Fitted q-distribution
81Bool_t LYZ1SUM = kTRUE; // Lee-Yang Zero (sum generating function), first pass over the data
82Bool_t LYZ1PROD = kTRUE; // Lee-Yang Zero (product generating function), first pass over the data
83Bool_t LYZ2SUM = kFALSE; // Lee-Yang Zero (sum generating function), second pass over the data
84Bool_t LYZ2PROD = kFALSE; // Lee-Yang Zero (product generating function), second pass over the data
85Bool_t LYZEP = kFALSE; // Lee-Yang Zero Event Plane
86Bool_t MH = kFALSE; // Mixed Harmonics (used for strong parity violation studies)
87Bool_t NL = kFALSE; // Nested Loops (neeed for debugging, only for developers)
e23cbd39 88
93417736 89// i) Define simple cuts for Reference Particle (RP) selection:
90Double_t ptMinRP = 0.0; // in GeV
91Double_t ptMaxRP = 10.0; // in GeV
92Double_t etaMinRP = -1.;
93Double_t etaMaxRP = 1.;
94Double_t phiMinRP = 0.0; // in degrees
95Double_t phiMaxRP = 360.0; // in degrees
96Bool_t bUseChargeRP = kFALSE; // if kFALSE, RPs with both sign of charges are taken
97Int_t chargeRP = 1; // +1 or -1
e23cbd39 98
93417736 99// j) Define simple cuts for Particle of Interest (POI) selection:
100Double_t ptMinPOI = 0.0; // in GeV
101Double_t ptMaxPOI = 10.0; // in GeV
102Double_t etaMinPOI = -1.; //
103Double_t etaMaxPOI = 1.;
104Double_t phiMinPOI = 0.0; // in degrees
105Double_t phiMaxPOI = 360.0; // in degrees
106Bool_t bUseChargePOI = kFALSE; // if kFALSE, POIs with both sign of charges are taken
107Int_t chargePOI = -1; // +1 or -1
7a01f4a7 108
93417736 109// k) Define the ranges for two subevents separated with eta gap (needed only for SP method):
110Double_t etaMinA = -0.8; // minimum eta of subevent A
111Double_t etaMaxA = -0.5; // maximum eta of subevent A
112Double_t etaMinB = 0.5; // minimum eta of subevent B
113Double_t etaMaxB = 0.8; // maximum eta of subevent B
524798bf 114
93417736 115// l) Enable/disable usage of particle weights:
116Bool_t usePhiWeights = kFALSE; // phi weights
117Bool_t usePtWeights = kFALSE; // pt weights
118Bool_t useEtaWeights = kFALSE; // eta weights
e23cbd39 119
4d603348 120enum anaModes {mLocal,mLocalSource,mLocalPAR};
93ff27bd 121// mLocal: Analyze data on your computer using aliroot
122// mLocalPAR: Analyze data on your computer using root + PAR files
123// mLocalSource: Analyze data on your computer using root + source files
124
93417736 125#include "TStopwatch.h"
126#include "TObjArray"
127#include "Riostream.h"
128#include "TFile.h"
129
130int runFlowAnalysisOnTheFly(Int_t mode=mLocal)
93ff27bd 131{
93417736 132 // Beging analysis 'on the fly'.
93ff27bd 133
93417736 134 // a) Formal necessities....;
135 // b) Initialize the flow event maker 'on the fly';
136 // c) If enabled, access particle weights from external file;
137 // d) Configure the flow analysis methods;
138 // e) Simple cuts for RPs;
139 // f) Simple cuts for POIs;
140 // g) Create and analyse events 'on the fly';
141 // h) Create the output file and directory structure for the final results of all methods;
142 // i) Calculate and store the final results of all methods.
93ff27bd 143
93417736 144 // a) Formal necessities....:
145 CheckUserSettings();
146 WelcomeMessage();
147 TStopwatch timer;
148 timer.Start();
93ff27bd 149 LoadLibraries(mode);
4d603348 150
93417736 151 // b) Initialize the flow event maker 'on the fly':
152 UInt_t uiSeed = 0; // if uiSeed is 0, the seed is determined uniquely in space and time via TUUID
153 if(bSameSeed){uiSeed = 44;}
154 AliFlowEventSimpleMakerOnTheFly* eventMakerOnTheFly = new AliFlowEventSimpleMakerOnTheFly(uiSeed);
155 eventMakerOnTheFly->SetMinMult(iMinMult);
156 eventMakerOnTheFly->SetMaxMult(iMaxMult);
157 eventMakerOnTheFly->SetMass(dMass);
158 eventMakerOnTheFly->SetTemperature(dTemperature);
159 eventMakerOnTheFly->SetV1(dV1);
160 eventMakerOnTheFly->SetV2(dV2);
161 eventMakerOnTheFly->SetV3(dV3);
162 eventMakerOnTheFly->SetV4(dV4);
dd1d2f85 163 eventMakerOnTheFly->SetV5(dV5);
164 eventMakerOnTheFly->SetV6(dV6);
93417736 165 if(bUniformFluctuationsV2)
4d603348 166 {
93417736 167 eventMakerOnTheFly->SetUniformFluctuationsV2(bUniformFluctuationsV2);
168 eventMakerOnTheFly->SetMinV2(dMinV2);
169 eventMakerOnTheFly->SetMaxV2(dMaxV2);
170 }
171 if(bPtDependentV2)
172 {
173 eventMakerOnTheFly->SetPtDependentV2(bPtDependentV2);
174 eventMakerOnTheFly->SetV2vsPtCutOff(dV2vsPtCutOff);
175 eventMakerOnTheFly->SetV2vsPtMax(dV2vsPtMax);
4d603348 176 }
93417736 177 eventMakerOnTheFly->SetSubeventEtaRange(etaMinA,etaMaxA,etaMinB,etaMaxB);
178 eventMakerOnTheFly->SetNTimes(nTimes);
179 if(!uniformAcceptance)
4d603348 180 {
93417736 181 eventMakerOnTheFly->SetUniformAcceptance(kFALSE);
182 eventMakerOnTheFly->SetFirstSectorPhiMin(phiMin1);
183 eventMakerOnTheFly->SetFirstSectorPhiMax(phiMax1);
184 eventMakerOnTheFly->SetFirstSectorProbability(p1);
185 eventMakerOnTheFly->SetSecondSectorPhiMin(phiMin2);
186 eventMakerOnTheFly->SetSecondSectorPhiMax(phiMax2);
187 eventMakerOnTheFly->SetSecondSectorProbability(p2);
188 }
4d603348 189 eventMakerOnTheFly->Init();
93417736 190
191 // c) If enabled, access particle weights from external file:
192 TFile *fileWithWeights = NULL;
193 TList *listWithWeights = NULL;
194 if(usePhiWeights||usePtWeights||useEtaWeights)
195 {
196 fileWithWeights = TFile::Open("weights.root","READ");
197 if(fileWithWeights)
198 {
199 listWithWeights = (TList*)fileWithWeights->Get("weights");
200 }
201 else
202 {
203 cout << " WARNING: the file <weights.root> with weights from the previous run was not found."<<endl;
204 break;
205 }
206 } // end of if(usePhiWeights||usePtWeights||useEtaWeights)
207
208 // d) Configure the flow analysis methods:
ce4a88f5 209 AliFlowAnalysisWithQCumulants *qc = NULL;
210 AliFlowAnalysisWithCumulants *gfc = NULL;
211 AliFlowAnalysisWithFittingQDistribution *fqd = NULL;
212 AliFlowAnalysisWithLeeYangZeros *lyz1sum = NULL;
213 AliFlowAnalysisWithLeeYangZeros *lyz1prod = NULL;
214 AliFlowAnalysisWithLeeYangZeros *lyz2sum = NULL;
215 AliFlowAnalysisWithLeeYangZeros *lyz2prod = NULL;
216 AliFlowAnalysisWithLYZEventPlane *lyzep = NULL;
217 AliFlowAnalysisWithScalarProduct *sp = NULL;
d66d46f7 218 AliFlowAnalysisWithMixedHarmonics *mh = NULL;
83bc3e95 219 AliFlowAnalysisWithNestedLoops *nl = NULL;
ce4a88f5 220 AliFlowAnalysisWithMCEventPlane *mcep = NULL;
93ff27bd 221 // MCEP = monte carlo event plane
93417736 222 if(MCEP)
d66d46f7 223 {
93417736 224 //AliFlowAnalysisWithMCEventPlane *mcep = new AliFlowAnalysisWithMCEventPlane();
225 mcep = new AliFlowAnalysisWithMCEventPlane();
226 mcep->SetHarmonic(2); // default is v2
227 mcep->Init();
228 } // end of if(MCEP)
229 // SP = Scalar Product
230 if(SP)
d66d46f7 231 {
93417736 232 sp = new AliFlowAnalysisWithScalarProduct();
233 if(listWithWeights){sp->SetWeightsList(listWithWeights);}
234 sp->SetUsePhiWeights(usePhiWeights);
235 sp->SetHarmonic(2);
236 sp->SetApplyCorrectionForNUA(kFALSE);
237 sp->Init();
238 } // end of if(SP)
ce4a88f5 239 // QC = Q-cumulants
93417736 240 if(QC)
241 {
242 qc = new AliFlowAnalysisWithQCumulants();
243 if(listWithWeights){qc->SetWeightsList(listWithWeights);}
244 if(usePhiWeights){qc->SetUsePhiWeights(usePhiWeights);}
245 if(usePtWeights){qc->SetUsePtWeights(usePtWeights);}
246 if(useEtaWeights){qc->SetUseEtaWeights(useEtaWeights);}
247 qc->SetHarmonic(2);
248 qc->SetCalculateDiffFlow(kTRUE);
249 qc->SetCalculate2DDiffFlow(kFALSE); // vs (pt,eta)
250 qc->SetApplyCorrectionForNUA(kFALSE);
251 qc->SetFillMultipleControlHistograms(kFALSE);
252 qc->SetMultiplicityWeight("combinations"); // default (other supported options are "unit" and "multiplicity")
253 qc->SetCalculateCumulantsVsM(kFALSE);
62e36168 254 qc->SetCalculateAllCorrelationsVsM(kFALSE); // calculate all correlations in mixed harmonics "vs M"
93417736 255 qc->SetnBinsMult(10000);
256 qc->SetMinMult(0);
257 qc->SetMaxMult(10000);
62e36168 258 qc->SetBookOnlyBasicCCH(kFALSE); // book only basic common control histograms
259 qc->SetCalculateDiffFlowVsEta(kTRUE); // if you set kFALSE only differential flow vs pt is calculated
dd1d2f85 260 qc->SetCalculateMixedHarmonics(kFALSE); // calculate all multi-partice mixed-harmonics correlators
93417736 261 qc->Init();
262 } // end of if(QC)
93ff27bd 263 // GFC = Generating Function Cumulants
93417736 264 if(GFC)
265 {
266 gfc = new AliFlowAnalysisWithCumulants();
267 if(listWithWeights){gfc->SetWeightsList(listWithWeights);}
268 if(usePhiWeights){gfc->SetUsePhiWeights(usePhiWeights);}
269 if(usePtWeights){gfc->SetUsePtWeights(usePtWeights);}
270 if(useEtaWeights){gfc->SetUseEtaWeights(useEtaWeights);}
271 // calculation vs multiplicity:
272 gfc->SetCalculateVsMultiplicity(kFALSE);
273 gfc->SetnBinsMult(10000);
274 gfc->SetMinMult(0);
275 gfc->SetMaxMult(10000);
276 gfc->SetHarmonic(2);
277 // tuning of interpolating parameters:
278 gfc->SetTuneParameters(kFALSE);
279 Double_t r0[10] = {1.8,1.9,2.0,2.1,2.2,2.3,2.4,2.5,2.6,2.7}; // up to 10 values allowed
280 for(Int_t r=0;r<10;r++){gfc->SetTuningR0(r0[r],r);}
281 gfc->Init();
282 } // end of if(GFC)
93ff27bd 283 // FQD = Fitting q-distribution
93417736 284 if(FQD)
285 {
286 fqd = new AliFlowAnalysisWithFittingQDistribution();
287 if(listWithWeights){fqd->SetWeightsList(listWithWeights);}
288 if(usePhiWeights){fqd->SetUsePhiWeights(usePhiWeights);}
289 fqd->SetHarmonic(2);
290 fqd->Init();
291 } // end of if(FQD)
93ff27bd 292 // LYZ1 = Lee-Yang Zeroes first run
93417736 293 if(LYZ1SUM)
294 {
295 lyz1sum = new AliFlowAnalysisWithLeeYangZeros();
296 lyz1sum->SetFirstRun(kTRUE);
297 lyz1sum->SetUseSum(kTRUE);
298 lyz1sum->Init();
299 } // end of if(LYZ1SUM)
300 if(LYZ1PROD)
301 {
302 lyz1prod = new AliFlowAnalysisWithLeeYangZeros();
303 lyz1prod->SetFirstRun(kTRUE);
304 lyz1prod->SetUseSum(kFALSE);
305 lyz1prod->Init();
306 } // end of if(LYZ1PROD)
93ff27bd 307 // LYZ2 = Lee-Yang Zeroes second run
93417736 308 if(LYZ2SUM)
309 {
310 lyz2sum = new AliFlowAnalysisWithLeeYangZeros();
311 // read the input file from the first run
312 TString inputFileNameLYZ2SUM = "outputLYZ1SUManalysis.root" ;
313 TFile* inputFileLYZ2SUM = new TFile(inputFileNameLYZ2SUM.Data(),"READ");
314 if(!inputFileLYZ2SUM || inputFileLYZ2SUM->IsZombie())
315 {
316 cerr<<" ERROR: To run LYZ2SUM you need the output file from LYZ1SUM. This file is not there! Please run LYZ1SUM first."<<endl;
317 break;
318 } else
319 {
c741f5d0 320 TList* inputListLYZ2SUM = (TList*)inputFileLYZ2SUM->Get("cobjLYZ1SUM");
93417736 321 if(!inputListLYZ2SUM){cout<<"Input list LYZ2SUM is NULL pointer!"<<endl;break;}
322 else
323 {
324 cout<<"LYZ2SUM input file/list read..."<<endl;
325 lyz2sum->SetFirstRunList(inputListLYZ2SUM);
326 lyz2sum->SetFirstRun(kFALSE);
327 lyz2sum->SetUseSum(kTRUE);
328 lyz2sum->Init();
93ff27bd 329 }
330 }
93417736 331 } // end of if(LYZ2SUM)
332 if(LYZ2PROD)
333 {
334 lyz2prod = new AliFlowAnalysisWithLeeYangZeros();
335 // read the input file from the first run
336 TString inputFileNameLYZ2PROD = "outputLYZ1PRODanalysis.root" ;
337 TFile* inputFileLYZ2PROD = new TFile(inputFileNameLYZ2PROD.Data(),"READ");
338 if(!inputFileLYZ2PROD || inputFileLYZ2PROD->IsZombie())
339 {
340 cerr<<" ERROR: To run LYZ2PROD you need the output file from LYZ1PROD. This file is not there! Please run LYZ1PROD first."<<endl;
341 break;
342 } else
343 {
c741f5d0 344 TList* inputListLYZ2PROD = (TList*)inputFileLYZ2PROD->Get("cobjLYZ1PROD");
93417736 345 if(!inputListLYZ2PROD){cout<<"Input list LYZ2PROD is NULL pointer!"<<endl;break;}
346 else
347 {
348 cout<<"LYZ2PROD input file/list read..."<<endl;
349 lyz2prod->SetFirstRunList(inputListLYZ2PROD);
350 lyz2prod->SetFirstRun(kFALSE);
351 lyz2prod->SetUseSum(kFALSE);
352 lyz2prod->Init();
c741f5d0 353 }
354 }
93417736 355 } // end of if(LYZ2PROD)
93ff27bd 356 // LYZEP = Lee-Yang Zeroes event plane
93417736 357 if(LYZEP)
358 {
359 AliFlowLYZEventPlane *ep = new AliFlowLYZEventPlane() ;
360 AliFlowAnalysisWithLYZEventPlane* lyzep = new AliFlowAnalysisWithLYZEventPlane();
93ff27bd 361 // read the input file from the second lyz run
c741f5d0 362 TString inputFileNameLYZEP = "outputLYZ2SUManalysis.root" ;
93ff27bd 363 TFile* inputFileLYZEP = new TFile(inputFileNameLYZEP.Data(),"READ");
364 if(!inputFileLYZEP || inputFileLYZEP->IsZombie()) {
5b40431d 365 cerr << " ERROR: To run LYZEP you need the output file from LYZ2SUM. This file is not there! Please run LYZ2SUM first." << endl ;
93ff27bd 366 break;
367 }
368 else {
c741f5d0 369 TList* inputListLYZEP = (TList*)inputFileLYZEP->Get("cobjLYZ2SUM");
370 if (!inputListLYZEP) {cout<<"Input list LYZEP is NULL pointer!"<<endl; break;}
93ff27bd 371 else {
372 cout<<"LYZEP input file/list read..."<<endl;
373 ep ->SetSecondRunList(inputListLYZEP);
374 lyzep->SetSecondRunList(inputListLYZEP);
375 ep ->Init();
376 lyzep->Init();
377 }
378 }
379 }
93417736 380 // Mixed Harmonics:
381 if(MH)
bdb5c432 382 {
93417736 383 mh = new AliFlowAnalysisWithMixedHarmonics();
384 mh->SetHarmonic(1); // integer n in expression cos[n(2phi1-phi2-phi3)] = v2n*vn^2
385 mh->SetMinMultiplicity(100);
386 mh->SetNoOfMultipicityBins(5);
387 mh->SetMultipicityBinWidth(200);
388 mh->Init();
389 } // end of if(MH)
390 // NL = Nested Loops:
391 if(NL)
e23cbd39 392 {
93417736 393 nl = new AliFlowAnalysisWithNestedLoops();
394 nl->Init();
395 } // end of if(NL)
e1911c19 396
93417736 397 // e) Simple cuts for RPs:
e1911c19 398 AliFlowTrackSimpleCuts *cutsRP = new AliFlowTrackSimpleCuts();
399 cutsRP->SetPtMax(ptMaxRP);
400 cutsRP->SetPtMin(ptMinRP);
401 cutsRP->SetEtaMax(etaMaxRP);
402 cutsRP->SetEtaMin(etaMinRP);
403 cutsRP->SetPhiMax(phiMaxRP*TMath::Pi()/180.);
404 cutsRP->SetPhiMin(phiMinRP*TMath::Pi()/180.);
93417736 405 if(bUseChargeRP){cutsRP->SetCharge(chargeRP);}
406
407 // f) Simple cuts for POIs:
e1911c19 408 AliFlowTrackSimpleCuts *cutsPOI = new AliFlowTrackSimpleCuts();
409 cutsPOI->SetPtMax(ptMaxPOI);
410 cutsPOI->SetPtMin(ptMinPOI);
411 cutsPOI->SetEtaMax(etaMaxPOI);
412 cutsPOI->SetEtaMin(etaMinPOI);
413 cutsPOI->SetPhiMax(phiMaxPOI*TMath::Pi()/180.);
414 cutsPOI->SetPhiMin(phiMinPOI*TMath::Pi()/180.);
93417736 415 if(bUseChargePOI){cutsPOI->SetCharge(chargePOI);}
e1911c19 416
93417736 417 // g) Create and analyse events 'on the fly':
418 for(Int_t i=0;i<iNevts;i++)
419 {
420 // Creating the event 'on the fly':
421 AliFlowEventSimple *event = eventMakerOnTheFly->CreateEventOnTheFly(cutsRP,cutsPOI);
422 // Passing the created event to flow analysis methods:
423 if(MCEP){mcep->Make(event);}
424 if(QC){qc->Make(event);}
425 if(GFC){gfc->Make(event);}
426 if(FQD){fqd->Make(event);}
427 if(LYZ1SUM){lyz1sum->Make(event);}
428 if(LYZ1PROD){lyz1prod->Make(event);}
429 if(LYZ2SUM){lyz2sum->Make(event);}
430 if(LYZ2PROD){lyz2prod->Make(event);}
431 if(LYZEP){lyzep->Make(event,ep);}
432 if(SP){sp->Make(event);}
433 if(MH){mh->Make(event);}
434 if(NL){nl->Make(event);}
435 delete event;
436 } // end of for(Int_t i=0;i<iNevts;i++)
93ff27bd 437
93417736 438 // h) Create the output file and directory structure for the final results of all methods:
ad87ae62 439 TString outputFileName = "AnalysisResults.root";
440 TFile *outputFile = new TFile(outputFileName.Data(),"RECREATE");
93417736 441 const Int_t nMethods = 12;
442 TString method[nMethods] = {"MCEP","SP","GFC","QC","FQD","LYZ1SUM","LYZ1PROD","LYZ2SUM","LYZ2PROD","LYZEP","MH","NL"};
ad87ae62 443 TDirectoryFile *dirFileFinal[nMethods] = {NULL};
444 TString fileName[nMethods];
445 for(Int_t i=0;i<nMethods;i++)
446 {
ad87ae62 447 fileName[i]+="output";
448 fileName[i]+=method[i].Data();
449 fileName[i]+="analysis";
450 dirFileFinal[i] = new TDirectoryFile(fileName[i].Data(),fileName[i].Data());
451 }
452
93417736 453 // i) Calculate and store the final results of all methods:
454 if(MCEP){mcep->Finish();mcep->WriteHistograms(dirFileFinal[0]);}
455 if(SP){sp->Finish();sp->WriteHistograms(dirFileFinal[1]);}
456 if(GFC){gfc->Finish();gfc->WriteHistograms(dirFileFinal[2]);}
457 if(QC){qc->Finish();qc->WriteHistograms(dirFileFinal[3]);}
458 if(FQD){fqd->Finish();fqd->WriteHistograms(dirFileFinal[4]);}
459 if(LYZ1SUM){lyz1sum->Finish();lyz1sum->WriteHistograms(dirFileFinal[5]);}
ad87ae62 460 if(LYZ1PROD){lyz1prod->Finish();lyz1prod->WriteHistograms(dirFileFinal[6]);}
93417736 461 if(LYZ2SUM){lyz2sum->Finish();lyz2sum->WriteHistograms(dirFileFinal[7]);}
ad87ae62 462 if(LYZ2PROD){lyz2prod->Finish();lyz2prod->WriteHistograms(dirFileFinal[8]);}
93417736 463 if(LYZEP){lyzep->Finish();lyzep->WriteHistograms(dirFileFinal[9]);}
464 if(MH){mh->Finish();mh->WriteHistograms(dirFileFinal[10]);}
465 if(NL){nl->Finish();nl->WriteHistograms(dirFileFinal[11]);}
93ff27bd 466
ad87ae62 467 outputFile->Close();
468 delete outputFile;
93ff27bd 469
470 cout<<endl;
471 cout<<endl;
472 cout<<" ---- LANDED SUCCESSFULLY ---- "<<endl;
473 cout<<endl;
474
475 timer.Stop();
476 cout << endl;
477 timer.Print();
93417736 478
479} // end of int runFlowAnalysisOnTheFly(Int_t mode=mLocal)
93ff27bd 480
481void SetupPar(char* pararchivename)
482{
483 //Load par files, create analysis libraries
484 //For testing, if par file already decompressed and modified
485 //classes then do not decompress.
486
487 TString cdir(Form("%s", gSystem->WorkingDirectory() )) ;
488 TString parpar(Form("%s.par", pararchivename)) ;
489 if ( gSystem->AccessPathName(parpar.Data()) ) {
490 gSystem->ChangeDirectory(gSystem->Getenv("ALICE_ROOT")) ;
491 TString processline(Form(".! make %s", parpar.Data())) ;
492 gROOT->ProcessLine(processline.Data()) ;
493 gSystem->ChangeDirectory(cdir) ;
494 processline = Form(".! mv /tmp/%s .", parpar.Data()) ;
495 gROOT->ProcessLine(processline.Data()) ;
496 }
497 if ( gSystem->AccessPathName(pararchivename) ) {
498 TString processline = Form(".! tar xvzf %s",parpar.Data()) ;
499 gROOT->ProcessLine(processline.Data());
500 }
501
502 TString ocwd = gSystem->WorkingDirectory();
503 gSystem->ChangeDirectory(pararchivename);
504
505 // check for BUILD.sh and execute
506 if (!gSystem->AccessPathName("PROOF-INF/BUILD.sh")) {
507 printf("*******************************\n");
508 printf("*** Building PAR archive ***\n");
509 cout<<pararchivename<<endl;
510 printf("*******************************\n");
511
512 if (gSystem->Exec("PROOF-INF/BUILD.sh")) {
513 Error("runProcess","Cannot Build the PAR Archive! - Abort!");
514 return -1;
515 }
516 }
517 // check for SETUP.C and execute
518 if (!gSystem->AccessPathName("PROOF-INF/SETUP.C")) {
519 printf("*******************************\n");
520 printf("*** Setup PAR archive ***\n");
521 cout<<pararchivename<<endl;
522 printf("*******************************\n");
523 gROOT->Macro("PROOF-INF/SETUP.C");
524 }
525
526 gSystem->ChangeDirectory(ocwd.Data());
527 printf("Current dir: %s\n", ocwd.Data());
528}
529
bdb5c432 530void CheckUserSettings()
531{
93417736 532 // Check if user settings make sense before taking off.
bdb5c432 533
93417736 534 if(iNevts <= 0)
bdb5c432 535 {
93417736 536 printf("\n WARNING: nEvts <= 0 !!!! Please check your settings before taking off.\n\n");
537 exit(0);
538 }
539 if(iMinMult < 0.)
540 {
541 printf("\n WARNING: iMinMult < 0 !!!! Please check your settings before taking off.\n\n");
542 exit(0);
bdb5c432 543 }
93417736 544 if(iMaxMult <= 0.)
bdb5c432 545 {
93417736 546 printf("\n WARNING: iMaxMult <= 0 !!!! Please check your settings before taking off.\n\n");
547 exit(0);
bdb5c432 548 }
93417736 549 if(iMinMult >= iMaxMult)
bdb5c432 550 {
93417736 551 printf("\n WARNING: iMinMult >= iMaxMult !!!! Please check your settings before taking off.\n\n");
552 exit(0);
bdb5c432 553 }
93417736 554 if(dMass < 0.)
bdb5c432 555 {
93417736 556 printf("\n WARNING: dMass < 0 !!!! Please check your settings before taking off.\n\n");
557 exit(0);
bdb5c432 558 }
93417736 559 if(dTemperature <= 1e-44)
560 {
561 printf("\n WARNING: dTemperature <= 0 !!!! Please check your settings before taking off.\n\n");
562 exit(0);
598d292e 563 }
93417736 564 if(TMath::Abs(dV1) > 0.5)
565 {
566 printf("\n WARNING: |dV1| > 0.5 !!!! Please check your settings before taking off.\n\n");
567 exit(0);
568 }
569 if(TMath::Abs(dV2) > 0.5)
570 {
571 printf("\n WARNING: |dV2| > 0.5 !!!! Please check your settings before taking off.\n\n");
572 exit(0);
573 }
574 if(TMath::Abs(dV3) > 0.5)
575 {
576 printf("\n WARNING: |dV3| > 0.5 !!!! Please check your settings before taking off.\n\n");
577 exit(0);
578 }
579 if(TMath::Abs(dV4) > 0.5)
bdb5c432 580 {
93417736 581 printf("\n WARNING: |dV4| > 0.5 !!!! Please check your settings before taking off.\n\n");
bdb5c432 582 exit(0);
583 }
dd1d2f85 584 if(TMath::Abs(dV5) > 0.5)
585 {
586 printf("\n WARNING: |dV5| > 0.5 !!!! Please check your settings before taking off.\n\n");
587 exit(0);
588 }
589 if(TMath::Abs(dV6) > 0.5)
590 {
591 printf("\n WARNING: |dV6| > 0.5 !!!! Please check your settings before taking off.\n\n");
592 exit(0);
593 }
93417736 594 if(LYZ1SUM && LYZ2SUM)
bdb5c432 595 {
93417736 596 cout<<" WARNING: You cannot run LYZ1 and LYZ2 at the same time! LYZ2 needs the output from LYZ1."<<endl;
bdb5c432 597 exit(0);
598 }
93417736 599 if(LYZ1PROD && LYZ2PROD)
bdb5c432 600 {
93417736 601 cout<<" WARNING: You cannot run LYZ1 and LYZ2 at the same time! LYZ2 needs the output from LYZ1."<<endl;
bdb5c432 602 exit(0);
603 }
93417736 604 if(LYZ2SUM && LYZEP)
bdb5c432 605 {
93417736 606 cout<<" WARNING: You cannot run LYZ2 and LYZEP at the same time! LYZEP needs the output from LYZ2."<<endl;
bdb5c432 607 exit(0);
608 }
93417736 609 if(LYZ1SUM && LYZEP)
1f89f68d 610 {
93417736 611 cout<<" WARNING: You cannot run LYZ1 and LYZEP at the same time! LYZEP needs the output from LYZ2."<<endl;
1f89f68d 612 exit(0);
613 }
93417736 614
615 if(!uniformAcceptance && phiMin1 > phiMax1)
616 {
617 cout<<" WARNING: You must have phiMin1 < phiMax1 !!!!"<<endl;
618 exit(0);
619 }
620 if(!uniformAcceptance && !((TMath::Abs(phiMin2) < 1.e-44) && (TMath::Abs(phiMax2) < 1.e-44) && (TMath::Abs(p2) < 1.e-44))
621 && (phiMin2 < phiMax1 || phiMin2 > phiMax2))
622 {
623 cout<<" WARNING: You must have phiMin2 > phiMax1 and phiMin2 < phiMax2 !!!!"<<endl;
624 exit(0);
625 }
626 if((phiMin1 < 0 || phiMin1 > 360) || (phiMax1 < 0 || phiMax1 > 360) ||
627 (phiMin2 < 0 || phiMin2 > 360) || (phiMax2 < 0 || phiMax2 > 360) )
628 {
629 cout<<" WARNING: You must take azimuthal angles from interval [0,360] !!!!"<<endl;
630 exit(0);
631 }
632 if((p1 < 0 || p1 > 1) || (p2 < 0 || p2 > 1))
633 {
634 cout<<" WARNING: you must take p1 and p2 from interval [0,1] !!!!"<<endl;
635 exit(0);
636 }
637 if(bPtDependentV2 && bUniformFluctuationsV2)
638 {
639 cout<<" WARNING: Uniform fluctuations not supported for pt denependent v2 !!!!"<<endl;
640 exit(0);
641 }
642
bdb5c432 643} // end of void CheckUserSettings()
644
93417736 645void WelcomeMessage()
646{
647 // Welcome.
648
649 cout<<endl;
650 cout<<endl;
651 cout<<" ---- ARE YOU READY TO FLY ? ---- "<<endl;
652 cout<<endl;
653
654 gSystem->Sleep(1544);
655
656 cout<<endl;
657 cout<<" ---- BEGIN FLOW ANALYSIS 'ON THE FLY' ---- "<<endl;
658 cout<<endl;
659 cout<<endl;
660
661 gSystem->Sleep(1544);
662
663} // end of void WelcomeMessage()
664
93ff27bd 665void LoadLibraries(const anaModes mode) {
666
667 //--------------------------------------
668 // Load the needed libraries most of them already loaded by aliroot
669 //--------------------------------------
ad87ae62 670 //gSystem->Load("libTree");
5d040cf3 671 gSystem->Load("libGeom");
672 gSystem->Load("libVMC");
673 gSystem->Load("libXMLIO");
674 gSystem->Load("libPhysics");
93ff27bd 675
676 //----------------------------------------------------------
677 // >>>>>>>>>>> Local mode <<<<<<<<<<<<<<
678 //----------------------------------------------------------
679 if (mode==mLocal) {
680 //--------------------------------------------------------
681 // If you want to use already compiled libraries
682 // in the aliroot distribution
683 //--------------------------------------------------------
684 gSystem->Load("libSTEERBase");
685 gSystem->Load("libESD");
686 gSystem->Load("libAOD");
687 gSystem->Load("libANALYSIS");
688 gSystem->Load("libANALYSISalice");
5d040cf3 689 gSystem->Load("libCORRFW");
690 cerr<<"libCORRFW loaded..."<<endl;
2e311896 691 gSystem->Load("libPWGflowBase");
692 cerr<<"libPWGflowBase loaded..."<<endl;
693 gSystem->Load("libPWGflowTasks");
694 cerr<<"libPWGflowTasks loaded..."<<endl;
93ff27bd 695 }
696
697 else if (mode == mLocalPAR) {
698 //--------------------------------------------------------
699 //If you want to use root and par files from aliroot
700 //--------------------------------------------------------
701 //If you want to use root and par files from aliroot
702 //--------------------------------------------------------
703 SetupPar("STEERBase");
704 SetupPar("ESD");
705 SetupPar("AOD");
706 SetupPar("ANALYSIS");
707 SetupPar("ANALYSISalice");
708 SetupPar("PWG2AOD");
709 SetupPar("CORRFW");
2e311896 710 SetupPar("PWGflowBase");
711 cerr<<"PWGflowBase.par loaded..."<<endl;
712 SetupPar("PWGflowTasks");
713 cerr<<"PWGflowTasks.par loaded..."<<endl;
93ff27bd 714 }
715
716 //---------------------------------------------------------
717 // <<<<<<<<<< Source mode >>>>>>>>>>>>
718 //---------------------------------------------------------
719 else if (mode==mLocalSource) {
720
721 // In root inline compile
722
93ff27bd 723 // Constants
2e311896 724 gROOT->LoadMacro("Base/AliFlowCommonConstants.cxx+");
725 gROOT->LoadMacro("Base/AliFlowLYZConstants.cxx+");
93ff27bd 726
727 // Flow event
2e311896 728 gROOT->LoadMacro("Base/AliFlowVector.cxx+");
729 gROOT->LoadMacro("Base/AliFlowTrackSimple.cxx+");
730 gROOT->LoadMacro("Base/AliFlowTrackSimpleCuts.cxx+");
731 gROOT->LoadMacro("Base/AliFlowEventSimple.cxx+");
93ff27bd 732
733 // Output histosgrams
2e311896 734 gROOT->LoadMacro("Base/AliFlowCommonHist.cxx+");
735 gROOT->LoadMacro("Base/AliFlowCommonHistResults.cxx+");
736 gROOT->LoadMacro("Base/AliFlowLYZHist1.cxx+");
737 gROOT->LoadMacro("Base/AliFlowLYZHist2.cxx+");
93ff27bd 738
739 // Functions needed for various methods
2e311896 740 gROOT->LoadMacro("Base/AliCumulantsFunctions.cxx+");
741 gROOT->LoadMacro("Base/AliFlowLYZEventPlane.cxx+");
93ff27bd 742
743 // Flow Analysis code for various methods
2e311896 744 gROOT->LoadMacro("Base/AliFlowAnalysisWithMCEventPlane.cxx+");
745 gROOT->LoadMacro("Base/AliFlowAnalysisWithScalarProduct.cxx+");
746 gROOT->LoadMacro("Base/AliFlowAnalysisWithLYZEventPlane.cxx+");
747 gROOT->LoadMacro("Base/AliFlowAnalysisWithLeeYangZeros.cxx+");
748 gROOT->LoadMacro("Base/AliFlowAnalysisWithCumulants.cxx+");
749 gROOT->LoadMacro("Base/AliFlowAnalysisWithQCumulants.cxx+");
750 gROOT->LoadMacro("Base/AliFlowAnalysisWithFittingQDistribution.cxx+");
751 gROOT->LoadMacro("Base/AliFlowAnalysisWithMixedHarmonics.cxx+");
752 gROOT->LoadMacro("Base/AliFlowAnalysisWithNestedLoops.cxx+");
93ff27bd 753
a47dc7e4 754 // Class to fill the FlowEvent on the fly (generate Monte Carlo events)
2e311896 755 gROOT->LoadMacro("Base/AliFlowEventSimpleMakerOnTheFly.cxx+");
93ff27bd 756
757 cout << "finished loading macros!" << endl;
758
759 }
760
761}
762
763