]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/FLOW/macros/redoFinish.C
71a634654b65c85e6d1dd9b9e14e4a8a2184ba55
[u/mrichter/AliRoot.git] / PWG2 / FLOW / macros / redoFinish.C
1 // Macro redoFinish.C is used after macro mergeOutput.C has been used for merging. 
2 // Before using this macro first read the explanation at the beginning of macro mergeOutput.C. 
3
4 enum libModes {mLocal,mLocalSource};
5
6 void redoFinish(TString type="ESD", Int_t mode=mLocal)
7 {
8  // type: type of analysis can be ESD, AOD, MC, ESDMC0, ESDMC1
9  //       (if type="" output files are from MC simulation (default))
10  // mode: if mode = mLocal: analyze data on your computer using aliroot
11  //       if mode = mLocalSource: analyze data on your computer using root + source files
12
13  // Name of merged, large statistics file obtained with macro mergeOutput.C: 
14  TString mergedFileName = "mergedAnalysisResults.root";
15  // Final output file name holding final results for large statistics sample:
16  TString outputFileName = "AnalysisResults.root"; 
17  
18  const Int_t nMethods = 10;
19  TString method[nMethods] = {"MCEP","SP","GFC","QC","FQD","LYZ1SUM","LYZ1PROD","LYZ2SUM","LYZ2PROD","LYZEP"};
20  
21  // Load needed libraries:                       
22  LoadLibrariesRF(mode);  
23    
24  // Accessing the merged, large statistics file obtained with macro mergeOutput.C.
25  // On this file the flow analysis will be redone, its content modified to store
26  // the correct final results and eventually renamed into TString outputFileName.
27  TString pwd(gSystem->pwd());
28  pwd+="/";
29  pwd+=mergedFileName.Data();
30  TFile *mergedFile = NULL;
31  if(gSystem->AccessPathName(pwd.Data(),kFileExists))
32  {
33   cout<<"WARNING: You do not have a merged output file:"<<endl;
34   cout<<"         "<<pwd.Data()<<endl;
35   cout<<endl;
36   cout<<"In order to get that file use macro mergeOutput.C first."<<endl;
37   cout<<endl;
38   exit(0);
39  } else 
40    {
41     // Create temporarily copy of "mergedAnalysisResults.root":
42     TSystemFile *fileTemp = new TSystemFile(mergedFileName.Data(),".");
43     fileTemp->Copy("mergedAnalysisResultsTemp.root");
44     delete fileTemp;
45     // Access merged file:
46     mergedFile = TFile::Open(pwd.Data(),"UPDATE");
47    }
48    
49  // Access from mergedFile the merged files for each method and from them the lists holding histograms:
50  TString fileName[nMethods]; 
51  TDirectoryFile *dirFile[nMethods] = {NULL}; 
52  TString listName[nMethods]; 
53  TList *list[nMethods] = {NULL};
54  for(Int_t i=0;i<nMethods;i++)
55  {
56   // Form a file name for each method:
57   fileName[i]+="output";
58   fileName[i]+=method[i].Data();
59   fileName[i]+="analysis";
60   fileName[i]+=type.Data();
61   // Access this file:
62   dirFile[i] = (TDirectoryFile*)mergedFile->FindObjectAny(fileName[i].Data());
63   // Form a list name for each method:
64   listName[i]+="cobj";
65   listName[i]+=method[i].Data();
66   // Access this list:
67   if(dirFile[i])
68   {
69    dirFile[i]->GetObject(listName[i].Data(),list[i]);
70   } else 
71     {
72      cout<<"WARNING: Couldn't find a file "<<fileName[i].Data()<<".root !!!!"<<endl;
73     }
74  } // End of for(Int_t i=0;i<nMethods;i++)
75
76  // Redo finish for each method (REMARK: this implementation can be dramatically improved!):
77  // MCEP:
78  for(Int_t i=0;i<nMethods;i++)
79  {
80   if(list[i] && strcmp(list[i]->GetName(),"cobjMCEP")==0)
81   {
82    AliFlowAnalysisWithMCEventPlane* mcep = new AliFlowAnalysisWithMCEventPlane();
83    mcep->GetOutputHistograms(list[i]);
84    mcep->Finish();
85    dirFile[i]->Add(list[i],kTRUE);
86    dirFile[i]->Write(dirFile[i]->GetName(),TObject::kSingleKey+TObject::kOverwrite);
87   } 
88   // SP:
89   else if(list[i] && strcmp(list[i]->GetName(),"cobjSP")==0)
90   {
91    AliFlowAnalysisWithScalarProduct* sp = new AliFlowAnalysisWithScalarProduct();
92    sp->GetOutputHistograms(list[i]);
93    sp->Finish();
94    dirFile[i]->Add(list[i],kTRUE);
95    dirFile[i]->Write(dirFile[i]->GetName(),TObject::kSingleKey+TObject::kOverwrite);
96   } 
97   // GFC:
98   else if(list[i] && strcmp(list[i]->GetName(),"cobjGFC")==0)
99   {
100    AliFlowAnalysisWithCumulants* gfc = new AliFlowAnalysisWithCumulants();
101    gfc->GetOutputHistograms(list[i]);
102    gfc->Finish();
103    dirFile[i]->Add(list[i],kTRUE);
104    dirFile[i]->Write(dirFile[i]->GetName(),TObject::kSingleKey+TObject::kOverwrite);
105   } 
106   // QC:
107   else if(list[i] && strcmp(list[i]->GetName(),"cobjQC")==0)
108   {
109    AliFlowAnalysisWithQCumulants* qc = new AliFlowAnalysisWithQCumulants();
110    qc->GetOutputHistograms(list[i]);
111    qc->Finish();
112    dirFile[i]->Add(list[i],kTRUE);
113    dirFile[i]->Write(dirFile[i]->GetName(),TObject::kSingleKey+TObject::kOverwrite);
114   } 
115   // FQD:
116   else if(list[i] && strcmp(list[i]->GetName(),"cobjFQD")==0)
117   {
118    AliFlowAnalysisWithFittingQDistribution* fqd = new AliFlowAnalysisWithFittingQDistribution();
119    fqd->GetOutputHistograms(list[i]);
120    fqd->Finish(kTRUE);
121    dirFile[i]->Add(list[i],kTRUE);
122    dirFile[i]->Write(dirFile[i]->GetName(),TObject::kSingleKey+TObject::kOverwrite);
123   }
124   // LYZ1SUM:
125   else if(list[i] && strcmp(list[i]->GetName(),"cobjLYZ1SUM")==0)
126   {
127    AliFlowAnalysisWithLeeYangZeros* lyz1sum = new AliFlowAnalysisWithLeeYangZeros();
128    lyz1sum->SetFirstRun(kTRUE);   
129    lyz1sum->SetUseSum(kTRUE);       
130    lyz1sum->GetOutputHistograms(list[i]);
131    lyz1sum->Finish();
132    dirFile[i]->Add(list[i],kTRUE);
133    dirFile[i]->Write(dirFile[i]->GetName(),TObject::kSingleKey+TObject::kOverwrite);
134   } 
135   // LYZ2SUM:
136   else if(list[i] && strcmp(list[i]->GetName(),"cobjLYZ2SUM")==0)
137   {
138    AliFlowAnalysisWithLeeYangZeros* lyz2sum = new AliFlowAnalysisWithLeeYangZeros();
139    lyz2sum->SetFirstRun(kFALSE);   
140    lyz2sum->SetUseSum(kTRUE);       
141    lyz2sum->GetOutputHistograms(list[i]);
142    lyz2sum->Finish();
143    dirFile[i]->Add(list[i],kTRUE);
144    dirFile[i]->Write(dirFile[i]->GetName(),TObject::kSingleKey+TObject::kOverwrite);
145   }
146   // LYZ1PROD:
147   else if(list[i] && strcmp(list[i]->GetName(),"cobjLYZ1PROD")==0)
148   {
149    AliFlowAnalysisWithLeeYangZeros* lyz1prod = new AliFlowAnalysisWithLeeYangZeros();
150    lyz1prod->SetFirstRun(kTRUE);   
151    lyz1prod->SetUseSum(kFALSE);       
152    lyz1prod->GetOutputHistograms(list[i]);
153    lyz1prod->Finish();
154    dirFile[i]->Add(list[i],kTRUE);
155    dirFile[i]->Write(dirFile[i]->GetName(),TObject::kSingleKey+TObject::kOverwrite);
156   }    
157   // LYZ2PROD:
158   else if(list[i] && strcmp(list[i]->GetName(),"cobjLYZ2PROD")==0)
159   {
160    AliFlowAnalysisWithLeeYangZeros* lyz2prod = new AliFlowAnalysisWithLeeYangZeros();
161    lyz2prod->SetFirstRun(kFALSE);   
162    lyz2prod->SetUseSum(kFALSE);       
163    lyz2prod->GetOutputHistograms(list[i]);
164    lyz2prod->Finish();
165    dirFile[i]->Add(list[i],kTRUE);
166    dirFile[i]->Write(dirFile[i]->GetName(),TObject::kSingleKey+TObject::kOverwrite);
167   }
168   // LYZEP:
169   else if(list[i] && strcmp(list[i]->GetName(),"cobjLYZEP")==0)
170   {
171    AliFlowAnalysisWithLYZEventPlane* lyzep = new AliFlowAnalysisWithLYZEventPlane();
172    lyzep->GetOutputHistograms(list[i]);
173    lyzep->Finish();
174    dirFile[i]->Add(list[i],kTRUE);
175    dirFile[i]->Write(dirFile[i]->GetName(),TObject::kSingleKey+TObject::kOverwrite);
176   }                
177  } // End of for(Int_t i=0;i<nMethods;i++)
178
179  // Close the final output file:
180  mergedFile->Close();
181  delete mergedFile;
182  
183  // Giving the final names:
184  TSystemFile *outputFileFinal = new TSystemFile(mergedFileName.Data(),".");
185  outputFileFinal->Rename(outputFileName.Data());
186  delete outputFileFinal; 
187  TSystemFile *mergedFileFinal = new TSystemFile("mergedAnalysisResultsTemp.root",".");
188  mergedFileFinal->Rename(mergedFileName.Data());
189  delete mergedFileFinal;
190      
191 } // End of void redoFinish(Int_t mode=mLocal)
192  
193 void LoadLibrariesRF(const libModes mode) {
194   
195   //--------------------------------------
196   // Load the needed libraries most of them already loaded by aliroot
197   //--------------------------------------
198   //gSystem->Load("libTree");
199   gSystem->Load("libGeom");
200   gSystem->Load("libVMC");
201   gSystem->Load("libXMLIO");
202   gSystem->Load("libPhysics");
203   
204   //----------------------------------------------------------
205   // >>>>>>>>>>> Local mode <<<<<<<<<<<<<< 
206   //----------------------------------------------------------
207   if (mode==mLocal) {
208     //--------------------------------------------------------
209     // If you want to use already compiled libraries 
210     // in the aliroot distribution
211     //--------------------------------------------------------
212
213   //==================================================================================  
214   //load needed libraries:
215   gSystem->AddIncludePath("-I$ROOTSYS/include");
216   //gSystem->Load("libTree");
217
218   // for AliRoot
219   gSystem->AddIncludePath("-I$ALICE_ROOT/include");
220   gSystem->Load("libANALYSIS");
221   gSystem->Load("libPWG2flowCommon");
222   //cerr<<"libPWG2flowCommon loaded ..."<<endl;
223   
224   }
225   
226   else if (mode==mLocalSource) {
227  
228     // In root inline compile
229   
230     // Constants  
231     gROOT->LoadMacro("AliFlowCommon/AliFlowCommonConstants.cxx+");
232     gROOT->LoadMacro("AliFlowCommon/AliFlowLYZConstants.cxx+");
233     gROOT->LoadMacro("AliFlowCommon/AliFlowCumuConstants.cxx+");
234     
235     // Flow event
236     gROOT->LoadMacro("AliFlowCommon/AliFlowVector.cxx+"); 
237     gROOT->LoadMacro("AliFlowCommon/AliFlowTrackSimple.cxx+");    
238     gROOT->LoadMacro("AliFlowCommon/AliFlowEventSimple.cxx+");
239     
240     // Cuts
241     gROOT->LoadMacro("AliFlowCommon/AliFlowTrackSimpleCuts.cxx+");    
242     
243     // Output histosgrams
244     gROOT->LoadMacro("AliFlowCommon/AliFlowCommonHist.cxx+");
245     gROOT->LoadMacro("AliFlowCommon/AliFlowCommonHistResults.cxx+");
246     gROOT->LoadMacro("AliFlowCommon/AliFlowLYZHist1.cxx+");
247     gROOT->LoadMacro("AliFlowCommon/AliFlowLYZHist2.cxx+");
248        
249     cout << "finished loading macros!" << endl;  
250     
251   } // end of else if (mode==mLocalSource) 
252   
253 } // end of void LoadLibrariesRF(const libModes mode)
254
255