]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/FLOW/macros/runFlowAnalysisOnTheFly.C
fix product method with <M> and better error messages macros
[u/mrichter/AliRoot.git] / PWG2 / FLOW / macros / runFlowAnalysisOnTheFly.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;
19 //--------------------------------------------------------------------------------------
20
21 // Weights 
22 // use weights for Q-vector:
23 Bool_t usePhiWeights = kFALSE; // phi weights (correction for non-uniform azimuthal acceptance)
24 Bool_t usePtWeights  = kFALSE; // pt weights 
25 Bool_t useEtaWeights = kFALSE; // eta weights
26
27 // Parameters for the simulation of events 'on the fly': 
28 Bool_t bSameSeed = kTRUE; // use always the same seed for random generators. 
29                            // usage of same seed (kTRUE) is relevant in two cases:
30                            // 1.) If you want to use LYZ method to calcualte differential flow;
31                            // 2.) If you want to use phi weights for GFC, QC and FQD
32                            
33 Bool_t bConstantHarmonics = kTRUE; // harmonics V1, V2, V4... are constant (kTRUE) or functions of pt and eta (kFALSE)
34
35 Int_t iLoops = 1; // number of times to use each track (to simulate nonflow)
36
37 Bool_t bMultDistrOfRPsIsGauss = kFALSE; // 1.) if kTRUE  = multiplicitiy of RPs is sampled e-b-e from Gaussian distribution with
38                                         //                 mean = iMultiplicityOfRP and spread = dMultiplicitySpreadOfRP
39                                         // 2.) if kFALSE = multiplicitiy of RPs is sampled e-b-e uniformly from 
40                                         //                 interval [iMinMultOfRP,iMaxMultOfRP]
41                                         // 3.) for a fixed multiplicity use Gaussian with zero spread or use uniform with iMinMult=iMaxMult
42                                         
43 Bool_t bV2DistrOfRPsIsGauss = kFALSE; // 1.) if kTRUE  = elliptic flow of RPs is sampled e-b-e from Gaussian distribution with
44                                       //                 mean = dV2RP and spread = dV2SpreadRP
45                                       // 2.) if kFALSE = elliptic flow of RPs is sampled e-b-e uniformly from 
46                                       //                 interval [dMinV2RP,dMaxV2RP]
47                                       // 3.) for a fixed elliptic flow use Gaussian with zero spread or use uniform with dMinV2RP=dMaxV2RP
48                                                                                    
49 Int_t iMultiplicityOfRP = 500;        // mean multiplicity of RPs (if sampled from Gaussian)
50 Double_t dMultiplicitySpreadOfRP = 0; // multiplicity spread of RPs (if sampled from Gaussian)
51 Int_t iMinMultOfRP = 400;             // minimal multiplicity of RPs (if sampled uniformly)
52 Int_t iMaxMultOfRP = 600;             // maximal multiplicity of RPs (if sampled uniformly)
53
54 Double_t dTemperatureOfRP = 0.44; // 'temperature' of RPs in GeV/c (increase this parameter to get more high pt RPs) 
55
56 //......................................................................................  
57 // if you use (pt,eta) dependent harmonics (bConstantHarmonics = kFALSE):
58 Double_t dPtCutOff = 2.0; // V2(pt) is linear up to pt = 2 GeV and for pt > 2 GeV it is constant: V2(pt) = dVRPMax
59 Double_t dV2RPMax = 0.20; // maximum value of V2(pt) for pt >= 2GeV
60 //...................................................................................... 
61
62 //......................................................................................  
63 // if you use constant harmonics (bConstantHarmonics = kTRUE) (i.e. no pt dependence):
64 Double_t dV2RP = 0.05;       // elliptic flow of RPs (if sampled from Gaussian)
65 Double_t dV2SpreadRP = 0.0;  // elliptic flow spread of RPs (if sampled from Gaussian)
66 Double_t dMinV2RP = 0.04;    // minimal elliptic flow of RPs (if sampled uniformly)
67 Double_t dMaxV2RP = 0.06;    // maximal elliptic flow of RPs (if sampled uniformly)
68
69 Double_t dV1RP = 0.0; // directed flow of RPs
70 Double_t dV1SpreadRP = 0.0; // directed flow spread of RPs
71
72 Double_t dV4RP = 0.0; // harmonic V4 of RPs (to be improved: name needed)
73 Double_t dV4SpreadRP = 0.0; // harmonic V4's spread of RPs (to be improved: name needed)
74 //......................................................................................  
75
76 enum anaModes {mLocal,mLocalSource,mLocalPAR};
77 // mLocal: Analyze data on your computer using aliroot
78 // mLocalPAR: Analyze data on your computer using root + PAR files
79 // mLocalSource: Analyze data on your computer using root + source files
80                                           
81 int runFlowAnalysisOnTheFly(Int_t mode=mLocal, Int_t nEvts=10)
82 {
83  TStopwatch timer;
84  timer.Start();
85  
86  if (LYZ1SUM && LYZ2SUM)  {cout<<"WARNING: you cannot run LYZ1 and LYZ2 at the same time! LYZ2 needs the output from LYZ1.  "<<endl; exit(); }
87  if (LYZ1PROD && LYZ2PROD)  {cout<<"WARNING: you cannot run LYZ1 and LYZ2 at the same time! LYZ2 needs the output from LYZ1.  "<<endl; exit(); }
88  if (LYZ2SUM && LYZEP) {cout<<"WARNING: you cannot run LYZ2 and LYZEP at the same time! LYZEP needs the output from LYZ2."<<endl; exit(); }
89  if (LYZ1SUM && LYZEP) {cout<<"WARNING: you cannot run LYZ1 and LYZEP at the same time! LYZEP needs the output from LYZ2."<<endl; exit(); }
90  
91  cout<<endl;
92  cout<<endl;
93  cout<<"      ---- ARE YOU READY TO FLY ? ----      "<<endl;
94  cout<<endl;
95  
96  cout<<endl;
97  cout<<" ---- BEGIN FLOW ANALYSIS 'ON THE FLY' ---- "<<endl;
98  cout<<endl;
99  cout<<endl;
100  
101  LoadLibraries(mode);
102
103  // Initialize the seed for random generator
104  UInt_t sseed = 0;
105  
106  if(bSameSeed) 
107  {
108   sseed = 44; // the default constant value for seed for random generators
109  } 
110
111  if(!bSameSeed)
112  {
113   TTimeStamp dt;
114   sseed = dt.GetNanoSec()/1000;
115  }
116  
117  //---------------------------------------------------------------------------------------
118  // If the weights are used: 
119  TFile *fileWithWeights = NULL;
120  TList *listWithWeights = NULL;
121   
122  if(usePhiWeights||usePtWeights||useEtaWeights) {
123    fileWithWeights = TFile::Open("weights.root","READ");
124    if(fileWithWeights) {
125      listWithWeights = (TList*)fileWithWeights->Get("weights");
126    }
127    else
128      {cout << " WARNING: the file <weights.root> with weights from the previous run was not found."<<endl;
129       break;
130      }    
131  }
132  
133  //---------------------------------------------------------------------------------------
134  // Initialize the flowevent maker
135  AliFlowEventSimpleMakerOnTheFly* eventMakerOnTheFly = new AliFlowEventSimpleMakerOnTheFly(sseed);
136  eventMakerOnTheFly->Init();
137   
138  //---------------------------------------------------------------------------------------
139  // Initialize all the flow methods:  
140  AliFlowAnalysisWithQCumulants    *qc       = NULL;
141  AliFlowAnalysisWithCumulants     *gfc      = NULL;
142  AliFittingQDistribution          *fqd      = NULL;
143  AliFlowAnalysisWithLeeYangZeros  *lyz1sum  = NULL;
144  AliFlowAnalysisWithLeeYangZeros  *lyz1prod = NULL;
145  AliFlowAnalysisWithLeeYangZeros  *lyz2sum  = NULL;
146  AliFlowAnalysisWithLeeYangZeros  *lyz2prod = NULL;
147  AliFlowAnalysisWithLYZEventPlane *lyzep    = NULL;
148  AliFlowAnalysisWithScalarProduct *sp       = NULL;
149  AliFlowAnalysisWithMCEventPlane  *mcep     = NULL;   
150
151  // MCEP = monte carlo event plane
152  if (MCEP) {
153    AliFlowAnalysisWithMCEventPlane *mcep = new AliFlowAnalysisWithMCEventPlane();
154    mcep->Init();
155  }
156
157   // QC = Q-cumulants  
158  if(QC) { 
159    AliFlowAnalysisWithQCumulants* qc = new AliFlowAnalysisWithQCumulants();
160    qc->Init();
161    if(listWithWeights) qc->SetWeightsList(listWithWeights);
162    if(usePhiWeights) qc->SetUsePhiWeights(usePhiWeights);
163    if(usePtWeights) qc->SetUsePtWeights(usePtWeights);
164    if(useEtaWeights) qc->SetUseEtaWeights(useEtaWeights);
165  }
166   
167  // GFC = Generating Function Cumulants 
168  if(GFC) {
169    AliFlowAnalysisWithCumulants* gfc = new AliFlowAnalysisWithCumulants();
170    gfc->Init();
171    if(listWithWeights) gfc->SetWeightsList(listWithWeights);
172    if(usePhiWeights) gfc->SetUsePhiWeights(usePhiWeights);
173    if(usePtWeights) gfc->SetUsePtWeights(usePtWeights);
174    if(useEtaWeights) gfc->SetUseEtaWeights(useEtaWeights);
175  }
176  
177  // FQD = Fitting q-distribution 
178  if(FQD) {
179    AliFittingQDistribution* fqd = new AliFittingQDistribution();
180    fqd->Init();
181    if(listWithWeights) fqd->SetWeightsList(listWithWeights);
182    if(usePhiWeights) fqd->SetUsePhiWeights(usePhiWeights);  
183  }
184  
185  // SP = Scalar Product 
186  if(SP) {
187    AliFlowAnalysisWithScalarProduct* sp = new AliFlowAnalysisWithScalarProduct();
188    sp->Init();
189  }
190
191  // LYZ1 = Lee-Yang Zeroes first run
192  if(LYZ1SUM) {
193    AliFlowAnalysisWithLeeYangZeros* lyz1sum = new AliFlowAnalysisWithLeeYangZeros();
194    lyz1sum->SetFirstRun(kTRUE);
195    lyz1sum->SetUseSum(kTRUE);
196    lyz1sum->Init();
197  }
198  if(LYZ1PROD) {
199    AliFlowAnalysisWithLeeYangZeros* lyz1prod = new AliFlowAnalysisWithLeeYangZeros();
200    lyz1prod->SetFirstRun(kTRUE);
201    lyz1prod->SetUseSum(kFALSE);
202    lyz1prod->Init();
203  }
204  // LYZ2 = Lee-Yang Zeroes second run
205  if(LYZ2SUM) {
206    AliFlowAnalysisWithLeeYangZeros* lyz2sum = new AliFlowAnalysisWithLeeYangZeros();
207    // read the input file from the first run 
208    TString inputFileNameLYZ2SUM = "outputLYZ1SUManalysis.root" ;
209    TFile* inputFileLYZ2SUM = new TFile(inputFileNameLYZ2SUM.Data(),"READ");
210    if(!inputFileLYZ2SUM || inputFileLYZ2SUM->IsZombie()) { 
211      cerr << " ERROR: To run LYZ2SUM you need the output file from LYZ1SUM. This file is not there! Please run LYZ1SUM first." << endl ;
212      break; 
213    }
214    else { 
215      TList* inputListLYZ2SUM = (TList*)inputFileLYZ2SUM->Get("cobjLYZ1SUM");  
216      if (!inputListLYZ2SUM) {cout<<"Input list LYZ2SUM is NULL pointer!"<<endl; break;}
217      else {
218        cout<<"LYZ2SUM input file/list read..."<<endl;
219        lyz2sum->SetFirstRunList(inputListLYZ2SUM);
220        lyz2sum->SetFirstRun(kFALSE);
221        lyz2sum->SetUseSum(kTRUE);
222        lyz2sum->Init();
223      }
224    }
225  }
226  if(LYZ2PROD) {
227    AliFlowAnalysisWithLeeYangZeros* lyz2prod = new AliFlowAnalysisWithLeeYangZeros();
228    // read the input file from the first run 
229    TString inputFileNameLYZ2PROD = "outputLYZ1PRODanalysis.root" ;
230    TFile* inputFileLYZ2PROD = new TFile(inputFileNameLYZ2PROD.Data(),"READ");
231    if(!inputFileLYZ2PROD || inputFileLYZ2PROD->IsZombie()) { 
232      cerr << " ERROR: To run LYZ2PROD you need the output file from LYZ1PROD. This file is not there! Please run LYZ1PROD first." << endl ;
233      break; 
234    }
235    else { 
236      TList* inputListLYZ2PROD = (TList*)inputFileLYZ2PROD->Get("cobjLYZ1PROD");  
237      if (!inputListLYZ2PROD) {cout<<"Input list LYZ2PROD is NULL pointer!"<<endl; break;}
238      else {
239        cout<<"LYZ2PROD input file/list read..."<<endl;
240        lyz2prod->SetFirstRunList(inputListLYZ2PROD);
241        lyz2prod->SetFirstRun(kFALSE);
242        lyz2prod->SetUseSum(kFALSE);
243        lyz2prod->Init();
244      }
245    }
246  }
247
248  // LYZEP = Lee-Yang Zeroes event plane
249  if(LYZEP) {
250    AliFlowLYZEventPlane* ep = new AliFlowLYZEventPlane() ;
251    AliFlowAnalysisWithLYZEventPlane* lyzep = new AliFlowAnalysisWithLYZEventPlane();
252    // read the input file from the second lyz run 
253    TString inputFileNameLYZEP = "outputLYZ2SUManalysis.root" ;
254    TFile* inputFileLYZEP = new TFile(inputFileNameLYZEP.Data(),"READ");
255    if(!inputFileLYZEP || inputFileLYZEP->IsZombie()) { 
256      cerr << " ERROR: To run LYZEP you need the output file from LYZ2SUM. This file is not there! Please run LYZ2SUM first." << endl ; 
257      break;
258    }
259    else { 
260      TList* inputListLYZEP = (TList*)inputFileLYZEP->Get("cobjLYZ2SUM");  
261      if (!inputListLYZEP) {cout<<"Input list LYZEP is NULL pointer!"<<endl; break;}
262      else {
263        cout<<"LYZEP input file/list read..."<<endl;
264        ep   ->SetSecondRunList(inputListLYZEP);
265        lyzep->SetSecondRunList(inputListLYZEP);
266        ep   ->Init();
267        lyzep->Init();
268      }
269    }
270  }
271  //---------------------------------------------------------------------------------------
272   
273  // set the global event parameters: 
274  eventMakerOnTheFly->SetNoOfLoops(iLoops);
275  
276  if(bMultDistrOfRPsIsGauss)
277  {
278   eventMakerOnTheFly->SetMultDistrOfRPsIsGauss(bMultDistrOfRPsIsGauss);
279   eventMakerOnTheFly->SetMultiplicityOfRP(iMultiplicityOfRP);
280   eventMakerOnTheFly->SetMultiplicitySpreadOfRP(dMultiplicitySpreadOfRP);
281  } else
282    {
283     eventMakerOnTheFly->SetMultDistrOfRPsIsGauss(bMultDistrOfRPsIsGauss);
284     eventMakerOnTheFly->SetMinMultOfRP(iMinMultOfRP);
285     eventMakerOnTheFly->SetMaxMultOfRP(iMaxMultOfRP); 
286    }
287   
288  eventMakerOnTheFly->SetTemperatureOfRP(dTemperatureOfRP);
289
290  eventMakerOnTheFly->SetV1RP(dV1RP);
291  eventMakerOnTheFly->SetV1SpreadRP(dV1SpreadRP);  
292  eventMakerOnTheFly->SetV4RP(dV4RP);
293  eventMakerOnTheFly->SetV4SpreadRP(dV4SpreadRP);  
294  
295  // constant harmonic V2:
296  if(bConstantHarmonics)
297  { 
298   eventMakerOnTheFly->SetUseConstantHarmonics(bConstantHarmonics);
299   if(bV2DistrOfRPsIsGauss)
300   {
301    eventMakerOnTheFly->SetV2DistrOfRPsIsGauss(bV2DistrOfRPsIsGauss);
302    eventMakerOnTheFly->SetV2RP(dV2RP);
303    eventMakerOnTheFly->SetV2SpreadRP(dV2SpreadRP);  
304   } else
305     {
306      eventMakerOnTheFly->SetV2DistrOfRPsIsGauss(bV2DistrOfRPsIsGauss);
307      eventMakerOnTheFly->SetMinV2RP(dMinV2RP);
308      eventMakerOnTheFly->SetMaxV2RP(dMaxV2RP);  
309     }
310  }
311  // (pt,eta) dependent harmonic V2:
312  if(!bConstantHarmonics)
313  {
314   eventMakerOnTheFly->SetUseConstantHarmonics(bConstantHarmonics);
315   eventMakerOnTheFly->SetV2RPMax(dV2RPMax);
316   eventMakerOnTheFly->SetPtCutOff(dPtCutOff);  
317  }
318        
319  //---------------------------------------------------------------------------------------  
320  // create and analyze events 'on the fly':
321
322  for(Int_t i=0;i<nEvts;i++) {   
323    // creating the event with above settings:
324    AliFlowEventSimple *event = eventMakerOnTheFly->CreateEventOnTheFly(); 
325    
326    // analyzing the created event 'on the fly':
327    // do flow analysis for various methods:
328    if(MCEP)    mcep->Make(event);
329    if(QC)      qc->Make(event);
330    if(GFC)     gfc->Make(event);
331    if(FQD)     fqd->Make(event);
332    if(LYZ1SUM) lyz1sum->Make(event);
333    if(LYZ1PROD)lyz1prod->Make(event);
334    if(LYZ2SUM) lyz2sum->Make(event);
335    if(LYZ2PROD)lyz2prod->Make(event);
336    if(LYZEP)   lyzep->Make(event,ep);
337    if(SP)      sp->Make(event);
338    
339    delete event;
340  } // end of for(Int_t i=0;i<nEvts;i++)
341  //---------------------------------------------------------------------------------------  
342
343
344
345  //---------------------------------------------------------------------------------------  
346  // calculating and storing the final results of flow analysis
347  if(MCEP)    {mcep->Finish();    mcep->WriteHistograms("outputMCEPanalysis.root");}
348  if(SP)      {sp->Finish();      sp->WriteHistograms("outputSPanalysis.root");}
349  if(QC)      {qc->Finish();      qc->WriteHistograms("outputQCanalysis.root");}
350  if(GFC)     {gfc->Finish();     gfc->WriteHistograms("outputGFCanalysis.root");}
351  if(FQD)     {fqd->Finish();     fqd->WriteHistograms("outputFQDanalysis.root");}
352  if(LYZ1SUM) {lyz1sum->Finish(); lyz1sum->WriteHistograms("outputLYZ1SUManalysis.root");}
353  if(LYZ1PROD){lyz1prod->Finish();lyz1prod->WriteHistograms("outputLYZ1PRODanalysis.root");}
354  if(LYZ2SUM) {lyz2sum->Finish(); lyz2sum->WriteHistograms("outputLYZ2SUManalysis.root");}
355  if(LYZ2PROD){lyz2prod->Finish();lyz2prod->WriteHistograms("outputLYZ2PRODanalysis.root");}
356  if(LYZEP)   {lyzep->Finish();   lyzep->WriteHistograms("outputLYZEPanalysis.root");}
357  //---------------------------------------------------------------------------------------  
358  
359  
360  
361  cout<<endl;
362  cout<<endl;
363  cout<<" ---- LANDED SUCCESSFULLY ---- "<<endl;
364  cout<<endl; 
365  
366  timer.Stop();
367  cout << endl;
368  timer.Print();
369 }
370
371 void SetupPar(char* pararchivename)
372 {
373   //Load par files, create analysis libraries
374   //For testing, if par file already decompressed and modified
375   //classes then do not decompress.
376  
377   TString cdir(Form("%s", gSystem->WorkingDirectory() )) ; 
378   TString parpar(Form("%s.par", pararchivename)) ; 
379   if ( gSystem->AccessPathName(parpar.Data()) ) {
380     gSystem->ChangeDirectory(gSystem->Getenv("ALICE_ROOT")) ;
381     TString processline(Form(".! make %s", parpar.Data())) ; 
382     gROOT->ProcessLine(processline.Data()) ;
383     gSystem->ChangeDirectory(cdir) ; 
384     processline = Form(".! mv /tmp/%s .", parpar.Data()) ;
385     gROOT->ProcessLine(processline.Data()) ;
386   } 
387   if ( gSystem->AccessPathName(pararchivename) ) {  
388     TString processline = Form(".! tar xvzf %s",parpar.Data()) ;
389     gROOT->ProcessLine(processline.Data());
390   }
391   
392   TString ocwd = gSystem->WorkingDirectory();
393   gSystem->ChangeDirectory(pararchivename);
394   
395   // check for BUILD.sh and execute
396   if (!gSystem->AccessPathName("PROOF-INF/BUILD.sh")) {
397     printf("*******************************\n");
398     printf("*** Building PAR archive    ***\n");
399     cout<<pararchivename<<endl;
400     printf("*******************************\n");
401     
402     if (gSystem->Exec("PROOF-INF/BUILD.sh")) {
403       Error("runProcess","Cannot Build the PAR Archive! - Abort!");
404       return -1;
405     }
406   }
407   // check for SETUP.C and execute
408   if (!gSystem->AccessPathName("PROOF-INF/SETUP.C")) {
409     printf("*******************************\n");
410     printf("*** Setup PAR archive       ***\n");
411     cout<<pararchivename<<endl;
412     printf("*******************************\n");
413     gROOT->Macro("PROOF-INF/SETUP.C");
414   }
415   
416   gSystem->ChangeDirectory(ocwd.Data());
417   printf("Current dir: %s\n", ocwd.Data());
418 }
419
420 void LoadLibraries(const anaModes mode) {
421   
422   //--------------------------------------
423   // Load the needed libraries most of them already loaded by aliroot
424   //--------------------------------------
425   gSystem->Load("libTree.so");
426   gSystem->Load("libGeom.so");
427   gSystem->Load("libVMC.so");
428   gSystem->Load("libXMLIO.so");
429   gSystem->Load("libPhysics.so");
430   
431   //----------------------------------------------------------
432   // >>>>>>>>>>> Local mode <<<<<<<<<<<<<< 
433   //----------------------------------------------------------
434   if (mode==mLocal) {
435     //--------------------------------------------------------
436     // If you want to use already compiled libraries 
437     // in the aliroot distribution
438     //--------------------------------------------------------
439     gSystem->Load("libSTEERBase");
440     gSystem->Load("libESD");
441     gSystem->Load("libAOD");
442     gSystem->Load("libANALYSIS");
443     gSystem->Load("libANALYSISalice");
444     gSystem->Load("libCORRFW.so");
445     cerr<<"libCORRFW.so loaded..."<<endl;
446     gSystem->Load("libPWG2flowCommon.so");
447     cerr<<"libPWG2flowCommon.so loaded..."<<endl;
448     gSystem->Load("libPWG2flowTasks.so");
449     cerr<<"libPWG2flowTasks.so loaded..."<<endl;
450   }
451   
452   else if (mode == mLocalPAR) {
453     //--------------------------------------------------------
454     //If you want to use root and par files from aliroot
455     //--------------------------------------------------------  
456      //If you want to use root and par files from aliroot
457     //--------------------------------------------------------  
458     SetupPar("STEERBase");
459     SetupPar("ESD");
460     SetupPar("AOD");
461     SetupPar("ANALYSIS");
462     SetupPar("ANALYSISalice");
463     SetupPar("PWG2AOD");
464     SetupPar("CORRFW");
465     SetupPar("PWG2flowCommon");
466     cerr<<"PWG2flowCommon.par loaded..."<<endl;
467     SetupPar("PWG2flowTasks");
468     cerr<<"PWG2flowTasks.par loaded..."<<endl;
469   }
470   
471   //---------------------------------------------------------
472   // <<<<<<<<<< Source mode >>>>>>>>>>>>
473   //---------------------------------------------------------
474   else if (mode==mLocalSource) {
475  
476     // In root inline compile
477
478    
479     // Constants  
480     gROOT->LoadMacro("AliFlowCommon/AliFlowCommonConstants.cxx+");
481     gROOT->LoadMacro("AliFlowCommon/AliFlowLYZConstants.cxx+");
482     gROOT->LoadMacro("AliFlowCommon/AliFlowCumuConstants.cxx+");
483     
484     // Flow event
485     gROOT->LoadMacro("AliFlowCommon/AliFlowVector.cxx+"); 
486     gROOT->LoadMacro("AliFlowCommon/AliFlowTrackSimple.cxx+");    
487     gROOT->LoadMacro("AliFlowCommon/AliFlowEventSimple.cxx+");
488     
489     // Cuts
490     gROOT->LoadMacro("AliFlowCommon/AliFlowTrackSimpleCuts.cxx+");    
491     
492     // Output histosgrams
493     gROOT->LoadMacro("AliFlowCommon/AliFlowCommonHist.cxx+");
494     gROOT->LoadMacro("AliFlowCommon/AliFlowCommonHistResults.cxx+");
495     gROOT->LoadMacro("AliFlowCommon/AliFlowLYZHist1.cxx+");
496     gROOT->LoadMacro("AliFlowCommon/AliFlowLYZHist2.cxx+");
497     
498     // Functions needed for various methods
499     gROOT->LoadMacro("AliFlowCommon/AliCumulantsFunctions.cxx+");
500     gROOT->LoadMacro("AliFlowCommon/AliFittingFunctionsForQDistribution.cxx+");
501     gROOT->LoadMacro("AliFlowCommon/AliFlowLYZEventPlane.cxx+");
502     
503     // Flow Analysis code for various methods
504     gROOT->LoadMacro("AliFlowCommon/AliFlowAnalysisWithMCEventPlane.cxx+"); 
505     gROOT->LoadMacro("AliFlowCommon/AliFlowAnalysisWithScalarProduct.cxx+");
506     gROOT->LoadMacro("AliFlowCommon/AliFlowAnalysisWithLYZEventPlane.cxx+");
507     gROOT->LoadMacro("AliFlowCommon/AliFlowAnalysisWithLeeYangZeros.cxx+");
508     gROOT->LoadMacro("AliFlowCommon/AliFlowAnalysisWithCumulants.cxx+");
509     gROOT->LoadMacro("AliFlowCommon/AliFlowAnalysisWithQCumulants.cxx+"); 
510     gROOT->LoadMacro("AliFlowCommon/AliFittingQDistribution.cxx+");
511     
512     // Class to fill the FlowEvent on the fly (generate Monte Carlo events)
513     gROOT->LoadMacro("AliFlowCommon/AliFlowEventSimpleMakerOnTheFly.cxx+");   
514     
515     cout << "finished loading macros!" << endl;  
516     
517   }  
518   
519 }
520
521