- workaround for copying and merging of ESDs: The HLTOUT contains ESDs for every...
[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           TTree* pTmpTree=AliHLTEsdManager::EmbedIntoTree(pESD);
121           if (pTmpTree) {
122             tgtesd->ReadFromTree(pTmpTree);
123             pTmpTree->GetEvent(0);
124             pTmpTree->GetUserInfo()->Clear();
125             delete pTmpTree;
126             HLTDebug("data block %s written to target ESD", AliHLTComponent::DataType2Text(dt).c_str());
127           } else {
128             iResult=-ENOMEM;
129           }
130         } else {
131         entry=Find(dt);
132         if (entry) {
133           entry->WriteESD(pESD, eventno);
134         } else {
135           HLTError("internal mismatch, can not create list entry");
136           iResult=-ENOMEM;
137         }
138         }
139       } else {
140         HLTWarning("data block %s is not of class type AliESDEvent, ignoring ...", AliHLTComponent::DataType2Text(dt).c_str());
141       }
142       if (pTree) {
143         // ESD has been created and must be cleaned up
144         delete pESD;
145         pESD=NULL;
146       }
147       delete pObj;
148       pObj=NULL;
149     } else {
150     }
151   }
152   return iResult;
153 }
154
155 int AliHLTEsdManager::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 AliHLTEsdManager::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 AliHLTEsdManager::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* AliHLTEsdManager::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 AliHLTEsdManager::AliHLTEsdListEntry::AliHLTEsdListEntry(AliHLTComponentDataType dt)
214   :
215   fName(),
216   fDirectory(),
217   fDt(dt)
218 {
219   // see header file for class documentation
220 }
221
222 AliHLTEsdManager::AliHLTEsdListEntry::~AliHLTEsdListEntry()
223 {
224   // see header file for class documentation
225 }
226
227 bool AliHLTEsdManager::AliHLTEsdListEntry::operator==(AliHLTComponentDataType dt) const
228 {
229   // see header file for class documentation
230   return fDt==dt;
231 }
232
233 int AliHLTEsdManager::AliHLTEsdListEntry::WriteESD(AliESDEvent* pSrcESD, int eventno)
234 {
235   // we need to copy the ESD, I did not find an approptiate
236   // method, the workaround is to save the ESD in a temporary
237   // tree, read the content back into the ESD structure
238   // used for filling.
239   // Unfortunately the following code crashes at the second event.
240   // The expert on the ESD (Christian Klein Boesig) does not have
241   // a solution either. It seems to be a problem in ROOT.
242   //  TTree* dummy=new TTree("dummy","dummy");
243   //  dummy->SetDirectory(0);
244   //  pESD->WriteToTree(dummy);
245   //  dummy->Fill();
246   //  dummy->GetUserInfo()->Add(pESD);
247   //  fpEsd->ReadFromTree(dummy);
248   //  dummy->GetEvent(0);
249   //  fpEsd->WriteToTree(fpTree);
250   //  fpTree->Fill();
251   //  dummy->GetUserInfo()->Clear();
252   //  delete dummy;
253   //
254   // The only way is via TChain, which is working on files only at the
255   // time of writing.
256   // We use temporary files for the new event to be copied into the
257   // existing tree.
258   //
259   int iResult=0;
260   if (fName.IsNull()) {
261     // this is the first event, create the file on disk and write ESD
262     TString origin;
263     origin.Insert(0, fDt.fOrigin, kAliHLTComponentDataTypefOriginSize);
264     origin.Remove(TString::kTrailing, ' ');
265     origin.ToUpper();
266     fName="";
267     if (!fDirectory.IsNull()) {
268       fName+=fDirectory; fName+="/";
269     }
270     fName+="AliHLT"; fName+=origin;
271     if (fDt!=kAliHLTDataTypeESDObject &&
272         fDt!=kAliHLTDataTypeESDTree) {
273
274       HLTWarning("non-standard ESD type %s", AliHLTComponent::DataType2Text(fDt).c_str());
275       TString id;
276       id.Insert(0, fDt.fID, kAliHLTComponentDataTypefIDsize);
277       id.Remove(TString::kTrailing, ' ');
278       id.ToUpper();
279       fName+="_"; fName+=id; fName+=".root";
280     } else {
281       fName+="ESDs.root";
282     }
283
284     if (!gSystem->AccessPathName(fName)) {
285       // file exists, delete
286       TString shellcmd="rm -f ";
287       shellcmd+=fName;
288       gSystem->Exec(shellcmd);
289     }
290   }
291
292   TChain chain("esdTree");
293   TList cleanup;
294   cleanup.SetOwner();
295
296   int nofCurrentEvents=0;
297   if (iResult>=0) {
298     if (!gSystem->AccessPathName(fName)) {
299       // these are the other events, use the target file and temporary files to merge
300       // with TChain
301       chain.Add(fName);
302
303       if (eventno>=0) {
304         TFile file(fName);
305         if (!file.IsZombie()) {
306           TTree* pSrcTree;
307           file.GetObject("esdTree", pSrcTree);
308           if (pSrcTree) {
309             nofCurrentEvents=pSrcTree->GetEntries();
310           }
311           file.Close();
312         }
313       }
314     }
315   }
316
317   // synchronize and add empty events
318   if (nofCurrentEvents<eventno) {
319     iResult=1; // indicate files to merge
320     TTree* pTgtTree=new TTree("esdTree", "Tree with HLT ESD objects");
321     if (pTgtTree) {
322       pTgtTree->SetDirectory(0);
323       AliESDEvent* pTmpESD=new AliESDEvent;
324       if (pTmpESD) {
325         TString tmpfilename;
326         FILE* pTmpFile=gSystem->TempFileName(tmpfilename);
327         if (pTmpFile) {
328           fclose(pTmpFile);
329           pTmpFile=NULL;
330           cleanup.Add(new TObjString(tmpfilename));
331           TFile emptyevents(tmpfilename, "RECREATE");
332           if (!emptyevents.IsZombie()) {
333             pTmpESD->CreateStdContent();
334             pTmpESD->WriteToTree(pTgtTree);
335             HLTDebug("adding %d empty events to file %s", eventno-nofCurrentEvents, fName.Data());
336             for (int i=nofCurrentEvents; i<eventno; i++) {
337               pTgtTree->Fill();
338             }
339             pTgtTree->GetUserInfo()->Add(pTmpESD);
340             emptyevents.cd();
341             pTgtTree->Write();
342             emptyevents.Close();
343             chain.Add(tmpfilename);
344             pTgtTree->GetUserInfo()->Clear();
345           }
346         }
347         delete pTmpESD;
348       } else {
349         iResult=-ENOMEM;
350       }
351       delete pTgtTree;
352     } else {
353       iResult=-ENOMEM;
354     }
355   }
356
357   if (iResult>=0 && pSrcESD) {
358     // add the new event to the chain
359     iResult=1; // indicate files to merge
360     TString tmpfilename=WriteTempFile(pSrcESD);
361     if (!tmpfilename.IsNull()) {
362       chain.Add(tmpfilename);
363       cleanup.Add(new TObjString(tmpfilename));
364     }
365   }
366
367   if (iResult>0) {
368     // build temporary file name for chain output
369     TString tgtName;
370     FILE* pTmpFile=gSystem->TempFileName(tgtName);
371     if (pTmpFile) {
372       fclose(pTmpFile);
373       pTmpFile=NULL;
374
375       chain.Merge(tgtName);
376       // rename the merged file to the original file
377       TString shellcmd="mv ";
378       shellcmd+=tgtName + " " + fName;
379       if (gSystem->Exec(shellcmd)==0) {
380         HLTDebug("renaming %s to %s", tgtName.Data(), fName.Data());
381       } else {
382         HLTError("can not rename temporary file %s to %s", tgtName.Data(), fName.Data());
383       }
384     } else {
385       HLTError("can not get temporary file name from system");
386       iResult=-EBADF;
387     }
388   }
389
390   // delete temporary files
391   // the list objects are cleaned up be the TList destructor as the
392   // list is owner
393   TIter entry(&cleanup);
394   while (TObject* pObj=entry.Next()) {
395     if (dynamic_cast<TObjString*>(pObj)) {
396       TString shellcmd="rm -f ";
397       shellcmd+=(dynamic_cast<TObjString*>(pObj))->GetString();
398       gSystem->Exec(shellcmd);
399     }
400   }
401
402   return iResult;
403 }
404
405 TString AliHLTEsdManager::AliHLTEsdListEntry::WriteTempFile(AliESDEvent* pESD) const
406 {
407   // see header file for class documentation
408   int iResult;
409   TString tmpfilename;
410   FILE* pTmpFile=gSystem->TempFileName(tmpfilename);
411   if (pTmpFile) {
412     fclose(pTmpFile);
413     pTmpFile=NULL;
414
415     TFile file(tmpfilename, "RECREATE");
416     if (!file.IsZombie()) {
417       TTree* pTree=AliHLTEsdManager::EmbedIntoTree(pESD);
418       if (pTree) {
419         file.cd();
420         if (pTree->Write()>0) {
421         } else {
422           HLTError("can not write esd tree to temporary file %s", tmpfilename.Data());
423         }
424
425         pTree->GetUserInfo()->Clear();
426         delete pTree;
427       } else {
428         iResult=-ENOMEM;
429       }
430       file.Close();
431     } else {
432       HLTError("can not open file %s", tmpfilename.Data());
433     }
434   } else {
435     HLTError("can not get temporary file name from system");
436     iResult=-EBADF;
437   }
438
439   if (iResult<0) {
440     if (gSystem->AccessPathName(tmpfilename)==0) {
441       TString shellcmd="rm -f ";
442       shellcmd+=tmpfilename;
443       gSystem->Exec(shellcmd);
444     }
445     tmpfilename="";
446   }
447   return tmpfilename;
448 }
449
450 void AliHLTEsdManager::AliHLTEsdListEntry::SetDirectory(const char* directory)
451 {
452   // see header file for class documentation
453   if (!directory) return;
454   if (!fName.IsNull()) {
455     HLTWarning("ESD entry already in writing mode (%s), ignoring directory", fName.Data());
456     return;
457   }
458   fDirectory=directory;
459 }
460
461 void AliHLTEsdManager::AliHLTEsdListEntry::Delete()
462 {
463   // see header file for class documentation
464   if (fName.IsNull()) return;
465   if (gSystem->AccessPathName(fName)!=0) return;
466
467   TString shellcmd="rm -f ";
468   shellcmd+=fName;
469   gSystem->Exec(shellcmd);
470   fName="";
471 }
472
473 const char* AliHLTEsdManager::AliHLTEsdListEntry::GetFileName() const
474 {
475   // see header file for class documentation
476   return fName.Data();
477 }