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