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