]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/TPCLib/tracking-ca/AliTPCtrackerCA.cxx
Making the TPCLib compatible with the trunk of Jan 2009
[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     const double kCLight = 0.000299792458;
100
101     float bz = AliTracker::GetBz() * kCLight;
102
103     float inRmin = fkParam->GetInnerRadiusLow();
104     //float inRmax = fkParam->GetInnerRadiusUp();
105     //float outRmin = fkParam->GetOuterRadiusLow();
106     float outRmax = fkParam->GetOuterRadiusUp();
107     float plusZmin = 0.0529937;
108     float plusZmax = 249.778;
109     float minusZmin = -249.645;
110     float minusZmax = -0.0799937;
111     float dalpha = 0.349066;
112     float alpha = 0.174533 + dalpha * iSlice;
113
114     bool zPlus = ( iSlice < 18 );
115     float zMin =  zPlus ? plusZmin : minusZmin;
116     float zMax =  zPlus ? plusZmax : minusZmax;
117     //TPCZmin = -249.645, ZMax = 249.778
118     //float rMin =  inRmin;
119     //float rMax =  outRmax;
120
121     float padPitch = 0.4;
122     float sigmaZ = 0.228808;
123
124     int nRows = fkParam->GetNRowLow() + fkParam->GetNRowUp();
125     float rowX[200];
126     for ( int irow = 0; irow < fkParam->GetNRowLow(); irow++ ) {
127       rowX[irow] = fkParam->GetPadRowRadiiLow( irow );
128     }
129     for ( int irow = 0; irow < fkParam->GetNRowUp(); irow++ ) {
130       rowX[fkParam->GetNRowLow()+irow] = fkParam->GetPadRowRadiiUp( irow );
131     }
132     AliHLTTPCCAParam param;
133     param.Initialize( iSlice, nRows, rowX, alpha, dalpha,
134                       inRmin, outRmax, zMin, zMax, padPitch, sigmaZ, bz );
135     param.SetHitPickUpFactor( 1. );
136     param.SetMaxTrackMatchDRow( 5 );
137     param.SetTrackConnectionFactor( 3.5 );
138
139     AliTPCClusterParam * clparam = AliTPCcalibDB::Instance()->GetClusterParam();
140     for ( int iRow = 0; iRow < nRows; iRow++ ) {
141       int    type = ( iRow < 63 ) ? 0 : ( ( iRow > 126 ) ? 1 : 2 );
142       for ( int iyz = 0; iyz < 2; iyz++ ) {
143         for ( int k = 0; k < 7; k++ ) {
144           //std::cout<<param.fParamS0Par[iyz][type][k]<<" "<<clparam->fParamS0Par[iyz][type][k] - param.fParamS0Par[iyz][type][k]<<std::endl;
145           param.SetParamS0Par( iyz, type, k, clparam->fParamS0Par[iyz][type][k] );
146         }
147       }
148     }
149     hlt.SliceTracker( iSlice ).Initialize( param );
150   }
151 }
152
153
154
155 int AliTPCtrackerCA::LoadClusters ( TTree * fromTree )
156 {
157   // load clusters to the local arrays
158   fNClusters = 0;
159   delete[] fClusters;
160   delete[] fClusterSliceRow;
161
162   AliHLTTPCCAStandaloneFramework &hlt = AliHLTTPCCAStandaloneFramework::Instance();
163
164   if ( fDoHLTPerformance ) {
165     AliHLTTPCCAPerformance::Instance().StartEvent();
166     if ( fDoHLTPerformanceClusters ) AliHLTTPCCAPerformance::Instance().SetDoClusterPulls( 1 );
167   }
168
169   if ( !fkParam ) return 1;
170
171   // load mc tracks
172   while ( fDoHLTPerformance ) {
173     if ( !gAlice ) break;
174 #ifndef HAVE_NOT_ALIRUNLOADER30859
175     AliRunLoader *rl = AliRunLoader::Instance();//gAlice->GetRunLoader();
176 #else
177     // the old way before rev 30859
178     AliRunLoader *rl = AliRunLoader::GetRunLoader();
179 #endif
180     if ( !rl ) break;
181     rl->LoadKinematics();
182     AliStack *stack = rl->Stack();
183     if ( !stack ) break;
184
185     AliHLTTPCCAPerformance::Instance().SetNMCTracks( stack->GetNtrack() );
186
187     for ( int itr = 0; itr < stack->GetNtrack(); itr++ ) {
188       TParticle *part = stack->Particle( itr );
189       AliHLTTPCCAPerformance::Instance().ReadMCTrack( itr, part );
190     }
191
192     { // check for MC tracks at the TPC entrance
193
194       bool *isTPC = 0;
195       isTPC = new bool [stack->GetNtrack()];
196       for ( int i = 0; i < stack->GetNtrack(); i++ ) isTPC[i] = 0;
197       rl->LoadTrackRefs();
198       TTree *mcTree = rl->TreeTR();
199       if ( !mcTree ) break;
200       TBranch *branch = mcTree->GetBranch( "TrackReferences" );
201       if ( !branch ) break;
202       TClonesArray tpcdummy( "AliTrackReference", 1000 ), *tpcRefs = &tpcdummy;
203       branch->SetAddress( &tpcRefs );
204       int nr = ( int )mcTree->GetEntries();
205       for ( int r = 0; r < nr; r++ ) {
206         mcTree->GetEvent( r );
207         int nref = tpcRefs->GetEntriesFast();
208         if ( !nref ) continue;
209         AliTrackReference *tpcRef = 0x0;
210         for ( int iref = 0; iref < nref; ++iref ) {
211           tpcRef = ( AliTrackReference* )tpcRefs->UncheckedAt( iref );
212           if ( tpcRef->DetectorId() == AliTrackReference::kTPC ) break;
213           tpcRef = 0x0;
214         }
215         if ( !tpcRef ) continue;
216
217         if ( isTPC[tpcRef->Label()] ) continue;
218
219         AliHLTTPCCAPerformance::Instance().ReadMCTPCTrack( tpcRef->Label(),
220             tpcRef->X(), tpcRef->Y(), tpcRef->Z(),
221             tpcRef->Px(), tpcRef->Py(), tpcRef->Pz() );
222         isTPC[tpcRef->Label()] = 1;
223         tpcRefs->Clear();
224       }
225       if ( isTPC ) delete[] isTPC;
226     }
227
228     while ( fDoHLTPerformanceClusters ) {
229       AliTPCLoader *tpcl = ( AliTPCLoader* ) rl->GetDetectorLoader( "TPC" );
230       if ( !tpcl ) break;
231       if ( tpcl->TreeH() == 0x0 ) {
232         if ( tpcl->LoadHits() ) break;
233       }
234       if ( tpcl->TreeH() == 0x0 ) break;
235
236       AliTPC *tpc = ( AliTPC* ) gAlice->GetDetector( "TPC" );
237       int nEnt = ( int )tpcl->TreeH()->GetEntries();
238       int nPoints = 0;
239       for ( int iEnt = 0; iEnt < nEnt; iEnt++ ) {
240         tpc->ResetHits();
241         tpcl->TreeH()->GetEvent( iEnt );
242         AliTPChit *phit = ( AliTPChit* )tpc->FirstHit( -1 );
243         for ( ; phit; phit = ( AliTPChit* )tpc->NextHit() ) nPoints++;
244       }
245       AliHLTTPCCAPerformance::Instance().SetNMCPoints( nPoints );
246
247       for ( int iEnt = 0; iEnt < nEnt; iEnt++ ) {
248         tpc->ResetHits();
249         tpcl->TreeH()->GetEvent( iEnt );
250         AliTPChit *phit = ( AliTPChit* )tpc->FirstHit( -1 );
251         for ( ; phit; phit = ( AliTPChit* )tpc->NextHit() ) {
252           AliHLTTPCCAPerformance::Instance().ReadMCPoint( phit->GetTrack(), phit->X(), phit->Y(), phit->Z(), phit->Time(), phit->fSector % 36 );
253         }
254       }
255       break;
256     }
257     break;
258   }
259
260   TBranch * br = fromTree->GetBranch( "Segment" );
261   if ( !br ) return 1;
262
263   AliTPCClustersRow *clrow = new AliTPCClustersRow;
264   clrow->SetClass( "AliTPCclusterMI" );
265   clrow->SetArray( 0 );
266   clrow->GetArray()->ExpandCreateFast( 10000 );
267
268   br->SetAddress( &clrow );
269
270   //
271   int nEnt = int( fromTree->GetEntries() );
272
273   fNClusters = 0;
274   for ( int i = 0; i < nEnt; i++ ) {
275     br->GetEntry( i );
276     int sec, row;
277     fkParam->AdjustSectorRow( clrow->GetID(), sec, row );
278     fNClusters += clrow->GetArray()->GetEntriesFast();
279   }
280
281   fClusters = new AliTPCclusterMI [fNClusters];
282   fClusterSliceRow = new unsigned int [fNClusters];
283
284   hlt.StartDataReading( fNClusters );
285
286   if ( fDoHLTPerformance ) AliHLTTPCCAPerformance::Instance().SetNHits( fNClusters );
287
288   int ind = 0;
289   for ( int i = 0; i < nEnt; i++ ) {
290     br->GetEntry( i );
291     int sec, row;
292     fkParam->AdjustSectorRow( clrow->GetID(), sec, row );
293     int nClu = clrow->GetArray()->GetEntriesFast();
294     float x = fkParam->GetPadRowRadii( sec, row );
295
296     if ( sec >= 36 ) {
297       sec = sec - 36;
298       row = row + fkParam->GetNRowLow();
299     }
300
301     unsigned int sliceRow = ( sec << 8 ) + row;
302
303     for ( int icl = 0; icl < nClu; icl++ ) {
304       int lab0 = -1;
305       int lab1 = -1;
306       int lab2 = -1;
307       AliTPCclusterMI* cluster = ( AliTPCclusterMI* )( clrow->GetArray()->At( icl ) );
308       if ( !cluster ) continue;
309       lab0 = cluster->GetLabel( 0 );
310       lab1 = cluster->GetLabel( 1 );
311       lab2 = cluster->GetLabel( 2 );
312
313       AliTPCTransform *transform = AliTPCcalibDB::Instance()->GetTransform() ;
314       if ( !transform ) {
315         AliFatal( "Tranformations not in calibDB" );
316       }
317       double xx[3] = {cluster->GetRow(), cluster->GetPad(), cluster->GetTimeBin()};
318       int id[1] = {cluster->GetDetector()};
319       transform->Transform( xx, id, 0, 1 );
320       //if (!AliTPCReconstructor::GetRecoParam()->GetBYMirror()){
321       //if (cluster->GetDetector()%36>17){
322       //xx[1]*=-1;
323       //}
324       //}
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 ( slice.Param().Alpha() != alpha ) {
415             if ( ! t0.Rotate(  slice.Param().Alpha() - alpha, .999 ) ) continue;
416             alpha = slice.Param().Alpha();
417           }
418           float x = slice.Row( row ).X();
419           if ( !t0.TransportToX( x, slice.Param().GetBz( t0 ), .999 ) ) continue;
420           float sy2, sz2;
421           slice.GetErrors2( row, t0, sy2, sz2 );
422           point.SetSigmaY( c->GetSigmaY2() / sy2 );
423           point.SetSigmaZ( c->GetSigmaZ2() / sz2 );
424           point.SetAngleY( TMath::Abs( t0.GetSinPhi() / t0.GetCosPhi() ) );
425           point.SetAngleZ( TMath::Abs( t0.GetDzDs() ) );
426         }
427       }
428       tTPC.CookdEdx( 0.02, 0.6 );
429
430       CookLabel( &tTPC, 0.1 );
431
432       if ( 1 ) { // correction like in off-line --- Adding systematic error
433
434         const double *param = AliTPCReconstructor::GetRecoParam()->GetSystematicError();
435         double covar[15];
436         for ( int i = 0; i < 15; i++ ) covar[i] = 0;
437         covar[0] = param[0] * param[0];
438         covar[2] = param[1] * param[1];
439         covar[5] = param[2] * param[2];
440         covar[9] = param[3] * param[3];
441         double facC =  AliTracker::GetBz() * kB2C;
442         covar[14] = param[4] * param[4] * facC * facC;
443         tTPC.AddCovariance( covar );
444       }
445
446       AliESDtrack tESD;
447       tESD.UpdateTrackParams( &( tTPC ), AliESDtrack::kTPCin );
448       //tESD.SetStatus( AliESDtrack::kTPCrefit );
449       //tESD.SetTPCPoints(tTPC.GetPoints());
450       int   ndedx = tTPC.GetNCDEDX( 0 );
451       float sdedx = tTPC.GetSDEDX( 0 );
452       float dedx  = tTPC.GetdEdx();
453       tESD.SetTPCsignal( dedx, sdedx, ndedx );
454       //tESD.myTPC = tTPC;
455
456       event->AddTrack( &tESD );
457     }
458   }
459
460
461   timer.Stop();
462
463   fStatCPUTime += timer.CpuTime();
464   fStatRealTime += timer.RealTime();
465
466   //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;
467
468
469   //cout<<"Do performance.."<<endl;
470
471   if ( fDoHLTPerformance ) AliHLTTPCCAPerformance::Instance().Performance();
472
473   if ( 0 ) {// Write Event
474     if ( fStatNEvents == 0 ) {
475       fstream geo;
476       geo.open( "CAEvents/settings.dat", ios::out );
477       if ( geo.is_open() ) {
478         hlt.WriteSettings( geo );
479       }
480       geo.close();
481     }
482
483     fstream hits;
484     char name[255];
485     sprintf( name, "CAEvents/%i.event.dat", fStatNEvents );
486     hits.open( name, ios::out );
487     if ( hits.is_open() ) {
488       hlt.WriteEvent( hits );
489       fstream tracks;
490       sprintf( name, "CAEvents/%i.tracks.dat", fStatNEvents );
491       tracks.open( name, ios::out );
492       hlt.WriteTracks( tracks );
493     }
494     hits.close();
495     if ( fDoHLTPerformance ) {
496       fstream mcevent, mcpoints;
497       char mcname[255];
498       sprintf( mcname, "CAEvents/%i.mcevent.dat", fStatNEvents );
499       mcevent.open( mcname, ios::out );
500       if ( mcevent.is_open() ) {
501         AliHLTTPCCAPerformance::Instance().WriteMCEvent( mcevent );
502       }
503       if ( 1 && fDoHLTPerformanceClusters ) {
504         sprintf( mcname, "CAEvents/%i.mcpoints.dat", fStatNEvents );
505         mcpoints.open( mcname, ios::out );
506         if ( mcpoints.is_open() ) {
507           AliHLTTPCCAPerformance::Instance().WriteMCPoints( mcpoints );
508         }
509         mcpoints.close();
510       }
511       mcevent.close();
512     }
513   }
514
515   fStatNEvents++;
516
517
518   //cout<<"End of AliTPCtrackerCA"<<endl;
519   return 0;
520 }
521
522
523 int AliTPCtrackerCA::RefitInward ( AliESDEvent *event )
524 {
525   //* forward propagation of ESD tracks
526
527   float xTPC = fkParam->GetInnerRadiusLow();
528   float dAlpha = fkParam->GetInnerAngle() / 180.*TMath::Pi();
529   float yMax = xTPC * TMath::Tan( dAlpha / 2. );
530
531   AliHLTTPCCAStandaloneFramework &hlt = AliHLTTPCCAStandaloneFramework::Instance();
532
533   int nentr = event->GetNumberOfTracks();
534
535   for ( int itr = 0; itr < nentr; itr++ ) {
536     AliESDtrack *esd = event->GetTrack( itr );
537     ULong_t status = esd->GetStatus();
538     if ( !( status&AliESDtrack::kTPCin ) ) continue;
539     AliHLTTPCCATrackParam t0;
540     AliHLTTPCCATrackConvertor::SetExtParam( t0, *esd );
541     AliHLTTPCCATrackParam t = t0;
542     float alpha = esd->GetAlpha();
543     //float dEdX=0;
544     AliHLTTPCCAMerger::AliHLTTPCCAClusterInfo infos[500];
545     int hits[500], hits1[500];
546     int nHits = esd->GetTPCclusters( hits );
547     for ( int i = 0; i < nHits; i++ ) {
548       hits1[i] = i;
549       int index = hits[i];
550       infos[i].SetISlice( fClusterSliceRow[index] >> 8 );
551       infos[i].SetX( fClusters[index].GetX() );
552       infos[i].SetY( fClusters[index].GetY() );
553       infos[i].SetZ( fClusters[index].GetZ() );
554     }
555
556     bool ok = hlt.Merger().FitTrack( t, alpha, t0, alpha, hits1, nHits, 0, infos );
557
558     if ( ok &&  nHits > 15 ) {
559       if ( t.TransportToXWithMaterial( xTPC, hlt.Merger().SliceParam().GetBz( t ) ) ) {
560         if ( t.GetY() > yMax ) {
561           if ( t.Rotate( dAlpha ) ) {
562             alpha += dAlpha;
563             t.TransportToXWithMaterial( xTPC, hlt.Merger().SliceParam().GetBz( t ) );
564           }
565         } else if ( t.GetY() < -yMax ) {
566           if ( t.Rotate( -dAlpha ) ) {
567             alpha += -dAlpha;
568             t.TransportToXWithMaterial( xTPC, hlt.Merger().SliceParam().GetBz( t ) );
569           }
570         }
571       }
572
573       AliTPCtrack tt( *esd );
574       if ( AliHLTTPCCATrackConvertor::GetExtParam( t, tt, alpha ) ) {
575         if ( t.X() > 50 ) esd->UpdateTrackParams( &tt, AliESDtrack::kTPCrefit );
576       }
577     }
578   }
579   return 0;
580 }
581
582 int AliTPCtrackerCA::PropagateBack( AliESDEvent *event )
583 {
584
585   //* backward propagation of ESD tracks
586
587   AliHLTTPCCAStandaloneFramework &hlt = AliHLTTPCCAStandaloneFramework::Instance();
588
589   int nentr = event->GetNumberOfTracks();
590
591   for ( int itr = 0; itr < nentr; itr++ ) {
592
593     AliESDtrack *esd = event->GetTrack( itr );
594     ULong_t status = esd->GetStatus();
595     if ( !( status&AliESDtrack::kTPCin ) ) continue;
596
597     AliHLTTPCCATrackParam t0;
598     AliHLTTPCCATrackConvertor::SetExtParam( t0, *esd  );
599     AliHLTTPCCATrackParam t = t0;
600     float alpha = esd->GetAlpha();
601     //float dEdX=0;
602     AliHLTTPCCAMerger::AliHLTTPCCAClusterInfo infos[500];
603     int hits[500], hits1[500];
604     int nHits = esd->GetTPCclusters( hits );
605     for ( int i = 0; i < nHits; i++ ) {
606       hits1[i] = i;
607       int index = hits[i];
608       infos[i].SetISlice( fClusterSliceRow[index] >> 8 );
609       infos[i].SetX( fClusters[index].GetX() );
610       infos[i].SetY( fClusters[index].GetY() );
611       infos[i].SetZ( fClusters[index].GetZ() );
612     }
613
614     bool ok = hlt.Merger().FitTrack( t, alpha, t0, alpha, hits1, nHits, 1, infos );
615
616     if ( ok &&  nHits > 15 ) {
617       AliTPCtrack tt( *esd );
618       if ( AliHLTTPCCATrackConvertor::GetExtParam( t, tt, alpha ) ) {
619         if ( t.X() > 50 ) esd->UpdateTrackParams( &tt, AliESDtrack::kTPCout );
620       }
621     }
622   }
623   return 0;
624 }