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