]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWG2/FLOW/macros/runFlowAnalysis.C
start of moving various methods into the flowevent
[u/mrichter/AliRoot.git] / PWG2 / FLOW / macros / runFlowAnalysis.C
CommitLineData
929098e4 1#include <Riostream.h>
2#include <TStopwatch.h>
3#include <TObjArray.h>
4#include <TFile.h>
5#include <TRandom3.h>
6#include <TTimeStamp.h>
7
8#include <AliFlowCommon/AliFlowCommonConstants.h>
9#include <AliFlowCommon/AliFlowLYZConstants.h>
10#include <AliFlowCommon/AliFlowCumuConstants.h>
11#include <AliFlowCommon/AliFlowVector.h>
12#include <AliFlowCommon/AliFlowTrackSimple.h>
13#include <AliFlowCommon/AliFlowEventSimple.h>
14#include <AliFlowCommon/AliFlowTrackSimpleCuts.h>
15#include <AliFlowCommon/AliFlowCommonHist.h>
16#include <AliFlowCommon/AliFlowCommonHistResults.h>
17#include <AliFlowCommon/AliFlowLYZHist1.h>
18#include <AliFlowCommon/AliFlowLYZHist2.h>
19#include <AliFlowCommon/AliCumulantsFunctions.h>
20#include <AliFlowCommon/AliFlowLYZEventPlane.h>
21#include <AliFlowCommon/AliFlowAnalysisWithMCEventPlane.h>
22#include <AliFlowCommon/AliFlowAnalysisWithScalarProduct.h>
23#include <AliFlowCommon/AliFlowAnalysisWithLYZEventPlane.h>
24#include <AliFlowCommon/AliFlowAnalysisWithLeeYangZeros.h>
25#include <AliFlowCommon/AliFlowAnalysisWithCumulants.h>
26#include <AliFlowCommon/AliFlowAnalysisWithQCumulants.h>
27#include <AliFlowCommon/AliFlowAnalysisWithFittingQDistribution.h>
28#include <AliFlowCommon/AliFlowAnalysisWithMixedHarmonics.h>
29#include <AliFlowCommon/AliFlowAnalysisWithNestedLoops.h>
30#include <AliFlowTasks/AliFlowEventSimpleMaker.h>
31#include <FlowEventMakers/FlowEventSimpleMaker.h>
334e3256 32
03a02aca 33//--------------------------------------------------------------------------------------
c852719e 34// Run flow analysis on local data with custom FlowEvent maker
03a02aca 35// RUN SETTINGS
334e3256 36//flow analysis method can be: (set to kTRUE or kFALSE)
c741f5d0 37Bool_t SP = kTRUE;
5b40431d 38Bool_t LYZ1SUM = kTRUE;
39Bool_t LYZ1PROD = kTRUE;
c741f5d0 40Bool_t LYZ2SUM = kFALSE;
41Bool_t LYZ2PROD = kFALSE;
5b40431d 42Bool_t LYZEP = kFALSE;
c741f5d0 43Bool_t GFC = kTRUE;
44Bool_t QC = kTRUE;
b5425e51 45Bool_t FQD = kTRUE;
ecac11c2 46Bool_t MH = kTRUE;
47Bool_t NL = kFALSE;
5b40431d 48Bool_t MCEP = kFALSE; //does not work yet 24/12/08
03a02aca 49//--------------------------------------------------------------------------------------
334e3256 50
20e52807 51// Weights
52// Use weights for Q vector
53Bool_t usePhiWeights = kFALSE; //Phi (correction for non-uniform azimuthal acceptance)
54Bool_t usePtWeights = kFALSE; //v'(pt) (differential flow in pt)
55Bool_t useEtaWeights = kFALSE; //v'(eta) (differential flow in eta)
56
03a02aca 57//--------------------------------------------------------------------------------------
58// CUT SETTINGS
334e3256 59//integrated selection
60Double_t ptMaxInt = 10.;
61Double_t ptMinInt = 0.;
62Double_t etaMaxInt = 1.;
63Double_t etaMinInt = -1.;
64Double_t phiMaxInt = 7.5;
65Double_t phiMinInt = 0.;
66Int_t PIDInt = 211;
67
68//differential selection
69Double_t ptMaxDiff = 10.;
70Double_t ptMinDiff = 0.;
71Double_t etaMaxDiff = 1.;
72Double_t etaMinDiff = -1.;
73Double_t phiMaxDiff = 7.5;
74Double_t phiMinDiff = 0.;
75Int_t PIDDiff = 211;
c741f5d0 76/*
03a02aca 77//--------------------------------------------------------------------------------------
33314629 78// FLOW SETTINGS (R.Rietkerk)
79Int_t nLoops=1; // Number of times to use the same particle (nonflow).
c741f5d0 80Double_t xEllipticFlowValue=0.1;// Add Elliptic Flow. Must be in range [0,1].
81Int_t nMultiplicityOfEvent=500; // Set Average Multiplicity.
33314629 82Double_t xSigmaFlow=0.00; // Add Elliptic Flow. Must be in range [0,1].
c741f5d0 83Int_t nSigmaMult=50; // Set Average Multiplicity.
33314629 84//--------------------------------------------------------------------------------------
c741f5d0 85*/
26c4cbb9 86
fdd8c18f 87enum anaModes {mLocal,mLocalSource,mLocalPAR,};
88//mLocal: Analyze data on your computer using aliroot
89//mLocalPAR: Analyze data on your computer using root + PAR files
90//mLocalSource: Analyze data on your computer using root + source files
91
929098e4 92void LoadLibraries(const anaModes mode);
93
26c4cbb9 94Int_t offset = 0;
95
929098e4 96int runFlowAnalysis(const anaModes mode=mLocal, Int_t aRuns = 100, const char*
c741f5d0 97 dir="/data/alice1/kolk/KineOnly3/")
98 // dir="/Users/snelling/alice_data/KineOnly3/")
33314629 99 // dir="/Users/snelling/alice_data/stoomboot/5b/")
334e3256 100{
101 TStopwatch timer;
102 timer.Start();
103
929098e4 104 if (LYZ1SUM && LYZ2SUM) {cout<<"WARNING: you cannot run LYZ1 and LYZ2 at the same time! LYZ2 needs the output from LYZ1."<<endl; exit(1); }
105 if (LYZ1PROD && LYZ2PROD) {cout<<"WARNING: you cannot run LYZ1 and LYZ2 at the same time! LYZ2 needs the output from LYZ1."<<endl; exit(1); }
106 if (LYZ2SUM && LYZEP) {cout<<"WARNING: you cannot run LYZ2 and LYZEP at the same time! LYZEP needs the output from LYZ2."<<endl; exit(1); }
107 if (LYZ1SUM && LYZEP) {cout<<"WARNING: you cannot run LYZ1 and LYZEP at the same time! LYZEP needs the output from LYZ2."<<endl; exit(1); }
334e3256 108
109
110 cout<<endl;
111 cout<<" ---- BEGIN ANALYSIS ---- "<<endl;
112 cout<<endl;
113
fdd8c18f 114 LoadLibraries(mode);
2a413443 115
33314629 116 TRandom3 random3Temp; //init for manual settings (R.Rietkerk)
117 TTimeStamp dt;
118 Int_t sseed = dt.GetNanoSec()/1000;
119 random3Temp.SetSeed(sseed);
120
fdd8c18f 121 if (mode == mLocal || mode == mLocalPAR) {
122 // AliFlow event in aliroot or with pars
123 AliFlowEventSimpleMaker* fEventMaker = new AliFlowEventSimpleMaker();
124 }
125 else if (mode == mLocalSource) {
126 // flow event in source mode
127 FlowEventSimpleMaker* fEventMaker = new FlowEventSimpleMaker();
128 }
129 else{
130 cout << "No supported running mode selected!" << endl;
131 break;
132 }
2a413443 133
2a413443 134
334e3256 135 //------------------------------------------------------------------------
e085f1a9 136 //cuts:
334e3256 137 AliFlowTrackSimpleCuts* cutsInt = new AliFlowTrackSimpleCuts();
138 cutsInt->SetPtMax(ptMaxInt);
139 cutsInt->SetPtMin(ptMinInt);
140 cutsInt->SetEtaMax(etaMaxInt);
141 cutsInt->SetEtaMin(etaMinInt);
142 cutsInt->SetPhiMax(phiMaxInt);
143 cutsInt->SetPhiMin(phiMinInt);
144 cutsInt->SetPID(PIDInt);
145
146 AliFlowTrackSimpleCuts* cutsDiff = new AliFlowTrackSimpleCuts();
147 cutsDiff->SetPtMax(ptMaxDiff);
148 cutsDiff->SetPtMin(ptMinDiff);
149 cutsDiff->SetEtaMax(etaMaxDiff);
150 cutsDiff->SetEtaMin(etaMinDiff);
151 cutsDiff->SetPhiMax(phiMaxDiff);
152 cutsDiff->SetPhiMin(phiMinDiff);
153 cutsDiff->SetPID(PIDDiff);
154
03a02aca 155 //if the weights are used:
156 TFile *fileWithWeights = NULL;
157 TList *listWithWeights = NULL;
334e3256 158
14db9c04 159 if(usePhiWeights||usePtWeights||useEtaWeights) {
160 fileWithWeights = TFile::Open("weights.root","READ");
161 if(fileWithWeights) {
162 listWithWeights = (TList*)fileWithWeights->Get("weights");
163 }
164 else
165 {cout << " WARNING: the file <weights.root> with weights from the previous run was not found."<<endl;
166 break;
167 }
e085f1a9 168 }
169
170 //flow methods:
ce4a88f5 171 AliFlowAnalysisWithQCumulants *qc = NULL;
172 AliFlowAnalysisWithCumulants *gfc = NULL;
173 AliFlowAnalysisWithFittingQDistribution *fqd = NULL;
174 AliFlowAnalysisWithLeeYangZeros *lyz1sum = NULL;
175 AliFlowAnalysisWithLeeYangZeros *lyz1prod = NULL;
176 AliFlowAnalysisWithLeeYangZeros *lyz2sum = NULL;
177 AliFlowAnalysisWithLeeYangZeros *lyz2prod = NULL;
178 AliFlowAnalysisWithLYZEventPlane *lyzep = NULL;
179 AliFlowAnalysisWithScalarProduct *sp = NULL;
ecac11c2 180 AliFlowAnalysisWithMCEventPlane *mcep = NULL;
181 AliFlowAnalysisWithMixedHarmonics *mh = NULL;
182 AliFlowAnalysisWithNestedLoops *nl = NULL;
deebd72f 183
334e3256 184 //MCEP = monte carlo event plane
deebd72f 185 if (MCEP) {
186 AliFlowAnalysisWithMCEventPlane *mcep = new AliFlowAnalysisWithMCEventPlane();
187 mcep->Init();
188 }
334e3256 189
190 //QC = Q-cumulants
deebd72f 191 if(QC) {
192 AliFlowAnalysisWithQCumulants* qc = new AliFlowAnalysisWithQCumulants();
03a02aca 193 if(listWithWeights) qc->SetWeightsList(listWithWeights);
194 if(usePhiWeights) qc->SetUsePhiWeights(usePhiWeights);
195 if(usePtWeights) qc->SetUsePtWeights(usePtWeights);
196 if(useEtaWeights) qc->SetUseEtaWeights(useEtaWeights);
ecac11c2 197 qc->Init();
deebd72f 198 }
334e3256 199
200 //GFC = Generating Function Cumulants
deebd72f 201 if(GFC) {
202 AliFlowAnalysisWithCumulants* gfc = new AliFlowAnalysisWithCumulants();
e5e75b58 203 if(listWithWeights) gfc->SetWeightsList(listWithWeights);
204 if(usePhiWeights) gfc->SetUsePhiWeights(usePhiWeights);
205 if(usePtWeights) gfc->SetUsePtWeights(usePtWeights);
206 if(useEtaWeights) gfc->SetUseEtaWeights(useEtaWeights);
ecac11c2 207 gfc->Init();
deebd72f 208 }
334e3256 209
210 //FQD = Fitting q-distribution
deebd72f 211 if(FQD) {
ce4a88f5 212 AliFlowAnalysisWithFittingQDistribution* fqd = new AliFlowAnalysisWithFittingQDistribution();
14db9c04 213 if(listWithWeights) fqd->SetWeightsList(listWithWeights);
214 if(usePhiWeights) fqd->SetUsePhiWeights(usePhiWeights);
ecac11c2 215 fqd->Init();
deebd72f 216 }
60ec66fb 217
218 //SP = Scalar Product
219 if(SP) {
220 AliFlowAnalysisWithScalarProduct* sp = new AliFlowAnalysisWithScalarProduct();
ecac11c2 221 if(usePhiWeights) sp->SetUsePhiWeights(usePhiWeights);
60ec66fb 222 sp->Init();
223 }
224
334e3256 225 //LYZ1 = Lee-Yang Zeroes first run
c741f5d0 226 if(LYZ1SUM) {
227 AliFlowAnalysisWithLeeYangZeros* lyz1sum = new AliFlowAnalysisWithLeeYangZeros();
228 lyz1sum->SetFirstRun(kTRUE);
229 lyz1sum->SetUseSum(kTRUE);
230 lyz1sum->Init();
231 }
232 if(LYZ1PROD) {
233 AliFlowAnalysisWithLeeYangZeros* lyz1prod = new AliFlowAnalysisWithLeeYangZeros();
234 lyz1prod->SetFirstRun(kTRUE);
235 lyz1prod->SetUseSum(kFALSE);
236 lyz1prod->Init();
deebd72f 237 }
334e3256 238 //LYZ2 = Lee-Yang Zeroes second run
c741f5d0 239 if(LYZ2SUM) {
240 AliFlowAnalysisWithLeeYangZeros* lyz2sum = new AliFlowAnalysisWithLeeYangZeros();
deebd72f 241 // read the input file from the first run
c741f5d0 242 TString inputFileNameLYZ2SUM = "outputLYZ1SUManalysis.root" ;
243 TFile* inputFileLYZ2SUM = new TFile(inputFileNameLYZ2SUM.Data(),"READ");
244 if(!inputFileLYZ2SUM || inputFileLYZ2SUM->IsZombie()) {
5b40431d 245 cerr << " ERROR: To run LYZ2SUM you need the output file from LYZ1SUM. This file is not there! Please run LYZ1SUM first." << endl ;
ce78b98c 246 break;
deebd72f 247 }
248 else {
c741f5d0 249 TList* inputListLYZ2SUM = (TList*)inputFileLYZ2SUM->Get("cobjLYZ1SUM");
250 if (!inputListLYZ2SUM) {cout<<"SUM Input list is NULL pointer!"<<endl; break;}
deebd72f 251 else {
c741f5d0 252 cout<<"LYZ2SUM input file/list read..."<<endl;
253 lyz2sum->SetFirstRunList(inputListLYZ2SUM);
254 lyz2sum->SetFirstRun(kFALSE);
255 lyz2sum->SetUseSum(kTRUE);
256 lyz2sum->Init();
257 }
258 }
259 }
260 if(LYZ2PROD) {
261 AliFlowAnalysisWithLeeYangZeros* lyz2prod = new AliFlowAnalysisWithLeeYangZeros();
262 // read the input file from the first run
263 TString inputFileNameLYZ2PROD = "outputLYZ1PRODanalysis.root" ;
264 TFile* inputFileLYZ2PROD = new TFile(inputFileNameLYZ2PROD.Data(),"READ");
265 if(!inputFileLYZ2PROD || inputFileLYZ2PROD->IsZombie()) {
5b40431d 266 cerr << " ERROR: To run LYZ2PROD you need the output file from LYZ1PROD. This file is not there! Please run LYZ1PROD first." << endl ;
c741f5d0 267 break;
268 }
269 else {
270 TList* inputListLYZ2PROD = (TList*)inputFileLYZ2PROD->Get("cobjLYZ1PROD");
271 if (!inputListLYZ2PROD) {cout<<"PROD Input list is NULL pointer!"<<endl; break;}
272 else {
273 cout<<"LYZ2PROD input file/list read..."<<endl;
274 lyz2prod->SetFirstRunList(inputListLYZ2PROD);
275 lyz2prod->SetFirstRun(kFALSE);
276 lyz2prod->SetUseSum(kTRUE);
277 lyz2prod->Init();
334e3256 278 }
279 }
deebd72f 280 }
334e3256 281 //LYZEP = Lee-Yang Zeroes event plane
deebd72f 282 if(LYZEP) {
283 AliFlowLYZEventPlane* ep = new AliFlowLYZEventPlane() ;
284 AliFlowAnalysisWithLYZEventPlane* lyzep = new AliFlowAnalysisWithLYZEventPlane();
285 // read the input file from the second lyz run
c741f5d0 286 TString inputFileNameLYZEP = "outputLYZ2SUManalysis.root" ;
deebd72f 287 TFile* inputFileLYZEP = new TFile(inputFileNameLYZEP.Data(),"READ");
288 if(!inputFileLYZEP || inputFileLYZEP->IsZombie()) {
5b40431d 289 cerr << " ERROR: To run LYZEP you need the output file from LYZ2SUM. This file is not there! Please run LYZ2SUM first." << endl ;
ce78b98c 290 break;
291 }
deebd72f 292 else {
c741f5d0 293 TList* inputListLYZEP = (TList*)inputFileLYZEP->Get("cobjLYZ2SUM");
ace38bad 294 if (!inputListLYZEP) {cout<<"Input list is NULL pointer!"<<endl; break;}
deebd72f 295 else {
296 cout<<"LYZEP input file/list read..."<<endl;
297 ep ->SetSecondRunList(inputListLYZEP);
298 lyzep->SetSecondRunList(inputListLYZEP);
299 ep ->Init();
300 lyzep->Init();
334e3256 301 }
302 }
deebd72f 303 }
ecac11c2 304 // MH = Mixed Harmonics:
305 if(MH) {
306 AliFlowAnalysisWithMixedHarmonics* mh = new AliFlowAnalysisWithMixedHarmonics();
307 if(listWithWeights) mh->SetWeightsList(listWithWeights);
308 //if(usePhiWeights) mh->SetUsePhiWeights(usePhiWeights); // to be improved (enabled)
309 //if(usePtWeights) mh->SetUsePtWeights(usePtWeights); // to be improved (enabled)
310 //if(useEtaWeights) mh->SetUseEtaWeights(useEtaWeights); // to be improved (enabled)
311 mh->Init();
312 }
313 // NL = Nested Loops:
314 if(NL) {
315 AliFlowAnalysisWithNestedLoops* nl = new AliFlowAnalysisWithNestedLoops();
316 if(listWithWeights) nl->SetWeightsList(listWithWeights);
317 //if(usePhiWeights) nl->SetUsePhiWeights(usePhiWeights); // to be improved (enabled)
318 //if(usePtWeights) nl->SetUsePtWeights(usePtWeights); // to be improved (enabled)
319 //if(useEtaWeights) nl->SetUseEtaWeights(useEtaWeights); // to be improved (enabled)
320 nl->Init();
321 }
322
334e3256 323 //------------------------------------------------------------------------
324
26c4cbb9 325
b25cc698 326 //standard code to read files in directory
334e3256 327 Int_t fCount = 0;
328 TString execDir(gSystem->pwd());
329 TString targetDir(dir);
330 TSystemDirectory* baseDir = new TSystemDirectory(".", dir);
331 TList* dirList = baseDir->GetListOfFiles();
332 if (!dirList) {
333 cout << endl << "No input files in: " << targetDir.Data() << endl;
334 break;
335 }
336 Int_t nDirs = dirList->GetEntries();
337 cout<<endl;
338 cout<<"Int_t nDirs = "<<nDirs<<endl;
339 gSystem->cd(execDir);
340
341 for(Int_t iDir=0;iDir<nDirs;++iDir)
342 {
343 TSystemFile* presentDir = (TSystemFile*)dirList->At(iDir);
344 if(!presentDir || !presentDir->IsDirectory() ||
345 strcmp(presentDir->GetName(), ".") == 0 ||
346 strcmp(presentDir->GetName(), "..") == 0)
347 {
348 cout << endl;
349 cout << "Directory (" << iDir << "): " << presentDir->GetName() <<
350 " - Skipping ... " << endl;
351 continue ;
352 }
353
354 if(offset > 0) { --offset ; continue ; }
355 if((aRuns > 0) && (fCount >= aRuns)) { break ; }
356
357 TString presentDirName(dir); // aDataDir
358 presentDirName += presentDir->GetName();
359 presentDirName += "/";
360 //cerr<<" presentDirName = "<<presentDirName<<endl;
361
362 TString fileName = presentDirName;
363 fileName += "galice.root";
364 Long_t *id, *size, *flags, *modtime;
365 if(gSystem->GetPathInfo(fileName.Data(),id,size,flags,modtime))
366 {
367 cout << " File : " << fileName << " does NOT exist ! - Skipping ... "
368 << endl;
369 continue;
370 }
deebd72f 371 // cout << endl ; cout << "Directory (" << iDir << "): " << presentDirName << " ... " << endl;
334e3256 372
373 //loop (simulations in the present dir)
374 TSystemDirectory* evtsDir = new TSystemDirectory(".",
375 presentDirName.Data());
376 TList* fileList = evtsDir->GetListOfFiles();
377 Int_t nFiles = fileList->GetEntries();
deebd72f 378 //cout<<" Int_t nFiles = "<<nFiles<<endl;
334e3256 379 gSystem->cd(execDir);
380 for(Int_t iFiles=0; iFiles<nFiles; ++iFiles)
381 {
382 TSystemFile* presentFile = (TSystemFile*) fileList->At(iFiles);
383 TString presentFileName(presentDirName);
384 presentFileName += presentFile->GetName();
385
386 if(!(presentFileName.Contains("Kinematics") &&
387 presentFileName.Contains("root"))) { continue ; }
388
deebd72f 389 //cout << " found: " << presentFileName.Data() << endl;
334e3256 390
391 TFile* kineFile = new TFile(presentFileName.Data(), "READ");
392 // kineFile->ls();
393 Int_t nEvts = kineFile->GetNkeys() ;
deebd72f 394 //cout << " . found: " << nEvts << " KineTree(s) in " << presentFileName.Data() << endl;
334e3256 395 TList* kineEventsList = (TList*)kineFile->GetListOfKeys();
396 TTree* kTree;
397 TIter next(kineEventsList);
398 TKey* key;
33314629 399
400 Double_t xRPAngle;
334e3256 401 //loop over the events
402 while( key=(TKey *)next() )
403 {
404 TDirectory* tDir = (TDirectory*)key->ReadObj();
405 if(!tDir) break;
406
407 TString evtDir(tDir->GetName());
deebd72f 408 //cout << " . . found: " << tDir->GetName() << endl;
334e3256 409
410 kTree = (TTree *)tDir->Get("TreeK");
411 if(!kTree) break;
412
413 Int_t nPart = kTree->GetEntries();
deebd72f 414 //cout << " . . . kTree " << fCount << " has " << nPart << " particles " << endl;
334e3256 415
416 //-----------------------------------------------------------
26c4cbb9 417 //fill and save the flow event
33314629 418
c741f5d0 419 /*
33314629 420 Int_t nNewMultOfEvent = random3Temp.Gaus(nMultiplicityOfEvent,nSigmaMult);
421 cout << "new multiplicity: " << nNewMultOfEvent << endl;
422 Double_t xNewFlowValue = random3Temp.Gaus(xEllipticFlowValue,xSigmaFlow);
4f678629 423 if ( (fCount % 100) == 0) {
424 cout << "new multiplicity: " << nNewMultOfEvent << endl;
425 cout << "new flow value: " << xNewFlowValue << endl;
426 }
33314629 427 fEventMaker->SetNoOfLoops(nLoops);
428 fEventMaker->SetEllipticFlowValue(xNewFlowValue);
429 fEventMaker->SetMultiplicityOfEvent(nNewMultOfEvent);
430 xRPAngle=TMath::TwoPi()*random3Temp.Rndm();
431 fEventMaker->SetMCReactionPlaneAngle(xRPAngle);
c741f5d0 432 */
33314629 433
26c4cbb9 434 AliFlowEventSimple *fEvent = fEventMaker->FillTracks(kTree, cutsInt, cutsDiff);
435
b25cc698 436 // do flow analysis for various methods
c741f5d0 437 if(MCEP) mcep->Make(fEvent);
438 if(QC) qc->Make(fEvent);
439 if(GFC) gfc->Make(fEvent);
440 if(FQD) fqd->Make(fEvent);
441 if(LYZ1SUM) lyz1sum->Make(fEvent);
442 if(LYZ1PROD)lyz1prod->Make(fEvent);
443 if(LYZ2SUM) lyz2sum->Make(fEvent);
444 if(LYZ2PROD)lyz2prod->Make(fEvent);
445 if(LYZEP) lyzep->Make(fEvent,ep);
ecac11c2 446 if(SP) sp->Make(fEvent);
447 if(MH) mh->Make(fEvent);
448 if(NL) nl->Make(fEvent);
334e3256 449 //-----------------------------------------------------------
334e3256 450 fCount++;
951817d5 451 //cout << "# " << fCount << " events processed" << endl;
334e3256 452 delete kTree;
453 delete fEvent;
454 }
455 delete kineFile ;
456 }
457 delete evtsDir ;
458 }
459
b5425e51 460 //---------------------------------------------------------------------------------------
461 // create a new file which will hold the final results of all methods:
462 TString outputFileName = "AnalysisResults.root";
463 TFile *outputFile = new TFile(outputFileName.Data(),"RECREATE");
464 // create a new file for each method wich will hold list with final results:
ecac11c2 465 const Int_t nMethods = 12;
466 TString method[nMethods] = {"MCEP","SP","GFC","QC","FQD","LYZ1SUM","LYZ1PROD","LYZ2SUM","LYZ2PROD","LYZEP","MH","NL"};
b5425e51 467 TDirectoryFile *dirFileFinal[nMethods] = {NULL};
468 TString fileNameMethod[nMethods];
469 for(Int_t i=0;i<nMethods;i++)
470 {
471 // form a file name for each method:
472 fileNameMethod[i]+="output";
473 fileNameMethod[i]+=method[i].Data();
474 fileNameMethod[i]+="analysis";
475 dirFileFinal[i] = new TDirectoryFile(fileNameMethod[i].Data(),fileNameMethod[i].Data());
476 }
477
478 // calculating and storing the final results of default flow analysis:
479 if(MCEP) {mcep->Finish(); mcep->WriteHistograms(dirFileFinal[0]);}
480 if(SP) {sp->Finish(); sp->WriteHistograms(dirFileFinal[1]);}
481 if(GFC) {gfc->Finish(); gfc->WriteHistograms(dirFileFinal[2]);}
482 if(QC) {qc->Finish(); qc->WriteHistograms(dirFileFinal[3]);}
483 if(FQD) {fqd->Finish(); fqd->WriteHistograms(dirFileFinal[4]);}
484 if(LYZ1SUM) {lyz1sum->Finish(); lyz1sum->WriteHistograms(dirFileFinal[5]);}
485 if(LYZ1PROD){lyz1prod->Finish();lyz1prod->WriteHistograms(dirFileFinal[6]);}
486 if(LYZ2SUM) {lyz2sum->Finish(); lyz2sum->WriteHistograms(dirFileFinal[7]);}
487 if(LYZ2PROD){lyz2prod->Finish();lyz2prod->WriteHistograms(dirFileFinal[8]);}
488 if(LYZEP) {lyzep->Finish(); lyzep->WriteHistograms(dirFileFinal[9]);}
ecac11c2 489 if(MH) {mh->Finish(); mh->WriteHistograms(dirFileFinal[10]);}
490 if(NL) {nl->Finish(); nl->WriteHistograms(dirFileFinal[11]);}
b5425e51 491 //---------------------------------------------------------------------------------------
492
493 outputFile->Close();
494 delete outputFile;
334e3256 495
496 cout << endl;
ace38bad 497 cout << " Fini ... " << endl;
334e3256 498 cout << endl;
499
500 timer.Stop();
501 cout << endl;
502 timer.Print();
503}
21e694dd 504
505void SetupPar(char* pararchivename)
506{
507 //Load par files, create analysis libraries
508 //For testing, if par file already decompressed and modified
509 //classes then do not decompress.
510
511 TString cdir(Form("%s", gSystem->WorkingDirectory() )) ;
512 TString parpar(Form("%s.par", pararchivename)) ;
513 if ( gSystem->AccessPathName(parpar.Data()) ) {
514 gSystem->ChangeDirectory(gSystem->Getenv("ALICE_ROOT")) ;
515 TString processline(Form(".! make %s", parpar.Data())) ;
516 gROOT->ProcessLine(processline.Data()) ;
517 gSystem->ChangeDirectory(cdir) ;
518 processline = Form(".! mv /tmp/%s .", parpar.Data()) ;
519 gROOT->ProcessLine(processline.Data()) ;
520 }
521 if ( gSystem->AccessPathName(pararchivename) ) {
522 TString processline = Form(".! tar xvzf %s",parpar.Data()) ;
523 gROOT->ProcessLine(processline.Data());
524 }
525
526 TString ocwd = gSystem->WorkingDirectory();
527 gSystem->ChangeDirectory(pararchivename);
528
529 // check for BUILD.sh and execute
530 if (!gSystem->AccessPathName("PROOF-INF/BUILD.sh")) {
531 printf("*******************************\n");
532 printf("*** Building PAR archive ***\n");
533 cout<<pararchivename<<endl;
534 printf("*******************************\n");
535
536 if (gSystem->Exec("PROOF-INF/BUILD.sh")) {
537 Error("runProcess","Cannot Build the PAR Archive! - Abort!");
538 return -1;
539 }
540 }
541 // check for SETUP.C and execute
542 if (!gSystem->AccessPathName("PROOF-INF/SETUP.C")) {
543 printf("*******************************\n");
544 printf("*** Setup PAR archive ***\n");
545 cout<<pararchivename<<endl;
546 printf("*******************************\n");
547 gROOT->Macro("PROOF-INF/SETUP.C");
548 }
549
550 gSystem->ChangeDirectory(ocwd.Data());
551 printf("Current dir: %s\n", ocwd.Data());
552}
553
fdd8c18f 554void LoadLibraries(const anaModes mode) {
555
556 //--------------------------------------
557 // Load the needed libraries most of them already loaded by aliroot
558 //--------------------------------------
b5425e51 559 //gSystem->Load("libTree");
5d040cf3 560 gSystem->Load("libGeom");
561 gSystem->Load("libVMC");
562 gSystem->Load("libXMLIO");
563 gSystem->Load("libPhysics");
fdd8c18f 564
565 //----------------------------------------------------------
566 // >>>>>>>>>>> Local mode <<<<<<<<<<<<<<
567 //----------------------------------------------------------
568 if (mode==mLocal) {
569 //--------------------------------------------------------
570 // If you want to use already compiled libraries
571 // in the aliroot distribution
572 //--------------------------------------------------------
573 gSystem->Load("libSTEERBase");
574 gSystem->Load("libESD");
575 gSystem->Load("libAOD");
576 gSystem->Load("libANALYSIS");
577 gSystem->Load("libANALYSISalice");
5d040cf3 578 gSystem->Load("libCORRFW");
579 cerr<<"libCORRFW loaded..."<<endl;
580 gSystem->Load("libPWG2flowCommon");
581 cerr<<"libPWG2flowCommon loaded..."<<endl;
582 gSystem->Load("libPWG2flowTasks");
583 cerr<<"libPWG2flowTasks loaded..."<<endl;
fdd8c18f 584 }
585
586 else if (mode == mLocalPAR) {
587 //--------------------------------------------------------
588 //If you want to use root and par files from aliroot
589 //--------------------------------------------------------
590 //If you want to use root and par files from aliroot
591 //--------------------------------------------------------
592 SetupPar("STEERBase");
593 SetupPar("ESD");
594 SetupPar("AOD");
595 SetupPar("ANALYSIS");
596 SetupPar("ANALYSISalice");
597 SetupPar("PWG2AOD");
598 SetupPar("CORRFW");
599 SetupPar("PWG2flowCommon");
600 cerr<<"PWG2flowCommon.par loaded..."<<endl;
601 SetupPar("PWG2flowTasks");
602 cerr<<"PWG2flowTasks.par loaded..."<<endl;
603 }
604
605 //---------------------------------------------------------
606 // <<<<<<<<<< Source mode >>>>>>>>>>>>
607 //---------------------------------------------------------
608 else if (mode==mLocalSource) {
609
610 // In root inline compile
611
612
613 // Constants
614 gROOT->LoadMacro("AliFlowCommon/AliFlowCommonConstants.cxx+");
615 gROOT->LoadMacro("AliFlowCommon/AliFlowLYZConstants.cxx+");
616 gROOT->LoadMacro("AliFlowCommon/AliFlowCumuConstants.cxx+");
617
618 // Flow event
619 gROOT->LoadMacro("AliFlowCommon/AliFlowVector.cxx+");
620 gROOT->LoadMacro("AliFlowCommon/AliFlowTrackSimple.cxx+");
621 gROOT->LoadMacro("AliFlowCommon/AliFlowEventSimple.cxx+");
622
623 // Cuts
624 gROOT->LoadMacro("AliFlowCommon/AliFlowTrackSimpleCuts.cxx+");
625
626 // Output histosgrams
627 gROOT->LoadMacro("AliFlowCommon/AliFlowCommonHist.cxx+");
628 gROOT->LoadMacro("AliFlowCommon/AliFlowCommonHistResults.cxx+");
629 gROOT->LoadMacro("AliFlowCommon/AliFlowLYZHist1.cxx+");
630 gROOT->LoadMacro("AliFlowCommon/AliFlowLYZHist2.cxx+");
631
632 // Functions needed for various methods
633 gROOT->LoadMacro("AliFlowCommon/AliCumulantsFunctions.cxx+");
fdd8c18f 634 gROOT->LoadMacro("AliFlowCommon/AliFlowLYZEventPlane.cxx+");
635
636 // Flow Analysis code for various methods
637 gROOT->LoadMacro("AliFlowCommon/AliFlowAnalysisWithMCEventPlane.cxx+");
638 gROOT->LoadMacro("AliFlowCommon/AliFlowAnalysisWithScalarProduct.cxx+");
639 gROOT->LoadMacro("AliFlowCommon/AliFlowAnalysisWithLYZEventPlane.cxx+");
640 gROOT->LoadMacro("AliFlowCommon/AliFlowAnalysisWithLeeYangZeros.cxx+");
641 gROOT->LoadMacro("AliFlowCommon/AliFlowAnalysisWithCumulants.cxx+");
642 gROOT->LoadMacro("AliFlowCommon/AliFlowAnalysisWithQCumulants.cxx+");
ecac11c2 643 gROOT->LoadMacro("AliFlowCommon/AliFlowAnalysisWithFittingQDistribution.cxx+");
644 gROOT->LoadMacro("AliFlowCommon/AliFlowAnalysisWithMixedHarmonics.cxx+");
645 gROOT->LoadMacro("AliFlowCommon/AliFlowAnalysisWithNestedLoops.cxx+");
fdd8c18f 646
647 // Class to fill the FlowEvent without aliroot dependence
648 // can be found in the directory FlowEventMakers
649 gROOT->LoadMacro("FlowEventMakers/FlowEventSimpleMaker.cxx+");
650
651 cout << "finished loading macros!" << endl;
652
653 }
654
655}
656
657