]>
Commit | Line | Data |
---|---|---|
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> | |
29 | ||
30 | #include "AliRawReaderChain.h" | |
31 | #include "AliRawEvent.h" | |
32 | ||
33 | ClassImp(AliRawReaderChain) | |
34 | ||
35 | AliRawReaderChain::AliRawReaderChain() : | |
36 | AliRawReaderRoot(), | |
37 | fChain(NULL) | |
38 | { | |
39 | // default constructor | |
40 | } | |
41 | ||
acce5036 | 42 | AliRawReaderChain::AliRawReaderChain(const char* listFileName) : |
6923e953 | 43 | AliRawReaderRoot(), |
44 | fChain(NULL) | |
45 | { | |
46 | // create raw-reader objects which takes as an input a root chain | |
47 | // from the file list found in 'listFileName' | |
48 | ||
49 | TFileCollection collection("RAW", | |
50 | "Collection with raw-data files", | |
acce5036 | 51 | listFileName); |
6923e953 | 52 | |
acce5036 | 53 | fChain = new TChain("RAW"); |
6923e953 | 54 | if (!fChain->AddFileInfoList((TCollection*)(collection.GetList()))) { |
55 | Error("AliRawReaderChain","Bad file list in collection, the chain is empty"); | |
56 | return; | |
57 | } | |
58 | ||
636c1780 | 59 | fChain->SetBranchStatus("*",1); |
6923e953 | 60 | |
61 | fEvent = new AliRawEvent; | |
62 | fChain->SetBranchAddress("rawevent", &fEvent); | |
63 | } | |
64 | ||
65 | AliRawReaderChain::AliRawReaderChain(TFileCollection *collection) : | |
66 | AliRawReaderRoot(), | |
67 | fChain(NULL) | |
68 | { | |
69 | // create raw-reader objects which takes as an input a root chain | |
70 | // from a root file collection | |
71 | ||
acce5036 | 72 | fChain = new TChain("RAW"); |
6923e953 | 73 | if (!fChain->AddFileInfoList((TCollection*)(collection->GetList()))) { |
74 | Error("AliRawReaderChain","Bad file list in collection, the chain is empty"); | |
75 | return; | |
76 | } | |
77 | ||
636c1780 | 78 | fChain->SetBranchStatus("*",1); |
79 | ||
80 | fEvent = new AliRawEvent; | |
81 | fChain->SetBranchAddress("rawevent", &fEvent); | |
82 | } | |
83 | ||
84 | AliRawReaderChain::AliRawReaderChain(TChain *chain) : | |
85 | AliRawReaderRoot(), | |
86 | fChain(chain) | |
87 | { | |
88 | // create raw-reader objects which takes as an input a root chain | |
89 | // from a root file collection | |
90 | ||
91 | fChain->SetBranchStatus("*",1); | |
6923e953 | 92 | |
93 | fEvent = new AliRawEvent; | |
94 | fChain->SetBranchAddress("rawevent", &fEvent); | |
95 | } | |
96 | ||
97 | AliRawReaderChain::AliRawReaderChain(const AliRawReaderChain& rawReader) : | |
98 | AliRawReaderRoot(rawReader), | |
99 | fChain(rawReader.fChain) | |
100 | { | |
101 | // copy constructor | |
102 | } | |
103 | ||
104 | AliRawReaderChain& AliRawReaderChain::operator = (const AliRawReaderChain& | |
105 | rawReader) | |
106 | { | |
107 | // assignment operator | |
108 | ||
109 | this->~AliRawReaderChain(); | |
110 | new(this) AliRawReaderChain(rawReader); | |
111 | return *this; | |
112 | } | |
113 | ||
114 | AliRawReaderChain::~AliRawReaderChain() | |
115 | { | |
116 | // delete objects and close root file | |
117 | ||
118 | if (fChain) { | |
119 | delete fChain; | |
120 | fChain = NULL; | |
121 | } | |
122 | } | |
123 | ||
124 | Bool_t AliRawReaderChain::NextEvent() | |
125 | { | |
126 | // go to the next event in the root file | |
127 | ||
128 | if (!fChain || !fChain->GetListOfFiles()->GetEntriesFast()) return kFALSE; | |
129 | ||
130 | do { | |
131 | delete fEvent; | |
132 | fEvent = new AliRawEvent; | |
636c1780 | 133 | TBranch *branch = fChain->GetBranch("rawevent"); |
134 | if (!branch) | |
135 | return kFALSE; | |
136 | if (branch->GetEntry(fEventIndex+1) <= 0) | |
6923e953 | 137 | return kFALSE; |
138 | fEventIndex++; | |
139 | } while (!IsEventSelected()); | |
140 | fEventNumber++; | |
141 | return Reset(); | |
142 | } | |
143 | ||
144 | Bool_t AliRawReaderChain::RewindEvents() | |
145 | { | |
146 | // go back to the beginning of the root file | |
147 | ||
148 | fEventIndex = -1; | |
149 | delete fEvent; | |
150 | fEvent = new AliRawEvent; | |
151 | fEventNumber = -1; | |
152 | return Reset(); | |
153 | } | |
636c1780 | 154 | |
155 | Bool_t AliRawReaderChain::GotoEvent(Int_t event) | |
156 | { | |
157 | // go to a particular event | |
158 | // Uses the absolute event index inside the | |
159 | // chain with raw data | |
160 | ||
161 | if (!fChain || !fChain->GetListOfFiles()->GetEntriesFast()) return kFALSE; | |
162 | ||
163 | delete fEvent; | |
164 | fEvent = new AliRawEvent; | |
165 | TBranch *branch = fChain->GetBranch("rawevent"); | |
166 | if (!branch) | |
167 | return kFALSE; | |
168 | if (branch->GetEntry(event) <= 0) | |
169 | return kFALSE; | |
170 | fEventIndex = event; | |
171 | fEventNumber++; | |
172 | return Reset(); | |
173 | } |