]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ANALYSIS/AliFileMerger.cxx
Forgot to increase class version
[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
32 /*
33   Usage:
34   // Libraries for all classes to be merged should be loaded before using the class
35   gSystem->Load("libANALYSIS");
36   gSystem->Load("libANALYSIScalib");
37   gSystem->Load("libTPCcalib"); 
38   TH1::AddDirectory(0);
39
40   //Example usage starting from the input data list in text file:
41   //
42   AliFileMerger merger;
43   merger.AddReject("esdFriend");
44   merger.IterTXT("calib.list","CalibObjects.root",kFALSE);
45   //
46
47 */
48 //////////////////////////////////////////////////////////////////////////
49  
50
51 #include <fstream>
52 #include "TSystem.h"
53 #include "TFile.h"
54 #include "TGrid.h"
55 #include "TGridResult.h"
56 #include "TObjString.h"
57 #include "TObjArray.h"
58 #include "TMethodCall.h"
59
60 #include "AliSysInfo.h"
61 #include "AliFileMerger.h"
62
63 ClassImp(AliFileMerger)
64
65 ////////////////////////////////////////////////////////////////////////
66
67 AliFileMerger::AliFileMerger():
68   TNamed(),
69   fRejectMask(0),
70   fAcceptMask(0)
71 {
72   //
73   // Default constructor
74   //
75 }
76
77 //______________________________________________________________________
78
79 AliFileMerger::AliFileMerger(const char* name):
80   TNamed(name,name),
81   fRejectMask(0),
82   fAcceptMask(0)
83 {
84   //
85   // 
86   //
87 }
88
89
90 void AliFileMerger::IterAlien(const char* outputDir, const char* outputFileName, const char* pattern){
91
92   //
93   // Merge the files coming out of the calibration job
94   // 
95   
96   TObjArray * mergeArray= new TObjArray;
97   
98   TString outputFile(outputFileName);
99   TString command;
100   // looking for files to be merged in the output directory
101   command = Form("find %s/ *%s", outputDir, pattern);
102   printf("command: %s\n", command.Data());
103   TGrid::Connect("alien://");
104   TGridResult *res = gGrid->Command(command);
105   if (!res) return;
106   TIter nextmap(res);
107   TMap *map = 0;
108   // loop over the results
109   while((map=(TMap*)nextmap())) {
110     // getting the turl
111     TObjString *objs = dynamic_cast<TObjString*>(map->GetValue("turl"));
112     if (!objs || !objs->GetString().Length()) {
113       // Nothing found - skip this output
114       delete res;
115       break;
116     } 
117     printf("looking for file %s\n",(objs->GetString()).Data());
118     TFile* currentFile=TFile::Open((objs->GetString()).Data());
119     Merge(currentFile, mergeArray);
120   }
121   Bool_t separate = kFALSE;
122   if (separate) 
123     StoreSeparateResults(mergeArray,outputFileName);
124   else
125     StoreResults(mergeArray,outputFileName);
126   delete mergeArray;
127 }
128
129
130
131 void AliFileMerger::IterTXT( const char * fileList,  const char* outputFileName, Bool_t separate){
132   
133   // Merge the files indicated in the list - fileList
134   // ASCII file option example: 
135   // find `pwd`/ | grep AliESDfriends_v1.root > calib.list
136   
137   TObjArray * mergeArray= new TObjArray;
138   
139   // Open the input stream
140   
141   ifstream in;
142   in.open(fileList);
143   // Read the input list of files 
144   TString objfile;
145   Int_t counter=0;
146   while(in.good()) {
147     in >> objfile;
148     if (!objfile.Contains("root")) continue; // protection    
149     printf("Open file:Counter\t%d\tMerging file %s\n",counter++,objfile.Data());
150     TFile currentFile(objfile.Data());
151     Merge(&currentFile, mergeArray);
152   }
153   if (separate) 
154     StoreSeparateResults(mergeArray, outputFileName);
155   else
156     StoreResults(mergeArray, outputFileName);
157   delete mergeArray;
158 }
159
160 void AliFileMerger::StoreResults(TObjArray * array, const char* outputFileName){
161   //
162   // Storing the results in one single file
163   //
164   TFile *f = new TFile(outputFileName,"recreate");
165   for (Int_t i=0; i<array->GetEntries(); i++){
166     TObject *object0 = array->At(i);
167     if (!object0) continue;
168     object0->Write();
169   }
170   f->Close();
171   delete f;
172 }
173
174
175 void AliFileMerger::StoreSeparateResults(TObjArray * array, const char* outputFileName){
176   //
177   // Store the results in separate files (one per object)
178   //
179   for (Int_t i=0; i<array->GetEntries(); i++){
180     TObject *object0 = array->At(i);
181     if (!object0) continue;
182     TFile *f = new TFile(Form("%s_%s.root",outputFileName,object0->GetName()),"recreate");
183     object0->Write();
184     f->Close();
185     delete f;
186   }
187 }
188
189
190
191 void AliFileMerger::Merge(TFile* fileIn, TObjArray * array){
192   //
193   // Merging procedure
194   //
195   static Int_t counter=-1;
196   counter++;
197   TObjArray *carray = new TObjArray;   //array of the objects inside current file
198   carray->SetOwner(kTRUE);
199   
200   // load all objects to  memory
201   
202   TList *farr = fileIn->GetListOfKeys();
203   if (!farr) return;
204   for (Int_t ical=0; ical<farr->GetEntries(); ical++){
205     if (!farr->At(ical)) continue;
206     TString name(farr->At(ical)->GetName());
207     if (!IsAccepted(name)) continue;                        // skip not accepted entries
208     TObject *obj = fileIn->Get(name.Data());
209     if (obj) carray->AddLast(obj);
210     AliSysInfo::AddStamp(name.Data(),1,ical,counter);  
211   }
212   
213   if (carray->GetEntries()==0) return;
214   TMethodCall callEnv;
215   
216   for (Int_t i=0; i<carray->GetEntries(); i++){
217     
218     TObjArray *templist = new TObjArray(1);
219     templist->SetOwner(kFALSE);
220     TObject *currentObject = carray->At(i);
221     if (!currentObject) continue;
222     printf("%s\n",currentObject->GetName());
223     callEnv.InitWithPrototype(currentObject->IsA(), "Merge", "TCollection*");
224     if (!callEnv.IsValid()) {continue;}
225     TString oname=currentObject->GetName();
226     TObject *mergedObject = array->FindObject(currentObject->GetName());
227     if (!mergedObject) {
228       array->AddLast(currentObject);
229       carray->RemoveAt(i);
230       continue;
231     }
232     templist->AddLast(currentObject);
233     callEnv.SetParam((Long_t) templist);
234     callEnv.Execute(mergedObject);
235     AliSysInfo::AddStamp(currentObject->GetName(),2,i,counter);  
236     delete templist;
237   }
238   delete carray;
239 }
240
241 Bool_t AliFileMerger::IsAccepted(TString name){
242   //
243   // Accept/reject logic
244   // name - name of the entry
245   //
246   //  if fAcceptMask specified   - entry has to be in list of selected
247   //  if fRejectMask speciefied  - entry with name speciief in the list are rejected 
248   //
249   Bool_t accept=kTRUE;
250   if (fAcceptMask){
251     //
252     accept=kFALSE;
253     for (Int_t iaccept=0; iaccept<fAcceptMask->GetEntries(); iaccept++){
254       if (name.Contains(fAcceptMask->At(iaccept)->GetName())) accept=kTRUE;   // entry was selected
255     }
256   }
257   if (!accept) return kFALSE;
258
259   if (fRejectMask){
260     //
261     for (Int_t ireject=0; ireject<fRejectMask->GetEntries(); ireject++){
262       if (name.Contains(fRejectMask->At(ireject)->GetName())) accept=kFALSE;   // entry was rejected
263     }
264   }
265   return accept;
266 }
267
268
269
270
271 void AliFileMerger::AddReject(const char *reject){
272   //
273   // add reject string to the list of entries to be rejected for merging
274   //
275   if (!fRejectMask) fRejectMask = new TObjArray;
276   fRejectMask->AddLast(new TObjString(reject));
277 }
278 void AliFileMerger::AddAccept(const char *accept){
279   //
280   // add reject string to the list of entries to be rejected for merging
281   //
282   if (!fAcceptMask) fAcceptMask = new TObjArray;
283   fAcceptMask->AddLast(new TObjString(accept));
284
285
286 }
287