c5123824 |
1 | // $Id$ |
2 | |
3 | //************************************************************************** |
4 | //* This file is property of and copyright by the ALICE HLT Project * |
5 | //* ALICE Experiment at CERN, All rights reserved. * |
6 | //* * |
7 | //* Primary Authors: Matthias Richter <Matthias.Richter@ift.uib.no> * |
8 | //* for The ALICE HLT Project. * |
9 | //* * |
10 | //* Permission to use, copy, modify and distribute this software and its * |
11 | //* documentation strictly for non-commercial purposes is hereby granted * |
12 | //* without fee, provided that the above copyright notice appears in all * |
13 | //* copies and that both the copyright notice and this permission notice * |
14 | //* appear in the supporting documentation. The authors make no claims * |
15 | //* about the suitability of this software for any purpose. It is * |
16 | //* provided "as is" without express or implied warranty. * |
17 | //************************************************************************** |
18 | |
19 | /** @file AliHLTOUTDigitReader.cxx |
20 | @author Matthias Richter |
21 | @date |
22 | @brief HLTOUT data wrapper for simulated AliRoot HLT digit data. |
23 | */ |
24 | |
25 | #include "AliHLTOUTDigitReader.h" |
26 | #include "AliRawDataHeader.h" |
27 | #include "AliDAQ.h" |
28 | #include "TTree.h" |
29 | #include "TFile.h" |
30 | #include "TArrayC.h" |
31 | #include "TBranch.h" |
32 | #include "TString.h" |
33 | #include <cassert> |
34 | |
35 | /** ROOT macro for the implementation of ROOT specific class methods */ |
36 | ClassImp(AliHLTOUTDigitReader) |
37 | |
032c5e5e |
38 | AliHLTOUTDigitReader::AliHLTOUTDigitReader(int event, AliHLTEsdManager* pEsdManager, const char* digitFile) |
c5123824 |
39 | : |
40 | AliHLTOUTHomerCollection(event, pEsdManager), |
032c5e5e |
41 | fDigitFileName(digitFile), |
c5123824 |
42 | fpDigitFile(NULL), |
43 | fpDigitTree(NULL), |
44 | fMinDDL(-1), |
45 | fMaxDDL(-1), |
46 | fppDigitArrays(NULL), |
47 | fpEquipments(NULL), |
48 | fNofDDLs(0), |
49 | fCurrent(-1) |
50 | { |
51 | // see header file for class documentation |
52 | // or |
53 | // refer to README to build package |
54 | // or |
55 | // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt |
56 | } |
57 | |
58 | AliHLTOUTDigitReader::~AliHLTOUTDigitReader() |
59 | { |
60 | // see header file for class documentation |
61 | CloseTree(); |
62 | |
63 | if (fpDigitFile) { |
64 | fpDigitFile->Close(); |
65 | delete fpDigitFile; |
66 | fpDigitFile=NULL; |
67 | } |
68 | } |
69 | |
70 | Bool_t AliHLTOUTDigitReader::ReadNextData(UChar_t*& data) |
71 | { |
72 | // see header file for class documentation |
73 | if (!fppDigitArrays && (!ReadArrays() || !fppDigitArrays)) |
74 | return kFALSE; |
75 | |
76 | if (fCurrent>=fNofDDLs) |
77 | return kFALSE; |
78 | |
79 | while (++fCurrent<fNofDDLs) { |
80 | if (fMinDDL>-1 && fMinDDL>fpEquipments[fCurrent]) continue; |
81 | if (fMaxDDL>-1 && fMaxDDL<fpEquipments[fCurrent]) continue; |
82 | if (fppDigitArrays[fCurrent]->GetSize()>(int)sizeof(AliRawDataHeader)) { |
83 | data=reinterpret_cast<UChar_t*>(fppDigitArrays[fCurrent]->GetArray()); |
84 | data+=sizeof(AliRawDataHeader); |
85 | return kTRUE; |
86 | } |
87 | } |
88 | return kFALSE; |
89 | } |
90 | |
91 | int AliHLTOUTDigitReader::Reset() |
92 | { |
93 | // see header file for class documentation |
94 | fCurrent=-1; |
95 | return 0; |
96 | } |
97 | |
98 | int AliHLTOUTDigitReader::GetDataSize() |
99 | { |
100 | // see header file for class documentation |
101 | if (fCurrent>=0 && fCurrent<fNofDDLs && fppDigitArrays) { |
102 | return fppDigitArrays[fCurrent]->GetSize()-sizeof(AliRawDataHeader); |
103 | } |
104 | return 0; |
105 | } |
106 | |
107 | const AliRawDataHeader* AliHLTOUTDigitReader::GetDataHeader() |
108 | { |
109 | // see header file for class documentation |
110 | if (fCurrent>=0 && fCurrent<fNofDDLs && fppDigitArrays) { |
111 | return reinterpret_cast<AliRawDataHeader*>(fppDigitArrays[fCurrent]->GetArray()); |
112 | } |
113 | return NULL; |
114 | } |
115 | |
116 | void AliHLTOUTDigitReader::SelectEquipment(int /*equipmentType*/, int minEquipmentId, int maxEquipmentId) |
117 | { |
118 | // see header file for class documentation |
119 | fMinDDL=minEquipmentId; |
120 | fMaxDDL=maxEquipmentId; |
121 | } |
122 | |
123 | int AliHLTOUTDigitReader::GetEquipmentId() |
124 | { |
125 | // see header file for class documentation |
126 | if (fCurrent>=0 && fCurrent<fNofDDLs && fpEquipments) { |
127 | return fpEquipments[fCurrent]; |
128 | } |
129 | return -1; |
130 | } |
131 | |
132 | bool AliHLTOUTDigitReader::ReadArrays() |
133 | { |
134 | // see header file for class documentation |
135 | bool result=false; |
136 | |
137 | if (GetCurrentEventNo()<0) { |
138 | HLTWarning("no event selected, no data available"); |
139 | return false; |
140 | } |
141 | |
7ae41d91 |
142 | // 2011-06-06 in order to support AliHLTReconstructor option 'ignore-hltout' for |
143 | // digits, the file name can be set to NULL, nothing done in that case |
144 | if (!fDigitFileName) return false; |
145 | |
c5123824 |
146 | if (!fpDigitFile) { |
032c5e5e |
147 | fpDigitFile=new TFile(fDigitFileName); |
c5123824 |
148 | } |
149 | if (!fpDigitFile) return false; |
150 | if (fpDigitFile->IsZombie()) return false; |
151 | |
152 | fpDigitFile->GetObject("rawhltout", fpDigitTree); |
153 | if (!fpDigitTree) { |
154 | return false; |
155 | } |
156 | |
157 | vector<int> equipments; |
158 | int iTotal=0; |
159 | TIter iter(fpDigitTree->GetListOfBranches()); |
160 | while (TBranch* br=(TBranch *)iter.Next()) { |
161 | iTotal++; |
162 | TString bname=br->GetName(); |
163 | if (bname.BeginsWith("HLT_") && bname.EndsWith(".ddl")) { |
164 | bname.ReplaceAll("HLT_", ""); |
165 | bname.ReplaceAll(".ddl", ""); |
166 | if (bname.IsDigit()) { |
167 | int equipment=bname.Atoi(); |
168 | int index=equipment-AliDAQ::DdlIDOffset("HLT"); |
169 | if (index>=0 && strcmp(br->GetName(), AliDAQ::DdlFileName("HLT", index))==0) { |
170 | equipments.push_back(equipment); |
171 | } else { |
172 | HLTWarning("equipment mismatch: skipping id %d (out of range of HLT equipments)", equipment); |
173 | } |
174 | } |
175 | } |
176 | } |
177 | HLTDebug("found %d HLT branche(s) out of %d", equipments.size(), iTotal); |
178 | |
179 | if (equipments.size()>0) { |
180 | fNofDDLs=equipments.size(); |
181 | fpEquipments=new int[fNofDDLs]; |
182 | fppDigitArrays=new TArrayC*[fNofDDLs]; |
183 | if (fpEquipments && fppDigitArrays) { |
184 | memcpy(fpEquipments, &equipments[0], fNofDDLs*sizeof(int)); |
185 | int i=0; |
186 | for (i=0; i<fNofDDLs; i++) { |
187 | fppDigitArrays[i]=NULL; |
188 | fpDigitTree->SetBranchAddress(AliDAQ::DdlFileName("HLT", fpEquipments[i]-AliDAQ::DdlIDOffset("HLT")), &fppDigitArrays[i]); |
189 | } |
190 | fpDigitTree->GetEvent(GetCurrentEventNo()); |
191 | for (i=0; i<fNofDDLs; i++) { |
192 | if (fppDigitArrays[i]) { |
193 | HLTDebug("branch %s: %d byte(s)", AliDAQ::DdlFileName("HLT", fpEquipments[i]-AliDAQ::DdlIDOffset("HLT")), fppDigitArrays[i]); |
194 | } |
195 | } |
196 | result=true; |
197 | } else { |
198 | fNofDDLs=0; |
199 | } |
200 | } |
201 | return result; |
202 | } |
203 | |
204 | int AliHLTOUTDigitReader::CloseTree() |
205 | { |
206 | // see header file for class documentation |
207 | int iResult=0; |
d92e81d0 |
208 | for (int i=0; i<fNofDDLs; i++) { |
209 | if (fppDigitArrays[i]) delete fppDigitArrays[i]; |
210 | fppDigitArrays[i]=NULL; |
211 | } |
c5123824 |
212 | if (fppDigitArrays) delete[] fppDigitArrays; |
213 | fppDigitArrays=NULL; |
214 | if (fpEquipments) delete[] fpEquipments; |
215 | fpEquipments=NULL; |
216 | if (fpDigitTree) delete fpDigitTree; |
217 | fpDigitTree=NULL; |
218 | fNofDDLs=0; |
219 | fCurrent=-1; |
220 | |
221 | return iResult; |
222 | } |
c1292031 |
223 | |
224 | void AliHLTOUTDigitReader::SetParam(TTree* /*pDigitTree*/, int event) |
225 | { |
226 | // see header file for class documentation |
227 | |
228 | // TODO: here we have the implemented to correct loading of |
229 | // the digit file from the run loader. At time of writing |
230 | // (Jun 2008) the AliLoader for the HLT is missing in the AliRoot |
231 | // framework. |
232 | fEvent=event; |
233 | } |
466d4e62 |
234 | |
235 | void AliHLTOUTDigitReader::SetParam(const char* filename, int event) |
236 | { |
237 | // see header file for class documentation |
238 | |
239 | if (filename && filename[0]!=0) { |
240 | fDigitFileName=filename; |
241 | } else { |
242 | HLTWarning("no valid digit file provided, using default file %s", fDigitFileName.Data()); |
243 | } |
244 | fEvent=event; |
245 | } |