]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG3/muon/mergeGridFiles.C
dc207fb1b989affeee61a3bdd7e6acc3bbde6a7d
[u/mrichter/AliRoot.git] / PWG3 / muon / mergeGridFiles.C
1 #if !defined(__CINT__) || defined(__MAKECINT__)
2
3 #include <Riostream.h>
4
5 // ROOT includes
6 #include "TString.h"
7 #include "TObjArray.h"
8 #include "TObjString.h"
9 #include "TFileMerger.h"
10 #include "TGrid.h"
11 #include "TList.h"
12 #include "TSystem.h"
13 #include "TGridResult.h"
14 #include "TMap.h"
15 #include "TFile.h"
16 #include "TROOT.h"
17 #endif
18
19 Bool_t AddFileList(TString, TString&, TString);
20 void ReadListFromFile(TString, TString&, TString);
21 Int_t GetLastStage(TGridResult*);
22
23 void mergeGridFiles(TString outFilename, TString inFileList, TString addPrefix = "alien://", Int_t nFilesPerStep = 100, Bool_t copyLocal = kTRUE, TString dirsToMerge = "")
24 {
25   TString fileList = "";
26   ReadListFromFile(inFileList, fileList, addPrefix);
27
28   if ( fileList.IsNull() ) {
29     printf("List of merging files is null: Nothing done!\n");
30     return;
31   }
32
33   TString logFilename = "toBeMerged.txt";
34
35   if ( ! gGrid ) copyLocal = kFALSE;
36
37   TString currList = fileList;
38   TString currOutput = "", mergedFiles = "";
39   TFileMerger* fileMerger = 0x0;
40   TObjArray* dirsToMergeArray = 0x0;
41   TList* sourceList = new TList();
42   if ( ! dirsToMerge.IsNull() ) {
43     dirsToMergeArray = dirsToMerge.Tokenize(" ");
44     dirsToMergeArray->SetOwner();
45   }
46
47   Bool_t showProgressBar = ! gROOT->IsBatch();
48
49   for ( Int_t istep = 0; istep < 100; istep++ ) {
50     TObjArray* array = currList.Tokenize(" ");
51     array->SetOwner();
52     currList = "";
53     Int_t nFiles = array->GetEntries();
54     Int_t subStep = -1;
55     for (Int_t ifile = 0; ifile < nFiles; ifile++ ) {
56       if ( ! fileMerger ) fileMerger = new TFileMerger(copyLocal);
57       TString currFilename = array->At(ifile)->GetName();
58       Bool_t isFileAdded = fileMerger->AddFile(currFilename.Data(), showProgressBar);
59       if ( ! isFileAdded ) continue;
60       //printf("New file %s\n", gROOT->GetListOfFiles()->Last()->GetName()); // REMEMBER TO CUT
61       if ( dirsToMergeArray ) {
62         if ( ! sourceList ) sourceList = new TList();
63         sourceList->Add(gROOT->GetListOfFiles()->Last());
64       }
65       
66       Int_t nFilesToMerge = fileMerger->GetMergeList()->GetEntries();
67       if ( nFilesToMerge % nFilesPerStep != 0 && ifile < nFiles - 1 ) 
68         continue;
69       // The following part is executed only at the end of each step
70       currOutput = outFilename;
71       if ( nFiles > nFilesPerStep ) {
72         subStep++;
73         currOutput.ReplaceAll(".root",Form("_%i_%i.root", istep, subStep));
74         AddFileList(currOutput, currList, "");
75       }
76       if ( dirsToMergeArray ) {
77         TFile* outFile = new TFile(currOutput.Data(), "recreate");
78         for ( Int_t idir=0; idir<dirsToMergeArray->GetEntries(); idir++ ) {
79           outFile->cd();
80           TDirectory* currDir = outFile->mkdir(dirsToMergeArray->At(idir)->GetName());
81           fileMerger->MergeRecursive(currDir, sourceList);
82         }
83         outFile->Close();
84       }
85       else {
86         fileMerger->OutputFile(currOutput.Data());
87         fileMerger->Merge();
88       }
89       printf("\nMerged in %s:\n", currOutput.Data());
90       mergedFiles = "";
91       for ( Int_t ientry=0; ientry<nFilesToMerge; ientry++ )
92         mergedFiles += Form("%s ", fileMerger->GetMergeList()->At(ientry)->GetName());
93       printf("%s\n\n", mergedFiles.Data());
94
95       // remove merged files
96       if ( istep > 0 )
97         gSystem->Exec(Form("rm %s", mergedFiles.Data()));
98
99       delete fileMerger;
100       fileMerger = 0x0;
101       if ( dirsToMergeArray ) {
102         delete sourceList;
103         sourceList = 0x0;
104       }
105
106       // Write log file to keep trace of files to be merged
107       ofstream logFile(logFilename.Data());
108       TString logString = "";
109       for ( Int_t jfile = ifile + 1; jfile < nFiles; jfile++ ) {
110         logString += Form("%s ", array->At(jfile)->GetName());
111       }
112       logString.Append(currList.Data());
113       logString.ReplaceAll(" ", "\n");
114       logFile << logString.Data() << endl;;
115       logFile.close();
116     } // loop on files
117
118
119     delete array;
120     printf("Step %i completed!\n", istep);
121
122     if ( nFiles <= nFilesPerStep ) break;
123   } // loop on steps
124
125   gSystem->Exec(Form("rm %s", logFilename.Data()));
126 }
127
128
129 //___________________________________________________
130 Bool_t AddFileList(TString filename, TString& fileList, TString addPrefix)
131 {
132   if ( filename.IsNull() || ! filename.Contains(".root") ) return kFALSE;
133
134   if ( ! addPrefix.IsNull() && ! filename.Contains(addPrefix.Data()) )
135     filename.Prepend(addPrefix.Data());
136
137   if ( filename.Contains("alien://") && ! gGrid )
138     TGrid::Connect("alien://");
139   
140   if ( ! fileList.IsNull() )
141     fileList.Append(" ");
142   fileList.Append(filename.Data());
143
144   return kTRUE;
145 }
146
147 //___________________________________________________
148 void ReadListFromFile(TString filename, TString& fileList, TString addPrefix)
149 {
150   ifstream inFile(filename.Data());
151   TString currLine = "";
152   if ( inFile.is_open() ) {
153     while ( ! inFile.eof() ) {
154       currLine.ReadLine(inFile,kTRUE); // Read line
155       AddFileList(currLine, fileList, addPrefix);
156     }
157     inFile.close();
158   }
159 }
160
161
162 //___________________________________________________
163 void completeProd(TString runListName="runList.txt", TString prodDir = "", TString baseDir="/alice/data/2010/LHC10h", TString outTaskFilename="QAresults.root", Int_t nFilesPerStep = 50, TString dirsToMerge = "MUON_QA MUON.TriggerEfficiencyMap", Bool_t mergeFast = kFALSE, Bool_t overwriteExisting = kFALSE)
164 {
165   TString outFilename = "completeFileList.txt";
166
167   // Get run list from file
168   ifstream inFile(runListName.Data());
169   TObjArray runList;
170   runList.SetOwner();
171   TString currRun;
172   Int_t minRun = 99999999;
173   Int_t maxRun = -1;
174   if (inFile.is_open()) {
175     while (! inFile.eof() ) {
176       currRun.ReadLine(inFile,kTRUE); // Read line
177       if ( currRun.IsNull() || ! currRun.IsDigit() ) continue;
178       Int_t currRunInt = currRun.Atoi();
179       minRun = TMath::Min(currRunInt, minRun);
180       maxRun = TMath::Max(currRunInt, maxRun);
181       runList.Add(new TObjString(Form("%d", currRunInt)));
182     }
183     inFile.close();
184   }
185
186   outFilename.ReplaceAll(".txt", Form("_%i_%i.txt", minRun, maxRun));
187
188   //TString filePattern[3] = {"", "Stage*/","*/"};
189   TString filePattern[2] = {"","*/"};
190
191   ofstream outFile(outFilename.Data());
192
193   // if ( outTaskFilename.Contains("QAresults.root") ) {
194   const Int_t kNlibs = 5; // 1
195   //TString loadLibs[kNlibs] = {"libPWG3base.so"};
196   //TString loadLibs[kNlibs] = {"libANALYSIS.so", "libANALYSISalice.so", "libTENDER.so", "libPWG1.so", "libPWG3base.so"};
197   TString loadLibs[kNlibs] = {"libANALYSIS.so", "libOADB.so", "libANALYSISalice.so", "libCORRFW.so", "libPWG3base.so"};
198   for ( Int_t ilib=0; ilib<kNlibs; ilib++ ) {
199     Int_t exitVal = gSystem->Load(loadLibs[ilib].Data());
200     if ( exitVal < 0 ) {
201       printf("Please run with aliroot if you're merging QA objects!\n");
202       return;
203     }
204   }
205   //}
206
207   if ( ! gGrid )
208     TGrid::Connect("alien://");
209
210   baseDir.ReplaceAll("alien://","");
211
212   TMap* map = 0x0;
213   TString stageName = "";
214   TString runsWithoutOut = "";
215   for ( Int_t irun=0; irun<runList.GetEntries(); irun++ ) {
216     TString currRunString = ((TObjString*)runList.At(irun))->GetString();
217
218     TString localOut = outTaskFilename;
219     localOut.ReplaceAll(".root", Form("_%s.root", currRunString.Data()));
220     if ( ! gSystem->AccessPathName(localOut.Data()) ) {
221       if ( overwriteExisting )
222         printf("Overwriting existing file %s\n", localOut.Data());
223       else {
224         printf("Warning: merged file %s already exist: do not overwrite\n", localOut.Data());
225         outFile << gSystem->pwd() << "/" << localOut.Data() << endl;
226         continue;
227       }
228     }
229
230     TString tmpFilename = Form("/tmp/mergeListRun%s.txt", currRunString.Data());
231     ofstream tmpFile(tmpFilename.Data());
232     TString mergeFilename = "";
233
234     Int_t nPatterns = ( mergeFast ) ? 1 : 2;
235     
236     for ( Int_t ipattern=0; ipattern<nPatterns; ipattern++ ) {
237       TString command = ( prodDir.Contains("private") ) ? Form("find %s/ *%s/%s%s", baseDir.Data(), currRunString.Data(), filePattern[ipattern].Data(), outTaskFilename.Data()) : Form("find %s/*%s %s/%s%s", baseDir.Data(), currRunString.Data(), prodDir.Data(), filePattern[ipattern].Data(), outTaskFilename.Data());
238
239       printf("%s\n", command.Data());
240
241       TGridResult* res = gGrid->Command(command);
242
243       if ( ! res || res->GetEntries() == 0 ) continue;
244
245       Int_t mergeStage = ( ipattern == 1 ) ? GetLastStage(res)  : -1;
246       stageName = Form("Stage_%i", mergeStage);
247
248       TIter nextmap(res);
249       while ( ( map = (TMap*)nextmap() ) ) {
250         // Loop 'find' results and get next LFN
251         TObjString *objs = dynamic_cast<TObjString*>(map->GetValue("turl"));
252         if (!objs || !objs->GetString().Length()) 
253           continue;
254
255         mergeFilename = objs->GetString();
256
257         if ( mergeStage > 0 && ! mergeFilename.Contains(stageName.Data()) ) continue;
258
259         tmpFile << mergeFilename.Data() << endl;
260       } // loop on grid lfns
261
262       delete res;
263
264       tmpFile.close();
265
266       // If the merging of specific directories is required
267       // run merging also if there is 1 file
268       // so that we get rid of other sub-directories
269       if ( ipattern == 1 || ! dirsToMerge.IsNull() ) {
270         mergeFilename = outTaskFilename;
271         mergeFilename.ReplaceAll(".root", Form("_%s.root", currRunString.Data()));
272         mergeGridFiles(mergeFilename, tmpFilename, "alien://", nFilesPerStep, kTRUE, dirsToMerge);
273       }
274
275       if ( ! mergeFilename.Contains("alien://") )
276         outFile << gSystem->pwd() << "/";
277       outFile << mergeFilename.Data() << endl;
278       gSystem->Exec(Form("rm %s", tmpFilename.Data()));
279       break;
280     } // loop on pattern
281     if ( mergeFilename.IsNull() ) runsWithoutOut += currRunString + " ";
282   } // loop on runs
283   if ( ! runsWithoutOut.IsNull() )
284     printf("\nNo output found in runs\n%s\n",runsWithoutOut.Data());
285   printf("\nOutput written in:\n%s\n", outFilename.Data());
286 }
287
288
289 //___________________________________________________
290 Int_t GetLastStage(TGridResult* res)
291 {
292   Int_t lastStage = 0;
293
294   TMap* map = 0x0;
295   TIter nextmap(res);
296   TString filename = "", currToken = "";
297   while ( ( map = (TMap*)nextmap() ) ) {
298     // Loop 'find' results and get next LFN
299     TObjString *objs = dynamic_cast<TObjString*>(map->GetValue("turl"));
300     if (!objs || !objs->GetString().Length()) 
301       continue;
302
303     filename = objs->GetString();
304     
305     if ( ! filename.Contains("Stage_") ) continue;
306
307     TObjArray* array = filename.Tokenize("/");
308     array->SetOwner();
309     for ( Int_t ientry=0; ientry<array->GetEntries(); ientry++ ) {
310       currToken = array->At(ientry)->GetName();
311       if ( currToken.Contains("Stage_") ) {
312         currToken.Remove(0,6);
313         Int_t currStage = currToken.Atoi();
314         lastStage = TMath::Max(currStage, lastStage);
315         break;
316       }
317     }
318     delete array;
319   } // loop on grid lfns
320
321   return lastStage;
322 }