]>
Commit | Line | Data |
---|---|---|
652cf9d2 | 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 testAliHLTEsdManager.C | |
20 | @author Matthias Richter | |
21 | @date | |
22 | @brief Test macro/program for the AliHLTEsdManager class | |
23 | */ | |
24 | ||
25 | #ifndef __CINT__ | |
26 | #include "TFile.h" | |
27 | #include "TDatime.h" | |
28 | #include "TRandom.h" | |
29 | #include "TArrayI.h" | |
30 | #include "TSystem.h" | |
31 | #include "TTree.h" | |
32 | #include "TFile.h" | |
33 | #include "AliESDEvent.h" | |
34 | #include "AliHLTDataTypes.h" | |
35 | #include "AliHLTComponent.h" | |
36 | #include "AliHLTEsdManager.h" | |
37 | #include "AliHLTMessage.h" | |
38 | #include "AliHLTSystem.h" | |
39 | #include <ostream> | |
40 | #endif //__CINT__ | |
41 | ||
42 | ///////////////////////////////////////////////////////////////// | |
43 | ///////////////////////////////////////////////////////////////// | |
44 | ///////////////////////////////////////////////////////////////// | |
45 | // | |
46 | // configuration of the test program | |
47 | // | |
48 | ||
49 | // printouts or not | |
50 | const bool bVerbose=false; | |
51 | ||
52 | ||
53 | ///////////////////////////////////////////////////////////////// | |
54 | ///////////////////////////////////////////////////////////////// | |
55 | ///////////////////////////////////////////////////////////////// | |
56 | // | |
57 | // forward declarations | |
58 | // | |
59 | class AliHLTEsdManager; | |
60 | int GetRandom(int min, int max); | |
521766bd | 61 | int CreateAndWriteESD(AliHLTEsdManager* manager, int eventno, AliHLTComponentDataType dt, AliESDEvent* pTgt); |
652cf9d2 | 62 | int CheckFields(const char* file, TArrayI* fields); |
63 | int CheckFields(TTree* pTree, TArrayI* fields, const char* file); | |
64 | ||
65 | ///////////////////////////////////////////////////////////////// | |
66 | ///////////////////////////////////////////////////////////////// | |
67 | ///////////////////////////////////////////////////////////////// | |
68 | // | |
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() | |
74 | { | |
75 | cout << "macro still not working with CINT, sorry!" << endl; | |
76 | return 0; | |
77 | ||
78 | int iResult=0; | |
79 | #ifdef __CINT__ | |
80 | if (gSystem->Load("libHLTrec.so")<0) { | |
521766bd | 81 | cerr << "error: error loading libHLTrec.so library" << endl; |
652cf9d2 | 82 | return -1; |
83 | } | |
84 | #endif | |
85 | ||
521766bd | 86 | AliHLTEsdManager* pManager=AliHLTEsdManager::New(); |
87 | if (!pManager) { | |
88 | cerr << "error: can not create manager instance" << endl; | |
89 | return -1; | |
90 | } | |
91 | pManager->SetDirectory(gSystem->TempDirectory()); | |
652cf9d2 | 92 | |
93 | int nofEvents=10; | |
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; | |
521766bd | 99 | CreateAndWriteESD(pManager, event, tpcesd, NULL); |
652cf9d2 | 100 | } |
101 | ||
521766bd | 102 | AliHLTEsdManager::Delete(pManager); |
652cf9d2 | 103 | return iResult; |
104 | } | |
105 | ||
106 | ///////////////////////////////////////////////////////////////// | |
107 | ///////////////////////////////////////////////////////////////// | |
108 | ///////////////////////////////////////////////////////////////// | |
109 | // | |
110 | int main(int /*argc*/, const char** /*argv*/) | |
111 | { | |
112 | int iResult=0; | |
113 | int nofEvents=10; | |
114 | AliHLTSystem gHLT; | |
115 | gHLT.SetGlobalLoggingLevel(kHLTLogDefault); | |
116 | ||
521766bd | 117 | AliHLTEsdManager* pManager=AliHLTEsdManager::New(); |
118 | if (!pManager) { | |
119 | cerr << "error: can not create manager instance" << endl; | |
120 | return -1; | |
121 | } | |
122 | pManager->SetDirectory(gSystem->TempDirectory()); | |
652cf9d2 | 123 | |
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, | |
131 | kAliHLTVoidDataType | |
132 | }; | |
133 | ||
134 | TTree* pMasterTree=new TTree("esdTree", "Tree with HLT ESD objects"); | |
521766bd | 135 | pMasterTree->SetDirectory(0); |
652cf9d2 | 136 | AliESDEvent* pMasterESD=new AliESDEvent; |
137 | pMasterESD->CreateStdContent(); | |
138 | pMasterESD->WriteToTree(pMasterTree); | |
139 | ||
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)); | |
146 | } | |
147 | AliESDEvent* pTgt=NULL; | |
148 | //if (type==0) pTgt=pMasterESD; | |
521766bd | 149 | int field=CreateAndWriteESD(pManager, event, types[type], pTgt); |
652cf9d2 | 150 | if (field>=0) { |
151 | (*randomFields[type])[event]=field; | |
152 | } else { | |
153 | iResult=-1; | |
154 | break; | |
155 | } | |
156 | } | |
157 | pMasterTree->Fill(); | |
158 | } | |
159 | ||
160 | if (iResult>=0) { | |
521766bd | 161 | pManager->PadESDs(nofEvents); |
652cf9d2 | 162 | } |
163 | ||
ba93bdb0 | 164 | vector<TString> filenames; |
652cf9d2 | 165 | for (int type=0; types[type]!=kAliHLTVoidDataType; type++) { |
521766bd | 166 | TString filename=pManager->GetFileNames(types[type]); |
ba93bdb0 | 167 | filenames.push_back(filename); |
168 | } | |
169 | // delete the manager instance to make sure the files are closed | |
170 | AliHLTEsdManager::Delete(pManager); | |
171 | pManager=NULL; | |
172 | for (int type=0; types[type]!=kAliHLTVoidDataType; type++) { | |
652cf9d2 | 173 | if (iResult>=0) { |
ba93bdb0 | 174 | iResult=CheckFields(filenames[type], dynamic_cast<TArrayI*>(randomFields[type])); |
652cf9d2 | 175 | } |
176 | TString shellcmd="rm -f "; | |
ba93bdb0 | 177 | shellcmd+=filenames[type]; |
652cf9d2 | 178 | gSystem->Exec(shellcmd); |
179 | } | |
180 | ||
181 | vector<TArrayI*>::iterator element; | |
182 | while ((element=randomFields.end())!=randomFields.begin()) { | |
183 | element--; | |
184 | if (*element) delete *element; | |
185 | randomFields.pop_back(); | |
186 | } | |
187 | ||
188 | delete pMasterESD; | |
189 | ||
190 | return iResult; | |
191 | } | |
192 | ||
193 | Bool_t seedSet=kFALSE; | |
194 | ||
195 | /** | |
196 | * Get a random number in the given range. | |
197 | */ | |
198 | int GetRandom(int min, int max) | |
199 | { | |
200 | if (max-min<2) return min; | |
201 | static TRandom rand; | |
202 | if (!seedSet) { | |
203 | TDatime dt; | |
204 | rand.SetSeed(dt.Get()); | |
205 | seedSet=kTRUE; | |
206 | } | |
207 | return rand.Integer(max-min); | |
208 | } | |
209 | ||
210 | /** | |
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 | |
215 | */ | |
521766bd | 216 | int CreateAndWriteESD(AliHLTEsdManager* pManager, int eventno, AliHLTComponentDataType dt, AliESDEvent* pTgt) |
652cf9d2 | 217 | { |
218 | int iResult=0; | |
219 | int magField=0; | |
521766bd | 220 | if (!pManager) { |
221 | cerr << "error: missing manager instance" << endl; | |
222 | return -1; | |
223 | } | |
652cf9d2 | 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); | |
233 | pTree->Fill(); | |
234 | pTree->GetUserInfo()->Add(pESD); | |
235 | AliHLTMessage msg(kMESS_OBJECT); | |
236 | msg.WriteObject(pTree); | |
237 | Int_t iMsgLength=msg.Length(); | |
238 | if (iMsgLength>0) { | |
239 | msg.SetLength(); // sets the length to the first (reserved) word | |
521766bd | 240 | iResult=pManager->WriteESD((AliHLTUInt8_t*)msg.Buffer(), iMsgLength, dt, pTgt, eventno); |
652cf9d2 | 241 | } |
242 | pTree->GetUserInfo()->Clear(); | |
243 | delete pTree; | |
244 | delete pESD; | |
245 | } else { | |
246 | message=": omitting block "; | |
247 | } | |
248 | if (iResult>=0) iResult=magField; | |
249 | if (bVerbose) cout << "event " << eventno << message << AliHLTComponent::DataType2Text(dt).c_str() << ": " << iResult << endl; | |
250 | return iResult; | |
251 | } | |
252 | ||
253 | /** | |
254 | * Read the ESD from the file and compare with the | |
255 | * random field values previously set to the ESDs | |
256 | */ | |
257 | int CheckFields(const char* file, TArrayI* fields) | |
258 | { | |
259 | if (!file || !fields) { | |
521766bd | 260 | cerr << "error: invalid parameters" << endl; |
652cf9d2 | 261 | return 0; |
262 | } | |
263 | TFile esdfile(file); | |
264 | if (!esdfile.IsZombie()) { | |
265 | TTree* pTree=NULL; | |
266 | esdfile.GetObject("esdTree", pTree); | |
267 | if (pTree) { | |
268 | int res=CheckFields(pTree, fields, file); | |
269 | if (res<0) return res; | |
270 | } else { | |
521766bd | 271 | cerr << "error: can not find esdTree in file " << file << endl; |
652cf9d2 | 272 | } |
273 | } else { | |
521766bd | 274 | cerr << "error: can not open file " << file << endl; |
652cf9d2 | 275 | return -1; |
276 | } | |
277 | cout << "checking: " << file << " ok" << endl; | |
278 | return 0; | |
279 | } | |
280 | ||
281 | /** | |
282 | * Compare ESD from tree with the | |
283 | * random field values previously set to the ESDs | |
284 | */ | |
285 | int CheckFields(TTree* pTree, TArrayI* fields, const char* file) | |
286 | { | |
287 | if (fields->GetSize()!=pTree->GetEntries()) { | |
521766bd | 288 | cerr << "error: event number mismatch in file " << file << " : expected " << fields->GetSize() << " found " << pTree->GetEntries() << endl; |
652cf9d2 | 289 | return -1; |
290 | } | |
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()) { | |
521766bd | 296 | cerr << "error: magnetic field mismatch in file " << file << " event " << event << ": expected " << fields->At(event) << " found " << pESD->GetMagneticField() << endl; |
652cf9d2 | 297 | return -1; |
298 | } | |
299 | } | |
300 | ||
301 | delete pESD; | |
302 | return 0; | |
303 | } |