]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ANALYSIS/AliFileMerger.cxx
DQM flag reassigned due to new PMT's calibration
[u/mrichter/AliRoot.git] / ANALYSIS / AliFileMerger.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /* $Id$ */
17
18 //////////////////////////////////////////////////////////////////////////
19 //  marian.ivanov@cern.ch
20 //  Utilities for file merging.
21 //  Additional functionality on top of the standard TFileMerger:
22 //
23 //  1. Possibility to Set the reject/accept list.
24 //     1.a)  Only entries selected in accept list are merged. By default all entries are selected
25 //           use AddAccept 0 to specify your desired entry
26 //     1.b)  Entries selected in reject list are not merged. By default the reject list is empty.
27 //
28 //  2. syswatch.log is created diring mergin procedure. 
29 //     Memeory consumption - for reading and for merging can be monitored
30
31 //  RS: Changed merger to respect the structure of files being merged (directories, collections...)
32 //      Additional option: SetNoTrees (default false) to not merge any tree
33 //      The code mostly taken from root's hadd.cxx
34 /*
35   Usage:
36   // Libraries for all classes to be merged should be loaded before using the class
37   gSystem->Load("libANALYSIS");
38   gSystem->Load("libANALYSIScalib");
39   gSystem->Load("libTPCcalib"); 
40   TH1::AddDirectory(0);
41
42   //Example usage starting from the input data list in text file:
43   //
44   AliFileMerger merger;
45   merger.AddReject("esdFriend");
46   merger.IterTXT("calib.list","CalibObjects.root",kFALSE);
47   //
48
49 */
50 //////////////////////////////////////////////////////////////////////////
51  
52
53 #include <fstream>
54 #include <THashList.h>
55 #include <TChain.h>
56 #include <TKey.h>
57 #include <TH1.h>
58 #include <THStack.h>
59 #include "TSystem.h"
60 #include "TFile.h"
61 #include "TGrid.h"
62 #include "TGridResult.h"
63 #include "TObjString.h"
64 #include "TObjArray.h"
65 #include "TMethodCall.h"
66 #include "Riostream.h"
67 #include "AliSysInfo.h"
68 #include "AliFileMerger.h"
69
70 ClassImp(AliFileMerger)
71
72 ////////////////////////////////////////////////////////////////////////
73
74 AliFileMerger::AliFileMerger():
75   TNamed(),
76   fRejectMask(0),
77   fAcceptMask(0),
78   fNoTrees(kFALSE)
79 {
80   //
81   // Default constructor
82   //
83 }
84
85 //______________________________________________________________________
86
87 AliFileMerger::AliFileMerger(const char* name):
88   TNamed(name,name),
89   fRejectMask(0),
90   fAcceptMask(0),
91   fNoTrees(kFALSE)
92 {
93   //
94   // 
95   //
96 }
97
98
99 void AliFileMerger::IterAlien(const char* outputDir, const char* outputFileName, const char* pattern, Bool_t dontOverwrite){
100
101   //
102   // Merge the files coming out of the calibration job
103   // 
104   TString outputFile(outputFileName);
105   gSystem->ExpandPathName(outputFile);
106   TString command;
107   // looking for files to be merged in the output directory
108   command = Form("find %s/ *%s", outputDir, pattern);
109   printf("command: %s\n", command.Data());
110   TGrid::Connect("alien://");
111   TGridResult *res = gGrid->Command(command);
112   if (!res) return;
113   TIter nextmap(res);
114   TMap *map = 0;
115   // loop over the results
116   TList sourcelist;  
117   sourcelist.SetOwner(kTRUE);
118   //
119   while((map=(TMap*)nextmap())) {
120     // getting the turl
121     TObjString *objs = dynamic_cast<TObjString*>(map->GetValue("turl"));
122     if (!objs || !objs->GetString().Length()) {
123       // Nothing found - skip this output
124       delete res;
125       break;
126     } 
127     printf("looking for file %s\n",(objs->GetString()).Data());
128     AddFile(&sourcelist, (objs->GetString()).Data());;
129   }
130   printf("List of files to be merged:\n");
131   sourcelist.Print();
132   //  
133   TFile* target = TFile::Open( outputFile.Data(), (dontOverwrite ? "CREATE":"RECREATE") );
134   if (!target || target->IsZombie()) {
135     cerr << "Error opening target file (does " << outputFileName << " exist?)." << endl;
136     cerr << "Use force = kTRUE to re-creation of output file." << endl;
137     return;
138   }  
139   MergeRootfile( target, &sourcelist);
140   delete res;
141   delete target;
142 }
143
144 void AliFileMerger::IterTXT( const char * fileList,  const char* outputFileName, Bool_t dontOverwrite){
145   
146   // Merge the files indicated in the list - fileList
147   // ASCII file option example: 
148   // find `pwd`/ | grep AliESDfriends_v1.root > calib.list
149   
150   // Open the input stream
151   ifstream in;
152   in.open(fileList);
153   // Read the input list of files 
154   TString objfile;
155   Int_t counter=0;
156   TList sourcelist;  
157   sourcelist.SetOwner(kTRUE);
158   while(in.good()) {
159     in >> objfile;
160     if (!objfile.Contains(".root")) continue; // protection
161     gSystem->ExpandPathName(objfile);
162     printf("Add file:Counter\t%d\tMerging file %s\n",counter++,objfile.Data());
163     AddFile(&sourcelist, objfile.Data());
164   }
165   //
166   printf("List of files to be merged:\n");
167   sourcelist.Print();
168   //
169   TString outputFile(outputFileName);
170   gSystem->ExpandPathName(outputFile);
171   TFile* target = TFile::Open( outputFile.Data(), (dontOverwrite ? "CREATE":"RECREATE") );
172   if (!target || target->IsZombie()) {
173     cerr << "Error opening target file (does " << outputFileName << " exist?)." << endl;
174     cerr << "Use force = kTRUE to re-creation of output file." << endl;
175     return;
176   }
177
178   MergeRootfile( target, &sourcelist);
179   delete target;
180 }
181
182 void AliFileMerger::StoreResults(TObjArray * array, const char* outputFileName){
183   //
184   // Storing the results in one single file
185   //
186   TFile *f = new TFile(outputFileName,"recreate");
187   for (Int_t i=0; i<array->GetEntries(); i++){
188     TObject *object0 = array->At(i);
189     if (!object0) continue;
190     object0->Write();
191   }
192   f->Close();
193   delete f;
194 }
195
196
197 void AliFileMerger::StoreSeparateResults(TObjArray * array, const char* outputFileName){
198   //
199   // Store the results in separate files (one per object)
200   //
201   for (Int_t i=0; i<array->GetEntries(); i++){
202     TObject *object0 = array->At(i);
203     if (!object0) continue;
204     TFile *f = new TFile(Form("%s_%s.root",outputFileName,object0->GetName()),"recreate");
205     object0->Write();
206     f->Close();
207     delete f;
208   }
209 }
210
211 void AliFileMerger::Merge(TFile* fileIn, TObjArray * array){
212   //
213   // Merging procedure
214   //
215   if (!array) return;
216   static Int_t counter=-1;
217   counter++;
218   TObjArray *carray = new TObjArray;   //array of the objects inside current file
219   carray->SetOwner(kTRUE);
220   
221   // load all objects to  memory
222   
223   TList *farr = fileIn->GetListOfKeys();
224   if (!farr) { 
225     delete carray;
226     return;
227   }
228   for (Int_t ical=0; ical<farr->GetEntries(); ical++){
229     if (!farr->At(ical)) continue;
230     TString name(farr->At(ical)->GetName());
231     if (!IsAccepted(name)) continue;                        // skip not accepted entries
232     TObject *obj = fileIn->Get(name.Data());
233     if (obj) carray->AddLast(obj);
234     AliSysInfo::AddStamp(name.Data(),1,ical,counter);  
235   }
236   
237   if (carray->GetEntries()==0)  { 
238     delete carray;
239     return;
240   }
241   TMethodCall callEnv;
242   Int_t entries =carray->GetEntriesFast();
243   for (Int_t i=0; i<entries; i++){
244     
245     TObjArray *templist = new TObjArray(1);
246     templist->SetOwner(kFALSE);
247     TObject *currentObject = carray->At(i);
248     if (!currentObject) { 
249       delete templist;
250       continue;
251     }
252     printf("%s\n",currentObject->GetName());
253     callEnv.InitWithPrototype(currentObject->IsA(), "Merge", "TCollection*");
254     if (!callEnv.IsValid()) {
255       delete templist; 
256       continue;
257     }
258     TString oname=currentObject->GetName();
259     TObject *mergedObject = array->FindObject(currentObject->GetName());
260     if (!mergedObject) {
261       array->AddLast(currentObject);
262       carray->RemoveAt(i);
263       delete templist; 
264       continue;
265     }
266     templist->AddLast(currentObject);
267     callEnv.SetParam((Long_t) templist);
268     callEnv.Execute(mergedObject);
269     AliSysInfo::AddStamp(currentObject->GetName(),2,i,counter);  
270     delete templist;
271   }
272   carray->Delete();
273   delete carray;
274 }
275
276 Bool_t AliFileMerger::IsAccepted(TString name){
277   //
278   // Accept/reject logic
279   // name - name of the entry
280   //
281   //  if fAcceptMask specified   - entry has to be in list of selected
282   //  if fRejectMask speciefied  - entry with name speciief in the list are rejected 
283   //
284   Bool_t accept=kTRUE;
285   if (fAcceptMask){
286     //
287     accept=kFALSE;
288     for (Int_t iaccept=0; iaccept<fAcceptMask->GetEntries(); iaccept++){
289       if (name.Contains(fAcceptMask->At(iaccept)->GetName())) accept=kTRUE;   // entry was selected
290     }
291   }
292   if (!accept) return kFALSE;
293
294   if (fRejectMask){
295     //
296     for (Int_t ireject=0; ireject<fRejectMask->GetEntries(); ireject++){
297       if (name.Contains(fRejectMask->At(ireject)->GetName())) accept=kFALSE;   // entry was rejected
298     }
299   }
300   return accept;
301 }
302
303
304
305
306 void AliFileMerger::AddReject(const char *reject){
307   //
308   // add reject string to the list of entries to be rejected for merging
309   //
310   if (!fRejectMask) fRejectMask = new TObjArray;
311   fRejectMask->AddLast(new TObjString(reject));
312 }
313 void AliFileMerger::AddAccept(const char *accept){
314   //
315   // add reject string to the list of entries to be rejected for merging
316   //
317   if (!fAcceptMask) fAcceptMask = new TObjArray;
318   fAcceptMask->AddLast(new TObjString(accept));
319
320
321 }
322
323 //___________________________________________________________________________
324 int AliFileMerger::MergeRootfile( TDirectory *target, TList *sourcelist)
325 {
326   // Merge all objects in a directory
327   // modified version of root's hadd.cxx
328   Int_t counterF = -1;
329   int status = 0;
330   cout << "Target path: " << target->GetPath() << endl;
331   TString path( (char*)strstr( target->GetPath(), ":" ) );
332   path.Remove( 0, 2 );
333   //
334   // find 1st valid file
335   TDirectory *first_source = (TDirectory*)sourcelist->First();
336   //
337   Int_t nguess = sourcelist->GetSize()+1000;
338   THashList allNames(nguess);
339   ((THashList*)target->GetList())->Rehash(nguess);
340   ((THashList*)target->GetListOfKeys())->Rehash(nguess);
341   TList listH;
342   TString listHargs;
343   listHargs.Form("((TCollection*)0x%lx)", (ULong_t)&listH);
344   //
345   while(first_source) {
346     counterF++;
347     TDirectory *current_sourcedir = first_source->GetDirectory(path);
348     if (!current_sourcedir) {
349       first_source = (TDirectory*)sourcelist->After(first_source);
350       continue;
351     }
352     // loop over all keys in this directory
353     TChain *globChain = 0;
354     TIter nextkey( current_sourcedir->GetListOfKeys() );
355     TKey *key, *oldkey=0;
356     //gain time, do not add the objects in the list in memory
357     TH1::AddDirectory(kFALSE);
358     //
359     int counterK = 0;
360     //
361     while ( (key = (TKey*)nextkey())) {
362       if (current_sourcedir == target) break;
363       //
364       // check if we don't reject this name
365       TString nameK(key->GetName());
366       if (!IsAccepted(nameK)) {
367         if (!counterF) printf("Object %s is in rejection list, skipping...\n",nameK.Data());
368         continue;
369       }
370       //
371       //keep only the highest cycle number for each key
372       if (oldkey && !strcmp(oldkey->GetName(),key->GetName())) continue;
373       if (!strcmp(key->GetClassName(),"TProcessID")) {key->ReadObj(); continue;}
374       if (allNames.FindObject(key->GetName())) continue;
375       TClass *cl = TClass::GetClass(key->GetClassName());
376       if (!cl || !cl->InheritsFrom(TObject::Class())) {
377         cout << "Cannot merge object type, name: "
378              << key->GetName() << " title: " << key->GetTitle() << endl;
379         continue;
380       }
381       allNames.Add(new TObjString(key->GetName()));
382       AliSysInfo::AddStamp(nameK.Data(),1,counterK++,counterF-1); 
383       // read object from first source file
384       //current_sourcedir->cd();
385       TObject *obj = key->ReadObj();
386       //printf("keyname=%s, obj=%x\n",key->GetName(),obj);
387       
388       if ( obj->IsA()->InheritsFrom( TTree::Class() ) ) {
389         
390         // loop over all source files create a chain of Trees "globChain"
391         if (!fNoTrees) { // 
392           TString obj_name;
393           if (path.Length()) {
394             obj_name = path + "/" + obj->GetName();
395           } else {
396             obj_name = obj->GetName();
397           }
398           globChain = new TChain(obj_name);
399           globChain->Add(first_source->GetName());
400           TFile *nextsource = (TFile*)sourcelist->After( first_source );
401           while ( nextsource ) {
402             //do not add to the list a file that does not contain this Tree
403             TFile *curf = TFile::Open(nextsource->GetName());
404             if (curf) {
405               Bool_t mustAdd = kFALSE;
406               if (curf->FindKey(obj_name)) {
407                 mustAdd = kTRUE;
408               } else {
409                 //we could be more clever here. No need to import the object
410                 //we are missing a function in TDirectory
411                 TObject *aobj = curf->Get(obj_name);
412                 if (aobj) { mustAdd = kTRUE; delete aobj;}
413               }
414               if (mustAdd) {
415                 globChain->Add(nextsource->GetName());
416               }
417             }
418             delete curf;
419             nextsource = (TFile*)sourcelist->After( nextsource );
420           }
421         }
422       } else if ( obj->IsA()->InheritsFrom( TDirectory::Class() ) ) {
423         // it's a subdirectory
424         
425         cout << "Found subdirectory " << obj->GetName() << endl;
426         // create a new subdir of same name and title in the target file
427         target->cd();
428         TDirectory *newdir = target->mkdir( obj->GetName(), obj->GetTitle() );
429         
430         // newdir is now the starting point of another round of merging
431         // newdir still knows its depth within the target file via
432         // GetPath(), so we can still figure out where we are in the recursion
433         status = MergeRootfile( newdir, sourcelist);
434         if (status) return status;
435         
436       } else if ( obj->InheritsFrom(TObject::Class())
437                   && obj->IsA()->GetMethodWithPrototype("Merge", "TCollection*") ) {
438         // object implements Merge(TCollection*)
439         
440         // loop over all source files and merge same-name object
441         TFile *nextsource = (TFile*)sourcelist->After( first_source );
442         while ( nextsource ) {
443           // make sure we are at the correct directory level by cd'ing to path
444           TDirectory *ndir = nextsource->GetDirectory(path);
445           if (ndir) {
446             ndir->cd();
447             TKey *key2 = (TKey*)gDirectory->GetListOfKeys()->FindObject(key->GetName());
448             if (key2) {
449               TObject *hobj = key2->ReadObj();
450               hobj->ResetBit(kMustCleanup);
451               listH.Add(hobj);
452               Int_t error = 0;
453               obj->Execute("Merge", listHargs.Data(), &error);
454               if (error) {
455                 cerr << "Error calling Merge() on " << obj->GetName()
456                      << " with the corresponding object in " << nextsource->GetName() << endl;
457               }
458               listH.Delete();
459             }
460           }
461           nextsource = (TFile*)sourcelist->After( nextsource );
462         }
463       } else if ( obj->IsA()->InheritsFrom( THStack::Class() ) ) {
464         THStack *hstack1 = (THStack*) obj;
465         TList* l = new TList();
466         
467         // loop over all source files and merge the histos of the
468         // corresponding THStacks with the one pointed to by "hstack1"
469         TFile *nextsource = (TFile*)sourcelist->After( first_source );
470         while ( nextsource ) {
471           // make sure we are at the correct directory level by cd'ing to path
472           TDirectory *ndir = nextsource->GetDirectory(path);
473           if (ndir) {
474             ndir->cd();
475             TKey *key2 = (TKey*)gDirectory->GetListOfKeys()->FindObject(hstack1->GetName());
476             if (key2) {
477               THStack *hstack2 = (THStack*) key2->ReadObj();
478               l->Add(hstack2->GetHists()->Clone());
479               delete hstack2;
480             }
481           }
482           
483           nextsource = (TFile*)sourcelist->After( nextsource );
484         }
485         hstack1->GetHists()->Merge(l);
486         l->Delete();
487       } else {
488         // object is of no type that we can merge
489         cout << "Cannot merge object type, name: "
490              << obj->GetName() << " title: " << obj->GetTitle() << endl;
491         
492         // loop over all source files and write similar objects directly to the output file
493         TFile *nextsource = (TFile*)sourcelist->After( first_source );
494         while ( nextsource ) {
495           // make sure we are at the correct directory level by cd'ing to path
496           TDirectory *ndir = nextsource->GetDirectory(path);
497           if (ndir) {
498             ndir->cd();
499             TKey *key2 = (TKey*)gDirectory->GetListOfKeys()->FindObject(key->GetName());
500             if (key2) {
501               TObject *nobj = key2->ReadObj();
502               nobj->ResetBit(kMustCleanup);
503               int nbytes1 = target->WriteTObject(nobj, key2->GetName(), "SingleKey" );
504               if (nbytes1 <= 0) status = -1;
505               delete nobj;
506             }
507           }
508           nextsource = (TFile*)sourcelist->After( nextsource );
509         }
510       }
511       
512       // now write the merged histogram (which is "in" obj) to the target file
513       // note that this will just store obj in the current directory level,
514       // which is not persistent until the complete directory itself is stored
515       // by "target->Write()" below
516       target->cd();
517       
518       //!!if the object is a tree, it is stored in globChain...
519       if(obj->IsA()->InheritsFrom( TDirectory::Class() )) {
520         //printf("cas d'une directory\n");
521       } else if(obj->IsA()->InheritsFrom( TTree::Class() )) {
522         if (!fNoTrees) {
523           globChain->ls();
524           globChain->Merge(target->GetFile(),0,"keep fast");
525           delete globChain;
526         }
527       } else {
528         int nbytes2 = obj->Write( key->GetName(), TObject::kSingleKey );
529         if (nbytes2 <= 0) status = -1;
530       }
531       oldkey = key;
532       delete obj;
533     } // while ( ( TKey *key = (TKey*)nextkey() ) )
534     first_source = (TDirectory*)sourcelist->After(first_source);
535   }
536   // save modifications to target file
537   target->SaveSelf(kTRUE);
538   //
539   return status;
540 }
541
542 //___________________________________________________________________________
543 int AliFileMerger::AddFile(TList* sourcelist, std::string entry)
544 {
545   // add a new file to the list of files
546   //  static int count(0);
547   if( entry.empty() ) return 0;
548   size_t j =entry.find_first_not_of(' ');
549   if( j==std::string::npos ) return 0;
550   entry = entry.substr(j);
551   if( entry.substr(0,1)=="@") {
552     std::ifstream indirect_file(entry.substr(1).c_str() );
553     if( ! indirect_file.is_open() ) {
554       std::cerr<< "Could not open indirect file " << entry.substr(1) << std::endl;
555       return 1;
556     }
557     while( indirect_file ){
558       std::string line;
559       std::getline(indirect_file, line);
560       if( AddFile(sourcelist, line)!=0 )return 1;;
561     }
562     return 0;
563   }
564   //  cout << "Source file " << (++count) << ": " << entry << endl;
565   
566   TFile* source = TFile::Open( entry.c_str());
567   if( source==0 ) {
568     cout << "Failed to open " << entry << " will skip" << endl;
569     return 0;
570   }
571   sourcelist->Add(source);
572   return 0;
573 }