]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/global/AliHLTGlobalEsdConverterComponent.cxx
- fixing binning
[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             if (fESD->GetList()->FindObject(id)) {
253               HLTDebug("removing object %s", id);
254               fESD->GetList()->Remove(fESD->GetList()->FindObject(id));
255             } else {
256               HLTWarning("failed to remove object '%s' from ESD", id);
257             }
258           }
259           // id="AliESDRun";           if (fESD->GetList()->FindObject(id)) fESD->GetList()->Remove(fESD->GetList()->FindObject(id));
260           // id="AliESDHeader";        if (fESD->GetList()->FindObject(id)) fESD->GetList()->Remove(fESD->GetList()->FindObject(id));
261           // id="AliESDZDC";           if (fESD->GetList()->FindObject(id)) fESD->GetList()->Remove(fESD->GetList()->FindObject(id));
262           // id="AliESDFMD";           if (fESD->GetList()->FindObject(id)) fESD->GetList()->Remove(fESD->GetList()->FindObject(id));
263           // id="AliESDVZERO";         if (fESD->GetList()->FindObject(id)) fESD->GetList()->Remove(fESD->GetList()->FindObject(id));
264           // id="AliESDTZERO";         if (fESD->GetList()->FindObject(id)) fESD->GetList()->Remove(fESD->GetList()->FindObject(id));
265           // id="TPCVertex";           if (fESD->GetList()->FindObject(id)) fESD->GetList()->Remove(fESD->GetList()->FindObject(id));
266           // id="SPDVertex";           if (fESD->GetList()->FindObject(id)) fESD->GetList()->Remove(fESD->GetList()->FindObject(id));
267           // id="PrimaryVertex";       if (fESD->GetList()->FindObject(id)) fESD->GetList()->Remove(fESD->GetList()->FindObject(id));
268           // id="AliMultiplicity";     if (fESD->GetList()->FindObject(id)) fESD->GetList()->Remove(fESD->GetList()->FindObject(id));
269           // id="PHOSTrigger";         if (fESD->GetList()->FindObject(id)) fESD->GetList()->Remove(fESD->GetList()->FindObject(id));
270           // id="EMCALTrigger";        if (fESD->GetList()->FindObject(id)) fESD->GetList()->Remove(fESD->GetList()->FindObject(id));
271           // id="SPDPileupVertices";   if (fESD->GetList()->FindObject(id)) fESD->GetList()->Remove(fESD->GetList()->FindObject(id));
272           // id="TrkPileupVertices";   if (fESD->GetList()->FindObject(id)) fESD->GetList()->Remove(fESD->GetList()->FindObject(id));
273           // id="Cascades";            if (fESD->GetList()->FindObject(id)) fESD->GetList()->Remove(fESD->GetList()->FindObject(id));
274           // id="Kinks";               if (fESD->GetList()->FindObject(id)) fESD->GetList()->Remove(fESD->GetList()->FindObject(id));
275           // id="AliRawDataErrorLogs"; if (fESD->GetList()->FindObject(id)) fESD->GetList()->Remove(fESD->GetList()->FindObject(id));
276           // id="AliESDACORDE";        if (fESD->GetList()->FindObject(id)) fESD->GetList()->Remove(fESD->GetList()->FindObject(id));
277           fESD->GetStdContent();
278           delete pTokens;
279         }
280       }
281     } else {
282       iResult=-ENOMEM;
283     }
284
285     SetupCTPData();
286   }
287
288   fBenchmark.SetTimer(0,"total");
289
290   return iResult;
291 }
292
293 int AliHLTGlobalEsdConverterComponent::DoDeinit()
294 {
295   // see header file for class documentation
296   if (fESD) delete fESD;
297   fESD=NULL;
298
299   return 0;
300 }
301
302 int AliHLTGlobalEsdConverterComponent::DoEvent(const AliHLTComponentEventData& /*evtData*/, 
303                                                AliHLTComponentTriggerData& trigData)
304 {
305   // see header file for class documentation
306   int iResult=0;
307   if (!fESD) return -ENODEV;
308
309   if (IsDataEvent()) fBenchmark.StartNewEvent();
310   fBenchmark.Start(0);
311
312   AliESDEvent* pESD = fESD;
313
314   pESD->Reset(); 
315   pESD->SetMagneticField(fSolenoidBz);
316   pESD->SetRunNumber(GetRunNo());
317   pESD->SetPeriodNumber(GetPeriodNumber());
318   pESD->SetOrbitNumber(GetOrbitNumber());
319   pESD->SetBunchCrossNumber(GetBunchCrossNumber());
320   pESD->SetTimeStamp(GetTimeStamp());
321
322   const AliHLTCTPData* pCTPData=CTPData();
323   if (pCTPData) {
324     AliHLTUInt64_t mask=pCTPData->ActiveTriggers(trigData);
325     for (int index=0; index<gkNCTPTriggerClasses; index++) {
326       if ((mask&((AliHLTUInt64_t)0x1<<index)) == 0) continue;
327       pESD->SetTriggerClass(pCTPData->Name(index), index);
328     }
329     pESD->SetTriggerMask(mask);
330   }
331
332   TTree* pTree = NULL;
333   if (fWriteTree)
334     pTree = new TTree("esdTree", "Tree with HLT ESD objects");
335  
336   if (pTree) {
337     pTree->SetDirectory(0);
338   }
339
340   if ((iResult=ProcessBlocks(pTree, pESD))>=0) {
341     // TODO: set the specification correctly
342     if (pTree) {
343       // the esd structure is written to the user info and is
344       // needed in te ReadFromTree method to read all objects correctly
345       pTree->GetUserInfo()->Add(pESD);
346       pESD->WriteToTree(pTree);
347       iResult=PushBack(pTree, kAliHLTDataTypeESDTree|kAliHLTDataOriginOut, 0);
348     } else {
349       iResult=PushBack(pESD, kAliHLTDataTypeESDObject|kAliHLTDataOriginOut, 0);
350     }
351     fBenchmark.AddOutput(GetLastObjectSize());
352   }
353   if (pTree) {
354     // clear user info list to prevent objects from being deleted
355     pTree->GetUserInfo()->Clear();
356     delete pTree;
357   }
358
359   fBenchmark.Stop(0);
360   HLTInfo( fBenchmark.GetStatistics() );
361
362   return iResult;
363 }
364
365 int AliHLTGlobalEsdConverterComponent::ProcessBlocks(TTree* pTree, AliESDEvent* pESD)
366 {
367   // see header file for class documentation
368
369   int iResult=0;
370   int iAddedDataBlocks=0;
371
372   // Barrel tracking
373
374   // in the first attempt this component reads the TPC tracks and updates in the
375   // second step from the ITS tracks
376
377   // first read MC information (if present)
378   std::map<int,int> mcLabels;
379
380   for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeTrackMC|kAliHLTDataOriginTPC);
381        pBlock!=NULL; pBlock=GetNextInputBlock()) {
382
383     fBenchmark.AddInput(pBlock->fSize);
384
385     AliHLTTrackMCData* dataPtr = reinterpret_cast<AliHLTTrackMCData*>( pBlock->fPtr );
386     if (sizeof(AliHLTTrackMCData)+dataPtr->fCount*sizeof(AliHLTTrackMCLabel)==pBlock->fSize) {
387       for( unsigned int il=0; il<dataPtr->fCount; il++ ){
388         AliHLTTrackMCLabel &lab = dataPtr->fLabels[il];
389         mcLabels[lab.fTrackID] = lab.fMCLabel;
390       }
391     } else {
392       HLTWarning("data mismatch in block %s (0x%08x): count %d, size %d -> ignoring track MC information", 
393                  DataType2Text(pBlock->fDataType).c_str(), pBlock->fSpecification, 
394                  dataPtr->fCount, pBlock->fSize);
395     }
396   }
397
398   // read dEdx information (if present)
399
400   AliHLTFloat32_t *dEdxTPC = 0; 
401   Int_t ndEdxTPC = 0;
402   for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypedEdx|kAliHLTDataOriginTPC);
403        pBlock!=NULL; pBlock=GetNextInputBlock()) {
404     fBenchmark.AddInput(pBlock->fSize);
405     dEdxTPC = reinterpret_cast<AliHLTFloat32_t*>( pBlock->fPtr );
406     ndEdxTPC = pBlock->fSize / sizeof(AliHLTFloat32_t);
407     break;
408   }
409
410   // convert the TPC tracks to ESD tracks
411   for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeTrack|kAliHLTDataOriginTPC);
412        pBlock!=NULL; pBlock=GetNextInputBlock()) {
413     fBenchmark.AddInput(pBlock->fSize);
414     vector<AliHLTGlobalBarrelTrack> tracks;
415     if ((iResult=AliHLTGlobalBarrelTrack::ConvertTrackDataArray(reinterpret_cast<const AliHLTTracksData*>(pBlock->fPtr), pBlock->fSize, tracks))>=0) {
416       for (vector<AliHLTGlobalBarrelTrack>::iterator element=tracks.begin();
417            element!=tracks.end(); element++) {
418         Float_t points[4] = {
419           element->GetX(),
420           element->GetY(),
421           element->GetLastPointX(),
422           element->GetLastPointY()
423         };
424
425         Int_t mcLabel = -1;
426         if( mcLabels.find(element->TrackID())!=mcLabels.end() )
427           mcLabel = mcLabels[element->TrackID()];
428         element->SetLabel( mcLabel );
429
430         AliESDtrack iotrack;
431
432         // for the moment, the number of clusters are not set when processing the
433         // kTPCin update, only at kTPCout
434         // there ar emainly three parameters updated for kTPCout
435         //   number of clusters
436         //   chi2
437         //   pid signal
438         // The first one can be updated already at that stage here, while the two others
439         // eventually require to update from the ITS tracks before. The exact scheme
440         // needs to be checked 
441         iotrack.SetID( element->TrackID() );
442         iotrack.UpdateTrackParams(&(*element),AliESDtrack::kTPCin);
443         {
444           AliHLTGlobalBarrelTrack outPar(*element);       
445           outPar.AliExternalTrackParam::PropagateTo( element->GetLastPointX(), fSolenoidBz );
446           iotrack.UpdateTrackParams(&outPar,AliESDtrack::kTPCout);
447         }
448         iotrack.SetTPCPoints(points);
449         if( element->TrackID()<ndEdxTPC ){
450           iotrack.SetTPCsignal( dEdxTPC[element->TrackID()], 0, 0 ); 
451           //AliTPCseed s;
452           //s.Set( element->GetX(), element->GetAlpha(),
453           //element->GetParameter(), element->GetCovariance() );
454           //s.SetdEdx( dEdxTPC[element->TrackID()] );
455           //s.CookPID();
456           //iotrack.SetTPCpid(s.TPCrPIDs() );
457         } else {
458           if( dEdxTPC ) HLTWarning("Wrong number of dEdx TPC labels");
459         }
460         pESD->AddTrack(&iotrack);
461         if (fVerbosity>0) element->Print();
462       }
463       HLTInfo("converted %d track(s) to AliESDtrack and added to ESD", tracks.size());
464       iAddedDataBlocks++;
465     } else if (iResult<0) {
466       HLTError("can not extract tracks from data block of type %s (specification %08x) of size %d: error %d", 
467                DataType2Text(pBlock->fDataType).c_str(), pBlock->fSpecification, pBlock->fSize, iResult);
468     }
469   }
470
471
472   // Get ITS SPD vertex
473   for( const AliHLTComponentBlockData *i= GetFirstInputBlock(kAliHLTDataTypeESDVertex|kAliHLTDataOriginITS); i!=NULL; i=GetNextInputBlock() ){
474     fBenchmark.AddInput(i->fSize);
475   }
476
477   for ( const TObject *iter = GetFirstInputObject(kAliHLTDataTypeESDVertex|kAliHLTDataOriginITS); iter != NULL; iter = GetNextInputObject() ) {
478     AliESDVertex *vtx = dynamic_cast<AliESDVertex*>(const_cast<TObject*>( iter ) );
479     pESD->SetPrimaryVertexSPD( vtx );
480   }
481
482   // now update ESD tracks with the ITSOut info
483   for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeTrack|kAliHLTDataOriginITSOut);
484        pBlock!=NULL; pBlock=GetNextInputBlock()) {
485     fBenchmark.AddInput(pBlock->fSize);
486     vector<AliHLTGlobalBarrelTrack> tracks;
487     if ((iResult=AliHLTGlobalBarrelTrack::ConvertTrackDataArray(reinterpret_cast<const AliHLTTracksData*>(pBlock->fPtr), pBlock->fSize, tracks))>0) {
488       for (vector<AliHLTGlobalBarrelTrack>::iterator element=tracks.begin();
489            element!=tracks.end(); element++) {
490         int tpcID=element->TrackID();
491         // the ITS tracker assigns the TPC track used as seed for a certain track to
492         // the trackID
493         if( tpcID<0 || tpcID>=pESD->GetNumberOfTracks()) continue;
494         AliESDtrack *tESD = pESD->GetTrack( tpcID );
495         if( tESD ) tESD->UpdateTrackParams( &(*element), AliESDtrack::kITSout );
496       }
497     }
498   }
499
500   // now update ESD tracks with the ITS info
501   for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeTrack|kAliHLTDataOriginITS);
502        pBlock!=NULL; pBlock=GetNextInputBlock()) {
503     fBenchmark.AddInput(pBlock->fSize);
504     vector<AliHLTGlobalBarrelTrack> tracks;
505     if ((iResult=AliHLTGlobalBarrelTrack::ConvertTrackDataArray(reinterpret_cast<const AliHLTTracksData*>(pBlock->fPtr), pBlock->fSize, tracks))>0) {
506       for (vector<AliHLTGlobalBarrelTrack>::iterator element=tracks.begin();
507            element!=tracks.end(); element++) {
508         int tpcID=element->TrackID();
509         // the ITS tracker assigns the TPC track used as seed for a certain track to
510         // the trackID
511         if( tpcID<0 || tpcID>=pESD->GetNumberOfTracks()) continue;
512         AliESDtrack *tESD = pESD->GetTrack( tpcID );
513         if( tESD ) tESD->UpdateTrackParams( &(*element), AliESDtrack::kITSin );
514       }
515     }
516   }
517
518   // update with  vertices and vertex-fitted tracks
519
520   for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeGlobalVertexer);
521        pBlock!=NULL; pBlock=GetNextInputBlock()) {
522     fBenchmark.AddInput(pBlock->fSize);   
523     AliHLTGlobalVertexerComponent::FillESD( pESD, reinterpret_cast<AliHLTGlobalVertexerComponent::AliHLTGlobalVertexerData* >(pBlock->fPtr) );
524   }
525
526
527   // convert the HLT TRD tracks to ESD tracks                        
528   for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeTrack | kAliHLTDataOriginTRD);
529        pBlock!=NULL; pBlock=GetNextInputBlock()) {
530     fBenchmark.AddInput(pBlock->fSize);
531     vector<AliHLTGlobalBarrelTrack> tracks;
532     if ((iResult=AliHLTGlobalBarrelTrack::ConvertTrackDataArray(reinterpret_cast<const AliHLTTracksData*>(pBlock->fPtr), pBlock->fSize, tracks))>0) {
533       for (vector<AliHLTGlobalBarrelTrack>::iterator element=tracks.begin();
534            element!=tracks.end(); element++) {
535         
536         Double_t TRDpid[AliPID::kSPECIES], eProb(0.2), restProb((1-eProb)/(AliPID::kSPECIES-1)); //eprob(element->GetTRDpid...);
537         for(Int_t i=0; i<AliPID::kSPECIES; i++){
538           switch(i){
539           case AliPID::kElectron: TRDpid[AliPID::kElectron]=eProb; break;
540           default: TRDpid[i]=restProb; break;
541           }
542         }
543         
544         AliESDtrack iotrack;
545         iotrack.UpdateTrackParams(&(*element),AliESDtrack::kTRDin);
546         iotrack.SetTRDpid(TRDpid);
547         
548         pESD->AddTrack(&iotrack);
549         if (fVerbosity>0) element->Print();
550       }
551       HLTInfo("converted %d track(s) to AliESDtrack and added to ESD", tracks.size());
552       iAddedDataBlocks++;
553     } else if (iResult<0) {
554       HLTError("can not extract tracks from data block of type %s (specification %08x) of size %d: error %d", 
555                DataType2Text(pBlock->fDataType).c_str(), pBlock->fSpecification, pBlock->fSize, iResult);
556     }
557   }
558   for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeCaloCluster | kAliHLTDataOriginAny); pBlock!=NULL; pBlock=GetNextInputBlock()) 
559     {
560       fBenchmark.AddInput(pBlock->fSize);
561       AliHLTCaloClusterHeaderStruct *caloClusterHeaderPtr = reinterpret_cast<AliHLTCaloClusterHeaderStruct*>(pBlock->fPtr);
562
563       HLTDebug("%d HLT clusters from spec: 0x%X", caloClusterHeaderPtr->fNClusters, pBlock->fSpecification);
564
565       //AliHLTCaloClusterReader reader;
566       //reader.SetMemory(caloClusterHeaderPtr);
567
568       AliHLTESDCaloClusterMaker clusterMaker;
569
570       int nClusters = clusterMaker.FillESD(pESD, caloClusterHeaderPtr);
571      
572       HLTInfo("converted %d cluster(s) to AliESDCaloCluster and added to ESD", nClusters);
573       iAddedDataBlocks++;
574     }
575   
576   // Add tracks from MUON.
577   for( const AliHLTComponentBlockData *i= GetFirstInputBlock(kAliHLTAnyDataType | kAliHLTDataOriginMUON); i!=NULL; i=GetNextInputBlock() ){
578     fBenchmark.AddInput(i->fSize);
579   }
580
581   for (const TObject* obj = GetFirstInputObject(kAliHLTAnyDataType | kAliHLTDataOriginMUON);
582        obj != NULL;
583        obj = GetNextInputObject()
584       )
585   {
586     const TClonesArray* tracklist = NULL;
587     if (obj->IsA() == AliESDEvent::Class())
588     {
589       const AliESDEvent* event = static_cast<const AliESDEvent*>(obj);
590       HLTDebug("Received a MUON ESD with specification: 0x%X", GetSpecification(obj));
591       if (event->GetList() == NULL) continue;
592       tracklist = dynamic_cast<const TClonesArray*>(event->GetList()->FindObject("MuonTracks"));
593       if (tracklist == NULL) continue;
594     }
595     else if (obj->IsA() == TClonesArray::Class())
596     {
597       tracklist = static_cast<const TClonesArray*>(obj);
598       HLTDebug("Received a MUON TClonesArray of tracks with specification: 0x%X", GetSpecification(obj));
599     }
600     else
601     {
602       // Cannot handle this object type.
603       continue;
604     }
605     HLTDebug("Received %d MUON tracks.", tracklist->GetEntriesFast());
606     if (tracklist->GetEntriesFast() > 0)
607     {
608       const AliESDMuonTrack* track = dynamic_cast<const AliESDMuonTrack*>(tracklist->UncheckedAt(0));
609       if (track == NULL)
610       {
611         HLTError(Form("%s from MUON does not contain AliESDMuonTrack objects.", obj->ClassName()));
612         continue;
613       }
614     }
615     for (Int_t i = 0; i < tracklist->GetEntriesFast(); ++i)
616     {
617       const AliESDMuonTrack* track = static_cast<const AliESDMuonTrack*>(tracklist->UncheckedAt(i));
618       pESD->AddMuonTrack(track);
619     }
620   }
621   
622   if (iAddedDataBlocks>0 && pTree) {
623     pTree->Fill();
624   }
625   
626   if (iResult>=0) iResult=iAddedDataBlocks;
627   return iResult;
628 }