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