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 const double kCLight = 0.000299792458;
101 float bz = AliTracker::GetBz() * kCLight;
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;
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;
121 float padPitch = 0.4;
122 float sigmaZ = 0.228808;
124 int nRows = fkParam->GetNRowLow() + fkParam->GetNRowUp();
126 for ( int irow = 0; irow < fkParam->GetNRowLow(); irow++ ) {
127 rowX[irow] = fkParam->GetPadRowRadiiLow( irow );
129 for ( int irow = 0; irow < fkParam->GetNRowUp(); irow++ ) {
130 rowX[fkParam->GetNRowLow()+irow] = fkParam->GetPadRowRadiiUp( irow );
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 );
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] );
149 hlt.SliceTracker( iSlice ).Initialize( param );
155 int AliTPCtrackerCA::LoadClusters ( TTree * fromTree )
157 // load clusters to the local arrays
160 delete[] fClusterSliceRow;
162 AliHLTTPCCAStandaloneFramework &hlt = AliHLTTPCCAStandaloneFramework::Instance();
164 if ( fDoHLTPerformance ) {
165 AliHLTTPCCAPerformance::Instance().StartEvent();
166 if ( fDoHLTPerformanceClusters ) AliHLTTPCCAPerformance::Instance().SetDoClusterPulls( 1 );
169 if ( !fkParam ) return 1;
172 while ( fDoHLTPerformance ) {
173 if ( !gAlice ) break;
174 #ifndef HAVE_NOT_ALIRUNLOADER30859
175 AliRunLoader *rl = AliRunLoader::Instance();//gAlice->GetRunLoader();
177 // the old way before rev 30859
178 AliRunLoader *rl = AliRunLoader::GetRunLoader();
181 rl->LoadKinematics();
182 AliStack *stack = rl->Stack();
185 AliHLTTPCCAPerformance::Instance().SetNMCTracks( stack->GetNtrack() );
187 for ( int itr = 0; itr < stack->GetNtrack(); itr++ ) {
188 TParticle *part = stack->Particle( itr );
189 AliHLTTPCCAPerformance::Instance().ReadMCTrack( itr, part );
192 { // check for MC tracks at the TPC entrance
195 isTPC = new bool [stack->GetNtrack()];
196 for ( int i = 0; i < stack->GetNtrack(); i++ ) isTPC[i] = 0;
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;
215 if ( !tpcRef ) continue;
217 if ( isTPC[tpcRef->Label()] ) continue;
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;
225 if ( isTPC ) delete[] isTPC;
228 while ( fDoHLTPerformanceClusters ) {
229 AliTPCLoader *tpcl = ( AliTPCLoader* ) rl->GetDetectorLoader( "TPC" );
231 if ( tpcl->TreeH() == 0x0 ) {
232 if ( tpcl->LoadHits() ) break;
234 if ( tpcl->TreeH() == 0x0 ) break;
236 AliTPC *tpc = ( AliTPC* ) gAlice->GetDetector( "TPC" );
237 int nEnt = ( int )tpcl->TreeH()->GetEntries();
239 for ( int iEnt = 0; iEnt < nEnt; iEnt++ ) {
241 tpcl->TreeH()->GetEvent( iEnt );
242 AliTPChit *phit = ( AliTPChit* )tpc->FirstHit( -1 );
243 for ( ; phit; phit = ( AliTPChit* )tpc->NextHit() ) nPoints++;
245 AliHLTTPCCAPerformance::Instance().SetNMCPoints( nPoints );
247 for ( int iEnt = 0; iEnt < nEnt; iEnt++ ) {
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 );
260 TBranch * br = fromTree->GetBranch( "Segment" );
263 AliTPCClustersRow *clrow = new AliTPCClustersRow;
264 clrow->SetClass( "AliTPCclusterMI" );
265 clrow->SetArray( 0 );
266 clrow->GetArray()->ExpandCreateFast( 10000 );
268 br->SetAddress( &clrow );
271 int nEnt = int( fromTree->GetEntries() );
274 for ( int i = 0; i < nEnt; i++ ) {
277 fkParam->AdjustSectorRow( clrow->GetID(), sec, row );
278 fNClusters += clrow->GetArray()->GetEntriesFast();
281 fClusters = new AliTPCclusterMI [fNClusters];
282 fClusterSliceRow = new unsigned int [fNClusters];
284 hlt.StartDataReading( fNClusters );
286 if ( fDoHLTPerformance ) AliHLTTPCCAPerformance::Instance().SetNHits( fNClusters );
289 for ( int i = 0; i < nEnt; i++ ) {
292 fkParam->AdjustSectorRow( clrow->GetID(), sec, row );
293 int nClu = clrow->GetArray()->GetEntriesFast();
294 float x = fkParam->GetPadRowRadii( sec, row );
298 row = row + fkParam->GetNRowLow();
301 unsigned int sliceRow = ( sec << 8 ) + row;
303 for ( int icl = 0; icl < nClu; icl++ ) {
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 );
313 AliTPCTransform *transform = AliTPCcalibDB::Instance()->GetTransform() ;
315 AliFatal( "Tranformations not in calibDB" );
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){
326 cluster->SetX( xx[0] );
327 cluster->SetY( xx[1] );
328 cluster->SetZ( xx[2] );
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 );
335 // chack Loading of Geo matrices from GeoManager - TEMPORARY FIX
337 cluster->SetX( posC[0] );
338 cluster->SetY( posC[1] );
339 cluster->SetZ( posC[2] );
342 float y = cluster->GetY();
343 float z = cluster->GetZ();
347 fClusters[index] = *cluster;
349 fClusterSliceRow[index] = sliceRow;
351 hlt.ReadCluster( index, sec, row, x, y, z, cluster->GetQ() );
353 if ( fDoHLTPerformance ) AliHLTTPCCAPerformance::Instance().ReadHitLabel( index, lab0, lab1, lab2 );
358 hlt.FinishDataReading();
360 //AliHLTTPCCAPerformance::Instance().SmearClustersMC();
365 AliCluster * AliTPCtrackerCA::GetCluster( int index ) const
367 return &( fClusters[index] );
370 int AliTPCtrackerCA::Clusters2Tracks( AliESDEvent *event )
373 //cout<<"Start of AliTPCtrackerCA"<<endl;
376 AliHLTTPCCAStandaloneFramework &hlt = AliHLTTPCCAStandaloneFramework::Instance();
382 const AliHLTTPCCAMergerOutput &hltOut = *( hlt.Merger().Output() );
384 for ( int itr = 0; itr < hltOut.NTracks(); itr++ ) {
388 const AliHLTTPCCAMergedTrack &tCA = hltOut.Track( itr );
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();
397 firstHit = nhits - 160;
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;
410 tTPC.SetClusterPointer( row, c );
411 AliTPCTrackerPoint &point = *( tTPC.GetTrackPoint( row ) );
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();
418 float x = slice.Row( row ).X();
419 if ( !t0.TransportToX( x, slice.Param().GetBz( t0 ), .999 ) ) continue;
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() ) );
428 tTPC.CookdEdx( 0.02, 0.6 );
430 CookLabel( &tTPC, 0.1 );
432 if ( 1 ) { // correction like in off-line --- Adding systematic error
434 const double *param = AliTPCReconstructor::GetRecoParam()->GetSystematicError();
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 );
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 );
456 event->AddTrack( &tESD );
463 fStatCPUTime += timer.CpuTime();
464 fStatRealTime += timer.RealTime();
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;
469 //cout<<"Do performance.."<<endl;
471 if ( fDoHLTPerformance ) AliHLTTPCCAPerformance::Instance().Performance();
473 if ( 0 ) {// Write Event
474 if ( fStatNEvents == 0 ) {
476 geo.open( "CAEvents/settings.dat", ios::out );
477 if ( geo.is_open() ) {
478 hlt.WriteSettings( geo );
485 sprintf( name, "CAEvents/%i.event.dat", fStatNEvents );
486 hits.open( name, ios::out );
487 if ( hits.is_open() ) {
488 hlt.WriteEvent( hits );
490 sprintf( name, "CAEvents/%i.tracks.dat", fStatNEvents );
491 tracks.open( name, ios::out );
492 hlt.WriteTracks( tracks );
495 if ( fDoHLTPerformance ) {
496 fstream mcevent, mcpoints;
498 sprintf( mcname, "CAEvents/%i.mcevent.dat", fStatNEvents );
499 mcevent.open( mcname, ios::out );
500 if ( mcevent.is_open() ) {
501 AliHLTTPCCAPerformance::Instance().WriteMCEvent( mcevent );
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 );
518 //cout<<"End of AliTPCtrackerCA"<<endl;
523 int AliTPCtrackerCA::RefitInward ( AliESDEvent *event )
525 //* forward propagation of ESD tracks
527 float xTPC = fkParam->GetInnerRadiusLow();
528 float dAlpha = fkParam->GetInnerAngle() / 180.*TMath::Pi();
529 float yMax = xTPC * TMath::Tan( dAlpha / 2. );
531 AliHLTTPCCAStandaloneFramework &hlt = AliHLTTPCCAStandaloneFramework::Instance();
533 int nentr = event->GetNumberOfTracks();
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();
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++ ) {
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() );
556 bool ok = hlt.Merger().FitTrack( t, alpha, t0, alpha, hits1, nHits, 0, infos );
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 ) ) {
563 t.TransportToXWithMaterial( xTPC, hlt.Merger().SliceParam().GetBz( t ) );
565 } else if ( t.GetY() < -yMax ) {
566 if ( t.Rotate( -dAlpha ) ) {
568 t.TransportToXWithMaterial( xTPC, hlt.Merger().SliceParam().GetBz( t ) );
573 AliTPCtrack tt( *esd );
574 if ( AliHLTTPCCATrackConvertor::GetExtParam( t, tt, alpha ) ) {
575 if ( t.X() > 50 ) esd->UpdateTrackParams( &tt, AliESDtrack::kTPCrefit );
582 int AliTPCtrackerCA::PropagateBack( AliESDEvent *event )
585 //* backward propagation of ESD tracks
587 AliHLTTPCCAStandaloneFramework &hlt = AliHLTTPCCAStandaloneFramework::Instance();
589 int nentr = event->GetNumberOfTracks();
591 for ( int itr = 0; itr < nentr; itr++ ) {
593 AliESDtrack *esd = event->GetTrack( itr );
594 ULong_t status = esd->GetStatus();
595 if ( !( status&AliESDtrack::kTPCin ) ) continue;
597 AliHLTTPCCATrackParam t0;
598 AliHLTTPCCATrackConvertor::SetExtParam( t0, *esd );
599 AliHLTTPCCATrackParam t = t0;
600 float alpha = esd->GetAlpha();
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++ ) {
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() );
614 bool ok = hlt.Merger().FitTrack( t, alpha, t0, alpha, hits1, nHits, 1, infos );
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 );