]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/rec/AliHLTOUTDigitReader.cxx
f4e58f105093602f2b2539cb25deeb3757f467a0
[u/mrichter/AliRoot.git] / HLT / rec / AliHLTOUTDigitReader.cxx
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
38 AliHLTOUTDigitReader::AliHLTOUTDigitReader(int event, AliHLTEsdManager* pEsdManager, const char* digitFile)
39   :
40   AliHLTOUTHomerCollection(event, pEsdManager),
41   fDigitFileName(digitFile),
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
142   if (!fpDigitFile) {
143     fpDigitFile=new TFile(fDigitFileName);
144   }
145   if (!fpDigitFile) return false;
146   if (fpDigitFile->IsZombie()) return false;
147
148   fpDigitFile->GetObject("rawhltout", fpDigitTree);
149   if (!fpDigitTree) {
150     return false;
151   }
152
153   vector<int> equipments;
154   int iTotal=0;
155   TIter iter(fpDigitTree->GetListOfBranches());
156   while (TBranch* br=(TBranch *)iter.Next()) {
157     iTotal++;
158     TString bname=br->GetName();
159     if (bname.BeginsWith("HLT_") && bname.EndsWith(".ddl")) {
160       bname.ReplaceAll("HLT_", "");
161       bname.ReplaceAll(".ddl", "");
162       if (bname.IsDigit()) {
163         int equipment=bname.Atoi();
164         int index=equipment-AliDAQ::DdlIDOffset("HLT");
165         if (index>=0 && strcmp(br->GetName(), AliDAQ::DdlFileName("HLT", index))==0) {
166           equipments.push_back(equipment);
167         } else {
168           HLTWarning("equipment mismatch: skipping id %d (out of range of HLT equipments)", equipment);
169         }
170       }
171     }
172   }
173   HLTDebug("found %d HLT branche(s) out of %d", equipments.size(), iTotal);
174
175   if (equipments.size()>0) {
176     fNofDDLs=equipments.size();
177     fpEquipments=new int[fNofDDLs];
178     fppDigitArrays=new TArrayC*[fNofDDLs];
179     if (fpEquipments && fppDigitArrays) {
180       memcpy(fpEquipments, &equipments[0], fNofDDLs*sizeof(int));
181       int i=0;
182       for (i=0; i<fNofDDLs; i++) {
183         fppDigitArrays[i]=NULL;
184         fpDigitTree->SetBranchAddress(AliDAQ::DdlFileName("HLT", fpEquipments[i]-AliDAQ::DdlIDOffset("HLT")), &fppDigitArrays[i]);
185       }
186       fpDigitTree->GetEvent(GetCurrentEventNo());
187       for (i=0; i<fNofDDLs; i++) {
188         if (fppDigitArrays[i]) {
189           HLTDebug("branch %s: %d byte(s)", AliDAQ::DdlFileName("HLT", fpEquipments[i]-AliDAQ::DdlIDOffset("HLT")), fppDigitArrays[i]);
190         }
191       }
192       result=true;
193     } else {
194       fNofDDLs=0;
195     }
196   }
197   return result;
198 }
199
200 int AliHLTOUTDigitReader::CloseTree()
201 {
202   // see header file for class documentation
203   int iResult=0;
204   if (fppDigitArrays) delete[] fppDigitArrays;
205   fppDigitArrays=NULL;
206   if (fpEquipments) delete[] fpEquipments;
207   fpEquipments=NULL;
208   if (fpDigitTree) delete fpDigitTree;
209   fpDigitTree=NULL;
210   fNofDDLs=0;
211   fCurrent=-1;
212
213   return iResult;
214 }
215
216 void AliHLTOUTDigitReader::SetParam(TTree* /*pDigitTree*/, int event)
217 {
218   // see header file for class documentation
219
220   // TODO: here we have the implemented to correct loading of
221   // the digit file from the run loader. At time of writing
222   // (Jun 2008) the AliLoader for the HLT is missing in the AliRoot
223   // framework.
224   fEvent=event;
225 }