]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGCF/FLOW/Documentation/examples/manual/runFlowOnTheFlyExample.C
Merge branch 'feature-movesplit'
[u/mrichter/AliRoot.git] / PWGCF / FLOW / Documentation / examples / manual / runFlowOnTheFlyExample.C
1 #include "TStopwatch.h"
2 #include "Riostream.h"
3 #include "TFile.h"
4
5
6 int runFlowOnTheFlyExample(Int_t nEvts=2000, Int_t mult=1000, Float_t v2=0.05, Int_t iseed=7669)
7 {
8   TStopwatch timer;
9   timer.Start();
10   
11   // Load the needed libraries for root (in AliRoot already loaded)
12   gSystem->Load("libGeom");
13   gSystem->Load("libVMC");
14   gSystem->Load("libXMLIO");
15   gSystem->Load("libPhysics");
16   gSystem->Load("libPWGflowBase");
17
18   fMyTRandom3 = new TRandom3(iseed);   
19   gRandom->SetSeed(fMyTRandom3->Integer(65539));
20
21   // Initialize the flow methods for default analysis:
22   AliFlowAnalysisWithMCEventPlane *mcep = new AliFlowAnalysisWithMCEventPlane();
23   //mcep->SetHarmonic(2); // default is v2
24   mcep->Init();
25
26   AliFlowAnalysisWithQCumulants* qc = new AliFlowAnalysisWithQCumulants();
27   // qc->SetHarmonic(2); // default is v2
28   qc->Init();
29   
30
31   // set cuts for the Reference Particles and Particles Of Interest:
32   AliFlowTrackSimpleCuts *cutsRP = new AliFlowTrackSimpleCuts();
33   //  cutsRP->SetPtMax(ptMaxRP);
34   AliFlowTrackSimpleCuts *cutsPOI = new AliFlowTrackSimpleCuts();
35   cutsPOI->SetPtMin(0.2);
36   cutsPOI->SetPtMax(2.0);
37
38   Printf("starting the main event loop..");
39   // create and analyze events 'on the fly':
40   for(Int_t i=0; i<nEvts; i++)
41     {
42       // creating the event with above settings:
43       AliFlowEventSimple* event = new AliFlowEventSimple(mult,AliFlowEventSimple::kGenerate);
44        event->AddV2(v2);
45       //event->TagTracks(cutsRP, cutsPOI);
46       event->TagPOI(cutsPOI);
47       event->TagRP(cutsRP);
48       // event->Print();
49       
50       // do flow analysis for various methods:
51       mcep->Make(event);
52       qc->Make(event);
53       cout <<"Event: " << i+1 << "\r"; cout.flush();
54       delete event;
55     } // end of for(Int_t i=0;i<nEvts;i++)
56   
57  
58   // calculate the final results
59   mcep->Finish();
60   qc->Finish();
61  
62   // open a new file which will hold the final results of all methods:
63   TString outputFileName = "AnalysisResults.root";
64   TFile *outputFile = new TFile(outputFileName.Data(),"RECREATE");
65   const Int_t nMethods = 2;
66   TString method[nMethods] = {"MCEP","QC"};
67   TDirectoryFile *dirFileFinal[nMethods] = {NULL};
68   TString fileName[nMethods];
69   for(Int_t i=0; i<nMethods; i++)
70     {
71       // form a file name for each method:
72       fileName[i]+="output";
73       fileName[i]+=method[i].Data();
74       fileName[i]+="analysis";
75       dirFileFinal[i] = new TDirectoryFile(fileName[i].Data(),fileName[i].Data());
76     }
77   
78   // store the final results
79   mcep->WriteHistograms(dirFileFinal[0]);
80   qc->WriteHistograms(dirFileFinal[1]);
81   
82   outputFile->Close();
83   
84   delete outputFile;
85   cout<<endl; cout<<" ---- Fini ---- "<<endl; cout<<endl;
86   
87   timer.Stop(); timer.Print();
88 }
89
90
91
92