]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/rec/AliHLTOUTDigitReader.cxx
- HLT simulation writes digit data in addition to raw data
[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)
39   :
40   AliHLTOUTHomerCollection(event, pEsdManager),
41   fpDigitFile(NULL),
42   fpDigitTree(NULL),
43   fMinDDL(-1),
44   fMaxDDL(-1),
45   fppDigitArrays(NULL),
46   fpEquipments(NULL),
47   fNofDDLs(0),
48   fCurrent(-1)
49 {
50   // see header file for class documentation
51   // or
52   // refer to README to build package
53   // or
54   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
55 }
56
57 AliHLTOUTDigitReader::~AliHLTOUTDigitReader()
58 {
59   // see header file for class documentation
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   // see header file for class documentation
72   if (!fppDigitArrays && (!ReadArrays() || !fppDigitArrays))
73     return kFALSE;
74
75   if (fCurrent>=fNofDDLs)
76     return kFALSE;
77
78   while (++fCurrent<fNofDDLs) {
79     if (fMinDDL>-1 && fMinDDL>fpEquipments[fCurrent]) continue;
80     if (fMaxDDL>-1 && fMaxDDL<fpEquipments[fCurrent]) continue;
81     if (fppDigitArrays[fCurrent]->GetSize()>(int)sizeof(AliRawDataHeader)) {
82       data=reinterpret_cast<UChar_t*>(fppDigitArrays[fCurrent]->GetArray());
83       data+=sizeof(AliRawDataHeader);
84       return kTRUE;
85     }
86   }
87   return kFALSE;
88 }
89
90 int AliHLTOUTDigitReader::Reset()
91 {
92   // see header file for class documentation
93   fCurrent=-1;
94   return 0;
95 }
96
97 int AliHLTOUTDigitReader::GetDataSize()
98 {
99   // see header file for class documentation
100   if (fCurrent>=0 && fCurrent<fNofDDLs && fppDigitArrays) {
101     return fppDigitArrays[fCurrent]->GetSize()-sizeof(AliRawDataHeader);
102   }
103   return 0;
104 }
105
106 const AliRawDataHeader* AliHLTOUTDigitReader::GetDataHeader()
107 {
108   // see header file for class documentation
109   if (fCurrent>=0 && fCurrent<fNofDDLs && fppDigitArrays) {
110     return reinterpret_cast<AliRawDataHeader*>(fppDigitArrays[fCurrent]->GetArray());
111   }
112   return NULL;
113 }
114
115 void AliHLTOUTDigitReader::SelectEquipment(int /*equipmentType*/, int minEquipmentId, int maxEquipmentId)
116 {
117   // see header file for class documentation
118   fMinDDL=minEquipmentId;
119   fMaxDDL=maxEquipmentId;
120 }
121
122 int AliHLTOUTDigitReader::GetEquipmentId()
123 {
124   // see header file for class documentation
125   if (fCurrent>=0 && fCurrent<fNofDDLs && fpEquipments) {
126     return fpEquipments[fCurrent];
127   }
128   return -1;
129 }
130
131 bool AliHLTOUTDigitReader::ReadArrays()
132 {
133   // see header file for class documentation
134   bool result=false;
135
136   if (GetCurrentEventNo()<0) {
137     HLTWarning("no event selected, no data available");
138     return false;
139   }
140
141   if (!fpDigitFile) {
142     fpDigitFile=new TFile("HLT.Digits.root");
143   }
144   if (!fpDigitFile) return false;
145   if (fpDigitFile->IsZombie()) return false;
146
147   fpDigitFile->GetObject("rawhltout", fpDigitTree);
148   if (!fpDigitTree) {
149     return false;
150   }
151
152   vector<int> equipments;
153   int iTotal=0;
154   TIter iter(fpDigitTree->GetListOfBranches());
155   while (TBranch* br=(TBranch *)iter.Next()) {
156     iTotal++;
157     TString bname=br->GetName();
158     if (bname.BeginsWith("HLT_") && bname.EndsWith(".ddl")) {
159       bname.ReplaceAll("HLT_", "");
160       bname.ReplaceAll(".ddl", "");
161       if (bname.IsDigit()) {
162         int equipment=bname.Atoi();
163         int index=equipment-AliDAQ::DdlIDOffset("HLT");
164         if (index>=0 && strcmp(br->GetName(), AliDAQ::DdlFileName("HLT", index))==0) {
165           equipments.push_back(equipment);
166         } else {
167           HLTWarning("equipment mismatch: skipping id %d (out of range of HLT equipments)", equipment);
168         }
169       }
170     }
171   }
172   HLTDebug("found %d HLT branche(s) out of %d", equipments.size(), iTotal);
173
174   if (equipments.size()>0) {
175     fNofDDLs=equipments.size();
176     fpEquipments=new int[fNofDDLs];
177     fppDigitArrays=new TArrayC*[fNofDDLs];
178     if (fpEquipments && fppDigitArrays) {
179       memcpy(fpEquipments, &equipments[0], fNofDDLs*sizeof(int));
180       int i=0;
181       for (i=0; i<fNofDDLs; i++) {
182         fppDigitArrays[i]=NULL;
183         fpDigitTree->SetBranchAddress(AliDAQ::DdlFileName("HLT", fpEquipments[i]-AliDAQ::DdlIDOffset("HLT")), &fppDigitArrays[i]);
184       }
185       fpDigitTree->GetEvent(GetCurrentEventNo());
186       for (i=0; i<fNofDDLs; i++) {
187         if (fppDigitArrays[i]) {
188           HLTDebug("branch %s: %d byte(s)", AliDAQ::DdlFileName("HLT", fpEquipments[i]-AliDAQ::DdlIDOffset("HLT")), fppDigitArrays[i]);
189         }
190       }
191       result=true;
192     } else {
193       fNofDDLs=0;
194     }
195   }
196   return result;
197 }
198
199 int AliHLTOUTDigitReader::CloseTree()
200 {
201   // see header file for class documentation
202   int iResult=0;
203   if (fppDigitArrays) delete[] fppDigitArrays;
204   fppDigitArrays=NULL;
205   if (fpEquipments) delete[] fpEquipments;
206   fpEquipments=NULL;
207   if (fpDigitTree) delete fpDigitTree;
208   fpDigitTree=NULL;
209   fNofDDLs=0;
210   fCurrent=-1;
211
212   return iResult;
213 }