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