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