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