]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/BASE/test/testAliHLTEventDDLBackwardCompatibility.C
Adding tests for backward compatibility of AliHLTReadoutList class.
[u/mrichter/AliRoot.git] / HLT / BASE / test / testAliHLTEventDDLBackwardCompatibility.C
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: Artur Szostak <artursz@iafrica.com>                   *
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   testAliHLTEventDDLBackwardCompatibility.C
20 /// @author Artur Szostak <artursz@iafrica.com>
21 /// @date   25 June 2010
22 /// @brief  Test program for backward compatibility of the AliHLTEventDDL structure.
23 ///
24
25 #if defined(__CINT__) && (! defined(__MAKECINT__))
26 #error This macro must be compiled. Try running as testAliHLTEventDDLBackwardCompatibility.C++.
27 #endif
28
29 #if !defined(__CINT__) || defined(__MAKECINT__)
30 #include "AliHLTDataTypes.h"
31 #include "AliHLTReadoutList.h"
32 #include "AliHLTDAQ.h"
33 #include "TRandom3.h"
34 #include "TString.h"
35 #include "TFile.h"
36 #include "TTree.h"
37 #include "Riostream.h"
38 #include <vector>
39 #include <algorithm>
40 #endif
41
42 /**
43  * Tests to see if the AliHLTReadoutList class handles both AliHLTEventDDLV0
44  * and AliHLTEventDDLV1 correctly.
45  */
46 bool CheckReadoutListConvertedCorrectly()
47 {
48         // Initialise the different versions of the structures so that every detector
49         // has only its first DDL set.
50         union
51         {
52                 AliHLTEventDDL eventddlV0;
53                 AliHLTEventDDLV0 bitsV0;
54         };
55         bitsV0.fCount = gkAliHLTDDLListSizeV0;
56         if (gkAliHLTDDLListSizeV0 != 30)
57         {
58                 cerr << "ERROR: gkAliHLTDDLListSizeV0 has a value of " << gkAliHLTDDLListSizeV0
59                         << " but expected a value of 30." << endl;
60                 return false;
61         }
62         bitsV0.fList[0] = 0x00000001; // kITSSPD
63         bitsV0.fList[1] = 0x00000001; // kITSSDD
64         bitsV0.fList[2] = 0x00000001; // kITSSSD
65         bitsV0.fList[3] = 0x00000001; // kTPC
66         bitsV0.fList[4] = 0x00000000; // kTPC
67         bitsV0.fList[5] = 0x00000000; // kTPC
68         bitsV0.fList[6] = 0x00000000; // kTPC
69         bitsV0.fList[7] = 0x00000000; // kTPC
70         bitsV0.fList[8] = 0x00000000; // kTPC
71         bitsV0.fList[9] = 0x00000000; // kTPC
72         bitsV0.fList[10] = 0x00000000; // kTPC
73         bitsV0.fList[11] = 0x00000001; // kTRD
74         bitsV0.fList[12] = 0x00000001; // kTOF
75         bitsV0.fList[13] = 0x00000000; // kTOF
76         bitsV0.fList[14] = 0x00000000; // kTOF
77         bitsV0.fList[15] = 0x00000001; // kHMPID
78         bitsV0.fList[16] = 0x00000001; // kPHOS
79         bitsV0.fList[17] = 0x00000001; // kCPV
80         bitsV0.fList[18] = 0x00000001; // kPMD
81         bitsV0.fList[19] = 0x00000001; // kMUONTRK
82         bitsV0.fList[20] = 0x00000001; // kMUONTRG
83         bitsV0.fList[21] = 0x00000001; // kFMD
84         bitsV0.fList[22] = 0x00000001; // kT0
85         bitsV0.fList[23] = 0x00000001; // kV0
86         bitsV0.fList[24] = 0x00000001; // kZDC
87         bitsV0.fList[25] = 0x00000001; // kACORDE
88         bitsV0.fList[26] = 0x00000001; // kTRG
89         bitsV0.fList[27] = 0x00000001; // kEMCAL
90         bitsV0.fList[28] = 0x00000001; // kDAQTEST
91         bitsV0.fList[29] = 0x00000001; // kHLT
92         
93         union
94         {
95                 AliHLTEventDDL eventddlV1;
96                 AliHLTEventDDLV1 bitsV1;
97         };
98         bitsV1.fCount = gkAliHLTDDLListSizeV1;
99         if (gkAliHLTDDLListSizeV1 != 31)
100         {
101                 cerr << "ERROR: gkAliHLTDDLListSizeV1 has a value of " << gkAliHLTDDLListSizeV1
102                         << " but expected a value of 31." << endl;
103                 return false;
104         }
105         bitsV1.fList[0] = 0x00000001; // kITSSPD
106         bitsV1.fList[1] = 0x00000001; // kITSSDD
107         bitsV1.fList[2] = 0x00000001; // kITSSSD
108         bitsV1.fList[3] = 0x00000001; // kTPC
109         bitsV1.fList[4] = 0x00000000; // kTPC
110         bitsV1.fList[5] = 0x00000000; // kTPC
111         bitsV1.fList[6] = 0x00000000; // kTPC
112         bitsV1.fList[7] = 0x00000000; // kTPC
113         bitsV1.fList[8] = 0x00000000; // kTPC
114         bitsV1.fList[9] = 0x00000000; // kTPC
115         bitsV1.fList[10] = 0x00000000; // kTPC
116         bitsV1.fList[11] = 0x00000001; // kTRD
117         bitsV1.fList[12] = 0x00000001; // kTOF
118         bitsV1.fList[13] = 0x00000000; // kTOF
119         bitsV1.fList[14] = 0x00000000; // kTOF
120         bitsV1.fList[15] = 0x00000001; // kHMPID
121         bitsV1.fList[16] = 0x00000001; // kPHOS
122         bitsV1.fList[17] = 0x00000001; // kCPV
123         bitsV1.fList[18] = 0x00000001; // kPMD
124         bitsV1.fList[19] = 0x00000001; // kMUONTRK
125         bitsV1.fList[20] = 0x00000001; // kMUONTRG
126         bitsV1.fList[21] = 0x00000001; // kFMD
127         bitsV1.fList[22] = 0x00000001; // kT0
128         bitsV1.fList[23] = 0x00000001; // kV0
129         bitsV1.fList[24] = 0x00000001; // kZDC
130         bitsV1.fList[25] = 0x00000001; // kACORDE
131         bitsV1.fList[26] = 0x00000001; // kTRG
132         bitsV1.fList[27] = 0x00000001; // kEMCAL
133         bitsV1.fList[28] = 0x00000000; // kEMCAL
134         bitsV1.fList[29] = 0x00000001; // kDAQTEST
135         bitsV1.fList[30] = 0x00000001; // kHLT
136         
137         AliHLTReadoutList rlV0(eventddlV0);
138         AliHLTReadoutList rlV1(eventddlV1);
139         
140         // Check that for both readout list versions only the first DDLs are
141         // enabled as expected.
142         for (Int_t i = 0; i < AliHLTDAQ::NumberOfDetectors(); ++i)
143         for (Int_t j = 0; j < AliHLTDAQ::NumberOfDdls(i); ++j)
144         {
145                 Int_t ddlid = AliHLTDAQ::DdlIDOffset(i) | (j & 0xFF);
146                 if (j == 0)
147                 {
148                         if (rlV0.IsDDLDisabled(ddlid))
149                         {
150                                 cerr << "ERROR: The first DDL for detector " << AliHLTDAQ::DetectorName(i)
151                                         << " was not enabled for readout list initialised from AliHLTEventDDLV0."
152                                         << endl;
153                                 return false;
154                         }
155                         if (rlV1.IsDDLDisabled(ddlid))
156                         {
157                                 cerr << "ERROR: The first DDL for detector " << AliHLTDAQ::DetectorName(i)
158                                         << " was not enabled for readout list initialised from AliHLTEventDDLV1."
159                                         << endl;
160                                 return false;
161                         }
162                 }
163                 else
164                 {
165                         if (rlV0.IsDDLEnabled(ddlid))
166                         {
167                                 cerr << "ERROR: DDL " << ddlid << " for detector " << AliHLTDAQ::DetectorName(i)
168                                         << " was marked enabled for readout list initialised from AliHLTEventDDLV0."
169                                         << endl;
170                                 return false;
171                         }
172                         if (rlV1.IsDDLEnabled(ddlid))
173                         {
174                                 cerr << "ERROR: DDL " << ddlid << " for detector " << AliHLTDAQ::DetectorName(i)
175                                         << " was marked enabled for readout list initialised from AliHLTEventDDLV1."
176                                         << endl;
177                                 return false;
178                         }
179                 }
180         }
181         
182         // Check that the internal buffers are identical.
183         if (rlV0.BufferSize() != rlV1.BufferSize())
184         {
185                 cerr << "ERROR: Buffer sizes for readout lists are different: rlV0.BufferSize() = "
186                         << rlV0.BufferSize() << ", rlV1.BufferSize() = " << rlV1.BufferSize() << endl;
187                 return false;
188         }
189         if (memcmp(rlV0.Buffer(), rlV1.Buffer(), rlV0.BufferSize()) != 0)
190         {
191                 cerr << "ERROR: Buffers for the two readout list versions are different." << endl;
192                 return false;
193         }
194         
195         return true;
196 }
197
198 /**
199  * Tests to see if reading old AliHLTReadoutList versions from root files works.
200  * \param filename  The name of the file generated by GenerateReadoutListFile.C,
201  *      which contains the readout list objects to test.
202  */
203 bool CheckReadingOldFormat(const char* filename = "oldAliHLTReadoutListFormat.root")
204 {
205         TFile file(filename, "READ");
206         
207         // Load the readout list objects.
208         AliHLTReadoutList* rl[5] = {NULL, NULL, NULL, NULL, NULL};
209         for (int i = 0; i < 5; ++i)
210         {
211                 char name[1024];
212                 sprintf(name, "readoutlist%d", i+1);
213                 rl[i] = dynamic_cast<AliHLTReadoutList*>(file.Get(name));
214                 if (rl[i] == NULL)
215                 {
216                         cerr << "ERROR: Could not get object '" << name
217                                 << "' from file '" << filename << "'." << endl;
218                         return false;
219                 }
220         }
221         
222         // Now load the tree and see that the objects are the same as the readout
223         // list objects stored directly to the TFile.
224         const char* treename = "rltree";
225         TTree* tree = dynamic_cast<TTree*>(file.Get(treename));
226         if (tree == NULL)
227         {
228                 cerr << "ERROR: Could not get TTree '" << treename
229                         << "' from file '" << filename << "'." << endl;
230                 return false;
231         }
232         AliHLTReadoutList* r = new AliHLTReadoutList;
233         tree->SetBranchAddress("readoutlist", &r);
234         for (int i = 0; i < 5; ++i)
235         {
236                 tree->GetEvent(i);
237                 if (r->BufferSize() != rl[i]->BufferSize())
238                 {
239                         cerr << "ERROR: readoutlist" << i+1
240                                 << " and the one from the TTree have different sizes."
241                                 << endl;
242                         return false;
243                 }
244                 if (memcmp(r->Buffer(), rl[i]->Buffer(), r->BufferSize()) != 0)
245                 {
246                         cerr << "ERROR: readoutlist" << i+1
247                                 << " and the one from the TTree are different."
248                                 << endl;
249                         return false;
250                 }
251         }
252         
253         // Now check each readout list individually.
254         typedef AliHLTReadoutList RL;
255         Int_t alwaysoff = RL::kITSSPD
256                 | RL::kITSSDD
257                 | RL::kITSSSD
258                 | RL::kTPC
259                 | RL::kTRD
260                 | RL::kTOF
261                 | RL::kHMPID
262                 | RL::kPHOS
263                 | RL::kCPV
264                 | RL::kPMD
265                 | RL::kMUONTRK
266                 | RL::kMUONTRG
267                 | RL::kFMD
268                 | RL::kT0
269                 | RL::kV0
270                 | RL::kZDC
271                 | RL::kACORDE;
272         
273         // We will need to try and set the missing EMCAL DDL bits.
274         // Otherwise the DetectorEnabled and DetectorDisabled methods will not
275         // give the expected answers.
276         for (int i = 24; i < AliHLTDAQ::NumberOfDdls("EMCAL"); ++i)
277         {
278                 for (int j = 1; j <= 3; ++j)
279                 {
280                         Int_t ddlid = AliHLTDAQ::DdlIDOffset("EMCAL") | (i & 0xFF);
281                         rl[j]->EnableDDLBit(ddlid);
282                 }
283         }
284         
285         int rlnum = 0;
286         if (not (rl[rlnum]->DetectorEnabled(RL::kTRG) and
287                  rl[rlnum]->DetectorDisabled(alwaysoff | RL::kEMCAL | RL::kDAQTEST | RL::kHLT)
288            ))
289         {
290                 rl[0]->Print();
291                 cerr << "ERROR: readoutlist" << rlnum+1 << " does not have the correct bits set." << endl;
292                 rl[rlnum]->Print();
293                 return false;
294         }
295         rlnum = 1;
296         if (not (rl[rlnum]->DetectorEnabled(RL::kTRG | RL::kEMCAL) and
297                  rl[rlnum]->DetectorDisabled(alwaysoff | RL::kDAQTEST | RL::kHLT)
298            ))
299         {
300                 cerr << "ERROR: readoutlist" << rlnum+1 << " does not have the correct bits set." << endl;
301                 rl[rlnum]->Print();
302                 return false;
303         }
304         rlnum = 2;
305         if (not (rl[rlnum]->DetectorEnabled(RL::kTRG | RL::kEMCAL | RL::kDAQTEST) and
306                  rl[rlnum]->DetectorDisabled(alwaysoff | RL::kHLT)
307            ))
308         {
309                 cerr << "ERROR: readoutlist" << rlnum+1 << " does not have the correct bits set." << endl;
310                 rl[rlnum]->Print();
311                 return false;
312         }
313         rlnum = 3;
314         if (not (rl[rlnum]->DetectorEnabled(RL::kTRG | RL::kEMCAL | RL::kDAQTEST | RL::kHLT) and
315                  rl[rlnum]->DetectorDisabled(alwaysoff)
316            ))
317         {
318                 cerr << "ERROR: readoutlist" << rlnum+1 << " does not have the correct bits set." << endl;
319                 rl[rlnum]->Print();
320                 return false;
321         }
322         rlnum = 4;
323         if (not (rl[rlnum]->DetectorEnabled(RL::kTRG | RL::kDAQTEST | RL::kHLT) and
324                  rl[rlnum]->DetectorDisabled(alwaysoff | RL::kEMCAL)
325            ))
326         {
327                 cerr << "ERROR: readoutlist" << rlnum+1 << " does not have the correct bits set." << endl;
328                 rl[rlnum]->Print();
329                 return false;
330         }
331         
332         return true;
333 }
334
335 /**
336  * Runs the unit tests for backward compatibility of AliHLTEventDDL succeeded.
337  * \returns true if the tests were passed and false otherwise.
338  */
339 bool testAliHLTEventDDLBackwardCompatibility()
340 {
341         if (not CheckReadoutListConvertedCorrectly()) return false;
342         if (not CheckReadingOldFormat()) return false;
343         return true;
344 }
345
346 #ifndef __MAKECINT__
347
348 int main(int /*argc*/, const char** /*argv*/)
349 {
350         bool resultOk = testAliHLTEventDDLBackwardCompatibility();
351         if (not resultOk) return 1;
352         return 0;
353 }
354
355 #endif // __MAKECINT__