]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGPP/MUON/lite/mergeGridFiles.C
Better use of TFileMerger to remove intermediate files. Add pass name in final output...
[u/mrichter/AliRoot.git] / PWGPP / MUON / lite / 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
40   Bool_t showProgressBar = ! gROOT->IsBatch();
41
42   for ( Int_t istep = 0; istep < 100; istep++ ) {
43     TFileMerger fileMerger(copyLocal);
44     Int_t mergeType = ( TFileMerger::kRegular | TFileMerger::kAll );
45     if ( ! dirsToMerge.IsNull() ) {
46       fileMerger.AddObjectNames(dirsToMerge.Data());
47       mergeType |= TFileMerger::kOnlyListed;
48     }
49     TObjArray* array = currList.Tokenize(" ");
50     currList = "";
51     Int_t nFiles = array->GetEntries();
52     Int_t subStep = -1;
53     for (Int_t ifile = 0; ifile < nFiles; ifile++ ) {
54       TString currFilename = array->At(ifile)->GetName();
55       Bool_t isFileAdded = fileMerger.AddFile(currFilename.Data(), showProgressBar);
56       if ( ! isFileAdded ) continue;
57       
58       Int_t nFilesToMerge = fileMerger.GetMergeList()->GetEntries();
59       if ( nFilesToMerge % nFilesPerStep != 0 && ifile < nFiles - 1 ) 
60         continue;
61       // The following part is executed only at the end of each step
62       currOutput = outFilename;
63       if ( nFiles > nFilesPerStep ) {
64         subStep++;
65         currOutput.ReplaceAll(".root",Form("_%i_%i.root", istep, subStep));
66         AddFileList(currOutput, currList, "");
67       }
68       fileMerger.OutputFile(currOutput.Data(),kTRUE,1); // needed when merging single files for specific directories
69       fileMerger.PartialMerge(mergeType);
70       printf("\nMerged in %s:\n", currOutput.Data());
71       mergedFiles = "";
72       for ( Int_t ientry=0; ientry<nFilesToMerge; ientry++ )
73         mergedFiles += Form("%s ", fileMerger.GetMergeList()->At(ientry)->GetName());
74       printf("%s\n\n", mergedFiles.Data());
75
76       // remove merged files
77       if ( istep > 0 )
78         gSystem->Exec(Form("rm %s", mergedFiles.Data()));
79
80       // Write log file to keep trace of files to be merged
81       ofstream logFile(logFilename.Data());
82       TString logString = "";
83       for ( Int_t jfile = ifile + 1; jfile < nFiles; jfile++ ) {
84         logString += Form("%s ", array->At(jfile)->GetName());
85       }
86       logString.Append(currList.Data());
87       logString.ReplaceAll(" ", "\n");
88       logFile << logString.Data() << endl;;
89       logFile.close();
90     } // loop on files
91
92     delete array;
93     printf("Step %i completed!\n", istep);
94
95     if ( nFiles <= nFilesPerStep ) break;
96   } // loop on steps
97
98   gSystem->Exec(Form("rm %s", logFilename.Data()));
99 }
100
101
102 //___________________________________________________
103 Bool_t AddFileList(TString filename, TString& fileList, TString addPrefix)
104 {
105   if ( filename.IsNull() || ! filename.Contains(".root") ) return kFALSE;
106
107   if ( ! addPrefix.IsNull() && ! filename.Contains(addPrefix.Data()) )
108     filename.Prepend(addPrefix.Data());
109
110   if ( filename.Contains("alien://") && ! gGrid )
111     TGrid::Connect("alien://");
112   
113   if ( ! fileList.IsNull() )
114     fileList.Append(" ");
115   fileList.Append(filename.Data());
116
117   return kTRUE;
118 }
119
120 //___________________________________________________
121 void ReadListFromFile(TString filename, TString& fileList, TString addPrefix)
122 {
123   ifstream inFile(filename.Data());
124   TString currLine = "";
125   if ( inFile.is_open() ) {
126     while ( ! inFile.eof() ) {
127       currLine.ReadLine(inFile,kTRUE); // Read line
128       AddFileList(currLine, fileList, addPrefix);
129     }
130     inFile.close();
131   }
132 }
133
134
135 //___________________________________________________
136 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 MTR_ChamberEffMap", Bool_t mergeFast = kFALSE, Bool_t overwriteExisting = kFALSE)
137 {
138   TString outFilename = "completeFileList.txt";
139
140   // Get run list from file
141   ifstream inFile(runListName.Data());
142   TObjArray runList;
143   runList.SetOwner();
144   TString currRun;
145   Int_t minRun = 99999999;
146   Int_t maxRun = -1;
147   if (inFile.is_open()) {
148     while (! inFile.eof() ) {
149       currRun.ReadLine(inFile,kTRUE); // Read line
150       if ( currRun.IsNull() || ! currRun.IsDigit() ) continue;
151       Int_t currRunInt = currRun.Atoi();
152       minRun = TMath::Min(currRunInt, minRun);
153       maxRun = TMath::Max(currRunInt, maxRun);
154       runList.Add(new TObjString(Form("%d", currRunInt)));
155     }
156     inFile.close();
157   }
158
159   outFilename.ReplaceAll(".txt", Form("_%i_%i.txt", minRun, maxRun));
160
161   //TString filePattern[3] = {"", "Stage*/","*/"};
162   TString filePattern[2] = {"","*/"};
163
164   ofstream outFile(outFilename.Data());
165
166   const Int_t kNlibs = 5; // 1
167   TString loadLibs[kNlibs] = {"libANALYSIS.so", "libOADB.so", "libANALYSISalice.so", "libCORRFW.so", "libPWGmuon.so"};
168   for ( Int_t ilib=0; ilib<kNlibs; ilib++ ) {
169     Int_t exitVal = gSystem->Load(loadLibs[ilib].Data());
170     if ( exitVal < 0 ) {
171       printf("Please run with aliroot if you're merging QA objects!\n");
172       return;
173     }
174   }
175
176   if ( ! gGrid )
177     TGrid::Connect("alien://");
178
179   baseDir.ReplaceAll("alien://","");
180
181   TMap* map = 0x0;
182   TString stageName = "";
183   TString runsWithoutOut = "";
184   for ( Int_t irun=0; irun<runList.GetEntries(); irun++ ) {
185     TString currRunString = ((TObjString*)runList.At(irun))->GetString();
186
187     TString localOut = outTaskFilename;
188     localOut.ReplaceAll(".root", Form("_%s.root", currRunString.Data()));
189     if ( ! gSystem->AccessPathName(localOut.Data()) ) {
190       if ( overwriteExisting )
191         printf("Overwriting existing file %s\n", localOut.Data());
192       else {
193         printf("Warning: merged file %s already exist: do not overwrite\n", localOut.Data());
194         outFile << gSystem->pwd() << "/" << localOut.Data() << endl;
195         continue;
196       }
197     }
198
199     TString tmpFilename = Form("/tmp/mergeListRun%s.txt", currRunString.Data());
200     ofstream tmpFile(tmpFilename.Data());
201     TString mergeFilename = "";
202
203     Int_t nPatterns = ( mergeFast ) ? 1 : 2;
204     
205     for ( Int_t ipattern=0; ipattern<nPatterns; ipattern++ ) {
206       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());
207
208       printf("%s\n", command.Data());
209
210       TGridResult* res = gGrid->Command(command);
211
212       if ( ! res || res->GetEntries() == 0 ) continue;
213
214       Int_t mergeStage = ( ipattern == 1 ) ? GetLastStage(res)  : -1;
215       stageName = Form("Stage_%i", mergeStage);
216
217       TIter nextmap(res);
218       while ( ( map = (TMap*)nextmap() ) ) {
219         // Loop 'find' results and get next LFN
220         TObjString *objs = dynamic_cast<TObjString*>(map->GetValue("turl"));
221         if (!objs || !objs->GetString().Length()) 
222           continue;
223
224         mergeFilename = objs->GetString();
225
226         if ( mergeStage > 0 && ! mergeFilename.Contains(stageName.Data()) ) continue;
227         if ( mergeFilename.Contains(".resubmit") ) continue; // Avoid double counting for runs where a resubmission was performed
228
229         tmpFile << mergeFilename.Data() << endl;
230       } // loop on grid lfns
231
232       delete res;
233
234       tmpFile.close();
235
236       // If the merging of specific directories is required
237       // run merging also if there is 1 file
238       // so that we get rid of other sub-directories
239       if ( ipattern == 1 || ! dirsToMerge.IsNull() ) {
240         mergeFilename = outTaskFilename;
241         mergeFilename.ReplaceAll(".root", Form("_%s.root", currRunString.Data()));
242         mergeGridFiles(mergeFilename, tmpFilename, "alien://", nFilesPerStep, kTRUE, dirsToMerge);
243       }
244
245       if ( ! mergeFilename.Contains("alien://") )
246         outFile << gSystem->pwd() << "/";
247       outFile << mergeFilename.Data() << endl;
248       gSystem->Exec(Form("rm %s", tmpFilename.Data()));
249       break;
250     } // loop on pattern
251     if ( mergeFilename.IsNull() ) runsWithoutOut += currRunString + " ";
252   } // loop on runs
253   if ( ! runsWithoutOut.IsNull() )
254     printf("\nNo output found in runs\n%s\n",runsWithoutOut.Data());
255   printf("\nOutput written in:\n%s\n", outFilename.Data());
256 }
257
258
259 //___________________________________________________
260 Int_t GetLastStage(TGridResult* res)
261 {
262   Int_t lastStage = 0;
263
264   TMap* map = 0x0;
265   TIter nextmap(res);
266   TString filename = "", currToken = "";
267   while ( ( map = (TMap*)nextmap() ) ) {
268     // Loop 'find' results and get next LFN
269     TObjString *objs = dynamic_cast<TObjString*>(map->GetValue("turl"));
270     if (!objs || !objs->GetString().Length()) 
271       continue;
272
273     filename = objs->GetString();
274     
275     if ( ! filename.Contains("Stage_") ) continue;
276
277     TObjArray* array = filename.Tokenize("/");
278     array->SetOwner();
279     for ( Int_t ientry=0; ientry<array->GetEntries(); ientry++ ) {
280       currToken = array->At(ientry)->GetName();
281       if ( currToken.Contains("Stage_") ) {
282         currToken.Remove(0,6);
283         Int_t currStage = currToken.Atoi();
284         lastStage = TMath::Max(currStage, lastStage);
285         break;
286       }
287     }
288     delete array;
289   } // loop on grid lfns
290
291   return lastStage;
292 }