]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/rec/AliHLTOUTDigitReader.cxx
Update master to aliroot
[u/mrichter/AliRoot.git] / HLT / rec / AliHLTOUTDigitReader.cxx
1 // $Id$
2
3 //**************************************************************************
4 //* This file is property of and copyright by the                          * 
5 //* ALICE Experiment at CERN, All rights reserved.                         *
6 //*                                                                        *
7 //* Primary Authors: Matthias Richter <Matthias.Richter@ift.uib.no>        *
8 //*                                                                        *
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 //**************************************************************************
17
18 /// @file   AliHLTOUTDigitReader.cxx
19 /// @author Matthias Richter
20 /// @date   
21 /// @brief  HLTOUT data wrapper for simulated AliRoot HLT digit data.
22 ///
23
24 #include "AliHLTOUTDigitReader.h"
25 #include "AliHLTCDHWrapper.h"
26 #include "AliDAQ.h"
27 #include "TTree.h"
28 #include "TFile.h"
29 #include "TArrayC.h"
30 #include "TBranch.h"
31 #include "TString.h"
32 #include <cassert>
33
34 /** ROOT macro for the implementation of ROOT specific class methods */
35 ClassImp(AliHLTOUTDigitReader)
36
37 AliHLTOUTDigitReader::AliHLTOUTDigitReader(int event, AliHLTEsdManager* pEsdManager, const char* digitFile)
38   :
39   AliHLTOUTHomerCollection(event, pEsdManager),
40   fDigitFileName(digitFile),
41   fpDigitFile(NULL),
42   fpDigitTree(NULL),
43   fMinDDL(-1),
44   fMaxDDL(-1),
45   fppDigitArrays(NULL),
46   fpEquipments(NULL),
47   fNofDDLs(0),
48   fCurrentLink(-1)
49 {
50   // constructor
51   //
52   // HLTOUT data wrapper for simulated AliRoot HLT digit data
53   // 
54   // see header file for class documentation
55 }
56
57 AliHLTOUTDigitReader::~AliHLTOUTDigitReader()
58 {
59   // destructor
60   CloseTree();
61
62   if (fpDigitFile) {
63     fpDigitFile->Close();
64     delete fpDigitFile;
65     fpDigitFile=NULL;
66   }
67 }
68
69 Bool_t AliHLTOUTDigitReader::ReadNextData(UChar_t*& data)
70 {
71   // overloaded from AliHLTOUTHomerCollection: switch to next DDL
72   if (!fppDigitArrays && (!ReadArrays() || !fppDigitArrays))
73     return kFALSE;
74
75   if (fCurrentLink>=fNofDDLs)
76     return kFALSE;
77
78   while (++fCurrentLink<fNofDDLs) {
79     if (fMinDDL>-1 && fMinDDL>fpEquipments[fCurrentLink]) continue;
80     if (fMaxDDL>-1 && fMaxDDL<fpEquipments[fCurrentLink]) continue;
81     AliHLTCDHWrapper cdh((void*)fppDigitArrays[fCurrentLink]->GetArray());
82     UInt_t headerSize=cdh.GetHeaderSize();
83     if (fppDigitArrays[fCurrentLink]->GetSize()>(int)headerSize) {
84       data=reinterpret_cast<UChar_t*>(fppDigitArrays[fCurrentLink]->GetArray());
85       data+=headerSize;
86       return kTRUE;
87     }
88   }
89   return kFALSE;
90 }
91
92 int AliHLTOUTDigitReader::Reset()
93 {
94   // overloaded from AliHLTOUTHomerCollection: reset DDL position
95   fCurrentLink=-1;
96   return 0;
97 }
98
99 int AliHLTOUTDigitReader::GetDataSize()
100 {
101   // overloaded from AliHLTOUTHomerCollection: get size of current DDL
102   if (fCurrentLink>=0 && fCurrentLink<fNofDDLs && fppDigitArrays) {
103     AliHLTCDHWrapper cdh((void*)fppDigitArrays[fCurrentLink]->GetArray());
104     return fppDigitArrays[fCurrentLink]->GetSize()-cdh.GetHeaderSize();
105   }
106   return 0;
107 }
108
109 AliHLTCDHWrapper AliHLTOUTDigitReader::GetDataHeader()
110 {
111   // overloaded from AliHLTOUTHomerCollection: get data header of current DDL
112   if (fCurrentLink>=0 && fCurrentLink<fNofDDLs && fppDigitArrays) {
113     return AliHLTCDHWrapper(fppDigitArrays[fCurrentLink]->GetArray());
114   }
115   return NULL;
116 }
117
118 void AliHLTOUTDigitReader::SelectEquipment(int /*equipmentType*/, int minEquipmentId, int maxEquipmentId)
119 {
120   // overloaded from AliHLTOUTHomerCollection: select equipment range
121   fMinDDL=minEquipmentId;
122   fMaxDDL=maxEquipmentId;
123 }
124
125 int AliHLTOUTDigitReader::GetEquipmentId()
126 {
127   // overloaded from AliHLTOUTHomerCollection: get id of current DDL
128   if (fCurrentLink>=0 && fCurrentLink<fNofDDLs && fpEquipments) {
129     return fpEquipments[fCurrentLink];
130   }
131   return -1;
132 }
133
134 bool AliHLTOUTDigitReader::ReadArrays()
135 {
136   // Read the data from the root file and HLTOUT raw tree.
137   // Retrieve the number of branches and allocate arrays acording
138   // to that. After initialization of the arrays and variables, the
139   // event fEnvent is read.
140   bool result=false;
141
142   if (GetCurrentEventNo()<0) {
143     HLTWarning("no event selected, no data available");
144     return false;
145   }
146
147   // 2011-06-06 in order to support AliHLTReconstructor option 'ignore-hltout' for
148   // digits, the file name can be set to NULL, nothing done in that case 
149   if (!fDigitFileName) return false;
150
151   if (!fpDigitFile) {
152     fpDigitFile=new TFile(fDigitFileName);
153   }
154   if (!fpDigitFile) return false;
155   if (fpDigitFile->IsZombie()) return false;
156
157   fpDigitFile->GetObject("rawhltout", fpDigitTree);
158   if (!fpDigitTree) {
159     return false;
160   }
161
162   vector<int> equipments;
163   int iTotal=0;
164   TIter iter(fpDigitTree->GetListOfBranches());
165   while (TBranch* br=(TBranch *)iter.Next()) {
166     iTotal++;
167     TString bname=br->GetName();
168     if (bname.BeginsWith("HLT_") && bname.EndsWith(".ddl")) {
169       bname.ReplaceAll("HLT_", "");
170       bname.ReplaceAll(".ddl", "");
171       if (bname.IsDigit()) {
172         int equipment=bname.Atoi();
173         int index=equipment-AliDAQ::DdlIDOffset("HLT");
174         if (index>=0 && strcmp(br->GetName(), AliDAQ::DdlFileName("HLT", index))==0) {
175           equipments.push_back(equipment);
176         } else {
177           HLTWarning("equipment mismatch: skipping id %d (out of range of HLT equipments)", equipment);
178         }
179       }
180     }
181   }
182   HLTDebug("found %d HLT branche(s) out of %d", equipments.size(), iTotal);
183
184   if (equipments.size()>0) {
185     fNofDDLs=equipments.size();
186     fpEquipments=new int[fNofDDLs];
187     fppDigitArrays=new TArrayC*[fNofDDLs];
188     if (fpEquipments && fppDigitArrays) {
189       memcpy(fpEquipments, &equipments[0], fNofDDLs*sizeof(int));
190       int i=0;
191       for (i=0; i<fNofDDLs; i++) {
192         fppDigitArrays[i]=NULL;
193         fpDigitTree->SetBranchAddress(AliDAQ::DdlFileName("HLT", fpEquipments[i]-AliDAQ::DdlIDOffset("HLT")), &fppDigitArrays[i]);
194       }
195       fpDigitTree->GetEvent(GetCurrentEventNo());
196       for (i=0; i<fNofDDLs; i++) {
197         if (fppDigitArrays[i]) {
198           HLTDebug("branch %s: %d byte(s)", AliDAQ::DdlFileName("HLT", fpEquipments[i]-AliDAQ::DdlIDOffset("HLT")), fppDigitArrays[i]);
199         }
200       }
201       result=true;
202     } else {
203       fNofDDLs=0;
204     }
205   }
206   return result;
207 }
208
209 int AliHLTOUTDigitReader::CloseTree()
210 {
211   // Cleanup tree and data arrays.
212   int iResult=0;
213   for (int i=0; i<fNofDDLs; i++) {
214     if (fppDigitArrays[i]) delete fppDigitArrays[i];
215     fppDigitArrays[i]=NULL;
216   }
217   if (fppDigitArrays) delete[] fppDigitArrays;
218   fppDigitArrays=NULL;
219   if (fpEquipments) delete[] fpEquipments;
220   fpEquipments=NULL;
221   if (fpDigitTree) delete fpDigitTree;
222   fpDigitTree=NULL;
223   fNofDDLs=0;
224   fCurrentLink=-1;
225
226   return iResult;
227 }
228
229 void AliHLTOUTDigitReader::SetParam(TTree* /*pDigitTree*/, int event)
230 {
231   // set parameter for this HLTOUT instance
232   // The function is for internal use only in conjunction with the
233   // AliHLTOUT::New() functions.
234
235   // TODO: here we have the implemented to correct loading of
236   // the digit file from the run loader. At time of writing
237   // (Jun 2008) the AliLoader for the HLT is missing in the AliRoot
238   // framework.
239   fEvent=event;
240 }
241
242 void AliHLTOUTDigitReader::SetParam(const char* filename, int event)
243 {
244   // set parameter for this HLTOUT instance
245   // The function is for internal use only in conjunction with the
246   // AliHLTOUT::New() functions.
247
248   if (filename && filename[0]!=0) {
249     fDigitFileName=filename;
250   } else {
251     HLTWarning("no valid digit file provided, using default file %s", fDigitFileName.Data());
252   }
253   fEvent=event;
254 }