]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGCF/FLOW/macros/redoFinish.C
updates to check K0 charged kaon dif
[u/mrichter/AliRoot.git] / PWGCF / FLOW / macros / redoFinish.C
1 // Macro redoFinish.C is typically used after the merging macros (mergeOutput.C or
2 // mergeOutputOnGrid.C) have been used to produce the merged, large statistics
3 // file of flow analysis. Results stored in merged file are WRONG because after
4 // merging the results from small statistics files are trivially summed up in all
5 // histograms. This is taken into account and corrected for with macro redoFinish.C.
6 // Another typical use of the macro redoFinish.C is to repeat the call to Finish()
7 // in all classes, but with different values of some settings which might modify
8 // the final results (Example: redo the Finish() and apply correction for detector
9 // effects in QC code because by default this correction is switched off).
10
11 // Name of the merged, large statistics file obtained with the merging macros:
12 TString mergedFileName = "mergedAnalysisResults.root";
13 // Final output file name holding correct final results for large statistics sample:
14 TString outputFileName = "AnalysisResults.root";
15
16 Bool_t bApplyCorrectionForNUA = kFALSE; // apply correction for non-uniform acceptance
17 Bool_t bApplyCorrectionForNUAVsM = kFALSE; // apply correction for non-uniform acceptance in each multiplicity bin independently
18 Bool_t bPropagateErrorAlsoFromNIT = kFALSE; // propagate error also from non-isotropic terms
19 Bool_t bMinimumBiasReferenceFlow = kTRUE; // store in CRH for reference flow the result obtained wihout rebinning in multiplicity (kTRUE)
20
21 void redoFinish()
22 {
23   LoadLibraries();
24
25   TString mergedFileFullPathName(gSystem->pwd());
26   mergedFileFullPathName+="/";
27   mergedFileFullPathName+=mergedFileName.Data();
28   TFile *mergedFile = NULL;
29   if(gSystem->AccessPathName(mergedFileFullPathName.Data(),kFileExists))
30   {
31     cout<<endl;
32     cout<<" WARNING: Couldn't find a file: "<<mergedFileName.Data()<<endl;
33     cout<<"          in directory "<<gSystem->pwd()<<" !!!!"<<endl;
34     cout<<endl;
35     exit(0);
36   }
37   else
38   {
39     // Create temporarily copy of <mergedFileName> if neccessary:
40     if(!(mergedFileName == outputFileName))
41     {
42       TSystemFile *fileTemp = new TSystemFile(mergedFileFullPathName.Data(),".");
43       fileTemp->Copy("mergedAnalysisResultsTemp.root");
44       delete fileTemp;
45     }
46     // Access <mergedFileName>:
47     mergedFile = TFile::Open(mergedFileFullPathName.Data(),"UPDATE");
48   }
49
50   // Access from <mergedFileName> the merged TDirectoryFile's for each method and from them the lists holding histograms:
51   
52   TList* mergedFileKeys = mergedFile->GetListOfKeys();
53   for(Int_t i=0; i<mergedFileKeys->GetEntries(); i++)
54   {
55     TDirectory* directory = dynamic_cast<TDirectory*>(mergedFile->Get(mergedFileKeys->At(i)->GetName()));
56     if (!directory) continue;
57
58     TList* listTemp = directory->GetListOfKeys();
59     if(!listTemp) continue;
60     for (Int_t icent=0; icent<listTemp->GetEntries(); icent++)
61     {
62       TList* list = dynamic_cast<TList*>(directory->Get(listTemp->At(icent)->GetName()));
63       if (!list) continue;
64
65       ////////////////////
66       if(TString(list->GetName()).Contains("MCEP"))
67       {
68         AliFlowAnalysisWithMCEventPlane* mcep = new AliFlowAnalysisWithMCEventPlane();
69         mcep->GetOutputHistograms(list);
70         mcep->Finish();
71         directory->Add(list,kTRUE);
72         directory->Write(directory->GetName(),TObject::kSingleKey+TObject::kWriteDelete);
73       }
74       // SP:
75       else if(TString(list->GetName()).Contains("SP"))
76       {
77         AliFlowAnalysisWithScalarProduct* sp = new AliFlowAnalysisWithScalarProduct();
78         sp->GetOutputHistograms(list);
79         //sp->GetHistProFlags()->SetBinContent(1,(Int_t)bApplyCorrectionForNUA);
80         sp->Finish();
81         directory->Add(list,kTRUE);
82         directory->Write(directory->GetName(),TObject::kSingleKey+TObject::kWriteDelete);
83       }
84       // GFC:
85       else if(TString(list->GetName()).Contains("GFC"))
86       {
87         AliFlowAnalysisWithCumulants* gfc = new AliFlowAnalysisWithCumulants();
88         gfc->GetOutputHistograms(list);
89         gfc->Finish();
90         directory->Add(list,kTRUE);
91         directory->Write(directory->GetName(),TObject::kSingleKey+TObject::kWriteDelete);
92       }
93       // QC:
94       else if(TString(list->GetName()).Contains("QC"))
95       {
96         AliFlowAnalysisWithQCumulants* qc = new AliFlowAnalysisWithQCumulants();
97         qc->GetOutputHistograms(list);
98         qc->GetIntFlowFlags()->SetBinContent(3,(Int_t)bApplyCorrectionForNUA);
99         qc->GetIntFlowFlags()->SetBinContent(8,(Int_t)bApplyCorrectionForNUAVsM);
100         qc->GetIntFlowFlags()->SetBinContent(9,(Int_t)bPropagateErrorAlsoFromNIT);
101         qc->GetIntFlowFlags()->SetBinContent(11,(Int_t)bMinimumBiasReferenceFlow);
102         qc->Finish();
103         directory->Add(list,kTRUE);
104         directory->Write(directory->GetName(),TObject::kSingleKey+TObject::kWriteDelete);
105       }
106       // FQD:
107       else if(TString(list->GetName()).Contains("FQD"))
108       {
109         AliFlowAnalysisWithFittingQDistribution* fqd = new AliFlowAnalysisWithFittingQDistribution();
110         fqd->GetOutputHistograms(list);
111         fqd->Finish(kTRUE);
112         directory->Add(list,kTRUE);
113         directory->Write(directory->GetName(),TObject::kSingleKey+TObject::kWriteDelete);
114       }
115       // LYZ1SUM:
116       else if(TString(list->GetName()).Contains("LYZ1SUM"))
117       {
118         AliFlowAnalysisWithLeeYangZeros* lyz1sum = new AliFlowAnalysisWithLeeYangZeros();
119         lyz1sum->SetFirstRun(kTRUE);
120         lyz1sum->SetUseSum(kTRUE);
121         lyz1sum->GetOutputHistograms(list);
122         lyz1sum->Finish();
123         directory->Add(list,kTRUE);
124         directory->Write(directory->GetName(),TObject::kSingleKey+TObject::kWriteDelete);
125       }
126       // LYZ2SUM:
127       else if(TString(list->GetName()).Contains("LYZ2SUM"))
128       {
129         AliFlowAnalysisWithLeeYangZeros* lyz2sum = new AliFlowAnalysisWithLeeYangZeros();
130         lyz2sum->SetFirstRun(kFALSE);
131         lyz2sum->SetUseSum(kTRUE);
132         lyz2sum->GetOutputHistograms(list);
133         lyz2sum->Finish();
134         directory->Add(list,kTRUE);
135         directory->Write(directory->GetName(),TObject::kSingleKey+TObject::kWriteDelete);
136       }
137       // LYZ1PROD:
138       else if(TString(list->GetName()).Contains("LYZ1PROD"))
139       {
140         AliFlowAnalysisWithLeeYangZeros* lyz1prod = new AliFlowAnalysisWithLeeYangZeros();
141         lyz1prod->SetFirstRun(kTRUE);
142         lyz1prod->SetUseSum(kFALSE);
143         lyz1prod->GetOutputHistograms(list);
144         lyz1prod->Finish();
145         directory->Add(list,kTRUE);
146         directory->Write(directory->GetName(),TObject::kSingleKey+TObject::kWriteDelete);
147       }
148       // LYZ2PROD:
149       else if(TString(list->GetName()).Contains("LYZ2PROD"))
150       {
151         AliFlowAnalysisWithLeeYangZeros* lyz2prod = new AliFlowAnalysisWithLeeYangZeros();
152         lyz2prod->SetFirstRun(kFALSE);
153         lyz2prod->SetUseSum(kFALSE);
154         lyz2prod->GetOutputHistograms(list);
155         lyz2prod->Finish();
156         directory->Add(list,kTRUE);
157         directory->Write(directory->GetName(),TObject::kSingleKey+TObject::kWriteDelete);
158       }
159       // LYZEP:
160       else if(TString(list->GetName()).Contains("LYZEP"))
161       {
162         AliFlowAnalysisWithLYZEventPlane* lyzep = new AliFlowAnalysisWithLYZEventPlane();
163         lyzep->GetOutputHistograms(list);
164         lyzep->Finish();
165         directory->Add(list,kTRUE);
166         directory->Write(directory->GetName(),TObject::kSingleKey+TObject::kWriteDelete);
167       }
168       // MH:
169       else if(TString(list->GetName()).Contains("MH"))
170       {
171         AliFlowAnalysisWithMixedHarmonics* mh = new AliFlowAnalysisWithMixedHarmonics();
172         mh->GetOutputHistograms(list);
173         mh->GetAnalysisSettings()->SetBinContent(1,(Int_t)bApplyCorrectionForNUA);
174         mh->Finish();
175         directory->Add(list,kTRUE);
176         directory->Write(directory->GetName(),TObject::kSingleKey+TObject::kWriteDelete);
177       }
178       // NL:
179       else if(TString(list->GetName()).Contains("NL"))
180       {
181         AliFlowAnalysisWithNestedLoops* nl = new AliFlowAnalysisWithNestedLoops();
182         nl->GetOutputHistograms(list);
183         nl->Finish();
184         directory->Add(list,kTRUE);
185         directory->Write(directory->GetName(),TObject::kSingleKey+TObject::kWriteDelete);
186       }
187     }//for icent
188   }//
189
190     // Close the final output file:
191     delete mergedFile;
192
193     // Giving the final names if neccessary:
194     if(!(mergedFileName == outputFileName))
195     {
196       TSystemFile *outputFileFinal = new TSystemFile(mergedFileName.Data(),".");
197       outputFileFinal->Rename(outputFileName.Data());
198       delete outputFileFinal;
199       TSystemFile *mergedFileFinal = new TSystemFile("mergedAnalysisResultsTemp.root",".");
200       mergedFileFinal->Rename(mergedFileName.Data());
201       delete mergedFileFinal;
202     } // end of if(!(mergedFileName == outputFileName))
203
204     cout<<endl;
205
206 } // end of void redoFinish(Int_t mode=mLocal)
207
208 void LoadLibraries()
209 {
210   //--------------------------------------
211   // Load the needed libraries most of them already loaded by aliroot
212   //--------------------------------------
213   //gSystem->Load("libTree");
214   gSystem->Load("libGeom");
215   gSystem->Load("libVMC");
216   gSystem->Load("libXMLIO");
217   gSystem->Load("libPhysics");
218
219   // for AliRoot
220   gSystem->Load("libANALYSIS");
221   gSystem->Load("libANALYSISalice");
222   gSystem->Load("libPWGflowBase");
223   gSystem->Load("libPWGflowTasks");
224 } // end of void LoadLibrariesRF(const libModes mode)
225
226