]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/TPCLib/tracking-ca/AliTPCtrackerCA.cxx
Temporary protection if one runs raw->sdigits for the real data.
[u/mrichter/AliRoot.git] / HLT / TPCLib / tracking-ca / AliTPCtrackerCA.cxx
1 // $Id$
2 // **************************************************************************
3 // This file is property of and copyright by the ALICE HLT Project          *
4 // ALICE Experiment at CERN, All rights reserved.                           *
5 //                                                                          *
6 // Primary Authors: Sergey Gorbunov <sergey.gorbunov@kip.uni-heidelberg.de> *
7 //                  Ivan Kisel <kisel@kip.uni-heidelberg.de>                *
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
20 #include "AliTPCtrackerCA.h"
21
22 #include "TTree.h"
23 #include "Riostream.h"
24 //#include "AliCluster.h"
25 #include "AliTPCClustersRow.h"
26 #include "AliTPCParam.h"
27 #include "AliTPCClusterParam.h"
28
29 #include "AliRun.h"
30 #include "AliRunLoader.h"
31 #include "AliStack.h"
32
33 #include "AliHLTTPCCAPerformance.h"
34 #include "AliHLTTPCCAParam.h"
35 #include "AliHLTTPCCATrackConvertor.h"
36 #include "AliHLTTPCCATracker.h"
37 #include "AliHLTTPCCAStandaloneFramework.h"
38
39 #include "TMath.h"
40 #include "AliTPCLoader.h"
41 #include "AliTPC.h"
42 #include "AliTPCclusterMI.h"
43 #include "AliTPCTransform.h"
44 #include "AliTPCcalibDB.h"
45 #include "AliTPCtrack.h"
46 #include "AliTPCseed.h"
47 #include "AliESDtrack.h"
48 #include "AliESDEvent.h"
49 #include "AliTrackReference.h"
50 #include "TStopwatch.h"
51 #include "AliTPCReconstructor.h"
52 #include "AliHLTTPCCAMergerOutput.h"
53
54 //#include <fstream.h>
55
56 ClassImp( AliTPCtrackerCA )
57
58 AliTPCtrackerCA::AliTPCtrackerCA()
59     : AliTracker(), fkParam( 0 ), fClusters( 0 ), fClusterSliceRow( 0 ), fNClusters( 0 ), fDoHLTPerformance( 0 ), fDoHLTPerformanceClusters( 0 ), fStatCPUTime( 0 ), fStatRealTime( 0 ), fStatNEvents( 0 )
60 {
61   //* default constructor
62 }
63
64 AliTPCtrackerCA::AliTPCtrackerCA( const AliTPCtrackerCA & ):
65     AliTracker(), fkParam( 0 ), fClusters( 0 ), fClusterSliceRow( 0 ), fNClusters( 0 ), fDoHLTPerformance( 0 ), fDoHLTPerformanceClusters( 0 ), fStatCPUTime( 0 ), fStatRealTime( 0 ), fStatNEvents( 0 )
66 {
67   //* dummy
68 }
69
70 const AliTPCtrackerCA & AliTPCtrackerCA::operator=( const AliTPCtrackerCA& ) const
71 {
72   //* dummy
73   return *this;
74 }
75
76
77 AliTPCtrackerCA::~AliTPCtrackerCA()
78 {
79   //* destructor
80   delete[] fClusters;
81   delete[] fClusterSliceRow;
82 }
83
84 //#include "AliHLTTPCCADisplay.h"
85
86 AliTPCtrackerCA::AliTPCtrackerCA( const AliTPCParam *par ):
87     AliTracker(), fkParam( par ), fClusters( 0 ), fClusterSliceRow( 0 ), fNClusters( 0 ), fDoHLTPerformance( 0 ), fDoHLTPerformanceClusters( 0 ), fStatCPUTime( 0 ), fStatRealTime( 0 ), fStatNEvents( 0 )
88 {
89   //* constructor
90
91   fDoHLTPerformance = 0;
92   fDoHLTPerformanceClusters = 0;
93
94   AliHLTTPCCAStandaloneFramework &hlt = AliHLTTPCCAStandaloneFramework::Instance();
95
96
97   for ( int iSlice = 0; iSlice < hlt.NSlices(); iSlice++ ) {
98
99     float bz = AliTracker::GetBz();
100
101     float inRmin = fkParam->GetInnerRadiusLow();
102     //float inRmax = fkParam->GetInnerRadiusUp();
103     //float outRmin = fkParam->GetOuterRadiusLow();
104     float outRmax = fkParam->GetOuterRadiusUp();
105     float plusZmin = 0.0529937;
106     float plusZmax = 249.778;
107     float minusZmin = -249.645;
108     float minusZmax = -0.0799937;
109     float dalpha = 0.349066;
110     float alpha = 0.174533 + dalpha * iSlice;
111
112     bool zPlus = ( iSlice < 18 );
113     float zMin =  zPlus ? plusZmin : minusZmin;
114     float zMax =  zPlus ? plusZmax : minusZmax;
115     //TPCZmin = -249.645, ZMax = 249.778
116     //float rMin =  inRmin;
117     //float rMax =  outRmax;
118
119     float padPitch = 0.4;
120     float sigmaZ = 0.228808;
121
122     int nRows = fkParam->GetNRowLow() + fkParam->GetNRowUp();
123     float rowX[200];
124     for ( int irow = 0; irow < fkParam->GetNRowLow(); irow++ ) {
125       rowX[irow] = fkParam->GetPadRowRadiiLow( irow );
126     }
127     for ( int irow = 0; irow < fkParam->GetNRowUp(); irow++ ) {
128       rowX[fkParam->GetNRowLow()+irow] = fkParam->GetPadRowRadiiUp( irow );
129     }
130     AliHLTTPCCAParam param;
131     param.Initialize( iSlice, nRows, rowX, alpha, dalpha,
132                       inRmin, outRmax, zMin, zMax, padPitch, sigmaZ, bz );
133     param.SetHitPickUpFactor( 1. );
134     param.SetMaxTrackMatchDRow( 5 );
135     param.SetTrackConnectionFactor( 3.5 );
136
137     AliTPCClusterParam * clparam = AliTPCcalibDB::Instance()->GetClusterParam();
138     for ( int iRow = 0; iRow < nRows; iRow++ ) {
139       int    type = ( iRow < 63 ) ? 0 : ( ( iRow > 126 ) ? 1 : 2 );
140       for ( int iyz = 0; iyz < 2; iyz++ ) {
141         for ( int k = 0; k < 7; k++ ) {
142           //std::cout<<param.fParamS0Par[iyz][type][k]<<" "<<clparam->fParamS0Par[iyz][type][k] - param.fParamS0Par[iyz][type][k]<<std::endl;
143 #ifndef HAVE_NOT_ALITPCCLUSTERPARAM_r40128
144           param.SetParamS0Par( iyz, type, k, clparam->ParamS0Par(iyz, type, k));
145 #else
146           param.SetParamS0Par( iyz, type, k, clparam->fParamS0Par[iyz][type][k] );
147 #endif //HAVE_NOT_ALITPCCLUSTERPARAM_r40128
148         }
149       }
150     }
151     //hlt.SliceTracker( iSlice ).Initialize( param );
152         hlt.InitializeSliceParam(iSlice, param);
153   }
154 }
155
156
157
158 int AliTPCtrackerCA::LoadClusters ( TTree * fromTree )
159 {
160   // load clusters to the local arrays
161   fNClusters = 0;
162   delete[] fClusters;
163   delete[] fClusterSliceRow;
164
165   AliHLTTPCCAStandaloneFramework &hlt = AliHLTTPCCAStandaloneFramework::Instance();
166
167   if ( fDoHLTPerformance ) {
168     AliHLTTPCCAPerformance::Instance().StartEvent();
169     if ( fDoHLTPerformanceClusters ) AliHLTTPCCAPerformance::Instance().SetDoClusterPulls( 1 );
170   }
171
172   if ( !fkParam ) return 1;
173
174   // load mc tracks
175   while ( fDoHLTPerformance ) {
176     if ( !gAlice ) break;
177 #ifndef HAVE_NOT_ALIRUNLOADER30859
178     AliRunLoader *rl = AliRunLoader::Instance();//gAlice->GetRunLoader();
179 #else
180     // the old way before rev 30859
181     AliRunLoader *rl = AliRunLoader::GetRunLoader();
182 #endif
183     if ( !rl ) break;
184     rl->LoadKinematics();
185     AliStack *stack = rl->Stack();
186     if ( !stack ) break;
187
188     AliHLTTPCCAPerformance::Instance().SetNMCTracks( stack->GetNtrack() );
189
190     for ( int itr = 0; itr < stack->GetNtrack(); itr++ ) {
191       TParticle *part = stack->Particle( itr );
192       AliHLTTPCCAPerformance::Instance().ReadMCTrack( itr, part );
193     }
194
195     { // check for MC tracks at the TPC entrance
196
197       bool *isTPC = 0;
198       isTPC = new bool [stack->GetNtrack()];
199       for ( int i = 0; i < stack->GetNtrack(); i++ ) isTPC[i] = 0;
200       rl->LoadTrackRefs();
201       TTree *mcTree = rl->TreeTR();
202       if ( !mcTree ) break;
203       TBranch *branch = mcTree->GetBranch( "TrackReferences" );
204       if ( !branch ) break;
205       TClonesArray tpcdummy( "AliTrackReference", 1000 ), *tpcRefs = &tpcdummy;
206       branch->SetAddress( &tpcRefs );
207       int nr = ( int )mcTree->GetEntries();
208       for ( int r = 0; r < nr; r++ ) {
209         mcTree->GetEvent( r );
210         int nref = tpcRefs->GetEntriesFast();
211         if ( !nref ) continue;
212         AliTrackReference *tpcRef = 0x0;
213         for ( int iref = 0; iref < nref; ++iref ) {
214           tpcRef = ( AliTrackReference* )tpcRefs->UncheckedAt( iref );
215           if ( tpcRef->DetectorId() == AliTrackReference::kTPC ) break;
216           tpcRef = 0x0;
217         }
218         if ( !tpcRef ) continue;
219
220         if ( isTPC[tpcRef->Label()] ) continue;
221
222         AliHLTTPCCAPerformance::Instance().ReadMCTPCTrack( tpcRef->Label(),
223             tpcRef->X(), tpcRef->Y(), tpcRef->Z(),
224             tpcRef->Px(), tpcRef->Py(), tpcRef->Pz() );
225         isTPC[tpcRef->Label()] = 1;
226         tpcRefs->Clear();
227       }
228       if ( isTPC ) delete[] isTPC;
229     }
230
231     while ( fDoHLTPerformanceClusters ) {
232       AliTPCLoader *tpcl = ( AliTPCLoader* ) rl->GetDetectorLoader( "TPC" );
233       if ( !tpcl ) break;
234       if ( tpcl->TreeH() == 0x0 ) {
235         if ( tpcl->LoadHits() ) break;
236       }
237       if ( tpcl->TreeH() == 0x0 ) break;
238
239       AliTPC *tpc = ( AliTPC* ) gAlice->GetDetector( "TPC" );
240       int nEnt = ( int )tpcl->TreeH()->GetEntries();
241       int nPoints = 0;
242       for ( int iEnt = 0; iEnt < nEnt; iEnt++ ) {
243         tpc->ResetHits();
244         tpcl->TreeH()->GetEvent( iEnt );
245         AliTPChit *phit = ( AliTPChit* )tpc->FirstHit( -1 );
246         for ( ; phit; phit = ( AliTPChit* )tpc->NextHit() ) nPoints++;
247       }
248       AliHLTTPCCAPerformance::Instance().SetNMCPoints( nPoints );
249
250       for ( int iEnt = 0; iEnt < nEnt; iEnt++ ) {
251         tpc->ResetHits();
252         tpcl->TreeH()->GetEvent( iEnt );
253         AliTPChit *phit = ( AliTPChit* )tpc->FirstHit( -1 );
254         for ( ; phit; phit = ( AliTPChit* )tpc->NextHit() ) {
255           AliHLTTPCCAPerformance::Instance().ReadMCPoint( phit->GetTrack(), phit->X(), phit->Y(), phit->Z(), phit->Time(), phit->fSector % 36 );
256         }
257       }
258       break;
259     }
260     break;
261   }
262
263   TBranch * br = fromTree->GetBranch( "Segment" );
264   if ( !br ) return 1;
265
266   AliTPCClustersRow *clrow = new AliTPCClustersRow;
267   clrow->SetClass( "AliTPCclusterMI" );
268   clrow->SetArray( 0 );
269   clrow->GetArray()->ExpandCreateFast( 10000 );
270
271   br->SetAddress( &clrow );
272
273   //
274   int nEnt = int( fromTree->GetEntries() );
275
276   fNClusters = 0;
277   for ( int i = 0; i < nEnt; i++ ) {
278     br->GetEntry( i );
279     int sec, row;
280     fkParam->AdjustSectorRow( clrow->GetID(), sec, row );
281     fNClusters += clrow->GetArray()->GetEntriesFast();
282   }
283
284   fClusters = new AliTPCclusterMI [fNClusters];
285   fClusterSliceRow = new unsigned int [fNClusters];
286
287   hlt.StartDataReading( fNClusters );
288
289   if ( fDoHLTPerformance ) AliHLTTPCCAPerformance::Instance().SetNHits( fNClusters );
290
291   int ind = 0;
292   for ( int i = 0; i < nEnt; i++ ) {
293     br->GetEntry( i );
294     int sec, row;
295     fkParam->AdjustSectorRow( clrow->GetID(), sec, row );
296     int nClu = clrow->GetArray()->GetEntriesFast();
297     float x = fkParam->GetPadRowRadii( sec, row );
298
299     if ( sec >= 36 ) {
300       sec = sec - 36;
301       row = row + fkParam->GetNRowLow();
302     }
303
304     unsigned int sliceRow = ( sec << 8 ) + row;
305
306     for ( int icl = 0; icl < nClu; icl++ ) {
307       int lab0 = -1;
308       int lab1 = -1;
309       int lab2 = -1;
310       AliTPCclusterMI* cluster = ( AliTPCclusterMI* )( clrow->GetArray()->At( icl ) );
311       if ( !cluster ) continue;
312       lab0 = cluster->GetLabel( 0 );
313       lab1 = cluster->GetLabel( 1 );
314       lab2 = cluster->GetLabel( 2 );
315
316       AliTPCTransform *transform = AliTPCcalibDB::Instance()->GetTransform() ;
317       transform->SetCurrentRecoParam((AliTPCRecoParam*)AliTPCReconstructor::GetRecoParam());
318
319       if ( !transform ) {
320         AliFatal( "Tranformations not in calibDB" );
321       }
322       double xx[3] = {cluster->GetRow(), cluster->GetPad(), cluster->GetTimeBin()};
323       int id[1] = {cluster->GetDetector()};
324       transform->Transform( xx, id, 0, 1 );
325
326       cluster->SetX( xx[0] );
327       cluster->SetY( xx[1] );
328       cluster->SetZ( xx[2] );
329
330       TGeoHMatrix  *mat = fkParam->GetClusterMatrix( cluster->GetDetector() );
331       double pos[3] = {cluster->GetX(), cluster->GetY(), cluster->GetZ()};
332       double posC[3] = {cluster->GetX(), cluster->GetY(), cluster->GetZ()};
333       if ( mat ) mat->LocalToMaster( pos, posC );
334       else {
335         // chack Loading of Geo matrices from GeoManager - TEMPORARY FIX
336       }
337       cluster->SetX( posC[0] );
338       cluster->SetY( posC[1] );
339       cluster->SetZ( posC[2] );
340
341       x = cluster->GetX();
342       float y = cluster->GetY();
343       float z = cluster->GetZ();
344
345
346       int index = ind++;
347       fClusters[index] = *cluster;
348
349       fClusterSliceRow[index] = sliceRow;
350
351       hlt.ReadCluster( index, sec, row, x, y, z, cluster->GetQ() );
352
353       if ( fDoHLTPerformance ) AliHLTTPCCAPerformance::Instance().ReadHitLabel( index, lab0, lab1, lab2 );
354     }
355   }
356   delete clrow;
357
358   hlt.FinishDataReading();
359
360   //AliHLTTPCCAPerformance::Instance().SmearClustersMC();
361
362   return 0;
363 }
364
365 AliCluster * AliTPCtrackerCA::GetCluster( int index ) const
366 {
367   return &( fClusters[index] );
368 }
369
370 int AliTPCtrackerCA::Clusters2Tracks( AliESDEvent *event )
371 {
372   // reconstruction
373   //cout<<"Start of AliTPCtrackerCA"<<endl;
374   TStopwatch timer;
375
376   AliHLTTPCCAStandaloneFramework &hlt = AliHLTTPCCAStandaloneFramework::Instance();
377
378   hlt.ProcessEvent();
379
380   if ( event ) {
381
382     const AliHLTTPCCAMergerOutput &hltOut = *( hlt.Merger().Output() );
383
384     for ( int itr = 0; itr < hltOut.NTracks(); itr++ ) {
385
386       AliTPCseed tTPC;
387
388       const AliHLTTPCCAMergedTrack &tCA = hltOut.Track( itr );
389
390       AliHLTTPCCATrackParam par = tCA.InnerParam();
391       AliHLTTPCCATrackConvertor::GetExtParam( par, tTPC, tCA.InnerAlpha() );
392       tTPC.SetMass( 0.13957 );
393       if ( TMath::Abs( tTPC.GetSigned1Pt() ) > 1. / 0.02 ) continue;
394       int nhits = tCA.NClusters();
395       int firstHit = 0;
396       if ( nhits > 160 ) {
397         firstHit = nhits - 160;
398         nhits = 160;
399       }
400       tTPC.SetNumberOfClusters( nhits );
401       float alpha = tCA.InnerAlpha();
402       AliHLTTPCCATrackParam t0 = par;
403       for ( int ih = 0; ih < nhits; ih++ ) {
404         int index = hltOut.ClusterId( tCA.FirstClusterRef() + firstHit + ih );
405         tTPC.SetClusterIndex( ih, index );
406         AliTPCclusterMI *c = &( fClusters[index] );
407         int iSlice = fClusterSliceRow[index] >> 8;
408         int row = fClusterSliceRow[index] & 0xff;
409
410         tTPC.SetClusterPointer( row, c );
411         AliTPCTrackerPoint &point = *( tTPC.GetTrackPoint( row ) );
412         {
413           //AliHLTTPCCATracker &slice = hlt.SliceTracker( iSlice );
414           if ( hlt.Param(iSlice).Alpha() != alpha ) {
415             if ( ! t0.Rotate(  hlt.Param(iSlice).Alpha() - alpha, .999 ) ) continue;
416             alpha = hlt.Param(iSlice).Alpha();
417           }
418           float x = hlt.Row(iSlice, row).X();
419           if ( !t0.TransportToX( x, hlt.Param(iSlice).GetBz( t0 ), .999 ) ) continue;
420           float sy2, sz2;
421           //slice.GetErrors2( row, t0, sy2, sz2 );
422                   hlt.Param(iSlice).GetClusterErrors2( row, t0.GetZ(), t0.SinPhi(), t0.GetCosPhi(), t0.DzDs(), sy2, sz2 );
423           point.SetSigmaY( c->GetSigmaY2() / sy2 );
424           point.SetSigmaZ( c->GetSigmaZ2() / sz2 );
425           point.SetAngleY( TMath::Abs( t0.GetSinPhi() / t0.GetCosPhi() ) );
426           point.SetAngleZ( TMath::Abs( t0.GetDzDs() ) );
427         }
428       }
429       tTPC.CookdEdx( 0.02, 0.6 );
430
431       CookLabel( &tTPC, 0.1 );
432
433       if ( 1 ) { // correction like in off-line --- Adding systematic error
434
435         const double *param = AliTPCReconstructor::GetRecoParam()->GetSystematicError();
436         double covar[15];
437         for ( int i = 0; i < 15; i++ ) covar[i] = 0;
438         covar[0] = param[0] * param[0];
439         covar[2] = param[1] * param[1];
440         covar[5] = param[2] * param[2];
441         covar[9] = param[3] * param[3];
442         double facC =  AliTracker::GetBz() * kB2C;
443         covar[14] = param[4] * param[4] * facC * facC;
444         tTPC.AddCovariance( covar );
445       }
446
447       AliESDtrack tESD;
448       tESD.UpdateTrackParams( &( tTPC ), AliESDtrack::kTPCin );
449       //tESD.SetStatus( AliESDtrack::kTPCrefit );
450       //tESD.SetTPCPoints(tTPC.GetPoints());
451       int   ndedx = tTPC.GetNCDEDX( 0 );
452       float sdedx = tTPC.GetSDEDX( 0 );
453       float dedx  = tTPC.GetdEdx();
454       tESD.SetTPCsignal( dedx, sdedx, ndedx );
455       //tESD.myTPC = tTPC;
456
457       event->AddTrack( &tESD );
458     }
459   }
460
461
462   timer.Stop();
463
464   fStatCPUTime += timer.CpuTime();
465   fStatRealTime += timer.RealTime();
466
467   //cout << "\n\nCA tracker speed: cpu = " << fStatCPUTime / ( fStatNEvents + 1 )*1.e3 << " [ms/ev], real = " << fStatRealTime / ( fStatNEvents + 1 )*1.e3 << " [ms/ev], n calls = " << ( fStatNEvents + 1 ) << endl << endl;
468
469
470   //cout<<"Do performance.."<<endl;
471
472   if ( fDoHLTPerformance ) AliHLTTPCCAPerformance::Instance().Performance();
473
474   if ( 0 ) {// Write Event
475     if ( fStatNEvents == 0 ) {
476       fstream geo;
477       geo.open( "CAEvents/settings.dat", ios::out );
478       if ( geo.is_open() ) {
479         hlt.WriteSettings( geo );
480       }
481       geo.close();
482     }
483
484     fstream hits;
485     char name[255];
486     sprintf( name, "CAEvents/%i.event.dat", fStatNEvents );
487     hits.open( name, ios::out );
488     if ( hits.is_open() ) {
489       hlt.WriteEvent( hits );
490       fstream tracks;
491       sprintf( name, "CAEvents/%i.tracks.dat", fStatNEvents );
492       tracks.open( name, ios::out );
493       hlt.WriteTracks( tracks );
494     }
495     hits.close();
496     if ( fDoHLTPerformance ) {
497       fstream mcevent, mcpoints;
498       char mcname[255];
499       sprintf( mcname, "CAEvents/%i.mcevent.dat", fStatNEvents );
500       mcevent.open( mcname, ios::out );
501       if ( mcevent.is_open() ) {
502         AliHLTTPCCAPerformance::Instance().WriteMCEvent( mcevent );
503       }
504       if ( 1 && fDoHLTPerformanceClusters ) {
505         sprintf( mcname, "CAEvents/%i.mcpoints.dat", fStatNEvents );
506         mcpoints.open( mcname, ios::out );
507         if ( mcpoints.is_open() ) {
508           AliHLTTPCCAPerformance::Instance().WriteMCPoints( mcpoints );
509         }
510         mcpoints.close();
511       }
512       mcevent.close();
513     }
514   }
515
516   fStatNEvents++;
517
518
519   //cout<<"End of AliTPCtrackerCA"<<endl;
520   return 0;
521 }
522
523
524 int AliTPCtrackerCA::RefitInward ( AliESDEvent *event )
525 {
526   //* forward propagation of ESD tracks
527
528   float xTPC = fkParam->GetInnerRadiusLow();
529   float dAlpha = fkParam->GetInnerAngle() / 180.*TMath::Pi();
530   float yMax = xTPC * TMath::Tan( dAlpha / 2. );
531
532   AliHLTTPCCAStandaloneFramework &hlt = AliHLTTPCCAStandaloneFramework::Instance();
533
534   int nentr = event->GetNumberOfTracks();
535
536   for ( int itr = 0; itr < nentr; itr++ ) {
537     AliESDtrack *esd = event->GetTrack( itr );
538     ULong_t status = esd->GetStatus();
539     if ( !( status&AliESDtrack::kTPCin ) ) continue;
540     AliHLTTPCCATrackParam t0;
541     AliHLTTPCCATrackConvertor::SetExtParam( t0, *esd );
542     AliHLTTPCCATrackParam t = t0;
543     float alpha = esd->GetAlpha();
544     //float dEdX=0;
545     AliHLTTPCCAMerger::AliHLTTPCCAClusterInfo infos[500];
546     int hits[500], hits1[500];
547     int nHits = esd->GetTPCclusters( hits );
548     for ( int i = 0; i < nHits; i++ ) {
549       hits1[i] = i;
550       int index = hits[i];
551       infos[i].SetISlice( fClusterSliceRow[index] >> 8 );
552       infos[i].SetIRow( fClusterSliceRow[index] & 0xff );
553       infos[i].SetId( index );
554       infos[i].SetX( fClusters[index].GetX() );
555       infos[i].SetY( fClusters[index].GetY() );
556       infos[i].SetZ( fClusters[index].GetZ() );
557     }
558
559     bool ok = hlt.Merger().FitTrack( t, alpha, t0, alpha, hits1, nHits, 0, 1,infos );
560
561     if ( ok &&  nHits > 15 ) {
562       if ( t.TransportToXWithMaterial( xTPC, hlt.Merger().SliceParam().GetBz( t ) ) ) {
563         if ( t.GetY() > yMax ) {
564           if ( t.Rotate( dAlpha ) ) {
565             alpha += dAlpha;
566             t.TransportToXWithMaterial( xTPC, hlt.Merger().SliceParam().GetBz( t ) );
567           }
568         } else if ( t.GetY() < -yMax ) {
569           if ( t.Rotate( -dAlpha ) ) {
570             alpha += -dAlpha;
571             t.TransportToXWithMaterial( xTPC, hlt.Merger().SliceParam().GetBz( t ) );
572           }
573         }
574       }
575
576       AliTPCtrack tt( *esd );
577       if ( AliHLTTPCCATrackConvertor::GetExtParam( t, tt, alpha ) ) {
578         if ( t.X() > 50 ) esd->UpdateTrackParams( &tt, AliESDtrack::kTPCrefit );
579       }
580     }
581   }
582   return 0;
583 }
584
585 int AliTPCtrackerCA::PropagateBack( AliESDEvent *event )
586 {
587
588   //* backward propagation of ESD tracks
589
590   AliHLTTPCCAStandaloneFramework &hlt = AliHLTTPCCAStandaloneFramework::Instance();
591
592   int nentr = event->GetNumberOfTracks();
593
594   for ( int itr = 0; itr < nentr; itr++ ) {
595
596     AliESDtrack *esd = event->GetTrack( itr );
597     ULong_t status = esd->GetStatus();
598     if ( !( status&AliESDtrack::kTPCin ) ) continue;
599
600     AliHLTTPCCATrackParam t0;
601     AliHLTTPCCATrackConvertor::SetExtParam( t0, *esd  );
602     AliHLTTPCCATrackParam t = t0;
603     float alpha = esd->GetAlpha();
604     //float dEdX=0;
605     AliHLTTPCCAMerger::AliHLTTPCCAClusterInfo infos[500];
606     int hits[500], hits1[500];
607     int nHits = esd->GetTPCclusters( hits );
608     for ( int i = 0; i < nHits; i++ ) {
609       hits1[i] = i;
610       int index = hits[i];
611       infos[i].SetISlice( fClusterSliceRow[index] >> 8 );
612       infos[i].SetIRow( fClusterSliceRow[index] & 0xff );
613       infos[i].SetId( index );
614       infos[i].SetX( fClusters[index].GetX() );
615       infos[i].SetY( fClusters[index].GetY() );
616       infos[i].SetZ( fClusters[index].GetZ() );
617     }
618
619     bool ok = hlt.Merger().FitTrack( t, alpha, t0, alpha, hits1, nHits, 1, 1, infos );
620
621     if ( ok &&  nHits > 15 ) {
622       AliTPCtrack tt( *esd );
623       if ( AliHLTTPCCATrackConvertor::GetExtParam( t, tt, alpha ) ) {
624         if ( t.X() > 50 ) esd->UpdateTrackParams( &tt, AliESDtrack::kTPCout );
625       }
626     }
627   }
628   return 0;
629 }