]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/global/AliHLTGlobalEsdConverterComponent.cxx
Cosmetics
[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
358   // in the first attempt this component reads the TPC tracks and updates in the
359   // second step from the ITS tracks
360
361   // first read MC information (if present)
362   std::map<int,int> mcLabels;
363   std::map<int,int> mcLabelsITS;
364
365   for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeTrackMC|kAliHLTDataOriginTPC);
366        pBlock!=NULL; pBlock=GetNextInputBlock()) {
367
368     fBenchmark.AddInput(pBlock->fSize);
369
370     AliHLTTrackMCData* dataPtr = reinterpret_cast<AliHLTTrackMCData*>( pBlock->fPtr );
371     if (sizeof(AliHLTTrackMCData)+dataPtr->fCount*sizeof(AliHLTTrackMCLabel)==pBlock->fSize) {
372       for( unsigned int il=0; il<dataPtr->fCount; il++ ){
373         AliHLTTrackMCLabel &lab = dataPtr->fLabels[il];
374         mcLabels[lab.fTrackID] = lab.fMCLabel;
375       }
376     } else {
377       HLTWarning("data mismatch in block %s (0x%08x): count %d, size %d -> ignoring track MC information", 
378                  DataType2Text(pBlock->fDataType).c_str(), pBlock->fSpecification, 
379                  dataPtr->fCount, pBlock->fSize);
380     }
381   }
382  
383   for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeTrackMC|kAliHLTDataOriginITS);
384        pBlock!=NULL; pBlock=GetNextInputBlock()) {
385
386     fBenchmark.AddInput(pBlock->fSize);
387
388     AliHLTTrackMCData* dataPtr = reinterpret_cast<AliHLTTrackMCData*>( pBlock->fPtr );
389     if (sizeof(AliHLTTrackMCData)+dataPtr->fCount*sizeof(AliHLTTrackMCLabel)==pBlock->fSize) {
390       for( unsigned int il=0; il<dataPtr->fCount; il++ ){
391         AliHLTTrackMCLabel &lab = dataPtr->fLabels[il];
392         mcLabelsITS[lab.fTrackID] = lab.fMCLabel;
393       }
394     } else {
395       HLTWarning("data mismatch in block %s (0x%08x): count %d, size %d -> ignoring track MC information", 
396                  DataType2Text(pBlock->fDataType).c_str(), pBlock->fSpecification, 
397                  dataPtr->fCount, pBlock->fSize);
398     }
399   }
400
401
402   // read dEdx information (if present)
403
404   AliHLTFloat32_t *dEdxTPC = 0; 
405   Int_t ndEdxTPC = 0;
406   for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypedEdx|kAliHLTDataOriginTPC);
407        pBlock!=NULL; pBlock=GetNextInputBlock()) {
408     fBenchmark.AddInput(pBlock->fSize);
409     dEdxTPC = reinterpret_cast<AliHLTFloat32_t*>( pBlock->fPtr );
410     ndEdxTPC = pBlock->fSize / sizeof(AliHLTFloat32_t);
411     break;
412   }
413
414   // convert the TPC tracks to ESD tracks
415   for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeTrack|kAliHLTDataOriginTPC);
416        pBlock!=NULL; pBlock=GetNextInputBlock()) {
417     fBenchmark.AddInput(pBlock->fSize);
418     vector<AliHLTGlobalBarrelTrack> tracks;
419     if ((iResult=AliHLTGlobalBarrelTrack::ConvertTrackDataArray(reinterpret_cast<const AliHLTTracksData*>(pBlock->fPtr), pBlock->fSize, tracks))>=0) {
420       for (vector<AliHLTGlobalBarrelTrack>::iterator element=tracks.begin();
421            element!=tracks.end(); element++) {
422         Float_t points[4] = {
423           element->GetX(),
424           element->GetY(),
425           element->GetLastPointX(),
426           element->GetLastPointY()
427         };
428
429         Int_t mcLabel = -1;
430         if( mcLabels.find(element->TrackID())!=mcLabels.end() )
431           mcLabel = mcLabels[element->TrackID()];
432         element->SetLabel( mcLabel );
433
434         AliESDtrack iotrack;
435
436         // for the moment, the number of clusters are not set when processing the
437         // kTPCin update, only at kTPCout
438         // there ar emainly three parameters updated for kTPCout
439         //   number of clusters
440         //   chi2
441         //   pid signal
442         // The first one can be updated already at that stage here, while the two others
443         // eventually require to update from the ITS tracks before. The exact scheme
444         // needs to be checked
445         // 2010-05-12 TODO: the outer parameter is set when updating with kTPCout but
446         // also when updated with kITSout. So the value from here is overwritten further
447         // down. Comes along with the necessity to check the full sequence.
448         iotrack.SetID( element->TrackID() );
449         {
450           AliHLTGlobalBarrelTrack outPar(*element);       
451           outPar.AliExternalTrackParam::PropagateTo( element->GetLastPointX(), fSolenoidBz );
452           outPar.SetLabel(element->GetLabel());
453           iotrack.UpdateTrackParams(&outPar,AliESDtrack::kTPCout);
454         }
455         iotrack.UpdateTrackParams(&(*element),AliESDtrack::kTPCin);
456         iotrack.SetTPCPoints(points);
457         if( element->TrackID()<ndEdxTPC ){
458           iotrack.SetTPCsignal( dEdxTPC[element->TrackID()], 0, 0 ); 
459           //AliTPCseed s;
460           //s.Set( element->GetX(), element->GetAlpha(),
461           //element->GetParameter(), element->GetCovariance() );
462           //s.SetdEdx( dEdxTPC[element->TrackID()] );
463           //s.CookPID();
464           //iotrack.SetTPCpid(s.TPCrPIDs() );
465         } else {
466           if( dEdxTPC ) HLTWarning("Wrong number of dEdx TPC labels");
467         }
468         iotrack.SetLabel(mcLabel);
469         pESD->AddTrack(&iotrack);
470         if (fVerbosity>0) element->Print();
471       }
472       HLTInfo("converted %d track(s) to AliESDtrack and added to ESD", tracks.size());
473       iAddedDataBlocks++;
474     } else if (iResult<0) {
475       HLTError("can not extract tracks from data block of type %s (specification %08x) of size %d: error %d", 
476                DataType2Text(pBlock->fDataType).c_str(), pBlock->fSpecification, pBlock->fSize, iResult);
477     }
478   }
479
480
481   // Get ITS SPD vertex
482   for( const AliHLTComponentBlockData *i= GetFirstInputBlock(kAliHLTDataTypeESDVertex|kAliHLTDataOriginITS); i!=NULL; i=GetNextInputBlock() ){
483     fBenchmark.AddInput(i->fSize);
484   }
485
486   for ( const TObject *iter = GetFirstInputObject(kAliHLTDataTypeESDVertex|kAliHLTDataOriginITS); iter != NULL; iter = GetNextInputObject() ) {
487     AliESDVertex *vtx = dynamic_cast<AliESDVertex*>(const_cast<TObject*>( iter ) );
488     pESD->SetPrimaryVertexSPD( vtx );
489   }
490
491   // now update ESD tracks with the ITSOut info
492   for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeTrack|kAliHLTDataOriginITSOut);
493        pBlock!=NULL; pBlock=GetNextInputBlock()) {
494     fBenchmark.AddInput(pBlock->fSize);
495     vector<AliHLTGlobalBarrelTrack> tracks;
496     if ((iResult=AliHLTGlobalBarrelTrack::ConvertTrackDataArray(reinterpret_cast<const AliHLTTracksData*>(pBlock->fPtr), pBlock->fSize, tracks))>0) {
497       for (vector<AliHLTGlobalBarrelTrack>::iterator element=tracks.begin();
498            element!=tracks.end(); element++) {
499         int tpcID=element->TrackID();
500         // the ITS tracker assigns the TPC track used as seed for a certain track to
501         // the trackID
502         if( tpcID<0 || tpcID>=pESD->GetNumberOfTracks()) continue;
503         AliESDtrack *tESD = pESD->GetTrack( tpcID );
504         element->SetLabel(tESD->GetLabel());
505         if( tESD ) tESD->UpdateTrackParams( &(*element), AliESDtrack::kITSout );
506       }
507     }
508   }
509
510   // now update ESD tracks with the ITS info
511   for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeTrack|kAliHLTDataOriginITS);
512        pBlock!=NULL; pBlock=GetNextInputBlock()) {
513     fBenchmark.AddInput(pBlock->fSize);
514     vector<AliHLTGlobalBarrelTrack> tracks;
515     if ((iResult=AliHLTGlobalBarrelTrack::ConvertTrackDataArray(reinterpret_cast<const AliHLTTracksData*>(pBlock->fPtr), pBlock->fSize, tracks))>0) {
516       for (vector<AliHLTGlobalBarrelTrack>::iterator element=tracks.begin();
517            element!=tracks.end(); element++) {
518         int tpcID=element->TrackID();
519         // the ITS tracker assigns the TPC track used as seed for a certain track to
520         // the trackID
521         if( tpcID<0 || tpcID>=pESD->GetNumberOfTracks()) continue;
522         Int_t mcLabel = -1;
523         if( mcLabelsITS.find(element->TrackID())!=mcLabelsITS.end() )
524           mcLabel = mcLabelsITS[element->TrackID()];
525         AliESDtrack *tESD = pESD->GetTrack( tpcID );
526         if (!tESD) continue;
527         if (mcLabel>=0 && tESD->GetLabel()>=0 && mcLabel!=tESD->GetLabel()) {
528           // label should be equal at this point
529           HLTWarning("mismatch in TPC and ITS MC label: %d vs. %d, ignoring ITS label", tESD->GetLabel(), mcLabel);
530           mcLabel=tESD->GetLabel();
531         } else if (mcLabel<0 && tESD->GetLabel()>=0) {
532           // keep the TPC label
533           mcLabel=tESD->GetLabel();
534         }
535         element->SetLabel( mcLabel );
536         tESD->UpdateTrackParams( &(*element), AliESDtrack::kITSin );
537
538         // TODO: add a proper refit
539         //tESD->UpdateTrackParams( &(*element), AliESDtrack::kTPCrefit );
540       }
541     }
542   }
543
544   // update with  vertices and vertex-fitted tracks
545
546   for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeGlobalVertexer);
547        pBlock!=NULL; pBlock=GetNextInputBlock()) {
548     fBenchmark.AddInput(pBlock->fSize);   
549     AliHLTGlobalVertexerComponent::FillESD( pESD, reinterpret_cast<AliHLTGlobalVertexerComponent::AliHLTGlobalVertexerData* >(pBlock->fPtr) );
550   }
551
552   // loop over all tracks and set the TPC refit flag by updating with the
553   // original TPC inner parameter if not yet set
554   // TODO: replace this by a proper refit
555   // code is comented for the moment as it does not fully solve the problems with
556   // the display
557   /*
558   for (int i=0; i<pESD->GetNumberOfTracks(); i++) {
559     if (!pESD->GetTrack(i) || 
560         !pESD->GetTrack(i)->GetTPCInnerParam() ||
561         pESD->GetTrack(i)->IsOn(AliESDtrack::kTPCrefit)) continue;
562     AliESDtrack* tESD=pESD->GetTrack(i);
563     AliHLTGlobalBarrelTrack inner(*tESD->GetTPCInnerParam());
564     inner.SetLabel(tESD->GetLabel());
565     tESD->UpdateTrackParams(&inner, AliESDtrack::kTPCrefit);
566   }
567   */
568
569   // convert the HLT TRD tracks to ESD tracks                        
570   for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeTrack | kAliHLTDataOriginTRD);
571        pBlock!=NULL; pBlock=GetNextInputBlock()) {
572     fBenchmark.AddInput(pBlock->fSize);
573     vector<AliHLTGlobalBarrelTrack> tracks;
574     if ((iResult=AliHLTGlobalBarrelTrack::ConvertTrackDataArray(reinterpret_cast<const AliHLTTracksData*>(pBlock->fPtr), pBlock->fSize, tracks))>0) {
575       for (vector<AliHLTGlobalBarrelTrack>::iterator element=tracks.begin();
576            element!=tracks.end(); element++) {
577         
578         Double_t TRDpid[AliPID::kSPECIES], eProb(0.2), restProb((1-eProb)/(AliPID::kSPECIES-1)); //eprob(element->GetTRDpid...);
579         for(Int_t i=0; i<AliPID::kSPECIES; i++){
580           switch(i){
581           case AliPID::kElectron: TRDpid[AliPID::kElectron]=eProb; break;
582           default: TRDpid[i]=restProb; break;
583           }
584         }
585         
586         AliESDtrack iotrack;
587         iotrack.UpdateTrackParams(&(*element),AliESDtrack::kTRDout);
588         iotrack.SetTRDpid(TRDpid);
589         
590         pESD->AddTrack(&iotrack);
591         if (fVerbosity>0) element->Print();
592       }
593       HLTInfo("converted %d track(s) to AliESDtrack and added to ESD", tracks.size());
594       iAddedDataBlocks++;
595     } else if (iResult<0) {
596       HLTError("can not extract tracks from data block of type %s (specification %08x) of size %d: error %d", 
597                DataType2Text(pBlock->fDataType).c_str(), pBlock->fSpecification, pBlock->fSize, iResult);
598     }
599   }
600   for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeCaloCluster | kAliHLTDataOriginAny); pBlock!=NULL; pBlock=GetNextInputBlock()) 
601     {
602       fBenchmark.AddInput(pBlock->fSize);
603       AliHLTCaloClusterHeaderStruct *caloClusterHeaderPtr = reinterpret_cast<AliHLTCaloClusterHeaderStruct*>(pBlock->fPtr);
604
605       HLTDebug("%d HLT clusters from spec: 0x%X", caloClusterHeaderPtr->fNClusters, pBlock->fSpecification);
606
607       //AliHLTCaloClusterReader reader;
608       //reader.SetMemory(caloClusterHeaderPtr);
609
610       AliHLTESDCaloClusterMaker clusterMaker;
611
612       int nClusters = clusterMaker.FillESD(pESD, caloClusterHeaderPtr);
613      
614       HLTInfo("converted %d cluster(s) to AliESDCaloCluster and added to ESD", nClusters);
615       iAddedDataBlocks++;
616     }
617   
618   // Add tracks from MUON.
619   for( const AliHLTComponentBlockData *i= GetFirstInputBlock(kAliHLTAnyDataType | kAliHLTDataOriginMUON); i!=NULL; i=GetNextInputBlock() ){
620     fBenchmark.AddInput(i->fSize);
621   }
622
623   for (const TObject* obj = GetFirstInputObject(kAliHLTAnyDataType | kAliHLTDataOriginMUON);
624        obj != NULL;
625        obj = GetNextInputObject()
626       )
627   {
628     const TClonesArray* tracklist = NULL;
629     if (obj->IsA() == AliESDEvent::Class())
630     {
631       const AliESDEvent* event = static_cast<const AliESDEvent*>(obj);
632       HLTDebug("Received a MUON ESD with specification: 0x%X", GetSpecification(obj));
633       if (event->GetList() == NULL) continue;
634       tracklist = dynamic_cast<const TClonesArray*>(event->GetList()->FindObject("MuonTracks"));
635       if (tracklist == NULL) continue;
636     }
637     else if (obj->IsA() == TClonesArray::Class())
638     {
639       tracklist = static_cast<const TClonesArray*>(obj);
640       HLTDebug("Received a MUON TClonesArray of tracks with specification: 0x%X", GetSpecification(obj));
641     }
642     else
643     {
644       // Cannot handle this object type.
645       continue;
646     }
647     HLTDebug("Received %d MUON tracks.", tracklist->GetEntriesFast());
648     if (tracklist->GetEntriesFast() > 0)
649     {
650       const AliESDMuonTrack* track = dynamic_cast<const AliESDMuonTrack*>(tracklist->UncheckedAt(0));
651       if (track == NULL)
652       {
653         HLTError(Form("%s from MUON does not contain AliESDMuonTrack objects.", obj->ClassName()));
654         continue;
655       }
656     }
657     for (Int_t i = 0; i < tracklist->GetEntriesFast(); ++i)
658     {
659       const AliESDMuonTrack* track = static_cast<const AliESDMuonTrack*>(tracklist->UncheckedAt(i));
660       pESD->AddMuonTrack(track);
661     }
662   }
663   
664   if (iAddedDataBlocks>0 && pTree) {
665     pTree->Fill();
666   }
667   
668   if (iResult>=0) iResult=iAddedDataBlocks;
669   return iResult;
670 }