3 //**************************************************************************
4 //* This file is property of and copyright by the *
5 //* ALICE Experiment at CERN, All rights reserved. *
7 //* Primary Authors: Matthias Richter <Matthias.Richter@ift.uib.no> *
9 //* Permission to use, copy, modify and distribute this software and its *
10 //* documentation strictly for non-commercial purposes is hereby granted *
11 //* without fee, provided that the above copyright notice appears in all *
12 //* copies and that both the copyright notice and this permission notice *
13 //* appear in the supporting documentation. The authors make no claims *
14 //* about the suitability of this software for any purpose. It is *
15 //* provided "as is" without express or implied warranty. *
16 //**************************************************************************
18 /// @file AliHLTOUTDigitReader.cxx
19 /// @author Matthias Richter
21 /// @brief HLTOUT data wrapper for simulated AliRoot HLT digit data.
24 #include "AliHLTOUTDigitReader.h"
25 #include "AliRawDataHeader.h"
34 /** ROOT macro for the implementation of ROOT specific class methods */
35 ClassImp(AliHLTOUTDigitReader)
37 AliHLTOUTDigitReader::AliHLTOUTDigitReader(int event, AliHLTEsdManager* pEsdManager, const char* digitFile)
39 AliHLTOUTHomerCollection(event, pEsdManager),
40 fDigitFileName(digitFile),
52 // HLTOUT data wrapper for simulated AliRoot HLT digit data
54 // see header file for class documentation
57 AliHLTOUTDigitReader::~AliHLTOUTDigitReader()
69 Bool_t AliHLTOUTDigitReader::ReadNextData(UChar_t*& data)
71 // overloaded from AliHLTOUTHomerCollection: switch to next DDL
72 if (!fppDigitArrays && (!ReadArrays() || !fppDigitArrays))
75 if (fCurrentLink>=fNofDDLs)
78 while (++fCurrentLink<fNofDDLs) {
79 if (fMinDDL>-1 && fMinDDL>fpEquipments[fCurrentLink]) continue;
80 if (fMaxDDL>-1 && fMaxDDL<fpEquipments[fCurrentLink]) continue;
81 if (fppDigitArrays[fCurrentLink]->GetSize()>(int)sizeof(AliRawDataHeader)) {
82 data=reinterpret_cast<UChar_t*>(fppDigitArrays[fCurrentLink]->GetArray());
83 data+=sizeof(AliRawDataHeader);
90 int AliHLTOUTDigitReader::Reset()
92 // overloaded from AliHLTOUTHomerCollection: reset DDL position
97 int AliHLTOUTDigitReader::GetDataSize()
99 // overloaded from AliHLTOUTHomerCollection: get size of current DDL
100 if (fCurrentLink>=0 && fCurrentLink<fNofDDLs && fppDigitArrays) {
101 return fppDigitArrays[fCurrentLink]->GetSize()-sizeof(AliRawDataHeader);
106 const AliRawDataHeader* AliHLTOUTDigitReader::GetDataHeader()
108 // overloaded from AliHLTOUTHomerCollection: get data header of current DDL
109 if (fCurrentLink>=0 && fCurrentLink<fNofDDLs && fppDigitArrays) {
110 return reinterpret_cast<AliRawDataHeader*>(fppDigitArrays[fCurrentLink]->GetArray());
115 void AliHLTOUTDigitReader::SelectEquipment(int /*equipmentType*/, int minEquipmentId, int maxEquipmentId)
117 // overloaded from AliHLTOUTHomerCollection: select equipment range
118 fMinDDL=minEquipmentId;
119 fMaxDDL=maxEquipmentId;
122 int AliHLTOUTDigitReader::GetEquipmentId()
124 // overloaded from AliHLTOUTHomerCollection: get id of current DDL
125 if (fCurrentLink>=0 && fCurrentLink<fNofDDLs && fpEquipments) {
126 return fpEquipments[fCurrentLink];
131 bool AliHLTOUTDigitReader::ReadArrays()
133 // Read the data from the root file and HLTOUT raw tree.
134 // Retrieve the number of branches and allocate arrays acording
135 // to that. After initialization of the arrays and variables, the
136 // event fEnvent is read.
139 if (GetCurrentEventNo()<0) {
140 HLTWarning("no event selected, no data available");
144 // 2011-06-06 in order to support AliHLTReconstructor option 'ignore-hltout' for
145 // digits, the file name can be set to NULL, nothing done in that case
146 if (!fDigitFileName) return false;
149 fpDigitFile=new TFile(fDigitFileName);
151 if (!fpDigitFile) return false;
152 if (fpDigitFile->IsZombie()) return false;
154 fpDigitFile->GetObject("rawhltout", fpDigitTree);
159 vector<int> equipments;
161 TIter iter(fpDigitTree->GetListOfBranches());
162 while (TBranch* br=(TBranch *)iter.Next()) {
164 TString bname=br->GetName();
165 if (bname.BeginsWith("HLT_") && bname.EndsWith(".ddl")) {
166 bname.ReplaceAll("HLT_", "");
167 bname.ReplaceAll(".ddl", "");
168 if (bname.IsDigit()) {
169 int equipment=bname.Atoi();
170 int index=equipment-AliDAQ::DdlIDOffset("HLT");
171 if (index>=0 && strcmp(br->GetName(), AliDAQ::DdlFileName("HLT", index))==0) {
172 equipments.push_back(equipment);
174 HLTWarning("equipment mismatch: skipping id %d (out of range of HLT equipments)", equipment);
179 HLTDebug("found %d HLT branche(s) out of %d", equipments.size(), iTotal);
181 if (equipments.size()>0) {
182 fNofDDLs=equipments.size();
183 fpEquipments=new int[fNofDDLs];
184 fppDigitArrays=new TArrayC*[fNofDDLs];
185 if (fpEquipments && fppDigitArrays) {
186 memcpy(fpEquipments, &equipments[0], fNofDDLs*sizeof(int));
188 for (i=0; i<fNofDDLs; i++) {
189 fppDigitArrays[i]=NULL;
190 fpDigitTree->SetBranchAddress(AliDAQ::DdlFileName("HLT", fpEquipments[i]-AliDAQ::DdlIDOffset("HLT")), &fppDigitArrays[i]);
192 fpDigitTree->GetEvent(GetCurrentEventNo());
193 for (i=0; i<fNofDDLs; i++) {
194 if (fppDigitArrays[i]) {
195 HLTDebug("branch %s: %d byte(s)", AliDAQ::DdlFileName("HLT", fpEquipments[i]-AliDAQ::DdlIDOffset("HLT")), fppDigitArrays[i]);
206 int AliHLTOUTDigitReader::CloseTree()
208 // Cleanup tree and data arrays.
210 for (int i=0; i<fNofDDLs; i++) {
211 if (fppDigitArrays[i]) delete fppDigitArrays[i];
212 fppDigitArrays[i]=NULL;
214 if (fppDigitArrays) delete[] fppDigitArrays;
216 if (fpEquipments) delete[] fpEquipments;
218 if (fpDigitTree) delete fpDigitTree;
226 void AliHLTOUTDigitReader::SetParam(TTree* /*pDigitTree*/, int event)
228 // set parameter for this HLTOUT instance
229 // The function is for internal use only in conjunction with the
230 // AliHLTOUT::New() functions.
232 // TODO: here we have the implemented to correct loading of
233 // the digit file from the run loader. At time of writing
234 // (Jun 2008) the AliLoader for the HLT is missing in the AliRoot
239 void AliHLTOUTDigitReader::SetParam(const char* filename, int event)
241 // set parameter for this HLTOUT instance
242 // The function is for internal use only in conjunction with the
243 // AliHLTOUT::New() functions.
245 if (filename && filename[0]!=0) {
246 fDigitFileName=filename;
248 HLTWarning("no valid digit file provided, using default file %s", fDigitFileName.Data());