]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/global/AliHLTGlobalFlatEsdConverterComponent.cxx
flat friends update
[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 "AliFlatESDTrigger.h"
35 #include "AliFlatESDVertex.h"
36 #include "AliFlatESDV0.h"
37 #include "AliFlatESDTrack.h"
38 #include "AliFlatExternalTrackParam.h"
39 #include "AliExternalTrackParam.h"
40
41 #include "AliHLTGlobalBarrelTrack.h"
42 #include "AliHLTExternalTrackParam.h"
43 #include "AliHLTTrackMCLabel.h"
44 #include "AliHLTCTPData.h"
45 #include "AliHLTErrorGuard.h"
46 #include "AliESDEvent.h"
47 #include "AliESDtrack.h"
48 #include "AliESDMuonTrack.h"
49 #include "AliESDMuonCluster.h"
50 #include "AliCDBEntry.h"
51 #include "AliCDBManager.h"
52 #include "AliPID.h"
53 #include "TTree.h"
54 #include "TList.h"
55 #include "TClonesArray.h"
56 //#include "AliHLTESDCaloClusterMaker.h"
57 //#include "AliHLTCaloClusterDataStruct.h"
58 //#include "AliHLTCaloClusterReader.h"
59 //#include "AliESDCaloCluster.h"
60 //#include "AliESDVZERO.h"
61 #include "AliHLTGlobalVertexerComponent.h"
62 #include "AliHLTVertexFinderBase.h"
63 #include "AliHLTTPCSpacePointData.h"
64 #include "AliHLTTPCClusterDataFormat.h"
65 #include "AliHLTTPCDefinitions.h"
66 #include "AliHLTTPCClusterMCData.h"
67 #include "AliHLTTPCTransform.h"
68
69 /** ROOT macro for the implementation of ROOT specific class methods */
70 ClassImp(AliHLTGlobalFlatEsdConverterComponent)
71
72 AliHLTGlobalFlatEsdConverterComponent::AliHLTGlobalFlatEsdConverterComponent()
73   : AliHLTProcessor()
74   , fVerbosity(0)    
75   , fBenchmark("FlatEsdConverter")
76   , fProduceFriend(1)
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 AliHLTGlobalFlatEsdConverterComponent::~AliHLTGlobalFlatEsdConverterComponent()
86 {
87   // see header file for class documentation
88 }
89
90 int AliHLTGlobalFlatEsdConverterComponent::Configure(const char* arguments)
91 {
92   // see header file for class documentation
93   int iResult=0;
94   if (!arguments) return iResult;
95
96   TString allArgs=arguments;
97   TString argument;
98   int bMissingParam=0;
99
100   TObjArray* pTokens=allArgs.Tokenize(" ");
101   if (pTokens) {
102     for (int i=0; i<pTokens->GetEntries() && iResult>=0; i++) {
103       argument=((TObjString*)pTokens->At(i))->String(); 
104       if (argument.IsNull()) continue;      
105       HLTError("unknown argument %s", argument.Data());
106       iResult=-EINVAL;
107       break;
108     }  
109     delete pTokens;
110   }
111   if (bMissingParam) {
112     HLTError("missing parameter for argument %s", argument.Data());
113     iResult=-EINVAL;
114   }
115
116   return iResult;
117 }
118
119 int AliHLTGlobalFlatEsdConverterComponent::Reconfigure(const char* cdbEntry, const char* chainId)
120 {
121   // see header file for class documentation
122   int iResult=0;
123   const char* path=NULL;
124   const char* defaultNotify="";
125   if (cdbEntry) {
126     path=cdbEntry;
127     defaultNotify=" (default)";
128   }
129   if (path) {
130     HLTInfo("reconfigure from entry %s%s, chain id %s", path, defaultNotify,(chainId!=NULL && chainId[0]!=0)?chainId:"<none>");
131     AliCDBEntry *pEntry = AliCDBManager::Instance()->Get(path/*,GetRunNo()*/);
132     if (pEntry) {
133       TObjString* pString=dynamic_cast<TObjString*>(pEntry->GetObject());
134       if (pString) {
135         HLTInfo("received configuration object string: \'%s\'", pString->String().Data());
136         iResult=Configure(pString->String().Data());
137       } else {
138         HLTError("configuration object \"%s\" has wrong type, required TObjString", path);
139       }
140     } else {
141       HLTError("can not fetch object \"%s\" from CDB", path);
142     }
143   }
144   
145   return iResult;
146 }
147
148 void AliHLTGlobalFlatEsdConverterComponent::GetInputDataTypes(AliHLTComponentDataTypeList& list)
149 {
150   // see header file for class documentation
151   list.push_back(kAliHLTDataTypeTrack);
152   list.push_back(kAliHLTDataTypeTrackMC);
153   list.push_back(kAliHLTDataTypeCaloCluster);
154   list.push_back(kAliHLTDataTypedEdx );
155   list.push_back(kAliHLTDataTypeESDVertex );
156   list.push_back(kAliHLTDataTypeESDObject);
157   list.push_back(kAliHLTDataTypeTObject);
158   list.push_back(kAliHLTDataTypeGlobalVertexer);
159   list.push_back(kAliHLTDataTypeV0Finder); // array of track ids for V0s
160   list.push_back(kAliHLTDataTypeKFVertex); // KFVertex object from vertexer
161   list.push_back(kAliHLTDataTypePrimaryFinder); // array of track ids for prim vertex
162   list.push_back(kAliHLTDataTypeESDContent);
163   list.push_back(AliHLTTPCDefinitions::fgkClustersDataType| kAliHLTDataOriginTPC);
164   list.push_back(AliHLTTPCDefinitions::fgkAliHLTDataTypeClusterMCInfo| kAliHLTDataOriginTPC);
165 }
166
167 AliHLTComponentDataType AliHLTGlobalFlatEsdConverterComponent::GetOutputDataType()
168 {
169   // see header file for class documentation
170   return kAliHLTDataTypeFlatESD|kAliHLTDataOriginOut;
171 }
172
173 void AliHLTGlobalFlatEsdConverterComponent::GetOutputDataSize(unsigned long& constBase, double& inputMultiplier)
174 {
175   // see header file for class documentation
176   constBase=2000000;
177   inputMultiplier=10.0;
178 }
179
180 int AliHLTGlobalFlatEsdConverterComponent::DoInit(int argc, const char** argv)
181 {
182   // see header file for class documentation
183   int iResult=0;
184   TString argument="";
185   int bMissingParam=0;
186
187   iResult=Reconfigure(NULL, NULL);
188   TString allArgs = "";
189   for ( int i = 0; i < argc; i++ ) {
190     if ( !allArgs.IsNull() ) allArgs += " ";
191     allArgs += argv[i];
192   }
193
194   TObjArray* pTokens=allArgs.Tokenize(" ");
195   if (pTokens) {
196     for (int i=0; i<pTokens->GetEntries() && iResult>=0; i++) {
197       argument=((TObjString*)pTokens->At(i))->String(); 
198       if (argument.IsNull()) continue;
199
200       // -noclusters
201       if (argument.CompareTo("-nofriend")==0) {
202         fProduceFriend=0;       
203         // -clusters
204       } else if (argument.CompareTo("-friend")==0) {
205         fProduceFriend=1;
206       } else {
207         HLTError("unknown argument %s", argument.Data());
208         iResult=-EINVAL;
209         break;
210       }
211     }
212   }
213   if (bMissingParam) {
214     HLTError("missing parameter for argument %s", argument.Data());
215     iResult=-EINVAL;
216   }
217
218   if (iResult>=0) {
219     SetupCTPData();
220   }
221
222   fBenchmark.SetTimer(0,"total");
223
224   return iResult;
225 }
226
227 int AliHLTGlobalFlatEsdConverterComponent::DoDeinit()
228 {
229   // see header file for class documentation
230
231   return 0;
232 }
233
234 int AliHLTGlobalFlatEsdConverterComponent::DoEvent( const AliHLTComponentEventData& /*evtData*/,
235                                                     const AliHLTComponentBlockData* /*blocks*/, 
236                                                     AliHLTComponentTriggerData& trigData,
237                                                     AliHLTUInt8_t* outputPtr, 
238                                                     AliHLTUInt32_t& size,
239                                                     AliHLTComponentBlockDataList& outputBlocks )
240 {
241   // see header file for class documentation
242
243   int iResult=0;
244
245   if (!IsDataEvent()) return iResult;
246
247   fBenchmark.StartNewEvent();
248   fBenchmark.Start(0);
249
250   size_t maxOutputSize = size;
251   size = 0;
252
253   // Read part of the input  data to local arrays
254   
255   // 1) first, read MC information if present
256
257   std::map<int,int> mcLabelsTPC;
258
259   for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeTrackMC|kAliHLTDataOriginTPC);
260        pBlock!=NULL; pBlock=GetNextInputBlock()) {
261     fBenchmark.AddInput(pBlock->fSize);
262     AliHLTTrackMCData* dataPtr = reinterpret_cast<AliHLTTrackMCData*>( pBlock->fPtr );
263     if (sizeof(AliHLTTrackMCData)+dataPtr->fCount*sizeof(AliHLTTrackMCLabel)==pBlock->fSize) {
264       for( unsigned int il=0; il<dataPtr->fCount; il++ ){
265         AliHLTTrackMCLabel &lab = dataPtr->fLabels[il];
266         mcLabelsTPC[lab.fTrackID] = lab.fMCLabel;
267       }
268     } else {
269       HLTWarning("data mismatch in block %s (0x%08x): count %d, size %d -> ignoring track MC information", 
270                  DataType2Text(pBlock->fDataType).c_str(), pBlock->fSpecification, 
271                  dataPtr->fCount, pBlock->fSize);
272     }
273   }
274  
275   // 2) read dEdx information (if present)
276
277   AliHLTFloat32_t *dEdxTPC = 0; 
278   Int_t ndEdxTPC = 0;
279   for (const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypedEdx|kAliHLTDataOriginTPC);
280        pBlock!=NULL; pBlock=NULL/*GetNextInputBlock() there is only one block*/) {
281     fBenchmark.AddInput(pBlock->fSize);
282     dEdxTPC = reinterpret_cast<AliHLTFloat32_t*>( pBlock->fPtr );
283     ndEdxTPC = pBlock->fSize / (3*sizeof(AliHLTFloat32_t));
284   }
285
286   // 3) read  TPC tracks, ITS refitted tracks, ITS OUT tracks
287
288   vector<AliHLTGlobalBarrelTrack> tracksTPC;
289   vector<AliHLTGlobalBarrelTrack> tracksTPCOut;
290   vector<AliHLTGlobalBarrelTrack> tracksITS;
291   vector<AliHLTGlobalBarrelTrack> tracksITSOut;
292
293   if( iResult>=0 ){
294     const AliHLTComponentBlockData* pBlock=GetFirstInputBlock(kAliHLTDataTypeTrack|kAliHLTDataOriginTPC);
295     if( pBlock ){
296       fBenchmark.AddInput(pBlock->fSize);
297       iResult = AliHLTGlobalBarrelTrack::ConvertTrackDataArray(reinterpret_cast<const AliHLTTracksData*>(pBlock->fPtr), pBlock->fSize, tracksTPC);
298     }  
299   
300     if( iResult>=0 ) { 
301       pBlock=GetFirstInputBlock(kAliHLTDataTypeTrack|kAliHLTDataOriginITS);
302       if( pBlock ){
303         fBenchmark.AddInput(pBlock->fSize);
304         iResult=AliHLTGlobalBarrelTrack::ConvertTrackDataArray(reinterpret_cast<const AliHLTTracksData*>(pBlock->fPtr), pBlock->fSize, tracksITS);
305       }
306     }
307     if( iResult>=0 ) { 
308       pBlock=GetFirstInputBlock(kAliHLTDataTypeTrack|kAliHLTDataOriginITSOut);
309       if( pBlock ){
310         fBenchmark.AddInput(pBlock->fSize);
311         iResult=AliHLTGlobalBarrelTrack::ConvertTrackDataArray(reinterpret_cast<const AliHLTTracksData*>(pBlock->fPtr), pBlock->fSize, tracksITSOut);
312       } 
313     }   
314     if( iResult<0 ){
315       HLTError("can not extract tracks from data block of type %s (specification %08x) of size %d: error %d", 
316                DataType2Text(pBlock->fDataType).c_str(), pBlock->fSpecification, pBlock->fSize, iResult);     
317     }
318   }
319
320   HLTWarning("converted %d TPC %d ITS %d ITSout track(s) to GlobalBarrelTrack", tracksTPC.size(), tracksITS.size(), tracksITSOut.size() );
321   
322   // Set TPC MC labels to tracks
323   for( UInt_t itr=0; itr < tracksTPC.size(); itr++) {
324     AliHLTGlobalBarrelTrack &track = tracksTPC[itr];
325     std::map<int,int>::const_iterator lab = mcLabelsTPC.find( track.TrackID() );
326     if( lab!=mcLabelsTPC.end() ) track.SetLabel( lab->second );
327     else track.SetLabel( -1 );
328   }
329
330   // Create TPC Out tracks - just propagate to the outermost TPC cluster
331   for( UInt_t itr=0; itr < tracksTPC.size(); itr++) {
332     tracksTPCOut.push_back( tracksTPC[itr] );
333     AliHLTGlobalBarrelTrack &track = tracksTPCOut.back();
334     const Int_t N=10; // number of steps.
335     const Float_t xRange = track.GetLastPointX() - track.GetX();
336     const Float_t xStep = xRange / N ;
337     for(int i = 1; i <= N; ++i) {
338       if(!track.AliExternalTrackParam::PropagateTo(track.GetX() + xStep, GetBz() )) break;
339     }
340   }
341
342
343   // ---------------------------------------------
344   //
345   // Start to fill the flat ESD structure
346   //
347
348   Int_t err = 0;
349   AliFlatESDEvent *flatEsd = reinterpret_cast<AliFlatESDEvent*>(outputPtr); 
350
351   do{ // single loop for easy break in case of output buffer overflow
352
353     size_t freeSpace = maxOutputSize;
354
355     err = ( freeSpace < sizeof( AliFlatESDEvent ) );    
356     if( err ) break;
357
358     new (flatEsd) AliFlatESDEvent;    
359  
360     freeSpace -= flatEsd->GetSize();
361   
362     // fill run info
363     {
364       flatEsd->SetMagneticField( GetBz() );
365       flatEsd->SetPeriodNumber( GetPeriodNumber() );
366       flatEsd->SetRunNumber( GetRunNo() );
367       flatEsd->SetOrbitNumber( GetOrbitNumber() );
368       flatEsd->SetBunchCrossNumber( GetBunchCrossNumber() );
369       flatEsd->SetTimeStamp( GetTimeStamp() );
370       //flatEsd->SetEventSpecie( GetEventSpecie() ); !!SG!! to do
371     }
372
373     // Fill trigger information  
374     {
375       const AliHLTCTPData* pCTPData=CTPData();
376       if (pCTPData) {
377         size_t triggerSize = 0;
378         int nTriggers = 0;
379         AliFlatESDTrigger *trigger = flatEsd->SetTriggersStart();
380         AliHLTTriggerMask_t mask = pCTPData->ActiveTriggers(trigData);
381         for (int index=0; index<gkNCTPTriggerClasses; index++) {
382           if ((mask&(AliHLTTriggerMask_t(0x1)<<index)) == 0) continue;
383           const char* name = pCTPData->Name(index);
384           if( name && name[0]!='\0' ){
385             err = trigger->SetTriggerClass( name, index, freeSpace );
386             if( err != 0 ) break;
387             nTriggers++;
388             freeSpace -= trigger->GetSize();
389             triggerSize += trigger->GetSize();
390             trigger = trigger->GetNextTriggerNonConst();
391           }
392         }    
393         flatEsd->SetTriggersEnd( nTriggers, triggerSize );
394         //first 50 triggers
395         AliHLTTriggerMask_t mask50;
396         mask50.set(); // set all bits
397         mask50 >>= 50; // shift 50 right
398         flatEsd->SetTriggerMask((mask&mask50).to_ulong());
399         //next 50
400         flatEsd->SetTriggerMaskNext50((mask>>50).to_ulong());
401       }
402     }
403  
404     if( err ) break;
405     
406     const AliESDVertex *primaryVertex = 0;
407
408     { // fill primary vertex Tracks
409
410       const AliESDVertex *primaryVertexTracks = dynamic_cast<const AliESDVertex*>( GetFirstInputObject( kAliHLTDataTypeESDVertex|kAliHLTDataOriginOut ) );               
411       primaryVertex = primaryVertexTracks;
412       err = flatEsd->SetPrimaryVertexTracks( primaryVertexTracks, freeSpace );
413       freeSpace = maxOutputSize - flatEsd->GetSize();
414     }
415
416     if( err ) break;
417
418     { // fill primary vertex SPD
419     
420       const AliESDVertex *primaryVertexSPD = dynamic_cast<const AliESDVertex*>( GetFirstInputObject( kAliHLTDataTypeESDVertex|kAliHLTDataOriginITS ) );
421       if( !primaryVertex ) primaryVertex = primaryVertexSPD;  
422       err = flatEsd->SetPrimaryVertexSPD( primaryVertexSPD, freeSpace );
423       freeSpace = maxOutputSize - flatEsd->GetSize();     
424     }  
425     
426     if( err ) break;
427
428     { // Fill track information to the flat ESD structure
429      
430       size_t trackSize = 0;
431       int nTracks = 0;
432       Long64_t *table = NULL;
433       AliFlatESDTrack *flatTrack = NULL;
434       err = flatEsd->SetTracksStart( flatTrack, table, tracksTPC.size(), freeSpace );
435       freeSpace = maxOutputSize - flatEsd->GetSize();
436     
437       if( err ) break;
438
439       for( UInt_t tpcIter=0, itsIter = 0, itsOutIter = 0; tpcIter < tracksTPC.size(); tpcIter++) {
440
441         // ----------------------------
442         // -- read track information
443
444         // TPC track parameters
445         
446         AliHLTGlobalBarrelTrack *tpcTrack = &(tracksTPC[tpcIter]);
447         AliHLTGlobalBarrelTrack *tpcOutTrack = &(tracksTPCOut[tpcIter]);
448
449         // ITS track parameters
450         
451         AliHLTGlobalBarrelTrack *itsRefit=0;
452         AliHLTGlobalBarrelTrack *itsOut=0;
453         
454         // ITS Refit track
455           
456         for(; itsIter< tracksITS.size() && tracksITS[itsIter].TrackID()<(int) tpcIter; itsIter++ );
457         
458         if( itsIter< tracksITS.size() && tracksITS[itsIter].TrackID() == (int) tpcIter ){
459           itsRefit = &(tracksITS[itsIter]);
460           itsIter++;
461         }
462         
463         // ITS Out track
464           
465         for(; itsOutIter< tracksITSOut.size() && tracksITSOut[itsOutIter].TrackID()<(int) tpcIter; itsOutIter++ );
466         
467         if( itsOutIter< tracksITSOut.size() && tracksITSOut[itsOutIter].TrackID() == (int) tpcIter ){
468           itsOut = &(tracksITSOut[itsOutIter]);
469           itsOutIter++;
470         }       
471         
472         // 
473
474         if (fVerbosity>0) tpcTrack->Print();  
475
476         Float_t tpcDeDx[3]={0,0,0};
477         
478         if( ndEdxTPC>0 ){       
479           if( tpcTrack->TrackID() < ndEdxTPC ){
480             AliHLTFloat32_t *val = &(dEdxTPC[3*tpcTrack->TrackID()]);
481             tpcDeDx[0] = val[0];
482             tpcDeDx[1] = val[1];
483             tpcDeDx[2] = val[2];
484           } else {
485             HLTWarning("Wrong number of dEdx TPC labels");
486           }
487         }
488     
489         // vertex-constrained parameters for TPC tracks 
490
491         const AliExternalTrackParam *tpcConstrained=0;       
492         AliESDtrack esdTrack;
493         if( primaryVertex ){
494           esdTrack.UpdateTrackParams(&(*tpcTrack),AliESDtrack::kTPCrefit);        
495           esdTrack.RelateToVertex( primaryVertex, GetBz(), 1000 );      
496           tpcConstrained = esdTrack.GetConstrainedParam();      
497         }
498
499         UInt_t nClustersTPC = tpcTrack->GetNumberOfPoints();
500         UInt_t nClustersITS = itsRefit ?itsRefit->GetNumberOfPoints() :0;
501
502         
503         // ----------------------------
504         // -- fill flat track structure
505
506         table[tpcIter] = trackSize;
507         err = ( freeSpace < flatTrack->EstimateSize() );
508         if( err ) break;
509         
510         new (flatTrack) AliFlatESDTrack;       
511         
512         flatTrack->SetExternalTrackParam( itsRefit, itsRefit, tpcTrack, tpcOutTrack, tpcConstrained, itsOut );
513         flatTrack->SetNumberOfTPCClusters( nClustersTPC );
514         flatTrack->SetNumberOfITSClusters( nClustersITS );
515         trackSize += flatTrack->GetSize();
516         freeSpace -= flatTrack->GetSize();
517         nTracks++;
518         flatTrack = flatTrack->GetNextTrackNonConst();    
519       }
520       flatEsd->SetTracksEnd( nTracks, trackSize );
521     }       
522
523     if( err ) break;
524
525     // Fill v0's
526   
527     {    
528       size_t v0size = 0;
529       int nV0s = 0; 
530       AliFlatESDV0 *flatV0 = flatEsd->SetV0sStart();
531
532       const AliHLTComponentBlockData* pBlock = GetFirstInputBlock(kAliHLTDataTypeGlobalVertexer|kAliHLTDataOriginOut);
533       if ( pBlock && pBlock->fSize && pBlock->fPtr) {
534         fBenchmark.AddInput(pBlock->fSize);
535         const AliHLTGlobalVertexerComponent::AliHLTGlobalVertexerData *data = reinterpret_cast<AliHLTGlobalVertexerComponent::AliHLTGlobalVertexerData*>(pBlock->fPtr);
536         const int* v0s = data->fTrackIndices + data->fNPrimTracks;
537         for (int i = 0; i < data->fNV0s; ++i) {
538           if(  freeSpace < flatV0->GetSize() ) { err = -1; break; }
539           new (flatV0) AliFlatESDV0;
540           flatV0->SetNegTrackID( v0s[2 * i] );
541           flatV0->SetPosTrackID( v0s[2 * i + 1] );
542           nV0s++;
543           v0size += flatV0->GetSize();
544           freeSpace -= flatV0->GetSize(); 
545           flatV0 = flatV0->GetNextV0NonConst();   
546         }
547       } else {
548         HLTWarning("xxx No V0 data block");
549       }
550
551       flatEsd->SetV0sEnd( nV0s, v0size );
552       cout<<"\nxxxx Found "<<nV0s<<" V0's\n"<<endl;
553     }
554     
555     if( err ) break;
556
557   }while(0);
558
559
560   if( err ){
561     HLTWarning( "Output buffer size %d exceeded, flat ESD event is not stored", maxOutputSize );
562   } else {
563
564     // set up the output block description
565     
566     AliHLTComponentBlockData outBlock;
567     FillBlockData( outBlock );
568     outBlock.fOffset = size;
569     outBlock.fSize = flatEsd->GetSize();
570     outBlock.fDataType = kAliHLTDataTypeFlatESD|kAliHLTDataOriginOut;
571     outputBlocks.push_back( outBlock );
572     fBenchmark.AddOutput(outBlock.fSize);
573     size += outBlock.fSize;
574   }
575
576   
577   // ---------------------------------------------
578   //
579   // Fill the flat ESD friend structure
580   //
581   
582   while( !err && fProduceFriend ){ // single loop for easy break in case of output buffer overflow
583
584     // ---------- Access to clusters --------------------
585
586     const AliHLTTPCClusterData  *partitionClusters[fkNPartition];  //! arrays of cluster data for each TPC partition
587     Int_t                        partitionNClusters[fkNPartition]; //! number of clusters for each TPC partition
588
589     {
590       for(Int_t i=0; i<fkNPartition; i++){
591         partitionClusters[i]  = 0;
592         partitionNClusters[i] = 0;    
593       }
594       for(const AliHLTComponentBlockData *iter = GetFirstInputBlock(AliHLTTPCDefinitions::fgkClustersDataType); iter != NULL; iter = GetNextInputBlock()){      
595         if(iter->fDataType != AliHLTTPCDefinitions::fgkClustersDataType) continue;
596         Int_t slice     = AliHLTTPCDefinitions::GetMinSliceNr(iter->fSpecification);
597         Int_t partition = AliHLTTPCDefinitions::GetMinPatchNr(iter->fSpecification);    
598         Int_t slicepartition = slice*6+partition;      
599         if(slicepartition<0 || slicepartition > fkNPartition){
600           HLTWarning("Wrong header of TPC cluster data, slice %d, partition %d", slice, partition );
601           continue;
602         }
603         const AliHLTTPCClusterData * clusterData = ( AliHLTTPCClusterData* )( iter->fPtr );
604         if( clusterData ) {
605           partitionClusters[slicepartition] = clusterData;
606           partitionNClusters[slicepartition] = clusterData->fSpacePointCnt;
607         }
608       } // end of loop over blocks of clusters    
609     }
610
611
612     AliFlatESDFriend *flatFriend = reinterpret_cast<AliFlatESDFriend*>(outputPtr + size);     
613
614     size_t freeSpaceTotal = maxOutputSize - size;
615     size_t freeSpace = freeSpaceTotal;
616
617     err = ( freeSpace < sizeof( AliFlatESDEvent ) );    
618     if( err ) break;
619
620     new (flatFriend) AliFlatESDFriend;
621       
622     freeSpace = freeSpaceTotal - flatFriend->GetSize();
623   
624     // fill event info
625     {
626       //flatFriend->SetSkipBit( 0 ); // SG!!
627       for( Int_t iSlice=0; iSlice<36; iSlice++ ){
628         int iSector = iSlice;
629         int nclu = 0;
630         for( Int_t iPartition=0; iPartition<3; iPartition++){       
631           int slicepartition = iSlice*6+iPartition;
632           nclu+= partitionNClusters[slicepartition];
633         }
634         flatFriend->SetNclustersTPC( iSector, nclu );
635         iSector = 36+iSlice;
636         nclu = 0;
637         for( Int_t iPartition=3; iPartition<6; iPartition++){       
638           int slicepartition = iSlice*6+iPartition;
639           nclu+= partitionNClusters[slicepartition];
640         }
641         flatFriend->SetNclustersTPC( iSector, nclu );
642         //SG!!!flatFriend->SetNclustersTPCused( iSector, esdFriend->GetNclustersTPCused(iSector) );
643       }
644     }
645       
646     { // Fill track friends information to the flat ESD friend structure
647      
648       size_t trackSize = 0;
649       int nTracks = 0;
650       int nTrackEntries = 0;
651       Long64_t *table = NULL;
652       AliFlatESDFriendTrack *flatTrack = NULL;
653       err = flatFriend->SetTracksStart( flatTrack, table, tracksTPC.size(), freeSpace );
654       if( err ) break;
655       freeSpace = freeSpaceTotal - flatFriend->GetSize();
656
657       for( UInt_t tpcIter=0, itsIter = 0, itsOutIter = 0; tpcIter < tracksTPC.size(); tpcIter++) {       
658         
659         // TPC track parameters
660         
661         AliHLTGlobalBarrelTrack *tpcTrack = &(tracksTPC[tpcIter]);
662         AliHLTGlobalBarrelTrack *tpcOutTrack = &(tracksTPCOut[tpcIter]);
663
664         // ITS track parameters
665         
666         AliHLTGlobalBarrelTrack *itsRefit=0;
667         AliHLTGlobalBarrelTrack *itsOut=0;
668         
669         // ITS Refit track
670           
671         for(; itsIter< tracksITS.size() && tracksITS[itsIter].TrackID()<(int) tpcIter; itsIter++ );
672         
673         if( itsIter< tracksITS.size() && tracksITS[itsIter].TrackID() == (int) tpcIter ){
674           itsRefit = &(tracksITS[itsIter]);
675           itsIter++;
676         }
677         
678         // ITS Out track
679           
680         for(; itsOutIter< tracksITSOut.size() && tracksITSOut[itsOutIter].TrackID()<(int) tpcIter; itsOutIter++ );
681           
682         if( itsOutIter< tracksITSOut.size() && tracksITSOut[itsOutIter].TrackID() == (int) tpcIter ){
683           itsOut = &(tracksITSOut[itsOutIter]);
684           itsOutIter++;
685         }
686         
687         // fill track parameters
688                   
689         table[tpcIter] = trackSize;
690         err = ( freeSpace < flatTrack->EstimateSize() );          
691         if( err ) break;
692         new (flatTrack) AliFlatESDFriendTrack;
693                   
694         freeSpace = freeSpaceTotal - flatFriend->GetSize();
695
696         flatTrack->SetSkipBit( 0 );
697         flatTrack->SetTrackParamTPCOut( tpcOutTrack );
698         flatTrack->SetTrackParamITSOut( itsOut );
699         // flatTrack->SetTrackParamTRDIn( track->GetTRDIn() );
700
701         // fill TPC seed
702
703         AliFlatTPCseed* seed = flatTrack->SetTPCseedStart();
704         new( seed ) AliFlatTPCseed;
705
706         seed->SetLabel( tpcTrack->GetLabel() );
707         seed->SetExternalTrackParam( tpcTrack );
708         
709         // clusters 
710
711         UInt_t nClusters = tpcTrack->GetNumberOfPoints();       
712         const UInt_t*clusterIDs = tpcTrack->GetPoints();
713         for(UInt_t ic=0; ic<nClusters; ic++){    
714           UInt_t id      = clusterIDs[ic];           
715           int iSlice = AliHLTTPCSpacePointData::GetSlice(id);
716           int iPartition = AliHLTTPCSpacePointData::GetPatch(id);
717           int iCluster = AliHLTTPCSpacePointData::GetNumber(id);
718           
719           if(iSlice<0 || iSlice>36 || iPartition<0 || iPartition>5){
720             HLTError("Corrupted TPC cluster Id: slice %d, partition %d, cluster %d", iSlice, iPartition, iCluster);
721             continue;
722           }
723           
724           const AliHLTTPCClusterData * clusterData = partitionClusters[iSlice*6 + iPartition];
725           if(!clusterData ){
726             HLTError("Clusters are missed for slice %d, partition %d", iSlice, iPartition );
727             continue;
728           }
729             
730           if(iCluster >= partitionNClusters[iSlice*6 + iPartition]){
731             HLTError("TPC slice %d, partition %d: ClusterID==%d >= N Cluaters==%d ", iSlice, iPartition,iCluster, partitionNClusters[iSlice*6 + iPartition] );
732             continue;
733           }
734
735             
736           const AliHLTTPCSpacePointData *chlt = &( clusterData->fSpacePoints[iCluster] );
737           AliTPCclusterMI cl;
738           cl.SetX(chlt->fX);
739           cl.SetY(chlt->fY);
740           cl.SetZ(chlt->fZ);
741           cl.SetSigmaY2(chlt->fSigmaY2);
742           cl.SetSigmaYZ( 0 );
743           cl.SetSigmaZ2(chlt->fSigmaZ2);
744           cl.SetQ( chlt->fCharge );
745           cl.SetMax( chlt->fQMax );
746           Int_t sector, row;
747           AliHLTTPCTransform::Slice2Sector(iSlice,chlt->fPadRow, sector, row);
748           cl.SetDetector( sector );
749           cl.SetRow( row );
750           //Float_t padtime[3] = {0,chlt->fY,chlt->fZ};
751           //AliHLTTPCTransform::Local2Raw( padtime, sector, row);
752           //cl.SetPad( (Int_t) padtime[1] );
753           //cl.SetTimeBin( (Int_t) padtime[2] );
754           
755                   
756           tpcTrack->Propagate( TMath::DegToRad()*(sector%18*20.+10.), cl.GetX(), GetBz() );
757           Double_t angle2 = tpcTrack->GetSnp()*tpcTrack->GetSnp();
758           angle2 = (angle2<1) ?TMath::Sqrt(angle2/(1-angle2)) :10.; 
759           AliTPCTrackerPoint point;
760           point.SetAngleY( angle2 );
761           point.SetAngleZ( tpcTrack->GetTgl() );
762
763           seed->AddCluster(&cl, &point ); 
764         } // end of associated cluster loop
765                 
766         
767         flatTrack->SetTPCseedEnd( seed->GetSize() );    
768           
769         trackSize += flatTrack->GetSize();
770         freeSpace -= flatTrack->GetSize();
771         nTrackEntries++;
772         nTracks++;
773         flatTrack = flatTrack->GetNextTrackNonConst();  
774         
775       } // fill tracks
776       
777
778       flatFriend->SetTracksEnd( nTracks, nTrackEntries, trackSize );
779       
780     }      
781
782     if( err ) break;
783
784     { // set up the output block description
785     
786       AliHLTComponentBlockData outBlock;
787       FillBlockData( outBlock );
788       outBlock.fOffset = size;
789       outBlock.fSize = flatFriend->GetSize();
790       outBlock.fDataType = kAliHLTDataTypeFlatESDFriend|kAliHLTDataOriginOut;
791       outputBlocks.push_back( outBlock );
792       fBenchmark.AddOutput(outBlock.fSize);      
793       size += outBlock.fSize;
794     }
795     
796     break;
797   }
798
799   if( err ){
800     HLTWarning( "Output buffer size %d exceeded, flat ESD friend event is not stored", maxOutputSize );
801     return -ENOSPC;
802   } 
803   
804
805   fBenchmark.Stop(0);
806   HLTWarning( fBenchmark.GetStatistics() );
807   if( err ){
808     HLTWarning( "Output buffer size %d exceeded, flat ESD event is not stored", maxOutputSize );
809     return -ENOSPC;
810   } 
811   
812   return 0;
813 }
814