]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ANALYSIS/AliFileMerger.cxx
Fixes for coverity warnings
[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     if(!currentFile) continue; // protection
120     Merge(currentFile, mergeArray);
121
122     if(currentFile) delete currentFile;
123   }
124   Bool_t separate = kFALSE;
125   if (separate) {
126     StoreSeparateResults(mergeArray,outputFileName);
127   }
128   else {
129     StoreResults(mergeArray,outputFileName);
130   }
131   delete mergeArray;
132   delete res;
133
134 }
135
136
137
138 void AliFileMerger::IterTXT( const char * fileList,  const char* outputFileName, Bool_t separate){
139   
140   // Merge the files indicated in the list - fileList
141   // ASCII file option example: 
142   // find `pwd`/ | grep AliESDfriends_v1.root > calib.list
143   
144   TObjArray * mergeArray= new TObjArray;
145   
146   // Open the input stream
147   
148   ifstream in;
149   in.open(fileList);
150   // Read the input list of files 
151   TString objfile;
152   Int_t counter=0;
153   while(in.good()) {
154     in >> objfile;
155     if (!objfile.Contains("root")) continue; // protection    
156     printf("Open file:Counter\t%d\tMerging file %s\n",counter++,objfile.Data());
157     TFile currentFile(objfile.Data());
158     Merge(&currentFile, mergeArray);
159   }
160   if (separate) { 
161     StoreSeparateResults(mergeArray, outputFileName);
162   }
163   else {
164     StoreResults(mergeArray, outputFileName);
165   }
166
167   delete mergeArray;
168 }
169
170 void AliFileMerger::StoreResults(TObjArray * array, const char* outputFileName){
171   //
172   // Storing the results in one single file
173   //
174   TFile *f = new TFile(outputFileName,"recreate");
175   for (Int_t i=0; i<array->GetEntries(); i++){
176     TObject *object0 = array->At(i);
177     if (!object0) continue;
178     object0->Write();
179   }
180   f->Close();
181   delete f;
182 }
183
184
185 void AliFileMerger::StoreSeparateResults(TObjArray * array, const char* outputFileName){
186   //
187   // Store the results in separate files (one per object)
188   //
189   for (Int_t i=0; i<array->GetEntries(); i++){
190     TObject *object0 = array->At(i);
191     if (!object0) continue;
192     TFile *f = new TFile(Form("%s_%s.root",outputFileName,object0->GetName()),"recreate");
193     object0->Write();
194     f->Close();
195     delete f;
196   }
197 }
198
199
200
201 void AliFileMerger::Merge(TFile* fileIn, TObjArray * array){
202   //
203   // Merging procedure
204   //
205   static Int_t counter=-1;
206   counter++;
207   TObjArray *carray = new TObjArray;   //array of the objects inside current file
208   carray->SetOwner(kTRUE);
209   
210   // load all objects to  memory
211   
212   TList *farr = fileIn->GetListOfKeys();
213   if (!farr) { 
214     delete carray;
215     return;
216   }
217   for (Int_t ical=0; ical<farr->GetEntries(); ical++){
218     if (!farr->At(ical)) continue;
219     TString name(farr->At(ical)->GetName());
220     if (!IsAccepted(name)) continue;                        // skip not accepted entries
221     TObject *obj = fileIn->Get(name.Data());
222     if (obj) carray->AddLast(obj);
223     AliSysInfo::AddStamp(name.Data(),1,ical,counter);  
224   }
225   
226   if (carray->GetEntries()==0)  { 
227     delete carray;
228     return;
229   }
230   TMethodCall callEnv;
231   
232   for (Int_t i=0; i<carray->GetEntries(); i++){
233     
234     TObjArray *templist = new TObjArray(1);
235     templist->SetOwner(kFALSE);
236     TObject *currentObject = carray->At(i);
237     if (!currentObject) { 
238       delete templist;
239       continue;
240     }
241     printf("%s\n",currentObject->GetName());
242     callEnv.InitWithPrototype(currentObject->IsA(), "Merge", "TCollection*");
243     if (!callEnv.IsValid()) {
244       delete templist; 
245       continue;
246     }
247     TString oname=currentObject->GetName();
248     TObject *mergedObject = array->FindObject(currentObject->GetName());
249     if (!mergedObject) {
250       array->AddLast(currentObject);
251       carray->RemoveAt(i);
252       delete templist; 
253       continue;
254     }
255     templist->AddLast(currentObject);
256     callEnv.SetParam((Long_t) templist);
257     callEnv.Execute(mergedObject);
258     AliSysInfo::AddStamp(currentObject->GetName(),2,i,counter);  
259     delete templist;
260   }
261   delete carray;
262 }
263
264 Bool_t AliFileMerger::IsAccepted(TString name){
265   //
266   // Accept/reject logic
267   // name - name of the entry
268   //
269   //  if fAcceptMask specified   - entry has to be in list of selected
270   //  if fRejectMask speciefied  - entry with name speciief in the list are rejected 
271   //
272   Bool_t accept=kTRUE;
273   if (fAcceptMask){
274     //
275     accept=kFALSE;
276     for (Int_t iaccept=0; iaccept<fAcceptMask->GetEntries(); iaccept++){
277       if (name.Contains(fAcceptMask->At(iaccept)->GetName())) accept=kTRUE;   // entry was selected
278     }
279   }
280   if (!accept) return kFALSE;
281
282   if (fRejectMask){
283     //
284     for (Int_t ireject=0; ireject<fRejectMask->GetEntries(); ireject++){
285       if (name.Contains(fRejectMask->At(ireject)->GetName())) accept=kFALSE;   // entry was rejected
286     }
287   }
288   return accept;
289 }
290
291
292
293
294 void AliFileMerger::AddReject(const char *reject){
295   //
296   // add reject string to the list of entries to be rejected for merging
297   //
298   if (!fRejectMask) fRejectMask = new TObjArray;
299   fRejectMask->AddLast(new TObjString(reject));
300 }
301 void AliFileMerger::AddAccept(const char *accept){
302   //
303   // add reject string to the list of entries to be rejected for merging
304   //
305   if (!fAcceptMask) fAcceptMask = new TObjArray;
306   fAcceptMask->AddLast(new TObjString(accept));
307
308
309 }
310