]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/global/AliHLTGlobalFlatEsdConverterComponent.cxx
added placement new after reinterpret_cast in order to get vtable; added empty contru...
[u/mrichter/AliRoot.git] / HLT / global / AliHLTGlobalFlatEsdConverterComponent.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   AliHLTGlobalFlatEsdConverterComponent.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 "AliHLTGlobalFlatEsdConverterComponent.h"
33 #include "AliFlatESDEvent.h"
34 #include "AliFlatESDTrack.h"
35 #include "AliFlatExternalTrackParam.h"
36 #include "AliExternalTrackParam.h"
37
38 #include "AliHLTGlobalBarrelTrack.h"
39 #include "AliHLTExternalTrackParam.h"
40 #include "AliHLTTrackMCLabel.h"
41 #include "AliHLTCTPData.h"
42 #include "AliHLTErrorGuard.h"
43 #include "AliESDEvent.h"
44 #include "AliESDtrack.h"
45 #include "AliESDMuonTrack.h"
46 #include "AliESDMuonCluster.h"
47 #include "AliCDBEntry.h"
48 #include "AliCDBManager.h"
49 #include "AliPID.h"
50 #include "TTree.h"
51 #include "TList.h"
52 #include "TClonesArray.h"
53 //#include "AliHLTESDCaloClusterMaker.h"
54 //#include "AliHLTCaloClusterDataStruct.h"
55 //#include "AliHLTCaloClusterReader.h"
56 //#include "AliESDCaloCluster.h"
57 //#include "AliESDVZERO.h"
58 #include "AliHLTGlobalVertexerComponent.h"
59 #include "AliHLTVertexFinderBase.h"
60 #include "AliHLTTPCSpacePointData.h"
61 #include "AliHLTTPCClusterDataFormat.h"
62 #include "AliHLTTPCDefinitions.h"
63 #include "AliHLTTPCClusterMCData.h"
64 #include "AliHLTTPCTransform.h"
65
66 #include "AliSysInfo.h"
67
68 /** ROOT macro for the implementation of ROOT specific class methods */
69 ClassImp(AliHLTGlobalFlatEsdConverterComponent)
70
71 AliHLTGlobalFlatEsdConverterComponent::AliHLTGlobalFlatEsdConverterComponent()
72   : AliHLTProcessor()
73   , fWriteClusters(0)
74   , fVerbosity(0)  
75   , fSolenoidBz(-5.00668)
76   , fBenchmark("FlatEsdConverter")
77 {
78   // see header file for class documentation
79   // or
80   // refer to README to build package
81   // or
82   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
83
84 }
85
86 AliHLTGlobalFlatEsdConverterComponent::~AliHLTGlobalFlatEsdConverterComponent()
87 {
88   // see header file for class documentation
89 }
90
91 int AliHLTGlobalFlatEsdConverterComponent::Configure(const char* arguments)
92 {
93   // see header file for class documentation
94   int iResult=0;
95   if (!arguments) return iResult;
96
97   TString allArgs=arguments;
98   TString argument;
99   int bMissingParam=0;
100
101   TObjArray* pTokens=allArgs.Tokenize(" ");
102   if (pTokens) {
103     for (int i=0; i<pTokens->GetEntries() && iResult>=0; i++) {
104       argument=((TObjString*)pTokens->At(i))->String(); 
105       if (argument.IsNull()) continue;      
106       HLTError("unknown argument %s", argument.Data());
107       iResult=-EINVAL;
108       break;
109     }  
110     delete pTokens;
111   }
112   if (bMissingParam) {
113     HLTError("missing parameter for argument %s", argument.Data());
114     iResult=-EINVAL;
115   }
116
117   return iResult;
118 }
119
120 int AliHLTGlobalFlatEsdConverterComponent::Reconfigure(const char* cdbEntry, const char* chainId)
121 {
122   // see header file for class documentation
123   int iResult=0;
124   const char* path=NULL;
125   const char* defaultNotify="";
126   if (cdbEntry) {
127     path=cdbEntry;
128     defaultNotify=" (default)";
129   }
130   if (path) {
131     HLTInfo("reconfigure from entry %s%s, chain id %s", path, defaultNotify,(chainId!=NULL && chainId[0]!=0)?chainId:"<none>");
132     AliCDBEntry *pEntry = AliCDBManager::Instance()->Get(path/*,GetRunNo()*/);
133     if (pEntry) {
134       TObjString* pString=dynamic_cast<TObjString*>(pEntry->GetObject());
135       if (pString) {
136         HLTInfo("received configuration object string: \'%s\'", pString->String().Data());
137         iResult=Configure(pString->String().Data());
138       } else {
139         HLTError("configuration object \"%s\" has wrong type, required TObjString", path);
140       }
141     } else {
142       HLTError("can not fetch object \"%s\" from CDB", path);
143     }
144   }
145   
146   return iResult;
147 }
148
149 void AliHLTGlobalFlatEsdConverterComponent::GetInputDataTypes(AliHLTComponentDataTypeList& list)
150 {
151   // see header file for class documentation
152   list.push_back(kAliHLTDataTypeTrack);
153   list.push_back(kAliHLTDataTypeTrackMC);
154   list.push_back(kAliHLTDataTypeCaloCluster);
155   list.push_back(kAliHLTDataTypedEdx );
156   list.push_back(kAliHLTDataTypeESDVertex );
157   list.push_back(kAliHLTDataTypeESDObject);
158   list.push_back(kAliHLTDataTypeTObject);
159   list.push_back(kAliHLTDataTypeGlobalVertexer);
160   list.push_back(kAliHLTDataTypeV0Finder); // array of track ids for V0s
161   list.push_back(kAliHLTDataTypeKFVertex); // KFVertex object from vertexer
162   list.push_back(kAliHLTDataTypePrimaryFinder); // array of track ids for prim vertex
163   list.push_back(kAliHLTDataTypeESDContent);
164   list.push_back(AliHLTTPCDefinitions::fgkClustersDataType| kAliHLTDataOriginTPC);
165   list.push_back(AliHLTTPCDefinitions::fgkAliHLTDataTypeClusterMCInfo| kAliHLTDataOriginTPC);
166 }
167
168 AliHLTComponentDataType AliHLTGlobalFlatEsdConverterComponent::GetOutputDataType()
169 {
170   // see header file for class documentation
171   return kAliHLTDataTypeFlatESD|kAliHLTDataOriginOut;
172 }
173
174 void AliHLTGlobalFlatEsdConverterComponent::GetOutputDataSize(unsigned long& constBase, double& inputMultiplier)
175 {
176   // see header file for class documentation
177   constBase=2000000;
178   inputMultiplier=10.0;
179 }
180
181 int AliHLTGlobalFlatEsdConverterComponent::DoInit(int argc, const char** argv)
182 {
183   // see header file for class documentation
184   int iResult=0;
185   TString argument="";
186   int bMissingParam=0;
187
188   // default list of skiped ESD objects
189   TString skipObjects=
190     // "AliESDRun,"
191     // "AliESDHeader,"
192     // "AliESDZDC,"
193     "AliESDFMD,"
194     // "AliESDVZERO,"
195     // "AliESDTZERO,"
196     // "TPCVertex,"
197     // "SPDVertex,"
198     // "PrimaryVertex,"
199     // "AliMultiplicity,"
200     // "PHOSTrigger,"
201     // "EMCALTrigger,"
202     // "SPDPileupVertices,"
203     // "TrkPileupVertices,"
204     "Cascades,"
205     "Kinks,"
206     "AliRawDataErrorLogs,"
207     "AliESDACORDE";
208
209   iResult=Reconfigure(NULL, NULL);
210   TString allArgs = "";
211   for ( int i = 0; i < argc; i++ ) {
212     if ( !allArgs.IsNull() ) allArgs += " ";
213     allArgs += argv[i];
214   }
215
216   TObjArray* pTokens=allArgs.Tokenize(" ");
217   if (pTokens) {
218     for (int i=0; i<pTokens->GetEntries() && iResult>=0; i++) {
219       argument=((TObjString*)pTokens->At(i))->String(); 
220       if (argument.IsNull()) continue;
221
222       // -noclusters
223       if (argument.CompareTo("-noclusters")==0) {
224         fWriteClusters=0;       
225         // -clusters
226       } else if (argument.CompareTo("-clusters")==0) {
227         fWriteClusters=1;
228       } else if (argument.Contains("-skipobject=")) {
229         argument.ReplaceAll("-skipobject=", "");
230         skipObjects=argument;
231       } else {
232         HLTError("unknown argument %s", argument.Data());
233         iResult=-EINVAL;
234         break;
235       }
236     }
237   }
238   if (bMissingParam) {
239     HLTError("missing parameter for argument %s", argument.Data());
240     iResult=-EINVAL;
241   }
242
243   fSolenoidBz=GetBz();
244
245   if (iResult>=0) {
246     SetupCTPData();
247   }
248
249   fBenchmark.SetTimer(0,"total");
250
251   return iResult;
252 }
253
254 int AliHLTGlobalFlatEsdConverterComponent::DoDeinit()
255 {
256   // see header file for class documentation
257
258   return 0;
259 }
260
261 int AliHLTGlobalFlatEsdConverterComponent::DoEvent( const AliHLTComponentEventData& /*evtData*/,
262                                                     const AliHLTComponentBlockData* /*blocks*/, 
263                                                     AliHLTComponentTriggerData& /*trigData*/,
264                                                     AliHLTUInt8_t* outputPtr, 
265                                                     AliHLTUInt32_t& size,
266                                                     AliHLTComponentBlockDataList& outputBlocks)
267 {
268
269   // see header file for class documentation
270
271
272 AliSysInfo::AddStamp("DoEvent.Start");
273
274
275   int iResult=0;
276   bool benchmark = true;
277
278   if (!IsDataEvent()) return iResult;
279
280   fBenchmark.StartNewEvent();
281   fBenchmark.Start(0);
282   
283   
284
285   size_t maxOutputSize = size;
286   size = 0;
287
288   AliFlatESDEvent *flatEsd = reinterpret_cast<AliFlatESDEvent*>(outputPtr); 
289   new (flatEsd) AliFlatESDEvent;    
290
291   /*
292   pESD->Reset(); 
293   pESD->SetMagneticField(fSolenoidBz);
294   pESD->SetRunNumber(GetRunNo());
295   pESD->SetPeriodNumber(GetPeriodNumber());
296   pESD->SetOrbitNumber(GetOrbitNumber());
297   pESD->SetBunchCrossNumber(GetBunchCrossNumber());
298   pESD->SetTimeStamp(GetTimeStamp());
299   
300   const AliHLTCTPData* pCTPData=CTPData();
301   if (pCTPData) {
302     AliHLTUInt64_t mask=pCTPData->ActiveTriggers(trigData);
303     for (int index=0; index<gkNCTPTriggerClasses; index++) {
304       if ((mask&((AliHLTUInt64_t)0x1<<index)) == 0) continue;
305       pESD->SetTriggerClass(pCTPData->Name(index), index);
306     }
307     pESD->SetTriggerMask(mask);
308   }
309   */
310
311   // Barrel tracking
312   // tracks are based on the TPC tracks, and only updated from the ITS information
313   // Sequence:
314   // 1) extract MC information for TPC and ITS from specific data blocks and store in
315   //    intermediate vector arrays
316   // 2) extract TPC tracks, update with MC labels if available, the track parameters
317   //    are estimated at the first cluster position
318   // 2.1) propagate to last cluster position and update kTPCout, sets also outer param (fOp)
319   // 2.2) update kTPCin, sets also inner param (fIp) and TPC inner param (fTPCInner)
320   // 2.3) update kTPCrefit using the same parameters at the first cluster position
321   //      HLT has strictly spoking no refit, but we want the flag to be set
322   //      can be changed to be done after all the individual barrel detector parameters
323   //      have been updated by looping over the tracks again
324   // 3) extract ITS tracks, the tracks are actually TPC tracks updated from the ITS
325   //    tracking information
326   // 3.1) TODO 2010-07-12: handle ITS standalone tracks by updating kITSout before kITSin
327   // 3.2) update with kITSin
328   //    TODO 2010-07-12 find out if the kITSrefit has to be set as well
329   // 4) extract TRD tracks and add to ESD
330   //    TODO 2010-07-12 at the moment there is no matching or merging of TPC and TRD tracks
331   // 5) Add Trigger Detectors 
332   //    VZERO, ZDC
333
334   // 1) first read MC information (if present)
335
336   std::map<int,int> mcLabelsTPC;
337   std::map<int,int> mcLabelsITS;
338
339   for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeTrackMC|kAliHLTDataOriginTPC);
340        pBlock!=NULL; pBlock=GetNextInputBlock()) {
341     fBenchmark.AddInput(pBlock->fSize);
342     AliHLTTrackMCData* dataPtr = reinterpret_cast<AliHLTTrackMCData*>( pBlock->fPtr );
343     if (sizeof(AliHLTTrackMCData)+dataPtr->fCount*sizeof(AliHLTTrackMCLabel)==pBlock->fSize) {
344       for( unsigned int il=0; il<dataPtr->fCount; il++ ){
345         AliHLTTrackMCLabel &lab = dataPtr->fLabels[il];
346         mcLabelsTPC[lab.fTrackID] = lab.fMCLabel;
347       }
348     } else {
349       HLTWarning("data mismatch in block %s (0x%08x): count %d, size %d -> ignoring track MC information", 
350                  DataType2Text(pBlock->fDataType).c_str(), pBlock->fSpecification, 
351                  dataPtr->fCount, pBlock->fSize);
352     }
353   }
354  
355   for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeTrackMC|kAliHLTDataOriginITS);
356        pBlock!=NULL; pBlock=GetNextInputBlock()) {
357     fBenchmark.AddInput(pBlock->fSize);
358     AliHLTTrackMCData* dataPtr = reinterpret_cast<AliHLTTrackMCData*>( pBlock->fPtr );
359     if (sizeof(AliHLTTrackMCData)+dataPtr->fCount*sizeof(AliHLTTrackMCLabel)==pBlock->fSize) {
360       for( unsigned int il=0; il<dataPtr->fCount; il++ ){
361         AliHLTTrackMCLabel &lab = dataPtr->fLabels[il];
362         mcLabelsITS[lab.fTrackID] = lab.fMCLabel;
363       }
364     } else {
365       HLTWarning("data mismatch in block %s (0x%08x): count %d, size %d -> ignoring track MC information", 
366                  DataType2Text(pBlock->fDataType).c_str(), pBlock->fSpecification, 
367                  dataPtr->fCount, pBlock->fSize);
368     }
369   }
370
371   // 2) read dEdx information (if present)
372
373   AliHLTFloat32_t *dEdxTPC = 0; 
374   Int_t ndEdxTPC = 0;
375   for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypedEdx|kAliHLTDataOriginTPC);
376        pBlock!=NULL; pBlock=NULL/*GetNextInputBlock() there is only one block*/) {
377     fBenchmark.AddInput(pBlock->fSize);
378     dEdxTPC = reinterpret_cast<AliHLTFloat32_t*>( pBlock->fPtr );
379     ndEdxTPC = pBlock->fSize / (3*sizeof(AliHLTFloat32_t));
380   }
381
382   // 3) read  TPC tracks 
383
384   vector<AliHLTGlobalBarrelTrack> tracksTPC;
385
386   { // there is only one block of TPC data expected
387     const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeTrack|kAliHLTDataOriginTPC);
388     if( pBlock ){
389       fBenchmark.AddInput(pBlock->fSize);
390       iResult=AliHLTGlobalBarrelTrack::ConvertTrackDataArray(reinterpret_cast<const AliHLTTracksData*>(pBlock->fPtr), pBlock->fSize, tracksTPC);
391     }
392     if( iResult>=0 ){
393             HLTWarning("converted %d track(s) to AliESDtrack and added to ESD", tracksTPC.size());
394     } else if (iResult<0) {
395             HLTError("can not extract tracks from data block of type %s (specification %08x) of size %d: error %d", 
396                      DataType2Text(pBlock->fDataType).c_str(), pBlock->fSpecification, pBlock->fSize, iResult);
397     }
398   }
399  
400   // 4) read ITS refitted tracks
401
402   vector<AliHLTGlobalBarrelTrack> tracksITS;
403   vector<AliHLTGlobalBarrelTrack> tracksITSOut;
404
405   if( iResult>=0 ) { // there is only one block of ITS tracks expected
406      const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeTrack|kAliHLTDataOriginITS);
407      if( pBlock ){
408        fBenchmark.AddInput(pBlock->fSize);
409        iResult=AliHLTGlobalBarrelTrack::ConvertTrackDataArray(reinterpret_cast<const AliHLTTracksData*>(pBlock->fPtr), pBlock->fSize, tracksITS);
410      }
411   } 
412
413   if( iResult>=0 ) { // there is only one block of ITS tracks expected
414     const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeTrack|kAliHLTDataOriginITSOut);
415     if( pBlock ){
416       fBenchmark.AddInput(pBlock->fSize);
417       iResult=AliHLTGlobalBarrelTrack::ConvertTrackDataArray(reinterpret_cast<const AliHLTTracksData*>(pBlock->fPtr), pBlock->fSize, tracksITS);
418     }
419   } 
420
421    
422   // read TPC clusters
423   const UInt_t kNSlices = 36;
424   const UInt_t kNPatches = 6;
425
426   const AliHLTTPCClusterData *clustersTPC[kNSlices][kNPatches];
427   const AliHLTTPCClusterMCLabel *clustersTPCMC[kNSlices][kNPatches];
428   for( UInt_t i=0; i<kNSlices; i++){
429      for( UInt_t j=0; j<kNPatches; j++){
430       clustersTPC[i][j] = 0;
431       clustersTPCMC[i][j] =0;
432      }
433   }
434
435   fWriteClusters = 1;
436
437   if( fWriteClusters ){
438
439     for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(AliHLTTPCDefinitions::fgkClustersDataType| kAliHLTDataOriginTPC);
440          pBlock!=NULL; pBlock=GetNextInputBlock()) {
441       fBenchmark.AddInput(pBlock->fSize);
442       UInt_t slice     = AliHLTTPCDefinitions::GetMinSliceNr(*pBlock); 
443       UInt_t patch  = AliHLTTPCDefinitions::GetMinPatchNr(*pBlock);
444       if( slice >= kNSlices || patch>= kNPatches ){
445         HLTWarning("Wrong slice / patch number of cluster block");
446                           continue;
447       } 
448       clustersTPC[slice][patch] = reinterpret_cast <const AliHLTTPCClusterData*> ( pBlock->fPtr );
449     }
450     
451     for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(AliHLTTPCDefinitions::fgkAliHLTDataTypeClusterMCInfo| kAliHLTDataOriginTPC);
452          pBlock!=NULL; pBlock=GetNextInputBlock()) {
453       fBenchmark.AddInput(pBlock->fSize);
454       UInt_t slice     = AliHLTTPCDefinitions::GetMinSliceNr(*pBlock); 
455       UInt_t patch  = AliHLTTPCDefinitions::GetMinPatchNr(*pBlock);
456       if( slice >= kNSlices || patch>= kNPatches ){
457         HLTWarning("Wrong slice / patch number of cluster MC block");
458         continue;
459       } 
460       clustersTPCMC[slice][patch] = reinterpret_cast <const AliHLTTPCClusterMCLabel*> ( pBlock->fPtr );
461     }
462
463   }
464
465   // Fill vertex information to the flat ESD
466
467   const AliESDVertex *primaryVertex = 0;
468   {
469     const AliESDVertex *primaryVertexSPD = dynamic_cast<const AliESDVertex*>( GetFirstInputObject( kAliHLTDataTypeESDVertex|kAliHLTDataOriginITS ) );
470     const AliESDVertex *primaryVertexTracks = dynamic_cast<const AliESDVertex*>( GetFirstInputObject( kAliHLTDataTypeESDVertex|kAliHLTDataOriginOut ) );
471     
472   //  cout<<endl<<" Primary vertex Tracks: "<<primaryVertexTracks<<", SPD: "<< primaryVertexSPD <<endl<<endl;
473
474     flatEsd->FillPrimaryVertices( primaryVertexSPD, primaryVertexTracks );
475     
476     primaryVertex = primaryVertexTracks;
477     if( !primaryVertex ) primaryVertex = primaryVertexSPD;  
478   }
479
480   // Fill the track information to the flat ESD structure
481   {
482     UInt_t itsIter = 0;
483     UInt_t itsOutIter = 0;
484     for( UInt_t tpcIter=0; tpcIter < tracksTPC.size(); tpcIter++) {
485        
486        const AliHLTGlobalBarrelTrack *tpcTrack = &(tracksTPC[tpcIter]);
487        /*
488        Float_t points[4] = {
489           static_cast<Float_t>(tpcTrack->GetX()),
490           static_cast<Float_t>(tpcTrack->GetY()),
491           static_cast<Float_t>(tpcTrack->GetLastPointX()),
492           static_cast<Float_t>(tpcTrack->GetLastPointY())
493         };
494        */
495         Int_t tpcLabel = -1;
496         if( mcLabelsTPC.find(tpcTrack->TrackID())!=mcLabelsTPC.end() )
497           tpcLabel = mcLabelsTPC[tpcTrack->TrackID()];
498
499         //tpcTrack->SetLabel( tpcLabel );
500         // iotrack.SetID( tpcTrack->TrackID() );
501
502         // set kTPCout - just propagate to the outermost TPC cluster
503         
504         AliHLTGlobalBarrelTrack outPar(*tpcTrack);
505         {
506           //outPar.AliExternalTrackParam::PropagateTo( tpcTrack->GetLastPointX(), fSolenoidBz );
507           const Int_t N=10; // number of steps.
508           const Float_t xRange = tpcTrack->GetLastPointX() - tpcTrack->GetX();
509           const Float_t xStep = xRange / N ;
510           for(int i = 1; i <= N; ++i) {
511             if(!outPar.AliExternalTrackParam::PropagateTo(tpcTrack->GetX() + xStep * i, fSolenoidBz)) break;
512           }
513         }
514         
515         //iotrack.SetTPCPoints(points);
516         /*
517         if( tpcTrack->TrackID()<ndEdxTPC ){
518           AliHLTFloat32_t *val = &(dEdxTPC[3*tpcTrack->TrackID()]);
519           iotrack.SetTPCsignal( val[0], val[1], (UChar_t) val[2] ); 
520           //AliTPCseed s;
521           //s.Set( tpcTrack->GetX(), tpcTrack->GetAlpha(),
522           //tpcTrack->GetParameter(), tpcTrack->GetCovariance() );
523           //s.SetdEdx( val[0] );
524           //s.CookPID();
525           //iotrack.SetTPCpid(s.TPCrPIDs() );
526         } else {
527           if( dEdxTPC ) HLTWarning("Wrong number of dEdx TPC labels");
528         }
529         */
530         //iotrack.SetLabel(mcLabel);
531
532         // ITS track
533
534         AliHLTGlobalBarrelTrack *itsRefit=0;
535         Int_t itsLabel = -1;
536         
537         for(; itsIter< tracksITS.size() && tracksITS[itsIter].TrackID()<(int) tpcIter; itsIter++ );
538
539         if( itsIter< tracksITS.size() && tracksITS[itsIter].TrackID() == (int) tpcIter ){
540           itsRefit = &(tracksITS[itsIter]);
541           if( mcLabelsITS.find(tpcIter)!=mcLabelsITS.end() ) itsLabel = mcLabelsITS[tpcIter];
542           itsIter++;
543         }
544
545         // ITS Out track
546
547         AliHLTGlobalBarrelTrack *itsOut=0;
548         
549         for(; itsOutIter< tracksITSOut.size() && tracksITSOut[itsOutIter].TrackID()<(int) tpcIter; itsOutIter++ );
550
551         if( itsOutIter< tracksITSOut.size() && tracksITSOut[itsOutIter].TrackID() == (int) tpcIter ){
552           itsOut = &(tracksITSOut[itsOutIter]);
553           itsOutIter++;
554         }
555
556         // Fill DCA parameters for TPC tracks
557         AliESDtrack cP;
558         
559         if( primaryVertex ){
560             cP.UpdateTrackParams( (itsRefit ?itsRefit :tpcTrack), AliESDtrack::kTPCin );
561             cP.RelateToVertex( primaryVertex, fSolenoidBz, 1000 );      
562         }
563         
564         AliFlatESDTrack *flatTrack = flatEsd->GetNextTrackPointer();
565         //cout<<"flatTrack: "<<flatTrack<<endl;
566
567         //cout<<"GetNumberOfTPCClusters before: "<<flatTrack->GetNumberOfTPCClusters()<<endl;
568
569         UInt_t nClustersTPC = tpcTrack->GetNumberOfPoints();
570         UInt_t nClustersITS = itsRefit ?itsRefit->GetNumberOfPoints() :0;
571
572         flatTrack->SetNumberOfITSClusters( nClustersITS );
573         //flatTrack->SetNumberOfTPCClusters(0);
574
575         //cout<<"GetNumberOfTPCClusters: "<<flatTrack->GetNumberOfTPCClusters()<<endl;
576
577         if( flatEsd->GetSize() + flatTrack->EstimateSize( kTRUE, nClustersTPC ) >= maxOutputSize ){
578                 cout<<endl<<endl<<"NOT ENOUGH MEMORY!!!!"<<endl<<endl;
579            iResult=-ENOMEM;
580            break;
581         }
582         
583         flatTrack->FillExternalTrackParam( itsRefit,  NULL, tpcTrack, &outPar, cP.GetConstrainedParam(), itsOut);
584         
585         if( fWriteClusters && tpcTrack->GetPoints() ){
586           const UInt_t* clusterIDs = tpcTrack->GetPoints();
587            for( UInt_t i=0; i<nClustersTPC; i++ ){
588               UInt_t id = clusterIDs[i];
589               UInt_t iSlice = AliHLTTPCSpacePointData::GetSlice(id);
590               UInt_t iPatch = AliHLTTPCSpacePointData::GetPatch(id);
591               UInt_t iCluster = AliHLTTPCSpacePointData::GetNumber(id);
592               if( iSlice >= kNSlices || iPatch>= kNPatches ){
593                  HLTWarning("Wrong slice / patch number of TPC cluster");
594                  continue;
595               } 
596               const AliHLTTPCClusterData *clusterBlock = clustersTPC[iSlice][iPatch];
597               if( !clusterBlock ){
598                  HLTWarning("no cluster block found for slice %d, patch %d",iSlice,iPatch);
599                  continue;
600               } 
601               if( iCluster>= clusterBlock->fSpacePointCnt ){
602                       HLTWarning("no cluster block found for slice %d, patch %d, cluster %d",iSlice,iPatch,iCluster);
603                  continue;
604               } 
605               const AliHLTTPCSpacePointData &cIn = clusterBlock->fSpacePoints[iCluster];
606
607               AliFlatTPCCluster *c= flatTrack->GetTPCCluster( flatTrack->GetNumberOfTPCClusters() );;
608               c->fX = cIn.GetX();
609               c->fY = cIn.GetY();
610               c->fZ = cIn.GetZ();
611               c->fPadRow  = cIn.GetPadRow() + AliHLTTPCTransform::GetFirstRow(iPatch);
612               c->fSigmaY2 = cIn.GetSigmaY2();
613               c->fSigmaZ2 = cIn.GetSigmaZ2();
614               c->fCharge  = cIn.GetCharge();
615               c->fQMax    = cIn.GetQMax();
616               flatTrack->StoreLastTPCCluster();
617            }
618         }
619
620                 //      cout<<"number of tpc clusters: "<<flatTrack->GetNumberOfTPCClusters()<<endl;
621                         //cout<<"number of its clusters: "<<flatTrack->GetNumberOfITSClusters()<<endl;
622         
623         flatEsd->StoreLastTrack();
624         
625         if (fVerbosity>0) tpcTrack->Print();
626     }    
627   }
628
629   // Fill v0's
630   
631     int nV0s =0;
632   {    
633     const AliHLTComponentBlockData* pP = GetFirstInputBlock(kAliHLTDataTypeGlobalVertexer|kAliHLTDataOriginOut);
634     if (pP && pP->fSize && pP->fPtr) {
635                 fBenchmark.AddInput(pP->fSize);
636       const AliHLTGlobalVertexerComponent::AliHLTGlobalVertexerData *data = reinterpret_cast<AliHLTGlobalVertexerComponent::AliHLTGlobalVertexerData*>(pP->fPtr);
637       const int* v0s = data->fTrackIndices + data->fNPrimTracks;
638       nV0s = data->fNV0s;
639       for (int i = 0; i < nV0s; ++i) {
640         AliFlatESDV0 *v0 = flatEsd->GetNextV0Pointer();
641         v0->fNegTrackID = v0s[2 * i];
642         v0->fPosTrackID = v0s[2 * i + 1];
643         flatEsd->StoreLastV0();
644       }
645     } else {
646       HLTWarning(" No V0 data block");
647     }
648   }
649   
650   // Get ITS SPD vertex
651   for( const AliHLTComponentBlockData *i= GetFirstInputBlock(kAliHLTDataTypeESDVertex|kAliHLTDataOriginITS); i!=NULL; i=GetNextInputBlock() ){
652     fBenchmark.AddInput(i->fSize);
653   }
654   // Get Track vertex
655   for( const AliHLTComponentBlockData *i= GetFirstInputBlock(kAliHLTDataTypeESDVertex|kAliHLTDataOriginOut); i!=NULL; i=GetNextInputBlock() ){
656     fBenchmark.AddInput(i->fSize);
657   }
658
659
660
661   /*
662   for ( const TObject *iter = GetFirstInputObject(kAliHLTDataTypeESDVertex|kAliHLTDataOriginITS); iter != NULL; iter = GetNextInputObject() ) {
663     AliESDVertex *vtx = dynamic_cast<AliESDVertex*>(const_cast<TObject*>( iter ) );
664     pESD->SetPrimaryVertexSPD( vtx );
665   }
666   
667
668   // update with  vertices and vertex-fitted tracks
669   // output of the GlobalVertexerComponent
670   for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeGlobalVertexer);
671        pBlock!=NULL; pBlock=GetNextInputBlock()) {
672     fBenchmark.AddInput(pBlock->fSize);   
673     AliHLTGlobalVertexerComponent::FillESD( pESD, reinterpret_cast<AliHLTGlobalVertexerComponent::AliHLTGlobalVertexerData* >(pBlock->fPtr) );
674   }
675
676   // update with  vertices and vertex-fitted tracks
677   // output of PrimaryVertexer and V0Finder components
678   TObject* pBase = (TObject*)GetFirstInputObject(kAliHLTDataTypeKFVertex | kAliHLTDataOriginOut);
679   if (pBase) {
680     AliKFVertex* kfVertex = dynamic_cast<AliKFVertex *>(pBase);
681     if (kfVertex) {
682       const AliHLTComponentBlockData* pP = GetFirstInputBlock(kAliHLTDataTypePrimaryFinder | kAliHLTDataOriginOut);
683       if (pP && pP->fSize && pP->fPtr) {
684         const AliHLTComponentBlockData* pV0 = GetFirstInputBlock(kAliHLTDataTypeV0Finder | kAliHLTDataOriginOut);
685         if (pV0 && pV0->fPtr && pInputESD && pInputESD->GetNumberOfV0s()>0) {
686           const int* v0s = static_cast<const int*>(pV0->fPtr);
687           HLTWarning("V0 array already filled from the input esd block, additional filling from V0 block of %d entries might cause inconsistent content", v0s[0]);
688         }
689         AliHLTVertexFinderBase::FillESD(pESD, kfVertex, pP->fPtr, pV0?pV0->fPtr:NULL);
690       } else
691         HLTWarning("Problem with primary finder's data block");
692     } else {
693       HLTWarning("primary vertex block of wrong type, expecting AliKFVertex instead of %s", pBase->GetName());
694     }
695   } else {
696     // throw an error if there is a V0 data block which can not be handled without
697     // the AliKFVertex object
698     if (GetFirstInputBlock(kAliHLTDataTypeV0Finder | kAliHLTDataOriginOut)!=NULL) {
699       ALIHLTERRORGUARD(1, "missing AliKFVertex object ignoring V0 data block of type %s",
700                        DataType2Text(kAliHLTDataTypeV0Finder|kAliHLTDataOriginOut).c_str());
701     }
702   }
703   */
704
705
706   // loop over all tracks and set the TPC refit flag by updating with the
707   // original TPC inner parameter if not yet set
708   // TODO: replace this by a proper refit
709   // code is comented for the moment as it does not fully solve the problems with
710   // the display
711   // - would set the main parameters to the TPC inner wall again, or
712   // - changes the inner param if the parameters are propagated, so we loose the track
713   //   reference point for the display
714   // with the current sequence we have the latter case as the DCA operations above
715   // change the TPC inner parameters
716   /*
717   for (int i=0; i<pESD->GetNumberOfTracks(); i++) {
718     if (!pESD->GetTrack(i) || 
719         !pESD->GetTrack(i)->GetTPCInnerParam() ||
720         pESD->GetTrack(i)->IsOn(AliESDtrack::kTPCrefit)) continue;
721     AliESDtrack* tESD=pESD->GetTrack(i);
722     AliHLTGlobalBarrelTrack inner(*tESD->GetTPCInnerParam());
723     inner.SetLabel(tESD->GetLabel());
724     tESD->UpdateTrackParams(&inner, AliESDtrack::kTPCrefit);
725   }
726   */
727  
728   if (iResult>=0) {            
729  
730     AliHLTComponentBlockData outBlock;
731     FillBlockData( outBlock );
732     outBlock.fOffset = size;
733     outBlock.fSize = flatEsd->GetSize();
734     outBlock.fDataType = kAliHLTDataTypeFlatESD|kAliHLTDataOriginOut;
735     outBlock.fSpecification = AliHLTTPCDefinitions::EncodeDataSpecification( 0, 35, 0, 5 );
736
737     outputBlocks.push_back( outBlock );
738
739     fBenchmark.AddOutput(outBlock.fSize);
740       
741     size += outBlock.fSize;
742   }
743
744   fBenchmark.Stop(0);
745   HLTWarning( fBenchmark.GetStatistics() );
746   
747   
748   
749   if(benchmark){
750         
751         Double_t statistics[10]; 
752         TString names[10];
753         fBenchmark.GetStatisticsData(statistics, names);
754         //  statistics[5] = tracksTPC.size();
755         //  statistics[7] = nV0s;
756           
757         //  FillBenchmarkHistos( statistics, names);
758           fBenchmark.Reset();
759
760         AliSysInfo::AddStamp("DoEvent.Stop", (int)(statistics[1]), (int)(statistics[2]) );
761   
762   }
763   return iResult;
764
765
766
767
768 }
769
770