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: Matthias Richter <Matthias.Richter@ift.uib.no> *
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 testAliHLTEsdManager.C
20 @author Matthias Richter
22 @brief Test macro/program for the AliHLTEsdManager class
33 #include "AliESDEvent.h"
34 #include "AliHLTDataTypes.h"
35 #include "AliHLTComponent.h"
36 #include "AliHLTEsdManager.h"
37 #include "AliHLTMessage.h"
38 #include "AliHLTSystem.h"
42 /////////////////////////////////////////////////////////////////
43 /////////////////////////////////////////////////////////////////
44 /////////////////////////////////////////////////////////////////
46 // configuration of the test program
50 const bool bVerbose=false;
53 /////////////////////////////////////////////////////////////////
54 /////////////////////////////////////////////////////////////////
55 /////////////////////////////////////////////////////////////////
57 // forward declarations
59 class AliHLTEsdManager;
60 int GetRandom(int min, int max);
61 int CreateAndWriteESD(AliHLTEsdManager* manager, int eventno, AliHLTComponentDataType dt, AliESDEvent* pTgt);
62 int CheckFields(const char* file, TArrayI* fields);
63 int CheckFields(TTree* pTree, TArrayI* fields, const char* file);
65 /////////////////////////////////////////////////////////////////
66 /////////////////////////////////////////////////////////////////
67 /////////////////////////////////////////////////////////////////
69 // compiled version and root macro are not equivalent for this unit test
70 // The macro just tests with one data type
71 // This is mainly due to some restrictions in CINT which can not handle
72 // the array of data types correctly
73 int testAliHLTEsdManager()
75 cout << "macro still not working with CINT, sorry!" << endl;
80 if (gSystem->Load("libHLTrec.so")<0) {
81 cerr << "error: error loading libHLTrec.so library" << endl;
86 AliHLTEsdManager* pManager=AliHLTEsdManager::New();
88 cerr << "error: can not create manager instance" << endl;
91 pManager->SetDirectory(gSystem->TempDirectory());
94 AliHLTComponentDataType tpcesd;
95 AliHLTComponent::SetDataType(tpcesd, "ESD_TREE", "TPC ");
96 cout << AliHLTComponent::DataType2Text(tpcesd).c_str() << endl;
97 for (int event=0; event<nofEvents && iResult>=0; event++) {
98 cout << AliHLTComponent::DataType2Text(tpcesd).c_str() << endl;
99 CreateAndWriteESD(pManager, event, tpcesd, NULL);
102 AliHLTEsdManager::Delete(pManager);
106 /////////////////////////////////////////////////////////////////
107 /////////////////////////////////////////////////////////////////
108 /////////////////////////////////////////////////////////////////
110 int main(int /*argc*/, const char** /*argv*/)
115 gHLT.SetGlobalLoggingLevel(kHLTLogDefault);
117 AliHLTEsdManager* pManager=AliHLTEsdManager::New();
119 cerr << "error: can not create manager instance" << endl;
122 pManager->SetDirectory(gSystem->TempDirectory());
124 AliHLTComponentDataType types[] = {
125 // first entry is special, ESD is written to the global target ESD
126 //kAliHLTDataTypeESDTree|kAliHLTDataOriginTPC,
127 kAliHLTDataTypeESDTree|kAliHLTDataOriginTPC,
128 kAliHLTDataTypeESDTree|kAliHLTDataOriginPHOS,
129 kAliHLTDataTypeESDTree|kAliHLTDataOriginTRD,
130 kAliHLTDataTypeDDLRaw|kAliHLTDataOriginTRD,
134 TTree* pMasterTree=new TTree("esdTree", "Tree with HLT ESD objects");
135 pMasterTree->SetDirectory(0);
136 AliESDEvent* pMasterESD=new AliESDEvent;
137 pMasterESD->CreateStdContent();
138 pMasterESD->WriteToTree(pMasterTree);
140 vector<TArrayI*> randomFields;
141 for (int event=0; event<nofEvents && iResult>=0; event++) {
142 pMasterESD->ResetStdContent();
143 for (unsigned int type=0; types[type]!=kAliHLTVoidDataType && iResult>=0; type++) {
144 if (randomFields.size()<=type) {
145 randomFields.push_back(new TArrayI(nofEvents));
147 AliESDEvent* pTgt=NULL;
148 //if (type==0) pTgt=pMasterESD;
149 int field=CreateAndWriteESD(pManager, event, types[type], pTgt);
151 (*randomFields[type])[event]=field;
161 pManager->PadESDs(nofEvents);
164 vector<TString> filenames;
165 for (int type=0; types[type]!=kAliHLTVoidDataType; type++) {
166 TString filename=pManager->GetFileNames(types[type]);
167 filenames.push_back(filename);
169 // delete the manager instance to make sure the files are closed
170 AliHLTEsdManager::Delete(pManager);
172 for (int type=0; types[type]!=kAliHLTVoidDataType; type++) {
174 iResult=CheckFields(filenames[type], dynamic_cast<TArrayI*>(randomFields[type]));
176 TString shellcmd="rm -f ";
177 shellcmd+=filenames[type];
178 gSystem->Exec(shellcmd);
181 vector<TArrayI*>::iterator element;
182 while ((element=randomFields.end())!=randomFields.begin()) {
184 if (*element) delete *element;
185 randomFields.pop_back();
193 Bool_t seedSet=kFALSE;
196 * Get a random number in the given range.
198 int GetRandom(int min, int max)
200 if (max-min<2) return min;
204 rand.SetSeed(dt.Get());
207 return rand.Integer(max-min);
211 * Creates a dummy ESD object and sets the magnetic field to a random number.
212 * The ESD is streamed via AliHLTMessage and processed by the AliHLTEsdmanager.
213 * @return 0 no ESD created, >0 random number of the magnetic field,
214 * neg error if failed
216 int CreateAndWriteESD(AliHLTEsdManager* pManager, int eventno, AliHLTComponentDataType dt, AliESDEvent* pTgt)
221 cerr << "error: missing manager instance" << endl;
224 const char* message="";
225 if ((GetRandom(0,10)%3)==0) {
226 message=": adding ESD for block ";
227 TTree* pTree=new TTree;
228 AliESDEvent* pESD=new AliESDEvent;
229 pESD->CreateStdContent();
230 magField=GetRandom(1, 1000);
231 pESD->SetMagneticField(magField);
232 pESD->WriteToTree(pTree);
234 pTree->GetUserInfo()->Add(pESD);
235 AliHLTMessage msg(kMESS_OBJECT);
236 msg.WriteObject(pTree);
237 Int_t iMsgLength=msg.Length();
239 msg.SetLength(); // sets the length to the first (reserved) word
240 iResult=pManager->WriteESD((AliHLTUInt8_t*)msg.Buffer(), iMsgLength, dt, pTgt, eventno);
242 pTree->GetUserInfo()->Clear();
246 message=": omitting block ";
248 if (iResult>=0) iResult=magField;
249 if (bVerbose) cout << "event " << eventno << message << AliHLTComponent::DataType2Text(dt).c_str() << ": " << iResult << endl;
254 * Read the ESD from the file and compare with the
255 * random field values previously set to the ESDs
257 int CheckFields(const char* file, TArrayI* fields)
259 if (!file || !fields) {
260 cerr << "error: invalid parameters" << endl;
264 if (!esdfile.IsZombie()) {
266 esdfile.GetObject("esdTree", pTree);
268 int res=CheckFields(pTree, fields, file);
269 if (res<0) return res;
271 cerr << "error: can not find esdTree in file " << file << endl;
274 cerr << "error: can not open file " << file << endl;
277 cout << "checking: " << file << " ok" << endl;
282 * Compare ESD from tree with the
283 * random field values previously set to the ESDs
285 int CheckFields(TTree* pTree, TArrayI* fields, const char* file)
287 if (fields->GetSize()!=pTree->GetEntries()) {
288 cerr << "error: event number mismatch in file " << file << " : expected " << fields->GetSize() << " found " << pTree->GetEntries() << endl;
291 AliESDEvent* pESD=new AliESDEvent;
292 pESD->ReadFromTree(pTree);
293 for (int event=0; event<pTree->GetEntries(); event++) {
294 pTree->GetEvent(event);
295 if (fields->At(event)!=pESD->GetMagneticField()) {
296 cerr << "error: magnetic field mismatch in file " << file << " event " << event << ": expected " << fields->At(event) << " found " << pESD->GetMagneticField() << endl;