]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/global/AliHLTGlobalEsdConverterComponent.cxx
general code audit, added documentation about AliESDtrack update sequence
[u/mrichter/AliRoot.git] / HLT / global / AliHLTGlobalEsdConverterComponent.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   AliHLTGlobalEsdConverterComponent.cxx
20 //  @author Matthias Richter
21 //  @date   
22 //  @brief  Global ESD converter component.
23 // 
24
25 // see header file for class documentation
26 // or
27 // refer to README to build package
28 // or
29 // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
30
31 #include <cassert>
32 #include "AliHLTGlobalEsdConverterComponent.h"
33 #include "AliHLTGlobalBarrelTrack.h"
34 #include "AliHLTExternalTrackParam.h"
35 #include "AliHLTTrackMCLabel.h"
36 #include "AliHLTCTPData.h"
37 #include "AliESDEvent.h"
38 #include "AliESDtrack.h"
39 #include "AliESDMuonTrack.h"
40 #include "AliCDBEntry.h"
41 #include "AliCDBManager.h"
42 #include "AliPID.h"
43 #include "TTree.h"
44 #include "TList.h"
45 #include "TClonesArray.h"
46 #include "AliHLTESDCaloClusterMaker.h"
47 #include "AliHLTCaloClusterDataStruct.h"
48 #include "AliHLTCaloClusterReader.h"
49 #include "AliESDCaloCluster.h"
50 #include "AliHLTGlobalVertexerComponent.h"
51
52 /** ROOT macro for the implementation of ROOT specific class methods */
53 ClassImp(AliHLTGlobalEsdConverterComponent)
54
55 AliHLTGlobalEsdConverterComponent::AliHLTGlobalEsdConverterComponent()
56   : AliHLTProcessor()
57   , fWriteTree(0)
58   , fVerbosity(0)
59   , fESD(NULL)
60   , fSolenoidBz(-5.00668)
61   , fBenchmark("EsdConverter")
62 {
63   // see header file for class documentation
64   // or
65   // refer to README to build package
66   // or
67   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
68 }
69
70 AliHLTGlobalEsdConverterComponent::~AliHLTGlobalEsdConverterComponent()
71 {
72   // see header file for class documentation
73   if (fESD) delete fESD;
74   fESD=NULL;
75 }
76
77 int AliHLTGlobalEsdConverterComponent::Configure(const char* arguments)
78 {
79   // see header file for class documentation
80   int iResult=0;
81   if (!arguments) return iResult;
82
83   TString allArgs=arguments;
84   TString argument;
85   int bMissingParam=0;
86
87   TObjArray* pTokens=allArgs.Tokenize(" ");
88   if (pTokens) {
89     for (int i=0; i<pTokens->GetEntries() && iResult>=0; i++) {
90       argument=((TObjString*)pTokens->At(i))->GetString();      
91       if (argument.IsNull()) continue;
92       
93       if (argument.CompareTo("-solenoidBz")==0) {
94         if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
95         HLTWarning("argument -solenoidBz is deprecated, magnetic field set up globally (%f)", GetBz());
96         continue;
97       } else {
98         HLTError("unknown argument %s", argument.Data());
99         iResult=-EINVAL;
100         break;
101       }
102     }
103     delete pTokens;
104   }
105   if (bMissingParam) {
106     HLTError("missing parameter for argument %s", argument.Data());
107     iResult=-EINVAL;
108   }
109
110   return iResult;
111 }
112
113 int AliHLTGlobalEsdConverterComponent::Reconfigure(const char* cdbEntry, const char* chainId)
114 {
115   // see header file for class documentation
116   int iResult=0;
117   const char* path=NULL;
118   const char* defaultNotify="";
119   if (cdbEntry) {
120     path=cdbEntry;
121     defaultNotify=" (default)";
122   }
123   if (path) {
124     HLTInfo("reconfigure from entry %s%s, chain id %s", path, defaultNotify,(chainId!=NULL && chainId[0]!=0)?chainId:"<none>");
125     AliCDBEntry *pEntry = AliCDBManager::Instance()->Get(path/*,GetRunNo()*/);
126     if (pEntry) {
127       TObjString* pString=dynamic_cast<TObjString*>(pEntry->GetObject());
128       if (pString) {
129         HLTInfo("received configuration object string: \'%s\'", pString->GetString().Data());
130         iResult=Configure(pString->GetString().Data());
131       } else {
132         HLTError("configuration object \"%s\" has wrong type, required TObjString", path);
133       }
134     } else {
135       HLTError("can not fetch object \"%s\" from CDB", path);
136     }
137   }
138   
139   return iResult;
140 }
141
142 void AliHLTGlobalEsdConverterComponent::GetInputDataTypes(AliHLTComponentDataTypeList& list)
143 {
144   // see header file for class documentation
145   list.push_back(kAliHLTDataTypeTrack);
146   list.push_back(kAliHLTDataTypeTrackMC);
147   list.push_back(kAliHLTDataTypeCaloCluster);
148   list.push_back(kAliHLTDataTypedEdx );
149   list.push_back(kAliHLTDataTypeESDVertex );
150   list.push_back(kAliHLTDataTypeESDObject);
151   list.push_back(kAliHLTDataTypeTObject);
152   list.push_back(kAliHLTDataTypeGlobalVertexer);
153 }
154
155 AliHLTComponentDataType AliHLTGlobalEsdConverterComponent::GetOutputDataType()
156 {
157   // see header file for class documentation
158   return kAliHLTDataTypeESDObject|kAliHLTDataOriginOut;
159 }
160
161 void AliHLTGlobalEsdConverterComponent::GetOutputDataSize(unsigned long& constBase, double& inputMultiplier)
162 {
163   // see header file for class documentation
164   constBase=2000000;
165   inputMultiplier=10.0;
166 }
167
168 int AliHLTGlobalEsdConverterComponent::DoInit(int argc, const char** argv)
169 {
170   // see header file for class documentation
171   int iResult=0;
172   TString argument="";
173   int bMissingParam=0;
174
175   // default list of skiped ESD objects
176   TString skipObjects=
177     // "AliESDRun,"
178     // "AliESDHeader,"
179     "AliESDZDC,"
180     "AliESDFMD,"
181     // "AliESDVZERO,"
182     // "AliESDTZERO,"
183     // "TPCVertex,"
184     // "SPDVertex,"
185     // "PrimaryVertex,"
186     // "AliMultiplicity,"
187     // "PHOSTrigger,"
188     // "EMCALTrigger,"
189     // "SPDPileupVertices,"
190     // "TrkPileupVertices,"
191     "Cascades,"
192     "Kinks,"
193     "AliRawDataErrorLogs,"
194     "AliESDACORDE";
195
196   iResult=Reconfigure(NULL, NULL);
197   TString allArgs = "";
198   for ( int i = 0; i < argc; i++ ) {
199     if ( !allArgs.IsNull() ) allArgs += " ";
200     allArgs += argv[i];
201   }
202
203   TObjArray* pTokens=allArgs.Tokenize(" ");
204   if (pTokens) {
205     for (int i=0; i<pTokens->GetEntries() && iResult>=0; i++) {
206       argument=((TObjString*)pTokens->At(i))->GetString();      
207       if (argument.IsNull()) continue;
208
209       // -notree
210       if (argument.CompareTo("-notree")==0) {
211         fWriteTree=0;
212         
213         // -tree
214       } else if (argument.CompareTo("-tree")==0) {
215         fWriteTree=1;
216       } else if (argument.CompareTo("-solenoidBz")==0) {
217         if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
218         HLTInfo("Magnetic Field set to: %s", ((TObjString*)pTokens->At(i))->GetString().Data());
219         HLTWarning("argument '-solenoidBz' is deprecated, solenoid field initiaized from CDB settings");
220         continue;
221       } else if (argument.Contains("-skipobject=")) {
222         argument.ReplaceAll("-skipobject=", "");
223         skipObjects=argument;
224       } else {
225         HLTError("unknown argument %s", argument.Data());
226         iResult=-EINVAL;
227         break;
228       }
229     }
230   }
231   if (bMissingParam) {
232     HLTError("missing parameter for argument %s", argument.Data());
233     iResult=-EINVAL;
234   }
235
236   fSolenoidBz=GetBz();
237
238   if (iResult>=0) {
239     fESD = new AliESDEvent;
240     if (fESD) {
241       fESD->CreateStdContent();
242
243       // remove some of the objects which are not needed
244       if (fESD->GetList() && !skipObjects.IsNull()) {
245         pTokens=skipObjects.Tokenize(",");
246         if (pTokens) {
247           const char* id=NULL;
248           TIter next(pTokens);
249           TObject* pObject=NULL;
250           while ((pObject=next())!=NULL) {
251             id=((TObjString*)pObject)->GetString().Data();
252             TObject* pObj=fESD->GetList()->FindObject(id);
253             if (pObj) {
254               HLTDebug("removing object %s", id);
255               fESD->GetList()->Remove(pObj);
256               delete pObj;
257             } else {
258               HLTWarning("failed to remove object '%s' from ESD", id);
259             }
260           }
261           fESD->GetStdContent();
262           delete pTokens;
263         }
264       }
265     } else {
266       iResult=-ENOMEM;
267     }
268
269     SetupCTPData();
270   }
271
272   fBenchmark.SetTimer(0,"total");
273
274   return iResult;
275 }
276
277 int AliHLTGlobalEsdConverterComponent::DoDeinit()
278 {
279   // see header file for class documentation
280   if (fESD) delete fESD;
281   fESD=NULL;
282
283   return 0;
284 }
285
286 int AliHLTGlobalEsdConverterComponent::DoEvent(const AliHLTComponentEventData& /*evtData*/, 
287                                                AliHLTComponentTriggerData& trigData)
288 {
289   // see header file for class documentation
290   int iResult=0;
291   if (!fESD) return -ENODEV;
292
293   if (IsDataEvent()) fBenchmark.StartNewEvent();
294   fBenchmark.Start(0);
295
296   AliESDEvent* pESD = fESD;
297
298   pESD->Reset(); 
299   pESD->SetMagneticField(fSolenoidBz);
300   pESD->SetRunNumber(GetRunNo());
301   pESD->SetPeriodNumber(GetPeriodNumber());
302   pESD->SetOrbitNumber(GetOrbitNumber());
303   pESD->SetBunchCrossNumber(GetBunchCrossNumber());
304   pESD->SetTimeStamp(GetTimeStamp());
305
306   const AliHLTCTPData* pCTPData=CTPData();
307   if (pCTPData) {
308     AliHLTUInt64_t mask=pCTPData->ActiveTriggers(trigData);
309     for (int index=0; index<gkNCTPTriggerClasses; index++) {
310       if ((mask&((AliHLTUInt64_t)0x1<<index)) == 0) continue;
311       pESD->SetTriggerClass(pCTPData->Name(index), index);
312     }
313     pESD->SetTriggerMask(mask);
314   }
315
316   TTree* pTree = NULL;
317   if (fWriteTree)
318     pTree = new TTree("esdTree", "Tree with HLT ESD objects");
319  
320   if (pTree) {
321     pTree->SetDirectory(0);
322   }
323
324   if ((iResult=ProcessBlocks(pTree, pESD))>=0) {
325     // TODO: set the specification correctly
326     if (pTree) {
327       // the esd structure is written to the user info and is
328       // needed in te ReadFromTree method to read all objects correctly
329       pTree->GetUserInfo()->Add(pESD);
330       pESD->WriteToTree(pTree);
331       iResult=PushBack(pTree, kAliHLTDataTypeESDTree|kAliHLTDataOriginOut, 0);
332     } else {
333       iResult=PushBack(pESD, kAliHLTDataTypeESDObject|kAliHLTDataOriginOut, 0);
334     }
335     fBenchmark.AddOutput(GetLastObjectSize());
336   }
337   if (pTree) {
338     // clear user info list to prevent objects from being deleted
339     pTree->GetUserInfo()->Clear();
340     delete pTree;
341   }
342
343   fBenchmark.Stop(0);
344   HLTInfo( fBenchmark.GetStatistics() );
345
346   return iResult;
347 }
348
349 int AliHLTGlobalEsdConverterComponent::ProcessBlocks(TTree* pTree, AliESDEvent* pESD)
350 {
351   // see header file for class documentation
352
353   int iResult=0;
354   int iAddedDataBlocks=0;
355
356   // Barrel tracking
357   // tracks are based on the TPC tracks, and only updated from the ITS information
358   // Sequence:
359   // 1) extract MC information for TPC and ITS from specific data blocks and store in
360   //    intermediate vector arrays
361   // 2) extract TPC tracks, update with MC labels if available, the track parameters
362   //    are estimated at the first cluster position
363   // 2.1) propagate to last cluster position and update kTPCout, sets also outer param (fOp)
364   // 2.2) update kTPCin, sets also inner param (fIp) and TPC inner param (fTPCInner)
365   // 2.3) update kTPCrefit using the same parameters at the first cluster position
366   //      HLT has strictly spoking no refit, but we want the flag to be set
367   //      can be changed to be done after all the individual barrel detector parameters
368   //      have been updated by looping over the tracks again
369   // 3) extract ITS tracks, the tracks are actually TPC tracks updated from the ITS
370   //    tracking information
371   // 3.1) TODO 2010-07-12: handle ITS standalone tracks by updating kITSout before kITSin
372   // 3.2) update with kITSin
373   //    TODO 2010-07-12 find out if the kITSrefit has to be set as well
374   // 4) extract TRD tracks and add to ESD
375   //    TODO 2010-07-12 at the moment there is no matching or merging of TPC and TRD tracks
376
377   // 1) first read MC information (if present)
378   std::map<int,int> mcLabelsTPC;
379   std::map<int,int> mcLabelsITS;
380
381   for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeTrackMC|kAliHLTDataOriginTPC);
382        pBlock!=NULL; pBlock=GetNextInputBlock()) {
383
384     fBenchmark.AddInput(pBlock->fSize);
385
386     AliHLTTrackMCData* dataPtr = reinterpret_cast<AliHLTTrackMCData*>( pBlock->fPtr );
387     if (sizeof(AliHLTTrackMCData)+dataPtr->fCount*sizeof(AliHLTTrackMCLabel)==pBlock->fSize) {
388       for( unsigned int il=0; il<dataPtr->fCount; il++ ){
389         AliHLTTrackMCLabel &lab = dataPtr->fLabels[il];
390         mcLabelsTPC[lab.fTrackID] = lab.fMCLabel;
391       }
392     } else {
393       HLTWarning("data mismatch in block %s (0x%08x): count %d, size %d -> ignoring track MC information", 
394                  DataType2Text(pBlock->fDataType).c_str(), pBlock->fSpecification, 
395                  dataPtr->fCount, pBlock->fSize);
396     }
397   }
398  
399   for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeTrackMC|kAliHLTDataOriginITS);
400        pBlock!=NULL; pBlock=GetNextInputBlock()) {
401
402     fBenchmark.AddInput(pBlock->fSize);
403
404     AliHLTTrackMCData* dataPtr = reinterpret_cast<AliHLTTrackMCData*>( pBlock->fPtr );
405     if (sizeof(AliHLTTrackMCData)+dataPtr->fCount*sizeof(AliHLTTrackMCLabel)==pBlock->fSize) {
406       for( unsigned int il=0; il<dataPtr->fCount; il++ ){
407         AliHLTTrackMCLabel &lab = dataPtr->fLabels[il];
408         mcLabelsITS[lab.fTrackID] = lab.fMCLabel;
409       }
410     } else {
411       HLTWarning("data mismatch in block %s (0x%08x): count %d, size %d -> ignoring track MC information", 
412                  DataType2Text(pBlock->fDataType).c_str(), pBlock->fSpecification, 
413                  dataPtr->fCount, pBlock->fSize);
414     }
415   }
416
417
418   // read dEdx information (if present)
419   // TODO 2010-07-12 this needs to be verified
420
421   AliHLTFloat32_t *dEdxTPC = 0; 
422   Int_t ndEdxTPC = 0;
423   for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypedEdx|kAliHLTDataOriginTPC);
424        pBlock!=NULL; pBlock=GetNextInputBlock()) {
425     fBenchmark.AddInput(pBlock->fSize);
426     dEdxTPC = reinterpret_cast<AliHLTFloat32_t*>( pBlock->fPtr );
427     ndEdxTPC = pBlock->fSize / sizeof(AliHLTFloat32_t);
428     break;
429   }
430
431   // 2) convert the TPC tracks to ESD tracks
432   for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeTrack|kAliHLTDataOriginTPC);
433        pBlock!=NULL; pBlock=GetNextInputBlock()) {
434     fBenchmark.AddInput(pBlock->fSize);
435     vector<AliHLTGlobalBarrelTrack> tracks;
436     if ((iResult=AliHLTGlobalBarrelTrack::ConvertTrackDataArray(reinterpret_cast<const AliHLTTracksData*>(pBlock->fPtr), pBlock->fSize, tracks))>=0) {
437       for (vector<AliHLTGlobalBarrelTrack>::iterator element=tracks.begin();
438            element!=tracks.end(); element++) {
439         Float_t points[4] = {
440           element->GetX(),
441           element->GetY(),
442           element->GetLastPointX(),
443           element->GetLastPointY()
444         };
445
446         Int_t mcLabel = -1;
447         if( mcLabelsTPC.find(element->TrackID())!=mcLabelsTPC.end() )
448           mcLabel = mcLabelsTPC[element->TrackID()];
449         element->SetLabel( mcLabel );
450
451         AliESDtrack iotrack;
452
453         // for the moment, the number of clusters are not set when processing the
454         // kTPCin update, only at kTPCout
455         // there ar emainly three parameters updated for kTPCout
456         //   number of clusters
457         //   chi2
458         //   pid signal
459         // The first one can be updated already at that stage here, while the two others
460         // eventually require to update from the ITS tracks before. The exact scheme
461         // needs to be checked
462         iotrack.SetID( element->TrackID() );
463
464         // 2.1 set kTPCout
465         // TODO 2010-07-12 disabled for the moment because of bug
466         // https://savannah.cern.ch/bugs/index.php?69875
467         // the propagation does not work, there is an offset in z
468         // propagate to last cluster and update parameters with flag kTPCout
469         // Note: updating with flag kITSout further down will overwrite the outer
470         // parameter again so this can be done only for ITS standalone tracks, meaning
471         // tracks in the ITS not associated with any TPC track
472         // HLT does not provide such standalone tracking
473         {
474           AliHLTGlobalBarrelTrack outPar(*element);       
475           outPar.AliExternalTrackParam::PropagateTo( element->GetLastPointX(), fSolenoidBz );
476           outPar.SetLabel(element->GetLabel());
477           //iotrack.UpdateTrackParams(&outPar,AliESDtrack::kTPCout);
478         }
479         // 2.2 TPC tracking estimates parameters at first cluster
480         iotrack.UpdateTrackParams(&(*element),AliESDtrack::kTPCin);
481
482         // 2.3 use the same parameters also for kTPCrefit, there is no proper refit in HLT
483         // maybe this can be done later, however also the inner parameter is changed which
484         // is used as a reference point in the display. That point should be at the first
485         // TPC cluster
486         iotrack.UpdateTrackParams(&(*element),AliESDtrack::kTPCrefit);
487         iotrack.SetTPCPoints(points);
488         if( element->TrackID()<ndEdxTPC ){
489           iotrack.SetTPCsignal( dEdxTPC[element->TrackID()], 0, 0 ); 
490           //AliTPCseed s;
491           //s.Set( element->GetX(), element->GetAlpha(),
492           //element->GetParameter(), element->GetCovariance() );
493           //s.SetdEdx( dEdxTPC[element->TrackID()] );
494           //s.CookPID();
495           //iotrack.SetTPCpid(s.TPCrPIDs() );
496         } else {
497           if( dEdxTPC ) HLTWarning("Wrong number of dEdx TPC labels");
498         }
499         iotrack.SetLabel(mcLabel);
500         pESD->AddTrack(&iotrack);
501         if (fVerbosity>0) element->Print();
502       }
503       HLTInfo("converted %d track(s) to AliESDtrack and added to ESD", tracks.size());
504       iAddedDataBlocks++;
505     } else if (iResult<0) {
506       HLTError("can not extract tracks from data block of type %s (specification %08x) of size %d: error %d", 
507                DataType2Text(pBlock->fDataType).c_str(), pBlock->fSpecification, pBlock->fSize, iResult);
508     }
509   }
510
511
512   // Get ITS SPD vertex
513   for( const AliHLTComponentBlockData *i= GetFirstInputBlock(kAliHLTDataTypeESDVertex|kAliHLTDataOriginITS); i!=NULL; i=GetNextInputBlock() ){
514     fBenchmark.AddInput(i->fSize);
515   }
516
517   for ( const TObject *iter = GetFirstInputObject(kAliHLTDataTypeESDVertex|kAliHLTDataOriginITS); iter != NULL; iter = GetNextInputObject() ) {
518     AliESDVertex *vtx = dynamic_cast<AliESDVertex*>(const_cast<TObject*>( iter ) );
519     pESD->SetPrimaryVertexSPD( vtx );
520   }
521
522   // 3.1. now update ESD tracks with the ITSOut info
523   // updating track parameters with flag kITSout will overwrite parameter set above with flag kTPCout
524   // TODO 2010-07-12 there are some issues with this updating sequence, for the moment update with
525   // flag kITSout is disabled, would require
526   // a) to update kTPCout after kITSout
527   // b) update only for ITS standalone tracks. The HLT ITS tracker does not provide standalone
528   //    tracking, every track is related to a TPC track
529   // Furthermore there seems to be a bug as the data blocks describe track parameters around the
530   // vertex, not at the outer ITS
531   // bug https://savannah.cern.ch/bugs/index.php?69872
532   for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeTrack|kAliHLTDataOriginITSOut);
533        pBlock!=NULL; pBlock=GetNextInputBlock()) {
534     fBenchmark.AddInput(pBlock->fSize);
535     vector<AliHLTGlobalBarrelTrack> tracks;
536     if ((iResult=AliHLTGlobalBarrelTrack::ConvertTrackDataArray(reinterpret_cast<const AliHLTTracksData*>(pBlock->fPtr), pBlock->fSize, tracks))>0) {
537       for (vector<AliHLTGlobalBarrelTrack>::iterator element=tracks.begin();
538            element!=tracks.end(); element++) {
539         int tpcID=element->TrackID();
540         // the ITS tracker assigns the TPC track used as seed for a certain track to
541         // the trackID
542         if( tpcID<0 || tpcID>=pESD->GetNumberOfTracks()) continue;
543         AliESDtrack *tESD = pESD->GetTrack( tpcID );
544         element->SetLabel(tESD->GetLabel());
545         // 2010-07-12 disabled, see above, bugfix https://savannah.cern.ch/bugs/index.php?69872
546         //if( tESD ) tESD->UpdateTrackParams( &(*element), AliESDtrack::kITSout );
547       }
548     }
549   }
550
551   // 3.2. now update ESD tracks with the ITS info
552   for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeTrack|kAliHLTDataOriginITS);
553        pBlock!=NULL; pBlock=GetNextInputBlock()) {
554     fBenchmark.AddInput(pBlock->fSize);
555     vector<AliHLTGlobalBarrelTrack> tracks;
556     if ((iResult=AliHLTGlobalBarrelTrack::ConvertTrackDataArray(reinterpret_cast<const AliHLTTracksData*>(pBlock->fPtr), pBlock->fSize, tracks))>0) {
557       for (vector<AliHLTGlobalBarrelTrack>::iterator element=tracks.begin();
558            element!=tracks.end(); element++) {
559         int tpcID=element->TrackID();
560         // the ITS tracker assigns the TPC track used as seed for a certain track to
561         // the trackID
562         if( tpcID<0 || tpcID>=pESD->GetNumberOfTracks()) continue;
563         Int_t mcLabel = -1;
564         if( mcLabelsITS.find(element->TrackID())!=mcLabelsITS.end() )
565           mcLabel = mcLabelsITS[element->TrackID()];
566         AliESDtrack *tESD = pESD->GetTrack( tpcID );
567         if (!tESD) continue;
568         // the labels for the TPC and ITS tracking params can be different, e.g.
569         // there can be a decay. The ITS label should then be the better one, the
570         // TPC label is saved in a member of AliESDtrack
571         if (mcLabel>=0) {
572           // upadte only if the ITS label is available, otherwise keep TPC label
573           element->SetLabel( mcLabel );
574         } else {
575           // bugfix https://savannah.cern.ch/bugs/?69713
576           element->SetLabel( tESD->GetLabel() );          
577         }
578         tESD->UpdateTrackParams( &(*element), AliESDtrack::kITSin );
579
580         // TODO: add a proper refit
581         //tESD->UpdateTrackParams( &(*element), AliESDtrack::kTPCrefit );
582       }
583     }
584   }
585
586   // update with  vertices and vertex-fitted tracks
587
588   for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeGlobalVertexer);
589        pBlock!=NULL; pBlock=GetNextInputBlock()) {
590     fBenchmark.AddInput(pBlock->fSize);   
591     AliHLTGlobalVertexerComponent::FillESD( pESD, reinterpret_cast<AliHLTGlobalVertexerComponent::AliHLTGlobalVertexerData* >(pBlock->fPtr) );
592   }
593
594   // Fill DCA parameters for TPC tracks
595   // TODO 2010-07-12 this propagates also the TPC inner param to beamline
596   // sounds not very reasonable
597   // https://savannah.cern.ch/bugs/index.php?69873
598   for (int i=0; i<pESD->GetNumberOfTracks(); i++) {
599     if (!pESD->GetTrack(i) || 
600         !pESD->GetTrack(i)->GetTPCInnerParam() ) continue;
601     pESD->GetTrack(i)->RelateToVertexTPC(pESD->GetPrimaryVertexTracks(), fSolenoidBz, 1000 );    
602   }
603
604   // loop over all tracks and set the TPC refit flag by updating with the
605   // original TPC inner parameter if not yet set
606   // TODO: replace this by a proper refit
607   // code is comented for the moment as it does not fully solve the problems with
608   // the display
609   // - would set the main parameters to the TPC inner wall again, or
610   // - changes the inner param if the parameters are propagated, so we loose the track
611   //   reference point for the display
612   // with the current sequence we have the latter case as the DCA operations above
613   // change the TPC inner parameters
614   /*
615   for (int i=0; i<pESD->GetNumberOfTracks(); i++) {
616     if (!pESD->GetTrack(i) || 
617         !pESD->GetTrack(i)->GetTPCInnerParam() ||
618         pESD->GetTrack(i)->IsOn(AliESDtrack::kTPCrefit)) continue;
619     AliESDtrack* tESD=pESD->GetTrack(i);
620     AliHLTGlobalBarrelTrack inner(*tESD->GetTPCInnerParam());
621     inner.SetLabel(tESD->GetLabel());
622     tESD->UpdateTrackParams(&inner, AliESDtrack::kTPCrefit);
623   }
624   */
625
626   // 4. convert the HLT TRD tracks to ESD tracks                        
627   for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeTrack | kAliHLTDataOriginTRD);
628        pBlock!=NULL; pBlock=GetNextInputBlock()) {
629     fBenchmark.AddInput(pBlock->fSize);
630     vector<AliHLTGlobalBarrelTrack> tracks;
631     if ((iResult=AliHLTGlobalBarrelTrack::ConvertTrackDataArray(reinterpret_cast<const AliHLTTracksData*>(pBlock->fPtr), pBlock->fSize, tracks))>0) {
632       for (vector<AliHLTGlobalBarrelTrack>::iterator element=tracks.begin();
633            element!=tracks.end(); element++) {
634         
635         Double_t TRDpid[AliPID::kSPECIES], eProb(0.2), restProb((1-eProb)/(AliPID::kSPECIES-1)); //eprob(element->GetTRDpid...);
636         for(Int_t i=0; i<AliPID::kSPECIES; i++){
637           switch(i){
638           case AliPID::kElectron: TRDpid[AliPID::kElectron]=eProb; break;
639           default: TRDpid[i]=restProb; break;
640           }
641         }
642         
643         AliESDtrack iotrack;
644         iotrack.UpdateTrackParams(&(*element),AliESDtrack::kTRDout);
645         iotrack.SetTRDpid(TRDpid);
646         
647         pESD->AddTrack(&iotrack);
648         if (fVerbosity>0) element->Print();
649       }
650       HLTInfo("converted %d track(s) to AliESDtrack and added to ESD", tracks.size());
651       iAddedDataBlocks++;
652     } else if (iResult<0) {
653       HLTError("can not extract tracks from data block of type %s (specification %08x) of size %d: error %d", 
654                DataType2Text(pBlock->fDataType).c_str(), pBlock->fSpecification, pBlock->fSize, iResult);
655     }
656   }
657   for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeCaloCluster | kAliHLTDataOriginAny); pBlock!=NULL; pBlock=GetNextInputBlock()) 
658     {
659       fBenchmark.AddInput(pBlock->fSize);
660       AliHLTCaloClusterHeaderStruct *caloClusterHeaderPtr = reinterpret_cast<AliHLTCaloClusterHeaderStruct*>(pBlock->fPtr);
661
662       HLTDebug("%d HLT clusters from spec: 0x%X", caloClusterHeaderPtr->fNClusters, pBlock->fSpecification);
663
664       //AliHLTCaloClusterReader reader;
665       //reader.SetMemory(caloClusterHeaderPtr);
666
667       AliHLTESDCaloClusterMaker clusterMaker;
668
669       int nClusters = clusterMaker.FillESD(pESD, caloClusterHeaderPtr);
670      
671       HLTInfo("converted %d cluster(s) to AliESDCaloCluster and added to ESD", nClusters);
672       iAddedDataBlocks++;
673     }
674   
675   // Add tracks from MUON.
676   for( const AliHLTComponentBlockData *i= GetFirstInputBlock(kAliHLTAnyDataType | kAliHLTDataOriginMUON); i!=NULL; i=GetNextInputBlock() ){
677     fBenchmark.AddInput(i->fSize);
678   }
679
680   for (const TObject* obj = GetFirstInputObject(kAliHLTAnyDataType | kAliHLTDataOriginMUON);
681        obj != NULL;
682        obj = GetNextInputObject()
683       )
684   {
685     const TClonesArray* tracklist = NULL;
686     if (obj->IsA() == AliESDEvent::Class())
687     {
688       const AliESDEvent* event = static_cast<const AliESDEvent*>(obj);
689       HLTDebug("Received a MUON ESD with specification: 0x%X", GetSpecification(obj));
690       if (event->GetList() == NULL) continue;
691       tracklist = dynamic_cast<const TClonesArray*>(event->GetList()->FindObject("MuonTracks"));
692       if (tracklist == NULL) continue;
693     }
694     else if (obj->IsA() == TClonesArray::Class())
695     {
696       tracklist = static_cast<const TClonesArray*>(obj);
697       HLTDebug("Received a MUON TClonesArray of tracks with specification: 0x%X", GetSpecification(obj));
698     }
699     else
700     {
701       // Cannot handle this object type.
702       continue;
703     }
704     HLTDebug("Received %d MUON tracks.", tracklist->GetEntriesFast());
705     if (tracklist->GetEntriesFast() > 0)
706     {
707       const AliESDMuonTrack* track = dynamic_cast<const AliESDMuonTrack*>(tracklist->UncheckedAt(0));
708       if (track == NULL)
709       {
710         HLTError(Form("%s from MUON does not contain AliESDMuonTrack objects.", obj->ClassName()));
711         continue;
712       }
713     }
714     for (Int_t i = 0; i < tracklist->GetEntriesFast(); ++i)
715     {
716       const AliESDMuonTrack* track = static_cast<const AliESDMuonTrack*>(tracklist->UncheckedAt(i));
717       pESD->AddMuonTrack(track);
718     }
719   }
720   
721   if (iAddedDataBlocks>0 && pTree) {
722     pTree->Fill();
723   }
724   
725   if (iResult>=0) iResult=iAddedDataBlocks;
726   return iResult;
727 }