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