]> git.uio.no Git - u/mrichter/AliRoot.git/blame - RAW/AliRawReaderChain.cxx
https://savannah.cern.ch/bugs/index.php?98544
[u/mrichter/AliRoot.git] / RAW / AliRawReaderChain.cxx
CommitLineData
6923e953 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///////////////////////////////////////////////////////////////////////////////
17///
18/// This is a class for reading raw data from a root chain.
19/// There are two constructors available - one from a text file containing the
20/// list of root raw-data files to be processed and one directly from
21/// TFileCollection.
22///
23/// cvetan.cheshkov@cern.ch 29/07/2008
24///
25///////////////////////////////////////////////////////////////////////////////
26
27#include <TChain.h>
28#include <TFileCollection.h>
3d456d99 29#include <TEntryList.h>
30#include "TGridCollection.h"
31#include <TPluginManager.h>
32#include <TROOT.h>
33#include <TSystem.h>
34#include <TFile.h>
35#include <TKey.h>
bb0edcaf 36#include <TGrid.h>
37#include <TGridResult.h>
6923e953 38
39#include "AliRawReaderChain.h"
33314186 40#include "AliRawVEvent.h"
9dafb98e 41#include "AliLog.h"
6923e953 42
43ClassImp(AliRawReaderChain)
44
9dafb98e 45TString AliRawReaderChain::fgSearchPath = "/alice/data";
46
6923e953 47AliRawReaderChain::AliRawReaderChain() :
48 AliRawReaderRoot(),
49 fChain(NULL)
50{
51 // default constructor
52}
53
3d456d99 54AliRawReaderChain::AliRawReaderChain(const char* fileName) :
6923e953 55 AliRawReaderRoot(),
56 fChain(NULL)
57{
58// create raw-reader objects which takes as an input a root chain
3d456d99 59// either from the file list found in 'fileName' (IsCollection = true)
60// or from entry list found in 'filename' (IsCollection = false)
61// The entry-list syntax follows root convetion: filename.root/listname
6923e953 62
acce5036 63 fChain = new TChain("RAW");
3d456d99 64
65 TString fileNameStr = fileName;
66 if (fileNameStr.EndsWith(".xml")) {
67
68 TGridCollection *collection = NULL;
69 TPluginManager* pluginManager = gROOT->GetPluginManager();
70 TPluginHandler* pluginHandler = pluginManager->FindHandler("TGridCollection", "alice");
71 if (!pluginHandler) {
72 pluginManager->AddHandler("TGridCollection", "alice",
73 "AliXMLCollection", "ANALYSISalice", "AliXMLCollection(const char*)");
74 pluginHandler = pluginManager->FindHandler("TGridCollection", "alice");
75 }
76 gSystem->Load("libANALYSIS");
77 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
78 collection = (TGridCollection*)pluginHandler->ExecPlugin(1,fileNameStr.Data());
79 }
80 else {
81 fIsValid = kFALSE;
82 return;
83 }
84 collection->Reset();
85 Bool_t elistsExist = kFALSE;
86 TEntryList *elist = new TEntryList();
87 while (collection->Next()) {
88 fChain->Add(collection->GetTURL(""));
89 TEntryList *list = (TEntryList *)collection->GetEntryList("");
90 if (list) {
91 list->SetTreeName("RAW");
92 list->SetFileName(collection->GetTURL(""));
93 elist->Add(list);
94 elistsExist = kTRUE;
95 }
96 }
97 if (elistsExist) {
98 fChain->SetEntryList(elist,"ne");
99 }
100 else {
101 Info("AliRawReaderChain", "no entry lists found in %s. Using all entries", fileNameStr.Data());
102 delete elist;
103 }
104 }
105 else if (fileNameStr.EndsWith(".root")) {
106
107 TDirectory* dir = gDirectory;
108 TFile *listFile = TFile::Open(fileNameStr.Data());
109 dir->cd();
110 if (!listFile || !listFile->IsOpen()) {
111 Error("AliRawReaderChain", "could not open file %s", fileNameStr.Data());
112 fIsValid = kFALSE;
113 return;
114 }
115
116 TEntryList *elist = NULL;
117 TKey *key = NULL;
118 TIter nextkey(listFile->GetListOfKeys());
119 while ((key=(TKey*)nextkey())){
120 if (strcmp("TEntryList", key->GetClassName())==0){
121 elist = (TEntryList*)key->ReadObj();
122 }
123 }
124 if (!elist) {
125 Error("AliRawReaderChain", "no TEntryList found in %s", fileNameStr.Data());
126 fIsValid = kFALSE;
127 return;
128 }
129
130 TEntryList *templist = NULL;
131 TList *elists = elist->GetLists();
132 TIter next(elists);
133 while((templist = (TEntryList*)next())){
134 Info("AliRawReaderChain", "%s added to the chain", templist->GetFileName());
135 fChain->Add(templist->GetFileName());
136 }
137 fChain->SetEntryList(elist,"ne");
138 }
139 else {
140
141 TFileCollection collection("RAW",
142 "Collection with raw-data files",
143 fileNameStr.Data());
144
145 if (!fChain->AddFileInfoList((TCollection*)(collection.GetList()))) {
146 Error("AliRawReaderChain","Bad file list in collection, the chain is empty");
147 fIsValid = kFALSE;
148 return;
149 }
6923e953 150 }
151
636c1780 152 fChain->SetBranchStatus("*",1);
4ebdfbcc 153 fChain->SetBranchAddress("rawevent",&fEvent,&fBranch);
6923e953 154}
155
156AliRawReaderChain::AliRawReaderChain(TFileCollection *collection) :
157 AliRawReaderRoot(),
158 fChain(NULL)
159{
160// create raw-reader objects which takes as an input a root chain
161// from a root file collection
162
acce5036 163 fChain = new TChain("RAW");
6923e953 164 if (!fChain->AddFileInfoList((TCollection*)(collection->GetList()))) {
165 Error("AliRawReaderChain","Bad file list in collection, the chain is empty");
a97af23d 166 fIsValid = kFALSE;
6923e953 167 return;
168 }
169
636c1780 170 fChain->SetBranchStatus("*",1);
4ebdfbcc 171 fChain->SetBranchAddress("rawevent",&fEvent,&fBranch);
636c1780 172}
173
174AliRawReaderChain::AliRawReaderChain(TChain *chain) :
175 AliRawReaderRoot(),
176 fChain(chain)
177{
178// create raw-reader objects which takes as an input a root chain
179// from a root file collection
180
68d8e599 181 if (!fChain) {
182 fIsValid = kFALSE;
183 return;
184 }
a97af23d 185
636c1780 186 fChain->SetBranchStatus("*",1);
4ebdfbcc 187 fChain->SetBranchAddress("rawevent",&fEvent,&fBranch);
6923e953 188}
189
3d456d99 190AliRawReaderChain::AliRawReaderChain(TEntryList *elist) :
191 AliRawReaderRoot(),
192 fChain(NULL)
193{
194// create raw-reader objects which takes as an input a root chain
195// from a root file collection
196
68d8e599 197 if (!elist) {
198 fIsValid = kFALSE;
199 return;
200 }
3d456d99 201
202 fChain = new TChain("RAW");
203
204 TEntryList *templist = NULL;
205 TList *elists = elist->GetLists();
206 TIter next(elists);
207 while((templist = (TEntryList*)next())){
208 fChain->Add(templist->GetFileName());
209 }
210 fChain->SetEntryList(elist,"ne");
211
212 fChain->SetBranchStatus("*",1);
213 fChain->SetBranchAddress("rawevent",&fEvent,&fBranch);
214}
215
bb0edcaf 216AliRawReaderChain::AliRawReaderChain(Int_t runNumber) :
217 AliRawReaderRoot(),
218 fChain(NULL)
219{
220// create raw-reader objects which takes as an input a root chain
221// with the raw-data files for a given run
222// It queries alien FC in order to do that and therefore
223// it needs alien API to be enabled
224
225 if (runNumber <= 0) {
226 Error("AliRawReaderChain","Bad run number:%d",runNumber);
227 fIsValid = kFALSE;
228 }
229
230 if (!gGrid) TGrid::Connect("alien://");
231 if (!gGrid) {
232 fIsValid = kFALSE;
233 return;
234 }
235
9dafb98e 236 if (fgSearchPath.IsNull()) fgSearchPath = "/alice/data";
237 TGridResult *res = gGrid->Query(fgSearchPath.Data(),Form("%09d/raw/*%09d*.root",runNumber,runNumber));
bb0edcaf 238 Int_t nFiles = res->GetEntries();
239 if (!nFiles) {
240 Error("AliRawReaderChain","No raw-data files found for run %d",runNumber);
241 fIsValid = kFALSE;
242 delete res;
243 return;
244 }
245
246 fChain = new TChain("RAW");
247 for (Int_t i = 0; i < nFiles; i++) {
248 TString filename = res->GetKey(i, "turl");
249 if(filename == "") continue;
250 fChain->Add(filename.Data());
251 }
252 delete res;
253
254 fChain->SetBranchStatus("*",1);
255 fChain->SetBranchAddress("rawevent",&fEvent,&fBranch);
256}
257
258
6923e953 259AliRawReaderChain::AliRawReaderChain(const AliRawReaderChain& rawReader) :
260 AliRawReaderRoot(rawReader),
261 fChain(rawReader.fChain)
262{
263// copy constructor
264}
265
266AliRawReaderChain& AliRawReaderChain::operator = (const AliRawReaderChain&
267 rawReader)
268{
269// assignment operator
270
271 this->~AliRawReaderChain();
272 new(this) AliRawReaderChain(rawReader);
273 return *this;
274}
275
276AliRawReaderChain::~AliRawReaderChain()
277{
278// delete objects and close root file
279
280 if (fChain) {
281 delete fChain;
282 fChain = NULL;
283 }
284}
285
286Bool_t AliRawReaderChain::NextEvent()
287{
288// go to the next event in the root file
289
290 if (!fChain || !fChain->GetListOfFiles()->GetEntriesFast()) return kFALSE;
291
292 do {
293 delete fEvent;
33314186 294 fEvent = NULL;
295 fEventHeader = NULL;
4ebdfbcc 296 Long64_t treeEntry = fChain->LoadTree(fEventIndex+1);
297 if (!fBranch)
636c1780 298 return kFALSE;
4ebdfbcc 299 if (fBranch->GetEntry(treeEntry) <= 0)
6923e953 300 return kFALSE;
33314186 301 fEventHeader = fEvent->GetHeader();
6923e953 302 fEventIndex++;
303 } while (!IsEventSelected());
304 fEventNumber++;
305 return Reset();
306}
307
308Bool_t AliRawReaderChain::RewindEvents()
309{
310// go back to the beginning of the root file
311
312 fEventIndex = -1;
313 delete fEvent;
33314186 314 fEvent = NULL;
315 fEventHeader = NULL;
6923e953 316 fEventNumber = -1;
317 return Reset();
318}
636c1780 319
320Bool_t AliRawReaderChain::GotoEvent(Int_t event)
321{
322 // go to a particular event
323 // Uses the absolute event index inside the
324 // chain with raw data
325
326 if (!fChain || !fChain->GetListOfFiles()->GetEntriesFast()) return kFALSE;
327
328 delete fEvent;
33314186 329 fEvent = NULL;
330 fEventHeader = NULL;
4ebdfbcc 331 Long64_t treeEntry = fChain->LoadTree(event);
332 if (!fBranch)
636c1780 333 return kFALSE;
4ebdfbcc 334 if (fBranch->GetEntry(treeEntry) <= 0)
636c1780 335 return kFALSE;
33314186 336 fEventHeader = fEvent->GetHeader();
636c1780 337 fEventIndex = event;
338 fEventNumber++;
339 return Reset();
340}
25e82ff5 341
342Int_t AliRawReaderChain::GetNumberOfEvents() const
343{
344 // Get the total number of events in the chain
345 // of raw-data files
346
347 if (!fChain) return -1;
348
349 return fChain->GetEntries();
350}
9dafb98e 351
352void AliRawReaderChain::SetSearchPath(const char* path)
353{
354 // set alien query search path
355 AliInfoGeneral("SetSearchPath",Form("Setting search path to \"%s\" (was \"%s\")",path,fgSearchPath.Data()));
356 fgSearchPath = path;
357}