]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/rec/AliHLTEsdManagerImplementation.cxx
added abstract AliEsdManager to libHLTbase and implementation to libHLTrec
[u/mrichter/AliRoot.git] / HLT / rec / AliHLTEsdManagerImplementation.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   AliHLTEsdManagerImplementation.cxx
20     @author Matthias Richter
21     @date   
22     @brief  Manager for merging and writing of HLT ESDs
23 */
24
25 #include "AliHLTEsdManagerImplementation.h"
26 #include "AliHLTComponent.h"
27 #include "AliESDEvent.h"
28 #include "AliHLTMessage.h"
29 #include "AliESDEvent.h"
30 #include "TFile.h"
31 #include "TTree.h"
32 #include "TClass.h"
33 #include "TObject.h"
34 #include "TObjectTable.h"
35 #include "TSystem.h"
36 #include "TChain.h"
37 #include "TList.h"
38
39 /** ROOT macro for the implementation of ROOT specific class methods */
40 ClassImp(AliHLTEsdManagerImplementation)
41
42 AliHLTEsdManagerImplementation::AliHLTEsdManagerImplementation()
43   :
44   fESDs(),
45   fDirectory()
46 {
47   // see header file for class documentation
48   // or
49   // refer to README to build package
50   // or
51   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
52 }
53
54 AliHLTEsdManagerImplementation::~AliHLTEsdManagerImplementation()
55 {
56   // see header file for class documentation
57   for (unsigned int i=0; i<fESDs.size(); i++) {
58     if (fESDs[i]) {
59       delete fESDs[i];
60     }
61     fESDs[i]=NULL;
62   }
63 }
64
65 AliHLTEsdManagerImplementation::AliHLTEsdListEntry* AliHLTEsdManagerImplementation::Find(AliHLTComponentDataType dt) const
66 {
67   // see header file for class documentation
68   AliHLTEsdListEntry* pEntry=NULL;
69   for (unsigned int i=0; i<fESDs.size(); i++) {
70     if (fESDs[i] && *(fESDs[i])==dt) {
71       pEntry=const_cast<AliHLTEsdListEntry*>(fESDs[i]);
72     }
73   }
74   return pEntry;
75 }
76
77 int AliHLTEsdManagerImplementation::WriteESD(const AliHLTUInt8_t* pBuffer, AliHLTUInt32_t size,
78                                AliHLTComponentDataType dt, AliESDEvent* tgtesd, int eventno)
79 {
80   // see header file for class documentation
81   if (!pBuffer && size<=0) return -EINVAL;
82   int iResult=0;
83   AliHLTUInt32_t firstWord=*((AliHLTUInt32_t*)pBuffer);
84   if (firstWord==size-sizeof(AliHLTUInt32_t)) {
85     HLTDebug("create object from block %s size %d", AliHLTComponent::DataType2Text(dt).c_str(), size);
86     AliHLTMessage msg(const_cast<AliHLTUInt8_t*>(pBuffer), size);
87     TClass* objclass=msg.GetClass();
88     TObject* pObj=msg.ReadObject(objclass);
89     if (pObj && objclass) {
90       HLTDebug("object %p type %s created", pObj, objclass->GetName());
91       AliESDEvent* pESD=dynamic_cast<AliESDEvent*>(pObj);
92       TTree* pTree=NULL;
93       if (!pESD) {
94         pTree=dynamic_cast<TTree*>(pObj);
95         if (pTree) {
96           pESD=new AliESDEvent;
97           pESD->CreateStdContent();
98           if (pTree->GetEntries()>0) {
99             if (pTree->GetEntries()>1) {
100               HLTWarning("only one entry allowed for ESD embedded into tree, data block %s contains tree with %d entires, taking first entry",
101                          AliHLTComponent::DataType2Text(dt).c_str(), pTree->GetEntries());
102             }
103             pESD->ReadFromTree(pTree);
104             pTree->GetEvent(0);
105           }
106         } else {
107           HLTWarning("tree of data block %s has no events, skipping data block", AliHLTComponent::DataType2Text(dt).c_str());
108         }
109       }
110       if (pESD) {
111         AliHLTEsdListEntry* entry=Find(dt);
112         if (!entry) {
113           AliHLTEsdListEntry* newEntry=new AliHLTEsdListEntry(dt);
114           if (!fDirectory.IsNull()) {
115             newEntry->SetDirectory(fDirectory);
116           }
117           fESDs.push_back(newEntry);
118         }
119         if (tgtesd) {
120 #ifndef HAVE_NOT_ESD_COPY
121           *tgtesd=*pESD;
122 #else //HAVE_NOT_ESD_COPY
123           static bool warningPrinted=false;
124           if (!warningPrinted) {
125             HLTWarning("old version of AliESDEvent does not provide assignment operator, skip merging to global hltEsd");
126           }
127           warningPrinted=true;
128 #endif //HAVE_NOT_ESD_COPY
129         } else {
130         entry=Find(dt);
131         if (entry) {
132           entry->WriteESD(pESD, eventno);
133         } else {
134           HLTError("internal mismatch, can not create list entry");
135           iResult=-ENOMEM;
136         }
137         }
138       } else {
139         HLTWarning("data block %s is not of class type AliESDEvent, ignoring ...", AliHLTComponent::DataType2Text(dt).c_str());
140       }
141       if (pTree) {
142         // ESD has been created and must be cleaned up
143         pESD->Reset();
144         delete pESD;
145         pESD=NULL;
146       }
147       delete pObj;
148       pObj=NULL;
149     } else {
150     }
151   }
152   return iResult;
153 }
154
155 int AliHLTEsdManagerImplementation::PadESDs(int eventno)
156 {
157   // see header file for class documentation
158   int iResult=0;
159   for (unsigned int i=0; i<fESDs.size(); i++) {
160     if (fESDs[i]) {
161       int res=fESDs[i]->WriteESD(NULL, eventno);
162       if (res<0 && iResult>=0) iResult=res;
163     }
164   }
165   return iResult;
166 }
167
168 void AliHLTEsdManagerImplementation::SetDirectory(const char* directory)
169 {
170   // see header file for class documentation
171   if (!directory) return;
172   fDirectory=directory;
173   for (unsigned int i=0; i<fESDs.size(); i++) {
174     if (fESDs[i]) {
175       fESDs[i]->SetDirectory(directory);
176     }
177   }
178 }
179
180 TString AliHLTEsdManagerImplementation::GetFileNames(AliHLTComponentDataType dt) const
181 {
182   TString result;
183   for (unsigned int i=0; i<fESDs.size(); i++) {
184     if (fESDs[i] && *(fESDs[i])==dt) {
185       if (!result.IsNull()) result+=" ";
186       result+=fESDs[i]->GetFileName();
187     }
188   }
189   return result;
190 }
191
192 TTree* AliHLTEsdManagerImplementation::EmbedIntoTree(AliESDEvent* pESD, const char* name, const char* title)
193 {
194   // see header file for class documentation
195   int iResult=0;
196   TTree* pTree=new TTree(name, title);
197   if (pTree) {
198     pESD->WriteToTree(pTree);
199     pTree->Fill();
200     pTree->GetUserInfo()->Add(pESD);
201   } else {
202     iResult=-ENOMEM;
203   }
204
205   if (iResult<0) {
206     pTree->GetUserInfo()->Clear();
207     delete pTree;
208   }
209
210   return pTree;
211 }
212
213 AliHLTEsdManagerImplementation::AliHLTEsdListEntry::AliHLTEsdListEntry(AliHLTComponentDataType dt)
214   :
215   fName(),
216   fDirectory(),
217   fDt(dt),
218   fpFile(NULL),
219   fpTree(NULL),
220   fpEsd(NULL)
221 {
222   // see header file for class documentation
223 }
224
225 AliHLTEsdManagerImplementation::AliHLTEsdListEntry::~AliHLTEsdListEntry()
226 {
227   // see header file for class documentation
228   if (fpEsd) delete fpEsd;
229   fpEsd=NULL;
230
231   if (fpTree) delete fpTree;
232   fpTree=NULL;
233
234   if (fpFile) delete fpFile;
235   fpFile=NULL;
236 }
237
238 bool AliHLTEsdManagerImplementation::AliHLTEsdListEntry::operator==(AliHLTComponentDataType dt) const
239 {
240   // see header file for class documentation
241   return fDt==dt;
242 }
243
244 int AliHLTEsdManagerImplementation::AliHLTEsdListEntry::WriteESD(AliESDEvent* pSrcESD, int eventno)
245 {
246   // see header file for class documentation
247   int iResult=0;
248 #ifndef HAVE_NOT_ESD_COPY
249   if (fName.IsNull()) {
250     // this is the first event, create the file name
251     TString origin;
252     origin.Insert(0, fDt.fOrigin, kAliHLTComponentDataTypefOriginSize);
253     origin.Remove(TString::kTrailing, ' ');
254     origin.ToUpper();
255     fName="";
256     if (!fDirectory.IsNull()) {
257       fName+=fDirectory; fName+="/";
258     }
259     fName+="AliHLT"; fName+=origin;
260     if (fDt!=kAliHLTDataTypeESDObject &&
261         fDt!=kAliHLTDataTypeESDTree) {
262
263       HLTWarning("non-standard ESD type %s", AliHLTComponent::DataType2Text(fDt).c_str());
264       TString id;
265       id.Insert(0, fDt.fID, kAliHLTComponentDataTypefIDsize);
266       id.Remove(TString::kTrailing, ' ');
267       id.ToUpper();
268       fName+="_"; fName+=id; fName+=".root";
269     } else {
270       fName+="ESDs.root";
271     }
272
273     fpFile=new TFile(fName, "RECREATE");
274     fpTree=new TTree("esdTree", "Tree with HLT ESD objects");
275     fpEsd=new AliESDEvent;
276     if (fpEsd) {
277       fpEsd->CreateStdContent();
278       if (fpTree) {
279         fpEsd->WriteToTree(fpTree);
280       }
281     }
282   }
283
284   if (fpFile && fpTree && fpEsd) {
285     // synchronize and add empty events
286     fpEsd->Reset();
287     int nofCurrentEvents=fpTree->GetEntries();
288     if (nofCurrentEvents<eventno) {
289       iResult=1; // indicate tree to be written
290       HLTDebug("adding %d empty events to file %s", eventno-nofCurrentEvents, fName.Data());
291       for (int i=nofCurrentEvents; i<eventno; i++) {
292         fpTree->Fill();
293       }
294     }
295
296     if (iResult>=0 && pSrcESD) {
297       *fpEsd=*pSrcESD;
298       fpTree->Fill();
299       iResult=1; // indicate tree to be written
300     }
301
302     if (iResult>0) {
303       fpFile->cd();
304       fpTree->GetUserInfo()->Add(fpEsd);
305       fpTree->Write(fpTree->GetName(),TObject::kOverwrite);
306       fpTree->GetUserInfo()->Clear();
307     }
308   }
309 #else //HAVE_NOT_ESD_COPY
310   // this is the old workaround, necessary for older AliRoot versions
311   // version<=v4-12-Release
312
313   // we need to copy the ESD, I did not find an approptiate
314   // method, the workaround is to save the ESD in a temporary
315   // tree, read the content back into the ESD structure
316   // used for filling.
317   // Unfortunately the following code crashes at the second event.
318   // The expert on the ESD (Christian Klein Boesig) does not have
319   // a solution either. It seems to be a problem in ROOT.
320   //  TTree* dummy=new TTree("dummy","dummy");
321   //  dummy->SetDirectory(0);
322   //  pESD->WriteToTree(dummy);
323   //  dummy->Fill();
324   //  dummy->GetUserInfo()->Add(pESD);
325   //  fpEsd->ReadFromTree(dummy);
326   //  dummy->GetEvent(0);
327   //  fpEsd->WriteToTree(fpTree);
328   //  fpTree->Fill();
329   //  dummy->GetUserInfo()->Clear();
330   //  delete dummy;
331   //
332   // The only way is via TChain, which is working on files only at the
333   // time of writing.
334   // We use temporary files for the new event to be copied into the
335   // existing tree.
336   //
337   if (fName.IsNull()) {
338     // this is the first event, create the file on disk and write ESD
339     TString origin;
340     origin.Insert(0, fDt.fOrigin, kAliHLTComponentDataTypefOriginSize);
341     origin.Remove(TString::kTrailing, ' ');
342     origin.ToUpper();
343     fName="";
344     if (!fDirectory.IsNull()) {
345       fName+=fDirectory; fName+="/";
346     }
347     fName+="AliHLT"; fName+=origin;
348     if (fDt!=kAliHLTDataTypeESDObject &&
349         fDt!=kAliHLTDataTypeESDTree) {
350
351       HLTWarning("non-standard ESD type %s", AliHLTComponent::DataType2Text(fDt).c_str());
352       TString id;
353       id.Insert(0, fDt.fID, kAliHLTComponentDataTypefIDsize);
354       id.Remove(TString::kTrailing, ' ');
355       id.ToUpper();
356       fName+="_"; fName+=id; fName+=".root";
357     } else {
358       fName+="ESDs.root";
359     }
360
361     if (!gSystem->AccessPathName(fName)) {
362       // file exists, delete
363       TString shellcmd="rm -f ";
364       shellcmd+=fName;
365       gSystem->Exec(shellcmd);
366     }
367   }
368
369   TChain chain("esdTree");
370   TList cleanup;
371   cleanup.SetOwner();
372
373   int nofCurrentEvents=0;
374   if (iResult>=0) {
375     if (!gSystem->AccessPathName(fName)) {
376       // these are the other events, use the target file and temporary files to merge
377       // with TChain
378       chain.Add(fName);
379
380       if (eventno>=0) {
381         TFile file(fName);
382         if (!file.IsZombie()) {
383           TTree* pSrcTree;
384           file.GetObject("esdTree", pSrcTree);
385           if (pSrcTree) {
386             nofCurrentEvents=pSrcTree->GetEntries();
387           }
388           file.Close();
389         }
390       }
391     }
392   }
393
394   // synchronize and add empty events
395   if (nofCurrentEvents<eventno) {
396     iResult=1; // indicate files to merge
397     TTree* pTgtTree=new TTree("esdTree", "Tree with HLT ESD objects");
398     if (pTgtTree) {
399       pTgtTree->SetDirectory(0);
400       AliESDEvent* pTmpESD=new AliESDEvent;
401       if (pTmpESD) {
402         TString tmpfilename;
403         FILE* pTmpFile=gSystem->TempFileName(tmpfilename);
404         if (pTmpFile) {
405           fclose(pTmpFile);
406           pTmpFile=NULL;
407           cleanup.Add(new TObjString(tmpfilename));
408           TFile emptyevents(tmpfilename, "RECREATE");
409           if (!emptyevents.IsZombie()) {
410             pTmpESD->CreateStdContent();
411             pTmpESD->WriteToTree(pTgtTree);
412             HLTDebug("adding %d empty events to file %s", eventno-nofCurrentEvents, fName.Data());
413             for (int i=nofCurrentEvents; i<eventno; i++) {
414               pTgtTree->Fill();
415             }
416             pTgtTree->GetUserInfo()->Add(pTmpESD);
417             emptyevents.cd();
418             pTgtTree->Write();
419             emptyevents.Close();
420             chain.Add(tmpfilename);
421             pTgtTree->GetUserInfo()->Clear();
422           }
423         }
424         delete pTmpESD;
425       } else {
426         iResult=-ENOMEM;
427       }
428       delete pTgtTree;
429     } else {
430       iResult=-ENOMEM;
431     }
432   }
433
434   if (iResult>=0 && pSrcESD) {
435     // add the new event to the chain
436     iResult=1; // indicate files to merge
437     TString tmpfilename=WriteTempFile(pSrcESD);
438     if (!tmpfilename.IsNull()) {
439       chain.Add(tmpfilename);
440       cleanup.Add(new TObjString(tmpfilename));
441     }
442   }
443
444   if (iResult>0) {
445     // build temporary file name for chain output
446     TString tgtName;
447     FILE* pTmpFile=gSystem->TempFileName(tgtName);
448     if (pTmpFile) {
449       fclose(pTmpFile);
450       pTmpFile=NULL;
451
452       // there have been problems with the memory consumption when using
453       // TChain::Merge
454       // but using a separate loop soemtimes crashes in AliESDEvent::ReadFromTree
455       // since this is for backward compatiblity only, we take the TChain::Merge
456       chain.Merge(tgtName);
457 //       TFile tgtFile(tgtName, "RECREATE");
458 //       TTree* pTgtTree=new TTree("esdTree", "Tree with HLT ESD objects");
459 //       AliESDEvent* pTgtESD=new AliESDEvent;
460 //       if (pTgtTree && pTgtESD) {
461 //      pTgtESD->ReadFromTree(&chain);
462 //      pTgtESD->WriteToTree(pTgtTree);
463
464 //      int nofEvents=chain.GetEntries();
465 //      for (int event=0; event<nofEvents; event++) {
466 //        chain.GetEntry(event);
467 //        pTgtTree->Fill();
468 //      }
469
470 //      pTgtTree->GetUserInfo()->Add(pTgtESD);
471 //      tgtFile.cd();
472 //      pTgtTree->Write();
473 //      pTgtTree->GetUserInfo()->Clear();
474 //       } else {
475 //      iResult=-ENOMEM;
476 //       }
477
478 //       if (pTgtTree) delete pTgtTree;
479 //       if (pTgtESD) delete pTgtESD;
480 //       tgtFile.Close();
481
482       // rename the merged file to the original file
483       TString shellcmd="mv ";
484       shellcmd+=tgtName + " " + fName;
485       if (gSystem->Exec(shellcmd)==0) {
486         HLTDebug("renaming %s to %s", tgtName.Data(), fName.Data());
487       } else {
488         HLTError("can not rename temporary file %s to %s", tgtName.Data(), fName.Data());
489       }
490     } else {
491       HLTError("can not get temporary file name from system");
492       iResult=-EBADF;
493     }
494   }
495
496   // delete temporary files
497   // the list objects are cleaned up by the TList destructor as the
498   // list is owner
499   TIter entry(&cleanup);
500   while (TObject* pObj=entry.Next()) {
501     if (dynamic_cast<TObjString*>(pObj)) {
502       TString shellcmd="rm -f ";
503       shellcmd+=(dynamic_cast<TObjString*>(pObj))->GetString();
504       gSystem->Exec(shellcmd);
505     }
506   }
507 #endif //HAVE_NOT_ESD_COPY
508
509   return iResult;
510 }
511
512 TString AliHLTEsdManagerImplementation::AliHLTEsdListEntry::WriteTempFile(AliESDEvent* pESD) const
513 {
514   // see header file for class documentation
515   int iResult=0;
516   TString tmpfilename;
517   FILE* pTmpFile=gSystem->TempFileName(tmpfilename);
518   if (pTmpFile) {
519     fclose(pTmpFile);
520     pTmpFile=NULL;
521
522     TFile file(tmpfilename, "RECREATE");
523     if (!file.IsZombie()) {
524       TTree* pTree=AliHLTEsdManagerImplementation::EmbedIntoTree(pESD);
525       if (pTree) {
526         file.cd();
527         if (pTree->Write()>0) {
528         } else {
529           HLTError("can not write esd tree to temporary file %s", tmpfilename.Data());
530         }
531
532         pTree->GetUserInfo()->Clear();
533         delete pTree;
534       } else {
535         iResult=-ENOMEM;
536       }
537       file.Close();
538     } else {
539       HLTError("can not open file %s", tmpfilename.Data());
540       iResult=-EBADF;
541     }
542   } else {
543     HLTError("can not get temporary file name from system");
544     iResult=-EBADF;
545   }
546
547   if (iResult<0) {
548     if (gSystem->AccessPathName(tmpfilename)==0) {
549       TString shellcmd="rm -f ";
550       shellcmd+=tmpfilename;
551       gSystem->Exec(shellcmd);
552     }
553     tmpfilename="";
554   }
555   return tmpfilename;
556 }
557
558 void AliHLTEsdManagerImplementation::AliHLTEsdListEntry::SetDirectory(const char* directory)
559 {
560   // see header file for class documentation
561   if (!directory) return;
562   if (!fName.IsNull()) {
563     HLTWarning("ESD entry already in writing mode (%s), ignoring directory", fName.Data());
564     return;
565   }
566   fDirectory=directory;
567 }
568
569 void AliHLTEsdManagerImplementation::AliHLTEsdListEntry::Delete()
570 {
571   // see header file for class documentation
572   if (fName.IsNull()) return;
573   if (gSystem->AccessPathName(fName)!=0) return;
574
575   TString shellcmd="rm -f ";
576   shellcmd+=fName;
577   gSystem->Exec(shellcmd);
578   fName="";
579 }
580
581 const char* AliHLTEsdManagerImplementation::AliHLTEsdListEntry::GetFileName() const
582 {
583   // see header file for class documentation
584   return fName.Data();
585 }