3 /**************************************************************************
4 * This file is property of and copyright by the ALICE HLT Project *
5 * ALICE Experiment at CERN, All rights reserved. *
7 * Primary Authors: Artur Szostak <artursz@iafrica.com> *
8 * for The ALICE HLT Project. *
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 **************************************************************************/
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.
25 #if defined(__CINT__) && (! defined(__MAKECINT__))
26 #error This macro must be compiled. Try running as testAliHLTEventDDLBackwardCompatibility.C++.
29 #if !defined(__CINT__) || defined(__MAKECINT__)
30 #include "AliHLTDataTypes.h"
31 #include "AliHLTReadoutList.h"
32 #include "AliHLTDAQ.h"
33 #include "AliHLTComponent.h"
34 #include "AliHLTCDHWrapper.h"
39 #include "Riostream.h"
45 * Tests to see if the AliHLTReadoutList class handles both AliHLTEventDDLV0
46 * and AliHLTEventDDLV1 correctly.
48 bool CheckReadoutListConvertedCorrectly()
50 // Initialise the different versions of the structures so that every detector
51 // has only its first DDL set.
54 AliHLTEventDDL eventddlV0;
55 AliHLTEventDDLV0 bitsV0;
57 bitsV0.fCount = gkAliHLTDDLListSizeV0;
58 if (gkAliHLTDDLListSizeV0 != 30)
60 cerr << "ERROR: gkAliHLTDDLListSizeV0 has a value of " << gkAliHLTDDLListSizeV0
61 << " but expected a value of 30." << endl;
64 bitsV0.fList[0] = 0x00000001; // kITSSPD
65 bitsV0.fList[1] = 0x00000001; // kITSSDD
66 bitsV0.fList[2] = 0x00000001; // kITSSSD
67 bitsV0.fList[3] = 0x00000001; // kTPC
68 bitsV0.fList[4] = 0x00000000; // kTPC
69 bitsV0.fList[5] = 0x00000000; // kTPC
70 bitsV0.fList[6] = 0x00000000; // kTPC
71 bitsV0.fList[7] = 0x00000000; // kTPC
72 bitsV0.fList[8] = 0x00000000; // kTPC
73 bitsV0.fList[9] = 0x00000000; // kTPC
74 bitsV0.fList[10] = 0x00000000; // kTPC
75 bitsV0.fList[11] = 0x00000001; // kTRD
76 bitsV0.fList[12] = 0x00000001; // kTOF
77 bitsV0.fList[13] = 0x00000000; // kTOF
78 bitsV0.fList[14] = 0x00000000; // kTOF
79 bitsV0.fList[15] = 0x00000001; // kHMPID
80 bitsV0.fList[16] = 0x00000001; // kPHOS
81 bitsV0.fList[17] = 0x00000001; // kCPV
82 bitsV0.fList[18] = 0x00000001; // kPMD
83 bitsV0.fList[19] = 0x00000001; // kMUONTRK
84 bitsV0.fList[20] = 0x00000001; // kMUONTRG
85 bitsV0.fList[21] = 0x00000001; // kFMD
86 bitsV0.fList[22] = 0x00000001; // kT0
87 bitsV0.fList[23] = 0x00000001; // kV0
88 bitsV0.fList[24] = 0x00000001; // kZDC
89 bitsV0.fList[25] = 0x00000001; // kACORDE
90 bitsV0.fList[26] = 0x00000001; // kTRG
91 bitsV0.fList[27] = 0x00000001; // kEMCAL
92 bitsV0.fList[28] = 0x00000001; // kDAQTEST
93 bitsV0.fList[29] = 0x00000001; // kHLT
97 AliHLTEventDDL eventddlV1;
98 AliHLTEventDDLV1 bitsV1;
100 bitsV1.fCount = gkAliHLTDDLListSizeV1;
101 if (gkAliHLTDDLListSizeV1 != 31)
103 cerr << "ERROR: gkAliHLTDDLListSizeV1 has a value of " << gkAliHLTDDLListSizeV1
104 << " but expected a value of 31." << endl;
107 bitsV1.fList[0] = 0x00000001; // kITSSPD
108 bitsV1.fList[1] = 0x00000001; // kITSSDD
109 bitsV1.fList[2] = 0x00000001; // kITSSSD
110 bitsV1.fList[3] = 0x00000001; // kTPC
111 bitsV1.fList[4] = 0x00000000; // kTPC
112 bitsV1.fList[5] = 0x00000000; // kTPC
113 bitsV1.fList[6] = 0x00000000; // kTPC
114 bitsV1.fList[7] = 0x00000000; // kTPC
115 bitsV1.fList[8] = 0x00000000; // kTPC
116 bitsV1.fList[9] = 0x00000000; // kTPC
117 bitsV1.fList[10] = 0x00000000; // kTPC
118 bitsV1.fList[11] = 0x00000001; // kTRD
119 bitsV1.fList[12] = 0x00000001; // kTOF
120 bitsV1.fList[13] = 0x00000000; // kTOF
121 bitsV1.fList[14] = 0x00000000; // kTOF
122 bitsV1.fList[15] = 0x00000001; // kHMPID
123 bitsV1.fList[16] = 0x00000001; // kPHOS
124 bitsV1.fList[17] = 0x00000001; // kCPV
125 bitsV1.fList[18] = 0x00000001; // kPMD
126 bitsV1.fList[19] = 0x00000001; // kMUONTRK
127 bitsV1.fList[20] = 0x00000001; // kMUONTRG
128 bitsV1.fList[21] = 0x00000001; // kFMD
129 bitsV1.fList[22] = 0x00000001; // kT0
130 bitsV1.fList[23] = 0x00000001; // kV0
131 bitsV1.fList[24] = 0x00000001; // kZDC
132 bitsV1.fList[25] = 0x00000001; // kACORDE
133 bitsV1.fList[26] = 0x00000001; // kTRG
134 bitsV1.fList[27] = 0x00000001; // kEMCAL
135 bitsV1.fList[28] = 0x00000000; // kEMCAL
136 bitsV1.fList[29] = 0x00000001; // kDAQTEST
137 bitsV1.fList[30] = 0x00000001; // kHLT
141 AliHLTEventDDL eventddlV2;
142 AliHLTEventDDLV2 bitsV2;
144 bitsV2.fCount = gkAliHLTDDLListSizeV2;
145 if (gkAliHLTDDLListSizeV2 != 32)
147 cerr << "ERROR: gkAliHLTDDLListSizeV2 has a value of " << gkAliHLTDDLListSizeV2
148 << " but expected a value of 32." << endl;
151 bitsV2.fList[0] = 0x00000001; // kITSSPD
152 bitsV2.fList[1] = 0x00000001; // kITSSDD
153 bitsV2.fList[2] = 0x00000001; // kITSSSD
154 bitsV2.fList[3] = 0x00000001; // kTPC
155 bitsV2.fList[4] = 0x00000000; // kTPC
156 bitsV2.fList[5] = 0x00000000; // kTPC
157 bitsV2.fList[6] = 0x00000000; // kTPC
158 bitsV2.fList[7] = 0x00000000; // kTPC
159 bitsV2.fList[8] = 0x00000000; // kTPC
160 bitsV2.fList[9] = 0x00000000; // kTPC
161 bitsV2.fList[10] = 0x00000000; // kTPC
162 bitsV2.fList[11] = 0x00000001; // kTRD
163 bitsV2.fList[12] = 0x00000001; // kTOF
164 bitsV2.fList[13] = 0x00000000; // kTOF
165 bitsV2.fList[14] = 0x00000000; // kTOF
166 bitsV2.fList[15] = 0x00000001; // kHMPID
167 bitsV2.fList[16] = 0x00000001; // kPHOS
168 bitsV2.fList[17] = 0x00000001; // kCPV
169 bitsV2.fList[18] = 0x00000001; // kPMD
170 bitsV2.fList[19] = 0x00000001; // kMUONTRK
171 bitsV2.fList[20] = 0x00000001; // kMUONTRG
172 bitsV2.fList[21] = 0x00000001; // kFMD
173 bitsV2.fList[22] = 0x00000001; // kT0
174 bitsV2.fList[23] = 0x00000001; // kV0
175 bitsV2.fList[24] = 0x00000001; // kZDC
176 bitsV2.fList[25] = 0x00000001; // kACORDE
177 bitsV2.fList[26] = 0x00000001; // kTRG
178 bitsV2.fList[27] = 0x00000001; // kEMCAL
179 bitsV2.fList[28] = 0x00000000; // kEMCAL
180 bitsV2.fList[29] = 0x00000001; // kDAQTEST
181 bitsV2.fList[30] = 0x00000000; // kAD
182 bitsV2.fList[31] = 0x00000001; // kHLT
184 AliHLTReadoutList rlV0(eventddlV0);
185 AliHLTReadoutList rlV1(eventddlV1);
186 AliHLTReadoutList rlV2(eventddlV2);
188 printf("eventddlV0.fCount: %i\n",eventddlV0.fCount);
189 printf("eventddlV1.fCount: %i\n",eventddlV1.fCount);
190 printf("eventddlV2.fCount: %i\n",eventddlV2.fCount);
191 printf("rlV0.BufferSize(): %i\n",rlV0.BufferSize());
192 printf("rlV1.BufferSize(): %i\n",rlV1.BufferSize());
193 printf("rlV2.BufferSize(): %i\n",rlV2.BufferSize());
194 printf("sizeof AliHLTEventDDLV0: %i\n", sizeof(AliHLTEventDDLV0));
195 printf("sizeof AliHLTEventDDLV1: %i\n", sizeof(AliHLTEventDDLV1));
196 printf("sizeof AliHLTEventDDLV2: %i\n", sizeof(AliHLTEventDDLV2));
197 printf("sizeof AliHLTEventDDL: %i\n", sizeof(AliHLTEventDDL));
198 printf("sizeof rlV0: %i, sizeof eventddlV0: %i\n",sizeof(rlV2), sizeof(eventddlV2));
199 printf("sizeof rlV1: %i, sizeof eventddlV1: %i\n",sizeof(rlV2), sizeof(eventddlV2));
200 printf("sizeof rlV2: %i, sizeof eventddlV2: %i\n",sizeof(rlV2), sizeof(eventddlV2));
202 // Check that for both readout list versions only the first DDLs are
203 // enabled as expected.
204 for (Int_t i = 0; i < AliHLTDAQ::NumberOfDetectors(); ++i)
205 for (Int_t j = 0; j < AliHLTDAQ::NumberOfDdls(i); ++j)
207 Int_t ddlid = AliHLTDAQ::DdlIDOffset(i) | (j & 0xFF);
210 if (rlV0.IsDDLDisabled(ddlid) && i!=20)
212 cerr << "ERROR: The first DDL for detector " << AliHLTDAQ::DetectorName(i)
213 << " was not enabled for readout list initialised from AliHLTEventDDLV0."
215 printf("i=%i, j=%i, ddlid=%i, det=%s, enabled: %i\n",i,j,ddlid,AliHLTDAQ::DetectorName(i), (rlV0.IsDDLDisabled(ddlid))?0:1);
218 if (rlV1.IsDDLDisabled(ddlid) && i!=20)
220 cerr << "ERROR: The first DDL for detector " << AliHLTDAQ::DetectorName(i)
221 << " was not enabled for readout list initialised from AliHLTEventDDLV1."
223 printf("i=%i, j=%i, ddlid=%i, det=%s, enabled: %i\n",i,j,ddlid,AliHLTDAQ::DetectorName(i), (rlV1.IsDDLDisabled(ddlid))?0:1);
226 if (rlV2.IsDDLDisabled(ddlid))
228 cerr << "ERROR: The first DDL for detector " << AliHLTDAQ::DetectorName(i)
229 << " was not enabled for readout list initialised from AliHLTEventDDLV2."
231 printf("i=%i, j=%i, ddlid=%i, det=%s, enabled: %i\n",i,j,ddlid,AliHLTDAQ::DetectorName(i), (rlV2.IsDDLDisabled(ddlid))?0:1);
237 if (rlV0.IsDDLEnabled(ddlid))
239 cerr << "ERROR: DDL " << ddlid << " for detector " << AliHLTDAQ::DetectorName(i)
240 << " was marked enabled for readout list initialised from AliHLTEventDDLV0."
244 if (rlV1.IsDDLEnabled(ddlid))
246 cerr << "ERROR: DDL " << ddlid << " for detector " << AliHLTDAQ::DetectorName(i)
247 << " was marked enabled for readout list initialised from AliHLTEventDDLV1."
251 if (rlV2.IsDDLEnabled(ddlid))
253 cerr << "ERROR: DDL " << ddlid << " for detector " << AliHLTDAQ::DetectorName(i)
254 << " was marked enabled for readout list initialised from AliHLTEventDDLV2."
261 // Check that the internal buffers are identical.
262 if (rlV0.BufferSize() != rlV1.BufferSize())
264 cerr << "ERROR: Buffer sizes for readout lists are different: rlV0.BufferSize() = "
265 << rlV0.BufferSize() << ", rlV1.BufferSize() = " << rlV1.BufferSize() << endl;
268 if (memcmp(rlV0.Buffer(), rlV1.Buffer(), rlV0.BufferSize()) != 0)
270 cerr << "ERROR: Buffers for the two readout list versions are different." << endl;
274 if (rlV1.BufferSize() != rlV2.BufferSize())
276 cerr << "ERROR: Buffer sizes for readout lists are different: rlV1.BufferSize() = "
277 << rlV1.BufferSize() << ", rlV2.BufferSize() = " << rlV2.BufferSize() << endl;
280 if (memcmp(rlV1.Buffer(), rlV2.Buffer(), rlV1.BufferSize()) != 0)
282 cerr << "ERROR: Buffers for the two readout list versions are different." << endl;
290 * Tests to see if AliHLTComponent::ExtractTriggerData recognises the old format
291 * of AliHLTEventTriggerData and handles it correctly.
293 bool CheckHandlingOfOldAliHLTEventTriggerData()
295 // Prepare old structure format.
296 struct AliHLTEventTriggerDataV0
298 AliHLTUInt8_t fAttributes[gkAliHLTBlockDAttributeCount];
299 AliHLTUInt64_t fHLTStatus;
300 AliHLTUInt32_t fCommonHeaderWordCnt;
301 AliHLTUInt32_t fCommonHeader[gkAliHLTCommonHeaderCount];
302 AliHLTEventDDLV0 fReadoutListV0;
304 AliHLTEventTriggerDataV0 eventTrigData;
305 memset(&eventTrigData, 0x0, sizeof(eventTrigData));
306 eventTrigData.fCommonHeaderWordCnt = 8;
307 eventTrigData.fHLTStatus = 0x123;
308 eventTrigData.fReadoutListV0.fCount = gkAliHLTDDLListSizeV0;
309 eventTrigData.fReadoutListV0.fList[0] = 0x00000001; // kITSSPD
310 eventTrigData.fReadoutListV0.fList[1] = 0x00000001; // kITSSDD
311 eventTrigData.fReadoutListV0.fList[2] = 0x00000001; // kITSSSD
312 eventTrigData.fReadoutListV0.fList[3] = 0x00000001; // kTPC
313 eventTrigData.fReadoutListV0.fList[4] = 0x00000000; // kTPC
314 eventTrigData.fReadoutListV0.fList[5] = 0x00000000; // kTPC
315 eventTrigData.fReadoutListV0.fList[6] = 0x00000000; // kTPC
316 eventTrigData.fReadoutListV0.fList[7] = 0x00000000; // kTPC
317 eventTrigData.fReadoutListV0.fList[8] = 0x00000000; // kTPC
318 eventTrigData.fReadoutListV0.fList[9] = 0x00000000; // kTPC
319 eventTrigData.fReadoutListV0.fList[10] = 0x00000000; // kTPC
320 eventTrigData.fReadoutListV0.fList[11] = 0x00000001; // kTRD
321 eventTrigData.fReadoutListV0.fList[12] = 0x00000001; // kTOF
322 eventTrigData.fReadoutListV0.fList[13] = 0x00000000; // kTOF
323 eventTrigData.fReadoutListV0.fList[14] = 0x00000000; // kTOF
324 eventTrigData.fReadoutListV0.fList[15] = 0x00000001; // kHMPID
325 eventTrigData.fReadoutListV0.fList[16] = 0x00000001; // kPHOS
326 eventTrigData.fReadoutListV0.fList[17] = 0x00000001; // kCPV
327 eventTrigData.fReadoutListV0.fList[18] = 0x00000001; // kPMD
328 eventTrigData.fReadoutListV0.fList[19] = 0x00000001; // kMUONTRK
329 eventTrigData.fReadoutListV0.fList[20] = 0x00000001; // kMUONTRG
330 eventTrigData.fReadoutListV0.fList[21] = 0x00000001; // kFMD
331 eventTrigData.fReadoutListV0.fList[22] = 0x00000001; // kT0
332 eventTrigData.fReadoutListV0.fList[23] = 0x00000001; // kV0
333 eventTrigData.fReadoutListV0.fList[24] = 0x00000001; // kZDC
334 eventTrigData.fReadoutListV0.fList[25] = 0x00000001; // kACORDE
335 eventTrigData.fReadoutListV0.fList[26] = 0x00000001; // kTRG
336 eventTrigData.fReadoutListV0.fList[27] = 0x00000001; // kEMCAL
337 eventTrigData.fReadoutListV0.fList[28] = 0x00000001; // kDAQTEST
338 eventTrigData.fReadoutListV0.fList[29] = 0x00000001; // kHLT
340 AliHLTComponentTriggerData trigData = {
341 sizeof(AliHLTComponentTriggerData),
342 sizeof(AliHLTEventTriggerDataV0),
346 AliHLTReadoutList readoutlistExpected;
347 for (Int_t i = 0; i < AliHLTDAQ::NumberOfDetectors(); ++i)
349 Int_t ddlid = AliHLTDAQ::DdlIDOffset(i);
350 readoutlistExpected.EnableDDLBit(ddlid);
353 const AliHLTUInt8_t (*attribs)[gkAliHLTBlockDAttributeCount];
354 AliHLTUInt64_t status = 0x0;
355 AliHLTCDHWrapper* const cdh = NULL;
356 AliHLTReadoutList readoutlist;
358 int result = AliHLTComponent::ExtractTriggerData(trigData, &attribs, &status, cdh, &readoutlist, true);
361 cerr << "ERROR: The method AliHLTComponent::ExtractTriggerData"
362 " fails for the old structure format." << endl;
365 if (attribs != &eventTrigData.fAttributes)
367 cerr << "ERROR: The method AliHLTComponent::ExtractTriggerData"
368 " fails to locate the attributes structure correctly." << endl;
373 cerr << "ERROR: The method AliHLTComponent::ExtractTriggerData"
374 " fails to locate the HLT status word correctly." << endl;
377 if ((const void*)cdh != (void*)&eventTrigData.fCommonHeader)
379 cerr << "ERROR: The method AliHLTComponent::ExtractTriggerData"
380 " fails to locate the Common Data Header (CDH) structure correctly." << endl;
383 if (memcmp(readoutlist.Buffer(), readoutlistExpected.Buffer(), readoutlist.BufferSize()) != 0)
385 cerr << "ERROR: The method AliHLTComponent::ExtractTriggerData"
386 " fails to extract the readout list correctly." << endl;
394 * Tests to see if reading old AliHLTReadoutList versions from root files works.
395 * \param filename The name of the file generated by GenerateReadoutListFile.C,
396 * which contains the readout list objects to test.
398 bool CheckReadingOldFormat(const char* filename = "$ALICE_ROOT/HLT/BASE/test/oldAliHLTReadoutListFormat.root")
400 TFile file(filename, "READ");
402 // Load the readout list objects.
403 AliHLTReadoutList* rl[5] = {NULL, NULL, NULL, NULL, NULL};
404 for (int i = 0; i < 5; ++i)
407 sprintf(name, "readoutlist%d", i+1);
408 rl[i] = dynamic_cast<AliHLTReadoutList*>(file.Get(name));
411 cerr << "ERROR: Could not get object '" << name
412 << "' from file '" << filename << "'." << endl;
417 // Now load the tree and see that the objects are the same as the readout
418 // list objects stored directly to the TFile.
419 const char* treename = "rltree";
420 TTree* tree = dynamic_cast<TTree*>(file.Get(treename));
423 cerr << "ERROR: Could not get TTree '" << treename
424 << "' from file '" << filename << "'." << endl;
427 AliHLTReadoutList* r = new AliHLTReadoutList;
428 tree->SetBranchAddress("readoutlist", &r);
429 for (int i = 0; i < 5; ++i)
432 if (r->BufferSize() != rl[i]->BufferSize())
434 cerr << "ERROR: readoutlist" << i+1
435 << " and the one from the TTree have different sizes."
439 #if 1 // ROOT_SVN_REVISION < 9999 //FIXME: after fixed https://savannah.cern.ch/bugs/?69241
440 r->SetDDLBit(9999999, kTRUE); // Triggers a reformating of the internal structure to the new version.
442 if (memcmp(r->Buffer(), rl[i]->Buffer(), r->BufferSize()) != 0)
444 cerr << "ERROR: readoutlist" << i+1
445 << " and the one from the TTree are different."
453 // Now check each readout list individually.
454 typedef AliHLTReadoutList RL;
455 Int_t alwaysoff = RL::kITSSPD
474 // We will need to try and set the missing EMCAL DDL bits.
475 // Otherwise the DetectorEnabled and DetectorDisabled methods will not
476 // give the expected answers.
477 for (int i = 24; i < AliHLTDAQ::NumberOfDdls("EMCAL"); ++i)
479 for (int j = 1; j <= 3; ++j)
481 Int_t ddlid = AliHLTDAQ::DdlIDOffset("EMCAL") | (i & 0xFF);
482 rl[j]->EnableDDLBit(ddlid);
487 if (not (rl[rlnum]->DetectorEnabled(RL::kTRG) and
488 rl[rlnum]->DetectorDisabled(alwaysoff | RL::kEMCAL | RL::kDAQTEST | RL::kHLT)
492 cerr << "ERROR: readoutlist" << rlnum+1 << " does not have the correct bits set." << endl;
497 if (not (rl[rlnum]->DetectorEnabled(RL::kTRG | RL::kEMCAL) and
498 rl[rlnum]->DetectorDisabled(alwaysoff | RL::kDAQTEST | RL::kHLT)
501 cerr << "ERROR: readoutlist" << rlnum+1 << " does not have the correct bits set." << endl;
506 if (not (rl[rlnum]->DetectorEnabled(RL::kTRG | RL::kEMCAL | RL::kDAQTEST) and
507 rl[rlnum]->DetectorDisabled(alwaysoff | RL::kHLT)
510 cerr << "ERROR: readoutlist" << rlnum+1 << " does not have the correct bits set." << endl;
515 if (not (rl[rlnum]->DetectorEnabled(RL::kTRG | RL::kEMCAL | RL::kDAQTEST | RL::kHLT) and
516 rl[rlnum]->DetectorDisabled(alwaysoff)
519 cerr << "ERROR: readoutlist" << rlnum+1 << " does not have the correct bits set." << endl;
520 printf("%i %i %i %i %i\n",rl[rlnum]->DetectorEnabled(RL::kTRG) ,rl[rlnum]->DetectorEnabled(RL::kEMCAL) ,rl[rlnum]->DetectorEnabled(RL::kDAQTEST) ,rl[rlnum]->DetectorEnabled(RL::kHLT), rl[rlnum]->DetectorEnabled(alwaysoff));
521 rl[rlnum]->Print("justlist");
525 if (not (rl[rlnum]->DetectorEnabled(RL::kTRG | RL::kDAQTEST | RL::kHLT) and
526 rl[rlnum]->DetectorDisabled(alwaysoff | RL::kEMCAL)
529 cerr << "ERROR: readoutlist" << rlnum+1 << " does not have the correct bits set." << endl;
530 rl[rlnum]->Print("justlist");
538 * Runs the unit tests for backward compatibility of AliHLTEventDDL succeeded.
539 * \returns true if the tests were passed and false otherwise.
541 bool testAliHLTEventDDLBackwardCompatibility()
543 if (not CheckReadoutListConvertedCorrectly()) return false;
544 if (not CheckHandlingOfOldAliHLTEventTriggerData()) return false;
545 if (not CheckReadingOldFormat()) return false;
551 int main(int /*argc*/, const char** /*argv*/)
553 bool resultOk = testAliHLTEventDDLBackwardCompatibility();
554 if (not resultOk) return 1;
558 #endif // __MAKECINT__