2 // **************************************************************************
3 // This file is property of and copyright by the ALICE HLT Project *
4 // ALICE Experiment at CERN, All rights reserved. *
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. *
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. *
18 //***************************************************************************
20 #include "AliTPCtrackerCA.h"
23 #include "Riostream.h"
24 //#include "AliCluster.h"
25 #include "AliTPCClustersRow.h"
26 #include "AliTPCParam.h"
27 #include "AliTPCClusterParam.h"
30 #include "AliRunLoader.h"
33 #include "AliHLTTPCCAPerformance.h"
34 #include "AliHLTTPCCAParam.h"
35 #include "AliHLTTPCCATrackConvertor.h"
36 #include "AliHLTTPCCATracker.h"
37 #include "AliHLTTPCCAStandaloneFramework.h"
40 #include "AliTPCLoader.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"
54 //#include <fstream.h>
56 ClassImp( AliTPCtrackerCA )
58 AliTPCtrackerCA::AliTPCtrackerCA()
59 : AliTracker(), fkParam( 0 ), fClusters( 0 ), fClusterSliceRow( 0 ), fNClusters( 0 ), fDoHLTPerformance( 0 ), fDoHLTPerformanceClusters( 0 ), fStatCPUTime( 0 ), fStatRealTime( 0 ), fStatNEvents( 0 )
61 //* default constructor
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 )
70 const AliTPCtrackerCA & AliTPCtrackerCA::operator=( const AliTPCtrackerCA& ) const
77 AliTPCtrackerCA::~AliTPCtrackerCA()
81 delete[] fClusterSliceRow;
84 //#include "AliHLTTPCCADisplay.h"
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 )
91 fDoHLTPerformance = 0;
92 fDoHLTPerformanceClusters = 0;
94 AliHLTTPCCAStandaloneFramework &hlt = AliHLTTPCCAStandaloneFramework::Instance();
97 for ( int iSlice = 0; iSlice < hlt.NSlices(); iSlice++ ) {
99 float bz = AliTracker::GetBz();
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;
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;
119 float padPitch = 0.4;
120 float sigmaZ = 0.228808;
122 int nRows = fkParam->GetNRowLow() + fkParam->GetNRowUp();
124 for ( int irow = 0; irow < fkParam->GetNRowLow(); irow++ ) {
125 rowX[irow] = fkParam->GetPadRowRadiiLow( irow );
127 for ( int irow = 0; irow < fkParam->GetNRowUp(); irow++ ) {
128 rowX[fkParam->GetNRowLow()+irow] = fkParam->GetPadRowRadiiUp( irow );
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 );
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 param.SetParamS0Par( iyz, type, k, clparam->fParamS0Par[iyz][type][k] );
147 //hlt.SliceTracker( iSlice ).Initialize( param );
148 hlt.InitializeSliceParam(iSlice, param);
154 int AliTPCtrackerCA::LoadClusters ( TTree * fromTree )
156 // load clusters to the local arrays
159 delete[] fClusterSliceRow;
161 AliHLTTPCCAStandaloneFramework &hlt = AliHLTTPCCAStandaloneFramework::Instance();
163 if ( fDoHLTPerformance ) {
164 AliHLTTPCCAPerformance::Instance().StartEvent();
165 if ( fDoHLTPerformanceClusters ) AliHLTTPCCAPerformance::Instance().SetDoClusterPulls( 1 );
168 if ( !fkParam ) return 1;
171 while ( fDoHLTPerformance ) {
172 if ( !gAlice ) break;
173 #ifndef HAVE_NOT_ALIRUNLOADER30859
174 AliRunLoader *rl = AliRunLoader::Instance();//gAlice->GetRunLoader();
176 // the old way before rev 30859
177 AliRunLoader *rl = AliRunLoader::GetRunLoader();
180 rl->LoadKinematics();
181 AliStack *stack = rl->Stack();
184 AliHLTTPCCAPerformance::Instance().SetNMCTracks( stack->GetNtrack() );
186 for ( int itr = 0; itr < stack->GetNtrack(); itr++ ) {
187 TParticle *part = stack->Particle( itr );
188 AliHLTTPCCAPerformance::Instance().ReadMCTrack( itr, part );
191 { // check for MC tracks at the TPC entrance
194 isTPC = new bool [stack->GetNtrack()];
195 for ( int i = 0; i < stack->GetNtrack(); i++ ) isTPC[i] = 0;
197 TTree *mcTree = rl->TreeTR();
198 if ( !mcTree ) break;
199 TBranch *branch = mcTree->GetBranch( "TrackReferences" );
200 if ( !branch ) break;
201 TClonesArray tpcdummy( "AliTrackReference", 1000 ), *tpcRefs = &tpcdummy;
202 branch->SetAddress( &tpcRefs );
203 int nr = ( int )mcTree->GetEntries();
204 for ( int r = 0; r < nr; r++ ) {
205 mcTree->GetEvent( r );
206 int nref = tpcRefs->GetEntriesFast();
207 if ( !nref ) continue;
208 AliTrackReference *tpcRef = 0x0;
209 for ( int iref = 0; iref < nref; ++iref ) {
210 tpcRef = ( AliTrackReference* )tpcRefs->UncheckedAt( iref );
211 if ( tpcRef->DetectorId() == AliTrackReference::kTPC ) break;
214 if ( !tpcRef ) continue;
216 if ( isTPC[tpcRef->Label()] ) continue;
218 AliHLTTPCCAPerformance::Instance().ReadMCTPCTrack( tpcRef->Label(),
219 tpcRef->X(), tpcRef->Y(), tpcRef->Z(),
220 tpcRef->Px(), tpcRef->Py(), tpcRef->Pz() );
221 isTPC[tpcRef->Label()] = 1;
224 if ( isTPC ) delete[] isTPC;
227 while ( fDoHLTPerformanceClusters ) {
228 AliTPCLoader *tpcl = ( AliTPCLoader* ) rl->GetDetectorLoader( "TPC" );
230 if ( tpcl->TreeH() == 0x0 ) {
231 if ( tpcl->LoadHits() ) break;
233 if ( tpcl->TreeH() == 0x0 ) break;
235 AliTPC *tpc = ( AliTPC* ) gAlice->GetDetector( "TPC" );
236 int nEnt = ( int )tpcl->TreeH()->GetEntries();
238 for ( int iEnt = 0; iEnt < nEnt; iEnt++ ) {
240 tpcl->TreeH()->GetEvent( iEnt );
241 AliTPChit *phit = ( AliTPChit* )tpc->FirstHit( -1 );
242 for ( ; phit; phit = ( AliTPChit* )tpc->NextHit() ) nPoints++;
244 AliHLTTPCCAPerformance::Instance().SetNMCPoints( nPoints );
246 for ( int iEnt = 0; iEnt < nEnt; iEnt++ ) {
248 tpcl->TreeH()->GetEvent( iEnt );
249 AliTPChit *phit = ( AliTPChit* )tpc->FirstHit( -1 );
250 for ( ; phit; phit = ( AliTPChit* )tpc->NextHit() ) {
251 AliHLTTPCCAPerformance::Instance().ReadMCPoint( phit->GetTrack(), phit->X(), phit->Y(), phit->Z(), phit->Time(), phit->fSector % 36 );
259 TBranch * br = fromTree->GetBranch( "Segment" );
262 AliTPCClustersRow *clrow = new AliTPCClustersRow;
263 clrow->SetClass( "AliTPCclusterMI" );
264 clrow->SetArray( 0 );
265 clrow->GetArray()->ExpandCreateFast( 10000 );
267 br->SetAddress( &clrow );
270 int nEnt = int( fromTree->GetEntries() );
273 for ( int i = 0; i < nEnt; i++ ) {
276 fkParam->AdjustSectorRow( clrow->GetID(), sec, row );
277 fNClusters += clrow->GetArray()->GetEntriesFast();
280 fClusters = new AliTPCclusterMI [fNClusters];
281 fClusterSliceRow = new unsigned int [fNClusters];
283 hlt.StartDataReading( fNClusters );
285 if ( fDoHLTPerformance ) AliHLTTPCCAPerformance::Instance().SetNHits( fNClusters );
288 for ( int i = 0; i < nEnt; i++ ) {
291 fkParam->AdjustSectorRow( clrow->GetID(), sec, row );
292 int nClu = clrow->GetArray()->GetEntriesFast();
293 float x = fkParam->GetPadRowRadii( sec, row );
297 row = row + fkParam->GetNRowLow();
300 unsigned int sliceRow = ( sec << 8 ) + row;
302 for ( int icl = 0; icl < nClu; icl++ ) {
306 AliTPCclusterMI* cluster = ( AliTPCclusterMI* )( clrow->GetArray()->At( icl ) );
307 if ( !cluster ) continue;
308 lab0 = cluster->GetLabel( 0 );
309 lab1 = cluster->GetLabel( 1 );
310 lab2 = cluster->GetLabel( 2 );
312 AliTPCTransform *transform = AliTPCcalibDB::Instance()->GetTransform() ;
313 transform->SetCurrentRecoParam((AliTPCRecoParam*)AliTPCReconstructor::GetRecoParam());
316 AliFatal( "Tranformations not in calibDB" );
318 double xx[3] = {cluster->GetRow(), cluster->GetPad(), cluster->GetTimeBin()};
319 int id[1] = {cluster->GetDetector()};
320 transform->Transform( xx, id, 0, 1 );
322 cluster->SetX( xx[0] );
323 cluster->SetY( xx[1] );
324 cluster->SetZ( xx[2] );
326 TGeoHMatrix *mat = fkParam->GetClusterMatrix( cluster->GetDetector() );
327 double pos[3] = {cluster->GetX(), cluster->GetY(), cluster->GetZ()};
328 double posC[3] = {cluster->GetX(), cluster->GetY(), cluster->GetZ()};
329 if ( mat ) mat->LocalToMaster( pos, posC );
331 // chack Loading of Geo matrices from GeoManager - TEMPORARY FIX
333 cluster->SetX( posC[0] );
334 cluster->SetY( posC[1] );
335 cluster->SetZ( posC[2] );
338 float y = cluster->GetY();
339 float z = cluster->GetZ();
343 fClusters[index] = *cluster;
345 fClusterSliceRow[index] = sliceRow;
347 hlt.ReadCluster( index, sec, row, x, y, z, cluster->GetQ() );
349 if ( fDoHLTPerformance ) AliHLTTPCCAPerformance::Instance().ReadHitLabel( index, lab0, lab1, lab2 );
354 hlt.FinishDataReading();
356 //AliHLTTPCCAPerformance::Instance().SmearClustersMC();
361 AliCluster * AliTPCtrackerCA::GetCluster( int index ) const
363 return &( fClusters[index] );
366 int AliTPCtrackerCA::Clusters2Tracks( AliESDEvent *event )
369 //cout<<"Start of AliTPCtrackerCA"<<endl;
372 AliHLTTPCCAStandaloneFramework &hlt = AliHLTTPCCAStandaloneFramework::Instance();
378 const AliHLTTPCCAMergerOutput &hltOut = *( hlt.Merger().Output() );
380 for ( int itr = 0; itr < hltOut.NTracks(); itr++ ) {
384 const AliHLTTPCCAMergedTrack &tCA = hltOut.Track( itr );
386 AliHLTTPCCATrackParam par = tCA.InnerParam();
387 AliHLTTPCCATrackConvertor::GetExtParam( par, tTPC, tCA.InnerAlpha() );
388 tTPC.SetMass( 0.13957 );
389 if ( TMath::Abs( tTPC.GetSigned1Pt() ) > 1. / 0.02 ) continue;
390 int nhits = tCA.NClusters();
393 firstHit = nhits - 160;
396 tTPC.SetNumberOfClusters( nhits );
397 float alpha = tCA.InnerAlpha();
398 AliHLTTPCCATrackParam t0 = par;
399 for ( int ih = 0; ih < nhits; ih++ ) {
400 int index = hltOut.ClusterId( tCA.FirstClusterRef() + firstHit + ih );
401 tTPC.SetClusterIndex( ih, index );
402 AliTPCclusterMI *c = &( fClusters[index] );
403 int iSlice = fClusterSliceRow[index] >> 8;
404 int row = fClusterSliceRow[index] & 0xff;
406 tTPC.SetClusterPointer( row, c );
407 AliTPCTrackerPoint &point = *( tTPC.GetTrackPoint( row ) );
409 //AliHLTTPCCATracker &slice = hlt.SliceTracker( iSlice );
410 if ( hlt.Param(iSlice).Alpha() != alpha ) {
411 if ( ! t0.Rotate( hlt.Param(iSlice).Alpha() - alpha, .999 ) ) continue;
412 alpha = hlt.Param(iSlice).Alpha();
414 float x = hlt.Row(iSlice, row).X();
415 if ( !t0.TransportToX( x, hlt.Param(iSlice).GetBz( t0 ), .999 ) ) continue;
417 //slice.GetErrors2( row, t0, sy2, sz2 );
418 hlt.Param(iSlice).GetClusterErrors2( row, t0.GetZ(), t0.SinPhi(), t0.GetCosPhi(), t0.DzDs(), sy2, sz2 );
419 point.SetSigmaY( c->GetSigmaY2() / sy2 );
420 point.SetSigmaZ( c->GetSigmaZ2() / sz2 );
421 point.SetAngleY( TMath::Abs( t0.GetSinPhi() / t0.GetCosPhi() ) );
422 point.SetAngleZ( TMath::Abs( t0.GetDzDs() ) );
425 tTPC.CookdEdx( 0.02, 0.6 );
427 CookLabel( &tTPC, 0.1 );
429 if ( 1 ) { // correction like in off-line --- Adding systematic error
431 const double *param = AliTPCReconstructor::GetRecoParam()->GetSystematicError();
433 for ( int i = 0; i < 15; i++ ) covar[i] = 0;
434 covar[0] = param[0] * param[0];
435 covar[2] = param[1] * param[1];
436 covar[5] = param[2] * param[2];
437 covar[9] = param[3] * param[3];
438 double facC = AliTracker::GetBz() * kB2C;
439 covar[14] = param[4] * param[4] * facC * facC;
440 tTPC.AddCovariance( covar );
444 tESD.UpdateTrackParams( &( tTPC ), AliESDtrack::kTPCin );
445 //tESD.SetStatus( AliESDtrack::kTPCrefit );
446 //tESD.SetTPCPoints(tTPC.GetPoints());
447 int ndedx = tTPC.GetNCDEDX( 0 );
448 float sdedx = tTPC.GetSDEDX( 0 );
449 float dedx = tTPC.GetdEdx();
450 tESD.SetTPCsignal( dedx, sdedx, ndedx );
453 event->AddTrack( &tESD );
460 fStatCPUTime += timer.CpuTime();
461 fStatRealTime += timer.RealTime();
463 //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;
466 //cout<<"Do performance.."<<endl;
468 if ( fDoHLTPerformance ) AliHLTTPCCAPerformance::Instance().Performance();
470 if ( 0 ) {// Write Event
471 if ( fStatNEvents == 0 ) {
473 geo.open( "CAEvents/settings.dat", ios::out );
474 if ( geo.is_open() ) {
475 hlt.WriteSettings( geo );
482 sprintf( name, "CAEvents/%i.event.dat", fStatNEvents );
483 hits.open( name, ios::out );
484 if ( hits.is_open() ) {
485 hlt.WriteEvent( hits );
487 sprintf( name, "CAEvents/%i.tracks.dat", fStatNEvents );
488 tracks.open( name, ios::out );
489 hlt.WriteTracks( tracks );
492 if ( fDoHLTPerformance ) {
493 fstream mcevent, mcpoints;
495 sprintf( mcname, "CAEvents/%i.mcevent.dat", fStatNEvents );
496 mcevent.open( mcname, ios::out );
497 if ( mcevent.is_open() ) {
498 AliHLTTPCCAPerformance::Instance().WriteMCEvent( mcevent );
500 if ( 1 && fDoHLTPerformanceClusters ) {
501 sprintf( mcname, "CAEvents/%i.mcpoints.dat", fStatNEvents );
502 mcpoints.open( mcname, ios::out );
503 if ( mcpoints.is_open() ) {
504 AliHLTTPCCAPerformance::Instance().WriteMCPoints( mcpoints );
515 //cout<<"End of AliTPCtrackerCA"<<endl;
520 int AliTPCtrackerCA::RefitInward ( AliESDEvent *event )
522 //* forward propagation of ESD tracks
524 float xTPC = fkParam->GetInnerRadiusLow();
525 float dAlpha = fkParam->GetInnerAngle() / 180.*TMath::Pi();
526 float yMax = xTPC * TMath::Tan( dAlpha / 2. );
528 AliHLTTPCCAStandaloneFramework &hlt = AliHLTTPCCAStandaloneFramework::Instance();
530 int nentr = event->GetNumberOfTracks();
532 for ( int itr = 0; itr < nentr; itr++ ) {
533 AliESDtrack *esd = event->GetTrack( itr );
534 ULong_t status = esd->GetStatus();
535 if ( !( status&AliESDtrack::kTPCin ) ) continue;
536 AliHLTTPCCATrackParam t0;
537 AliHLTTPCCATrackConvertor::SetExtParam( t0, *esd );
538 AliHLTTPCCATrackParam t = t0;
539 float alpha = esd->GetAlpha();
541 AliHLTTPCCAMerger::AliHLTTPCCAClusterInfo infos[500];
542 int hits[500], hits1[500];
543 int nHits = esd->GetTPCclusters( hits );
544 for ( int i = 0; i < nHits; i++ ) {
547 infos[i].SetISlice( fClusterSliceRow[index] >> 8 );
548 infos[i].SetIRow( fClusterSliceRow[index] & 0xff );
549 infos[i].SetId( index );
550 infos[i].SetX( fClusters[index].GetX() );
551 infos[i].SetY( fClusters[index].GetY() );
552 infos[i].SetZ( fClusters[index].GetZ() );
555 bool ok = hlt.Merger().FitTrack( t, alpha, t0, alpha, hits1, nHits, 0, 1,infos );
557 if ( ok && nHits > 15 ) {
558 if ( t.TransportToXWithMaterial( xTPC, hlt.Merger().SliceParam().GetBz( t ) ) ) {
559 if ( t.GetY() > yMax ) {
560 if ( t.Rotate( dAlpha ) ) {
562 t.TransportToXWithMaterial( xTPC, hlt.Merger().SliceParam().GetBz( t ) );
564 } else if ( t.GetY() < -yMax ) {
565 if ( t.Rotate( -dAlpha ) ) {
567 t.TransportToXWithMaterial( xTPC, hlt.Merger().SliceParam().GetBz( t ) );
572 AliTPCtrack tt( *esd );
573 if ( AliHLTTPCCATrackConvertor::GetExtParam( t, tt, alpha ) ) {
574 if ( t.X() > 50 ) esd->UpdateTrackParams( &tt, AliESDtrack::kTPCrefit );
581 int AliTPCtrackerCA::PropagateBack( AliESDEvent *event )
584 //* backward propagation of ESD tracks
586 AliHLTTPCCAStandaloneFramework &hlt = AliHLTTPCCAStandaloneFramework::Instance();
588 int nentr = event->GetNumberOfTracks();
590 for ( int itr = 0; itr < nentr; itr++ ) {
592 AliESDtrack *esd = event->GetTrack( itr );
593 ULong_t status = esd->GetStatus();
594 if ( !( status&AliESDtrack::kTPCin ) ) continue;
596 AliHLTTPCCATrackParam t0;
597 AliHLTTPCCATrackConvertor::SetExtParam( t0, *esd );
598 AliHLTTPCCATrackParam t = t0;
599 float alpha = esd->GetAlpha();
601 AliHLTTPCCAMerger::AliHLTTPCCAClusterInfo infos[500];
602 int hits[500], hits1[500];
603 int nHits = esd->GetTPCclusters( hits );
604 for ( int i = 0; i < nHits; i++ ) {
607 infos[i].SetISlice( fClusterSliceRow[index] >> 8 );
608 infos[i].SetIRow( fClusterSliceRow[index] & 0xff );
609 infos[i].SetId( index );
610 infos[i].SetX( fClusters[index].GetX() );
611 infos[i].SetY( fClusters[index].GetY() );
612 infos[i].SetZ( fClusters[index].GetZ() );
615 bool ok = hlt.Merger().FitTrack( t, alpha, t0, alpha, hits1, nHits, 1, 1, infos );
617 if ( ok && nHits > 15 ) {
618 AliTPCtrack tt( *esd );
619 if ( AliHLTTPCCATrackConvertor::GetExtParam( t, tt, alpha ) ) {
620 if ( t.X() > 50 ) esd->UpdateTrackParams( &tt, AliESDtrack::kTPCout );