1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 //////////////////////////////////////////////////////////////////////////
19 // marian.ivanov@cern.ch
20 // Utilities for file merging.
21 // Additional functionality on top of the standard TFileMerger:
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.
28 // 2. syswatch.log is created diring mergin procedure.
29 // Memeory consumption - for reading and for merging can be monitored
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");
40 //Example usage starting from the input data list in text file:
43 merger.AddReject("esdFriend");
44 merger.IterTXT("calib.list","CalibObjects.root",kFALSE);
48 //////////////////////////////////////////////////////////////////////////
55 #include "TGridResult.h"
56 #include "TObjString.h"
57 #include "TObjArray.h"
58 #include "TMethodCall.h"
60 #include "AliSysInfo.h"
61 #include "AliFileMerger.h"
63 ClassImp(AliFileMerger)
65 ////////////////////////////////////////////////////////////////////////
67 AliFileMerger::AliFileMerger():
73 // Default constructor
77 //______________________________________________________________________
79 AliFileMerger::AliFileMerger(const char* name):
90 void AliFileMerger::IterAlien(const char* outputDir, const char* outputFileName, const char* pattern){
93 // Merge the files coming out of the calibration job
96 TObjArray * mergeArray= new TObjArray;
98 TString outputFile(outputFileName);
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);
108 // loop over the results
109 while((map=(TMap*)nextmap())) {
111 TObjString *objs = dynamic_cast<TObjString*>(map->GetValue("turl"));
112 if (!objs || !objs->GetString().Length()) {
113 // Nothing found - skip this output
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);
122 if(currentFile) delete currentFile;
125 // StoreSeparateResults(mergeArray,outputFileName);
126 StoreResults(mergeArray,outputFileName);
135 void AliFileMerger::IterTXT( const char * fileList, const char* outputFileName, Bool_t separate){
137 // Merge the files indicated in the list - fileList
138 // ASCII file option example:
139 // find `pwd`/ | grep AliESDfriends_v1.root > calib.list
141 TObjArray * mergeArray= new TObjArray;
143 // Open the input stream
147 // Read the input list of files
152 if (!objfile.Contains("root")) continue; // protection
153 printf("Open file:Counter\t%d\tMerging file %s\n",counter++,objfile.Data());
154 TFile currentFile(objfile.Data());
155 Merge(¤tFile, mergeArray);
158 StoreSeparateResults(mergeArray, outputFileName);
161 StoreResults(mergeArray, outputFileName);
167 void AliFileMerger::StoreResults(TObjArray * array, const char* outputFileName){
169 // Storing the results in one single file
171 TFile *f = new TFile(outputFileName,"recreate");
172 for (Int_t i=0; i<array->GetEntries(); i++){
173 TObject *object0 = array->At(i);
174 if (!object0) continue;
182 void AliFileMerger::StoreSeparateResults(TObjArray * array, const char* outputFileName){
184 // Store the results in separate files (one per object)
186 for (Int_t i=0; i<array->GetEntries(); i++){
187 TObject *object0 = array->At(i);
188 if (!object0) continue;
189 TFile *f = new TFile(Form("%s_%s.root",outputFileName,object0->GetName()),"recreate");
198 void AliFileMerger::Merge(TFile* fileIn, TObjArray * array){
203 static Int_t counter=-1;
205 TObjArray *carray = new TObjArray; //array of the objects inside current file
206 carray->SetOwner(kTRUE);
208 // load all objects to memory
210 TList *farr = fileIn->GetListOfKeys();
215 for (Int_t ical=0; ical<farr->GetEntries(); ical++){
216 if (!farr->At(ical)) continue;
217 TString name(farr->At(ical)->GetName());
218 if (!IsAccepted(name)) continue; // skip not accepted entries
219 TObject *obj = fileIn->Get(name.Data());
220 if (obj) carray->AddLast(obj);
221 AliSysInfo::AddStamp(name.Data(),1,ical,counter);
224 if (carray->GetEntries()==0) {
229 Int_t entries =carray->GetEntriesFast();
230 for (Int_t i=0; i<entries; i++){
232 TObjArray *templist = new TObjArray(1);
233 templist->SetOwner(kFALSE);
234 TObject *currentObject = carray->At(i);
235 if (!currentObject) {
239 printf("%s\n",currentObject->GetName());
240 callEnv.InitWithPrototype(currentObject->IsA(), "Merge", "TCollection*");
241 if (!callEnv.IsValid()) {
245 TString oname=currentObject->GetName();
246 TObject *mergedObject = array->FindObject(currentObject->GetName());
248 array->AddLast(currentObject);
253 templist->AddLast(currentObject);
254 callEnv.SetParam((Long_t) templist);
255 callEnv.Execute(mergedObject);
256 AliSysInfo::AddStamp(currentObject->GetName(),2,i,counter);
263 Bool_t AliFileMerger::IsAccepted(TString name){
265 // Accept/reject logic
266 // name - name of the entry
268 // if fAcceptMask specified - entry has to be in list of selected
269 // if fRejectMask speciefied - entry with name speciief in the list are rejected
275 for (Int_t iaccept=0; iaccept<fAcceptMask->GetEntries(); iaccept++){
276 if (name.Contains(fAcceptMask->At(iaccept)->GetName())) accept=kTRUE; // entry was selected
279 if (!accept) return kFALSE;
283 for (Int_t ireject=0; ireject<fRejectMask->GetEntries(); ireject++){
284 if (name.Contains(fRejectMask->At(ireject)->GetName())) accept=kFALSE; // entry was rejected
293 void AliFileMerger::AddReject(const char *reject){
295 // add reject string to the list of entries to be rejected for merging
297 if (!fRejectMask) fRejectMask = new TObjArray;
298 fRejectMask->AddLast(new TObjString(reject));
300 void AliFileMerger::AddAccept(const char *accept){
302 // add reject string to the list of entries to be rejected for merging
304 if (!fAcceptMask) fAcceptMask = new TObjArray;
305 fAcceptMask->AddLast(new TObjString(accept));