]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/rec/AliHLTEsdManagerImplementation.cxx
adding methods for handling ESDEvent objects outside liHLTbase
[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 ///         This is an implementation of the abstract interface AliHLTEsdManager
24
25 #include "AliHLTEsdManagerImplementation.h"
26 #include "AliHLTComponent.h"
27 #include "AliESDEvent.h"
28 #include "AliHLTMessage.h"
29 #include "AliESDEvent.h"
30 #include "AliESDtrack.h"
31 #include "AliESDFMD.h"
32 #include "AliESDVZERO.h"
33 #include "AliESDTZERO.h"
34 #include "AliESDCaloCells.h"
35 #include "AliMultiplicity.h"
36 #include "AliESDACORDE.h"
37 #include "TFile.h"
38 #include "TTree.h"
39 #include "TClass.h"
40 #include "TObject.h"
41 #include "TList.h"
42
43 const float kAliESDVZERODefaultTime = -1024.;
44 const float kAliESDVZERODefaultTimeGlitch = 1e-6;
45 const float kAliESDZDCDefaultEMEnergy = 0.;
46 const float kAliESDZDCDefaultEMEnergyGlitch = 1e-6;
47
48 /** ROOT macro for the implementation of ROOT specific class methods */
49 ClassImp(AliHLTEsdManagerImplementation)
50
51 AliHLTEsdManagerImplementation::AliHLTEsdManagerImplementation()
52   :
53   fESDs()
54   , fDirectory()
55   , fTreeName()
56   , fWriteLocal(false)
57 {
58   // see header file for class documentation
59   // or
60   // refer to README to build package
61   // or
62   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
63
64   CheckClassConditions();
65 }
66
67 AliHLTEsdManagerImplementation::~AliHLTEsdManagerImplementation()
68 {
69   // see header file for class documentation
70   for (unsigned int i=0; i<fESDs.size(); i++) {
71     if (fESDs[i]) {
72       delete fESDs[i];
73     }
74     fESDs[i]=NULL;
75   }
76 }
77
78 int AliHLTEsdManagerImplementation::SetOption(const char* option)
79 {
80   // see header file for class documentation
81   int iResult=0;
82   TString strOptions=option;
83   TObjArray* pTokens=strOptions.Tokenize(" ");
84   if (pTokens) {
85     if (pTokens->GetEntriesFast()>0) {
86       for (int n=0; n<pTokens->GetEntriesFast(); n++) {
87         TString data=((TObjString*)pTokens->At(n))->GetString();
88         if (data.IsNull()) continue;
89
90         if (data.CompareTo("-writelocal")==0) {
91           fWriteLocal=true;
92         } else if (data.Contains("-directory=")) {
93           data.ReplaceAll("-directory=", "");
94           SetDirectory(data.Data());
95         } else if (data.Contains("-treename=")) {
96           data.ReplaceAll("-treename=", "");
97           fTreeName=data;
98         } else {
99           HLTError("unknown argument %s", data.Data());
100           iResult=-EINVAL;
101           break;
102         }
103       }
104     }
105     delete pTokens;
106   }
107   return iResult;
108 }
109
110 AliHLTEsdManagerImplementation::AliHLTEsdListEntry* AliHLTEsdManagerImplementation::Find(AliHLTComponentDataType dt) const
111 {
112   // see header file for class documentation
113   AliHLTEsdListEntry* pEntry=NULL;
114   for (unsigned int i=0; i<fESDs.size(); i++) {
115     if (fESDs[i] && *(fESDs[i])==dt) {
116       pEntry=const_cast<AliHLTEsdListEntry*>(fESDs[i]);
117     }
118   }
119   return pEntry;
120 }
121
122 int AliHLTEsdManagerImplementation::WriteESD(const AliHLTUInt8_t* pBuffer, AliHLTUInt32_t size,
123                                AliHLTComponentDataType dt, AliESDEvent* tgtesd, int eventno)
124 {
125   // see header file for class documentation
126   if (!pBuffer && size<=0) return -EINVAL;
127   int iResult=0;
128   AliHLTUInt32_t firstWord=*((AliHLTUInt32_t*)pBuffer);
129   if (firstWord==size-sizeof(AliHLTUInt32_t)) {
130     HLTDebug("create object from block %s size %d", AliHLTComponent::DataType2Text(dt).c_str(), size);
131     AliHLTMessage msg(const_cast<AliHLTUInt8_t*>(pBuffer), size);
132     TClass* objclass=msg.GetClass();
133     TObject* pObj=msg.ReadObject(objclass);
134     if (pObj && objclass) {
135       HLTDebug("object %p type %s created", pObj, objclass->GetName());
136       AliESDEvent* pESD=dynamic_cast<AliESDEvent*>(pObj);
137       TTree* pTree=NULL;
138       if (!pESD) {
139         pTree=dynamic_cast<TTree*>(pObj);
140         if (pTree) {
141           pESD=new AliESDEvent;
142           pESD->CreateStdContent();
143           if (pTree->GetEntries()>0) {
144             if (pTree->GetEntries()>1) {
145               HLTWarning("only one entry allowed for ESD embedded into tree, data block %s contains tree with %d entires, taking first entry",
146                          AliHLTComponent::DataType2Text(dt).c_str(), pTree->GetEntries());
147             }
148             pESD->ReadFromTree(pTree);
149             pTree->GetEvent(0);
150           }
151         } else {
152           HLTWarning("tree of data block %s has no events, skipping data block", AliHLTComponent::DataType2Text(dt).c_str());
153         }
154       }
155       if (pESD) {
156         AliHLTEsdListEntry* entry=Find(dt);
157         if (!entry) {
158           if ((entry=new AliHLTEsdListEntry(dt))!=NULL) {
159             if (!fDirectory.IsNull()) {
160               entry->SetDirectory(fDirectory);
161             }
162             if (!fTreeName.IsNull()) {
163               entry->SetTreeName(fTreeName);
164             }
165             fESDs.push_back(entry);
166           }
167         }
168         if (entry) {
169           if (tgtesd) {
170             Merge(tgtesd, pESD);
171           }
172
173           // Matthias 2009-06-06: writing of individual ESD files for the different origins was a
174           // first attempt when functionality was missing in the AliRoot framework and remained as
175           // debugging feature. ESD merging is now implemented and data written to the hltEsd, so
176           // the feature is now disabled by default because it causes increasing memory consumption.
177           // Presumably not because of a memory leak but the way the internal TTree is used and kept
178           // in memory.
179           // Writing of local files can be optionally switched on as e.g. by the EsdCollector component.
180           if (fWriteLocal) entry->WriteESD(pESD, eventno);
181         } else {
182           HLTError("internal mismatch, can not create list entry");
183           iResult=-ENOMEM;
184         }
185       } else {
186         HLTWarning("data block %s is not of class type AliESDEvent, ignoring ...", AliHLTComponent::DataType2Text(dt).c_str());
187       }
188       if (pTree) {
189         // ESD has been created and must be cleaned up
190         pESD->Reset();
191         delete pESD;
192         pESD=NULL;
193       }
194       delete pObj;
195       pObj=NULL;
196     } else {
197     }
198   }
199   return iResult;
200 }
201
202 int AliHLTEsdManagerImplementation::PadESDs(int eventno)
203 {
204   // see header file for class documentation
205   int iResult=0;
206   for (unsigned int i=0; i<fESDs.size(); i++) {
207     if (fESDs[i]) {
208       int res=fESDs[i]->WriteESD(NULL, eventno);
209       if (res<0 && iResult>=0) iResult=res;
210     }
211   }
212   return iResult;
213 }
214
215 void AliHLTEsdManagerImplementation::SetDirectory(const char* directory)
216 {
217   // see header file for class documentation
218   if (!directory) return;
219   fDirectory=directory;
220   for (unsigned int i=0; i<fESDs.size(); i++) {
221     if (fESDs[i]) {
222       fESDs[i]->SetDirectory(directory);
223     }
224   }
225 }
226
227 TString AliHLTEsdManagerImplementation::GetFileNames(AliHLTComponentDataType dt) const
228 {
229   TString result;
230   for (unsigned int i=0; i<fESDs.size(); i++) {
231     if (fESDs[i] && *(fESDs[i])==dt) {
232       if (!result.IsNull()) result+=" ";
233       result+=fESDs[i]->GetFileName();
234     }
235   }
236   return result;
237 }
238
239 TTree* AliHLTEsdManagerImplementation::EmbedIntoTree(AliESDEvent* pESD, const char* name, const char* title)
240 {
241   // see header file for class documentation
242   int iResult=0;
243   TTree* pTree=new TTree(name, title);
244   if (pTree) {
245     pESD->WriteToTree(pTree);
246     pTree->Fill();
247     pTree->GetUserInfo()->Add(pESD);
248   } else {
249     iResult=-ENOMEM;
250   }
251
252   if (iResult<0) {
253     pTree->GetUserInfo()->Clear();
254     delete pTree;
255   }
256
257   return pTree;
258 }
259
260 AliHLTEsdManagerImplementation::AliHLTEsdListEntry::AliHLTEsdListEntry(AliHLTComponentDataType dt)
261   :
262   fName(),
263   fDirectory(),
264   fDt(dt),
265   fpFile(NULL),
266   fpTree(NULL),
267   fpEsd(NULL),
268   fPrefix(),
269   fTreeName("esdTree")
270 {
271   // see header file for class documentation
272 }
273
274 AliHLTEsdManagerImplementation::AliHLTEsdListEntry::~AliHLTEsdListEntry()
275 {
276   // see header file for class documentation
277   if (fpEsd) delete fpEsd;
278   fpEsd=NULL;
279
280   if (fpTree) delete fpTree;
281   fpTree=NULL;
282
283   if (fpFile) {
284     fpFile->Close();
285     delete fpFile;
286   }
287   fpFile=NULL;
288 }
289
290 bool AliHLTEsdManagerImplementation::AliHLTEsdListEntry::operator==(AliHLTComponentDataType dt) const
291 {
292   // see header file for class documentation
293   return fDt==dt;
294 }
295
296 int AliHLTEsdManagerImplementation::AliHLTEsdListEntry::WriteESD(AliESDEvent* pSrcESD, int eventno)
297 {
298   // see header file for class documentation
299   int iResult=0;
300
301   if (fName.IsNull()) {
302     // this is the first event, create the file name
303     fName="";
304     if (!fDirectory.IsNull()) {
305       fName+=fDirectory; fName+="/";
306     }
307     fName+="Ali"; fName+=GetPrefix();
308     if (fDt!=kAliHLTDataTypeESDObject &&
309         fDt!=kAliHLTDataTypeESDTree) {
310
311       HLTWarning("non-standard ESD type %s", AliHLTComponent::DataType2Text(fDt).c_str());
312       TString id;
313       id.Insert(0, fDt.fID, kAliHLTComponentDataTypefIDsize);
314       id.Remove(TString::kTrailing, ' ');
315       id.ToUpper();
316       fName+="_"; fName+=id; fName+=".root";
317     } else {
318       fName+="ESDs.root";
319     }
320
321     fpFile=new TFile(fName, "RECREATE");
322     fpTree=new TTree(fTreeName, "Tree with HLT ESD objects");
323     fpTree->SetDirectory(0);
324     fpEsd=new AliESDEvent;
325     if (fpEsd) {
326       fpEsd->CreateStdContent();
327       *fpEsd=*pSrcESD;
328       if (fpTree) {
329         fpEsd->WriteToTree(fpTree);
330       }
331     }
332   }
333
334   if (fpFile && fpTree && fpEsd) {
335     // synchronize and add empty events
336     fpEsd->Reset();
337     int nofCurrentEvents=fpTree->GetEntries();
338     if (nofCurrentEvents<eventno) {
339       iResult=1; // indicate tree to be written
340       HLTDebug("adding %d empty events to file %s", eventno-nofCurrentEvents, fName.Data());
341       for (int i=nofCurrentEvents; i<eventno; i++) {
342         fpTree->Fill();
343       }
344     }
345
346     if (iResult>=0 && pSrcESD) {
347       int nofObjects=fpEsd->GetList()->GetEntries();
348       *fpEsd=*pSrcESD;
349       if (nofObjects!=fpEsd->GetList()->GetEntries()) {
350         // The source ESD contains object not present in the target ESD
351         // before. Those objects will not be written to the tree since
352         // the branch layout has been created earlier.
353         // Create new tree with the additional branches, copy the entries
354         // of the current tree into the new tree, and continue.
355         TTree* pNewTree=new TTree(fTreeName, "Tree with HLT ESD objects");
356         pNewTree->SetDirectory(0);
357         AliESDEvent* readESD=new AliESDEvent;
358         readESD->CreateStdContent();
359         readESD->ReadFromTree(fpTree);
360         fpEsd->Reset();
361         fpEsd->WriteToTree(pNewTree);
362         HLTDebug("cloning tree with %d entries", fpTree->GetEntries());
363         for (int event=0; event<fpTree->GetEntries(); event++) {
364           fpTree->GetEntry(event);
365           *fpEsd=*readESD;
366           pNewTree->Fill();
367           fpEsd->Reset();
368         }
369         fpFile->Close();
370         delete fpFile;
371         delete readESD;
372         delete fpTree;
373         fpFile=new TFile(fName, "RECREATE");
374         fpTree=pNewTree;
375         *fpEsd=*pSrcESD;
376         HLTDebug("new ESD with %d objects", fpEsd->GetList()->GetEntries());
377       }
378       fpTree->Fill();
379       iResult=1; // indicate tree to be written
380     }
381
382     if (iResult>0) {
383       fpFile->cd();
384       fpTree->GetUserInfo()->Add(fpEsd);
385       fpTree->Write(fpTree->GetName(),TObject::kOverwrite);
386       fpTree->GetUserInfo()->Clear();
387     }
388   }
389
390   return iResult;
391 }
392
393 void AliHLTEsdManagerImplementation::AliHLTEsdListEntry::SetDirectory(const char* directory)
394 {
395   // see header file for class documentation
396   if (!directory) return;
397   if (!fName.IsNull()) {
398     HLTWarning("ESD entry already in writing mode (%s), ignoring directory", fName.Data());
399     return;
400   }
401   fDirectory=directory;
402 }
403
404 const char* AliHLTEsdManagerImplementation::AliHLTEsdListEntry::GetFileName() const
405 {
406   // see header file for class documentation
407   return fName.Data();
408 }
409
410 const char* AliHLTEsdManagerImplementation::AliHLTEsdListEntry::GetPrefix()
411 {
412   // see header file for class documentation
413   if (fPrefix.IsNull()) {
414     fPrefix.Insert(0, fDt.fOrigin, kAliHLTComponentDataTypefOriginSize);
415     fPrefix.Remove(TString::kTrailing, ' ');
416     fPrefix.ToUpper();
417     if (!fPrefix.Contains("HLT")) {
418       fPrefix.Insert(0, "HLT");
419     }
420   }
421   return fPrefix.Data();
422 }
423
424 int AliHLTEsdManagerImplementation::Merge(AliESDEvent* pTgt, AliESDEvent* pSrc) const
425 {
426   // see header file for class documentation
427   int iResult=0;
428   if (!pTgt || !pSrc) return -EINVAL;
429
430   TIter next(pSrc->GetList());
431   TObject* pSrcObject=NULL;
432   static int warningCount=0;
433   while ((pSrcObject=next())) {
434     if(!pSrcObject->InheritsFrom("TCollection")){
435       // simple objects
436       // for every type of object we have to find out whether it is empty or not
437       // objects are only copied if non-empty, otherwise valid content would be
438       // overridden by empty objects during the merging
439       bool copy=false;
440       TString name=pSrcObject->GetName();
441       if(pSrcObject->InheritsFrom("AliHLTTriggerDecision")){
442         copy=true;
443       } else if (pSrcObject->IsA()==AliESDRun::Class()) {
444         AliESDRun* pESDRun=dynamic_cast<AliESDRun*>(pSrcObject);
445         // zero might be a valid run no in simulation, so we check also whether the CTP trigger classes are set
446         copy=(pESDRun && (pESDRun->GetRunNumber()>0 || !pESDRun->GetActiveTriggerClasses().IsNull()));
447       } else if (pSrcObject->IsA()==AliESDHeader::Class()) {
448         AliESDHeader* pESDHeader=dynamic_cast<AliESDHeader*>(pSrcObject);
449         copy=(pESDHeader && pESDHeader->GetTriggerMask()!=0);
450       } else if (pSrcObject->IsA()==AliESDVertex::Class()) {
451         AliESDVertex* pESDVertex=dynamic_cast<AliESDVertex*>(pSrcObject);
452         copy=(pESDVertex && pESDVertex->GetNContributors()>0);
453       } else if (pSrcObject->IsA()==AliESDTZERO::Class()) {
454         AliESDTZERO* pESDTZero=dynamic_cast<AliESDTZERO*>(pSrcObject);
455         copy=(pESDTZero && (pESDTZero->GetT0zVertex()!=0.0 || pESDTZero->GetT0()!=0.0));
456       } else if (pSrcObject->IsA()==AliESDVZERO::Class()) {
457         AliESDVZERO* pESDVZERO=dynamic_cast<AliESDVZERO*>(pSrcObject);
458         // object contains data if one of the times exceeds the default value
459         copy=pESDVZERO &&
460           (pESDVZERO->GetV0ATime() > kAliESDVZERODefaultTime+kAliESDVZERODefaultTimeGlitch ||
461            pESDVZERO->GetV0CTime() > kAliESDVZERODefaultTime+kAliESDVZERODefaultTimeGlitch);
462       } else if (pSrcObject->IsA()==AliESDFMD::Class()) {
463         AliESDFMD* pESDFMD=dynamic_cast<AliESDFMD*>(pSrcObject);
464         copy=(pESDFMD && false); // have to find an easy valid condition
465       } else if (pSrcObject->IsA()==AliESDZDC::Class()) {
466         AliESDZDC* pESDZDC=dynamic_cast<AliESDZDC*>(pSrcObject);
467         // object contains data if the EM energies are different from the default value
468         copy=pESDZDC &&
469           (TMath::Abs(pESDZDC->GetZDCEMEnergy(0)-kAliESDZDCDefaultEMEnergy) > kAliESDZDCDefaultEMEnergyGlitch ||
470            TMath::Abs(pESDZDC->GetZDCEMEnergy(1)-kAliESDZDCDefaultEMEnergy) > kAliESDZDCDefaultEMEnergyGlitch);
471       } else if (pSrcObject->IsA()==AliMultiplicity::Class()) {
472         AliMultiplicity* pMultiplicity=dynamic_cast<AliMultiplicity*>(pSrcObject);
473         copy=(pMultiplicity && pMultiplicity->GetNumberOfTracklets()>0);
474       } else if (pSrcObject->IsA()==AliESDCaloTrigger::Class()) {
475         AliESDCaloTrigger* pESDCaloTrigger=dynamic_cast<AliESDCaloTrigger*>(pSrcObject);
476         copy=(pESDCaloTrigger && false); // have to find an easy valid condition
477       } else if (pSrcObject->IsA()==AliESDCaloCells::Class()) {
478         AliESDCaloCells* pESDCaloCells=dynamic_cast<AliESDCaloCells*>(pSrcObject);
479         copy=(pESDCaloCells && false); // have to find an easy valid condition
480       } else if (pSrcObject->IsA()==AliESDACORDE::Class()) {
481         AliESDACORDE* pESDACORDE=dynamic_cast<AliESDACORDE*>(pSrcObject);
482         copy=(pESDACORDE && false); // have to find an easy valid condition
483       } else if (!AliHLTESDEventHelper::IsStdContent(name)) {
484         // this is likely to be ok as long as it is not any object of the std content.
485         copy=true;
486       } else {
487         HLTError("no merging implemented for object %s, omitting", name.Data());
488       }
489       if (copy) {
490         //pSrcObject->Print();
491         TObject* pTgtObject=pTgt->FindListObject(name);
492         if (pTgtObject) {
493           pSrcObject->Copy(*pTgtObject);
494         } else {
495           pTgt->AddObject(pSrcObject->Clone());
496         }
497       }
498     } else if(pSrcObject->InheritsFrom("TClonesArray")){
499       TClonesArray* pTClA=dynamic_cast<TClonesArray*>(pSrcObject);
500       if (pTClA!=NULL && pTClA->GetEntriesFast()>0) {
501         TString name=pTClA->GetName();
502         TObject* pTgtObject=pTgt->GetList()->FindObject(name);
503         TClonesArray* pTgtArray=NULL;
504         if (pTgtObject!=NULL && pTgtObject->InheritsFrom("TClonesArray")){
505           pTgtArray=dynamic_cast<TClonesArray*>(pTgtObject);
506           if (pTgtArray) {
507             TString classType=pTClA->Class()->GetName();
508             if (classType.CompareTo(pTgtArray->Class()->GetName())==0) {
509               if (pTgtArray->GetEntries()==0) {
510                 pTgtArray->ExpandCreate(pTClA->GetEntries());
511                 for(int i=0; i<pTClA->GetEntriesFast(); ++i){
512                   (*pTClA)[i]->Copy(*((*pTgtArray)[i]));
513                 }
514               } else {
515                 if (warningCount++<10) {
516                   HLTWarning("TClonesArray \"%s\"  in target ESD %p is already filled with %d entries",
517                              name.Data(), pTgt, pTgtArray->GetEntries());
518                 }
519                 iResult=-EBUSY;
520               }
521             } else {
522               if (warningCount++<10) {
523                 HLTWarning("TClonesArray \"%s\" exists in target ESD %p, but describes incompatible class type %s instead of %s",
524                            name.Data(), pTgt, pTgtArray->GetClass()->GetName(), pTClA->GetClass()->GetName());
525               }
526               iResult=-EBUSY;
527             }
528           } else {
529             if (warningCount++<10) {
530               HLTError("internal error: dynamic cast failed for object %s %p", pTgtObject->GetName(), pTgtObject);
531             }
532             iResult=-EBUSY;
533           }
534         } else if (pTgtObject) {
535           if (warningCount++<10) {
536             HLTWarning("object \"%s\" does already exist in target ESD %p, but is %s rather than TClonesArray"
537                        " skipping data",
538                        name.Data(), pTgt, pTgtObject->Class()->GetName());
539           }
540           iResult=-EBUSY;
541         } else {
542           if (warningCount++<10) {
543             HLTWarning("object \"%s\" does not exist in target ESD %p, data can not be copied because it will be lost when filling the tree",
544                        name.Data(), pTgt);
545           }
546           iResult=-ENOENT;
547         }
548       }
549     }
550   }
551   return iResult;
552 }
553
554 bool AliHLTEsdManagerImplementation::AliHLTESDEventHelper::IsStdContent(const char* key)
555 {
556   // check if the key denotes a std object
557   TString needle=key;
558   for (int i=0; i<kESDListN; i++) {
559     if (needle.CompareTo(fgkESDListName[i])==0) return true;
560   }
561   return false;
562 }
563
564 int AliHLTEsdManagerImplementation::CheckClassConditions() const
565 {
566   // this is a helper method which checks if some thing in the default
567   // object initialization changes during the evolution of the source code
568   // which is necessary for checking the validity of data in the object
569
570   // check AliESDVZERO
571   AliESDVZERO vzerodummy;
572   if (TMath::Abs(vzerodummy.GetV0ATime()-kAliESDVZERODefaultTime) > kAliESDVZERODefaultTimeGlitch ||
573       TMath::Abs(vzerodummy.GetV0CTime()-kAliESDVZERODefaultTime) > kAliESDVZERODefaultTimeGlitch) {
574     HLTWarning("initialization of AliESDVZERO has changed, default time values now %f/%f, "
575                "revision of condition might be needed",
576                vzerodummy.GetV0ATime(), vzerodummy.GetV0CTime());
577   }
578
579   // check AliESDZDC
580   AliESDZDC zdcdummy;
581   if (TMath::Abs(zdcdummy.GetZDCEMEnergy(0)-kAliESDZDCDefaultEMEnergy) > kAliESDZDCDefaultEMEnergyGlitch ||
582       TMath::Abs(zdcdummy.GetZDCEMEnergy(1)-kAliESDZDCDefaultEMEnergy) > kAliESDZDCDefaultEMEnergyGlitch) {
583     HLTWarning("initialization of AliESDZDC has changed, default em energy values now %f/%f, "
584                "revision of condition might be needed",
585                zdcdummy.GetZDCEMEnergy(0), zdcdummy.GetZDCEMEnergy(1));
586   }
587
588   return 0;
589 }
590
591 AliESDEvent* AliHLTEsdManagerImplementation::CreateEsdEvent(bool bCreateStdContent) const
592 {
593   // create ESDEvent object and optionally initialize standard content
594   AliESDEvent* pESD=new AliESDEvent;
595   if (pESD && bCreateStdContent) pESD->CreateStdContent();
596   return pESD;
597 }
598
599 int AliHLTEsdManagerImplementation::AddObject(AliESDEvent* pESD, const TObject* pObject, const char* branchname) const
600 {
601   // add object to the list of ESD contributors
602   if (!pESD || !pObject) return -EINVAL;
603   TObject* pESDObject=pESD->FindListObject(branchname);
604   if (pESDObject) {
605     // copy the content to the already existing object
606     pObject->Copy(*pESDObject);
607   } else {
608     // add a new object
609     pESD->AddObject(pObject->Clone());
610   }
611   return 0;
612 }
613
614 int AliHLTEsdManagerImplementation::ResetEsdEvent(AliESDEvent* pESD) const
615 {
616   // reset the specified ESD object
617   if (!pESD) return -EINVAL;
618   pESD->Reset();
619   return 0;
620 }